Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
electromagnetic/TestEm7/src/PhysListEmStandardSS.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 //
26 /// \file electromagnetic/TestEm7/src/PhysListEmStandardSS.cc
27 /// \brief Implementation of the PhysListEmStandardSS class
28 //
29 // $Id: PhysListEmStandardSS.cc 68585 2013-04-01 23:35:07Z adotti $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "PhysListEmStandardSS.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 
57 #include "G4eIonisation.hh"
58 #include "G4eBremsstrahlung.hh"
59 #include "G4Generator2BS.hh"
60 #include "G4SeltzerBergerModel.hh"
63 
64 #include "G4eplusAnnihilation.hh"
65 #include "G4UAtomicDeexcitation.hh"
66 
67 #include "G4MuIonisation.hh"
68 #include "G4MuBremsstrahlung.hh"
69 #include "G4MuPairProduction.hh"
70 
71 #include "G4hIonisation.hh"
72 #include "G4ionIonisation.hh"
74 
75 #include "G4PhysicsListHelper.hh"
76 #include "G4BuilderType.hh"
77 
78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
79 
81  : G4VPhysicsConstructor(name)
82 {
85 }
86 
87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88 
90 {}
91 
92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93 
95 {
97 
98  // muon & hadron bremsstrahlung and pair production
101 
102  aParticleIterator->reset();
103  while( (*aParticleIterator)() ){
104  G4ParticleDefinition* particle = aParticleIterator->value();
105  G4String particleName = particle->GetParticleName();
106 
107  if (particleName == "gamma") {
108 
109  // Compton scattering
111  cs->SetEmModel(new G4KleinNishinaModel(),1);
112  G4VEmModel* theLowEPComptonModel = new G4LowEPComptonModel();
113  theLowEPComptonModel->SetHighEnergyLimit(20*MeV);
114  cs->AddEmModel(0, theLowEPComptonModel);
115  ph->RegisterProcess(cs, particle);
116 
117  // Photoelectric
119  G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
120  theLivermorePEModel->SetHighEnergyLimit(10*GeV);
121  pe->SetEmModel(theLivermorePEModel,1);
122  ph->RegisterProcess(pe, particle);
123 
124  // Gamma conversion
126  G4VEmModel* thePenelopeGCModel = new G4PenelopeGammaConversionModel();
127  thePenelopeGCModel->SetHighEnergyLimit(1*GeV);
128  gc->SetEmModel(thePenelopeGCModel,1);
129  ph->RegisterProcess(gc, particle);
130 
131  // Rayleigh scattering
132  ph->RegisterProcess(new G4RayleighScattering(), particle);
133 
134  } else if (particleName == "e-") {
135 
136  // ionisation
137  G4eIonisation* eIoni = new G4eIonisation();
138  eIoni->SetStepFunction(0.2, 100*um);
139  G4VEmModel* theIoniPenelope = new G4PenelopeIonisationModel();
140  theIoniPenelope->SetHighEnergyLimit(0.1*MeV);
141  eIoni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
142 
143  // bremsstrahlung
144  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
145 
146  ph->RegisterProcess(eIoni, particle);
147  ph->RegisterProcess(eBrem, particle);
148  ph->RegisterProcess(new G4CoulombScattering(), particle);
149 
150  } else if (particleName == "e+") {
151  // ionisation
152  G4eIonisation* eIoni = new G4eIonisation();
153  eIoni->SetStepFunction(0.2, 100*um);
154  G4VEmModel* theIoniPenelope = new G4PenelopeIonisationModel();
155  theIoniPenelope->SetHighEnergyLimit(0.1*MeV);
156  eIoni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
157 
158  // bremsstrahlung
159  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
160 
161  ph->RegisterProcess(eIoni, particle);
162  ph->RegisterProcess(eBrem, particle);
163  ph->RegisterProcess(new G4CoulombScattering(), particle);
164 
165  // annihilation at rest and in flight
166  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
167 
168  } else if (particleName == "mu+" ||
169  particleName == "mu-" ) {
170 
171  G4MuIonisation* muIoni = new G4MuIonisation();
172  muIoni->SetStepFunction(0.2, 50*um);
173 
174  ph->RegisterProcess(muIoni, particle);
175  ph->RegisterProcess(mub, particle);
176  ph->RegisterProcess(mup, particle);
177  ph->RegisterProcess(new G4CoulombScattering(), particle);
178 
179  } else if (particleName == "alpha" || particleName == "He3") {
180 
181  G4ionIonisation* ionIoni = new G4ionIonisation();
182  ionIoni->SetStepFunction(0.1, 10*um);
183  ph->RegisterProcess(ionIoni, particle);
184 
186  cs->SetBuildTableFlag(false);
187  ph->RegisterProcess(cs, particle);
188 
189  } else if (particleName == "GenericIon" ) {
190 
191  G4ionIonisation* ionIoni = new G4ionIonisation();
192  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
193  ionIoni->SetStepFunction(0.1, 1*um);
194  ph->RegisterProcess(ionIoni, particle);
195 
197  cs->SetBuildTableFlag(false);
198  ph->RegisterProcess(cs, particle);
199 
200  } else if ((!particle->IsShortLived()) &&
201  (particle->GetPDGCharge() != 0.0) &&
202  (particle->GetParticleName() != "chargedgeantino")) {
203  //all others charged particles except geantino
204 
205  ph->RegisterProcess(new G4hIonisation(), particle);
207  cs->SetBuildTableFlag(false);
208  ph->RegisterProcess(cs, particle);
209 
210  }
211  }
212 
213  // Em options
214  //
215  // Main options and setting parameters are shown here.
216  // Several of them have default values.
217  //
218  G4EmProcessOptions emOptions;
219 
220  //physics tables
221  //
222  emOptions.SetMinEnergy(10*eV);
223  emOptions.SetMaxEnergy(10*TeV);
224  emOptions.SetDEDXBinning(12*20);
225  emOptions.SetLambdaBinning(12*20);
226 
227  // scattering
228  emOptions.SetPolarAngleLimit(0.0);
229 
230  // Deexcitation
233  de->SetFluo(true);
234 }
235 
236 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
237 
static G4LossTableManager * Instance()
PhysListEmStandardSS(const G4String &name="standardSS")
void SetBuildTableFlag(G4bool val)
void SetMinEnergy(G4double val)
void SetStepFunction(G4double v1, G4double v2)
const XML_Char * name
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
void SetDEDXBinning(G4int val)
void SetEmModel(G4VEmModel *, G4int index=1)
void SetLambdaBinning(G4int val)
#define aParticleIterator
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
void SetMaxEnergy(G4double val)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
static G4PhysicsListHelper * GetPhysicsListHelper()
void SetEmModel(G4VEmModel *, G4int index=1)
G4double GetPDGCharge() const
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetPolarAngleLimit(G4double val)