Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
B03PhysicsList.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 //
26 /// \file biasing/B03/src/B03PhysicsList.cc
27 /// \brief Implementation of the B03PhysicsList class
28 //
29 //
30 // $Id: B03PhysicsList.cc 75089 2013-10-25 23:25:21Z dwright $
31 //
32 
33 #include "globals.hh"
34 #include <iomanip>
35 
36 #include "B03PhysicsList.hh"
37 
38 #include "G4ParticleDefinition.hh"
39 #include "G4ParticleWithCuts.hh"
40 #include "G4ProcessManager.hh"
41 #include "G4ProcessVector.hh"
42 #include "G4ParticleTypes.hh"
43 #include "G4ParticleTable.hh"
44 #include "G4BosonConstructor.hh"
45 #include "G4LeptonConstructor.hh"
46 #include "G4MesonConstructor.hh"
47 #include "G4BaryonConstructor.hh"
48 #include "G4IonConstructor.hh"
50 #include "G4Material.hh"
51 #include "G4MaterialTable.hh"
52 #include "G4SystemOfUnits.hh"
53 
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55 
57 {
58  fParaWorldName.clear();
59  SetVerboseLevel(1);
60 }
61 
62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
63 
65 {
66  fParaWorldName.clear();
67 }
68 
69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70 
72 {
73  // In this method, static member functions should be called
74  // for all particles which you want to use.
75  // This ensures that objects of these particle types will be
76  // created in the program.
77 
84 }
85 
86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87 
89 {
90  // Construct all bosons
91  G4BosonConstructor pConstructor;
92  pConstructor.ConstructParticle();
93 }
94 
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96 
98 {
99  // Construct all leptons
100  G4LeptonConstructor pConstructor;
101  pConstructor.ConstructParticle();
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105 
107 {
108  // Construct all mesons
109  G4MesonConstructor pConstructor;
110  pConstructor.ConstructParticle();
111 }
112 
113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
114 
116 {
117  // Construct all barions
118  G4BaryonConstructor pConstructor;
119  pConstructor.ConstructParticle();
120 }
121 
122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
123 
125 {
126  // Construct light ions
127  G4IonConstructor pConstructor;
128  pConstructor.ConstructParticle();
129 }
130 
131 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
132 
134 {
135  // Construct resonaces and quarks
136  G4ShortLivedConstructor pConstructor;
137  pConstructor.ConstructParticle();
138 }
139 
140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141 
143 {
147  ConstructEM();
149  ConstructHad();
151 }
152 
153 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
154 
155 #include "G4ComptonScattering.hh"
156 #include "G4GammaConversion.hh"
157 #include "G4PhotoElectricEffect.hh"
158 
159 #include "G4eMultipleScattering.hh"
160 #include "G4MuMultipleScattering.hh"
161 #include "G4hMultipleScattering.hh"
162 
163 #include "G4eIonisation.hh"
164 #include "G4eBremsstrahlung.hh"
165 #include "G4eplusAnnihilation.hh"
166 
167 #include "G4MuIonisation.hh"
168 #include "G4MuBremsstrahlung.hh"
169 #include "G4MuPairProduction.hh"
170 
171 #include "G4hIonisation.hh"
172 
174 {
175  theParticleIterator->reset();
176  while( (*theParticleIterator)() ){
177  G4ParticleDefinition* particle = theParticleIterator->value();
178  G4ProcessManager* pmanager = particle->GetProcessManager();
179  G4String particleName = particle->GetParticleName();
180 
181  if (particleName == "gamma") {
182  // gamma
183  // Construct processes for gamma
184  pmanager->AddDiscreteProcess(new G4GammaConversion());
185  pmanager->AddDiscreteProcess(new G4ComptonScattering());
186  pmanager->AddDiscreteProcess(new G4PhotoElectricEffect());
187 
188  } else if (particleName == "e-") {
189  //electron
190  // Construct processes for electron
191  pmanager->AddProcess(new G4eMultipleScattering(),-1,1,1);
192  pmanager->AddProcess(new G4eIonisation(),-1,2,2);
193  pmanager->AddProcess(new G4eBremsstrahlung(),-1,-1,3);
194 
195  } else if (particleName == "e+") {
196  //positron
197  // Construct processes for positron
198  pmanager->AddProcess(new G4eMultipleScattering(),-1,1,1);
199 
200  pmanager->AddProcess(new G4eIonisation(),-1,2,2);
201  pmanager->AddProcess(new G4eBremsstrahlung(),-1,-1,3);
202  pmanager->AddProcess(new G4eplusAnnihilation(),0,-1,4);
203 
204  } else if( particleName == "mu+" ||
205  particleName == "mu-" ) {
206  //muon
207  // Construct processes for muon+
208  pmanager->AddProcess(new G4MuMultipleScattering(),-1,1,1);
209  pmanager->AddProcess(new G4MuIonisation(),-1,2,2);
210  pmanager->AddProcess(new G4MuBremsstrahlung(),-1,-1,3);
211  pmanager->AddProcess(new G4MuPairProduction(),-1,-1,4);
212 
213  } else if( particleName == "GenericIon" ) {
214  pmanager->AddProcess(new G4hMultipleScattering(),-1,1,1);
215  pmanager->AddProcess(new G4hIonisation(),-1,2,2);
216  } else {
217  if ((particle->GetPDGCharge() != 0.0) &&
218  (particle->GetParticleName() != "chargedgeantino")&&
219  (!particle->IsShortLived()) ) {
220  // all others charged particles except geantino
221  pmanager->AddProcess(new G4hMultipleScattering(),-1,1,1);
222  pmanager->AddProcess(new G4hIonisation(),-1,2,2);
223  }
224  }
225  }
226 }
227 
228 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
229 
230 // Hadron Processes
231 
232 #include "G4HadronElasticProcess.hh"
233 #include "G4HadronFissionProcess.hh"
234 #include "G4HadronCaptureProcess.hh"
235 
261 
262 // Low-energy Models
263 
264 #include "G4HadronElastic.hh"
265 #include "G4NeutronRadCapture.hh"
266 #include "G4LFission.hh"
267 
268 // -- generator models
269 #include "G4TheoFSGenerator.hh"
270 #include "G4ExcitationHandler.hh"
271 #include "G4Evaporation.hh"
272 #include "G4CompetitiveFission.hh"
273 #include "G4FermiBreakUp.hh"
274 #include "G4StatMF.hh"
276 #include "G4Fancy3DNucleus.hh"
277 #include "G4CascadeInterface.hh"
278 #include "G4StringModel.hh"
279 #include "G4PreCompoundModel.hh"
280 #include "G4FTFModel.hh"
281 #include "G4QGSMFragmentation.hh"
283 #include "G4ExcitedStringDecay.hh"
284 #include "G4QMDReaction.hh"
286 
287 // Cross sections
288 #include "G4IonsShenCrossSection.hh"
289 #include "G4TripathiCrossSection.hh"
291 
292 //
293 // ConstructHad()
294 //
295 // Makes discrete physics processes for the hadrons
296 // The processes are: Elastic scattering, Inelastic scattering,
297 // Fission (for neutron only), and Capture (neutron).
298 //
299 
301 {
302  // this will be the model class for high energies
303  G4TheoFSGenerator* theTheoModel = new G4TheoFSGenerator;
304  G4TheoFSGenerator* antiBHighEnergyModel = new G4TheoFSGenerator;
305 
306  // all models for treatment of thermal nucleus
307  G4Evaporation* theEvaporation = new G4Evaporation;
308  G4FermiBreakUp* theFermiBreakUp = new G4FermiBreakUp;
309  G4StatMF* theMF = new G4StatMF;
310 
311  // Evaporation logic
312  G4ExcitationHandler* theHandler = new G4ExcitationHandler;
313  theHandler->SetEvaporation(theEvaporation);
314  theHandler->SetFermiModel(theFermiBreakUp);
315  theHandler->SetMultiFragmentation(theMF);
316  theHandler->SetMaxAandZForFermiBreakUp(12, 6);
317  theHandler->SetMinEForMultiFrag(3*MeV);
318 
319  // Pre equilibrium stage
320  G4PreCompoundModel* thePreEquilib = new G4PreCompoundModel(theHandler);
321 
322  // a no-cascade generator-precompound interaface
323  G4GeneratorPrecompoundInterface* theCascade =
325  theCascade->SetDeExcitation(thePreEquilib);
326 
327  // Bertini cascade
328  G4CascadeInterface* bertini = new G4CascadeInterface;
329  bertini->SetMaxEnergy(22*MeV);
330 
331  // here come the high energy parts
332  G4VPartonStringModel* theStringModel;
333  theStringModel = new G4FTFModel;
334  theTheoModel->SetTransport(theCascade);
335  theTheoModel->SetHighEnergyGenerator(theStringModel);
336  theTheoModel->SetMinEnergy(19*GeV);
337  theTheoModel->SetMaxEnergy(100*TeV);
338 
339  G4VLongitudinalStringDecay* theFragmentation = new G4QGSMFragmentation;
340  G4ExcitedStringDecay* theStringDecay =
341  new G4ExcitedStringDecay(theFragmentation);
342  theStringModel->SetFragmentationModel(theStringDecay);
343 
344  // high energy model for anti-baryons
345  antiBHighEnergyModel = new G4TheoFSGenerator("ANTI-FTFP");
346  G4FTFModel* antiBStringModel = new G4FTFModel;
347  G4ExcitedStringDecay* stringDecay =
349  antiBStringModel->SetFragmentationModel(stringDecay);
350 
351  G4GeneratorPrecompoundInterface* antiBCascade =
353  G4PreCompoundModel* preEquilib =
355  antiBCascade->SetDeExcitation(preEquilib);
356 
357  antiBHighEnergyModel->SetTransport(antiBCascade);
358  antiBHighEnergyModel->SetHighEnergyGenerator(antiBStringModel);
359  antiBHighEnergyModel->SetMinEnergy(0.0);
360  antiBHighEnergyModel->SetMaxEnergy(20*TeV);
361 
362  // Light ion models
364  binaryCascade->SetMinEnergy(0.0);
365  binaryCascade->SetMaxEnergy(110*MeV);
366 
367  G4QMDReaction* qmd = new G4QMDReaction;
368  qmd->SetMinEnergy(100*MeV);
369  qmd->SetMaxEnergy(10*GeV);
370 
374 
375  // Elastic process
376  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
377  G4HadronElastic* theElasticModel = new G4HadronElastic;
378  theElasticProcess->RegisterMe(theElasticModel);
379 
380  theParticleIterator->reset();
381  while ((*theParticleIterator)()) {
382  G4ParticleDefinition* particle = theParticleIterator->value();
383  G4ProcessManager* pmanager = particle->GetProcessManager();
384  G4String particleName = particle->GetParticleName();
385 
386  if (particleName == "pi+") {
387  pmanager->AddDiscreteProcess(theElasticProcess);
388  G4PionPlusInelasticProcess* theInelasticProcess =
389  new G4PionPlusInelasticProcess("inelastic");
390  theInelasticProcess->RegisterMe(bertini);
391  theInelasticProcess->RegisterMe(theTheoModel);
392  pmanager->AddDiscreteProcess(theInelasticProcess);
393  } else if (particleName == "pi-") {
394  pmanager->AddDiscreteProcess(theElasticProcess);
395  G4PionMinusInelasticProcess* theInelasticProcess =
396  new G4PionMinusInelasticProcess("inelastic");
397  theInelasticProcess->RegisterMe(bertini);
398  theInelasticProcess->RegisterMe(theTheoModel);
399  pmanager->AddDiscreteProcess(theInelasticProcess);
400  } else if (particleName == "kaon+") {
401  pmanager->AddDiscreteProcess(theElasticProcess);
402  G4KaonPlusInelasticProcess* theInelasticProcess =
403  new G4KaonPlusInelasticProcess("inelastic");
404  theInelasticProcess->RegisterMe(bertini);
405  theInelasticProcess->RegisterMe(theTheoModel);
406  pmanager->AddDiscreteProcess(theInelasticProcess);
407  }
408  else if (particleName == "kaon0S") {
409  pmanager->AddDiscreteProcess(theElasticProcess);
410  G4KaonZeroSInelasticProcess* theInelasticProcess =
411  new G4KaonZeroSInelasticProcess("inelastic");
412  theInelasticProcess->RegisterMe(bertini);
413  theInelasticProcess->RegisterMe(theTheoModel);
414  pmanager->AddDiscreteProcess(theInelasticProcess);
415  }
416  else if (particleName == "kaon0L") {
417  pmanager->AddDiscreteProcess(theElasticProcess);
418  G4KaonZeroLInelasticProcess* theInelasticProcess =
419  new G4KaonZeroLInelasticProcess("inelastic");
420  theInelasticProcess->RegisterMe(bertini);
421  theInelasticProcess->RegisterMe(theTheoModel);
422  pmanager->AddDiscreteProcess(theInelasticProcess);
423  }
424  else if (particleName == "kaon-") {
425  pmanager->AddDiscreteProcess(theElasticProcess);
426  G4KaonMinusInelasticProcess* theInelasticProcess =
427  new G4KaonMinusInelasticProcess("inelastic");
428  theInelasticProcess->RegisterMe(bertini);
429  theInelasticProcess->RegisterMe(theTheoModel);
430  pmanager->AddDiscreteProcess(theInelasticProcess);
431  }
432  else if (particleName == "proton") {
433  pmanager->AddDiscreteProcess(theElasticProcess);
434  G4ProtonInelasticProcess* theInelasticProcess =
435  new G4ProtonInelasticProcess("inelastic");
436  theInelasticProcess->RegisterMe(bertini);
437  theInelasticProcess->RegisterMe(theTheoModel);
438  pmanager->AddDiscreteProcess(theInelasticProcess);
439  }
440  else if (particleName == "anti_proton") {
441  pmanager->AddDiscreteProcess(theElasticProcess);
442  G4AntiProtonInelasticProcess* theInelasticProcess =
443  new G4AntiProtonInelasticProcess("inelastic");
444  theInelasticProcess->RegisterMe(antiBHighEnergyModel);
445  pmanager->AddDiscreteProcess(theInelasticProcess);
446 
447  } else if (particleName == "neutron") {
448  // elastic scattering
449  pmanager->AddDiscreteProcess(theElasticProcess);
450 
451  // inelastic scattering
452  G4NeutronInelasticProcess* theInelasticProcess =
453  new G4NeutronInelasticProcess("inelastic");
454  theInelasticProcess->RegisterMe(bertini);
455  theInelasticProcess->RegisterMe(theTheoModel);
456  pmanager->AddDiscreteProcess(theInelasticProcess);
457 
458  // fission
459  G4HadronFissionProcess* theFissionProcess = new G4HadronFissionProcess;
460  G4LFission* theFissionModel = new G4LFission;
461  theFissionProcess->RegisterMe(theFissionModel);
462  pmanager->AddDiscreteProcess(theFissionProcess);
463 
464  // capture
465  G4HadronCaptureProcess* theCaptureProcess = new G4HadronCaptureProcess;
466  G4NeutronRadCapture* theCaptureModel = new G4NeutronRadCapture;
467  theCaptureProcess->RegisterMe(theCaptureModel);
468  pmanager->AddDiscreteProcess(theCaptureProcess);
469 
470  } else if (particleName == "anti_neutron") {
471  pmanager->AddDiscreteProcess(theElasticProcess);
472  G4AntiNeutronInelasticProcess* theInelasticProcess =
473  new G4AntiNeutronInelasticProcess("inelastic");
474  theInelasticProcess->RegisterMe(antiBHighEnergyModel);
475  pmanager->AddDiscreteProcess(theInelasticProcess);
476 
477  } else if (particleName == "lambda") {
478  pmanager->AddDiscreteProcess(theElasticProcess);
479  G4LambdaInelasticProcess* theInelasticProcess =
480  new G4LambdaInelasticProcess("inelastic");
481  theInelasticProcess->RegisterMe(bertini);
482  theInelasticProcess->RegisterMe(theTheoModel);
483  pmanager->AddDiscreteProcess(theInelasticProcess);
484  }
485  else if (particleName == "anti_lambda") {
486  pmanager->AddDiscreteProcess(theElasticProcess);
487  G4AntiLambdaInelasticProcess* theInelasticProcess =
488  new G4AntiLambdaInelasticProcess("inelastic");
489  theInelasticProcess->RegisterMe(antiBHighEnergyModel);
490  pmanager->AddDiscreteProcess(theInelasticProcess);
491  }
492  else if (particleName == "sigma+") {
493  pmanager->AddDiscreteProcess(theElasticProcess);
494  G4SigmaPlusInelasticProcess* theInelasticProcess =
495  new G4SigmaPlusInelasticProcess("inelastic");
496  theInelasticProcess->RegisterMe(bertini);
497  theInelasticProcess->RegisterMe(theTheoModel);
498  pmanager->AddDiscreteProcess(theInelasticProcess);
499  }
500  else if (particleName == "sigma-") {
501  pmanager->AddDiscreteProcess(theElasticProcess);
502  G4SigmaMinusInelasticProcess* theInelasticProcess =
503  new G4SigmaMinusInelasticProcess("inelastic");
504  theInelasticProcess->RegisterMe(bertini);
505  theInelasticProcess->RegisterMe(theTheoModel);
506  pmanager->AddDiscreteProcess(theInelasticProcess);
507  }
508  else if (particleName == "anti_sigma+") {
509  pmanager->AddDiscreteProcess(theElasticProcess);
510  G4AntiSigmaPlusInelasticProcess* theInelasticProcess =
511  new G4AntiSigmaPlusInelasticProcess("inelastic");
512  theInelasticProcess->RegisterMe(antiBHighEnergyModel);
513  pmanager->AddDiscreteProcess(theInelasticProcess);
514  }
515  else if (particleName == "anti_sigma-") {
516  pmanager->AddDiscreteProcess(theElasticProcess);
517  G4AntiSigmaMinusInelasticProcess* theInelasticProcess =
518  new G4AntiSigmaMinusInelasticProcess("inelastic");
519  theInelasticProcess->RegisterMe(antiBHighEnergyModel);
520  pmanager->AddDiscreteProcess(theInelasticProcess);
521  }
522  else if (particleName == "xi0") {
523  pmanager->AddDiscreteProcess(theElasticProcess);
524  G4XiZeroInelasticProcess* theInelasticProcess =
525  new G4XiZeroInelasticProcess("inelastic");
526  theInelasticProcess->RegisterMe(bertini);
527  theInelasticProcess->RegisterMe(theTheoModel);
528  pmanager->AddDiscreteProcess(theInelasticProcess);
529  }
530  else if (particleName == "xi-") {
531  pmanager->AddDiscreteProcess(theElasticProcess);
532  G4XiMinusInelasticProcess* theInelasticProcess =
533  new G4XiMinusInelasticProcess("inelastic");
534  theInelasticProcess->RegisterMe(bertini);
535  theInelasticProcess->RegisterMe(theTheoModel);
536  pmanager->AddDiscreteProcess(theInelasticProcess);
537  }
538  else if (particleName == "anti_xi0") {
539  pmanager->AddDiscreteProcess(theElasticProcess);
540  G4AntiXiZeroInelasticProcess* theInelasticProcess =
541  new G4AntiXiZeroInelasticProcess("inelastic");
542  theInelasticProcess->RegisterMe(antiBHighEnergyModel);
543  pmanager->AddDiscreteProcess(theInelasticProcess);
544  }
545  else if (particleName == "anti_xi-") {
546  pmanager->AddDiscreteProcess(theElasticProcess);
547  G4AntiXiMinusInelasticProcess* theInelasticProcess =
548  new G4AntiXiMinusInelasticProcess("inelastic");
549  theInelasticProcess->RegisterMe(antiBHighEnergyModel);
550  pmanager->AddDiscreteProcess(theInelasticProcess);
551  }
552  else if (particleName == "deuteron") {
553  pmanager->AddDiscreteProcess(theElasticProcess);
554  G4DeuteronInelasticProcess* theInelasticProcess =
555  new G4DeuteronInelasticProcess("inelastic");
556  theInelasticProcess->RegisterMe(binaryCascade);
557  theInelasticProcess->RegisterMe(qmd);
558  theInelasticProcess->AddDataSet(shenXS);
559  theInelasticProcess->AddDataSet(tripXS);
560  theInelasticProcess->AddDataSet(tripLightXS);
561  pmanager->AddDiscreteProcess(theInelasticProcess);
562  }
563  else if (particleName == "triton") {
564  pmanager->AddDiscreteProcess(theElasticProcess);
565  G4TritonInelasticProcess* theInelasticProcess =
566  new G4TritonInelasticProcess("inelastic");
567  theInelasticProcess->RegisterMe(binaryCascade);
568  theInelasticProcess->RegisterMe(qmd);
569  theInelasticProcess->AddDataSet(shenXS);
570  theInelasticProcess->AddDataSet(tripXS);
571  theInelasticProcess->AddDataSet(tripLightXS);
572  pmanager->AddDiscreteProcess(theInelasticProcess);
573  }
574  else if (particleName == "alpha") {
575  pmanager->AddDiscreteProcess(theElasticProcess);
576  G4AlphaInelasticProcess* theInelasticProcess =
577  new G4AlphaInelasticProcess("inelastic");
578  theInelasticProcess->RegisterMe(binaryCascade);
579  theInelasticProcess->RegisterMe(qmd);
580  theInelasticProcess->AddDataSet(shenXS);
581  theInelasticProcess->AddDataSet(tripXS);
582  theInelasticProcess->AddDataSet(tripLightXS);
583  pmanager->AddDiscreteProcess(theInelasticProcess);
584 
585  } else if (particleName == "omega-") {
586  pmanager->AddDiscreteProcess(theElasticProcess);
587  G4OmegaMinusInelasticProcess* theInelasticProcess =
588  new G4OmegaMinusInelasticProcess("inelastic");
589  theInelasticProcess->RegisterMe(bertini);
590  theInelasticProcess->RegisterMe(theTheoModel);
591  pmanager->AddDiscreteProcess(theInelasticProcess);
592 
593  } else if (particleName == "anti_omega-") {
594  pmanager->AddDiscreteProcess(theElasticProcess);
595  G4AntiOmegaMinusInelasticProcess* theInelasticProcess =
596  new G4AntiOmegaMinusInelasticProcess("inelastic");
597  theInelasticProcess->RegisterMe(antiBHighEnergyModel);
598  pmanager->AddDiscreteProcess(theInelasticProcess);
599  }
600  }
601 }
602 
603 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
604 
606 {;}
607 
608 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
609 
610 #include "G4Decay.hh"
612 {
613  G4Decay* theDecayProcess = new G4Decay();
614  theParticleIterator->reset();
615  while( (*theParticleIterator)() ){
616  G4ParticleDefinition* particle = theParticleIterator->value();
617  G4ProcessManager* pmanager = particle->GetProcessManager();
618  if (theDecayProcess->IsApplicable(*particle)) {
619  pmanager ->AddProcess(theDecayProcess);
620  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
621  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
622  }
623  }
624 }
625 
626 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
627 
629 {
630  if (verboseLevel >0)
631  {
632  G4cout << "B03PhysicsList::SetCuts:";
633  G4cout << "CutLength : " << defaultCutValue/mm << " (mm)" << G4endl;
634  }
635  // "G4VUserPhysicsList::SetCutsWithDefault" method sets
636  // the default cut value for all particle types
638 }
639 
640 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
641 
644 
645  G4int npw = fParaWorldName.size();
646  for ( G4int i = 0; i < npw; i++){
647  G4ParallelWorldScoringProcess* theParallelWorldScoringProcess
648  = new G4ParallelWorldScoringProcess("ParaWorldScoringProc");
649  theParallelWorldScoringProcess->SetParallelWorld(fParaWorldName[i]);
650 
651  theParticleIterator->reset();
652  while( (*theParticleIterator)() ){
653  G4ParticleDefinition* particle = theParticleIterator->value();
654  if ( !particle->IsShortLived() ){
655  G4ProcessManager* pmanager = particle->GetProcessManager();
656  pmanager->AddProcess(theParallelWorldScoringProcess);
657  pmanager->SetProcessOrderingToLast(theParallelWorldScoringProcess
658  ,idxAtRest);
659  pmanager->SetProcessOrdering(theParallelWorldScoringProcess
660  ,idxAlongStep,1);
661  pmanager->SetProcessOrderingToLast(theParallelWorldScoringProcess
662  ,idxPostStep);
663  }
664  }
665  }
666 
667 }
668 
669 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
670 
671 #include "G4ImportanceProcess.hh"
672 #include "G4IStore.hh"
674 
675 
676  G4cout << " Preparing Importance Sampling with GhostWorld "
677  << fBiasWorldName << G4endl;
678  fGeomSampler->SetParallel(true); // parallelworld
679  G4IStore* iStore = G4IStore::GetInstance(fBiasWorldName);
680  fGeomSampler->SetWorld(iStore->GetParallelWorldVolumePointer());
681  // fGeomSampler->PrepareImportanceSampling(G4IStore::
682  // GetInstance(fBiasWorldName), 0);
683  static G4bool first = true;
684  if(first) {
685  fGeomSampler->PrepareImportanceSampling(iStore, 0);
686 
687  fGeomSampler->Configure();
688  G4cout << " GeomSampler Configured!!! " << G4endl;
689  first = false;
690  }
691 
692 #ifdef G4MULTITHREADED
693  fGeomSampler->AddProcess();
694 #else
695  G4cout << " Running in singlethreaded mode!!! " << G4endl;
696 #endif
697 
698 }
699 
700 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void SetWorld(const G4VPhysicalVolume *world)
void ConstructAllBosons()
void AddScoringProcess()
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4Decay.cc:89
void SetFragmentationModel(G4VStringFragmentation *aModel)
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
Definition of the B03PhysicsList class.
void SetMinEForMultiFrag(G4double anE)
void SetParallel(G4bool paraflag)
static void ConstructParticle()
virtual void SetCuts()
static void ConstructParticle()
G4ProcessManager * GetProcessManager() const
int G4int
Definition: G4Types.hh:78
virtual const G4VPhysicalVolume * GetParallelWorldVolumePointer() const
Definition: G4IStore.cc:94
const G4String & GetParticleName() const
void SetHighEnergyGenerator(G4VHighEnergyGenerator *const value)
static void ConstructParticle()
void RegisterMe(G4HadronicInteraction *a)
void AddBiasingProcess()
void SetMinEnergy(G4double anEnergy)
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
virtual void ConstructGeneral()
G4GLOB_DLL std::ostream G4cout
static void ConstructParticle()
virtual void ConstructProcess()
virtual void ConstructLeptHad()
bool G4bool
Definition: G4Types.hh:79
virtual void AddProcess()
void SetFermiModel(G4VFermiBreakUp *ptr)
void SetVerboseLevel(G4int value)
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void SetParallelWorld(G4String parallelWorldName)
void SetMultiFragmentation(G4VMultiFragmentation *ptr)
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
void ConstructAllLeptons()
static void ConstructParticle()
virtual void ConstructHad()
void ConstructAllShortLiveds()
void ConstructAllBaryons()
static G4IStore * GetInstance()
Definition: G4IStore.cc:211
virtual void ConstructEM()
void SetProcessOrderingToLast(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt)
void ConstructAllMesons()
void SetMaxAandZForFermiBreakUp(G4int anA, G4int aZ)
virtual void ConstructParticle()
void SetEvaporation(G4VEvaporation *ptr)
void SetMaxEnergy(const G4double anEnergy)
void SetDeExcitation(G4VPreCompoundModel *ptr)
#define G4endl
Definition: G4ios.hh:61
void SetTransport(G4VIntraNuclearTransportModel *const value)
virtual void Configure()
G4double GetPDGCharge() const
#define theParticleIterator
virtual void PrepareImportanceSampling(G4VIStore *istore, const G4VImportanceAlgorithm *ialg)
virtual ~B03PhysicsList()