Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PhysListEmLivermore.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
28 //
29 //
30 // $Id: PhysListEmLivermore.cc 100278 2016-10-17 08:34:03Z gcosmo $
31 //
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
34 
35 #include "PhysListEmLivermore.hh"
36 #include "G4ParticleDefinition.hh"
37 #include "G4ProcessManager.hh"
38 
39 // gamma
40 
41 #include "G4PhotoElectricEffect.hh"
43 
44 #include "G4ComptonScattering.hh"
46 
47 #include "G4GammaConversion.hh"
49 
50 #include "G4RayleighScattering.hh"
52 
53 // e-
54 
55 #include "G4eIonisation.hh"
58 
59 #include "G4eBremsstrahlung.hh"
61 
62 // e+
63 
64 #include "G4eplusAnnihilation.hh"
65 
66 // mu
67 
68 #include "G4MuIonisation.hh"
69 #include "G4MuBremsstrahlung.hh"
70 #include "G4MuPairProduction.hh"
71 
72 // hadrons, ions
73 
74 #include "G4hIonisation.hh"
75 #include "G4ionIonisation.hh"
76 
77 #include "G4LossTableManager.hh"
78 #include "G4UAtomicDeexcitation.hh"
79 
80 #include "G4SystemOfUnits.hh"
81 
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83 
85  : G4VPhysicsConstructor(name)
86 { }
87 
88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
89 
91 { }
92 
93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
94 
96 {
97  // Add standard EM Processes
98 
100  particleIterator->reset();
101  while( (*particleIterator)() ){
102  G4ParticleDefinition* particle = particleIterator->value();
103  G4ProcessManager* pmanager = particle->GetProcessManager();
104  G4String particleName = particle->GetParticleName();
105 
106  //Applicability range for Livermore models
107  //for higher energies, the Standard models are used
108  G4double highEnergyLimit = 1*GeV;
109 
110  if (particleName == "gamma") {
111  // gamma
112 
115  photModel = new G4LivermorePhotoElectricModel();
116  photModel->SetHighEnergyLimit(highEnergyLimit);
117  phot->AddEmModel(0, photModel);
118  pmanager->AddDiscreteProcess(phot);
119 
122  comptModel = new G4LivermoreComptonModel();
123  comptModel->SetHighEnergyLimit(highEnergyLimit);
124  compt->AddEmModel(0, comptModel);
125  pmanager->AddDiscreteProcess(compt);
126 
127  G4GammaConversion* conv = new G4GammaConversion();
129  convModel = new G4LivermoreGammaConversionModel();
130  convModel->SetHighEnergyLimit(highEnergyLimit);
131  conv->AddEmModel(0, convModel);
132  pmanager->AddDiscreteProcess(conv);
133 
136  raylModel = new G4LivermoreRayleighModel();
137  raylModel->SetHighEnergyLimit(highEnergyLimit);
138  rayl->AddEmModel(0, raylModel);
139  pmanager->AddDiscreteProcess(rayl);
140 
141  } else if (particleName == "e-") {
142  //electron
143 
144  G4eIonisation* eIoni = new G4eIonisation();
146  eIoniModel = new G4LivermoreIonisationModel();
147  eIoniModel->SetHighEnergyLimit(highEnergyLimit);
148  eIoni->AddEmModel(0, eIoniModel, new G4UniversalFluctuation() );
149  pmanager->AddProcess(eIoni, -1,-1, 1);
150 
151  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
153  eBremModel = new G4LivermoreBremsstrahlungModel();
154  eBremModel->SetHighEnergyLimit(highEnergyLimit);
155  eBrem->AddEmModel(0, eBremModel);
156  pmanager->AddProcess(eBrem, -1,-1, 2);
157 
158  } else if (particleName == "e+") {
159  //positron
160  pmanager->AddProcess(new G4eIonisation, -1,-1, 1);
161  pmanager->AddProcess(new G4eBremsstrahlung, -1,-1, 2);
162  pmanager->AddProcess(new G4eplusAnnihilation, 0,-1, 3);
163 
164  } else if( particleName == "mu+" ||
165  particleName == "mu-" ) {
166  //muon
167  pmanager->AddProcess(new G4MuIonisation, -1,-1, 1);
168  pmanager->AddProcess(new G4MuBremsstrahlung, -1,-1, 2);
169  pmanager->AddProcess(new G4MuPairProduction, -1,-1, 3);
170 
171  } else if( particleName == "alpha" || particleName == "GenericIon" ) {
172  pmanager->AddProcess(new G4ionIonisation, -1,-1, 1);
173 
174  } else if ((!particle->IsShortLived()) &&
175  (particle->GetPDGCharge() != 0.0) &&
176  (particle->GetParticleName() != "chargedgeantino")) {
177  //all others charged particles except geantino
178  pmanager->AddProcess(new G4hIonisation, -1,-1, 1);
179  }
180  }
181 
182  // Deexcitation
183  //
185  de->SetFluo(true);
186  de->SetAuger(false);
187  de->SetPIXE(false);
189 }
190 
191 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
192 
const XML_Char * name
Definition: expat.h:151
static G4LossTableManager * Instance()
virtual void ConstructProcess()
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
G4ProcessManager * GetProcessManager() const
static constexpr double GeV
Definition: G4SIunits.hh:217
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
void SetAtomDeexcitation(G4VAtomDeexcitation *)
PhysListEmLivermore(const G4String &name="livermore")