Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
B01PhysicsList.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 "globals.hh"
34 #include <iomanip>
35 
36 #include "B01PhysicsList.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 
55 {
56  paraWorldName.clear();
57  SetVerboseLevel(1);
58 
59 }
60 
62 {
63  paraWorldName.clear();
64 }
65 
67 {
68  // In this method, static member functions should be called
69  // for all particles which you want to use.
70  // This ensures that objects of these particle types will be
71  // created in the program.
72 
79 }
80 
82 {
83  // Construct all bosons
84  G4BosonConstructor pConstructor;
85  pConstructor.ConstructParticle();
86 }
87 
89 {
90  // Construct all leptons
91  G4LeptonConstructor pConstructor;
92  pConstructor.ConstructParticle();
93 }
94 
96 {
97  // Construct all mesons
98  G4MesonConstructor pConstructor;
99  pConstructor.ConstructParticle();
100 }
101 
103 {
104  // Construct all barions
105  G4BaryonConstructor pConstructor;
106  pConstructor.ConstructParticle();
107 }
108 
110 {
111  // Construct light ions
112  G4IonConstructor pConstructor;
113  pConstructor.ConstructParticle();
114 }
115 
117 {
118  // Construct resonaces and quarks
119  G4ShortLivedConstructor pConstructor;
120  pConstructor.ConstructParticle();
121 }
122 
124 {
127  ConstructEM();
129  ConstructHad();
131 }
132 
133 #include "G4ComptonScattering.hh"
134 #include "G4GammaConversion.hh"
135 #include "G4PhotoElectricEffect.hh"
136 
137 #include "G4eMultipleScattering.hh"
138 #include "G4MuMultipleScattering.hh"
139 #include "G4hMultipleScattering.hh"
140 
141 #include "G4eIonisation.hh"
142 #include "G4eBremsstrahlung.hh"
143 #include "G4eplusAnnihilation.hh"
144 
145 #include "G4MuIonisation.hh"
146 #include "G4MuBremsstrahlung.hh"
147 #include "G4MuPairProduction.hh"
148 
149 #include "G4hIonisation.hh"
150 
152 {
154  while( (*theParticleIterator)() ){
156  G4ProcessManager* pmanager = particle->GetProcessManager();
157  G4String particleName = particle->GetParticleName();
158 
159  if (particleName == "gamma") {
160  // gamma
161  // Construct processes for gamma
162  pmanager->AddDiscreteProcess(new G4GammaConversion());
163  pmanager->AddDiscreteProcess(new G4ComptonScattering());
164  pmanager->AddDiscreteProcess(new G4PhotoElectricEffect());
165 
166  } else if (particleName == "e-") {
167  //electron
168  // Construct processes for electron
169  pmanager->AddProcess(new G4eMultipleScattering(),-1,1,1);
170  pmanager->AddProcess(new G4eIonisation(),-1,2,2);
171  pmanager->AddProcess(new G4eBremsstrahlung(),-1,-1,3);
172 
173  } else if (particleName == "e+") {
174  //positron
175  // Construct processes for positron
176  pmanager->AddProcess(new G4eMultipleScattering(),-1,1,1);
177 
178  pmanager->AddProcess(new G4eIonisation(),-1,2,2);
179  pmanager->AddProcess(new G4eBremsstrahlung(),-1,-1,3);
180  pmanager->AddProcess(new G4eplusAnnihilation(),0,-1,4);
181 
182  } else if( particleName == "mu+" ||
183  particleName == "mu-" ) {
184  //muon
185  // Construct processes for muon+
186  pmanager->AddProcess(new G4MuMultipleScattering(),-1,1,1);
187  pmanager->AddProcess(new G4MuIonisation(),-1,2,2);
188  pmanager->AddProcess(new G4MuBremsstrahlung(),-1,-1,3);
189  pmanager->AddProcess(new G4MuPairProduction(),-1,-1,4);
190 
191  } else if( particleName == "GenericIon" ) {
192  pmanager->AddProcess(new G4hMultipleScattering(),-1,1,1);
193  pmanager->AddProcess(new G4hIonisation(),-1,2,2);
194  } else {
195  if ((particle->GetPDGCharge() != 0.0) &&
196  (particle->GetParticleName() != "chargedgeantino")&&
197  (!particle->IsShortLived()) ) {
198  // all others charged particles except geantino
199  pmanager->AddProcess(new G4hMultipleScattering(),-1,1,1);
200  pmanager->AddProcess(new G4hIonisation(),-1,2,2);
201  }
202  }
203  }
204 }
205 
206 // Hadron Processes
207 
208 #include "G4HadronElasticProcess.hh"
209 #include "G4HadronFissionProcess.hh"
210 #include "G4HadronCaptureProcess.hh"
211 
237 
238 // Low-energy Models
239 
240 #include "G4LElastic.hh"
241 #include "G4LFission.hh"
242 #include "G4LCapture.hh"
243 
244 #include "G4LEPionPlusInelastic.hh"
245 #include "G4LEPionMinusInelastic.hh"
246 #include "G4LEKaonPlusInelastic.hh"
247 #include "G4LEKaonZeroSInelastic.hh"
248 #include "G4LEKaonZeroLInelastic.hh"
249 #include "G4LEKaonMinusInelastic.hh"
250 #include "G4LEProtonInelastic.hh"
252 #include "G4LENeutronInelastic.hh"
254 #include "G4LELambdaInelastic.hh"
256 #include "G4LESigmaPlusInelastic.hh"
260 #include "G4LEXiZeroInelastic.hh"
261 #include "G4LEXiMinusInelastic.hh"
264 #include "G4LEDeuteronInelastic.hh"
265 #include "G4LETritonInelastic.hh"
266 #include "G4LEAlphaInelastic.hh"
269 
270 // -- generator models
271 #include "G4TheoFSGenerator.hh"
272 #include "G4ExcitationHandler.hh"
273 #include "G4Evaporation.hh"
274 #include "G4CompetitiveFission.hh"
275 #include "G4FermiBreakUp.hh"
276 #include "G4StatMF.hh"
278 #include "G4Fancy3DNucleus.hh"
279 #include "G4LEProtonInelastic.hh"
280 #include "G4StringModel.hh"
281 #include "G4PreCompoundModel.hh"
282 #include "G4FTFModel.hh"
283 #include "G4QGSMFragmentation.hh"
284 #include "G4ExcitedStringDecay.hh"
285 
286 //
287 // ConstructHad()
288 //
289 // Makes discrete physics processes for the hadrons, at present limited
290 // to those particles with GHEISHA interactions (INTRC > 0).
291 // The processes are: Elastic scattering, Inelastic scattering,
292 // Fission (for neutron only), and Capture (neutron).
293 //
294 // F.W.Jones 06-JUL-1998
295 //
296 
298 {
299  // this will be the model class for high energies
300  G4TheoFSGenerator * theTheoModel = new G4TheoFSGenerator;
301 
302  // all models for treatment of thermal nucleus
303  G4Evaporation * theEvaporation = new G4Evaporation;
304  G4FermiBreakUp * theFermiBreakUp = new G4FermiBreakUp;
305  G4StatMF * theMF = new G4StatMF;
306 
307  // Evaporation logic
308  G4ExcitationHandler * theHandler = new G4ExcitationHandler;
309  theHandler->SetEvaporation(theEvaporation);
310  theHandler->SetFermiModel(theFermiBreakUp);
311  theHandler->SetMultiFragmentation(theMF);
312  theHandler->SetMaxAandZForFermiBreakUp(12, 6);
313  theHandler->SetMinEForMultiFrag(3*MeV);
314 
315  // Pre equilibrium stage
316  G4PreCompoundModel * thePreEquilib = new G4PreCompoundModel(theHandler);
317 
318 
319  // a no-cascade generator-precompound interaface
321  theCascade->SetDeExcitation(thePreEquilib);
322 
323  // here come the high energy parts
324  // the string model; still not quite according to design - Explicite use of the forseen interfaces
325  // will be tested and documented in this program by beta-02 at latest.
326  G4VPartonStringModel * theStringModel;
327  theStringModel = new G4FTFModel;
328  theTheoModel->SetTransport(theCascade);
329  theTheoModel->SetHighEnergyGenerator(theStringModel);
330  theTheoModel->SetMinEnergy(19*GeV);
331  theTheoModel->SetMaxEnergy(100*TeV);
332 
333  G4VLongitudinalStringDecay * theFragmentation = new G4QGSMFragmentation;
334  G4ExcitedStringDecay * theStringDecay = new G4ExcitedStringDecay(theFragmentation);
335  theStringModel->SetFragmentationModel(theStringDecay);
336 
337 // done with the generator model (most of the above is also available as default)
338  G4HadronElasticProcess* theElasticProcess =
340  G4LElastic* theElasticModel = new G4LElastic;
341  theElasticProcess->RegisterMe(theElasticModel);
342  G4HadronElasticProcess* theElasticProcess1 =
345  while ((*theParticleIterator)()) {
347  G4ProcessManager* pmanager = particle->GetProcessManager();
348  G4String particleName = particle->GetParticleName();
349 
350  if (particleName == "pi+") {
351  pmanager->AddDiscreteProcess(theElasticProcess);
352  G4PionPlusInelasticProcess* theInelasticProcess =
353  new G4PionPlusInelasticProcess("inelastic");
354  G4LEPionPlusInelastic* theInelasticModel =
356  theInelasticProcess->RegisterMe(theInelasticModel);
357  theInelasticProcess->RegisterMe(theTheoModel);
358  pmanager->AddDiscreteProcess(theInelasticProcess);
359  }
360  else if (particleName == "pi-") {
361  pmanager->AddDiscreteProcess(theElasticProcess);
362  G4PionMinusInelasticProcess* theInelasticProcess =
363  new G4PionMinusInelasticProcess("inelastic");
364  G4LEPionMinusInelastic* theInelasticModel =
366  theInelasticProcess->RegisterMe(theInelasticModel);
367  theInelasticProcess->RegisterMe(theTheoModel);
368  pmanager->AddDiscreteProcess(theInelasticProcess);
369  }
370  else if (particleName == "kaon+") {
371  pmanager->AddDiscreteProcess(theElasticProcess);
372  G4KaonPlusInelasticProcess* theInelasticProcess =
373  new G4KaonPlusInelasticProcess("inelastic");
374  G4LEKaonPlusInelastic* theInelasticModel = new G4LEKaonPlusInelastic;
375  theInelasticProcess->RegisterMe(theInelasticModel);
376  theInelasticProcess->RegisterMe(theTheoModel);
377  pmanager->AddDiscreteProcess(theInelasticProcess);
378  }
379  else if (particleName == "kaon0S") {
380  pmanager->AddDiscreteProcess(theElasticProcess);
381  G4KaonZeroSInelasticProcess* theInelasticProcess =
382  new G4KaonZeroSInelasticProcess("inelastic");
383  G4LEKaonZeroSInelastic* theInelasticModel =
385  theInelasticProcess->RegisterMe(theInelasticModel);
386  theInelasticProcess->RegisterMe(theTheoModel);
387  pmanager->AddDiscreteProcess(theInelasticProcess);
388  }
389  else if (particleName == "kaon0L") {
390  pmanager->AddDiscreteProcess(theElasticProcess);
391  G4KaonZeroLInelasticProcess* theInelasticProcess =
392  new G4KaonZeroLInelasticProcess("inelastic");
393  G4LEKaonZeroLInelastic* theInelasticModel =
395  theInelasticProcess->RegisterMe(theInelasticModel);
396  theInelasticProcess->RegisterMe(theTheoModel);
397  pmanager->AddDiscreteProcess(theInelasticProcess);
398  }
399  else if (particleName == "kaon-") {
400  pmanager->AddDiscreteProcess(theElasticProcess);
401  G4KaonMinusInelasticProcess* theInelasticProcess =
402  new G4KaonMinusInelasticProcess("inelastic");
403  G4LEKaonMinusInelastic* theInelasticModel =
405  theInelasticProcess->RegisterMe(theInelasticModel);
406  theInelasticProcess->RegisterMe(theTheoModel);
407  pmanager->AddDiscreteProcess(theInelasticProcess);
408  }
409  else if (particleName == "proton") {
410  pmanager->AddDiscreteProcess(theElasticProcess);
411  G4ProtonInelasticProcess* theInelasticProcess =
412  new G4ProtonInelasticProcess("inelastic");
413  G4LEProtonInelastic* theInelasticModel = new G4LEProtonInelastic;
414  theInelasticProcess->RegisterMe(theInelasticModel);
415  theInelasticProcess->RegisterMe(theTheoModel);
416  pmanager->AddDiscreteProcess(theInelasticProcess);
417  }
418  else if (particleName == "anti_proton") {
419  pmanager->AddDiscreteProcess(theElasticProcess);
420  G4AntiProtonInelasticProcess* theInelasticProcess =
421  new G4AntiProtonInelasticProcess("inelastic");
422  G4LEAntiProtonInelastic* theInelasticModel =
424  theInelasticProcess->RegisterMe(theInelasticModel);
425  theInelasticProcess->RegisterMe(theTheoModel);
426  pmanager->AddDiscreteProcess(theInelasticProcess);
427  }
428  else if (particleName == "neutron") {
429 
430  // elastic scattering
431  G4LElastic* theElasticModel1 = new G4LElastic;
432  theElasticProcess1->RegisterMe(theElasticModel1);
433  pmanager->AddDiscreteProcess(theElasticProcess1);
434  // inelastic scattering
435  G4NeutronInelasticProcess* theInelasticProcess =
436  new G4NeutronInelasticProcess("inelastic");
437  G4LENeutronInelastic* theInelasticModel = new G4LENeutronInelastic;
438  theInelasticProcess->RegisterMe(theInelasticModel);
439  theInelasticProcess->RegisterMe(theTheoModel);
440  pmanager->AddDiscreteProcess(theInelasticProcess);
441  // fission
442  G4HadronFissionProcess* theFissionProcess =
444  G4LFission* theFissionModel = new G4LFission;
445  theFissionProcess->RegisterMe(theFissionModel);
446  pmanager->AddDiscreteProcess(theFissionProcess);
447  // capture
448  G4HadronCaptureProcess* theCaptureProcess =
450  G4LCapture* theCaptureModel = new G4LCapture;
451  theCaptureProcess->RegisterMe(theCaptureModel);
452  pmanager->AddDiscreteProcess(theCaptureProcess);
453  }
454  else if (particleName == "anti_neutron") {
455  pmanager->AddDiscreteProcess(theElasticProcess);
456  G4AntiNeutronInelasticProcess* theInelasticProcess =
457  new G4AntiNeutronInelasticProcess("inelastic");
458  G4LEAntiNeutronInelastic* theInelasticModel =
460  theInelasticProcess->RegisterMe(theInelasticModel);
461  theInelasticProcess->RegisterMe(theTheoModel);
462  pmanager->AddDiscreteProcess(theInelasticProcess);
463  }
464  else if (particleName == "lambda") {
465  pmanager->AddDiscreteProcess(theElasticProcess);
466  G4LambdaInelasticProcess* theInelasticProcess =
467  new G4LambdaInelasticProcess("inelastic");
468  G4LELambdaInelastic* theInelasticModel = new G4LELambdaInelastic;
469  theInelasticProcess->RegisterMe(theInelasticModel);
470  theInelasticProcess->RegisterMe(theTheoModel);
471  pmanager->AddDiscreteProcess(theInelasticProcess);
472  }
473  else if (particleName == "anti_lambda") {
474  pmanager->AddDiscreteProcess(theElasticProcess);
475  G4AntiLambdaInelasticProcess* theInelasticProcess =
476  new G4AntiLambdaInelasticProcess("inelastic");
477  G4LEAntiLambdaInelastic* theInelasticModel =
479  theInelasticProcess->RegisterMe(theInelasticModel);
480  theInelasticProcess->RegisterMe(theTheoModel);
481  pmanager->AddDiscreteProcess(theInelasticProcess);
482  }
483  else if (particleName == "sigma+") {
484  pmanager->AddDiscreteProcess(theElasticProcess);
485  G4SigmaPlusInelasticProcess* theInelasticProcess =
486  new G4SigmaPlusInelasticProcess("inelastic");
487  G4LESigmaPlusInelastic* theInelasticModel =
489  theInelasticProcess->RegisterMe(theInelasticModel);
490  theInelasticProcess->RegisterMe(theTheoModel);
491  pmanager->AddDiscreteProcess(theInelasticProcess);
492  }
493  else if (particleName == "sigma-") {
494  pmanager->AddDiscreteProcess(theElasticProcess);
495  G4SigmaMinusInelasticProcess* theInelasticProcess =
496  new G4SigmaMinusInelasticProcess("inelastic");
497  G4LESigmaMinusInelastic* theInelasticModel =
499  theInelasticProcess->RegisterMe(theInelasticModel);
500  theInelasticProcess->RegisterMe(theTheoModel);
501  pmanager->AddDiscreteProcess(theInelasticProcess);
502  }
503  else if (particleName == "anti_sigma+") {
504  pmanager->AddDiscreteProcess(theElasticProcess);
505  G4AntiSigmaPlusInelasticProcess* theInelasticProcess =
506  new G4AntiSigmaPlusInelasticProcess("inelastic");
507  G4LEAntiSigmaPlusInelastic* theInelasticModel =
509  theInelasticProcess->RegisterMe(theInelasticModel);
510  theInelasticProcess->RegisterMe(theTheoModel);
511  pmanager->AddDiscreteProcess(theInelasticProcess);
512  }
513  else if (particleName == "anti_sigma-") {
514  pmanager->AddDiscreteProcess(theElasticProcess);
515  G4AntiSigmaMinusInelasticProcess* theInelasticProcess =
516  new G4AntiSigmaMinusInelasticProcess("inelastic");
517  G4LEAntiSigmaMinusInelastic* theInelasticModel =
519  theInelasticProcess->RegisterMe(theInelasticModel);
520  theInelasticProcess->RegisterMe(theTheoModel);
521  pmanager->AddDiscreteProcess(theInelasticProcess);
522  }
523  else if (particleName == "xi0") {
524  pmanager->AddDiscreteProcess(theElasticProcess);
525  G4XiZeroInelasticProcess* theInelasticProcess =
526  new G4XiZeroInelasticProcess("inelastic");
527  G4LEXiZeroInelastic* theInelasticModel =
529  theInelasticProcess->RegisterMe(theInelasticModel);
530  theInelasticProcess->RegisterMe(theTheoModel);
531  pmanager->AddDiscreteProcess(theInelasticProcess);
532  }
533  else if (particleName == "xi-") {
534  pmanager->AddDiscreteProcess(theElasticProcess);
535  G4XiMinusInelasticProcess* theInelasticProcess =
536  new G4XiMinusInelasticProcess("inelastic");
537  G4LEXiMinusInelastic* theInelasticModel =
539  theInelasticProcess->RegisterMe(theInelasticModel);
540  theInelasticProcess->RegisterMe(theTheoModel);
541  pmanager->AddDiscreteProcess(theInelasticProcess);
542  }
543  else if (particleName == "anti_xi0") {
544  pmanager->AddDiscreteProcess(theElasticProcess);
545  G4AntiXiZeroInelasticProcess* theInelasticProcess =
546  new G4AntiXiZeroInelasticProcess("inelastic");
547  G4LEAntiXiZeroInelastic* theInelasticModel =
549  theInelasticProcess->RegisterMe(theInelasticModel);
550  theInelasticProcess->RegisterMe(theTheoModel);
551  pmanager->AddDiscreteProcess(theInelasticProcess);
552  }
553  else if (particleName == "anti_xi-") {
554  pmanager->AddDiscreteProcess(theElasticProcess);
555  G4AntiXiMinusInelasticProcess* theInelasticProcess =
556  new G4AntiXiMinusInelasticProcess("inelastic");
557  G4LEAntiXiMinusInelastic* theInelasticModel =
559  theInelasticProcess->RegisterMe(theInelasticModel);
560  theInelasticProcess->RegisterMe(theTheoModel);
561  pmanager->AddDiscreteProcess(theInelasticProcess);
562  }
563  else if (particleName == "deuteron") {
564  pmanager->AddDiscreteProcess(theElasticProcess);
565  G4DeuteronInelasticProcess* theInelasticProcess =
566  new G4DeuteronInelasticProcess("inelastic");
567  G4LEDeuteronInelastic* theInelasticModel =
569  theInelasticProcess->RegisterMe(theInelasticModel);
570  theInelasticProcess->RegisterMe(theTheoModel);
571  pmanager->AddDiscreteProcess(theInelasticProcess);
572  }
573  else if (particleName == "triton") {
574  pmanager->AddDiscreteProcess(theElasticProcess);
575  G4TritonInelasticProcess* theInelasticProcess =
576  new G4TritonInelasticProcess("inelastic");
577  G4LETritonInelastic* theInelasticModel =
579  theInelasticProcess->RegisterMe(theInelasticModel);
580  theInelasticProcess->RegisterMe(theTheoModel);
581  pmanager->AddDiscreteProcess(theInelasticProcess);
582  }
583  else if (particleName == "alpha") {
584  pmanager->AddDiscreteProcess(theElasticProcess);
585  G4AlphaInelasticProcess* theInelasticProcess =
586  new G4AlphaInelasticProcess("inelastic");
587  G4LEAlphaInelastic* theInelasticModel =
588  new G4LEAlphaInelastic;
589  theInelasticProcess->RegisterMe(theInelasticModel);
590  theInelasticProcess->RegisterMe(theTheoModel);
591  pmanager->AddDiscreteProcess(theInelasticProcess);
592  }
593  else if (particleName == "omega-") {
594  pmanager->AddDiscreteProcess(theElasticProcess);
595  G4OmegaMinusInelasticProcess* theInelasticProcess =
596  new G4OmegaMinusInelasticProcess("inelastic");
597  G4LEOmegaMinusInelastic* theInelasticModel =
599  theInelasticProcess->RegisterMe(theInelasticModel);
600  theInelasticProcess->RegisterMe(theTheoModel);
601  pmanager->AddDiscreteProcess(theInelasticProcess);
602  }
603  else if (particleName == "anti_omega-") {
604  pmanager->AddDiscreteProcess(theElasticProcess);
605  G4AntiOmegaMinusInelasticProcess* theInelasticProcess =
606  new G4AntiOmegaMinusInelasticProcess("inelastic");
607  G4LEAntiOmegaMinusInelastic* theInelasticModel =
609  theInelasticProcess->RegisterMe(theInelasticModel);
610  theInelasticProcess->RegisterMe(theTheoModel);
611  pmanager->AddDiscreteProcess(theInelasticProcess);
612  }
613  }
614 }
615 
617 {;}
618 
619 #include "G4Decay.hh"
621 {
622  G4Decay* theDecayProcess = new G4Decay();
624  while( (*theParticleIterator)() ){
626  G4ProcessManager* pmanager = particle->GetProcessManager();
627  if (theDecayProcess->IsApplicable(*particle)) {
628  pmanager ->AddProcess(theDecayProcess);
629  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
630  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
631  }
632  }
633 }
634 
636 {
637  if (verboseLevel >0)
638  {
639  G4cout << "B01PhysicsList::SetCuts:";
640  G4cout << "CutLength : " << defaultCutValue/mm << " (mm)" << G4endl;
641  }
642  // "G4VUserPhysicsList::SetCutsWithDefault" method sets
643  // the default cut value for all particle types
645 }
646 
649 
650  G4int npw = paraWorldName.size();
651  for ( G4int i = 0; i < npw; i++){
652  G4ParallelWorldScoringProcess* theParallelWorldScoringProcess
653  = new G4ParallelWorldScoringProcess("ParaWorldScoringProc");
654  theParallelWorldScoringProcess->SetParallelWorld(paraWorldName[i]);
655 
657  while( (*theParticleIterator)() ){
659  if ( !particle->IsShortLived() ){
660  G4ProcessManager* pmanager = particle->GetProcessManager();
661  pmanager->AddProcess(theParallelWorldScoringProcess);
662  pmanager->SetProcessOrderingToLast(theParallelWorldScoringProcess,idxAtRest);
663  pmanager->SetProcessOrdering(theParallelWorldScoringProcess,idxAlongStep,1);
664  pmanager->SetProcessOrderingToLast(theParallelWorldScoringProcess,idxPostStep);
665  }
666  }
667  }
668 
669 }