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 | // ------------------------------------------------------------------- |
---|
27 | // $Id: MicrobeamPhysicsList.cc,v 1.8 2009/04/30 10:23:57 sincerti Exp $ |
---|
28 | // ------------------------------------------------------------------- |
---|
29 | |
---|
30 | #include "G4ParticleDefinition.hh" |
---|
31 | #include "G4ProcessManager.hh" |
---|
32 | #include "G4ParticleTypes.hh" |
---|
33 | #include "G4StepLimiter.hh" |
---|
34 | #include "G4BaryonConstructor.hh" |
---|
35 | #include "G4IonConstructor.hh" |
---|
36 | #include "G4MesonConstructor.hh" |
---|
37 | |
---|
38 | #include "MicrobeamPhysicsList.hh" |
---|
39 | |
---|
40 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
41 | |
---|
42 | MicrobeamPhysicsList::MicrobeamPhysicsList(): G4VUserPhysicsList() |
---|
43 | { |
---|
44 | defaultCutValue = 0.01*micrometer; |
---|
45 | cutForGamma = defaultCutValue; |
---|
46 | cutForElectron = defaultCutValue; |
---|
47 | cutForPositron = defaultCutValue; |
---|
48 | SetVerboseLevel(1); |
---|
49 | } |
---|
50 | |
---|
51 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
52 | |
---|
53 | MicrobeamPhysicsList::~MicrobeamPhysicsList() |
---|
54 | {} |
---|
55 | |
---|
56 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
57 | |
---|
58 | void MicrobeamPhysicsList::ConstructParticle() |
---|
59 | { |
---|
60 | ConstructBosons(); |
---|
61 | ConstructLeptons(); |
---|
62 | ConstructBaryons(); |
---|
63 | } |
---|
64 | |
---|
65 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
66 | |
---|
67 | void MicrobeamPhysicsList::ConstructBosons() |
---|
68 | { |
---|
69 | // gamma |
---|
70 | G4Gamma::GammaDefinition(); |
---|
71 | |
---|
72 | // optical photon |
---|
73 | G4OpticalPhoton::OpticalPhotonDefinition(); |
---|
74 | } |
---|
75 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
76 | |
---|
77 | void MicrobeamPhysicsList::ConstructLeptons() |
---|
78 | { |
---|
79 | // leptons |
---|
80 | G4Electron::ElectronDefinition(); |
---|
81 | G4Positron::PositronDefinition(); |
---|
82 | } |
---|
83 | |
---|
84 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
85 | |
---|
86 | void MicrobeamPhysicsList::ConstructBaryons() |
---|
87 | { |
---|
88 | // baryons |
---|
89 | G4BaryonConstructor bConstructor; |
---|
90 | bConstructor.ConstructParticle(); |
---|
91 | |
---|
92 | G4IonConstructor iConstructor; |
---|
93 | iConstructor.ConstructParticle(); |
---|
94 | |
---|
95 | G4MesonConstructor mConstructor; |
---|
96 | mConstructor.ConstructParticle(); |
---|
97 | } |
---|
98 | |
---|
99 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
100 | |
---|
101 | void MicrobeamPhysicsList::ConstructProcess() |
---|
102 | { |
---|
103 | AddTransportation(); |
---|
104 | ConstructEM(); |
---|
105 | ConstructHad(); |
---|
106 | ConstructGeneral(); |
---|
107 | } |
---|
108 | |
---|
109 | // *** Processes and models |
---|
110 | |
---|
111 | // gamma |
---|
112 | |
---|
113 | #include "G4PhotoElectricEffect.hh" |
---|
114 | #include "G4LivermorePhotoElectricModel.hh" |
---|
115 | |
---|
116 | #include "G4ComptonScattering.hh" |
---|
117 | #include "G4LivermoreComptonModel.hh" |
---|
118 | |
---|
119 | #include "G4GammaConversion.hh" |
---|
120 | #include "G4LivermoreGammaConversionModel.hh" |
---|
121 | |
---|
122 | #include "G4RayleighScattering.hh" |
---|
123 | #include "G4LivermoreRayleighModel.hh" |
---|
124 | |
---|
125 | // e- |
---|
126 | |
---|
127 | #include "G4eMultipleScattering.hh" |
---|
128 | #include "G4UniversalFluctuation.hh" |
---|
129 | |
---|
130 | #include "G4eIonisation.hh" |
---|
131 | #include "G4LivermoreIonisationModel.hh" |
---|
132 | |
---|
133 | #include "G4eBremsstrahlung.hh" |
---|
134 | #include "G4LivermoreBremsstrahlungModel.hh" |
---|
135 | |
---|
136 | // e+ |
---|
137 | |
---|
138 | #include "G4eplusAnnihilation.hh" |
---|
139 | |
---|
140 | // mu |
---|
141 | |
---|
142 | #include "G4MuIonisation.hh" |
---|
143 | #include "G4MuBremsstrahlung.hh" |
---|
144 | #include "G4MuPairProduction.hh" |
---|
145 | |
---|
146 | // hadrons |
---|
147 | |
---|
148 | #include "G4hMultipleScattering.hh" |
---|
149 | #include "G4MscStepLimitType.hh" |
---|
150 | |
---|
151 | #include "G4hBremsstrahlung.hh" |
---|
152 | #include "G4hPairProduction.hh" |
---|
153 | |
---|
154 | #include "G4hIonisation.hh" |
---|
155 | #include "G4ionIonisation.hh" |
---|
156 | #include "G4IonParametrisedLossModel.hh" |
---|
157 | |
---|
158 | // |
---|
159 | |
---|
160 | #include "G4LossTableManager.hh" |
---|
161 | #include "G4EmProcessOptions.hh" |
---|
162 | |
---|
163 | |
---|
164 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
165 | |
---|
166 | void MicrobeamPhysicsList::ConstructEM() |
---|
167 | { |
---|
168 | theParticleIterator->reset(); |
---|
169 | |
---|
170 | while( (*theParticleIterator)() ){ |
---|
171 | |
---|
172 | G4ParticleDefinition* particle = theParticleIterator->value(); |
---|
173 | G4ProcessManager* pmanager = particle->GetProcessManager(); |
---|
174 | G4String particleName = particle->GetParticleName(); |
---|
175 | |
---|
176 | if (particleName == "gamma") { |
---|
177 | |
---|
178 | |
---|
179 | G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect(); |
---|
180 | G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel = new G4LivermorePhotoElectricModel(); |
---|
181 | thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel); |
---|
182 | pmanager->AddDiscreteProcess(thePhotoElectricEffect); |
---|
183 | |
---|
184 | G4ComptonScattering* theComptonScattering = new G4ComptonScattering(); |
---|
185 | G4LivermoreComptonModel* theLivermoreComptonModel = new G4LivermoreComptonModel(); |
---|
186 | theComptonScattering->AddEmModel(0, theLivermoreComptonModel); |
---|
187 | pmanager->AddDiscreteProcess(theComptonScattering); |
---|
188 | |
---|
189 | G4GammaConversion* theGammaConversion = new G4GammaConversion(); |
---|
190 | G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel = new G4LivermoreGammaConversionModel(); |
---|
191 | theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel); |
---|
192 | pmanager->AddDiscreteProcess(theGammaConversion); |
---|
193 | |
---|
194 | G4RayleighScattering* theRayleigh = new G4RayleighScattering(); |
---|
195 | G4LivermoreRayleighModel* theRayleighModel = new G4LivermoreRayleighModel(); |
---|
196 | theRayleigh->AddEmModel(0, theRayleighModel); |
---|
197 | pmanager->AddDiscreteProcess(theRayleigh); |
---|
198 | |
---|
199 | pmanager->AddProcess(new G4StepLimiter(), -1, -1, 5); |
---|
200 | |
---|
201 | } else if (particleName == "e-") { |
---|
202 | |
---|
203 | G4eMultipleScattering* msc = new G4eMultipleScattering(); |
---|
204 | msc->SetStepLimitType(fUseDistanceToBoundary); |
---|
205 | pmanager->AddProcess(msc, -1, 1, 1); |
---|
206 | |
---|
207 | // Ionisation |
---|
208 | G4eIonisation* eIoni = new G4eIonisation(); |
---|
209 | eIoni->AddEmModel(0, new G4LivermoreIonisationModel(), new G4UniversalFluctuation() ); |
---|
210 | eIoni->SetStepFunction(0.2, 100*um); // |
---|
211 | pmanager->AddProcess(eIoni, -1, 2, 2); |
---|
212 | |
---|
213 | // Bremsstrahlung |
---|
214 | G4eBremsstrahlung* eBrem = new G4eBremsstrahlung(); |
---|
215 | eBrem->AddEmModel(0, new G4LivermoreBremsstrahlungModel()); |
---|
216 | pmanager->AddProcess(eBrem, -1,-3, 3); |
---|
217 | |
---|
218 | pmanager->AddProcess(new G4StepLimiter(), -1, -1, 4); |
---|
219 | |
---|
220 | } else if (particleName == "e+") { |
---|
221 | |
---|
222 | // Identical to G4EmStandardPhysics_option3 |
---|
223 | |
---|
224 | G4eMultipleScattering* msc = new G4eMultipleScattering(); |
---|
225 | msc->SetStepLimitType(fUseDistanceToBoundary); |
---|
226 | pmanager->AddProcess(msc, -1, 1, 1); |
---|
227 | |
---|
228 | G4eIonisation* eIoni = new G4eIonisation(); |
---|
229 | eIoni->SetStepFunction(0.2, 100*um); |
---|
230 | pmanager->AddProcess(eIoni, -1, 2, 2); |
---|
231 | |
---|
232 | pmanager->AddProcess(new G4eBremsstrahlung, -1,-3, 3); |
---|
233 | |
---|
234 | pmanager->AddProcess(new G4eplusAnnihilation,0,-1, 4); |
---|
235 | |
---|
236 | pmanager->AddProcess(new G4StepLimiter(), -1, -1, 5); |
---|
237 | |
---|
238 | } else if (particleName == "GenericIon") { |
---|
239 | |
---|
240 | pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1); |
---|
241 | |
---|
242 | G4ionIonisation* ionIoni = new G4ionIonisation(); |
---|
243 | ionIoni->SetEmModel(new G4IonParametrisedLossModel()); |
---|
244 | ionIoni->SetStepFunction(0.1, 20*um); |
---|
245 | pmanager->AddProcess(ionIoni, -1, 2, 2); |
---|
246 | |
---|
247 | pmanager->AddProcess(new G4StepLimiter(), -1, -1, 3); |
---|
248 | |
---|
249 | } else if (particleName == "alpha" || |
---|
250 | particleName == "He3" ) { |
---|
251 | |
---|
252 | // Identical to G4EmStandardPhysics_option3 |
---|
253 | |
---|
254 | pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1); |
---|
255 | |
---|
256 | G4ionIonisation* ionIoni = new G4ionIonisation(); |
---|
257 | ionIoni->SetStepFunction(0.1, 20*um); |
---|
258 | pmanager->AddProcess(ionIoni, -1, 2, 2); |
---|
259 | |
---|
260 | pmanager->AddProcess(new G4StepLimiter(), -1, -1, 3); |
---|
261 | } |
---|
262 | |
---|
263 | // |
---|
264 | |
---|
265 | } |
---|
266 | } |
---|
267 | |
---|
268 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
269 | |
---|
270 | #include "G4HadronElasticProcess.hh" |
---|
271 | #include "G4LElastic.hh" |
---|
272 | |
---|
273 | #include "G4AlphaInelasticProcess.hh" |
---|
274 | #include "G4BinaryLightIonReaction.hh" |
---|
275 | #include "G4TripathiCrossSection.hh" |
---|
276 | #include "G4IonsShenCrossSection.hh" |
---|
277 | #include "G4LEAlphaInelastic.hh" |
---|
278 | |
---|
279 | void MicrobeamPhysicsList::ConstructHad() |
---|
280 | { |
---|
281 | |
---|
282 | G4HadronElasticProcess * theElasticProcess = new G4HadronElasticProcess; |
---|
283 | theElasticProcess->RegisterMe( new G4LElastic() ); |
---|
284 | |
---|
285 | theParticleIterator->reset(); |
---|
286 | while( (*theParticleIterator)() ) |
---|
287 | { |
---|
288 | G4ParticleDefinition* particle = theParticleIterator->value(); |
---|
289 | G4ProcessManager* pManager = particle->GetProcessManager(); |
---|
290 | |
---|
291 | if (particle->GetParticleName() == "alpha") |
---|
292 | { |
---|
293 | |
---|
294 | // INELASTIC SCATTERING |
---|
295 | // Binary Cascade |
---|
296 | G4BinaryLightIonReaction* theBC = new G4BinaryLightIonReaction(); |
---|
297 | theBC -> SetMinEnergy(80.*MeV); |
---|
298 | theBC -> SetMaxEnergy(40.*GeV); |
---|
299 | |
---|
300 | // TRIPATHI CROSS SECTION |
---|
301 | // Implementation of formulas in analogy to NASA technical paper 3621 by |
---|
302 | // Tripathi, et al. Cross-sections for ion ion scattering |
---|
303 | G4TripathiCrossSection* TripathiCrossSection = new G4TripathiCrossSection; |
---|
304 | |
---|
305 | // IONS SHEN CROSS SECTION |
---|
306 | // Implementation of formulas |
---|
307 | // Shen et al. Nuc. Phys. A 491 130 (1989) |
---|
308 | // Total Reaction Cross Section for Heavy-Ion Collisions |
---|
309 | G4IonsShenCrossSection* aShen = new G4IonsShenCrossSection; |
---|
310 | |
---|
311 | // Final state production model for Alpha inelastic scattering below 20 GeV |
---|
312 | G4LEAlphaInelastic* theAIModel = new G4LEAlphaInelastic; |
---|
313 | theAIModel -> SetMaxEnergy(100.*MeV); |
---|
314 | |
---|
315 | G4AlphaInelasticProcess * theIPalpha = new G4AlphaInelasticProcess; |
---|
316 | theIPalpha->AddDataSet(TripathiCrossSection); |
---|
317 | theIPalpha->AddDataSet(aShen); |
---|
318 | |
---|
319 | // Register the Alpha Inelastic and Binary Cascade Model |
---|
320 | theIPalpha->RegisterMe(theAIModel); |
---|
321 | theIPalpha->RegisterMe(theBC); |
---|
322 | |
---|
323 | // Activate the alpha inelastic scattering using the alpha inelastic and binary cascade model |
---|
324 | pManager -> AddDiscreteProcess(theIPalpha); |
---|
325 | |
---|
326 | // Activate the Hadron Elastic Process |
---|
327 | pManager -> AddDiscreteProcess(theElasticProcess); |
---|
328 | |
---|
329 | } |
---|
330 | } |
---|
331 | |
---|
332 | } |
---|
333 | |
---|
334 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
335 | |
---|
336 | void MicrobeamPhysicsList::ConstructGeneral() |
---|
337 | { } |
---|
338 | |
---|
339 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
340 | |
---|
341 | void MicrobeamPhysicsList::SetCuts() |
---|
342 | { |
---|
343 | if (verboseLevel >0){ |
---|
344 | G4cout << "MicrobeamPhysicsList::SetCuts:"; |
---|
345 | G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl; |
---|
346 | } |
---|
347 | |
---|
348 | // set cut values for gamma at first and for e- second and next for e+, |
---|
349 | // because some processes for e+/e- need cut values for gamma |
---|
350 | SetCutValue(cutForGamma, "gamma"); |
---|
351 | SetCutValue(cutForElectron, "e-"); |
---|
352 | SetCutValue(cutForPositron, "e+"); |
---|
353 | |
---|
354 | if (verboseLevel>0) DumpCutValuesTable(); |
---|
355 | |
---|
356 | } |
---|
357 | |
---|