source: trunk/source/processes/electromagnetic/lowenergy/test/hTest/src/hTestPhysicsList.cc @ 1199

Last change on this file since 1199 was 1199, checked in by garnier, 15 years ago

nvx fichiers dans CVS

File size: 12.5 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer                                           *
4// *                                                                  *
5// * The  Geant4 software  is  copyright of the Copyright Holders  of *
6// * the Geant4 Collaboration.  It is provided  under  the terms  and *
7// * conditions of the Geant4 Software License,  included in the file *
8// * LICENSE and available at  http://cern.ch/geant4/license .  These *
9// * include a list of copyright holders.                             *
10// *                                                                  *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work  make  any representation or  warranty, express or implied, *
14// * regarding  this  software system or assume any liability for its *
15// * use.  Please see the license in the file  LICENSE  and URL above *
16// * for the full disclaimer and the limitation of liability.         *
17// *                                                                  *
18// * This  code  implementation is the result of  the  scientific and *
19// * technical work of the GEANT4 collaboration.                      *
20// * By using,  copying,  modifying or  distributing the software (or *
21// * any work based  on the software)  you  agree  to acknowledge its *
22// * use  in  resulting  scientific  publications,  and indicate your *
23// * acceptance of all terms of the Geant4 Software license.          *
24// ********************************************************************
25//
26//---------------------------------------------------------------------------
27//
28//---------------------------------------------------------------------------
29//
30// ClassName:   hTestPhysicsList
31// 
32// Description: hTest PhysicsList for Geant4 tests
33//
34// Authors:   07.04.01  V.Ivanchenko
35//
36// Modified:
37//
38//----------------------------------------------------------------------------
39//
40
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
43//
44
45#include "hTestPhysicsList.hh"
46#include "hTestDetectorConstruction.hh"
47#include "hTestPhysicsListMessenger.hh"
48#include "hTestLowEPhysicsList.hh"
49#include "hTestStanPhysicsList.hh"
50#include "hTestHadronPhysicsList.hh"
51#include "hTestHadronPhysicsList1.hh"
52#include "hTestStepCut.hh"
53
54#include "G4ParticleDefinition.hh"
55#include "G4ParticleWithCuts.hh"
56#include "G4IonC12.hh"
57#include "G4IonAr40.hh"
58#include "G4IonFe56.hh"
59#include "G4ProcessManager.hh"
60#include "G4ProcessVector.hh"
61#include "G4ParticleTypes.hh"
62#include "G4ParticleTable.hh"
63#include "G4Material.hh"
64#include "G4EnergyLossTables.hh"
65#include "G4ios.hh"
66#include <iomanip>                
67#include "G4Decay.hh"
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70
71hTestPhysicsList::hTestPhysicsList(const hTestDetectorConstruction* p):
72  pDet(p),
73  physicsIsDefined(false)
74{
75  InitializeMe(); 
76}
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79
80void hTestPhysicsList::InitializeMe()
81{
82  verbose = pDet->GetVerbose();
83 
84  // default cuts
85  cutForGamma     = 1.0*cm;
86  cutForElectron  = 1.0*mm;
87  cutForProton    = 1.0*mm;
88  maxChargedStep  = DBL_MAX; 
89  lowEnergyLimit  = 10.0*eV;
90  highEnergyLimit  = 100.0*GeV;
91
92  theMessenger = new hTestPhysicsListMessenger(this);
93
94  emPhysics = G4String("");
95  hadronPhysics = G4String("");
96  decayPhysics = G4String("none");
97  SetEMPhysicsList(G4String("LowEnergy")); 
98  SetHadronPhysicsList(G4String("none"));
99
100}
101
102//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
103
104hTestPhysicsList::~hTestPhysicsList()
105{
106  delete theMessenger; 
107}
108
109//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
110
111void hTestPhysicsList::ConstructParticle()
112{
113  // In this method, static member functions should be called
114  // for all particles which you want to use.
115  // This ensures that objects of these particle types will be
116  // created in the program.
117
118  ConstructMyBosons();
119  ConstructMyLeptons();
120  ConstructMyBarions();
121  ConstructMyMesons();
122  ConstructMyIons();
123}
124
125//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
126
127void hTestPhysicsList::ConstructMyBosons()
128{
129  // pseudo-particles
130  G4Geantino::GeantinoDefinition();
131  G4ChargedGeantino::ChargedGeantinoDefinition();
132 
133  // gamma
134  G4Gamma::GammaDefinition();
135}
136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
137
138void hTestPhysicsList::ConstructMyLeptons()
139{
140  // leptons
141  G4Electron::ElectronDefinition();
142  G4Positron::PositronDefinition();
143  G4MuonPlus::MuonPlusDefinition();
144  G4MuonMinus::MuonMinusDefinition();
145
146  G4NeutrinoE::NeutrinoEDefinition();
147  G4AntiNeutrinoE::AntiNeutrinoEDefinition();
148  G4NeutrinoMu::NeutrinoMuDefinition();
149  G4AntiNeutrinoMu::AntiNeutrinoMuDefinition();
150}
151
152//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
153
154void hTestPhysicsList::ConstructMyMesons()
155{
156 //  mesons
157  G4PionPlus::PionPlusDefinition();
158  G4PionMinus::PionMinusDefinition();
159  G4PionZero::PionZeroDefinition();
160  G4KaonPlus::KaonPlusDefinition();
161  G4KaonMinus::KaonMinusDefinition();
162  G4Eta::EtaDefinition();
163  G4EtaPrime::EtaPrimeDefinition();
164  G4KaonZero::KaonZeroDefinition();
165  G4AntiKaonZero::AntiKaonZeroDefinition();
166  G4KaonZeroLong::KaonZeroLongDefinition();
167  G4KaonZeroShort::KaonZeroShortDefinition();
168}
169
170//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
171
172void hTestPhysicsList::ConstructMyBarions()
173{
174//  barions
175  G4Proton::ProtonDefinition();
176  G4AntiProton::AntiProtonDefinition();
177  G4Neutron::NeutronDefinition();
178  G4AntiNeutron::AntiNeutronDefinition();
179  G4Lambda::LambdaDefinition();
180  G4SigmaZero::SigmaZeroDefinition();
181  G4SigmaPlus::SigmaPlusDefinition();
182  G4SigmaMinus::SigmaMinusDefinition();
183}
184
185//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
186
187void hTestPhysicsList::ConstructMyIons()
188{
189//  Ions
190  G4Deuteron::DeuteronDefinition();
191  G4Triton::TritonDefinition();
192  G4Alpha::AlphaDefinition();
193  G4IonC12::IonC12Definition();
194  G4IonAr40::IonAr40Definition();
195  G4IonFe56::IonFe56Definition();
196  G4GenericIon::GenericIonDefinition();
197}
198
199//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
200
201void hTestPhysicsList::ConstructProcess()
202{
203  verbose = pDet->GetVerbose();
204  AddTransportation();
205  if(theEMList)  theEMList->ConstructEM();
206  if(theHadList) theHadList->ConstructHad();
207  ConstructDecay();
208  physicsIsDefined = true;
209}
210
211//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
212
213void hTestPhysicsList::ConstructDecay()
214{
215  // Add Decay Process
216  G4Decay* theDecayProcess = new G4Decay();
217  theParticleIterator->reset();
218  hTestStepCut* theStepCut = new hTestStepCut(G4String("user cut"),this);
219
220  while( (*theParticleIterator)() ){
221    G4ParticleDefinition* particle = theParticleIterator->value();
222    G4ProcessManager* pmanager = particle->GetProcessManager();
223
224    if (theDecayProcess->IsApplicable(*particle) && decayPhysics != "none") { 
225      pmanager ->AddProcess(theDecayProcess);
226
227      // set ordering for PostStepDoIt and AtRestDoIt
228      pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
229      pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
230    }
231    if (particle->GetPDGCharge()) { 
232      pmanager ->AddProcess(theStepCut);
233      pmanager ->SetProcessOrdering(theStepCut, idxPostStep);
234    }
235  }
236}
237
238//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
239
240void hTestPhysicsList::SetCuts()
241{
242
243  //special for low energy physics
244  G4Gamma   ::SetEnergyRange(lowEnergyLimit,highEnergyLimit);
245  G4Electron::SetEnergyRange(lowEnergyLimit,highEnergyLimit);
246  G4Positron::SetEnergyRange(lowEnergyLimit,highEnergyLimit);
247   
248  if (verbose > 0){
249    G4cout << "hTestPhysicsList::SetCuts: "
250           << "CutLength = " << maxChargedStep/mm << " mm" << G4endl;
251  } 
252  // set cut values for gamma at first and for e- second and next for e+,
253  // because some processes for e+/e- need cut values for gamma
254
255   G4cout << "hTestPhysicsList: Set cuts for all particles " << G4endl; 
256
257   SetCutValue(cutForGamma,"gamma");
258
259   SetCutValue(cutForElectron,"e-");
260   SetCutValue(cutForElectron,"e+");
261
262   SetCutValue(cutForProton,"mu-");
263   SetCutValue(cutForProton,"mu+");
264
265  // set cut values for proton and anti_proton before all other hadrons
266  // because some processes for hadrons need cut values for proton/anti_proton
267
268  SetCutValue(cutForProton, "proton");
269
270  SetCutValue(cutForProton, "anti_proton");
271
272  SetCutValueForOthers(cutForProton);
273
274  if (verbose > 1) {
275    G4cout << "hTestPhysicsList: Dump the table" << G4endl;
276  }
277  DumpCutValuesTable();
278}
279
280//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
281
282void hTestPhysicsList::SetGammaCut(G4double val)
283{
284  if(physicsIsDefined) ResetCuts();
285  cutForGamma = val;
286}
287
288//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
289
290void hTestPhysicsList::SetElectronCut(G4double val)
291{
292  if(physicsIsDefined) ResetCuts();
293  cutForElectron = val;
294}
295
296//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
297
298void hTestPhysicsList::SetProtonCut(G4double val)
299{
300  if(physicsIsDefined) ResetCuts();
301  cutForProton = val;
302}
303
304//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
305
306void hTestPhysicsList::SetElectronCutByEnergy(G4double energy)
307{
308  G4Material* currMat = pDet->GetAbsorberMaterial();
309  G4ParticleDefinition* part = G4Electron::ElectronDefinition();
310
311  // Get range from electron energy and set it as a cut
312  SetElectronCut(G4EnergyLossTables::GetRange(part,energy,currMat));
313}
314
315//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
316
317void hTestPhysicsList::SetLowEnergyLimit(G4double energy)
318{
319  if(physicsIsDefined) ResetCuts();
320  lowEnergyLimit = energy;
321  G4cout << "hTestPhysicsList: lowEnergyLimit = " 
322         << lowEnergyLimit/eV << " eV" << G4endl;
323}
324
325//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
326
327void hTestPhysicsList::SetHighEnergyLimit(G4double energy)
328{
329  if(physicsIsDefined) ResetCuts();
330  highEnergyLimit = energy;
331  G4cout << "hTestPhysicsList: highEnergyLimit = " 
332         << highEnergyLimit/GeV << " GeV" << G4endl;
333}
334
335//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
336
337void hTestPhysicsList::SetMaxStep(G4double step)
338{
339  maxChargedStep = step;
340  G4cout << "hTestPhysicsList: MaxChargedStep = " 
341         << maxChargedStep/mm << " mm" << G4endl;
342}
343
344//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
345
346void hTestPhysicsList::SetEMPhysicsList(const G4String& name)
347{
348  verbose = pDet->GetVerbose();
349  // Name is not changed
350  if(name == emPhysics) return;
351
352  // Define hadronic process class
353
354  if(name == "Standard") {
355    theEMList = new hTestStanPhysicsList();
356    theEMList->SetVerbose(verbose);
357
358  } else if (name == "LowEnergy") {
359    theEMList = new hTestLowEPhysicsList();
360    theEMList->SetVerbose(verbose);
361
362  } else if (name == "none") {
363    theEMList = 0;
364
365  } else {
366    G4cout << "HsPhysicsList: There are no EMPhysicsList called <"
367           << name << ">, so old one is used" << G4endl;
368    return;
369  }
370
371  emPhysics = name;   
372  G4cout << "HsPhysicsList: <" << name
373         << "> EMPhysicsList is set" << G4endl;
374  if(physicsIsDefined) ResetCuts(); 
375}
376
377//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
378
379void hTestPhysicsList::SetHadronPhysicsList(const G4String& name)
380{
381  verbose = pDet->GetVerbose();
382  // Name is not changed
383  if(name == hadronPhysics) return;
384
385  // Define hadronic process class
386  G4String stand = G4String("Standard");
387  G4String hado  = G4String("LowEnergy");
388  G4String frag  = G4String("IonFragmentation");
389  G4String none  = G4String("none");
390
391  if(name == "LEparametrised") {
392    theHadList = new hTestHadronPhysicsList();
393    theHadList->SetVerbose(verbose);
394
395  } else if (name == "CHIPS") {
396    theHadList = new hTestHadronPhysicsList1();
397    theHadList->SetVerbose(verbose);
398
399  } else if (name == "none") {
400    theHadList = 0;
401
402  } else {
403    G4cout << "HsPhysicsList: There are no HadronicPhysicsList called <"
404           << name << ">, so old one is used" << G4endl;
405    return;
406  }
407
408  hadronPhysics = name;   
409  G4cout << "HsPhysicsList: <" << name
410         << "> HadronicPhysicsList is set" << G4endl;
411  if(physicsIsDefined) ResetCuts(); 
412}
413
414//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
415
416
417
418
419
420
Note: See TracBrowser for help on using the repository browser.