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

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

tag geant4.9.4 beta 1 + modifs locales

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