source: trunk/examples/extended/biasing/B01/src/B01PhysicsList.cc @ 1333

Last change on this file since 1333 was 1230, checked in by garnier, 15 years ago

update to geant4.9.3

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