Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PhysListEmStandardNR.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 // $Id: PhysListEmStandardNR.cc 100284 2016-10-17 08:41:43Z gcosmo $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "PhysListEmStandardNR.hh"
35 
36 #include "G4SystemOfUnits.hh"
37 #include "G4ParticleDefinition.hh"
38 #include "G4LossTableManager.hh"
39 #include "G4EmProcessOptions.hh"
40 
41 #include "G4ComptonScattering.hh"
42 #include "G4GammaConversion.hh"
43 #include "G4PhotoElectricEffect.hh"
44 #include "G4RayleighScattering.hh"
45 #include "G4PEEffectFluoModel.hh"
46 #include "G4KleinNishinaModel.hh"
47 #include "G4LowEPComptonModel.hh"
50 
51 #include "G4eMultipleScattering.hh"
53 #include "G4hMultipleScattering.hh"
54 #include "G4CoulombScattering.hh"
56 #include "G4UrbanMscModel.hh"
57 
58 #include "G4eIonisation.hh"
59 #include "G4eBremsstrahlung.hh"
60 #include "G4Generator2BS.hh"
61 #include "G4SeltzerBergerModel.hh"
64 
65 #include "G4eplusAnnihilation.hh"
66 #include "G4UAtomicDeexcitation.hh"
67 
68 #include "G4MuIonisation.hh"
69 #include "G4MuBremsstrahlung.hh"
70 #include "G4MuPairProduction.hh"
71 
72 #include "G4hIonisation.hh"
73 #include "G4ionIonisation.hh"
75 
77 
78 #include "G4PhysicsListHelper.hh"
79 #include "G4BuilderType.hh"
80 
81 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
82 
84  : G4VPhysicsConstructor(name)
85 {
88 }
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
91 
93 {}
94 
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96 
98 {
100 
101  // muon & hadron bremsstrahlung and pair production
104 
106  G4double energyLimit = 100.*MeV;
107  nucr->SetMaxEnergyForScattering(energyLimit);
109  csm->SetActivationLowEnergyLimit(energyLimit);
110 
112  particleIterator->reset();
113  while( (*particleIterator)() ){
114  G4ParticleDefinition* particle = particleIterator->value();
115  G4String particleName = particle->GetParticleName();
116 
117  if (particleName == "gamma") {
118 
119  // Compton scattering
121  cs->SetEmModel(new G4KleinNishinaModel(),1);
122  ph->RegisterProcess(cs, particle);
123 
124  // Photoelectric
126  G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
127  theLivermorePEModel->SetHighEnergyLimit(10*GeV);
128  pe->SetEmModel(theLivermorePEModel,1);
129  ph->RegisterProcess(pe, particle);
130 
131  // Gamma conversion
133  G4VEmModel* thePenelopeGCModel = new G4PenelopeGammaConversionModel();
134  thePenelopeGCModel->SetHighEnergyLimit(1*GeV);
135  gc->SetEmModel(thePenelopeGCModel,1);
136  ph->RegisterProcess(gc, particle);
137 
138  // Rayleigh scattering
139  ph->RegisterProcess(new G4RayleighScattering(), particle);
140 
141  } else if (particleName == "e-") {
142 
143  // ionisation
144  G4eIonisation* eIoni = new G4eIonisation();
145  eIoni->SetStepFunction(0.2, 100*um);
146 
147  // bremsstrahlung
148  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
149 
150  ph->RegisterProcess(new G4eMultipleScattering(), particle);
151  ph->RegisterProcess(eIoni, particle);
152  ph->RegisterProcess(eBrem, particle);
153 
154  } else if (particleName == "e+") {
155  // ionisation
156  G4eIonisation* eIoni = new G4eIonisation();
157  eIoni->SetStepFunction(0.2, 100*um);
158 
159  // bremsstrahlung
160  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
161 
162  ph->RegisterProcess(new G4eMultipleScattering(), particle);
163  ph->RegisterProcess(eIoni, particle);
164  ph->RegisterProcess(eBrem, particle);
165 
166  // annihilation at rest and in flight
167  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
168 
169  } else if (particleName == "mu+" ||
170  particleName == "mu-" ) {
171 
172  G4MuIonisation* muIoni = new G4MuIonisation();
173  muIoni->SetStepFunction(0.2, 50*um);
174 
175  ph->RegisterProcess(muIoni, particle);
176  ph->RegisterProcess(mub, particle);
177  ph->RegisterProcess(mup, particle);
178  ph->RegisterProcess(new G4CoulombScattering(), particle);
179 
180  } else if (particleName == "alpha" || particleName == "He3") {
181 
184  model->SetActivationLowEnergyLimit(energyLimit);
185  msc->SetEmModel(model, 1);
186  ph->RegisterProcess(msc, particle);
187 
188  G4ionIonisation* ionIoni = new G4ionIonisation();
189  ionIoni->SetStepFunction(0.1, 10*um);
190  ph->RegisterProcess(ionIoni, particle);
191 
192  ph->RegisterProcess(nucr, particle);
193 
194  } else if (particleName == "GenericIon" ) {
195 
198  model->SetActivationLowEnergyLimit(energyLimit);
199  msc->SetEmModel(model, 1);
200  ph->RegisterProcess(msc, particle);
201 
202  G4ionIonisation* ionIoni = new G4ionIonisation();
203  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
204  ionIoni->SetStepFunction(0.1, 1*um);
205  ph->RegisterProcess(ionIoni, particle);
206 
207  ph->RegisterProcess(nucr, particle);
208 
209  } else if (particleName == "proton" ||
210  particleName == "deuteron" ||
211  particleName == "triton") {
212 
215  model->SetActivationLowEnergyLimit(energyLimit);
216  msc->SetEmModel(model, 1);
217  ph->RegisterProcess(msc, particle);
218 
219  G4hIonisation* hIoni = new G4hIonisation();
220  hIoni->SetStepFunction(0.05, 1*um);
221  ph->RegisterProcess(hIoni, particle);
222 
223  ph->RegisterProcess(nucr, particle);
224 
225  } else if ((!particle->IsShortLived()) &&
226  (particle->GetPDGCharge() != 0.0) &&
227  (particle->GetParticleName() != "chargedgeantino")) {
228  //all others charged particles except geantino
229 
230  ph->RegisterProcess(new G4hMultipleScattering(), particle);
231  ph->RegisterProcess(new G4hIonisation(), particle);
232  }
233  }
234 
235  // Em options
236  //
237  // Main options and setting parameters are shown here.
238  // Several of them have default values.
239  //
240  G4EmProcessOptions emOptions;
241 
242  //physics tables
243  //
244  emOptions.SetMinEnergy(10*eV);
245  emOptions.SetMaxEnergy(10*TeV);
246  emOptions.SetDEDXBinning(12*20);
247  emOptions.SetLambdaBinning(12*20);
248 
249  // scattering
250  emOptions.SetPolarAngleLimit(0.0);
251 
252  // Deexcitation
255  de->SetFluo(true);
256 }
257 
258 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
259 
const XML_Char * name
Definition: expat.h:151
static G4LossTableManager * Instance()
Definition of the G4ScreenedNuclearRecoil class.
void SetMinEnergy(G4double val)
void SetEmModel(G4VMscModel *, G4int index=1)
void SetStepFunction(G4double v1, G4double v2, G4bool lock=true)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
static constexpr double TeV
Definition: G4SIunits.hh:218
void SetDEDXBinning(G4int val)
virtual void ConstructProcess()
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
void SetEmModel(G4VEmModel *, G4int index=1)
static constexpr double um
Definition: G4SIunits.hh:113
void SetLambdaBinning(G4int val)
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static constexpr double eV
Definition: G4SIunits.hh:215
void SetMaxEnergyForScattering(G4double energy)
set the upper energy beyond which this process has no cross section
PhysListEmStandardNR(const G4String &name="standardNR")
void SetMaxEnergy(G4double val)
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
static G4PhysicsListHelper * GetPhysicsListHelper()
void SetEmModel(G4VEmModel *, G4int index=1)
static constexpr double GeV
Definition: G4SIunits.hh:217
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
static constexpr double MeV
Definition: G4SIunits.hh:214
Definition of the PhysListEmStandardNR class.
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
A process which handles screened Coulomb collisions between nuclei.
const XML_Char XML_Content * model
Definition: expat.h:151
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetPolarAngleLimit(G4double val)