Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
F01PhysicsList.cc
Go to the documentation of this file.
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 //
28 //
29 //
30 // $Id$
31 //
32 
33 #include "G4Timer.hh"
34 
35 #include "F01PhysicsList.hh"
38 
39 #include "G4ProcessManager.hh"
40 #include "G4ProcessVector.hh"
41 #include "G4ParticleTypes.hh"
42 #include "G4ParticleTable.hh"
43 #include "G4Material.hh"
44 #include "G4EnergyLossTables.hh"
45 #include "G4UnitsTable.hh"
46 #include "G4SystemOfUnits.hh"
47 #include "G4ios.hh"
48 #include <iomanip>
49 
50 
51 
53 //
54 //
55 
57 : G4VUserPhysicsList(), MaxChargedStep(DBL_MAX),
58  thePhotoElectricEffect(0), theComptonScattering(0),
59  theGammaConversion(0), theeminusIonisation(0),
60  theeminusBremsstrahlung(0), theeplusIonisation(0),
61  theeplusBremsstrahlung(0), theeplusAnnihilation(0),
62  theeminusStepCut(0),theeplusStepCut(0)
63 {
64  pDet = p;
65 
66  defaultCutValue = 1.000*mm ;
67 
68  cutForGamma = defaultCutValue ;
69  cutForElectron = defaultCutValue ;
70 
71  SetVerboseLevel(1);
72  physicsListMessenger = new F01PhysicsListMessenger(this);
73 }
74 
76 //
77 //
78 
80 {
81  delete physicsListMessenger;
82 }
83 
85 //
86 //
87 
89 {
90  // In this method, static member functions should be called
91  // for all particles which you want to use.
92  // This ensures that objects of these particle types will be
93  // created in the program.
94 
99 }
100 
102 //
103 //
104 
106 {
107  // gamma
108 
110 
111  // charged geantino
112 
114 
115 
116 }
117 
119 {
120  // leptons
121 
126 
131 }
132 
134 {
135  // mesons
136 
142 }
143 
144 
146 {
147 // barions
148 
151 }
152 
153 
155 //
156 //
157 
159 {
161  // AddParameterisation();
162 
163  ConstructEM();
165 }
166 
168 //
169 //
170 
171 #include "G4ComptonScattering.hh"
172 #include "G4GammaConversion.hh"
173 #include "G4PhotoElectricEffect.hh"
174 
175 #include "G4eMultipleScattering.hh"
176 #include "G4MuMultipleScattering.hh"
177 #include "G4hMultipleScattering.hh"
178 
179 #include "G4eIonisation.hh"
180 #include "G4eBremsstrahlung.hh"
181 #include "G4eplusAnnihilation.hh"
182 
183 #include "G4MuIonisation.hh"
184 #include "G4MuBremsstrahlung.hh"
185 #include "G4MuPairProduction.hh"
186 
187 #include "G4hIonisation.hh"
188 
189 #include "F01StepCut.hh"
190 
192 {
194 
195  while( (*theParticleIterator)() )
196  {
198  G4ProcessManager* pmanager = particle->GetProcessManager();
199  G4String particleName = particle->GetParticleName();
200 
201  if (particleName == "gamma")
202  {
203  // Construct processes for gamma
204 
205  thePhotoElectricEffect = new G4PhotoElectricEffect();
206  theComptonScattering = new G4ComptonScattering();
207  theGammaConversion = new G4GammaConversion();
208 
209  pmanager->AddDiscreteProcess(thePhotoElectricEffect);
210  pmanager->AddDiscreteProcess(theComptonScattering);
211 
212  pmanager->AddDiscreteProcess(theGammaConversion);
213 
214  }
215  else if (particleName == "e-")
216  {
217  // Construct processes for electron
218 
219  // theeminusMultipleScattering = new G4eMultipleScattering();
220  theeminusIonisation = new G4eIonisation();
221  theeminusBremsstrahlung = new G4eBremsstrahlung();
222 
223 
224  theeminusStepCut = new F01StepCut();
225 
226  // pmanager->AddProcess(theeminusMultipleScattering,-1,1,1);
227 
228  pmanager->AddProcess(theeminusIonisation,-1,2,2);
229 
230 
231  pmanager->AddProcess(theeminusBremsstrahlung,-1,-1,3);
232 
233 
234  pmanager->AddProcess(theeminusStepCut,-1,-1,4);
235  theeminusStepCut->SetMaxStep(MaxChargedStep) ;
236 
237  }
238  else if (particleName == "e+")
239  {
240  // Construct processes for positron
241 
242  // theeplusMultipleScattering = new G4eMultipleScattering();
243  theeplusIonisation = new G4eIonisation();
244  theeplusBremsstrahlung = new G4eBremsstrahlung();
245  // theeplusAnnihilation = new G4eplusAnnihilation();
246 
247 
248  theeplusStepCut = new F01StepCut();
249 
250  // pmanager->AddProcess(theeplusMultipleScattering,-1,1,1);
251  pmanager->AddProcess(theeplusIonisation,-1,2,2);
252  pmanager->AddProcess(theeplusBremsstrahlung,-1,-1,3);
253  // pmanager->AddProcess(theeplusAnnihilation,0,-1,4);
254 
255 
256  pmanager->AddProcess(theeplusStepCut,-1,-1,5);
257  theeplusStepCut->SetMaxStep(MaxChargedStep) ;
258 
259  }
260  else if( particleName == "mu+" ||
261  particleName == "mu-" )
262  {
263  // Construct processes for muon+
264 
265  F01StepCut* muonStepCut = new F01StepCut();
266 
267  G4MuIonisation* themuIonisation = new G4MuIonisation() ;
268  pmanager->AddProcess(new G4MuMultipleScattering(),-1,1,1);
269  pmanager->AddProcess(themuIonisation,-1,2,2);
270  pmanager->AddProcess(new G4MuBremsstrahlung(),-1,-1,3);
271  pmanager->AddProcess(new G4MuPairProduction(),-1,-1,4);
272 
273  pmanager->AddProcess( muonStepCut,-1,-1,3);
274  muonStepCut->SetMaxStep(MaxChargedStep) ;
275 
276  }
277  else if (
278  particleName == "proton"
279  || particleName == "antiproton"
280  || particleName == "pi+"
281  || particleName == "pi-"
282  || particleName == "kaon+"
283  || particleName == "kaon-"
284  )
285  {
286  F01StepCut* thehadronStepCut = new F01StepCut();
287 
288  G4hIonisation* thehIonisation = new G4hIonisation() ;
289  G4hMultipleScattering* thehMultipleScattering =
290  new G4hMultipleScattering() ;
291 
292 
293  pmanager->AddProcess(thehMultipleScattering,-1,1,1);
294  pmanager->AddProcess(thehIonisation,-1,2,2);
295 
296 
297  pmanager->AddProcess( thehadronStepCut,-1,-1,3);
298  thehadronStepCut->SetMaxStep(MaxChargedStep) ;
299  thehadronStepCut->SetMaxStep(10*mm) ;
300 
301  }
302  }
303 }
304 
305 
306 #include "G4Decay.hh"
307 
309 {
310  // Add Decay Process
311 
312  G4Decay* theDecayProcess = new G4Decay();
314 
315  while( (*theParticleIterator)() )
316  {
318  G4ProcessManager* pmanager = particle->GetProcessManager();
319 
320  if (theDecayProcess->IsApplicable(*particle))
321  {
322  pmanager ->AddProcess(theDecayProcess);
323 
324  // set ordering for PostStepDoIt and AtRestDoIt
325 
326  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
327  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
328  }
329  }
330 }
331 
332 
333 
335 
337 {
338  G4Timer theTimer ;
339  theTimer.Start() ;
340  if (verboseLevel >0)
341  {
342  G4cout << "F01PhysicsList::SetCuts:";
343  G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
344  }
345  // set cut values for gamma at first and for e- second and next for e+,
346  // because some processes for e+/e- need cut values for gamma
347 
348  SetCutValue(cutForGamma,"gamma");
349 
350  SetCutValue(cutForElectron,"e-");
351  SetCutValue(cutForElectron,"e+");
352 
354 
355  theTimer.Stop();
356  G4cout.precision(6);
357  G4cout << G4endl ;
358  G4cout << "total time(SetCuts)=" << theTimer.GetUserElapsed() << " s " <<G4endl;
359 
360 }
361 
363 
365 {
366  cutForGamma = val;
367 }
368 
370 
372 {
373  cutForElectron = val;
374 }
375 
376 
378 
380 {
381  MaxChargedStep = step ;
382  G4cout << " MaxChargedStep=" << MaxChargedStep << G4endl;
383  G4cout << G4endl;
384 }