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