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 | // $Id: G4EmDNAPhysics.cc,v 1.2 2009/11/01 13:21:13 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-03-cand-03 $ |
---|
28 | |
---|
29 | #include "G4EmDNAPhysics.hh" |
---|
30 | |
---|
31 | #include "G4ProcessManager.hh" |
---|
32 | #include "G4DNAGenericIonsManager.hh" |
---|
33 | |
---|
34 | // *** Processes and models for Geant4-DNA |
---|
35 | |
---|
36 | #include "G4DNAElastic.hh" |
---|
37 | #include "G4DNAChampionElasticModel.hh" |
---|
38 | #include "G4DNAScreenedRutherfordElasticModel.hh" |
---|
39 | |
---|
40 | #include "G4DNAExcitation.hh" |
---|
41 | #include "G4DNAEmfietzoglouExcitationModel.hh" |
---|
42 | #include "G4DNAMillerGreenExcitationModel.hh" |
---|
43 | #include "G4DNABornExcitationModel.hh" |
---|
44 | |
---|
45 | #include "G4DNAIonisation.hh" |
---|
46 | #include "G4DNABornIonisationModel.hh" |
---|
47 | #include "G4DNARuddIonisationModel.hh" |
---|
48 | |
---|
49 | #include "G4DNAChargeDecrease.hh" |
---|
50 | #include "G4DNADingfelderChargeDecreaseModel.hh" |
---|
51 | |
---|
52 | #include "G4DNAChargeIncrease.hh" |
---|
53 | #include "G4DNADingfelderChargeIncreaseModel.hh" |
---|
54 | |
---|
55 | // particles |
---|
56 | |
---|
57 | #include "G4Electron.hh" |
---|
58 | #include "G4Proton.hh" |
---|
59 | |
---|
60 | // Warning : the following is needed in order to use EM Physics builders |
---|
61 | // e+ |
---|
62 | #include "G4Positron.hh" |
---|
63 | #include "G4eMultipleScattering.hh" |
---|
64 | #include "G4eIonisation.hh" |
---|
65 | #include "G4eBremsstrahlung.hh" |
---|
66 | #include "G4eplusAnnihilation.hh" |
---|
67 | // gamma |
---|
68 | #include "G4Gamma.hh" |
---|
69 | #include "G4PhotoElectricEffect.hh" |
---|
70 | #include "G4LivermorePhotoElectricModel.hh" |
---|
71 | #include "G4ComptonScattering.hh" |
---|
72 | #include "G4LivermoreComptonModel.hh" |
---|
73 | #include "G4GammaConversion.hh" |
---|
74 | #include "G4LivermoreGammaConversionModel.hh" |
---|
75 | #include "G4RayleighScattering.hh" |
---|
76 | #include "G4LivermoreRayleighModel.hh" |
---|
77 | // end of warning |
---|
78 | |
---|
79 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
80 | |
---|
81 | G4EmDNAPhysics::G4EmDNAPhysics( |
---|
82 | G4int ver, const G4String& name) |
---|
83 | : G4VPhysicsConstructor(name), verbose(ver) |
---|
84 | {} |
---|
85 | |
---|
86 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
87 | |
---|
88 | G4EmDNAPhysics::~G4EmDNAPhysics() |
---|
89 | {} |
---|
90 | |
---|
91 | #include "G4GenericIon.hh" |
---|
92 | |
---|
93 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
94 | |
---|
95 | void G4EmDNAPhysics::ConstructParticle() |
---|
96 | { |
---|
97 | // bosons |
---|
98 | G4Gamma::Gamma(); |
---|
99 | |
---|
100 | // leptons |
---|
101 | G4Electron::Electron(); |
---|
102 | G4Positron::Positron(); |
---|
103 | |
---|
104 | // baryons |
---|
105 | G4Proton::Proton(); |
---|
106 | |
---|
107 | G4GenericIon::GenericIonDefinition(); |
---|
108 | |
---|
109 | G4DNAGenericIonsManager * genericIonsManager; |
---|
110 | genericIonsManager=G4DNAGenericIonsManager::Instance(); |
---|
111 | genericIonsManager->GetIon("alpha++"); |
---|
112 | genericIonsManager->GetIon("alpha+"); |
---|
113 | genericIonsManager->GetIon("helium"); |
---|
114 | genericIonsManager->GetIon("hydrogen"); |
---|
115 | |
---|
116 | } |
---|
117 | |
---|
118 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
119 | |
---|
120 | void G4EmDNAPhysics::ConstructProcess() |
---|
121 | { |
---|
122 | theParticleIterator->reset(); |
---|
123 | while( (*theParticleIterator)() ) |
---|
124 | { |
---|
125 | G4ParticleDefinition* particle = theParticleIterator->value(); |
---|
126 | G4ProcessManager* pmanager = particle->GetProcessManager(); |
---|
127 | G4String particleName = particle->GetParticleName(); |
---|
128 | |
---|
129 | if (particleName == "e-") { |
---|
130 | |
---|
131 | // *** Elastic scattering (two alternative models available) *** |
---|
132 | |
---|
133 | G4DNAElastic* theDNAElasticProcess = new G4DNAElastic(); |
---|
134 | theDNAElasticProcess->SetModel(new G4DNAScreenedRutherfordElasticModel()); |
---|
135 | |
---|
136 | // or alternative model |
---|
137 | // theDNAElasticProcess->SetModel(new G4DNAChampionElasticModel()); |
---|
138 | |
---|
139 | pmanager->AddDiscreteProcess(theDNAElasticProcess); |
---|
140 | |
---|
141 | // *** Excitation *** |
---|
142 | pmanager->AddDiscreteProcess(new G4DNAExcitation()); |
---|
143 | |
---|
144 | // *** Ionisation *** |
---|
145 | pmanager->AddDiscreteProcess(new G4DNAIonisation()); |
---|
146 | |
---|
147 | |
---|
148 | } else if ( particleName == "proton" ) { |
---|
149 | pmanager->AddDiscreteProcess(new G4DNAExcitation()); |
---|
150 | pmanager->AddDiscreteProcess(new G4DNAIonisation()); |
---|
151 | pmanager->AddDiscreteProcess(new G4DNAChargeDecrease()); |
---|
152 | |
---|
153 | } else if ( particleName == "hydrogen" ) { |
---|
154 | pmanager->AddDiscreteProcess(new G4DNAIonisation()); |
---|
155 | pmanager->AddDiscreteProcess(new G4DNAChargeIncrease()); |
---|
156 | |
---|
157 | } else if ( particleName == "alpha" ) { |
---|
158 | pmanager->AddDiscreteProcess(new G4DNAExcitation()); |
---|
159 | pmanager->AddDiscreteProcess(new G4DNAIonisation()); |
---|
160 | pmanager->AddDiscreteProcess(new G4DNAChargeDecrease()); |
---|
161 | |
---|
162 | } else if ( particleName == "alpha+" ) { |
---|
163 | pmanager->AddDiscreteProcess(new G4DNAExcitation()); |
---|
164 | pmanager->AddDiscreteProcess(new G4DNAIonisation()); |
---|
165 | pmanager->AddDiscreteProcess(new G4DNAChargeDecrease()); |
---|
166 | pmanager->AddDiscreteProcess(new G4DNAChargeIncrease()); |
---|
167 | |
---|
168 | } else if ( particleName == "helium" ) { |
---|
169 | pmanager->AddDiscreteProcess(new G4DNAExcitation()); |
---|
170 | pmanager->AddDiscreteProcess(new G4DNAIonisation()); |
---|
171 | pmanager->AddDiscreteProcess(new G4DNAChargeIncrease()); |
---|
172 | |
---|
173 | } |
---|
174 | |
---|
175 | // Warning : the following particles and processes are needed by EM Physics builders |
---|
176 | // They are taken from the default Livermore Physics list |
---|
177 | // These particles are currently not handled by Geant4-DNA |
---|
178 | |
---|
179 | // e+ |
---|
180 | |
---|
181 | else if (particleName == "e+") { |
---|
182 | |
---|
183 | // Identical to G4EmStandardPhysics_option3 |
---|
184 | |
---|
185 | G4eMultipleScattering* msc = new G4eMultipleScattering(); |
---|
186 | msc->SetStepLimitType(fUseDistanceToBoundary); |
---|
187 | pmanager->AddProcess(msc, -1, 1, 1); |
---|
188 | |
---|
189 | G4eIonisation* eIoni = new G4eIonisation(); |
---|
190 | eIoni->SetStepFunction(0.2, 100*um); |
---|
191 | |
---|
192 | pmanager->AddProcess(eIoni, -1, 2, 2); |
---|
193 | pmanager->AddProcess(new G4eBremsstrahlung, -1,-3, 3); |
---|
194 | pmanager->AddProcess(new G4eplusAnnihilation,0,-1, 4); |
---|
195 | |
---|
196 | } else if (particleName == "gamma") { |
---|
197 | |
---|
198 | G4double LivermoreHighEnergyLimit = GeV; |
---|
199 | |
---|
200 | G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect(); |
---|
201 | G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel = |
---|
202 | new G4LivermorePhotoElectricModel(); |
---|
203 | theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit); |
---|
204 | thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel); |
---|
205 | pmanager->AddDiscreteProcess(thePhotoElectricEffect); |
---|
206 | |
---|
207 | G4ComptonScattering* theComptonScattering = new G4ComptonScattering(); |
---|
208 | G4LivermoreComptonModel* theLivermoreComptonModel = |
---|
209 | new G4LivermoreComptonModel(); |
---|
210 | theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit); |
---|
211 | theComptonScattering->AddEmModel(0, theLivermoreComptonModel); |
---|
212 | pmanager->AddDiscreteProcess(theComptonScattering); |
---|
213 | |
---|
214 | G4GammaConversion* theGammaConversion = new G4GammaConversion(); |
---|
215 | G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel = |
---|
216 | new G4LivermoreGammaConversionModel(); |
---|
217 | theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit); |
---|
218 | theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel); |
---|
219 | pmanager->AddDiscreteProcess(theGammaConversion); |
---|
220 | |
---|
221 | G4RayleighScattering* theRayleigh = new G4RayleighScattering(); |
---|
222 | G4LivermoreRayleighModel* theRayleighModel = new G4LivermoreRayleighModel(); |
---|
223 | theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit); |
---|
224 | theRayleigh->AddEmModel(0, theRayleighModel); |
---|
225 | pmanager->AddDiscreteProcess(theRayleigh); |
---|
226 | } |
---|
227 | |
---|
228 | // Warning : end of particles and processes are needed by EM Physics builders |
---|
229 | |
---|
230 | } |
---|
231 | } |
---|
232 | |
---|
233 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|