Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
TestEm13/src/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 //
26 /// \file electromagnetic/TestEm13/src/PhysListEmLivermore.cc
27 /// \brief Implementation of the PhysListEmLivermore class
28 //
29 //
30 // $Id: PhysListEmLivermore.cc 68585 2013-04-01 23:35:07Z adotti $
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 "G4SystemOfUnits.hh"
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80 
82  : G4VPhysicsConstructor(name)
83 { }
84 
85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
86 
88 { }
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
91 
93 {
94  // Add standard EM Processes
95 
96  aParticleIterator->reset();
97  while( (*aParticleIterator)() ){
98  G4ParticleDefinition* particle = aParticleIterator->value();
99  G4ProcessManager* pmanager = particle->GetProcessManager();
100  G4String particleName = particle->GetParticleName();
101 
102  //Applicability range for Livermore models
103  //for higher energies, the Standard models are used
104  G4double highEnergyLimit = 1*GeV;
105 
106  if (particleName == "gamma") {
107  // gamma
108 
111  photModel = new G4LivermorePhotoElectricModel();
112  photModel->SetHighEnergyLimit(highEnergyLimit);
113  phot->AddEmModel(0, photModel);
114  pmanager->AddDiscreteProcess(phot);
115 
118  comptModel = new G4LivermoreComptonModel();
119  comptModel->SetHighEnergyLimit(highEnergyLimit);
120  compt->AddEmModel(0, comptModel);
121  pmanager->AddDiscreteProcess(compt);
122 
123  G4GammaConversion* conv = new G4GammaConversion();
125  convModel = new G4LivermoreGammaConversionModel();
126  convModel->SetHighEnergyLimit(highEnergyLimit);
127  conv->AddEmModel(0, convModel);
128  pmanager->AddDiscreteProcess(conv);
129 
132  raylModel = new G4LivermoreRayleighModel();
133  raylModel->SetHighEnergyLimit(highEnergyLimit);
134  rayl->AddEmModel(0, raylModel);
135  pmanager->AddDiscreteProcess(rayl);
136 
137  } else if (particleName == "e-") {
138  //electron
139 
140  G4eIonisation* eIoni = new G4eIonisation();
142  eIoniModel = new G4LivermoreIonisationModel();
143  eIoniModel->SetHighEnergyLimit(highEnergyLimit);
144  eIoni->AddEmModel(0, eIoniModel, new G4UniversalFluctuation() );
145  pmanager->AddProcess(eIoni, -1,-1, 1);
146 
147  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
149  eBremModel = new G4LivermoreBremsstrahlungModel();
150  eBremModel->SetHighEnergyLimit(highEnergyLimit);
151  eBrem->AddEmModel(0, eBremModel);
152  pmanager->AddProcess(eBrem, -1,-1, 2);
153 
154  } else if (particleName == "e+") {
155  //positron
156  pmanager->AddProcess(new G4eIonisation, -1,-1, 1);
157  pmanager->AddProcess(new G4eBremsstrahlung, -1,-1, 2);
158  pmanager->AddProcess(new G4eplusAnnihilation, 0,-1, 3);
159 
160  } else if( particleName == "mu+" ||
161  particleName == "mu-" ) {
162  //muon
163  pmanager->AddProcess(new G4MuIonisation, -1,-1, 1);
164  pmanager->AddProcess(new G4MuBremsstrahlung, -1,-1, 2);
165  pmanager->AddProcess(new G4MuPairProduction, -1,-1, 3);
166 
167  } else if( particleName == "alpha" || particleName == "GenericIon" ) {
168  pmanager->AddProcess(new G4ionIonisation, -1,-1, 1);
169 
170  } else if ((!particle->IsShortLived()) &&
171  (particle->GetPDGCharge() != 0.0) &&
172  (particle->GetParticleName() != "chargedgeantino")) {
173  //all others charged particles except geantino
174  pmanager->AddProcess(new G4hIonisation, -1,-1, 1);
175  }
176  }
177 }
178 
179 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
180 
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
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
PhysListEmLivermore(const G4String &name="livermore")