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 | // ------------------------------------------------------------- |
---|
28 | // GEANT 4 class implementation file |
---|
29 | // |
---|
30 | // History: based on object model of |
---|
31 | // 2nd December 1995, G.Cosmo |
---|
32 | // ---------- G4hLowEnergyIonisation physics process ------- |
---|
33 | // by Vladimir Ivanchenko, 14 July 1999 |
---|
34 | // was made on the base of G4hIonisation class |
---|
35 | // developed by Laszlo Urban |
---|
36 | // ************************************************************ |
---|
37 | // It is the extention of the ionisation process for the slow |
---|
38 | // charged hadrons. |
---|
39 | // ************************************************************ |
---|
40 | // 28 July 1999 V.Ivanchenko cleen up |
---|
41 | // 17 August 1999 G.Mancinelli added ICRU parametrisations for protons |
---|
42 | // 20 August 1999 G.Mancinelli added ICRU tables for alpha |
---|
43 | // 31 August 1999 V.Ivanchenko update and cleen up |
---|
44 | // 30 Sept. 1999 V.Ivanchenko minor upgrade |
---|
45 | // 12 Dec. 1999 S. Chauvie added Barkas correction |
---|
46 | // 19 Jan. 2000 V.Ivanchenko minor changing in Barkas corrections |
---|
47 | // 02 April 2000 S. Chauvie linearization of Barkas effect |
---|
48 | // 03 April 2000 V.Ivanchenko Nuclear Stopping power for antiprotons |
---|
49 | // 23 May 2000 MG Pia Clean up for QAO model |
---|
50 | // 24 May 2000 MG Pia Code properly indented to improve legibility |
---|
51 | // 17 July 2000 V.Ivanchenko Bug in scaling AlongStepDoIt method |
---|
52 | // 25 July 2000 V.Ivanchenko New design iteration |
---|
53 | // 17 August 2000 V.Ivanchenko Add ion fluctuation models |
---|
54 | // 18 August 2000 V.Ivanchenko Bug fixed in GetConstrain |
---|
55 | // 22 August 2000 V.Ivanchenko Insert paramStepLimit and |
---|
56 | // reorganise access to Barkas and Bloch terms |
---|
57 | // 04 Sept. 2000 V.Ivanchenko rename fluctuations |
---|
58 | // 05 Sept. 2000 V.Ivanchenko clean up |
---|
59 | // 03 Oct. 2000 V.Ivanchenko CodeWizard clean up |
---|
60 | // 03 Nov. 2000 V.Ivanchenko MinKineticEnergy=LowestKineticEnergy=10eV |
---|
61 | // 05 Nov. 2000 MG Pia - Removed const cast previously introduced to get |
---|
62 | // the code compiled (const G4Material* now introduced in |
---|
63 | // electromagnetic/utils utils-V02-00-03 tag) |
---|
64 | // (this is going back and forth, to cope with Michel's |
---|
65 | // utils tag not being accepted yet by system testing) |
---|
66 | // 21 Nov. 2000 V.Ivanchenko Fix a problem in fluctuations |
---|
67 | // 23 Nov. 2000 V.Ivanchenko Ion type fluctuations only for charge>0 |
---|
68 | // 10 May 2001 V.Ivanchenko Clean up againist Linux compilation with -Wall |
---|
69 | // 23 May 2001 V.Ivanchenko Minor fix in PostStepDoIt |
---|
70 | // 07 June 2001 V.Ivanchenko Clean up AntiProtonDEDX + add print out |
---|
71 | // 18 June 2001 V.Ivanchenko Cleanup print out |
---|
72 | // 18 Oct. 2001 V.Ivanchenko Add fluorescence |
---|
73 | // 30 Oct. 2001 V.Ivanchenko Add minGammaEnergy and minElectronEnergy |
---|
74 | // 07 Dec 2001 V.Ivanchenko Add SetFluorescence method |
---|
75 | // 15 Feb 2002 V.Ivanchenko Fix problem of Generic Ions |
---|
76 | // 25 Mar 2002 V.Ivanchenko Fix problem of fluorescence below threshold |
---|
77 | // 28 Mar 2002 V.Ivanchenko Set fluorescence off by default |
---|
78 | // 09 Apr 2002 V.Ivanchenko Fix table problem of GenericIons |
---|
79 | // 28 May 2002 V.Ivanchenko Remove flag fStopAndKill |
---|
80 | // 31 May 2002 V.Ivanchenko Add path of Fluo + Auger cuts to |
---|
81 | // AtomicDeexcitation |
---|
82 | // 03 Jun 2002 MGP Restore fStopAndKill |
---|
83 | // 10 Jun 2002 V.Ivanchenko Restore fStopButAlive |
---|
84 | // 12 Jun 2002 V.Ivanchenko Fix in fluctuations - if tmax<2*Ipot Gaussian |
---|
85 | // fluctuations enables |
---|
86 | // 20 Sept 2002 V.Ivanchenko Clean up energy ranges for models |
---|
87 | // 07 Oct 2002 V.Ivanchenko Clean up initialisation of fluorescence |
---|
88 | // 28 Oct 2002 V.Ivanchenko Optimal binning for dE/dx |
---|
89 | // 10 Dec 2002 V.Ivanchenko antiProtonLowEnergy -> 25 keV, QEG model below |
---|
90 | // 21 Jan 2003 V.Ivanchenko Cut per region |
---|
91 | // 10 Mar 2003 V.Ivanchenko Use SubTypes for ions |
---|
92 | // 12 Apr 2003 V.Ivanchenko Cut per region for fluo AlongStep |
---|
93 | // 18 Apr 2003 V.Ivanchenko finalRange redefinition |
---|
94 | // 26 Apr 2003 V.Ivanchenko fix for stepLimit |
---|
95 | // 30 Mar 2004 S.Saliceti add shellCS data member and expFlag variable, |
---|
96 | // atom total cross section for the Empiric Model |
---|
97 | // 28 May 2004 V.Ivanchenko fix for ionisation of antiprotons in complex materials |
---|
98 | // 30 Aug 2004 V.Ivanchenko use energy limit for parameterisation from model |
---|
99 | // 03 Oct 2005 V.Ivanchenko change logic of definition of high energy limit for |
---|
100 | // parametrised proton model: min(user value, model limit) |
---|
101 | // 26 Jan 2005 S. Chauvie added PrintInfoDefinition() for antiproton |
---|
102 | // 30 Sep 2009 A.Mantero Removed dependencies to old shell Ionisation XS models |
---|
103 | // 07 Jun 2010 Code Celaning for June beta Release |
---|
104 | // ----------------------------------------------------------------------- |
---|
105 | |
---|
106 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
107 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
108 | |
---|
109 | #include "G4hLowEnergyIonisation.hh" |
---|
110 | #include "globals.hh" |
---|
111 | #include "G4ios.hh" |
---|
112 | #include "Randomize.hh" |
---|
113 | #include "G4Poisson.hh" |
---|
114 | #include "G4UnitsTable.hh" |
---|
115 | #include "G4EnergyLossTables.hh" |
---|
116 | #include "G4Material.hh" |
---|
117 | #include "G4DynamicParticle.hh" |
---|
118 | #include "G4ParticleDefinition.hh" |
---|
119 | #include "G4AtomicDeexcitation.hh" |
---|
120 | #include "G4AtomicTransitionManager.hh" |
---|
121 | #include "G4ShellVacancy.hh" |
---|
122 | #include "G4VhShellCrossSection.hh" |
---|
123 | #include "G4VEMDataSet.hh" |
---|
124 | #include "G4EMDataSet.hh" |
---|
125 | #include "G4CompositeEMDataSet.hh" |
---|
126 | #include "G4Gamma.hh" |
---|
127 | #include "G4LogLogInterpolation.hh" |
---|
128 | #include "G4SemiLogInterpolation.hh" |
---|
129 | #include "G4ProcessManager.hh" |
---|
130 | #include "G4ProductionCutsTable.hh" |
---|
131 | #include "G4teoCrossSection.hh" |
---|
132 | #include "G4empCrossSection.hh" |
---|
133 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
134 | |
---|
135 | G4hLowEnergyIonisation::G4hLowEnergyIonisation(const G4String& processName) |
---|
136 | : G4hLowEnergyLoss(processName), |
---|
137 | theBetheBlochModel(0), |
---|
138 | theProtonModel(0), |
---|
139 | theAntiProtonModel(0), |
---|
140 | theIonEffChargeModel(0), |
---|
141 | theNuclearStoppingModel(0), |
---|
142 | theIonChuFluctuationModel(0), |
---|
143 | theIonYangFluctuationModel(0), |
---|
144 | theProtonTable("ICRU_R49p"), |
---|
145 | theAntiProtonTable("ICRU_R49p"), |
---|
146 | theNuclearTable("ICRU_R49"), |
---|
147 | nStopping(true), |
---|
148 | theBarkas(true), |
---|
149 | theMeanFreePathTable(0), |
---|
150 | paramStepLimit (0.005), |
---|
151 | shellVacancy(0), |
---|
152 | shellCS(0), |
---|
153 | theFluo(false) |
---|
154 | { |
---|
155 | InitializeMe(); |
---|
156 | } |
---|
157 | |
---|
158 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
159 | |
---|
160 | void G4hLowEnergyIonisation::InitializeMe() |
---|
161 | { |
---|
162 | LowestKineticEnergy = 10.0*eV ; |
---|
163 | HighestKineticEnergy = 100.0*GeV ; |
---|
164 | MinKineticEnergy = 10.0*eV ; |
---|
165 | TotBin = 360 ; |
---|
166 | protonLowEnergy = 1.*keV ; |
---|
167 | protonHighEnergy = 100.*MeV ; |
---|
168 | antiProtonLowEnergy = 25.*keV ; |
---|
169 | antiProtonHighEnergy = 2.*MeV ; |
---|
170 | minGammaEnergy = 25.*keV; |
---|
171 | minElectronEnergy = 25.*keV; |
---|
172 | verboseLevel = 0; |
---|
173 | |
---|
174 | shellCS = new G4teoCrossSection("analytical"); |
---|
175 | |
---|
176 | deexcitationManager.InitialiseForNewRun(); |
---|
177 | deexcitationManager.SetAugerActive(false); |
---|
178 | deexcitationManager.SetPIXEActive(true); |
---|
179 | |
---|
180 | } |
---|
181 | |
---|
182 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
183 | |
---|
184 | G4hLowEnergyIonisation::~G4hLowEnergyIonisation() |
---|
185 | { |
---|
186 | if (theMeanFreePathTable) { |
---|
187 | theMeanFreePathTable->clearAndDestroy(); |
---|
188 | delete theMeanFreePathTable; |
---|
189 | } |
---|
190 | if(theBetheBlochModel)delete theBetheBlochModel; |
---|
191 | if(theProtonModel)delete theProtonModel; |
---|
192 | if(theAntiProtonModel)delete theAntiProtonModel; |
---|
193 | if(theNuclearStoppingModel)delete theNuclearStoppingModel; |
---|
194 | if(theIonEffChargeModel)delete theIonEffChargeModel; |
---|
195 | if(theIonChuFluctuationModel)delete theIonChuFluctuationModel; |
---|
196 | if(theIonYangFluctuationModel)delete theIonYangFluctuationModel; |
---|
197 | if(shellVacancy) delete shellVacancy; |
---|
198 | if(shellCS) delete shellCS; |
---|
199 | cutForDelta.clear(); |
---|
200 | G4int length = zFluoDataVector.size(); |
---|
201 | if(length) { |
---|
202 | for(G4int i=0; i<length; i++) { |
---|
203 | delete zFluoDataVector[i]; |
---|
204 | } |
---|
205 | zFluoDataVector.clear(); |
---|
206 | } |
---|
207 | } |
---|
208 | |
---|
209 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
210 | |
---|
211 | void G4hLowEnergyIonisation::SetElectronicStoppingPowerModel( |
---|
212 | const G4ParticleDefinition* aParticle, |
---|
213 | const G4String& dedxTable) |
---|
214 | // This method defines the ionisation parametrisation method via its name |
---|
215 | { |
---|
216 | if(0 < aParticle->GetPDGCharge()) { |
---|
217 | SetProtonElectronicStoppingPowerModel(dedxTable) ; |
---|
218 | } else { |
---|
219 | SetAntiProtonElectronicStoppingPowerModel(dedxTable) ; |
---|
220 | } |
---|
221 | } |
---|
222 | |
---|
223 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
224 | |
---|
225 | void G4hLowEnergyIonisation::InitializeParametrisation() |
---|
226 | |
---|
227 | { |
---|
228 | // Define models for parametrisation of electronic energy losses |
---|
229 | theBetheBlochModel = new G4hBetheBlochModel("Bethe-Bloch") ; |
---|
230 | theProtonModel = new G4hParametrisedLossModel(theProtonTable) ; |
---|
231 | protonHighEnergy = std::min(protonHighEnergy,theProtonModel->HighEnergyLimit(0, 0)); |
---|
232 | theAntiProtonModel = new G4QAOLowEnergyLoss(theAntiProtonTable) ; |
---|
233 | theNuclearStoppingModel = new G4hNuclearStoppingModel(theNuclearTable) ; |
---|
234 | theIonEffChargeModel = new G4hIonEffChargeSquare("Ziegler1988") ; |
---|
235 | theIonChuFluctuationModel = new G4IonChuFluctuationModel("Chu") ; |
---|
236 | theIonYangFluctuationModel = new G4IonYangFluctuationModel("Yang") ; |
---|
237 | } |
---|
238 | |
---|
239 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
240 | |
---|
241 | void G4hLowEnergyIonisation::BuildPhysicsTable( |
---|
242 | const G4ParticleDefinition& aParticleType) |
---|
243 | |
---|
244 | // just call BuildLossTable+BuildLambdaTable |
---|
245 | { |
---|
246 | if(verboseLevel > 0) { |
---|
247 | G4cout << "G4hLowEnergyIonisation::BuildPhysicsTable for " |
---|
248 | << aParticleType.GetParticleName() |
---|
249 | << " mass(MeV)= " << aParticleType.GetPDGMass()/MeV |
---|
250 | << " charge= " << aParticleType.GetPDGCharge()/eplus |
---|
251 | << " type= " << aParticleType.GetParticleType() |
---|
252 | << G4endl; |
---|
253 | |
---|
254 | if(verboseLevel > 1) { |
---|
255 | G4ProcessVector* pv = aParticleType.GetProcessManager()->GetProcessList(); |
---|
256 | |
---|
257 | G4cout << " 0: " << (*pv)[0]->GetProcessName() << " " << (*pv)[0] |
---|
258 | << " 1: " << (*pv)[1]->GetProcessName() << " " << (*pv)[1] |
---|
259 | // << " 2: " << (*pv)[2]->GetProcessName() << " " << (*pv)[2] |
---|
260 | << G4endl; |
---|
261 | G4cout << "ionModel= " << theIonEffChargeModel |
---|
262 | << " MFPtable= " << theMeanFreePathTable |
---|
263 | << " iniMass= " << initialMass |
---|
264 | << G4endl; |
---|
265 | } |
---|
266 | } |
---|
267 | |
---|
268 | if(aParticleType.GetParticleType() == "nucleus" && |
---|
269 | aParticleType.GetParticleName() != "GenericIon" && |
---|
270 | aParticleType.GetParticleSubType() == "generic") |
---|
271 | { |
---|
272 | |
---|
273 | G4EnergyLossTables::Register(&aParticleType, |
---|
274 | theDEDXpTable, |
---|
275 | theRangepTable, |
---|
276 | theInverseRangepTable, |
---|
277 | theLabTimepTable, |
---|
278 | theProperTimepTable, |
---|
279 | LowestKineticEnergy, HighestKineticEnergy, |
---|
280 | proton_mass_c2/aParticleType.GetPDGMass(), |
---|
281 | TotBin); |
---|
282 | |
---|
283 | return; |
---|
284 | } |
---|
285 | |
---|
286 | if( !CutsWhereModified() && theLossTable) return; |
---|
287 | |
---|
288 | InitializeParametrisation() ; |
---|
289 | G4Proton* theProton = G4Proton::Proton(); |
---|
290 | G4AntiProton* theAntiProton = G4AntiProton::AntiProton(); |
---|
291 | |
---|
292 | charge = aParticleType.GetPDGCharge()/eplus; |
---|
293 | chargeSquare = charge*charge ; |
---|
294 | |
---|
295 | const G4ProductionCutsTable* theCoupleTable= |
---|
296 | G4ProductionCutsTable::GetProductionCutsTable(); |
---|
297 | size_t numOfCouples = theCoupleTable->GetTableSize(); |
---|
298 | |
---|
299 | cutForDelta.clear(); |
---|
300 | cutForGamma.clear(); |
---|
301 | |
---|
302 | for (size_t j=0; j<numOfCouples; j++) { |
---|
303 | |
---|
304 | // get material parameters needed for the energy loss calculation |
---|
305 | const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(j); |
---|
306 | const G4Material* material= couple->GetMaterial(); |
---|
307 | |
---|
308 | // the cut cannot be below lowest limit |
---|
309 | G4double tCut = (*(theCoupleTable->GetEnergyCutsVector(1)))[j]; |
---|
310 | if(tCut > HighestKineticEnergy) tCut = HighestKineticEnergy; |
---|
311 | |
---|
312 | G4double excEnergy = material->GetIonisation()->GetMeanExcitationEnergy(); |
---|
313 | |
---|
314 | tCut = std::max(tCut,excEnergy); |
---|
315 | cutForDelta.push_back(tCut); |
---|
316 | |
---|
317 | // the cut cannot be below lowest limit |
---|
318 | tCut = (*(theCoupleTable->GetEnergyCutsVector(0)))[j]; |
---|
319 | if(tCut > HighestKineticEnergy) tCut = HighestKineticEnergy; |
---|
320 | tCut = std::max(tCut,minGammaEnergy); |
---|
321 | cutForGamma.push_back(tCut); |
---|
322 | } |
---|
323 | |
---|
324 | if(verboseLevel > 0) { |
---|
325 | G4cout << "Cuts are defined " << G4endl; |
---|
326 | } |
---|
327 | |
---|
328 | if(0.0 < charge) |
---|
329 | { |
---|
330 | { |
---|
331 | BuildLossTable(*theProton) ; |
---|
332 | |
---|
333 | // The following vector has a fixed dimension (see src/G4hLowEnergyLoss.cc for more details) |
---|
334 | // It happended in the past that caused memory corruption errors. The problem is still pending, even if temporary solved |
---|
335 | // G4cout << "[NOTE]: __LINE__=" << __LINE__ << ", aParticleType=" << aParticleType.GetParticleName() << ", theProton=" << theProton << ", theLossTable=" << theLossTable << ", CounterOfpProcess=" << CounterOfpProcess << G4endl; |
---|
336 | |
---|
337 | RecorderOfpProcess[CounterOfpProcess] = theLossTable ; |
---|
338 | CounterOfpProcess++; |
---|
339 | } |
---|
340 | } else { |
---|
341 | { |
---|
342 | BuildLossTable(*theAntiProton) ; |
---|
343 | |
---|
344 | // The following vector has a fixed dimension (see src/G4hLowEnergyLoss.cc for more details) |
---|
345 | // It happended in the past that caused memory corruption errors. The problem is still pending, even if temporary solved |
---|
346 | // G4cout << "[NOTE]: __LINE__=" << __LINE__ << ", aParticleType=" << aParticleType.GetParticleName() << ", theAntiProton=" << theAntiProton << ", theLossTable=" << theLossTable << ", CounterOfpbarProcess=" << CounterOfpbarProcess << G4endl; |
---|
347 | |
---|
348 | RecorderOfpbarProcess[CounterOfpbarProcess] = theLossTable ; |
---|
349 | CounterOfpbarProcess++; |
---|
350 | } |
---|
351 | } |
---|
352 | |
---|
353 | if(verboseLevel > 0) { |
---|
354 | G4cout << "G4hLowEnergyIonisation::BuildPhysicsTable: " |
---|
355 | << "Loss table is built " |
---|
356 | // << theLossTable |
---|
357 | << G4endl; |
---|
358 | } |
---|
359 | |
---|
360 | BuildLambdaTable(aParticleType) ; |
---|
361 | BuildDataForFluorescence(aParticleType); |
---|
362 | |
---|
363 | if(verboseLevel > 1) { |
---|
364 | G4cout << (*theMeanFreePathTable) << G4endl; |
---|
365 | } |
---|
366 | |
---|
367 | if(verboseLevel > 0) { |
---|
368 | G4cout << "G4hLowEnergyIonisation::BuildPhysicsTable: " |
---|
369 | << "DEDX table will be built " |
---|
370 | // << theDEDXpTable << " " << theDEDXpbarTable |
---|
371 | // << " " << theRangepTable << " " << theRangepbarTable |
---|
372 | << G4endl; |
---|
373 | } |
---|
374 | |
---|
375 | BuildDEDXTable(aParticleType) ; |
---|
376 | |
---|
377 | if(verboseLevel > 1) { |
---|
378 | G4cout << (*theDEDXpTable) << G4endl; |
---|
379 | } |
---|
380 | |
---|
381 | if((&aParticleType == theProton) || (&aParticleType == theAntiProton)) PrintInfoDefinition() ; |
---|
382 | |
---|
383 | if(verboseLevel > 0) { |
---|
384 | G4cout << "G4hLowEnergyIonisation::BuildPhysicsTable: end for " |
---|
385 | << aParticleType.GetParticleName() << G4endl; |
---|
386 | } |
---|
387 | } |
---|
388 | |
---|
389 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
390 | |
---|
391 | void G4hLowEnergyIonisation::BuildLossTable( |
---|
392 | const G4ParticleDefinition& aParticleType) |
---|
393 | { |
---|
394 | |
---|
395 | // Initialisation |
---|
396 | G4double lowEdgeEnergy , ionloss, ionlossBB, paramB ; |
---|
397 | G4double lowEnergy, highEnergy; |
---|
398 | G4Proton* theProton = G4Proton::Proton(); |
---|
399 | |
---|
400 | if(aParticleType == *theProton) { |
---|
401 | lowEnergy = protonLowEnergy ; |
---|
402 | highEnergy = protonHighEnergy ; |
---|
403 | charge = 1.0 ; |
---|
404 | } else { |
---|
405 | lowEnergy = antiProtonLowEnergy ; |
---|
406 | highEnergy = antiProtonHighEnergy ; |
---|
407 | charge = -1.0 ; |
---|
408 | } |
---|
409 | chargeSquare = 1.0 ; |
---|
410 | |
---|
411 | const G4ProductionCutsTable* theCoupleTable= |
---|
412 | G4ProductionCutsTable::GetProductionCutsTable(); |
---|
413 | size_t numOfCouples = theCoupleTable->GetTableSize(); |
---|
414 | |
---|
415 | if ( theLossTable) { |
---|
416 | theLossTable->clearAndDestroy(); |
---|
417 | delete theLossTable; |
---|
418 | } |
---|
419 | |
---|
420 | theLossTable = new G4PhysicsTable(numOfCouples); |
---|
421 | |
---|
422 | // loop for materials |
---|
423 | for (size_t j=0; j<numOfCouples; j++) { |
---|
424 | |
---|
425 | // create physics vector and fill it |
---|
426 | G4PhysicsLogVector* aVector = new G4PhysicsLogVector(LowestKineticEnergy, |
---|
427 | HighestKineticEnergy, |
---|
428 | TotBin); |
---|
429 | |
---|
430 | // get material parameters needed for the energy loss calculation |
---|
431 | const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(j); |
---|
432 | const G4Material* material= couple->GetMaterial(); |
---|
433 | |
---|
434 | if ( charge > 0.0 ) { |
---|
435 | ionloss = ProtonParametrisedDEDX(couple,highEnergy) ; |
---|
436 | } else { |
---|
437 | ionloss = AntiProtonParametrisedDEDX(couple,highEnergy) ; |
---|
438 | } |
---|
439 | |
---|
440 | ionlossBB = theBetheBlochModel->TheValue(&aParticleType,material,highEnergy) ; |
---|
441 | ionlossBB -= DeltaRaysEnergy(couple,highEnergy,proton_mass_c2) ; |
---|
442 | |
---|
443 | |
---|
444 | paramB = ionloss/ionlossBB - 1.0 ; |
---|
445 | |
---|
446 | // now comes the loop for the kinetic energy values |
---|
447 | for (G4int i = 0 ; i < TotBin ; i++) { |
---|
448 | lowEdgeEnergy = aVector->GetLowEdgeEnergy(i) ; |
---|
449 | |
---|
450 | // low energy part for this material, parametrised energy loss formulae |
---|
451 | if ( lowEdgeEnergy < highEnergy ) { |
---|
452 | |
---|
453 | if ( charge > 0.0 ) { |
---|
454 | ionloss = ProtonParametrisedDEDX(couple,lowEdgeEnergy) ; |
---|
455 | } else { |
---|
456 | ionloss = AntiProtonParametrisedDEDX(couple,lowEdgeEnergy) ; |
---|
457 | } |
---|
458 | |
---|
459 | } else { |
---|
460 | |
---|
461 | // high energy part for this material, Bethe-Bloch formula |
---|
462 | ionloss = theBetheBlochModel->TheValue(theProton,material, |
---|
463 | lowEdgeEnergy) ; |
---|
464 | |
---|
465 | ionloss -= DeltaRaysEnergy(couple,lowEdgeEnergy,proton_mass_c2) ; |
---|
466 | |
---|
467 | ionloss *= (1.0 + paramB*highEnergy/lowEdgeEnergy) ; |
---|
468 | } |
---|
469 | |
---|
470 | // now put the loss into the vector |
---|
471 | if(verboseLevel > 1) { |
---|
472 | G4cout << "E(MeV)= " << lowEdgeEnergy/MeV |
---|
473 | << " dE/dx(MeV/mm)= " << ionloss*mm/MeV |
---|
474 | << " in " << material->GetName() << G4endl; |
---|
475 | } |
---|
476 | aVector->PutValue(i,ionloss) ; |
---|
477 | } |
---|
478 | // Insert vector for this material into the table |
---|
479 | theLossTable->insert(aVector) ; |
---|
480 | } |
---|
481 | } |
---|
482 | |
---|
483 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
484 | |
---|
485 | void G4hLowEnergyIonisation::BuildDataForFluorescence( |
---|
486 | const G4ParticleDefinition& aParticleType) |
---|
487 | { |
---|
488 | |
---|
489 | if(verboseLevel > 1) { |
---|
490 | G4cout << "G4hLowEnergyIonisation::BuildDataForFluorescence for " |
---|
491 | << aParticleType.GetParticleName() << " is started" << G4endl; |
---|
492 | } |
---|
493 | |
---|
494 | // fill data for fluorescence |
---|
495 | |
---|
496 | deexcitationManager.SetCutForSecondaryPhotons(minGammaEnergy); |
---|
497 | deexcitationManager.SetCutForAugerElectrons(minElectronEnergy); |
---|
498 | |
---|
499 | G4double mass = aParticleType.GetPDGMass(); |
---|
500 | const G4ProductionCutsTable* theCoupleTable= |
---|
501 | G4ProductionCutsTable::GetProductionCutsTable(); |
---|
502 | size_t numOfCouples = theCoupleTable->GetTableSize(); |
---|
503 | |
---|
504 | if (shellVacancy != 0) delete shellVacancy; |
---|
505 | shellVacancy = new G4ShellVacancy(); |
---|
506 | G4DataVector* ksi = 0; |
---|
507 | G4DataVector* ksi1 = 0; |
---|
508 | G4DataVector* energy = 0; |
---|
509 | G4DataVector* energy1 = 0; |
---|
510 | size_t binForFluo = TotBin/10; |
---|
511 | G4int length = zFluoDataVector.size(); |
---|
512 | if(length > 0) { |
---|
513 | for(G4int i=0; i<length; i++) { |
---|
514 | G4VEMDataSet* x = zFluoDataVector[i]; |
---|
515 | delete x; |
---|
516 | } |
---|
517 | zFluoDataVector.clear(); |
---|
518 | } |
---|
519 | |
---|
520 | G4PhysicsLogVector* bVector = new G4PhysicsLogVector(LowestKineticEnergy, |
---|
521 | HighestKineticEnergy, |
---|
522 | binForFluo); |
---|
523 | const G4AtomicTransitionManager* transitionManager = |
---|
524 | G4AtomicTransitionManager::Instance(); |
---|
525 | |
---|
526 | G4double bindingEnergy; |
---|
527 | // G4double x; |
---|
528 | // G4double y; |
---|
529 | |
---|
530 | // loop for materials |
---|
531 | for (size_t j=0; j<numOfCouples; j++) { |
---|
532 | |
---|
533 | // get material parameters needed for the energy loss calculation |
---|
534 | const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(j); |
---|
535 | const G4Material* material= couple->GetMaterial(); |
---|
536 | |
---|
537 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
538 | size_t NumberOfElements = material->GetNumberOfElements() ; |
---|
539 | const G4double* theAtomicNumDensityVector = |
---|
540 | material->GetAtomicNumDensityVector(); |
---|
541 | G4VDataSetAlgorithm* interp = new G4SemiLogInterpolation(); |
---|
542 | G4VEMDataSet* xsis = new G4CompositeEMDataSet(interp, 1., 1.); |
---|
543 | G4VDataSetAlgorithm* interp1 = new G4SemiLogInterpolation(); |
---|
544 | G4VEMDataSet* xsis1 = new G4CompositeEMDataSet(interp1, 1., 1.); |
---|
545 | |
---|
546 | G4double tCut = cutForDelta[j]; |
---|
547 | G4double elDensity = 1.; |
---|
548 | |
---|
549 | for (size_t iel=0; iel<NumberOfElements; iel++ ) { |
---|
550 | |
---|
551 | G4int Z = (G4int)((*theElementVector)[iel]->GetZ()); |
---|
552 | G4int nShells = transitionManager->NumberOfShells(Z); |
---|
553 | energy = new G4DataVector(); |
---|
554 | ksi = new G4DataVector(); |
---|
555 | energy1= new G4DataVector(); |
---|
556 | ksi1 = new G4DataVector(); |
---|
557 | //if(NumberOfElements > 1) |
---|
558 | elDensity = theAtomicNumDensityVector[iel]/((G4double)nShells); |
---|
559 | |
---|
560 | for (size_t j = 0; j<binForFluo; j++) { |
---|
561 | |
---|
562 | G4double tkin = bVector->GetLowEdgeEnergy(j); |
---|
563 | G4double gamma = tkin/mass + 1.; |
---|
564 | G4double beta2 = 1.0 - 1.0/(gamma*gamma); |
---|
565 | G4double r = electron_mass_c2/mass; |
---|
566 | G4double tmax = 2.*electron_mass_c2*(gamma*gamma - 1.)/(1. + 2.*gamma*r + r*r); |
---|
567 | G4double cross = 0.; |
---|
568 | G4double cross1 = 0.; |
---|
569 | G4double eAverage= 0.; |
---|
570 | G4double tmin = std::min(tCut,tmax); |
---|
571 | G4double rel; |
---|
572 | |
---|
573 | for (G4int n=0; n<nShells; n++) { |
---|
574 | |
---|
575 | bindingEnergy = transitionManager->Shell(Z, n)->BindingEnergy(); |
---|
576 | if (tmin > bindingEnergy) { |
---|
577 | rel = std::log(tmin/bindingEnergy); |
---|
578 | eAverage += rel - beta2*(tmin - bindingEnergy)/tmax; |
---|
579 | cross += 1.0/bindingEnergy - 1.0/tmin - beta2*rel/tmax; |
---|
580 | } |
---|
581 | if (tmax > tmin) { |
---|
582 | cross1 += 1.0/tmin - 1.0/tmax - beta2*std::log(tmax/tmin)/tmax; |
---|
583 | } |
---|
584 | } |
---|
585 | |
---|
586 | cross1 *= elDensity; |
---|
587 | energy1->push_back(tkin); |
---|
588 | ksi1->push_back(cross1); |
---|
589 | |
---|
590 | if(eAverage > 0.) cross /= eAverage; |
---|
591 | else cross = 0.; |
---|
592 | |
---|
593 | energy->push_back(tkin); |
---|
594 | ksi->push_back(cross); |
---|
595 | } |
---|
596 | G4VDataSetAlgorithm* algo = interp->Clone(); |
---|
597 | G4VEMDataSet* set = new G4EMDataSet(Z,energy,ksi,algo,1.,1.); |
---|
598 | xsis->AddComponent(set); |
---|
599 | G4VDataSetAlgorithm* algo1 = interp1->Clone(); |
---|
600 | G4VEMDataSet* set1 = new G4EMDataSet(Z,energy1,ksi1,algo1,1.,1.); |
---|
601 | xsis1->AddComponent(set1); |
---|
602 | } |
---|
603 | if(verboseLevel > 1) { |
---|
604 | G4cout << "### Shell inverse cross sections for " |
---|
605 | << material->GetName() << G4endl; |
---|
606 | xsis->PrintData(); |
---|
607 | G4cout << "### Atom cross sections for " |
---|
608 | << material->GetName() << G4endl; |
---|
609 | xsis1->PrintData(); |
---|
610 | } |
---|
611 | shellVacancy->AddXsiTable(xsis); |
---|
612 | zFluoDataVector.push_back(xsis1); |
---|
613 | } |
---|
614 | delete bVector; |
---|
615 | } |
---|
616 | |
---|
617 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
618 | |
---|
619 | void G4hLowEnergyIonisation::BuildLambdaTable( |
---|
620 | const G4ParticleDefinition& aParticleType) |
---|
621 | |
---|
622 | { |
---|
623 | // Build mean free path tables for the delta ray production process |
---|
624 | // tables are built for MATERIALS |
---|
625 | |
---|
626 | if(verboseLevel > 1) { |
---|
627 | G4cout << "G4hLowEnergyIonisation::BuildLambdaTable for " |
---|
628 | << aParticleType.GetParticleName() << " is started" << G4endl; |
---|
629 | } |
---|
630 | |
---|
631 | |
---|
632 | G4double lowEdgeEnergy, value; |
---|
633 | charge = aParticleType.GetPDGCharge()/eplus ; |
---|
634 | chargeSquare = charge*charge ; |
---|
635 | initialMass = aParticleType.GetPDGMass(); |
---|
636 | |
---|
637 | const G4ProductionCutsTable* theCoupleTable= |
---|
638 | G4ProductionCutsTable::GetProductionCutsTable(); |
---|
639 | size_t numOfCouples = theCoupleTable->GetTableSize(); |
---|
640 | |
---|
641 | |
---|
642 | if (theMeanFreePathTable) { |
---|
643 | theMeanFreePathTable->clearAndDestroy(); |
---|
644 | delete theMeanFreePathTable; |
---|
645 | } |
---|
646 | |
---|
647 | theMeanFreePathTable = new G4PhysicsTable(numOfCouples); |
---|
648 | |
---|
649 | // loop for materials |
---|
650 | |
---|
651 | for (size_t J=0 ; J < numOfCouples; J++) { |
---|
652 | |
---|
653 | //create physics vector then fill it .... |
---|
654 | G4PhysicsLogVector* aVector = new G4PhysicsLogVector(LowestKineticEnergy, |
---|
655 | HighestKineticEnergy, |
---|
656 | TotBin); |
---|
657 | |
---|
658 | // compute the (macroscopic) cross section first |
---|
659 | const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(J); |
---|
660 | const G4Material* material= couple->GetMaterial(); |
---|
661 | |
---|
662 | const G4ElementVector* theElementVector = |
---|
663 | material->GetElementVector() ; |
---|
664 | const G4double* theAtomicNumDensityVector = |
---|
665 | material->GetAtomicNumDensityVector(); |
---|
666 | const G4int NumberOfElements = material->GetNumberOfElements() ; |
---|
667 | |
---|
668 | // get the electron kinetic energy cut for the actual material, |
---|
669 | // it will be used in ComputeMicroscopicCrossSection |
---|
670 | // ( it is the SAME for ALL the ELEMENTS in THIS MATERIAL ) |
---|
671 | // ------------------------------------------------------ |
---|
672 | |
---|
673 | G4double deltaCut = cutForDelta[J]; |
---|
674 | |
---|
675 | for ( G4int i = 0 ; i < TotBin ; i++ ) { |
---|
676 | lowEdgeEnergy = aVector->GetLowEdgeEnergy(i) ; |
---|
677 | G4double sigma = 0.0 ; |
---|
678 | G4int Z; |
---|
679 | |
---|
680 | for (G4int iel=0; iel<NumberOfElements; iel++ ) { |
---|
681 | Z = (G4int) (*theElementVector)[iel]->GetZ(); |
---|
682 | totalCrossSectionMap [Z] = ComputeMicroscopicCrossSection( |
---|
683 | aParticleType, |
---|
684 | lowEdgeEnergy, |
---|
685 | Z, |
---|
686 | deltaCut ) ; |
---|
687 | sigma += theAtomicNumDensityVector[iel]*ComputeMicroscopicCrossSection( |
---|
688 | aParticleType, |
---|
689 | lowEdgeEnergy, |
---|
690 | Z, |
---|
691 | deltaCut ) ; |
---|
692 | |
---|
693 | } |
---|
694 | |
---|
695 | // mean free path = 1./macroscopic cross section |
---|
696 | |
---|
697 | value = sigma<=0 ? DBL_MAX : 1./sigma ; |
---|
698 | |
---|
699 | aVector->PutValue(i, value) ; |
---|
700 | } |
---|
701 | |
---|
702 | theMeanFreePathTable->insert(aVector); |
---|
703 | } |
---|
704 | |
---|
705 | } |
---|
706 | |
---|
707 | |
---|
708 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
709 | |
---|
710 | G4double G4hLowEnergyIonisation::ComputeMicroscopicCrossSection( |
---|
711 | const G4ParticleDefinition& aParticleType, |
---|
712 | G4double kineticEnergy, |
---|
713 | G4double atomicNumber, |
---|
714 | G4double deltaCutInEnergy) const |
---|
715 | { |
---|
716 | //****************************************************************** |
---|
717 | // cross section formula is OK for spin=0, 1/2, 1 only ! |
---|
718 | // ***************************************************************** |
---|
719 | |
---|
720 | // calculates the microscopic cross section in GEANT4 internal units |
---|
721 | // ( it is called for elements , AtomicNumber = z ) |
---|
722 | |
---|
723 | G4double energy, gamma, beta2, tmax, var; |
---|
724 | G4double totalCrossSection = 0.0 ; |
---|
725 | |
---|
726 | G4double particleMass = initialMass; |
---|
727 | |
---|
728 | // get particle data ................................... |
---|
729 | |
---|
730 | energy = kineticEnergy + particleMass; |
---|
731 | |
---|
732 | // some kinematics...................... |
---|
733 | |
---|
734 | gamma = energy/particleMass; |
---|
735 | beta2 = 1.0 - 1.0/(gamma*gamma); |
---|
736 | var = electron_mass_c2/particleMass; |
---|
737 | tmax = 2.*electron_mass_c2*(gamma*gamma - 1.)/(1. + 2.*gamma*var + var*var); |
---|
738 | |
---|
739 | // now you can calculate the total cross section |
---|
740 | |
---|
741 | if( tmax > deltaCutInEnergy ) { |
---|
742 | |
---|
743 | var=deltaCutInEnergy/tmax; |
---|
744 | totalCrossSection = (1.0 - var*(1.0 - beta2*std::log(var))) / deltaCutInEnergy ; |
---|
745 | G4double spin = aParticleType.GetPDGSpin() ; |
---|
746 | |
---|
747 | // +term for spin=1/2 particle |
---|
748 | if( 0.5 == spin ) |
---|
749 | totalCrossSection += 0.5 * (tmax - deltaCutInEnergy) / (energy*energy); |
---|
750 | |
---|
751 | // +term for spin=1 particle |
---|
752 | else if( 0.9 < spin ) |
---|
753 | totalCrossSection += -std::log(var)/(3.0*deltaCutInEnergy) + |
---|
754 | (tmax - deltaCutInEnergy) * ( (5.0+ 1.0/var)*0.25 / (energy*energy) - |
---|
755 | beta2 / (tmax * deltaCutInEnergy) ) / 3.0 ; |
---|
756 | |
---|
757 | totalCrossSection *= twopi_mc2_rcl2 * atomicNumber / beta2 ; |
---|
758 | } |
---|
759 | |
---|
760 | return totalCrossSection ; |
---|
761 | } |
---|
762 | |
---|
763 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
764 | |
---|
765 | G4double G4hLowEnergyIonisation::GetMeanFreePath(const G4Track& trackData, |
---|
766 | G4double, // previousStepSize |
---|
767 | enum G4ForceCondition* condition) |
---|
768 | { |
---|
769 | const G4DynamicParticle* aParticle = trackData.GetDynamicParticle(); |
---|
770 | const G4MaterialCutsCouple* couple = trackData.GetMaterialCutsCouple(); |
---|
771 | const G4Material* material = couple->GetMaterial(); |
---|
772 | G4double meanFreePath; |
---|
773 | G4bool isOutRange ; |
---|
774 | |
---|
775 | *condition = NotForced ; |
---|
776 | |
---|
777 | G4double kineticEnergy = (aParticle->GetKineticEnergy())*initialMass/(aParticle->GetMass()); |
---|
778 | charge = aParticle->GetCharge()/eplus; |
---|
779 | chargeSquare = theIonEffChargeModel->TheValue(aParticle, material); |
---|
780 | |
---|
781 | if(kineticEnergy < LowestKineticEnergy) meanFreePath = DBL_MAX; |
---|
782 | |
---|
783 | else { |
---|
784 | if(kineticEnergy > HighestKineticEnergy) |
---|
785 | kineticEnergy = HighestKineticEnergy; |
---|
786 | meanFreePath = (((*theMeanFreePathTable)(couple->GetIndex()))-> |
---|
787 | GetValue(kineticEnergy,isOutRange))/chargeSquare; |
---|
788 | } |
---|
789 | |
---|
790 | return meanFreePath ; |
---|
791 | } |
---|
792 | |
---|
793 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
794 | |
---|
795 | G4double G4hLowEnergyIonisation::GetConstraints( |
---|
796 | const G4DynamicParticle* particle, |
---|
797 | const G4MaterialCutsCouple* couple) |
---|
798 | { |
---|
799 | // returns the Step limit |
---|
800 | // dEdx is calculated as well as the range |
---|
801 | // based on Effective Charge Approach |
---|
802 | |
---|
803 | const G4Material* material = couple->GetMaterial(); |
---|
804 | G4Proton* theProton = G4Proton::Proton(); |
---|
805 | G4AntiProton* theAntiProton = G4AntiProton::AntiProton(); |
---|
806 | |
---|
807 | G4double stepLimit = 0.0 ; |
---|
808 | G4double dx, highEnergy; |
---|
809 | |
---|
810 | G4double massRatio = proton_mass_c2/(particle->GetMass()) ; |
---|
811 | G4double kineticEnergy = particle->GetKineticEnergy() ; |
---|
812 | |
---|
813 | // Scale the kinetic energy |
---|
814 | |
---|
815 | G4double tscaled = kineticEnergy*massRatio ; |
---|
816 | fBarkas = 0.0; |
---|
817 | |
---|
818 | if(charge > 0.0) { |
---|
819 | |
---|
820 | highEnergy = protonHighEnergy ; |
---|
821 | |
---|
822 | fRangeNow = G4EnergyLossTables::GetRange(theProton, tscaled, couple); |
---|
823 | dx = G4EnergyLossTables::GetRange(theProton, highEnergy, couple); |
---|
824 | fdEdx = G4EnergyLossTables::GetDEDX(theProton, tscaled, couple) |
---|
825 | * chargeSquare ; |
---|
826 | |
---|
827 | // Correction for positive ions |
---|
828 | if(theBarkas && tscaled > highEnergy) { |
---|
829 | fBarkas = BarkasTerm(material,tscaled)*std::sqrt(chargeSquare)*chargeSquare |
---|
830 | + BlochTerm(material,tscaled,chargeSquare); |
---|
831 | } |
---|
832 | // Antiprotons and negative hadrons |
---|
833 | } else { |
---|
834 | |
---|
835 | highEnergy = antiProtonHighEnergy ; |
---|
836 | fRangeNow = G4EnergyLossTables::GetRange(theAntiProton, tscaled, couple); |
---|
837 | dx = G4EnergyLossTables::GetRange(theAntiProton, highEnergy, couple); |
---|
838 | fdEdx = G4EnergyLossTables::GetDEDX(theAntiProton, tscaled, couple) |
---|
839 | * chargeSquare ; |
---|
840 | |
---|
841 | if(theBarkas && tscaled > highEnergy) { |
---|
842 | fBarkas = -BarkasTerm(material,tscaled)*std::sqrt(chargeSquare)*chargeSquare |
---|
843 | + BlochTerm(material,tscaled,chargeSquare); |
---|
844 | } |
---|
845 | } |
---|
846 | /* |
---|
847 | const G4Material* mat = couple->GetMaterial(); |
---|
848 | G4double fac = gram/(MeV*cm2*mat->GetDensity()); |
---|
849 | G4cout << particle->GetDefinition()->GetParticleName() |
---|
850 | << " in " << mat->GetName() |
---|
851 | << " E(MeV)= " << kineticEnergy/MeV |
---|
852 | << " dedx(MeV*cm^2/g)= " << fdEdx*fac |
---|
853 | << " barcas(MeV*cm^2/gram)= " << fBarkas*fac |
---|
854 | << " Q^2= " << chargeSquare |
---|
855 | << G4endl; |
---|
856 | */ |
---|
857 | // scaling back |
---|
858 | fRangeNow /= (chargeSquare*massRatio) ; |
---|
859 | dx /= (chargeSquare*massRatio) ; |
---|
860 | |
---|
861 | stepLimit = fRangeNow ; |
---|
862 | G4double r = std::min(finalRange, couple->GetProductionCuts() |
---|
863 | ->GetProductionCut(idxG4ElectronCut)); |
---|
864 | |
---|
865 | if (fRangeNow > r) { |
---|
866 | stepLimit = dRoverRange*fRangeNow + r*(1.0 - dRoverRange)*(2.0 - r/fRangeNow); |
---|
867 | if (stepLimit > fRangeNow) stepLimit = fRangeNow; |
---|
868 | } |
---|
869 | // compute the (random) Step limit in standard energy range |
---|
870 | if(tscaled > highEnergy ) { |
---|
871 | |
---|
872 | // add Barkas correction directly to dedx |
---|
873 | fdEdx += fBarkas; |
---|
874 | |
---|
875 | if(stepLimit > fRangeNow - dx*0.9) stepLimit = fRangeNow - dx*0.9 ; |
---|
876 | |
---|
877 | // Step limit in low energy range |
---|
878 | } else { |
---|
879 | G4double x = dx*paramStepLimit; |
---|
880 | if (stepLimit > x) stepLimit = x; |
---|
881 | } |
---|
882 | return stepLimit ; |
---|
883 | } |
---|
884 | |
---|
885 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
886 | |
---|
887 | G4VParticleChange* G4hLowEnergyIonisation::AlongStepDoIt( |
---|
888 | const G4Track& trackData, |
---|
889 | const G4Step& stepData) |
---|
890 | { |
---|
891 | // compute the energy loss after a step |
---|
892 | G4Proton* theProton = G4Proton::Proton(); |
---|
893 | G4AntiProton* theAntiProton = G4AntiProton::AntiProton(); |
---|
894 | G4double finalT = 0.0 ; |
---|
895 | |
---|
896 | aParticleChange.Initialize(trackData) ; |
---|
897 | |
---|
898 | const G4MaterialCutsCouple* couple = trackData.GetMaterialCutsCouple(); |
---|
899 | const G4Material* material = couple->GetMaterial(); |
---|
900 | |
---|
901 | // get the actual (true) Step length from stepData |
---|
902 | const G4double step = stepData.GetStepLength() ; |
---|
903 | |
---|
904 | const G4DynamicParticle* particle = trackData.GetDynamicParticle() ; |
---|
905 | |
---|
906 | G4double kineticEnergy = particle->GetKineticEnergy() ; |
---|
907 | G4double massRatio = proton_mass_c2/(particle->GetMass()) ; |
---|
908 | G4double tscaled= kineticEnergy*massRatio ; |
---|
909 | G4double eloss = 0.0 ; |
---|
910 | G4double nloss = 0.0 ; |
---|
911 | |
---|
912 | |
---|
913 | // very small particle energy |
---|
914 | if(kineticEnergy < MinKineticEnergy) { |
---|
915 | |
---|
916 | eloss = kineticEnergy ; |
---|
917 | |
---|
918 | // particle energy outside tabulated energy range |
---|
919 | } else if( kineticEnergy > HighestKineticEnergy) { |
---|
920 | eloss = step*fdEdx ; |
---|
921 | |
---|
922 | // big step |
---|
923 | } else if(step >= fRangeNow ) { |
---|
924 | eloss = kineticEnergy ; |
---|
925 | |
---|
926 | // tabulated range |
---|
927 | } else { |
---|
928 | |
---|
929 | // step longer than linear step limit |
---|
930 | if(step > linLossLimit*fRangeNow) { |
---|
931 | |
---|
932 | G4double rscaled= fRangeNow*massRatio*chargeSquare ; |
---|
933 | G4double sscaled= step *massRatio*chargeSquare ; |
---|
934 | |
---|
935 | if(charge > 0.0) { |
---|
936 | eloss = G4EnergyLossTables::GetPreciseEnergyFromRange( |
---|
937 | theProton,rscaled, couple) - |
---|
938 | G4EnergyLossTables::GetPreciseEnergyFromRange( |
---|
939 | theProton,rscaled-sscaled,couple) ; |
---|
940 | |
---|
941 | } else { |
---|
942 | eloss = G4EnergyLossTables::GetPreciseEnergyFromRange( |
---|
943 | theAntiProton,rscaled,couple) - |
---|
944 | G4EnergyLossTables::GetPreciseEnergyFromRange( |
---|
945 | theAntiProton,rscaled-sscaled,couple) ; |
---|
946 | } |
---|
947 | eloss /= massRatio ; |
---|
948 | |
---|
949 | // Barkas correction at big step |
---|
950 | eloss += fBarkas*step; |
---|
951 | |
---|
952 | // step shorter than linear step limit |
---|
953 | } else { |
---|
954 | eloss = step*fdEdx ; |
---|
955 | } |
---|
956 | if(nStopping && tscaled < protonHighEnergy) { |
---|
957 | nloss = (theNuclearStoppingModel->TheValue(particle, material))*step; |
---|
958 | } |
---|
959 | } |
---|
960 | |
---|
961 | if(eloss < 0.0) eloss = 0.0; |
---|
962 | |
---|
963 | finalT = kineticEnergy - eloss - nloss; |
---|
964 | |
---|
965 | if( EnlossFlucFlag && 0.0 < eloss && finalT > MinKineticEnergy) { |
---|
966 | |
---|
967 | // now the electron loss with fluctuation |
---|
968 | eloss = ElectronicLossFluctuation(particle, couple, eloss, step) ; |
---|
969 | if(eloss < 0.0) eloss = 0.0; |
---|
970 | finalT = kineticEnergy - eloss - nloss; |
---|
971 | } |
---|
972 | |
---|
973 | // stop particle if the kinetic energy <= MinKineticEnergy |
---|
974 | if (finalT*massRatio <= MinKineticEnergy ) { |
---|
975 | |
---|
976 | finalT = 0.0; |
---|
977 | if(!particle->GetDefinition()->GetProcessManager()-> |
---|
978 | GetAtRestProcessVector()->size()) |
---|
979 | aParticleChange.ProposeTrackStatus(fStopAndKill); |
---|
980 | else |
---|
981 | aParticleChange.ProposeTrackStatus(fStopButAlive); |
---|
982 | } |
---|
983 | |
---|
984 | aParticleChange.ProposeEnergy( finalT ); |
---|
985 | eloss = kineticEnergy-finalT; |
---|
986 | |
---|
987 | // Deexcitation only of ionised atoms |
---|
988 | G4double hMass = particle->GetMass(); |
---|
989 | std::vector<G4DynamicParticle*>* newpart = 0; |
---|
990 | G4DynamicParticle* part = 0; |
---|
991 | |
---|
992 | if(theFluo) newpart = DeexciteAtom(couple, kineticEnergy, hMass, eloss); |
---|
993 | |
---|
994 | if(newpart != 0) { |
---|
995 | |
---|
996 | // G4cout << "AlongStep DEEXCTATION!!!" << G4endl; //debug |
---|
997 | size_t nSecondaries = newpart->size(); |
---|
998 | aParticleChange.SetNumberOfSecondaries(nSecondaries); |
---|
999 | G4Track* newtrack = 0; |
---|
1000 | const G4StepPoint* preStep = stepData.GetPreStepPoint(); |
---|
1001 | const G4StepPoint* postStep = stepData.GetPostStepPoint(); |
---|
1002 | G4ThreeVector r = preStep->GetPosition(); |
---|
1003 | G4ThreeVector deltaR = postStep->GetPosition(); |
---|
1004 | deltaR -= r; |
---|
1005 | G4double t = preStep->GetGlobalTime(); |
---|
1006 | G4double deltaT = postStep->GetGlobalTime(); |
---|
1007 | deltaT -= t; |
---|
1008 | G4double time, q, e; |
---|
1009 | G4ThreeVector position; |
---|
1010 | |
---|
1011 | for(size_t i=0; i<nSecondaries; i++) { |
---|
1012 | |
---|
1013 | part = (*newpart)[i]; |
---|
1014 | if(part) { |
---|
1015 | |
---|
1016 | e = part->GetKineticEnergy(); |
---|
1017 | if(e <= eloss) { |
---|
1018 | |
---|
1019 | eloss -= e; |
---|
1020 | q = G4UniformRand(); |
---|
1021 | time = deltaT*q + t; |
---|
1022 | position = deltaR*q; |
---|
1023 | position += r; |
---|
1024 | newtrack = new G4Track(part, time, position); |
---|
1025 | aParticleChange.AddSecondary(newtrack); |
---|
1026 | |
---|
1027 | } else { |
---|
1028 | |
---|
1029 | delete part; |
---|
1030 | |
---|
1031 | } |
---|
1032 | } |
---|
1033 | } |
---|
1034 | delete newpart; |
---|
1035 | } |
---|
1036 | |
---|
1037 | aParticleChange.ProposeLocalEnergyDeposit(eloss); |
---|
1038 | return &aParticleChange ; |
---|
1039 | } |
---|
1040 | |
---|
1041 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1042 | |
---|
1043 | G4double G4hLowEnergyIonisation::ProtonParametrisedDEDX( |
---|
1044 | const G4MaterialCutsCouple* couple, |
---|
1045 | G4double kineticEnergy) const |
---|
1046 | { |
---|
1047 | const G4Material* material = couple->GetMaterial(); |
---|
1048 | G4Proton* theProton = G4Proton::Proton(); |
---|
1049 | G4double eloss = 0.0; |
---|
1050 | |
---|
1051 | // Free Electron Gas Model |
---|
1052 | if(kineticEnergy < protonLowEnergy) { |
---|
1053 | eloss = (theProtonModel->TheValue(theProton, material, protonLowEnergy)) |
---|
1054 | * std::sqrt(kineticEnergy/protonLowEnergy) ; |
---|
1055 | |
---|
1056 | // Parametrisation |
---|
1057 | } else { |
---|
1058 | eloss = theProtonModel->TheValue(theProton, material, kineticEnergy) ; |
---|
1059 | } |
---|
1060 | |
---|
1061 | // Delta rays energy |
---|
1062 | eloss -= DeltaRaysEnergy(couple,kineticEnergy,proton_mass_c2) ; |
---|
1063 | |
---|
1064 | if(verboseLevel > 2) { |
---|
1065 | G4cout << "p E(MeV)= " << kineticEnergy/MeV |
---|
1066 | << " dE/dx(MeV/mm)= " << eloss*mm/MeV |
---|
1067 | << " for " << material->GetName() |
---|
1068 | << " model: " << theProtonModel << G4endl; |
---|
1069 | } |
---|
1070 | |
---|
1071 | if(eloss < 0.0) eloss = 0.0 ; |
---|
1072 | |
---|
1073 | return eloss ; |
---|
1074 | } |
---|
1075 | |
---|
1076 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1077 | |
---|
1078 | G4double G4hLowEnergyIonisation::AntiProtonParametrisedDEDX( |
---|
1079 | const G4MaterialCutsCouple* couple, |
---|
1080 | G4double kineticEnergy) const |
---|
1081 | { |
---|
1082 | const G4Material* material = couple->GetMaterial(); |
---|
1083 | G4AntiProton* theAntiProton = G4AntiProton::AntiProton(); |
---|
1084 | G4double eloss = 0.0 ; |
---|
1085 | |
---|
1086 | // Antiproton model is used |
---|
1087 | if(theAntiProtonModel->IsInCharge(theAntiProton,material)) { |
---|
1088 | if(kineticEnergy < antiProtonLowEnergy) { |
---|
1089 | eloss = theAntiProtonModel->TheValue(theAntiProton,material,antiProtonLowEnergy) |
---|
1090 | * std::sqrt(kineticEnergy/antiProtonLowEnergy) ; |
---|
1091 | |
---|
1092 | // Parametrisation |
---|
1093 | } else { |
---|
1094 | eloss = theAntiProtonModel->TheValue(theAntiProton,material, |
---|
1095 | kineticEnergy); |
---|
1096 | } |
---|
1097 | |
---|
1098 | // The proton model is used + Barkas correction |
---|
1099 | } else { |
---|
1100 | if(kineticEnergy < protonLowEnergy) { |
---|
1101 | eloss = theProtonModel->TheValue(G4Proton::Proton(),material,protonLowEnergy) |
---|
1102 | * std::sqrt(kineticEnergy/protonLowEnergy) ; |
---|
1103 | |
---|
1104 | // Parametrisation |
---|
1105 | } else { |
---|
1106 | eloss = theProtonModel->TheValue(G4Proton::Proton(),material, |
---|
1107 | kineticEnergy); |
---|
1108 | } |
---|
1109 | //if(theBarkas) eloss -= 2.0*BarkasTerm(material, kineticEnergy); |
---|
1110 | } |
---|
1111 | |
---|
1112 | // Delta rays energy |
---|
1113 | eloss -= DeltaRaysEnergy(couple,kineticEnergy,proton_mass_c2) ; |
---|
1114 | |
---|
1115 | if(verboseLevel > 2) { |
---|
1116 | G4cout << "pbar E(MeV)= " << kineticEnergy/MeV |
---|
1117 | << " dE/dx(MeV/mm)= " << eloss*mm/MeV |
---|
1118 | << " for " << material->GetName() |
---|
1119 | << " model: " << theProtonModel << G4endl; |
---|
1120 | } |
---|
1121 | |
---|
1122 | if(eloss < 0.0) eloss = 0.0 ; |
---|
1123 | |
---|
1124 | return eloss ; |
---|
1125 | } |
---|
1126 | |
---|
1127 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1128 | |
---|
1129 | G4double G4hLowEnergyIonisation::DeltaRaysEnergy( |
---|
1130 | const G4MaterialCutsCouple* couple, |
---|
1131 | G4double kineticEnergy, |
---|
1132 | G4double particleMass) const |
---|
1133 | { |
---|
1134 | G4double dloss = 0.0 ; |
---|
1135 | |
---|
1136 | G4double deltaCutNow = cutForDelta[(couple->GetIndex())] ; |
---|
1137 | const G4Material* material = couple->GetMaterial(); |
---|
1138 | G4double electronDensity = material->GetElectronDensity(); |
---|
1139 | G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy(); |
---|
1140 | |
---|
1141 | G4double tau = kineticEnergy/particleMass ; |
---|
1142 | G4double rateMass = electron_mass_c2/particleMass ; |
---|
1143 | |
---|
1144 | // some local variables |
---|
1145 | |
---|
1146 | G4double gamma,bg2,beta2,tmax,x ; |
---|
1147 | |
---|
1148 | gamma = tau + 1.0 ; |
---|
1149 | bg2 = tau*(tau+2.0) ; |
---|
1150 | beta2 = bg2/(gamma*gamma) ; |
---|
1151 | tmax = 2.*electron_mass_c2*bg2/(1.0+2.0*gamma*rateMass+rateMass*rateMass) ; |
---|
1152 | |
---|
1153 | // Validity range for delta electron cross section |
---|
1154 | G4double deltaCut = std::max(deltaCutNow, eexc); |
---|
1155 | |
---|
1156 | if ( deltaCut < tmax) { |
---|
1157 | x = deltaCut / tmax ; |
---|
1158 | dloss = ( beta2 * (x - 1.0) - std::log(x) ) * twopi_mc2_rcl2 |
---|
1159 | * electronDensity / beta2 ; |
---|
1160 | } |
---|
1161 | return dloss ; |
---|
1162 | } |
---|
1163 | |
---|
1164 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1165 | |
---|
1166 | G4VParticleChange* G4hLowEnergyIonisation::PostStepDoIt( |
---|
1167 | const G4Track& trackData, |
---|
1168 | const G4Step& stepData) |
---|
1169 | { |
---|
1170 | // Units are expressed in GEANT4 internal units. |
---|
1171 | |
---|
1172 | G4double KineticEnergy,TotalEnergy,TotalMomentum,betasquare, |
---|
1173 | DeltaKineticEnergy,DeltaTotalMomentum,costheta,sintheta,phi, |
---|
1174 | dirx,diry,dirz,finalKineticEnergy,finalPx,finalPy,finalPz, |
---|
1175 | x,xc,grej,Psquare,Esquare,rate,finalMomentum ; |
---|
1176 | |
---|
1177 | aParticleChange.Initialize(trackData) ; |
---|
1178 | const G4MaterialCutsCouple* couple = trackData.GetMaterialCutsCouple(); |
---|
1179 | |
---|
1180 | const G4DynamicParticle* aParticle = trackData.GetDynamicParticle() ; |
---|
1181 | |
---|
1182 | // some kinematics |
---|
1183 | |
---|
1184 | ParticleMass=aParticle->GetDefinition()->GetPDGMass(); |
---|
1185 | KineticEnergy=aParticle->GetKineticEnergy(); |
---|
1186 | TotalEnergy=KineticEnergy + ParticleMass ; |
---|
1187 | Psquare=KineticEnergy*(TotalEnergy+ParticleMass) ; |
---|
1188 | Esquare=TotalEnergy*TotalEnergy; |
---|
1189 | betasquare=Psquare/Esquare; |
---|
1190 | G4ThreeVector ParticleDirection = aParticle->GetMomentumDirection() ; |
---|
1191 | |
---|
1192 | G4double gamma= KineticEnergy/ParticleMass + 1.; |
---|
1193 | G4double r = electron_mass_c2/ParticleMass; |
---|
1194 | G4double tmax = 2.*electron_mass_c2*(gamma*gamma - 1.)/(1. + 2.*gamma*r + r*r); |
---|
1195 | |
---|
1196 | // Validity range for delta electron cross section |
---|
1197 | G4double DeltaCut = cutForDelta[couple->GetIndex()]; |
---|
1198 | |
---|
1199 | // This should not be a case |
---|
1200 | if(DeltaCut >= tmax) |
---|
1201 | return G4VContinuousDiscreteProcess::PostStepDoIt(trackData,stepData); |
---|
1202 | |
---|
1203 | xc = DeltaCut / tmax; |
---|
1204 | rate = tmax / TotalEnergy; |
---|
1205 | rate = rate*rate ; |
---|
1206 | G4double spin = aParticle->GetDefinition()->GetPDGSpin() ; |
---|
1207 | |
---|
1208 | // sampling follows ... |
---|
1209 | do { |
---|
1210 | x=xc/(1.-(1.-xc)*G4UniformRand()); |
---|
1211 | |
---|
1212 | if(0.0 == spin) { |
---|
1213 | grej = 1.0 - betasquare * x ; |
---|
1214 | |
---|
1215 | } else if (0.5 == spin) { |
---|
1216 | grej = (1.0 - betasquare * x + 0.5*x*x*rate) / (1.0 + 0.5 * rate) ; |
---|
1217 | |
---|
1218 | } else { |
---|
1219 | grej = (1.0 - betasquare * x ) * (1.0 + x/ (3.0*xc)) + |
---|
1220 | x * x * rate * (1.0 + 0.5 * x / xc) / 3.0 / |
---|
1221 | (1.0 + 1.0/(3.0*xc) + rate *(1.0+ 0.5/xc) /3.0) ; |
---|
1222 | } |
---|
1223 | |
---|
1224 | } while( G4UniformRand() > grej ); |
---|
1225 | |
---|
1226 | |
---|
1227 | DeltaKineticEnergy = x * tmax; |
---|
1228 | |
---|
1229 | DeltaTotalMomentum = std::sqrt(DeltaKineticEnergy * (DeltaKineticEnergy + |
---|
1230 | 2. * electron_mass_c2 )) ; |
---|
1231 | TotalMomentum = std::sqrt(Psquare) ; |
---|
1232 | costheta = DeltaKineticEnergy * (TotalEnergy + electron_mass_c2) |
---|
1233 | /(DeltaTotalMomentum * TotalMomentum) ; |
---|
1234 | |
---|
1235 | // protection against costheta > 1 or < -1 --------------- |
---|
1236 | if ( costheta < -1. ) |
---|
1237 | costheta = -1. ; |
---|
1238 | if ( costheta > +1. ) |
---|
1239 | costheta = +1. ; |
---|
1240 | |
---|
1241 | // direction of the delta electron ........ |
---|
1242 | phi = twopi * G4UniformRand() ; |
---|
1243 | sintheta = std::sqrt(1. - costheta*costheta); |
---|
1244 | dirx = sintheta * std::cos(phi) ; |
---|
1245 | diry = sintheta * std::sin(phi) ; |
---|
1246 | dirz = costheta ; |
---|
1247 | |
---|
1248 | G4ThreeVector DeltaDirection(dirx,diry,dirz) ; |
---|
1249 | DeltaDirection.rotateUz(ParticleDirection) ; |
---|
1250 | |
---|
1251 | // create G4DynamicParticle object for delta ray |
---|
1252 | G4DynamicParticle *theDeltaRay = new G4DynamicParticle; |
---|
1253 | theDeltaRay->SetKineticEnergy( DeltaKineticEnergy ); |
---|
1254 | theDeltaRay->SetMomentumDirection(DeltaDirection.x(), |
---|
1255 | DeltaDirection.y(), |
---|
1256 | DeltaDirection.z()); |
---|
1257 | theDeltaRay->SetDefinition(G4Electron::Electron()); |
---|
1258 | |
---|
1259 | // fill aParticleChange |
---|
1260 | finalKineticEnergy = KineticEnergy - DeltaKineticEnergy ; |
---|
1261 | |
---|
1262 | // Generation of Fluorescence and Auger |
---|
1263 | size_t nSecondaries = 0; |
---|
1264 | size_t totalNumber = 1; |
---|
1265 | std::vector<G4DynamicParticle*>* secondaryVector = 0; |
---|
1266 | G4DynamicParticle* aSecondary = 0; |
---|
1267 | G4ParticleDefinition* type = 0; |
---|
1268 | |
---|
1269 | // Select atom and shell |
---|
1270 | G4int Z = SelectRandomAtom(couple, KineticEnergy); |
---|
1271 | |
---|
1272 | // G4cout << "Fluorescence is switched :" << theFluo << G4endl; |
---|
1273 | |
---|
1274 | // Fluorescence data start from element 6 |
---|
1275 | if(theFluo && Z > 5) { |
---|
1276 | |
---|
1277 | |
---|
1278 | |
---|
1279 | // Atom total cross section |
---|
1280 | shellCS->SetTotalCS(totalCrossSectionMap[Z]); |
---|
1281 | |
---|
1282 | G4int shell = shellCS->SelectRandomShell(Z, KineticEnergy,ParticleMass,DeltaKineticEnergy); |
---|
1283 | |
---|
1284 | if (shell!=-1) { |
---|
1285 | |
---|
1286 | const G4AtomicShell* atomicShell = |
---|
1287 | (G4AtomicTransitionManager::Instance())->Shell(Z, shell); |
---|
1288 | G4double bindingEnergy = atomicShell->BindingEnergy(); |
---|
1289 | |
---|
1290 | if(verboseLevel > 1) { |
---|
1291 | G4cout << "PostStep Z= " << Z << " shell= " << shell |
---|
1292 | << " bindingE(keV)= " << bindingEnergy/keV |
---|
1293 | << " finalE(keV)= " << finalKineticEnergy/keV |
---|
1294 | << G4endl; |
---|
1295 | } |
---|
1296 | |
---|
1297 | |
---|
1298 | |
---|
1299 | if (finalKineticEnergy >= bindingEnergy |
---|
1300 | && (bindingEnergy >= minGammaEnergy |
---|
1301 | || bindingEnergy >= minElectronEnergy) ) { |
---|
1302 | |
---|
1303 | // G4int shellId = atomicShell->ShellId(); |
---|
1304 | deexcitationManager.GenerateParticles(secondaryVector, atomicShell, Z, minGammaEnergy, minElectronEnergy); |
---|
1305 | |
---|
1306 | if (secondaryVector != 0) { |
---|
1307 | // debug G4cout << "DEEXCTATION!!!" << G4endl; //debug |
---|
1308 | nSecondaries = secondaryVector->size(); |
---|
1309 | for (size_t i = 0; i<nSecondaries; i++) { |
---|
1310 | |
---|
1311 | aSecondary = (*secondaryVector)[i]; |
---|
1312 | if (aSecondary) { |
---|
1313 | |
---|
1314 | G4double e = aSecondary->GetKineticEnergy(); |
---|
1315 | type = aSecondary->GetDefinition(); |
---|
1316 | if (e < finalKineticEnergy && |
---|
1317 | ((type == G4Gamma::Gamma() && e > minGammaEnergy ) || |
---|
1318 | (type == G4Electron::Electron() && e > minElectronEnergy ))) { |
---|
1319 | |
---|
1320 | finalKineticEnergy -= e; |
---|
1321 | totalNumber++; |
---|
1322 | |
---|
1323 | } else { |
---|
1324 | |
---|
1325 | delete aSecondary; |
---|
1326 | (*secondaryVector)[i] = 0; |
---|
1327 | } |
---|
1328 | } |
---|
1329 | } |
---|
1330 | } |
---|
1331 | } |
---|
1332 | } |
---|
1333 | } |
---|
1334 | |
---|
1335 | // Save delta-electrons |
---|
1336 | |
---|
1337 | G4double edep = 0.0; |
---|
1338 | |
---|
1339 | if (finalKineticEnergy > MinKineticEnergy) |
---|
1340 | { |
---|
1341 | finalPx = TotalMomentum*ParticleDirection.x() |
---|
1342 | - DeltaTotalMomentum*DeltaDirection.x(); |
---|
1343 | finalPy = TotalMomentum*ParticleDirection.y() |
---|
1344 | - DeltaTotalMomentum*DeltaDirection.y(); |
---|
1345 | finalPz = TotalMomentum*ParticleDirection.z() |
---|
1346 | - DeltaTotalMomentum*DeltaDirection.z(); |
---|
1347 | finalMomentum = |
---|
1348 | std::sqrt(finalPx*finalPx+finalPy*finalPy+finalPz*finalPz) ; |
---|
1349 | finalPx /= finalMomentum ; |
---|
1350 | finalPy /= finalMomentum ; |
---|
1351 | finalPz /= finalMomentum ; |
---|
1352 | |
---|
1353 | aParticleChange.ProposeMomentumDirection( finalPx,finalPy,finalPz ); |
---|
1354 | } |
---|
1355 | else |
---|
1356 | { |
---|
1357 | edep = finalKineticEnergy; |
---|
1358 | finalKineticEnergy = 0.; |
---|
1359 | aParticleChange.ProposeMomentumDirection(ParticleDirection.x(), |
---|
1360 | ParticleDirection.y(),ParticleDirection.z()); |
---|
1361 | if(!aParticle->GetDefinition()->GetProcessManager()-> |
---|
1362 | GetAtRestProcessVector()->size()) |
---|
1363 | aParticleChange.ProposeTrackStatus(fStopAndKill); |
---|
1364 | else |
---|
1365 | aParticleChange.ProposeTrackStatus(fStopButAlive); |
---|
1366 | } |
---|
1367 | |
---|
1368 | aParticleChange.ProposeEnergy( finalKineticEnergy ); |
---|
1369 | aParticleChange.ProposeLocalEnergyDeposit (edep); |
---|
1370 | aParticleChange.SetNumberOfSecondaries(totalNumber); |
---|
1371 | aParticleChange.AddSecondary(theDeltaRay); |
---|
1372 | |
---|
1373 | // Save Fluorescence and Auger |
---|
1374 | |
---|
1375 | if (secondaryVector) { |
---|
1376 | |
---|
1377 | for (size_t l = 0; l < nSecondaries; l++) { |
---|
1378 | |
---|
1379 | aSecondary = (*secondaryVector)[l]; |
---|
1380 | if(aSecondary) { |
---|
1381 | aParticleChange.AddSecondary(aSecondary); |
---|
1382 | } |
---|
1383 | } |
---|
1384 | delete secondaryVector; |
---|
1385 | } |
---|
1386 | |
---|
1387 | return G4VContinuousDiscreteProcess::PostStepDoIt(trackData,stepData); |
---|
1388 | } |
---|
1389 | |
---|
1390 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1391 | |
---|
1392 | |
---|
1393 | |
---|
1394 | void G4hLowEnergyIonisation::SelectShellIonisationCS(G4String val) { |
---|
1395 | |
---|
1396 | if (val == "analytical" ) { |
---|
1397 | if (shellCS) delete shellCS; |
---|
1398 | shellCS = new G4teoCrossSection(val); |
---|
1399 | } |
---|
1400 | else if (val == "empirical") { |
---|
1401 | if (shellCS) delete shellCS; |
---|
1402 | shellCS = new G4empCrossSection(); |
---|
1403 | } |
---|
1404 | } |
---|
1405 | |
---|
1406 | |
---|
1407 | |
---|
1408 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1409 | |
---|
1410 | |
---|
1411 | std::vector<G4DynamicParticle*>* G4hLowEnergyIonisation::DeexciteAtom(const G4MaterialCutsCouple* couple, |
---|
1412 | G4double incidentEnergy, |
---|
1413 | G4double hMass, |
---|
1414 | G4double eLoss) |
---|
1415 | { |
---|
1416 | |
---|
1417 | if (verboseLevel > 1) { |
---|
1418 | G4cout << "DeexciteAtom: cutForPhotons(keV)= " << minGammaEnergy/keV |
---|
1419 | << " cutForElectrons(keV)= " << minElectronEnergy/keV |
---|
1420 | << " eLoss(MeV)= " << eLoss |
---|
1421 | << G4endl; |
---|
1422 | } |
---|
1423 | |
---|
1424 | |
---|
1425 | |
---|
1426 | if(eLoss < minGammaEnergy && eLoss < minElectronEnergy) return 0; |
---|
1427 | |
---|
1428 | const G4Material* material = couple->GetMaterial(); |
---|
1429 | G4int index = couple->GetIndex(); |
---|
1430 | // G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy(); |
---|
1431 | G4double gamma = incidentEnergy/hMass + 1; |
---|
1432 | G4double beta2 = 1.0 - 1.0/(gamma*gamma); |
---|
1433 | G4double r = electron_mass_c2/hMass; |
---|
1434 | G4double tmax = 2.*electron_mass_c2*(gamma*gamma - 1.)/(1. + 2.*gamma*r + r*r); |
---|
1435 | G4double tcut = std::min(tmax,cutForDelta[index]); |
---|
1436 | const G4AtomicTransitionManager* transitionManager = |
---|
1437 | G4AtomicTransitionManager::Instance(); |
---|
1438 | |
---|
1439 | size_t nElements = material->GetNumberOfElements(); |
---|
1440 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
1441 | G4bool stop = true; |
---|
1442 | |
---|
1443 | for (size_t j=0; j<nElements; j++) { |
---|
1444 | |
---|
1445 | G4int Z = (G4int)((*theElementVector)[j]->GetZ()); |
---|
1446 | G4double maxE = transitionManager->Shell(Z, 0)->BindingEnergy(); |
---|
1447 | |
---|
1448 | if (Z > 5 && maxE < tcut && (maxE > minGammaEnergy || maxE > minElectronEnergy) ) { |
---|
1449 | stop = false; |
---|
1450 | break; |
---|
1451 | } |
---|
1452 | } |
---|
1453 | |
---|
1454 | if(stop) return 0; |
---|
1455 | |
---|
1456 | // create vector of tracks of secondary particles |
---|
1457 | |
---|
1458 | std::vector<G4DynamicParticle*>* partVector = |
---|
1459 | new std::vector<G4DynamicParticle*>; |
---|
1460 | std::vector<G4DynamicParticle*>* secVector = new std::vector<G4DynamicParticle*>; |
---|
1461 | G4DynamicParticle* aSecondary = 0; |
---|
1462 | G4ParticleDefinition* type = 0; |
---|
1463 | G4double e, tkin, grej; |
---|
1464 | G4ThreeVector position; |
---|
1465 | G4int shell; |
---|
1466 | |
---|
1467 | // sample secondaries |
---|
1468 | |
---|
1469 | G4double etot = 0.0; |
---|
1470 | std::vector<G4int> n = shellVacancy->GenerateNumberOfIonisations(couple, |
---|
1471 | incidentEnergy, eLoss); |
---|
1472 | |
---|
1473 | for (size_t i=0; i<nElements; i++) { |
---|
1474 | |
---|
1475 | size_t nVacancies = n[i]; |
---|
1476 | G4int Z = (G4int)((*theElementVector)[i]->GetZ()); |
---|
1477 | G4double maxE = transitionManager->Shell(Z, 0)->BindingEnergy(); |
---|
1478 | |
---|
1479 | if (nVacancies && Z > 5 && maxE < tcut && (maxE > minGammaEnergy || maxE > minElectronEnergy)) { |
---|
1480 | for(size_t j=0; j<nVacancies; j++) { |
---|
1481 | |
---|
1482 | // sampling follows |
---|
1483 | do { |
---|
1484 | tkin = tcut/(1.0 + (tcut/maxE - 1.0)*G4UniformRand()); |
---|
1485 | grej = 1.0 - beta2 * tkin/tmax; |
---|
1486 | |
---|
1487 | } while( G4UniformRand() > grej ); |
---|
1488 | |
---|
1489 | // Atom total cross section |
---|
1490 | shellCS->SetTotalCS(totalCrossSectionMap[Z]); |
---|
1491 | |
---|
1492 | shell = shellCS->SelectRandomShell(Z,incidentEnergy,hMass,tkin); |
---|
1493 | |
---|
1494 | // shellId = transitionManager->Shell(Z, shell)->ShellId(); |
---|
1495 | G4double maxE = transitionManager->Shell(Z, shell)->BindingEnergy(); |
---|
1496 | |
---|
1497 | if (maxE>minGammaEnergy || maxE>minElectronEnergy ) |
---|
1498 | { |
---|
1499 | deexcitationManager.GenerateParticles(secVector, transitionManager->Shell(Z, shell), Z, minGammaEnergy, minElectronEnergy); |
---|
1500 | } |
---|
1501 | |
---|
1502 | if (!(secVector->empty())) { |
---|
1503 | size_t secN = secVector->size(); |
---|
1504 | for (size_t l = 0; l<secN; l++) { |
---|
1505 | |
---|
1506 | aSecondary = (*secVector)[l]; |
---|
1507 | if(aSecondary) { |
---|
1508 | |
---|
1509 | e = aSecondary->GetKineticEnergy(); |
---|
1510 | type = aSecondary->GetDefinition(); |
---|
1511 | if ( etot + e <= eLoss && |
---|
1512 | ( (type == G4Gamma::Gamma() && e > minGammaEnergy ) || |
---|
1513 | (type == G4Electron::Electron() && e > minElectronEnergy) ) ) |
---|
1514 | { |
---|
1515 | etot += e; |
---|
1516 | partVector->push_back(aSecondary); |
---|
1517 | } |
---|
1518 | else |
---|
1519 | { |
---|
1520 | delete aSecondary; |
---|
1521 | } |
---|
1522 | aSecondary = 0; |
---|
1523 | } |
---|
1524 | (*secVector)[l] = 0; |
---|
1525 | delete (*secVector)[l]; |
---|
1526 | } |
---|
1527 | |
---|
1528 | // secVector = 0; |
---|
1529 | |
---|
1530 | } |
---|
1531 | } |
---|
1532 | } |
---|
1533 | } |
---|
1534 | |
---|
1535 | delete secVector; |
---|
1536 | |
---|
1537 | if(partVector->empty()) { |
---|
1538 | delete partVector; |
---|
1539 | return 0; |
---|
1540 | } |
---|
1541 | |
---|
1542 | return partVector; |
---|
1543 | } |
---|
1544 | |
---|
1545 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1546 | |
---|
1547 | G4int G4hLowEnergyIonisation::SelectRandomAtom(const G4MaterialCutsCouple* couple, |
---|
1548 | G4double kineticEnergy) const |
---|
1549 | { |
---|
1550 | const G4Material* material = couple->GetMaterial(); |
---|
1551 | G4int nElements = material->GetNumberOfElements(); |
---|
1552 | G4int Z = 0; |
---|
1553 | |
---|
1554 | if(nElements == 1) { |
---|
1555 | Z = (G4int)(material->GetZ()); |
---|
1556 | return Z; |
---|
1557 | } |
---|
1558 | |
---|
1559 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
1560 | std::vector<G4double> p; |
---|
1561 | G4int index = couple->GetIndex(); |
---|
1562 | |
---|
1563 | G4double norm = 0.0; |
---|
1564 | for (G4int j=0; j<nElements; j++) { |
---|
1565 | |
---|
1566 | const G4VEMDataSet* set = (zFluoDataVector[index])->GetComponent(j); |
---|
1567 | G4double cross = set->FindValue(kineticEnergy); |
---|
1568 | |
---|
1569 | p.push_back(cross); |
---|
1570 | norm += cross; |
---|
1571 | } |
---|
1572 | |
---|
1573 | if(norm == 0.0) return 0; |
---|
1574 | |
---|
1575 | G4double q = norm*G4UniformRand(); |
---|
1576 | |
---|
1577 | for (G4int i=0; i<nElements; i++) { |
---|
1578 | |
---|
1579 | if(p[i] > q) { |
---|
1580 | Z = (G4int)((*theElementVector)[i]->GetZ()); |
---|
1581 | break; |
---|
1582 | } |
---|
1583 | q -= p[i]; |
---|
1584 | } |
---|
1585 | |
---|
1586 | return Z; |
---|
1587 | } |
---|
1588 | |
---|
1589 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1590 | |
---|
1591 | G4double G4hLowEnergyIonisation::ComputeDEDX( |
---|
1592 | const G4ParticleDefinition* aParticle, |
---|
1593 | const G4MaterialCutsCouple* couple, |
---|
1594 | G4double kineticEnergy) |
---|
1595 | { |
---|
1596 | const G4Material* material = couple->GetMaterial(); |
---|
1597 | G4Proton* theProton = G4Proton::Proton(); |
---|
1598 | G4AntiProton* theAntiProton = G4AntiProton::AntiProton(); |
---|
1599 | G4double dedx = 0.0 ; |
---|
1600 | |
---|
1601 | G4double tscaled = kineticEnergy*proton_mass_c2/(aParticle->GetPDGMass()) ; |
---|
1602 | charge = aParticle->GetPDGCharge() ; |
---|
1603 | |
---|
1604 | if(charge>0.0) { |
---|
1605 | if(tscaled > protonHighEnergy) { |
---|
1606 | dedx=G4EnergyLossTables::GetDEDX(theProton,tscaled,couple) ; |
---|
1607 | |
---|
1608 | } else { |
---|
1609 | dedx=ProtonParametrisedDEDX(couple,tscaled) ; |
---|
1610 | } |
---|
1611 | |
---|
1612 | } else { |
---|
1613 | if(tscaled > antiProtonHighEnergy) { |
---|
1614 | dedx=G4EnergyLossTables::GetDEDX(theAntiProton,tscaled,couple); |
---|
1615 | |
---|
1616 | } else { |
---|
1617 | dedx=AntiProtonParametrisedDEDX(couple,tscaled) ; |
---|
1618 | } |
---|
1619 | } |
---|
1620 | dedx *= theIonEffChargeModel->TheValue(aParticle, material, kineticEnergy) ; |
---|
1621 | |
---|
1622 | return dedx ; |
---|
1623 | } |
---|
1624 | |
---|
1625 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1626 | |
---|
1627 | G4double G4hLowEnergyIonisation::BarkasTerm(const G4Material* material, |
---|
1628 | G4double kineticEnergy) const |
---|
1629 | //Function to compute the Barkas term for protons: |
---|
1630 | // |
---|
1631 | //Ref. Z_1^3 effect in the stopping power of matter for charged particles |
---|
1632 | // J.C Ashley and R.H.Ritchie |
---|
1633 | // Physical review B Vol.5 No.7 1 April 1972 pagg. 2393-2397 |
---|
1634 | // |
---|
1635 | { |
---|
1636 | static double FTable[47][2] = { |
---|
1637 | { 0.02, 21.5}, |
---|
1638 | { 0.03, 20.0}, |
---|
1639 | { 0.04, 18.0}, |
---|
1640 | { 0.05, 15.6}, |
---|
1641 | { 0.06, 15.0}, |
---|
1642 | { 0.07, 14.0}, |
---|
1643 | { 0.08, 13.5}, |
---|
1644 | { 0.09, 13.}, |
---|
1645 | { 0.1, 12.2}, |
---|
1646 | { 0.2, 9.25}, |
---|
1647 | { 0.3, 7.0}, |
---|
1648 | { 0.4, 6.0}, |
---|
1649 | { 0.5, 4.5}, |
---|
1650 | { 0.6, 3.5}, |
---|
1651 | { 0.7, 3.0}, |
---|
1652 | { 0.8, 2.5}, |
---|
1653 | { 0.9, 2.0}, |
---|
1654 | { 1.0, 1.7}, |
---|
1655 | { 1.2, 1.2}, |
---|
1656 | { 1.3, 1.0}, |
---|
1657 | { 1.4, 0.86}, |
---|
1658 | { 1.5, 0.7}, |
---|
1659 | { 1.6, 0.61}, |
---|
1660 | { 1.7, 0.52}, |
---|
1661 | { 1.8, 0.5}, |
---|
1662 | { 1.9, 0.43}, |
---|
1663 | { 2.0, 0.42}, |
---|
1664 | { 2.1, 0.3}, |
---|
1665 | { 2.4, 0.2}, |
---|
1666 | { 3.0, 0.13}, |
---|
1667 | { 3.08, 0.1}, |
---|
1668 | { 3.1, 0.09}, |
---|
1669 | { 3.3, 0.08}, |
---|
1670 | { 3.5, 0.07}, |
---|
1671 | { 3.8, 0.06}, |
---|
1672 | { 4.0, 0.051}, |
---|
1673 | { 4.1, 0.04}, |
---|
1674 | { 4.8, 0.03}, |
---|
1675 | { 5.0, 0.024}, |
---|
1676 | { 5.1, 0.02}, |
---|
1677 | { 6.0, 0.013}, |
---|
1678 | { 6.5, 0.01}, |
---|
1679 | { 7.0, 0.009}, |
---|
1680 | { 7.1, 0.008}, |
---|
1681 | { 8.0, 0.006}, |
---|
1682 | { 9.0, 0.0032}, |
---|
1683 | { 10.0, 0.0025} }; |
---|
1684 | |
---|
1685 | // Information on particle and material |
---|
1686 | G4double kinE = kineticEnergy ; |
---|
1687 | if(0.5*MeV > kinE) kinE = 0.5*MeV ; |
---|
1688 | G4double gamma = 1.0 + kinE / proton_mass_c2 ; |
---|
1689 | G4double beta2 = 1.0 - 1.0/(gamma*gamma) ; |
---|
1690 | if(0.0 >= beta2) return 0.0; |
---|
1691 | |
---|
1692 | G4double BarkasTerm = 0.0; |
---|
1693 | G4double AMaterial = 0.0; |
---|
1694 | G4double ZMaterial = 0.0; |
---|
1695 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
1696 | G4int numberOfElements = material->GetNumberOfElements(); |
---|
1697 | |
---|
1698 | for (G4int i = 0; i<numberOfElements; i++) { |
---|
1699 | |
---|
1700 | AMaterial = (*theElementVector)[i]->GetA()*mole/g; |
---|
1701 | ZMaterial = (*theElementVector)[i]->GetZ(); |
---|
1702 | |
---|
1703 | G4double X = 137.0 * 137.0 * beta2 / ZMaterial; |
---|
1704 | |
---|
1705 | // Variables to compute L_1 |
---|
1706 | G4double Eta0Chi = 0.8; |
---|
1707 | G4double EtaChi = Eta0Chi * ( 1.0 + 6.02*std::pow( ZMaterial,-1.19 ) ); |
---|
1708 | G4double W = ( EtaChi * std::pow( ZMaterial,1.0/6.0 ) ) / std::sqrt(X); |
---|
1709 | G4double FunctionOfW = FTable[46][1]*FTable[46][0]/W ; |
---|
1710 | |
---|
1711 | for(G4int j=0; j<47; j++) { |
---|
1712 | |
---|
1713 | if( W < FTable[j][0] ) { |
---|
1714 | |
---|
1715 | if(0 == j) { |
---|
1716 | FunctionOfW = FTable[0][1] ; |
---|
1717 | |
---|
1718 | } else { |
---|
1719 | FunctionOfW = (FTable[j][1] - FTable[j-1][1]) * (W - FTable[j-1][0]) |
---|
1720 | / (FTable[j][0] - FTable[j-1][0]) |
---|
1721 | + FTable[j-1][1] ; |
---|
1722 | } |
---|
1723 | |
---|
1724 | break; |
---|
1725 | } |
---|
1726 | |
---|
1727 | } |
---|
1728 | |
---|
1729 | BarkasTerm += FunctionOfW /( std::sqrt(ZMaterial * X) * X); |
---|
1730 | } |
---|
1731 | |
---|
1732 | BarkasTerm *= twopi_mc2_rcl2 * (material->GetElectronDensity()) / beta2 ; |
---|
1733 | |
---|
1734 | return BarkasTerm; |
---|
1735 | } |
---|
1736 | |
---|
1737 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1738 | |
---|
1739 | G4double G4hLowEnergyIonisation::BlochTerm(const G4Material* material, |
---|
1740 | G4double kineticEnergy, |
---|
1741 | G4double cSquare) const |
---|
1742 | //Function to compute the Bloch term for protons: |
---|
1743 | // |
---|
1744 | //Ref. Z_1^3 effect in the stopping power of matter for charged particles |
---|
1745 | // J.C Ashley and R.H.Ritchie |
---|
1746 | // Physical review B Vol.5 No.7 1 April 1972 pagg. 2393-2397 |
---|
1747 | // |
---|
1748 | { |
---|
1749 | G4double eloss = 0.0 ; |
---|
1750 | G4double gamma = 1.0 + kineticEnergy / proton_mass_c2 ; |
---|
1751 | G4double beta2 = 1.0 - 1.0/(gamma*gamma) ; |
---|
1752 | G4double y = cSquare / (137.0*137.0*beta2) ; |
---|
1753 | |
---|
1754 | if(y < 0.05) { |
---|
1755 | eloss = 1.202 ; |
---|
1756 | |
---|
1757 | } else { |
---|
1758 | eloss = 1.0 / (1.0 + y) ; |
---|
1759 | G4double de = eloss ; |
---|
1760 | |
---|
1761 | for(G4int i=2; de>eloss*0.01; i++) { |
---|
1762 | de = 1.0/( i * (i*i + y)) ; |
---|
1763 | eloss += de ; |
---|
1764 | } |
---|
1765 | } |
---|
1766 | eloss *= -1.0 * y * cSquare * twopi_mc2_rcl2 * |
---|
1767 | (material->GetElectronDensity()) / beta2 ; |
---|
1768 | |
---|
1769 | return eloss; |
---|
1770 | } |
---|
1771 | |
---|
1772 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1773 | |
---|
1774 | G4double G4hLowEnergyIonisation::ElectronicLossFluctuation( |
---|
1775 | const G4DynamicParticle* particle, |
---|
1776 | const G4MaterialCutsCouple* couple, |
---|
1777 | G4double meanLoss, |
---|
1778 | G4double step) const |
---|
1779 | // calculate actual loss from the mean loss |
---|
1780 | // The model used to get the fluctuation is essentially the same |
---|
1781 | // as in Glandz in Geant3. |
---|
1782 | { |
---|
1783 | // data members to speed up the fluctuation calculation |
---|
1784 | // G4int imat ; |
---|
1785 | // G4double f1Fluct,f2Fluct,e1Fluct,e2Fluct,rateFluct,ipotFluct; |
---|
1786 | // G4double e1LogFluct,e2LogFluct,ipotLogFluct; |
---|
1787 | |
---|
1788 | static const G4double minLoss = 1.*eV ; |
---|
1789 | static const G4double kappa = 10. ; |
---|
1790 | static const G4double theBohrBeta2 = 50.0 * keV/proton_mass_c2 ; |
---|
1791 | |
---|
1792 | const G4Material* material = couple->GetMaterial(); |
---|
1793 | G4int imaterial = couple->GetIndex() ; |
---|
1794 | G4double ipotFluct = material->GetIonisation()->GetMeanExcitationEnergy() ; |
---|
1795 | G4double electronDensity = material->GetElectronDensity() ; |
---|
1796 | G4double zeff = electronDensity/(material->GetTotNbOfAtomsPerVolume()) ; |
---|
1797 | |
---|
1798 | // get particle data |
---|
1799 | G4double tkin = particle->GetKineticEnergy(); |
---|
1800 | G4double particleMass = particle->GetMass() ; |
---|
1801 | G4double deltaCutInKineticEnergyNow = cutForDelta[imaterial]; |
---|
1802 | |
---|
1803 | // shortcut for very very small loss |
---|
1804 | if(meanLoss < minLoss) return meanLoss ; |
---|
1805 | |
---|
1806 | // Validity range for delta electron cross section |
---|
1807 | G4double threshold = std::max(deltaCutInKineticEnergyNow,ipotFluct); |
---|
1808 | G4double loss, siga; |
---|
1809 | |
---|
1810 | G4double rmass = electron_mass_c2/particleMass; |
---|
1811 | G4double tau = tkin/particleMass; |
---|
1812 | G4double tau1 = tau+1.0; |
---|
1813 | G4double tau2 = tau*(tau+2.); |
---|
1814 | G4double tmax = 2.*electron_mass_c2*tau2/(1.+2.*tau1*rmass+rmass*rmass); |
---|
1815 | |
---|
1816 | |
---|
1817 | if(tmax > threshold) tmax = threshold; |
---|
1818 | G4double beta2 = tau2/(tau1*tau1); |
---|
1819 | |
---|
1820 | // Gaussian fluctuation |
---|
1821 | if(meanLoss > kappa*tmax || tmax < kappa*ipotFluct ) |
---|
1822 | { |
---|
1823 | siga = tmax * (1.0-0.5*beta2) * step * twopi_mc2_rcl2 |
---|
1824 | * electronDensity / beta2 ; |
---|
1825 | |
---|
1826 | // High velocity or negatively charged particle |
---|
1827 | if( beta2 > 3.0*theBohrBeta2*zeff || charge < 0.0) { |
---|
1828 | siga = std::sqrt( siga * chargeSquare ) ; |
---|
1829 | |
---|
1830 | // Low velocity - additional ion charge fluctuations according to |
---|
1831 | // Q.Yang et al., NIM B61(1991)149-155. |
---|
1832 | } else { |
---|
1833 | G4double chu = theIonChuFluctuationModel->TheValue(particle, material); |
---|
1834 | G4double yang = theIonYangFluctuationModel->TheValue(particle, material); |
---|
1835 | siga = std::sqrt( siga * (chargeSquare * chu + yang)) ; |
---|
1836 | } |
---|
1837 | |
---|
1838 | do { |
---|
1839 | loss = G4RandGauss::shoot(meanLoss,siga); |
---|
1840 | } while (loss < 0. || loss > 2.0*meanLoss); |
---|
1841 | return loss; |
---|
1842 | } |
---|
1843 | |
---|
1844 | // Non Gaussian fluctuation |
---|
1845 | static const G4double probLim = 0.01 ; |
---|
1846 | static const G4double sumaLim = -std::log(probLim) ; |
---|
1847 | static const G4double alim = 10.; |
---|
1848 | |
---|
1849 | G4double suma,w1,w2,C,e0,lossc,w; |
---|
1850 | G4double a1,a2,a3; |
---|
1851 | G4int p1,p2,p3; |
---|
1852 | G4int nb; |
---|
1853 | G4double corrfac, na,alfa,rfac,namean,sa,alfa1,ea,sea; |
---|
1854 | G4double dp3; |
---|
1855 | |
---|
1856 | G4double f1Fluct = material->GetIonisation()->GetF1fluct(); |
---|
1857 | G4double f2Fluct = material->GetIonisation()->GetF2fluct(); |
---|
1858 | G4double e1Fluct = material->GetIonisation()->GetEnergy1fluct(); |
---|
1859 | G4double e2Fluct = material->GetIonisation()->GetEnergy2fluct(); |
---|
1860 | G4double e1LogFluct = material->GetIonisation()->GetLogEnergy1fluct(); |
---|
1861 | G4double e2LogFluct = material->GetIonisation()->GetLogEnergy2fluct(); |
---|
1862 | G4double rateFluct = material->GetIonisation()->GetRateionexcfluct(); |
---|
1863 | G4double ipotLogFluct= material->GetIonisation()->GetLogMeanExcEnergy(); |
---|
1864 | |
---|
1865 | w1 = tmax/ipotFluct; |
---|
1866 | w2 = std::log(2.*electron_mass_c2*tau2); |
---|
1867 | |
---|
1868 | C = meanLoss*(1.-rateFluct)/(w2-ipotLogFluct-beta2); |
---|
1869 | |
---|
1870 | a1 = C*f1Fluct*(w2-e1LogFluct-beta2)/e1Fluct; |
---|
1871 | a2 = C*f2Fluct*(w2-e2LogFluct-beta2)/e2Fluct; |
---|
1872 | a3 = rateFluct*meanLoss*(tmax-ipotFluct)/(ipotFluct*tmax*std::log(w1)); |
---|
1873 | if(a1 < 0.0) a1 = 0.0; |
---|
1874 | if(a2 < 0.0) a2 = 0.0; |
---|
1875 | if(a3 < 0.0) a3 = 0.0; |
---|
1876 | |
---|
1877 | suma = a1+a2+a3; |
---|
1878 | |
---|
1879 | loss = 0.; |
---|
1880 | |
---|
1881 | |
---|
1882 | if(suma < sumaLim) // very small Step |
---|
1883 | { |
---|
1884 | e0 = material->GetIonisation()->GetEnergy0fluct(); |
---|
1885 | |
---|
1886 | if(tmax == ipotFluct) |
---|
1887 | { |
---|
1888 | a3 = meanLoss/e0; |
---|
1889 | |
---|
1890 | if(a3>alim) |
---|
1891 | { |
---|
1892 | siga=std::sqrt(a3) ; |
---|
1893 | p3 = std::max(0,G4int(G4RandGauss::shoot(a3,siga)+0.5)); |
---|
1894 | } |
---|
1895 | else |
---|
1896 | p3 = G4Poisson(a3); |
---|
1897 | |
---|
1898 | loss = p3*e0 ; |
---|
1899 | |
---|
1900 | if(p3 > 0) |
---|
1901 | loss += (1.-2.*G4UniformRand())*e0 ; |
---|
1902 | |
---|
1903 | } |
---|
1904 | else |
---|
1905 | { |
---|
1906 | tmax = tmax-ipotFluct+e0 ; |
---|
1907 | a3 = meanLoss*(tmax-e0)/(tmax*e0*std::log(tmax/e0)); |
---|
1908 | |
---|
1909 | if(a3>alim) |
---|
1910 | { |
---|
1911 | siga=std::sqrt(a3) ; |
---|
1912 | p3 = std::max(0,int(G4RandGauss::shoot(a3,siga)+0.5)); |
---|
1913 | } |
---|
1914 | else |
---|
1915 | p3 = G4Poisson(a3); |
---|
1916 | |
---|
1917 | if(p3 > 0) |
---|
1918 | { |
---|
1919 | w = (tmax-e0)/tmax ; |
---|
1920 | if(p3 > nmaxCont2) |
---|
1921 | { |
---|
1922 | dp3 = G4float(p3) ; |
---|
1923 | corrfac = dp3/G4float(nmaxCont2) ; |
---|
1924 | p3 = nmaxCont2 ; |
---|
1925 | } |
---|
1926 | else |
---|
1927 | corrfac = 1. ; |
---|
1928 | |
---|
1929 | for(G4int i=0; i<p3; i++) loss += 1./(1.-w*G4UniformRand()) ; |
---|
1930 | loss *= e0*corrfac ; |
---|
1931 | } |
---|
1932 | } |
---|
1933 | } |
---|
1934 | |
---|
1935 | else // not so small Step |
---|
1936 | { |
---|
1937 | // excitation type 1 |
---|
1938 | if(a1>alim) |
---|
1939 | { |
---|
1940 | siga=std::sqrt(a1) ; |
---|
1941 | p1 = std::max(0,G4int(G4RandGauss::shoot(a1,siga)+0.5)); |
---|
1942 | } |
---|
1943 | else |
---|
1944 | p1 = G4Poisson(a1); |
---|
1945 | |
---|
1946 | // excitation type 2 |
---|
1947 | if(a2>alim) |
---|
1948 | { |
---|
1949 | siga=std::sqrt(a2) ; |
---|
1950 | p2 = std::max(0,G4int(G4RandGauss::shoot(a2,siga)+0.5)); |
---|
1951 | } |
---|
1952 | else |
---|
1953 | p2 = G4Poisson(a2); |
---|
1954 | |
---|
1955 | loss = p1*e1Fluct+p2*e2Fluct; |
---|
1956 | |
---|
1957 | // smearing to avoid unphysical peaks |
---|
1958 | if(p2 > 0) |
---|
1959 | loss += (1.-2.*G4UniformRand())*e2Fluct; |
---|
1960 | else if (loss>0.) |
---|
1961 | loss += (1.-2.*G4UniformRand())*e1Fluct; |
---|
1962 | |
---|
1963 | // ionisation ....................................... |
---|
1964 | if(a3 > 0.) |
---|
1965 | { |
---|
1966 | if(a3>alim) |
---|
1967 | { |
---|
1968 | siga=std::sqrt(a3) ; |
---|
1969 | p3 = std::max(0,G4int(G4RandGauss::shoot(a3,siga)+0.5)); |
---|
1970 | } |
---|
1971 | else |
---|
1972 | p3 = G4Poisson(a3); |
---|
1973 | |
---|
1974 | lossc = 0.; |
---|
1975 | if(p3 > 0) |
---|
1976 | { |
---|
1977 | na = 0.; |
---|
1978 | alfa = 1.; |
---|
1979 | if (p3 > nmaxCont2) |
---|
1980 | { |
---|
1981 | dp3 = G4float(p3); |
---|
1982 | rfac = dp3/(G4float(nmaxCont2)+dp3); |
---|
1983 | namean = G4float(p3)*rfac; |
---|
1984 | sa = G4float(nmaxCont1)*rfac; |
---|
1985 | na = G4RandGauss::shoot(namean,sa); |
---|
1986 | if (na > 0.) |
---|
1987 | { |
---|
1988 | alfa = w1*G4float(nmaxCont2+p3)/ |
---|
1989 | (w1*G4float(nmaxCont2)+G4float(p3)); |
---|
1990 | alfa1 = alfa*std::log(alfa)/(alfa-1.); |
---|
1991 | ea = na*ipotFluct*alfa1; |
---|
1992 | sea = ipotFluct*std::sqrt(na*(alfa-alfa1*alfa1)); |
---|
1993 | lossc += G4RandGauss::shoot(ea,sea); |
---|
1994 | } |
---|
1995 | } |
---|
1996 | |
---|
1997 | nb = G4int(G4float(p3)-na); |
---|
1998 | if (nb > 0) |
---|
1999 | { |
---|
2000 | w2 = alfa*ipotFluct; |
---|
2001 | w = (tmax-w2)/tmax; |
---|
2002 | for (G4int k=0; k<nb; k++) lossc += w2/(1.-w*G4UniformRand()); |
---|
2003 | } |
---|
2004 | } |
---|
2005 | loss += lossc; |
---|
2006 | } |
---|
2007 | } |
---|
2008 | |
---|
2009 | return loss ; |
---|
2010 | } |
---|
2011 | |
---|
2012 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
2013 | |
---|
2014 | void G4hLowEnergyIonisation::SetCutForSecondaryPhotons(G4double cut) |
---|
2015 | { |
---|
2016 | minGammaEnergy = cut; |
---|
2017 | } |
---|
2018 | |
---|
2019 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
2020 | |
---|
2021 | void G4hLowEnergyIonisation::SetCutForAugerElectrons(G4double cut) |
---|
2022 | { |
---|
2023 | minElectronEnergy = cut; |
---|
2024 | } |
---|
2025 | |
---|
2026 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
2027 | |
---|
2028 | void G4hLowEnergyIonisation::ActivateAugerElectronProduction(G4bool val) |
---|
2029 | { |
---|
2030 | deexcitationManager.SetAugerActive(val); |
---|
2031 | } |
---|
2032 | |
---|
2033 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
2034 | |
---|
2035 | void G4hLowEnergyIonisation::PrintInfoDefinition() const |
---|
2036 | { |
---|
2037 | G4String comments = " Knock-on electron cross sections . "; |
---|
2038 | comments += "\n Good description above the mean excitation energy.\n"; |
---|
2039 | comments += " Delta ray energy sampled from differential Xsection."; |
---|
2040 | |
---|
2041 | G4cout << G4endl << GetProcessName() << ": " << comments |
---|
2042 | << "\n PhysicsTables from " << LowestKineticEnergy / eV << " eV " |
---|
2043 | << " to " << HighestKineticEnergy / TeV << " TeV " |
---|
2044 | << " in " << TotBin << " bins." |
---|
2045 | << "\n Electronic stopping power model is " |
---|
2046 | << theProtonTable |
---|
2047 | << "\n from " << protonLowEnergy / keV << " keV " |
---|
2048 | << " to " << protonHighEnergy / MeV << " MeV " << "." << G4endl ; |
---|
2049 | G4cout << "\n Parametrisation model for antiprotons is " |
---|
2050 | << theAntiProtonTable |
---|
2051 | << "\n from " << antiProtonLowEnergy / keV << " keV " |
---|
2052 | << " to " << antiProtonHighEnergy / MeV << " MeV " << "." << G4endl ; |
---|
2053 | if(theBarkas){ |
---|
2054 | G4cout << " Parametrization of the Barkas effect is switched on." |
---|
2055 | << G4endl ; |
---|
2056 | } |
---|
2057 | if(nStopping) { |
---|
2058 | G4cout << " Nuclear stopping power model is " << theNuclearTable |
---|
2059 | << G4endl ; |
---|
2060 | } |
---|
2061 | |
---|
2062 | G4bool printHead = true; |
---|
2063 | |
---|
2064 | const G4ProductionCutsTable* theCoupleTable= |
---|
2065 | G4ProductionCutsTable::GetProductionCutsTable(); |
---|
2066 | size_t numOfCouples = theCoupleTable->GetTableSize(); |
---|
2067 | |
---|
2068 | // loop for materials |
---|
2069 | |
---|
2070 | for (size_t j=0 ; j < numOfCouples; j++) { |
---|
2071 | |
---|
2072 | const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(j); |
---|
2073 | const G4Material* material= couple->GetMaterial(); |
---|
2074 | G4double deltaCutNow = cutForDelta[(couple->GetIndex())] ; |
---|
2075 | G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy(); |
---|
2076 | |
---|
2077 | if(eexc > deltaCutNow) { |
---|
2078 | if(printHead) { |
---|
2079 | printHead = false ; |
---|
2080 | |
---|
2081 | G4cout << " material min.delta energy(keV) " << G4endl; |
---|
2082 | G4cout << G4endl; |
---|
2083 | } |
---|
2084 | |
---|
2085 | G4cout << std::setw(20) << material->GetName() |
---|
2086 | << std::setw(15) << eexc/keV << G4endl; |
---|
2087 | } |
---|
2088 | } |
---|
2089 | } |
---|