source: trunk/source/processes/hadronic/models/photolepton_hadron/muon_nuclear/test/src/PlhPhysicsList.cc @ 1199

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

nvx fichiers dans CVS

File size: 25.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: PlhPhysicsList.cc,v 1.2 2008/12/18 13:02:26 gunter Exp $
28// GEANT4 tag $Name: geant4-09-03-cand-01 $
29//
30#include "globals.hh"
31#include "PlhPhysicsList.hh"
32#include "G4ParticleDefinition.hh"
33#include "G4ParticleWithCuts.hh"
34#include "G4ProcessManager.hh"
35#include "G4ProcessVector.hh"
36#include "G4ParticleTypes.hh"
37#include "G4ParticleTable.hh"
38#include "G4Material.hh"
39#include "G4MaterialTable.hh"
40#include "G4ios.hh"
41#include <iomanip>                
42
43
44PlhPhysicsList::PlhPhysicsList():  G4VUserPhysicsList()
45{
46  SetVerboseLevel(1);
47}
48
49PlhPhysicsList::~PlhPhysicsList()
50{
51}
52
53void PlhPhysicsList::ConstructParticle()
54{
55  // In this method, static member functions should be called
56  // for all particles which you want to use.
57  // This ensures that objects of these particle types will be
58  // created in the program.
59
60  ConstructBosons();
61  ConstructLeptons();
62  ConstructMesons();
63  ConstructBarions();
64  ConstructIons();
65}
66
67void PlhPhysicsList::ConstructBosons()
68{
69  // pseudo-particles
70  G4Geantino::GeantinoDefinition();
71  G4ChargedGeantino::ChargedGeantinoDefinition();
72
73  // gamma
74  G4Gamma::GammaDefinition();
75
76  // optical photon
77  G4OpticalPhoton::OpticalPhotonDefinition();
78}
79
80void PlhPhysicsList::ConstructLeptons()
81{
82  // leptons
83  G4MuonPlus::MuonPlusDefinition();
84  G4MuonMinus::MuonMinusDefinition();
85  G4TauMinus::TauMinusDefinition();
86  G4TauPlus::TauPlusDefinition();
87  G4Electron::ElectronDefinition();
88  G4Positron::PositronDefinition();
89  G4NeutrinoTau::NeutrinoTauDefinition();
90  G4AntiNeutrinoTau::AntiNeutrinoTauDefinition();
91  G4NeutrinoMu::NeutrinoMuDefinition();
92  G4AntiNeutrinoMu::AntiNeutrinoMuDefinition();
93  G4NeutrinoE::NeutrinoEDefinition();
94  G4AntiNeutrinoE::AntiNeutrinoEDefinition();
95}
96
97void PlhPhysicsList::ConstructMesons()
98{
99 //  mesons
100  G4PionPlus::PionPlusDefinition();
101  G4PionMinus::PionMinusDefinition();
102  G4PionZero::PionZeroDefinition();
103  G4Eta::EtaDefinition();
104  G4EtaPrime::EtaPrimeDefinition();
105  G4RhoZero::RhoZeroDefinition();
106  G4KaonPlus::KaonPlusDefinition();
107  G4KaonMinus::KaonMinusDefinition();
108  G4KaonZero::KaonZeroDefinition();
109  G4AntiKaonZero::AntiKaonZeroDefinition();
110  G4KaonZeroLong::KaonZeroLongDefinition();
111  G4KaonZeroShort::KaonZeroShortDefinition();
112
113  G4DMesonPlus::DMesonPlusDefinition();
114  G4DMesonMinus::DMesonMinusDefinition();
115  G4DMesonZero::DMesonZeroDefinition();
116  G4AntiDMesonZero::AntiDMesonZeroDefinition();
117  G4DsMesonPlus::DsMesonPlusDefinition();
118  G4DsMesonMinus::DsMesonMinusDefinition();
119  G4JPsi::JPsiDefinition();
120
121  G4BMesonPlus::BMesonPlusDefinition();
122  G4BMesonMinus::BMesonMinusDefinition();
123  G4BMesonZero::BMesonZeroDefinition();
124  G4AntiBMesonZero::AntiBMesonZeroDefinition();
125  G4BsMesonZero::BsMesonZeroDefinition();
126  G4AntiBsMesonZero::AntiBsMesonZeroDefinition();
127}
128
129void PlhPhysicsList::ConstructBarions()
130{
131//  barions
132  G4Proton::ProtonDefinition();
133  G4AntiProton::AntiProtonDefinition();
134  G4Neutron::NeutronDefinition();
135  G4AntiNeutron::AntiNeutronDefinition();
136  G4Lambda::LambdaDefinition();
137  G4AntiLambda::AntiLambdaDefinition();
138  G4SigmaZero::SigmaZeroDefinition();
139  G4AntiSigmaZero::AntiSigmaZeroDefinition();
140  G4SigmaPlus::SigmaPlusDefinition();
141  G4AntiSigmaPlus::AntiSigmaPlusDefinition();
142  G4SigmaMinus::SigmaMinusDefinition();
143  G4AntiSigmaMinus::AntiSigmaMinusDefinition();
144  G4XiZero::XiZeroDefinition();
145  G4AntiXiZero::AntiXiZeroDefinition();
146  G4XiMinus::XiMinusDefinition();
147  G4AntiXiMinus::AntiXiMinusDefinition();
148  G4OmegaMinus::OmegaMinusDefinition();
149  G4AntiOmegaMinus::AntiOmegaMinusDefinition();
150
151  G4LambdacPlus::LambdacPlusDefinition();
152  G4SigmacPlusPlus::SigmacPlusPlusDefinition();
153  G4SigmacPlus::SigmacPlusDefinition();
154  G4SigmacZero::SigmacZeroDefinition();
155  G4XicPlus::XicPlusDefinition();
156  G4XicZero::XicZeroDefinition();
157  G4OmegacZero::OmegacZeroDefinition();
158  G4AntiLambdacPlus::AntiLambdacPlusDefinition();
159  G4AntiSigmacPlusPlus::AntiSigmacPlusPlusDefinition();
160  G4AntiSigmacPlus::AntiSigmacPlusDefinition();
161  G4AntiSigmacZero::AntiSigmacZeroDefinition();
162  G4AntiXicPlus::AntiXicPlusDefinition();
163  G4AntiXicZero::AntiXicZeroDefinition();
164  G4AntiOmegacZero::AntiOmegacZeroDefinition();
165}
166
167void PlhPhysicsList::ConstructIons()
168{
169  //  nuclei
170  G4Alpha::AlphaDefinition();
171  G4Deuteron::DeuteronDefinition();
172  G4Triton::TritonDefinition();
173  G4He3::He3Definition();
174}
175
176void PlhPhysicsList::ConstructProcess()
177{
178  AddTransportation();
179  ConstructEM();
180  ConstructLeptHad();
181  // ConstructHad();
182  // ConstructGeneral();
183}
184
185#include "G4ComptonScattering.hh"
186#include "G4GammaConversion.hh"
187#include "G4PhotoElectricEffect.hh"
188
189#include "G4MultipleScattering.hh"
190
191#include "G4eIonisation.hh"
192#include "G4eBremsstrahlung.hh"
193#include "G4eplusAnnihilation.hh"
194
195#include "G4MuIonisation.hh"
196#include "G4MuBremsstrahlung.hh"
197#include "G4MuPairProduction.hh"
198
199#include "G4hIonisation.hh"
200void PlhPhysicsList::ConstructEM()
201{
202  theParticleIterator->reset();
203  while( (*theParticleIterator)() ){
204    G4ParticleDefinition* particle = theParticleIterator->value();
205    G4ProcessManager* pmanager = particle->GetProcessManager();
206    G4String particleName = particle->GetParticleName();
207     
208    if (particleName == "gamma") {
209    // gamma
210      // Construct processes for gamma
211      pmanager->AddDiscreteProcess(new G4GammaConversion());
212      pmanager->AddDiscreteProcess(new G4ComptonScattering());     
213      pmanager->AddDiscreteProcess(new G4PhotoElectricEffect());
214
215    } else if (particleName == "e-") {
216    //electron
217      // Construct processes for electron
218      pmanager->AddProcess(new G4MultipleScattering(),-1,1,1);
219      pmanager->AddProcess(new G4eIonisation(),-1,2,2);
220      pmanager->AddProcess(new G4eBremsstrahlung(),-1,-1,3);
221 
222    } else if (particleName == "e+") {
223    //positron
224      // Construct processes for positron
225     pmanager->AddProcess(new G4MultipleScattering(),-1,1,1);
226     
227     pmanager->AddProcess(new G4eIonisation(),-1,2,2);
228     pmanager->AddProcess(new G4eBremsstrahlung(),-1,-1,3);     
229     pmanager->AddProcess(new G4eplusAnnihilation(),0,-1,4);
230 
231    } else if( particleName == "mu+" || 
232               particleName == "mu-"    ) {
233    //muon 
234     // Construct processes for muon+
235     pmanager->AddProcess(new G4MultipleScattering(),-1,1,1);
236     pmanager->AddProcess(new G4MuIonisation(),-1,2,2);
237     pmanager->AddProcess(new G4MuBremsstrahlung(),-1,-1,3);
238     pmanager->AddProcess(new G4MuPairProduction(),-1,-1,4);       
239     
240    } else { 
241      if ((particle->GetPDGCharge() != 0.0) && 
242          (particle->GetParticleName() != "chargedgeantino")) {
243     // all others charged particles except geantino
244       pmanager->AddProcess(new G4MultipleScattering(),-1,1,1);
245       pmanager->AddProcess(new G4hIonisation(),-1,2,2);       
246     }
247    }
248  }
249}
250
251// Hadron Processes
252
253#include "G4HadronElasticProcess.hh"
254#include "G4HadronFissionProcess.hh"
255#include "G4HadronCaptureProcess.hh"
256
257#include "G4PionPlusInelasticProcess.hh"
258#include "G4PionMinusInelasticProcess.hh"
259#include "G4KaonPlusInelasticProcess.hh"
260#include "G4KaonZeroSInelasticProcess.hh"
261#include "G4KaonZeroLInelasticProcess.hh"
262#include "G4KaonMinusInelasticProcess.hh"
263#include "G4ProtonInelasticProcess.hh"
264#include "G4AntiProtonInelasticProcess.hh"
265#include "G4NeutronInelasticProcess.hh"
266#include "G4AntiNeutronInelasticProcess.hh"
267#include "G4LambdaInelasticProcess.hh"
268#include "G4AntiLambdaInelasticProcess.hh"
269#include "G4SigmaPlusInelasticProcess.hh"
270#include "G4SigmaMinusInelasticProcess.hh"
271#include "G4AntiSigmaPlusInelasticProcess.hh"
272#include "G4AntiSigmaMinusInelasticProcess.hh"
273#include "G4XiZeroInelasticProcess.hh"
274#include "G4XiMinusInelasticProcess.hh"
275#include "G4AntiXiZeroInelasticProcess.hh"
276#include "G4AntiXiMinusInelasticProcess.hh"
277#include "G4DeuteronInelasticProcess.hh"
278#include "G4TritonInelasticProcess.hh"
279#include "G4AlphaInelasticProcess.hh"
280#include "G4OmegaMinusInelasticProcess.hh"
281#include "G4AntiOmegaMinusInelasticProcess.hh"
282
283// Low-energy Models
284
285#include "G4LElastic.hh"
286#include "G4LFission.hh"
287#include "G4LCapture.hh"
288
289#include "G4LEPionPlusInelastic.hh"
290#include "G4LEPionMinusInelastic.hh"
291#include "G4LEKaonPlusInelastic.hh"
292#include "G4LEKaonZeroSInelastic.hh"
293#include "G4LEKaonZeroLInelastic.hh"
294#include "G4LEKaonMinusInelastic.hh"
295#include "G4LEProtonInelastic.hh"
296#include "G4LEAntiProtonInelastic.hh"
297#include "G4LENeutronInelastic.hh"
298#include "G4LEAntiNeutronInelastic.hh"
299#include "G4LELambdaInelastic.hh"
300#include "G4LEAntiLambdaInelastic.hh"
301#include "G4LESigmaPlusInelastic.hh"
302#include "G4LESigmaMinusInelastic.hh"
303#include "G4LEAntiSigmaPlusInelastic.hh"
304#include "G4LEAntiSigmaMinusInelastic.hh"
305#include "G4LEXiZeroInelastic.hh"
306#include "G4LEXiMinusInelastic.hh"
307#include "G4LEAntiXiZeroInelastic.hh"
308#include "G4LEAntiXiMinusInelastic.hh"
309#include "G4LEDeuteronInelastic.hh"
310#include "G4LETritonInelastic.hh"
311#include "G4LEAlphaInelastic.hh"
312#include "G4LEOmegaMinusInelastic.hh"
313#include "G4LEAntiOmegaMinusInelastic.hh"
314
315//
316// ConstructHad()
317//
318// Makes discrete physics processes for the hadrons, at present limited
319// to those particles with GHEISHA interactions (INTRC > 0).
320// The processes are: Elastic scattering, Inelastic scattering,
321// Fission (for neutron only), and Capture (neutron).
322//
323// F.W.Jones  06-JUL-1998
324//
325
326void PlhPhysicsList::ConstructHad()
327{
328   G4HadronElasticProcess* theElasticProcess = 
329                                    new G4HadronElasticProcess;
330   G4LElastic* theElasticModel = new G4LElastic;
331   theElasticProcess->RegisterMe(theElasticModel);
332
333   theParticleIterator->reset();
334   while ((*theParticleIterator)()) {
335      G4ParticleDefinition* particle = theParticleIterator->value();
336      G4ProcessManager* pmanager = particle->GetProcessManager();
337      G4String particleName = particle->GetParticleName();
338     
339      if (particleName == "pi+") {
340         pmanager->AddDiscreteProcess(theElasticProcess);
341         G4PionPlusInelasticProcess* theInelasticProcess = 
342                                new G4PionPlusInelasticProcess("inelastic");
343         G4LEPionPlusInelastic* theInelasticModel = 
344                                new G4LEPionPlusInelastic;
345         theInelasticProcess->RegisterMe(theInelasticModel);
346         pmanager->AddDiscreteProcess(theInelasticProcess);
347      }
348      else if (particleName == "pi-") {
349         pmanager->AddDiscreteProcess(theElasticProcess);
350         G4PionMinusInelasticProcess* theInelasticProcess = 
351                                new G4PionMinusInelasticProcess("inelastic");
352         G4LEPionMinusInelastic* theInelasticModel = 
353                                new G4LEPionMinusInelastic;
354         theInelasticProcess->RegisterMe(theInelasticModel);
355         pmanager->AddDiscreteProcess(theInelasticProcess);
356      }
357      else if (particleName == "kaon+") {
358         pmanager->AddDiscreteProcess(theElasticProcess);
359         G4KaonPlusInelasticProcess* theInelasticProcess = 
360                                  new G4KaonPlusInelasticProcess("inelastic");
361         G4LEKaonPlusInelastic* theInelasticModel = new G4LEKaonPlusInelastic;
362         theInelasticProcess->RegisterMe(theInelasticModel);
363         pmanager->AddDiscreteProcess(theInelasticProcess);
364      }
365      else if (particleName == "kaon0S") {
366         pmanager->AddDiscreteProcess(theElasticProcess);
367         G4KaonZeroSInelasticProcess* theInelasticProcess = 
368                             new G4KaonZeroSInelasticProcess("inelastic");
369         G4LEKaonZeroSInelastic* theInelasticModel = 
370                             new G4LEKaonZeroSInelastic;
371         theInelasticProcess->RegisterMe(theInelasticModel);
372         pmanager->AddDiscreteProcess(theInelasticProcess);
373      }
374      else if (particleName == "kaon0L") {
375         pmanager->AddDiscreteProcess(theElasticProcess);
376         G4KaonZeroLInelasticProcess* theInelasticProcess = 
377                             new G4KaonZeroLInelasticProcess("inelastic");
378         G4LEKaonZeroLInelastic* theInelasticModel = 
379                             new G4LEKaonZeroLInelastic;
380         theInelasticProcess->RegisterMe(theInelasticModel);
381         pmanager->AddDiscreteProcess(theInelasticProcess);
382      }
383      else if (particleName == "kaon-") {
384         pmanager->AddDiscreteProcess(theElasticProcess);
385         G4KaonMinusInelasticProcess* theInelasticProcess = 
386                                 new G4KaonMinusInelasticProcess("inelastic");
387         G4LEKaonMinusInelastic* theInelasticModel = 
388                                 new G4LEKaonMinusInelastic;
389         theInelasticProcess->RegisterMe(theInelasticModel);
390         pmanager->AddDiscreteProcess(theInelasticProcess);
391      }
392      else if (particleName == "proton") {
393         pmanager->AddDiscreteProcess(theElasticProcess);
394         G4ProtonInelasticProcess* theInelasticProcess = 
395                                    new G4ProtonInelasticProcess("inelastic");
396         G4LEProtonInelastic* theInelasticModel = new G4LEProtonInelastic;
397         theInelasticProcess->RegisterMe(theInelasticModel);
398         pmanager->AddDiscreteProcess(theInelasticProcess);
399      }
400      else if (particleName == "anti_proton") {
401         pmanager->AddDiscreteProcess(theElasticProcess);
402         G4AntiProtonInelasticProcess* theInelasticProcess = 
403                                new G4AntiProtonInelasticProcess("inelastic");
404         G4LEAntiProtonInelastic* theInelasticModel = 
405                                new G4LEAntiProtonInelastic;
406         theInelasticProcess->RegisterMe(theInelasticModel);
407         pmanager->AddDiscreteProcess(theInelasticProcess);
408      }
409      else if (particleName == "neutron") {
410         pmanager->AddDiscreteProcess(theElasticProcess);
411         G4NeutronInelasticProcess* theInelasticProcess = 
412                                    new G4NeutronInelasticProcess("inelastic");
413         G4LENeutronInelastic* theInelasticModel = new G4LENeutronInelastic;
414         theInelasticProcess->RegisterMe(theInelasticModel);
415         pmanager->AddDiscreteProcess(theInelasticProcess);
416         G4HadronFissionProcess* theFissionProcess =
417                                    new G4HadronFissionProcess;
418         G4LFission* theFissionModel = new G4LFission;
419         theFissionProcess->RegisterMe(theFissionModel);
420         pmanager->AddDiscreteProcess(theFissionProcess);
421         G4HadronCaptureProcess* theCaptureProcess =
422                                    new G4HadronCaptureProcess;
423         G4LCapture* theCaptureModel = new G4LCapture;
424         theCaptureProcess->RegisterMe(theCaptureModel);
425         pmanager->AddDiscreteProcess(theCaptureProcess);
426      } 
427      else if (particleName == "anti_neutron") {
428         pmanager->AddDiscreteProcess(theElasticProcess);
429         G4AntiNeutronInelasticProcess* theInelasticProcess = 
430                               new G4AntiNeutronInelasticProcess("inelastic");
431         G4LEAntiNeutronInelastic* theInelasticModel = 
432                               new G4LEAntiNeutronInelastic;
433         theInelasticProcess->RegisterMe(theInelasticModel);
434         pmanager->AddDiscreteProcess(theInelasticProcess);
435      }
436      else if (particleName == "lambda") {
437         pmanager->AddDiscreteProcess(theElasticProcess);
438         G4LambdaInelasticProcess* theInelasticProcess = 
439                                    new G4LambdaInelasticProcess("inelastic");
440         G4LELambdaInelastic* theInelasticModel = new G4LELambdaInelastic;
441         theInelasticProcess->RegisterMe(theInelasticModel);
442         pmanager->AddDiscreteProcess(theInelasticProcess);
443      }
444      else if (particleName == "anti_lambda") {
445         pmanager->AddDiscreteProcess(theElasticProcess);
446         G4AntiLambdaInelasticProcess* theInelasticProcess = 
447                                new G4AntiLambdaInelasticProcess("inelastic");
448         G4LEAntiLambdaInelastic* theInelasticModel = 
449                                new G4LEAntiLambdaInelastic;
450         theInelasticProcess->RegisterMe(theInelasticModel);
451         pmanager->AddDiscreteProcess(theInelasticProcess);
452      }
453      else if (particleName == "sigma+") {
454         pmanager->AddDiscreteProcess(theElasticProcess);
455         G4SigmaPlusInelasticProcess* theInelasticProcess = 
456                                 new G4SigmaPlusInelasticProcess("inelastic");
457         G4LESigmaPlusInelastic* theInelasticModel = 
458                                 new G4LESigmaPlusInelastic;
459         theInelasticProcess->RegisterMe(theInelasticModel);
460         pmanager->AddDiscreteProcess(theInelasticProcess);
461      }
462      else if (particleName == "sigma-") {
463         pmanager->AddDiscreteProcess(theElasticProcess);
464         G4SigmaMinusInelasticProcess* theInelasticProcess = 
465                                 new G4SigmaMinusInelasticProcess("inelastic");
466         G4LESigmaMinusInelastic* theInelasticModel = 
467                                 new G4LESigmaMinusInelastic;
468         theInelasticProcess->RegisterMe(theInelasticModel);
469         pmanager->AddDiscreteProcess(theInelasticProcess);
470      }
471      else if (particleName == "anti_sigma+") {
472         pmanager->AddDiscreteProcess(theElasticProcess);
473         G4AntiSigmaPlusInelasticProcess* theInelasticProcess = 
474                             new G4AntiSigmaPlusInelasticProcess("inelastic");
475         G4LEAntiSigmaPlusInelastic* theInelasticModel = 
476                                 new G4LEAntiSigmaPlusInelastic;
477         theInelasticProcess->RegisterMe(theInelasticModel);
478         pmanager->AddDiscreteProcess(theInelasticProcess);
479      }
480      else if (particleName == "anti_sigma-") {
481         pmanager->AddDiscreteProcess(theElasticProcess);
482         G4AntiSigmaMinusInelasticProcess* theInelasticProcess = 
483                            new G4AntiSigmaMinusInelasticProcess("inelastic");
484         G4LEAntiSigmaMinusInelastic* theInelasticModel = 
485                                 new G4LEAntiSigmaMinusInelastic;
486         theInelasticProcess->RegisterMe(theInelasticModel);
487         pmanager->AddDiscreteProcess(theInelasticProcess);
488      }
489      else if (particleName == "xi0") {
490         pmanager->AddDiscreteProcess(theElasticProcess);
491         G4XiZeroInelasticProcess* theInelasticProcess = 
492                            new G4XiZeroInelasticProcess("inelastic");
493         G4LEXiZeroInelastic* theInelasticModel = 
494                                 new G4LEXiZeroInelastic;
495         theInelasticProcess->RegisterMe(theInelasticModel);
496         pmanager->AddDiscreteProcess(theInelasticProcess);
497      }
498      else if (particleName == "xi-") {
499         pmanager->AddDiscreteProcess(theElasticProcess);
500         G4XiMinusInelasticProcess* theInelasticProcess = 
501                            new G4XiMinusInelasticProcess("inelastic");
502         G4LEXiMinusInelastic* theInelasticModel = 
503                                 new G4LEXiMinusInelastic;
504         theInelasticProcess->RegisterMe(theInelasticModel);
505         pmanager->AddDiscreteProcess(theInelasticProcess);
506      }
507      else if (particleName == "anti_xi0") {
508         pmanager->AddDiscreteProcess(theElasticProcess);
509         G4AntiXiZeroInelasticProcess* theInelasticProcess = 
510                            new G4AntiXiZeroInelasticProcess("inelastic");
511         G4LEAntiXiZeroInelastic* theInelasticModel = 
512                                 new G4LEAntiXiZeroInelastic;
513         theInelasticProcess->RegisterMe(theInelasticModel);
514         pmanager->AddDiscreteProcess(theInelasticProcess);
515      }
516      else if (particleName == "anti_xi-") {
517         pmanager->AddDiscreteProcess(theElasticProcess);
518         G4AntiXiMinusInelasticProcess* theInelasticProcess = 
519                            new G4AntiXiMinusInelasticProcess("inelastic");
520         G4LEAntiXiMinusInelastic* theInelasticModel = 
521                                 new G4LEAntiXiMinusInelastic;
522         theInelasticProcess->RegisterMe(theInelasticModel);
523         pmanager->AddDiscreteProcess(theInelasticProcess);
524      }
525      else if (particleName == "deuteron") {
526         pmanager->AddDiscreteProcess(theElasticProcess);
527         G4DeuteronInelasticProcess* theInelasticProcess = 
528                            new G4DeuteronInelasticProcess("inelastic");
529         G4LEDeuteronInelastic* theInelasticModel = 
530                                 new G4LEDeuteronInelastic;
531         theInelasticProcess->RegisterMe(theInelasticModel);
532         pmanager->AddDiscreteProcess(theInelasticProcess);
533      }
534      else if (particleName == "triton") {
535         pmanager->AddDiscreteProcess(theElasticProcess);
536         G4TritonInelasticProcess* theInelasticProcess = 
537                            new G4TritonInelasticProcess("inelastic");
538         G4LETritonInelastic* theInelasticModel = 
539                                 new G4LETritonInelastic;
540         theInelasticProcess->RegisterMe(theInelasticModel);
541         pmanager->AddDiscreteProcess(theInelasticProcess);
542      }
543      else if (particleName == "alpha") {
544         pmanager->AddDiscreteProcess(theElasticProcess);
545         G4AlphaInelasticProcess* theInelasticProcess = 
546                            new G4AlphaInelasticProcess("inelastic");
547         G4LEAlphaInelastic* theInelasticModel = 
548                                 new G4LEAlphaInelastic;
549         theInelasticProcess->RegisterMe(theInelasticModel);
550         pmanager->AddDiscreteProcess(theInelasticProcess);
551      }
552      else if (particleName == "omega-") {
553         pmanager->AddDiscreteProcess(theElasticProcess);
554         G4OmegaMinusInelasticProcess* theInelasticProcess = 
555                            new G4OmegaMinusInelasticProcess("inelastic");
556         G4LEOmegaMinusInelastic* theInelasticModel = 
557                                 new G4LEOmegaMinusInelastic;
558         theInelasticProcess->RegisterMe(theInelasticModel);
559         pmanager->AddDiscreteProcess(theInelasticProcess);
560      }
561      else if (particleName == "anti_omega-") {
562         pmanager->AddDiscreteProcess(theElasticProcess);
563         G4AntiOmegaMinusInelasticProcess* theInelasticProcess = 
564                            new G4AntiOmegaMinusInelasticProcess("inelastic");
565         G4LEAntiOmegaMinusInelastic* theInelasticModel = 
566                                 new G4LEAntiOmegaMinusInelastic;
567         theInelasticProcess->RegisterMe(theInelasticModel);
568         pmanager->AddDiscreteProcess(theInelasticProcess);
569      }
570   }
571}
572
573#include "G4MuonNucleusProcess.hh"
574void PlhPhysicsList::ConstructLeptHad()
575{
576  G4MuonNucleusProcess* theMuonNucleusProcess = new G4MuonNucleusProcess();
577  theParticleIterator->reset();
578  while( (*theParticleIterator)() ){
579    G4ParticleDefinition* particle = theParticleIterator->value();
580    G4ProcessManager* pmanager = particle->GetProcessManager();
581    G4String particleName = particle->GetParticleName();
582    if( particleName == "mu+" || particleName == "mu-" ) {
583      pmanager->AddDiscreteProcess(theMuonNucleusProcess);
584    }
585  }
586}
587
588
589#include "G4Decay.hh"
590void PlhPhysicsList::ConstructGeneral()
591{
592  G4Decay* theDecayProcess = new G4Decay();
593  theParticleIterator->reset();
594  while( (*theParticleIterator)() ){
595    G4ParticleDefinition* particle = theParticleIterator->value();
596    G4ProcessManager* pmanager = particle->GetProcessManager();
597    if (theDecayProcess->IsApplicable(*particle)) { 
598      pmanager ->AddProcess(theDecayProcess);
599      pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
600      pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
601    }
602  }
603}
604
605void PlhPhysicsList::SetCuts(G4double cut)
606{
607  if (verboseLevel >0){
608    G4cout << "PlhPhysicsList::SetCuts:";
609    G4cout << "CutLength : " << cut/mm << " (mm)" << G4endl;
610  } 
611
612  // set cut values for gamma at first and for e- second and next for e+,
613  // because some processes for e+/e- need cut values for gamma
614  SetCutValue(cut, "gamma");
615  SetCutValue(cut, "e-");
616  SetCutValue(cut, "e+");
617 
618  // set cut values for proton and anti_proton before all other hadrons
619  // because some processes for hadrons need cut values for proton/anti_proton
620  SetCutValue(cut, "proton");
621  SetCutValue(cut, "anti_proton");
622 
623  SetCutValueForOthers(cut);
624
625  if (verboseLevel>1) {
626    DumpCutValuesTable();
627  }
628}
629
630
Note: See TracBrowser for help on using the repository browser.