Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
src/PhysListEmStandard_SS.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 medical/fanoCavity/src/PhysListEmStandard_SS.cc
27 /// \brief Implementation of the PhysListEmStandard_SS class
28 //
29 // $Id: PhysListEmStandard_SS.cc 68525 2013-04-01 21:16:50Z adotti $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "PhysListEmStandard_SS.hh"
35 #include "DetectorConstruction.hh"
36 
37 #include "G4ParticleDefinition.hh"
38 #include "G4ProcessManager.hh"
39 
40 #include "G4ComptonScattering.hh"
41 #include "MyKleinNishinaCompton.hh"
42 #include "G4GammaConversion.hh"
43 #include "G4PhotoElectricEffect.hh"
44 
45 #include "G4CoulombScattering.hh"
46 
47 #include "G4eIonisation.hh"
48 #include "MyMollerBhabhaModel.hh"
49 #include "G4eBremsstrahlung.hh"
50 #include "G4eplusAnnihilation.hh"
51 
52 #include "G4hIonisation.hh"
53 
54 #include "G4EmProcessOptions.hh"
55 #include "G4MscStepLimitType.hh"
56 
57 #include "G4SystemOfUnits.hh"
58 
59 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60 
63 : G4VPhysicsConstructor(name), fDetector(det)
64 {}
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67 
69 {}
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72 
74 {
75  // Add standard EM Processes
76  //
77 
78  aParticleIterator->reset();
79  while( (*aParticleIterator)() ){
80  G4ParticleDefinition* particle = aParticleIterator->value();
81  G4ProcessManager* pmanager = particle->GetProcessManager();
82  G4String particleName = particle->GetParticleName();
83 
84  if (particleName == "gamma") {
85  // gamma
86 
88  MyKleinNishinaCompton* comptonModel = new MyKleinNishinaCompton(fDetector);
89  comptonModel->SetCSFactor(1000.);
90  compton->SetEmModel(comptonModel );
91 
93  pmanager->AddDiscreteProcess(compton);
95 
96  } else if (particleName == "e-") {
97  //electron
98 
99  G4eIonisation* eIoni = new G4eIonisation();
100  eIoni->SetEmModel(new MyMollerBhabhaModel);
101 
102  pmanager->AddProcess(new G4CoulombScattering, -1, -1, 1);
103  pmanager->AddProcess(eIoni, -1, 1, 2);
104 /// pmanager->AddProcess(new G4eBremsstrahlung, -1, 2, 3);
105 
106  } else if (particleName == "e+") {
107  //positron
108 
109  G4eIonisation* pIoni = new G4eIonisation();
110  pIoni->SetEmModel(new MyMollerBhabhaModel);
111 
112  pmanager->AddProcess(new G4CoulombScattering, -1, -1, 1);
113  pmanager->AddProcess(pIoni, -1, 1, 2);
114 /// pmanager->AddProcess(new G4eBremsstrahlung, -1, 2, 3);
115  pmanager->AddProcess(new G4eplusAnnihilation, 0,- 1, 3);
116 
117  } else if( particleName == "proton" ) {
118  //proton
119  pmanager->AddProcess(new G4CoulombScattering, -1, -1, 1);
120  pmanager->AddProcess(new G4hIonisation, -1, 1, 2);
121  }
122  }
123 
124  // Em options
125  //
126  // Main options and setting parameters are shown here.
127  // Several of them have default values.
128  //
129  G4EmProcessOptions emOptions;
130 
131  //physics tables
132  //
133  emOptions.SetMinEnergy(100*eV); //default
134  emOptions.SetMaxEnergy(10*GeV); //default
135  emOptions.SetDEDXBinning(8*20); //default=8*7
136  emOptions.SetLambdaBinning(8*20); //default=8*7
137  emOptions.SetSplineFlag(true); //default
138 
139  //multiple coulomb scattering
140  //
141  emOptions.SetPolarAngleLimit(0.0);
142 
143  //energy loss
144  //
145  emOptions.SetStepFunction(0.2, 10*um); //default=(0.2, 1*mm)
146 
147  //build CSDA range
148  //
149  emOptions.SetBuildCSDARange(true); //default=false
150  emOptions.SetMaxEnergyForCSDARange(10*GeV);
151  emOptions.SetDEDXBinningForCSDARange(8*20); //default=8*7
152 
153  //ionization
154  //
155  emOptions.SetSubCutoff(false); //default
156 }
157 
158 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
159 
void SetCSFactor(G4double factor)
void SetSplineFlag(G4bool val)
void SetMinEnergy(G4double val)
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
const XML_Char * name
void SetStepFunction(G4double v1, G4double v2)
void SetDEDXBinningForCSDARange(G4int val)
Definition of the MyKleinNishinaCompton class.
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
void SetMaxEnergyForCSDARange(G4double val)
PhysListEmStandard_SS(const G4String &name, DetectorConstruction *det)
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)
void SetMaxEnergy(G4double val)
void SetEmModel(G4VEmModel *, G4int index=1)
void SetBuildCSDARange(G4bool val)
void SetSubCutoff(G4bool val, const G4Region *r=0)
void SetPolarAngleLimit(G4double val)