Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ExGflashPhysicsList.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 #include "ExGflashPhysicsList.hh"
30 
31 #include "globals.hh"
32 #include "G4ParticleDefinition.hh"
33 #include "G4ParticleWithCuts.hh"
34 #include "G4ProcessManager.hh"
35 #include "G4ProcessVector.hh"
36 #include "G4ParticleTypes.hh"
37 #include "G4ParticleTable.hh"
38 #include "G4Material.hh"
39 #include "G4MaterialTable.hh"
40 #include "G4Region.hh"
41 #include "G4RegionStore.hh"
42 #include "G4ios.hh"
43 #include <iomanip>
44 
46 
47 using namespace std;
48 
50 {
51  SetVerboseLevel(0);
52 }
53 
55 {
56 }
57 
59 {
60  // In this method, static member functions should be called
61  // for all particles which you want to use.
62  // This ensures that objects of these particle types will be
63  // created in the program.
64 
65  std::cout<<"start construct particle"<<std::endl;
70  ConstructIons();
71  std::cout<<"end construct particle"<<std::endl;
72 }
73 
75 {
76  // pseudo-particles
79 
80  // gamma
82 
83  // optical photon
85 }
86 
87 #include "G4LeptonConstructor.hh"
89 {
90  // Construct all leptons
91  G4LeptonConstructor pConstructor;
92  pConstructor.ConstructParticle();
93 }
94 
95 #include "G4MesonConstructor.hh"
97 {
98  // Construct all mesons
99  G4MesonConstructor pConstructor;
100  pConstructor.ConstructParticle();
101 }
102 
103 #include "G4BaryonConstructor.hh"
105 {
106  // Construct all barions
107  G4BaryonConstructor pConstructor;
108  pConstructor.ConstructParticle();
109 }
110 
111 #include "G4IonConstructor.hh"
113 {
114  // Construct light ions
115  G4IonConstructor pConstructor;
116  pConstructor.ConstructParticle();
117 }
118 
120 {
121  // std::cout<<"1111"<<std::endl;
123  // std::cout<<"2222"<<std::endl;
125  std::cout<<"AddParameterisation"<<std::endl;
126 
127  ConstructEM();
128  std::cout<<"ConstructEM"<<std::endl;
130  // std::cout<<"5555"<<std::endl;
131 }
132 
133 
135 {
137 }
138 
139 #include "G4ComptonScattering.hh"
140 #include "G4GammaConversion.hh"
141 #include "G4PhotoElectricEffect.hh"
142 
143 #include "G4eMultipleScattering.hh"
144 #include "G4MuMultipleScattering.hh"
145 #include "G4hMultipleScattering.hh"
146 
147 #include "G4eIonisation.hh"
148 #include "G4eBremsstrahlung.hh"
149 #include "G4eplusAnnihilation.hh"
150 
151 #include "G4UserLimits.hh"
152 
153 #include "G4MuIonisation.hh"
154 #include "G4MuBremsstrahlung.hh"
155 #include "G4MuPairProduction.hh"
156 
157 #include "G4hIonisation.hh"
159 {
160 
161  G4cout<<"Physics List constructor"<<G4endl;
162  SetCuts();
164  while( (*theParticleIterator)() ){
166  G4ProcessManager* pmanager = particle->GetProcessManager();
167  G4String particleName = particle->GetParticleName();
168 
169  if (particleName == "gamma") {
170  // gamma
171  // Construct processes for gamma
172  G4VProcess* theGammaConversion = new G4GammaConversion();
173  G4VProcess* theComptonScattering = new G4ComptonScattering();
174  G4VProcess* thePhotoElectricEffect = new G4PhotoElectricEffect();
175  // G4VProcess* thegammacut = new G4UserLimits();
176  // thegammacut->SetUserMinEkine(1.0*MeV);
177 
178  pmanager->AddDiscreteProcess(theGammaConversion);
179  pmanager->AddDiscreteProcess(theComptonScattering);
180  pmanager->AddDiscreteProcess(thePhotoElectricEffect);
181  // pmanager->AddProcess(thegammacut);
182  // G4cout <<"theGammaConversion" << theGammaConversion <<endl;
183  //G4cout <<"theComptonScattering" << theComptonScattering <<endl;
184  //G4cout <<"thePhotoElectricEffect" << thePhotoElectricEffect <<endl;
185 
186  } else if (particleName == "e-") {
187  //electron
188  // Construct processes for electron
189  G4VProcess* theeminusMultipleScattering = new G4eMultipleScattering();
190  G4VProcess* theeminusIonisation = new G4eIonisation();
191  G4VProcess* theeminusBremsstrahlung = new G4eBremsstrahlung();
192  // G4VProcess* theeminuscut = new G4UserLimits();
193  // theeminuscut->SetUserMinEkine(1.0*MeV);
194  // add processes
195  pmanager->AddProcess(theeminusMultipleScattering);
196  pmanager->AddProcess(theeminusIonisation);
197  pmanager->AddProcess(theeminusBremsstrahlung);
198 
199  // pmanager->AddProcess( theeminuscut);
200 
201 
202 
203  // set ordering for AlongStepDoIt
204  pmanager->SetProcessOrdering(theeminusMultipleScattering, idxAlongStep, 1);
205  pmanager->SetProcessOrdering(theeminusIonisation, idxAlongStep, 2);
206  // set ordering for PostStepDoIt
207  pmanager->SetProcessOrdering(theeminusMultipleScattering, idxPostStep, 1);
208  pmanager->SetProcessOrdering(theeminusIonisation, idxPostStep, 2);
209  pmanager->SetProcessOrdering(theeminusBremsstrahlung, idxPostStep, 3);
210 
211  //G4cout <<"theeminusMultipleScattering" << theeminusMultipleScattering <<endl;
212  //G4cout <<"theeminusIonisation" << theeminusIonisation <<endl;
213  //G4cout <<"theeminusBremsstrahlung" << theeminusBremsstrahlung <<endl;
214 
215  } else if (particleName == "e+") {
216  //positron
217  // Construct processes for positron
218  G4VProcess* theeplusMultipleScattering = new G4eMultipleScattering();
219  G4VProcess* theeplusIonisation = new G4eIonisation();
220  G4VProcess* theeplusBremsstrahlung = new G4eBremsstrahlung();
221  G4VProcess* theeplusAnnihilation = new G4eplusAnnihilation();
222  // add processes
223  pmanager->AddProcess(theeplusMultipleScattering);
224  pmanager->AddProcess(theeplusIonisation);
225  pmanager->AddProcess(theeplusBremsstrahlung);
226  pmanager->AddProcess(theeplusAnnihilation);
227  // set ordering for AtRestDoIt
228  pmanager->SetProcessOrderingToFirst(theeplusAnnihilation, idxAtRest);
229  // set ordering for AlongStepDoIt
230  pmanager->SetProcessOrdering(theeplusMultipleScattering, idxAlongStep, 1);
231  pmanager->SetProcessOrdering(theeplusIonisation, idxAlongStep, 2);
232  // set ordering for PostStepDoIt
233  pmanager->SetProcessOrdering(theeplusMultipleScattering, idxPostStep, 1);
234  pmanager->SetProcessOrdering(theeplusIonisation, idxPostStep, 2);
235  pmanager->SetProcessOrdering(theeplusBremsstrahlung, idxPostStep, 3);
236  pmanager->SetProcessOrdering(theeplusAnnihilation, idxPostStep, 4);
237 
238  //G4cout <<"theeplusMultipleScattering" << theeplusMultipleScattering <<endl;
239  //G4cout <<"theeplusIonisation" << theeplusIonisation <<endl;
240  //G4cout <<"theeplusBremsstrahlung" << theeplusBremsstrahlung <<endl;
241 
242  } else if( particleName == "mu+" ||
243  particleName == "mu-" ) {
244  //muon
245  // Construct processes for muon+
246  G4VProcess* aMultipleScattering = new G4MuMultipleScattering();
247  G4VProcess* aBremsstrahlung = new G4MuBremsstrahlung();
248  G4VProcess* aPairProduction = new G4MuPairProduction();
249  G4VProcess* anIonisation = new G4MuIonisation();
250  // add processes
251  pmanager->AddProcess(anIonisation);
252  pmanager->AddProcess(aMultipleScattering);
253  pmanager->AddProcess(aBremsstrahlung);
254  pmanager->AddProcess(aPairProduction);
255  // set ordering for AlongStepDoIt
256  pmanager->SetProcessOrdering(aMultipleScattering, idxAlongStep, 1);
257  pmanager->SetProcessOrdering(anIonisation, idxAlongStep, 2);
258  // set ordering for PostStepDoIt
259  pmanager->SetProcessOrdering(aMultipleScattering, idxPostStep, 1);
260  pmanager->SetProcessOrdering(anIonisation, idxPostStep, 2);
261  pmanager->SetProcessOrdering(aBremsstrahlung, idxPostStep, 3);
262  pmanager->SetProcessOrdering(aPairProduction, idxPostStep, 4);
263 
264  } else if ((!particle->IsShortLived()) &&
265  (particle->GetPDGCharge() != 0.0) &&
266  (particle->GetParticleName() != "chargedgeantino")) {
267  // all others charged particles except geantino
268  G4VProcess* aMultipleScattering = new G4hMultipleScattering();
269  G4VProcess* anIonisation = new G4hIonisation();
270  // add processes
271  pmanager->AddProcess(anIonisation);
272  pmanager->AddProcess(aMultipleScattering);
273  // set ordering for AlongStepDoIt
274  pmanager->SetProcessOrdering(aMultipleScattering, idxAlongStep, 1);
275  pmanager->SetProcessOrdering(anIonisation, idxAlongStep, 2);
276  // set ordering for PostStepDoIt
277  pmanager->SetProcessOrdering(aMultipleScattering, idxPostStep, 1);
278  pmanager->SetProcessOrdering(anIonisation, idxPostStep, 2);
279  }
280  }
281 }
282 
283 
284 #include "G4Decay.hh"
286 {
287  // Add Decay Process
288  G4Decay* theDecayProcess = new G4Decay();
289  //G4cout << "decay" <<theDecayProcess<<endl;
291  while( (*theParticleIterator)() ){
293  G4ProcessManager* pmanager = particle->GetProcessManager();
294  if (theDecayProcess->IsApplicable(*particle)) {
295  pmanager ->AddProcess(theDecayProcess);
296  // set ordering for PostStepDoIt and AtRestDoIt
297  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
298  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
299  }
300  }
301 }
302 
303 // WARNING: This methode is mandatory if U want to use GFLASH
305 {
307  theFastSimulationManagerProcess =
309  G4cout << "FastSimulationManagerProcess" <<G4endl;
311  //std::cout<<"---"<<std::endl;
312  while( (*theParticleIterator)() ){
313  //std::cout<<"+++"<<std::endl;
314 
316  // std::cout<<"--- particle "<<particle->GetParticleName()<<std::endl;
317  G4ProcessManager* pmanager = particle->GetProcessManager();
318  // The fast simulation process becomes a discrete process only since 9.0:
319  pmanager->AddDiscreteProcess(theFastSimulationManagerProcess);
320  }
321 }
322 
324 {
325  if (verboseLevel >1){
326  G4cout << "ExGflashPhysicsList::SetCuts:";
327  }
328  // " G4VUserPhysicsList::SetCutsWithDefault" method sets
329  // the default cut value for all particle types
332  // SetCutValue(100*mm, "gamma");
333 // SetCutValue(0*mm, "e-");
334 // SetCutValue(0*mm, "e+");
335 
336 // SetCutValue(62*mm, "gamma");
337 // SetCutValue(0.73*mm, "e-");
338 // SetCutValue(0.78*mm, "e+");
339 
340 
341 
342 
343 
345 // set cuts for region crystals with default Cuts
346  G4Region* region = G4RegionStore::GetInstance()->GetRegion("crystals");
347  region->SetProductionCuts(
348  G4ProductionCutsTable::GetProductionCutsTable()->GetDefaultProductionCuts());
349 }
350 
351 
352 
353 
354 
355 
356 
357 
358