1 | // |
---|
2 | // ******************************************************************** |
---|
3 | // * License and Disclaimer * |
---|
4 | // * * |
---|
5 | // * The Geant4 software is copyright of the Copyright Holders of * |
---|
6 | // * the Geant4 Collaboration. It is provided under the terms and * |
---|
7 | // * conditions of the Geant4 Software License, included in the file * |
---|
8 | // * LICENSE and available at http://cern.ch/geant4/license . These * |
---|
9 | // * include a list of copyright holders. * |
---|
10 | // * * |
---|
11 | // * Neither the authors of this software system, nor their employing * |
---|
12 | // * institutes,nor the agencies providing financial support for this * |
---|
13 | // * work make any representation or warranty, express or implied, * |
---|
14 | // * regarding this software system or assume any liability for its * |
---|
15 | // * use. Please see the license in the file LICENSE and URL above * |
---|
16 | // * for the full disclaimer and the limitation of liability. * |
---|
17 | // * * |
---|
18 | // * This code implementation is the result of the scientific and * |
---|
19 | // * technical work of the GEANT4 collaboration. * |
---|
20 | // * By using, copying, modifying or distributing the software (or * |
---|
21 | // * any work based on the software) you agree to acknowledge its * |
---|
22 | // * use in resulting scientific publications, and indicate your * |
---|
23 | // * acceptance of all terms of the Geant4 Software license. * |
---|
24 | // ******************************************************************** |
---|
25 | // |
---|
26 | // $Id: G4VEmModel.hh,v 1.72 2009/09/23 14:42:47 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-03 $ |
---|
28 | // |
---|
29 | // ------------------------------------------------------------------- |
---|
30 | // |
---|
31 | // GEANT4 Class header file |
---|
32 | // |
---|
33 | // |
---|
34 | // File name: G4VEmModel |
---|
35 | // |
---|
36 | // Author: Vladimir Ivanchenko |
---|
37 | // |
---|
38 | // Creation date: 03.01.2002 |
---|
39 | // |
---|
40 | // Modifications: |
---|
41 | // |
---|
42 | // 23-12-02 V.Ivanchenko change interface before move to cut per region |
---|
43 | // 24-01-03 Cut per region (V.Ivanchenko) |
---|
44 | // 13-02-03 Add name (V.Ivanchenko) |
---|
45 | // 25-02-03 Add sample theta and displacement (V.Ivanchenko) |
---|
46 | // 23-07-03 Replace G4Material by G4MaterialCutCouple in dE/dx and CrossSection |
---|
47 | // calculation (V.Ivanchenko) |
---|
48 | // 01-03-04 L.Urban signature changed in SampleCosineTheta |
---|
49 | // 23-04-04 L.urban signature of SampleCosineTheta changed back |
---|
50 | // 17-11-04 Add method CrossSectionPerAtom (V.Ivanchenko) |
---|
51 | // 14-03-05 Reduce number of pure virtual methods and make inline part |
---|
52 | // separate (V.Ivanchenko) |
---|
53 | // 24-03-05 Remove IsInCharge and add G4VParticleChange in the constructor (VI) |
---|
54 | // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko) |
---|
55 | // 15-04-05 optimize internal interface for msc (V.Ivanchenko) |
---|
56 | // 08-05-05 A -> N (V.Ivanchenko) |
---|
57 | // 25-07-05 Move constructor and destructor to the body (V.Ivanchenko) |
---|
58 | // 02-02-06 ComputeCrossSectionPerAtom: default value A=0. (mma) |
---|
59 | // 06-02-06 add method ComputeMeanFreePath() (mma) |
---|
60 | // 07-03-06 Optimize msc methods (V.Ivanchenko) |
---|
61 | // 29-06-06 Add member currentElement and Get/Set methods (V.Ivanchenko) |
---|
62 | // 29-10-07 Added SampleScattering (V.Ivanchenko) |
---|
63 | // 15-07-08 Reorder class members and improve comments (VI) |
---|
64 | // 21-07-08 Added vector of G4ElementSelector and methods to use it (VI) |
---|
65 | // 12-09-08 Added methods GetParticleCharge, GetChargeSquareRatio, |
---|
66 | // CorrectionsAlongStep, ActivateNuclearStopping (VI) |
---|
67 | // 16-02-09 Moved implementations of virtual methods to source (VI) |
---|
68 | // 07-04-09 Moved msc methods from G4VEmModel to G4VMscModel (VI) |
---|
69 | // |
---|
70 | // Class Description: |
---|
71 | // |
---|
72 | // Abstract interface to energy loss models |
---|
73 | |
---|
74 | // ------------------------------------------------------------------- |
---|
75 | // |
---|
76 | |
---|
77 | #ifndef G4VEmModel_h |
---|
78 | #define G4VEmModel_h 1 |
---|
79 | |
---|
80 | #include "globals.hh" |
---|
81 | #include "G4DynamicParticle.hh" |
---|
82 | #include "G4ParticleDefinition.hh" |
---|
83 | #include "G4MaterialCutsCouple.hh" |
---|
84 | #include "G4Material.hh" |
---|
85 | #include "G4Element.hh" |
---|
86 | #include "G4ElementVector.hh" |
---|
87 | #include "G4DataVector.hh" |
---|
88 | #include "G4VEmFluctuationModel.hh" |
---|
89 | #include "G4EmElementSelector.hh" |
---|
90 | #include "Randomize.hh" |
---|
91 | #include <vector> |
---|
92 | |
---|
93 | class G4PhysicsTable; |
---|
94 | class G4Region; |
---|
95 | class G4VParticleChange; |
---|
96 | class G4ParticleChangeForLoss; |
---|
97 | class G4ParticleChangeForGamma; |
---|
98 | class G4Track; |
---|
99 | |
---|
100 | class G4VEmModel |
---|
101 | { |
---|
102 | |
---|
103 | public: |
---|
104 | |
---|
105 | G4VEmModel(const G4String& nam); |
---|
106 | |
---|
107 | virtual ~G4VEmModel(); |
---|
108 | |
---|
109 | //------------------------------------------------------------------------ |
---|
110 | // Virtual methods to be implemented for any concrete model |
---|
111 | //------------------------------------------------------------------------ |
---|
112 | |
---|
113 | virtual void Initialise(const G4ParticleDefinition*, |
---|
114 | const G4DataVector&) = 0; |
---|
115 | |
---|
116 | virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*, |
---|
117 | const G4MaterialCutsCouple*, |
---|
118 | const G4DynamicParticle*, |
---|
119 | G4double tmin = 0.0, |
---|
120 | G4double tmax = DBL_MAX) = 0; |
---|
121 | |
---|
122 | //------------------------------------------------------------------------ |
---|
123 | // Methods with standard implementation; may be overwritten if needed |
---|
124 | //------------------------------------------------------------------------ |
---|
125 | |
---|
126 | // main method to compute dEdx |
---|
127 | virtual G4double ComputeDEDXPerVolume(const G4Material*, |
---|
128 | const G4ParticleDefinition*, |
---|
129 | G4double kineticEnergy, |
---|
130 | G4double cutEnergy = DBL_MAX); |
---|
131 | |
---|
132 | // main method to compute cross section per Volume |
---|
133 | virtual G4double CrossSectionPerVolume(const G4Material*, |
---|
134 | const G4ParticleDefinition*, |
---|
135 | G4double kineticEnergy, |
---|
136 | G4double cutEnergy = 0.0, |
---|
137 | G4double maxEnergy = DBL_MAX); |
---|
138 | |
---|
139 | // main method to compute cross section per atom |
---|
140 | virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*, |
---|
141 | G4double kinEnergy, |
---|
142 | G4double Z, |
---|
143 | G4double A = 0., /* amu */ |
---|
144 | G4double cutEnergy = 0.0, |
---|
145 | G4double maxEnergy = DBL_MAX); |
---|
146 | |
---|
147 | // min cut in kinetic energy allowed by the model |
---|
148 | virtual G4double MinEnergyCut(const G4ParticleDefinition*, |
---|
149 | const G4MaterialCutsCouple*); |
---|
150 | |
---|
151 | // Compute effective ion charge square |
---|
152 | virtual G4double GetChargeSquareRatio(const G4ParticleDefinition*, |
---|
153 | const G4Material*, |
---|
154 | G4double kineticEnergy); |
---|
155 | |
---|
156 | // Compute ion charge |
---|
157 | virtual G4double GetParticleCharge(const G4ParticleDefinition*, |
---|
158 | const G4Material*, |
---|
159 | G4double kineticEnergy); |
---|
160 | |
---|
161 | // add correction to energy loss and compute non-ionizing energy loss |
---|
162 | virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*, |
---|
163 | const G4DynamicParticle*, |
---|
164 | G4double& eloss, |
---|
165 | G4double& niel, |
---|
166 | G4double length); |
---|
167 | |
---|
168 | // sample PIXE deexcitation |
---|
169 | virtual void SampleDeexcitationAlongStep(const G4Material*, |
---|
170 | const G4Track&, |
---|
171 | G4double& eloss); |
---|
172 | |
---|
173 | // initilisation at run time for a given material |
---|
174 | virtual void SetupForMaterial(const G4ParticleDefinition*, |
---|
175 | const G4Material*, |
---|
176 | G4double kineticEnergy); |
---|
177 | |
---|
178 | // add a region for the model |
---|
179 | virtual void DefineForRegion(const G4Region*); |
---|
180 | |
---|
181 | protected: |
---|
182 | |
---|
183 | // initialisation of the ParticleChange for the model |
---|
184 | G4ParticleChangeForLoss* GetParticleChangeForLoss(); |
---|
185 | |
---|
186 | // initialisation of the ParticleChange for the model |
---|
187 | G4ParticleChangeForGamma* GetParticleChangeForGamma(); |
---|
188 | |
---|
189 | // kinematically allowed max kinetic energy of a secondary |
---|
190 | virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition*, |
---|
191 | G4double kineticEnergy); |
---|
192 | |
---|
193 | public: |
---|
194 | |
---|
195 | //------------------------------------------------------------------------ |
---|
196 | // Generic methods common to all models |
---|
197 | //------------------------------------------------------------------------ |
---|
198 | |
---|
199 | // should be called at initialisation to build element selectors |
---|
200 | void InitialiseElementSelectors(const G4ParticleDefinition*, |
---|
201 | const G4DataVector&); |
---|
202 | |
---|
203 | // dEdx per unit length |
---|
204 | inline G4double ComputeDEDX(const G4MaterialCutsCouple*, |
---|
205 | const G4ParticleDefinition*, |
---|
206 | G4double kineticEnergy, |
---|
207 | G4double cutEnergy = DBL_MAX); |
---|
208 | |
---|
209 | // cross section per volume |
---|
210 | inline G4double CrossSection(const G4MaterialCutsCouple*, |
---|
211 | const G4ParticleDefinition*, |
---|
212 | G4double kineticEnergy, |
---|
213 | G4double cutEnergy = 0.0, |
---|
214 | G4double maxEnergy = DBL_MAX); |
---|
215 | |
---|
216 | // compute mean free path via cross section per volume |
---|
217 | inline G4double ComputeMeanFreePath(const G4ParticleDefinition*, |
---|
218 | G4double kineticEnergy, |
---|
219 | const G4Material*, |
---|
220 | G4double cutEnergy = 0.0, |
---|
221 | G4double maxEnergy = DBL_MAX); |
---|
222 | |
---|
223 | // generic cross section per element |
---|
224 | inline G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*, |
---|
225 | const G4Element*, |
---|
226 | G4double kinEnergy, |
---|
227 | G4double cutEnergy = 0.0, |
---|
228 | G4double maxEnergy = DBL_MAX); |
---|
229 | |
---|
230 | // atom can be selected effitiantly if element selectors are initialised |
---|
231 | inline const G4Element* SelectRandomAtom(const G4MaterialCutsCouple*, |
---|
232 | const G4ParticleDefinition*, |
---|
233 | G4double kineticEnergy, |
---|
234 | G4double cutEnergy = 0.0, |
---|
235 | G4double maxEnergy = DBL_MAX); |
---|
236 | |
---|
237 | // to select atom cross section per volume is recomputed for each element |
---|
238 | inline const G4Element* SelectRandomAtom(const G4Material*, |
---|
239 | const G4ParticleDefinition*, |
---|
240 | G4double kineticEnergy, |
---|
241 | G4double cutEnergy = 0.0, |
---|
242 | G4double maxEnergy = DBL_MAX); |
---|
243 | |
---|
244 | // select isotope in order to have precise mass of the nucleus |
---|
245 | inline G4int SelectIsotopeNumber(const G4Element*); |
---|
246 | |
---|
247 | //------------------------------------------------------------------------ |
---|
248 | // Get/Set methods |
---|
249 | //------------------------------------------------------------------------ |
---|
250 | |
---|
251 | inline G4VEmFluctuationModel* GetModelOfFluctuations(); |
---|
252 | |
---|
253 | inline G4double HighEnergyLimit() const; |
---|
254 | |
---|
255 | inline G4double LowEnergyLimit() const; |
---|
256 | |
---|
257 | inline G4double PolarAngleLimit() const; |
---|
258 | |
---|
259 | inline G4double SecondaryThreshold() const; |
---|
260 | |
---|
261 | inline G4bool LPMFlag() const; |
---|
262 | |
---|
263 | inline G4bool DeexcitationFlag() const; |
---|
264 | |
---|
265 | inline void SetHighEnergyLimit(G4double); |
---|
266 | |
---|
267 | inline void SetLowEnergyLimit(G4double); |
---|
268 | |
---|
269 | inline void SetActivationHighEnergyLimit(G4double); |
---|
270 | |
---|
271 | inline void SetActivationLowEnergyLimit(G4double); |
---|
272 | |
---|
273 | inline G4bool IsActive(G4double kinEnergy); |
---|
274 | |
---|
275 | inline void SetPolarAngleLimit(G4double); |
---|
276 | |
---|
277 | inline void SetSecondaryThreshold(G4double); |
---|
278 | |
---|
279 | inline void SetLPMFlag(G4bool val); |
---|
280 | |
---|
281 | inline void SetDeexcitationFlag(G4bool val); |
---|
282 | |
---|
283 | inline void ActivateNuclearStopping(G4bool); |
---|
284 | |
---|
285 | inline G4double MaxSecondaryKinEnergy(const G4DynamicParticle* dynParticle); |
---|
286 | |
---|
287 | inline const G4String& GetName() const; |
---|
288 | |
---|
289 | inline void SetParticleChange(G4VParticleChange*, G4VEmFluctuationModel*); |
---|
290 | |
---|
291 | inline void SetCurrentCouple(const G4MaterialCutsCouple*); |
---|
292 | |
---|
293 | protected: |
---|
294 | |
---|
295 | inline const G4MaterialCutsCouple* CurrentCouple() const; |
---|
296 | |
---|
297 | inline void SetCurrentElement(const G4Element*); |
---|
298 | |
---|
299 | inline const G4Element* GetCurrentElement() const; |
---|
300 | |
---|
301 | private: |
---|
302 | |
---|
303 | // hide assignment operator |
---|
304 | G4VEmModel & operator=(const G4VEmModel &right); |
---|
305 | G4VEmModel(const G4VEmModel&); |
---|
306 | |
---|
307 | // ======== Parameters of the class fixed at construction ========= |
---|
308 | |
---|
309 | G4VEmFluctuationModel* fluc; |
---|
310 | const G4String name; |
---|
311 | |
---|
312 | // ======== Parameters of the class fixed at initialisation ======= |
---|
313 | |
---|
314 | G4double lowLimit; |
---|
315 | G4double highLimit; |
---|
316 | G4double eMinActive; |
---|
317 | G4double eMaxActive; |
---|
318 | G4double polarAngleLimit; |
---|
319 | G4double secondaryThreshold; |
---|
320 | G4bool theLPMflag; |
---|
321 | |
---|
322 | G4int nSelectors; |
---|
323 | std::vector<G4EmElementSelector*> elmSelectors; |
---|
324 | |
---|
325 | protected: |
---|
326 | |
---|
327 | G4VParticleChange* pParticleChange; |
---|
328 | G4bool nuclearStopping; |
---|
329 | |
---|
330 | // ======== Cashed values - may be state dependent ================ |
---|
331 | |
---|
332 | private: |
---|
333 | |
---|
334 | const G4MaterialCutsCouple* currentCouple; |
---|
335 | const G4Element* currentElement; |
---|
336 | |
---|
337 | G4int nsec; |
---|
338 | G4bool flagDeexcitation; |
---|
339 | std::vector<G4double> xsec; |
---|
340 | |
---|
341 | }; |
---|
342 | |
---|
343 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
344 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
345 | |
---|
346 | inline G4double G4VEmModel::ComputeDEDX(const G4MaterialCutsCouple* c, |
---|
347 | const G4ParticleDefinition* p, |
---|
348 | G4double kinEnergy, |
---|
349 | G4double cutEnergy) |
---|
350 | { |
---|
351 | currentCouple = c; |
---|
352 | return ComputeDEDXPerVolume(c->GetMaterial(),p,kinEnergy,cutEnergy); |
---|
353 | } |
---|
354 | |
---|
355 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
356 | |
---|
357 | inline G4double G4VEmModel::CrossSection(const G4MaterialCutsCouple* c, |
---|
358 | const G4ParticleDefinition* p, |
---|
359 | G4double kinEnergy, |
---|
360 | G4double cutEnergy, |
---|
361 | G4double maxEnergy) |
---|
362 | { |
---|
363 | currentCouple = c; |
---|
364 | return CrossSectionPerVolume(c->GetMaterial(),p,kinEnergy,cutEnergy,maxEnergy); |
---|
365 | } |
---|
366 | |
---|
367 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
368 | |
---|
369 | inline G4double G4VEmModel::ComputeMeanFreePath(const G4ParticleDefinition* p, |
---|
370 | G4double ekin, |
---|
371 | const G4Material* material, |
---|
372 | G4double emin, |
---|
373 | G4double emax) |
---|
374 | { |
---|
375 | G4double mfp = DBL_MAX; |
---|
376 | G4double cross = CrossSectionPerVolume(material,p,ekin,emin,emax); |
---|
377 | if (cross > DBL_MIN) mfp = 1./cross; |
---|
378 | return mfp; |
---|
379 | } |
---|
380 | |
---|
381 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
382 | |
---|
383 | inline G4double G4VEmModel::ComputeCrossSectionPerAtom( |
---|
384 | const G4ParticleDefinition* part, |
---|
385 | const G4Element* elm, |
---|
386 | G4double kinEnergy, |
---|
387 | G4double cutEnergy, |
---|
388 | G4double maxEnergy) |
---|
389 | { |
---|
390 | currentElement = elm; |
---|
391 | return ComputeCrossSectionPerAtom(part,kinEnergy,elm->GetZ(),elm->GetN(), |
---|
392 | cutEnergy,maxEnergy); |
---|
393 | } |
---|
394 | |
---|
395 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
396 | |
---|
397 | inline |
---|
398 | const G4Element* G4VEmModel::SelectRandomAtom(const G4MaterialCutsCouple* couple, |
---|
399 | const G4ParticleDefinition* p, |
---|
400 | G4double kinEnergy, |
---|
401 | G4double cutEnergy, |
---|
402 | G4double maxEnergy) |
---|
403 | { |
---|
404 | currentCouple = couple; |
---|
405 | if(nSelectors > 0) { |
---|
406 | currentElement = |
---|
407 | elmSelectors[couple->GetIndex()]->SelectRandomAtom(kinEnergy); |
---|
408 | } else { |
---|
409 | currentElement = SelectRandomAtom(couple->GetMaterial(),p,kinEnergy, |
---|
410 | cutEnergy,maxEnergy); |
---|
411 | } |
---|
412 | return currentElement; |
---|
413 | } |
---|
414 | |
---|
415 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
416 | |
---|
417 | inline |
---|
418 | const G4Element* G4VEmModel::SelectRandomAtom(const G4Material* material, |
---|
419 | const G4ParticleDefinition* pd, |
---|
420 | G4double kinEnergy, |
---|
421 | G4double tcut, |
---|
422 | G4double tmax) |
---|
423 | { |
---|
424 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
425 | G4int n = material->GetNumberOfElements() - 1; |
---|
426 | currentElement = (*theElementVector)[n]; |
---|
427 | if (n > 0) { |
---|
428 | G4double x = G4UniformRand()* |
---|
429 | G4VEmModel::CrossSectionPerVolume(material,pd,kinEnergy,tcut,tmax); |
---|
430 | for(G4int i=0; i<n; i++) { |
---|
431 | if (x <= xsec[i]) { |
---|
432 | currentElement = (*theElementVector)[i]; |
---|
433 | break; |
---|
434 | } |
---|
435 | } |
---|
436 | } |
---|
437 | return currentElement; |
---|
438 | } |
---|
439 | |
---|
440 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
441 | |
---|
442 | inline G4int G4VEmModel::SelectIsotopeNumber(const G4Element* elm) |
---|
443 | { |
---|
444 | currentElement = elm; |
---|
445 | G4int N = G4int(elm->GetN() + 0.5); |
---|
446 | G4int ni = elm->GetNumberOfIsotopes(); |
---|
447 | if(ni > 0) { |
---|
448 | G4int idx = 0; |
---|
449 | if(ni > 1) { |
---|
450 | G4double* ab = elm->GetRelativeAbundanceVector(); |
---|
451 | G4double x = G4UniformRand(); |
---|
452 | for(; idx<ni; idx++) { |
---|
453 | x -= ab[idx]; |
---|
454 | if (x <= 0.0) break; |
---|
455 | } |
---|
456 | if(idx >= ni) idx = ni - 1; |
---|
457 | } |
---|
458 | N = elm->GetIsotope(idx)->GetN(); |
---|
459 | } |
---|
460 | return N; |
---|
461 | } |
---|
462 | |
---|
463 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
464 | |
---|
465 | inline G4VEmFluctuationModel* G4VEmModel::GetModelOfFluctuations() |
---|
466 | { |
---|
467 | return fluc; |
---|
468 | } |
---|
469 | |
---|
470 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
471 | |
---|
472 | inline G4double G4VEmModel::HighEnergyLimit() const |
---|
473 | { |
---|
474 | return highLimit; |
---|
475 | } |
---|
476 | |
---|
477 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
478 | |
---|
479 | inline G4double G4VEmModel::LowEnergyLimit() const |
---|
480 | { |
---|
481 | return lowLimit; |
---|
482 | } |
---|
483 | |
---|
484 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
485 | |
---|
486 | inline G4double G4VEmModel::PolarAngleLimit() const |
---|
487 | { |
---|
488 | return polarAngleLimit; |
---|
489 | } |
---|
490 | |
---|
491 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
492 | |
---|
493 | inline G4double G4VEmModel::SecondaryThreshold() const |
---|
494 | { |
---|
495 | return secondaryThreshold; |
---|
496 | } |
---|
497 | |
---|
498 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
499 | |
---|
500 | inline G4bool G4VEmModel::LPMFlag() const |
---|
501 | { |
---|
502 | return theLPMflag; |
---|
503 | } |
---|
504 | |
---|
505 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
506 | |
---|
507 | inline G4bool G4VEmModel::DeexcitationFlag() const |
---|
508 | { |
---|
509 | return flagDeexcitation; |
---|
510 | } |
---|
511 | |
---|
512 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
513 | |
---|
514 | inline void G4VEmModel::SetHighEnergyLimit(G4double val) |
---|
515 | { |
---|
516 | highLimit = val; |
---|
517 | } |
---|
518 | |
---|
519 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
520 | |
---|
521 | inline void G4VEmModel::SetLowEnergyLimit(G4double val) |
---|
522 | { |
---|
523 | lowLimit = val; |
---|
524 | } |
---|
525 | |
---|
526 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
527 | |
---|
528 | inline void G4VEmModel::SetActivationHighEnergyLimit(G4double val) |
---|
529 | { |
---|
530 | eMaxActive = val; |
---|
531 | } |
---|
532 | |
---|
533 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
534 | |
---|
535 | inline void G4VEmModel::SetActivationLowEnergyLimit(G4double val) |
---|
536 | { |
---|
537 | eMinActive = val; |
---|
538 | } |
---|
539 | |
---|
540 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
541 | |
---|
542 | inline G4bool G4VEmModel::IsActive(G4double kinEnergy) |
---|
543 | { |
---|
544 | return (kinEnergy >= eMinActive && kinEnergy <= eMaxActive); |
---|
545 | } |
---|
546 | |
---|
547 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
548 | |
---|
549 | inline void G4VEmModel::SetPolarAngleLimit(G4double val) |
---|
550 | { |
---|
551 | polarAngleLimit = val; |
---|
552 | } |
---|
553 | |
---|
554 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
555 | |
---|
556 | inline void G4VEmModel::SetSecondaryThreshold(G4double val) |
---|
557 | { |
---|
558 | secondaryThreshold = val; |
---|
559 | } |
---|
560 | |
---|
561 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
562 | |
---|
563 | inline void G4VEmModel::SetLPMFlag(G4bool val) |
---|
564 | { |
---|
565 | theLPMflag = val; |
---|
566 | } |
---|
567 | |
---|
568 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
569 | |
---|
570 | inline void G4VEmModel::SetDeexcitationFlag(G4bool val) |
---|
571 | { |
---|
572 | flagDeexcitation = val; |
---|
573 | } |
---|
574 | |
---|
575 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
576 | |
---|
577 | inline void G4VEmModel::ActivateNuclearStopping(G4bool val) |
---|
578 | { |
---|
579 | nuclearStopping = val; |
---|
580 | } |
---|
581 | |
---|
582 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
583 | |
---|
584 | inline |
---|
585 | G4double G4VEmModel::MaxSecondaryKinEnergy(const G4DynamicParticle* dynPart) |
---|
586 | { |
---|
587 | return MaxSecondaryEnergy(dynPart->GetDefinition(), |
---|
588 | dynPart->GetKineticEnergy()); |
---|
589 | } |
---|
590 | |
---|
591 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
592 | |
---|
593 | inline const G4String& G4VEmModel::GetName() const |
---|
594 | { |
---|
595 | return name; |
---|
596 | } |
---|
597 | |
---|
598 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
599 | |
---|
600 | inline void G4VEmModel::SetParticleChange(G4VParticleChange* p, |
---|
601 | G4VEmFluctuationModel* f = 0) |
---|
602 | { |
---|
603 | if(p && pParticleChange != p) pParticleChange = p; |
---|
604 | fluc = f; |
---|
605 | } |
---|
606 | |
---|
607 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
608 | |
---|
609 | inline void G4VEmModel::SetCurrentCouple(const G4MaterialCutsCouple* p) |
---|
610 | { |
---|
611 | currentCouple = p; |
---|
612 | } |
---|
613 | |
---|
614 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
615 | |
---|
616 | inline const G4MaterialCutsCouple* G4VEmModel::CurrentCouple() const |
---|
617 | { |
---|
618 | return currentCouple; |
---|
619 | } |
---|
620 | |
---|
621 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
622 | |
---|
623 | inline void G4VEmModel::SetCurrentElement(const G4Element* elm) |
---|
624 | { |
---|
625 | currentElement = elm; |
---|
626 | } |
---|
627 | |
---|
628 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
629 | |
---|
630 | inline const G4Element* G4VEmModel::GetCurrentElement() const |
---|
631 | { |
---|
632 | return currentElement; |
---|
633 | } |
---|
634 | |
---|
635 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
636 | |
---|
637 | #endif |
---|
638 | |
---|