Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
TestEm14/src/PhysListEmPenelope.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/TestEm14/src/PhysListEmPenelope.cc
27 /// \brief Implementation of the PhysListEmPenelope class
28 //
29 //
30 // $Id: PhysListEmPenelope.cc 68585 2013-04-01 23:35:07Z adotti $
31 //
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
34 
35 #include "PhysListEmPenelope.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"
66 
67 // mu
68 
69 #include "G4MuIonisation.hh"
70 #include "G4MuBremsstrahlung.hh"
71 #include "G4MuPairProduction.hh"
72 
73 // hadrons, ions
74 
75 #include "G4hIonisation.hh"
76 #include "G4ionIonisation.hh"
77 
78 #include "G4EmProcessOptions.hh"
79 #include "G4LossTableManager.hh"
80 #include "G4UAtomicDeexcitation.hh"
81 
82 #include "G4SystemOfUnits.hh"
83 
84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85 
87  : G4VPhysicsConstructor(name)
88 { }
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
91 
93 { }
94 
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96 
98 {
99  // Add standard EM Processes
100 
101  aParticleIterator->reset();
102  while( (*aParticleIterator)() ){
103  G4ParticleDefinition* particle = aParticleIterator->value();
104  G4ProcessManager* pmanager = particle->GetProcessManager();
105  G4String particleName = particle->GetParticleName();
106 
107  //Applicability range for Penelope models
108  //for higher energies, the Standard models are used
109  G4double highEnergyLimit = 1*GeV;
110 
111  if (particleName == "gamma") {
112  // gamma
113 
116  photModel = new G4PenelopePhotoElectricModel();
117  photModel->SetHighEnergyLimit(highEnergyLimit);
118  phot->AddEmModel(0, photModel);
119  pmanager->AddDiscreteProcess(phot);
120 
123  comptModel = new G4PenelopeComptonModel();
124  comptModel->SetHighEnergyLimit(highEnergyLimit);
125  compt->AddEmModel(0, comptModel);
126  pmanager->AddDiscreteProcess(compt);
127 
128  G4GammaConversion* conv = new G4GammaConversion();
130  convModel = new G4PenelopeGammaConversionModel();
131  convModel->SetHighEnergyLimit(highEnergyLimit);
132  conv->AddEmModel(0, convModel);
133  pmanager->AddDiscreteProcess(conv);
134 
137  raylModel = new G4PenelopeRayleighModel();
138  raylModel->SetHighEnergyLimit(highEnergyLimit);
139  rayl->AddEmModel(0, raylModel);
140  pmanager->AddDiscreteProcess(rayl);
141 
142  } else if (particleName == "e-") {
143  //electron
144 
145  G4eIonisation* eIoni = new G4eIonisation();
147  eIoniModel = new G4PenelopeIonisationModel();
148  eIoniModel->SetHighEnergyLimit(highEnergyLimit);
149  eIoni->AddEmModel(0, eIoniModel, new G4UniversalFluctuation() );
150  pmanager->AddProcess(eIoni, -1,-1, 1);
151 
152  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
154  eBremModel = new G4PenelopeBremsstrahlungModel();
155  eBremModel->SetHighEnergyLimit(highEnergyLimit);
156  eBrem->AddEmModel(0, eBremModel);
157  pmanager->AddProcess(eBrem, -1,-1, 2);
158 
159  } else if (particleName == "e+") {
160  //positron
161  G4eIonisation* eIoni = new G4eIonisation();
163  eIoniModel = new G4PenelopeIonisationModel();
164  eIoniModel->SetHighEnergyLimit(highEnergyLimit);
165  eIoni->AddEmModel(0, eIoniModel, new G4UniversalFluctuation() );
166  pmanager->AddProcess(eIoni, -1,-1, 1);
167 
168  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
170  eBremModel = new G4PenelopeBremsstrahlungModel();
171  eBremModel->SetHighEnergyLimit(highEnergyLimit);
172  eBrem->AddEmModel(0, eBremModel);
173  pmanager->AddProcess(eBrem, -1,-1, 2);
174 
177  eAnniModel = new G4PenelopeAnnihilationModel();
178  eAnniModel->SetHighEnergyLimit(highEnergyLimit);
179  eAnni->AddEmModel(0, eAnniModel);
180  pmanager->AddProcess(eAnni, 0,-1, 3);
181 
182  } else if( particleName == "mu+" ||
183  particleName == "mu-" ) {
184  //muon
185  pmanager->AddProcess(new G4MuIonisation, -1,-1, 1);
186  pmanager->AddProcess(new G4MuBremsstrahlung, -1,-1, 2);
187  pmanager->AddProcess(new G4MuPairProduction, -1,-1, 3);
188 
189  } else if( particleName == "alpha" || particleName == "GenericIon" ) {
190  pmanager->AddProcess(new G4ionIonisation, -1,-1, 1);
191 
192  } else if ((!particle->IsShortLived()) &&
193  (particle->GetPDGCharge() != 0.0) &&
194  (particle->GetParticleName() != "chargedgeantino")) {
195  //all others charged particles except geantino
196  pmanager->AddProcess(new G4hIonisation, -1,-1, 1);
197  }
198  }
199 
200  // Deexcitation
201  //
203  de->SetFluo(true);
204  de->SetAuger(false);
205  de->SetPIXE(false);
207 }
208 
209 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
210 
static G4LossTableManager * Instance()
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
PhysListEmPenelope(const G4String &name="penelope")
const XML_Char * name
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
#define aParticleIterator
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
void SetAtomDeexcitation(G4VAtomDeexcitation *)