Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
electromagnetic/TestEm12/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/TestEm12/src/PhysListEmStandardSS.cc
27 /// \brief Implementation of the PhysListEmStandardSS class
28 //
29 // $Id: PhysListEmStandardSS.cc 77092 2013-11-21 10:50:40Z gcosmo $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "PhysListEmStandardSS.hh"
35 
36 #include "G4ParticleDefinition.hh"
37 #include "G4ProcessManager.hh"
38 
39 #include "G4ComptonScattering.hh"
40 #include "G4GammaConversion.hh"
41 #include "G4PhotoElectricEffect.hh"
42 
43 #include "G4CoulombScattering.hh"
46 
47 #include "G4eIonisation.hh"
48 #include "G4eBremsstrahlung.hh"
49 #include "G4eplusAnnihilation.hh"
50 
51 #include "G4MuIonisation.hh"
52 #include "G4MuBremsstrahlung.hh"
53 #include "G4MuPairProduction.hh"
54 
55 #include "G4hIonisation.hh"
56 #include "G4ionIonisation.hh"
57 
58 #include "G4EmProcessOptions.hh"
59 
60 #include "G4SystemOfUnits.hh"
61 
62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
63 
65  : G4VPhysicsConstructor(name)
66 {}
67 
68 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69 
71 {}
72 
73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
74 
76 {
77  // Add standard EM Processes
78 
79  aParticleIterator->reset();
80  while( (*aParticleIterator)() ){
81  G4ParticleDefinition* particle = aParticleIterator->value();
82  G4ProcessManager* pmanager = particle->GetProcessManager();
83  G4String particleName = particle->GetParticleName();
84 
85  if (particleName == "gamma") {
86  // gamma
90 
91  } else if (particleName == "e-") {
92  //electron
93  pmanager->AddProcess(new G4eIonisation, -1, 1, 1);
94  pmanager->AddProcess(new G4eBremsstrahlung, -1, 2, 2);
95 
99  model->SetLowEnergyThreshold(10*eV);
100  cs->SetEmModel(model, 1);
101  pmanager->AddDiscreteProcess(cs);
102 
103  } else if (particleName == "e+") {
104  //positron
105  pmanager->AddProcess(new G4eIonisation, -1, 1, 1);
106  pmanager->AddProcess(new G4eBremsstrahlung, -1, 2, 2);
107  pmanager->AddProcess(new G4eplusAnnihilation, 0,-1, 3);
108  pmanager->AddDiscreteProcess(new G4CoulombScattering);
109 
110  } else if (particleName == "mu+" ||
111  particleName == "mu-" ) {
112  //muon
113  pmanager->AddProcess(new G4MuIonisation, -1, 1, 1);
114  pmanager->AddProcess(new G4MuBremsstrahlung, -1, 2, 2);
115  pmanager->AddProcess(new G4MuPairProduction, -1, 3, 3);
116  pmanager->AddDiscreteProcess(new G4CoulombScattering);
117 
118  } else if (particleName == "alpha" || particleName == "He3") {
119  pmanager->AddProcess(new G4ionIonisation, -1, 1, 1);
121  //cs->AddEmModel(0, new G4IonCoulombScatteringModel());
122  cs->SetBuildTableFlag(false);
123  pmanager->AddDiscreteProcess(cs);
124 
125  } else if (particleName == "GenericIon" ) {
126  pmanager->AddProcess(new G4ionIonisation, -1, 1, 1);
128  //cs->AddEmModel(0, new G4IonCoulombScatteringModel());
129  cs->SetBuildTableFlag(false);
130  pmanager->AddDiscreteProcess(cs);
131 
132  } else if ((!particle->IsShortLived()) &&
133  (particle->GetPDGCharge() != 0.0) &&
134  (particle->GetParticleName() != "chargedgeantino")) {
135  //all others charged particles except geantino
137  pmanager->AddProcess(new G4hIonisation, -1, 1, 1);
138  cs->SetBuildTableFlag(false);
139  pmanager->AddDiscreteProcess(cs);
140  }
141  }
142 
143  // Em options
144  //
145  // Main options and setting parameters are shown here.
146  // Several of them have default values.
147  //
148  G4EmProcessOptions emOptions;
149 
150  //physics tables
151  //
152  emOptions.SetMinEnergy(100*eV); //default
153  emOptions.SetMaxEnergy(100*TeV); //default
154  emOptions.SetDEDXBinning(12*20); //default=12*7
155  emOptions.SetLambdaBinning(12*20); //default=12*7
156  // emOptions.SetSplineFlag(true); //default
157 
158  //energy loss
159  //
160  //emOptions.SetStepFunction(0.2, 100*um); //default=(0.2, 1*mm)
161  //emOptions.SetLinearLossLimit(1.e-2); //default
162 
163  //ionization
164  //
165  emOptions.SetSubCutoff(false); //default
166 
167  // scattering
168  emOptions.SetPolarAngleLimit(0.0);
169 }
170 
171 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
172 
PhysListEmStandardSS(const G4String &name="standardSS")
void SetBuildTableFlag(G4bool val)
void SetMinEnergy(G4double val)
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
const XML_Char * name
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
void SetDEDXBinning(G4int val)
void SetEmModel(G4VEmModel *, G4int index=1)
void SetLambdaBinning(G4int val)
#define aParticleIterator
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
const XML_Char XML_Content * model
void SetMaxEnergy(G4double val)
G4double GetPDGCharge() const
void SetSubCutoff(G4bool val, const G4Region *r=0)
void SetPolarAngleLimit(G4double val)