source: trunk/examples/extended/biasing/B02/src/B02PhysicsList.cc @ 1342

Last change on this file since 1342 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

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