source: trunk/source/processes/electromagnetic/utils/include/G4VEmModel.hh@ 1057

Last change on this file since 1057 was 1055, checked in by garnier, 17 years ago

maj sur la beta de geant 4.9.3

File size: 18.6 KB
RevLine 
[819]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//
[1055]26// $Id: G4VEmModel.hh,v 1.69 2009/05/26 15:00:49 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-beta-cand-01 $
[819]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)
[961]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)
[1055]67// 16-02-09 Moved implementations of virtual methods to source (VI)
68// 07-04-09 Moved msc methods from G4VEmModel to G4VMscModel (VI)
[819]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"
[961]89#include "G4EmElementSelector.hh"
[819]90#include "Randomize.hh"
[961]91#include <vector>
[819]92
93class G4PhysicsTable;
94class G4Region;
95class G4VParticleChange;
[1055]96class G4ParticleChangeForLoss;
97class G4ParticleChangeForGamma;
[819]98class G4Track;
99
100class G4VEmModel
101{
102
103public:
104
105 G4VEmModel(const G4String& nam);
106
107 virtual ~G4VEmModel();
108
109 //------------------------------------------------------------------------
[961]110 // Virtual methods to be implemented for any concrete model
[819]111 //------------------------------------------------------------------------
112
[961]113 virtual void Initialise(const G4ParticleDefinition*,
114 const G4DataVector&) = 0;
[819]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
[961]126 // main method to compute dEdx
[819]127 virtual G4double ComputeDEDXPerVolume(const G4Material*,
128 const G4ParticleDefinition*,
129 G4double kineticEnergy,
130 G4double cutEnergy = DBL_MAX);
131
[961]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);
[819]138
[961]139 // main method to compute cross section depending on atom
[819]140 virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
141 G4double kinEnergy,
142 G4double Z,
[961]143 G4double A = 0., /* amu */
[819]144 G4double cutEnergy = 0.0,
145 G4double maxEnergy = DBL_MAX);
[961]146
147 // min cut in kinetic energy allowed by the model
148 virtual G4double MinEnergyCut(const G4ParticleDefinition*,
149 const G4MaterialCutsCouple*);
[819]150
[961]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 ompute non-ionizing energy loss
162 virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*,
163 const G4DynamicParticle*,
164 G4double& eloss,
165 G4double& niel,
166 G4double length);
167
[1055]168 // sample PIXE deexcitation
169 virtual void SampleDeexcitationAlongStep(const G4Material*,
170 const G4Track&,
171 G4double& eloss);
172
173 // add region for the model
174 virtual void DefineForRegion(const G4Region*);
175
176 // initilisation at run time for given material
177 virtual void SetupForMaterial(const G4ParticleDefinition*,
178 const G4Material*,
179 G4double kineticEnergy);
180
[819]181protected:
182
[1055]183 // initialisation of the ParticleChange for the model
184 G4ParticleChangeForLoss* GetParticleChangeForLoss();
185
186 // initialisation of the ParticleChange for the model
187 G4ParticleChangeForGamma* GetParticleChangeForGamma();
188
[961]189 // kinematically allowed max kinetic energy of a secondary
[819]190 virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
191 G4double kineticEnergy);
192
193public:
194
195 //------------------------------------------------------------------------
196 // Generic methods common to all models
197 //------------------------------------------------------------------------
198
[961]199 // should be called at initialisation to build element selectors
200 void InitialiseElementSelectors(const G4ParticleDefinition*,
201 const G4DataVector&);
[819]202
[1055]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
[961]216 // compute mean free path via cross section per volume
[1055]217 inline G4double ComputeMeanFreePath(const G4ParticleDefinition*,
218 G4double kineticEnergy,
219 const G4Material*,
220 G4double cutEnergy = 0.0,
221 G4double maxEnergy = DBL_MAX);
[819]222
[961]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);
[819]229
[961]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);
[819]236
[1055]237 // to select atom cross section per volume is recomputed for each element
[961]238 inline const G4Element* SelectRandomAtom(const G4Material*,
239 const G4ParticleDefinition*,
240 G4double kineticEnergy,
241 G4double cutEnergy = 0.0,
242 G4double maxEnergy = DBL_MAX);
[819]243
[961]244 // select isotope in order to have precise mass of the nucleus
245 inline G4int SelectIsotopeNumber(const G4Element*);
[819]246
[961]247 //------------------------------------------------------------------------
248 // Get/Set methods
249 //------------------------------------------------------------------------
[819]250
[961]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
[1055]263 inline G4bool DeexcitationFlag() const;
264
[961]265 inline void SetHighEnergyLimit(G4double);
266
267 inline void SetLowEnergyLimit(G4double);
268
269 inline void SetPolarAngleLimit(G4double);
270
271 inline void SetSecondaryThreshold(G4double);
272
273 inline void SetLPMFlag(G4bool val);
274
[1055]275 inline void SetDeexcitationFlag(G4bool val);
276
[961]277 inline void ActivateNuclearStopping(G4bool);
278
279 inline G4double MaxSecondaryKinEnergy(const G4DynamicParticle* dynParticle);
280
281 inline const G4String& GetName() const;
282
283 inline void SetParticleChange(G4VParticleChange*, G4VEmFluctuationModel*);
284
[1055]285 inline void SetCurrentCouple(const G4MaterialCutsCouple*);
286
[819]287protected:
288
[1055]289 inline const G4MaterialCutsCouple* CurrentCouple() const;
[819]290
[961]291 inline void SetCurrentElement(const G4Element*);
[819]292
[1055]293 inline const G4Element* GetCurrentElement() const;
294
[819]295private:
296
297 // hide assignment operator
298 G4VEmModel & operator=(const G4VEmModel &right);
299 G4VEmModel(const G4VEmModel&);
300
[961]301 // ======== Parameters of the class fixed at construction =========
302
303 G4VEmFluctuationModel* fluc;
304 const G4String name;
305
306 // ======== Parameters of the class fixed at initialisation =======
307
[819]308 G4double lowLimit;
309 G4double highLimit;
[961]310 G4double polarAngleLimit;
311 G4double secondaryThreshold;
312 G4bool theLPMflag;
[819]313
[961]314 G4int nSelectors;
315 std::vector<G4EmElementSelector*> elmSelectors;
[819]316
317protected:
318
319 G4VParticleChange* pParticleChange;
[961]320 G4bool nuclearStopping;
321
322 // ======== Cashed values - may be state dependent ================
323
324private:
325
[1055]326 const G4MaterialCutsCouple* currentCouple;
327 const G4Element* currentElement;
328
[961]329 G4int nsec;
[1055]330 G4bool flagDeexcitation;
[961]331 std::vector<G4double> xsec;
332
[819]333};
334
335//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
336//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
337
[1055]338inline G4double G4VEmModel::ComputeDEDX(const G4MaterialCutsCouple* c,
339 const G4ParticleDefinition* p,
340 G4double kinEnergy,
341 G4double cutEnergy)
[819]342{
[1055]343 currentCouple = c;
344 return ComputeDEDXPerVolume(c->GetMaterial(),p,kinEnergy,cutEnergy);
[819]345}
346
347//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
348
[1055]349inline G4double G4VEmModel::CrossSection(const G4MaterialCutsCouple* c,
350 const G4ParticleDefinition* p,
351 G4double kinEnergy,
352 G4double cutEnergy,
353 G4double maxEnergy)
[819]354{
[1055]355 currentCouple = c;
356 return CrossSectionPerVolume(c->GetMaterial(),p,kinEnergy,cutEnergy,maxEnergy);
[819]357}
358
359//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
360
[1055]361inline G4double G4VEmModel::ComputeMeanFreePath(const G4ParticleDefinition* p,
362 G4double ekin,
363 const G4Material* material,
364 G4double emin,
365 G4double emax)
[819]366{
[1055]367 G4double mfp = DBL_MAX;
368 G4double cross = CrossSectionPerVolume(material,p,ekin,emin,emax);
369 if (cross > DBL_MIN) mfp = 1./cross;
370 return mfp;
[819]371}
372
373//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
374
[961]375inline G4double G4VEmModel::ComputeCrossSectionPerAtom(
376 const G4ParticleDefinition* part,
377 const G4Element* elm,
378 G4double kinEnergy,
379 G4double cutEnergy,
380 G4double maxEnergy)
[819]381{
[961]382 currentElement = elm;
383 return ComputeCrossSectionPerAtom(part,kinEnergy,elm->GetZ(),elm->GetN(),
384 cutEnergy,maxEnergy);
[819]385}
386
387//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
388
[961]389inline
390const G4Element* G4VEmModel::SelectRandomAtom(const G4MaterialCutsCouple* couple,
391 const G4ParticleDefinition* p,
392 G4double kinEnergy,
393 G4double cutEnergy,
394 G4double maxEnergy)
[819]395{
[1055]396 currentCouple = couple;
[961]397 if(nSelectors > 0) {
398 currentElement =
399 elmSelectors[couple->GetIndex()]->SelectRandomAtom(kinEnergy);
400 } else {
401 currentElement = SelectRandomAtom(couple->GetMaterial(),p,kinEnergy,
402 cutEnergy,maxEnergy);
403 }
404 return currentElement;
[819]405}
406
407//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
408
[961]409inline
410const G4Element* G4VEmModel::SelectRandomAtom(const G4Material* material,
411 const G4ParticleDefinition* pd,
412 G4double kinEnergy,
413 G4double tcut,
414 G4double tmax)
415{
416 const G4ElementVector* theElementVector = material->GetElementVector();
417 G4int n = material->GetNumberOfElements() - 1;
418 currentElement = (*theElementVector)[n];
419 if (n > 0) {
420 G4double x = G4UniformRand()*
421 G4VEmModel::CrossSectionPerVolume(material,pd,kinEnergy,tcut,tmax);
422 for(G4int i=0; i<n; i++) {
423 if (x <= xsec[i]) {
424 currentElement = (*theElementVector)[i];
425 break;
426 }
427 }
428 }
429 return currentElement;
430}
[819]431
[961]432//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
433
434inline G4int G4VEmModel::SelectIsotopeNumber(const G4Element* elm)
435{
436 G4int N = G4int(elm->GetN() + 0.5);
437 G4int ni = elm->GetNumberOfIsotopes();
438 if(ni > 0) {
439 G4int idx = 0;
440 if(ni > 1) {
441 G4double* ab = currentElement->GetRelativeAbundanceVector();
442 G4double x = G4UniformRand();
443 for(; idx<ni; idx++) {
444 x -= ab[idx];
445 if (x <= 0.0) break;
446 }
447 if(idx >= ni) idx = ni - 1;
448 }
449 N = elm->GetIsotope(idx)->GetN();
450 }
451 return N;
452}
453
454//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
455
[1055]456inline G4VEmFluctuationModel* G4VEmModel::GetModelOfFluctuations()
[819]457{
[1055]458 return fluc;
[819]459}
460
461//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
462
[1055]463inline G4double G4VEmModel::HighEnergyLimit() const
[819]464{
[1055]465 return highLimit;
[819]466}
467
468//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
469
[1055]470inline G4double G4VEmModel::LowEnergyLimit() const
[819]471{
[1055]472 return lowLimit;
[819]473}
474
475//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
476
[1055]477inline G4double G4VEmModel::PolarAngleLimit() const
[819]478{
[1055]479 return polarAngleLimit;
[819]480}
481
482//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
483
[1055]484inline G4double G4VEmModel::SecondaryThreshold() const
[819]485{
[1055]486 return secondaryThreshold;
[819]487}
488
[1055]489//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
490
491inline G4bool G4VEmModel::LPMFlag() const
492{
493 return theLPMflag;
494}
495
496//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
497
498inline G4bool G4VEmModel::DeexcitationFlag() const
499{
500 return flagDeexcitation;
501}
502
[819]503//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[1055]504
505inline void G4VEmModel::SetHighEnergyLimit(G4double val)
506{
507 highLimit = val;
508}
509
[819]510//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
511
[1055]512inline void G4VEmModel::SetLowEnergyLimit(G4double val)
513{
514 lowLimit = val;
515}
[819]516
517//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
518
[1055]519inline void G4VEmModel::SetPolarAngleLimit(G4double val)
[819]520{
[1055]521 polarAngleLimit = val;
[819]522}
523
524//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
525
[1055]526inline void G4VEmModel::SetSecondaryThreshold(G4double val)
[819]527{
[1055]528 secondaryThreshold = val;
[819]529}
530
[1055]531//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
532
533inline void G4VEmModel::SetLPMFlag(G4bool val)
534{
535 theLPMflag = val;
536}
537
538//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
539
540inline void G4VEmModel::SetDeexcitationFlag(G4bool val)
541{
542 flagDeexcitation = val;
543}
544
545//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
546
547inline void G4VEmModel::ActivateNuclearStopping(G4bool val)
548{
549 nuclearStopping = val;
550}
551
[819]552//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
553
[1055]554inline
555G4double G4VEmModel::MaxSecondaryKinEnergy(const G4DynamicParticle* dynPart)
[819]556{
[1055]557 return MaxSecondaryEnergy(dynPart->GetDefinition(),
558 dynPart->GetKineticEnergy());
[819]559}
560
561//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
562
[1055]563inline const G4String& G4VEmModel::GetName() const
564{
565 return name;
566}
[819]567
568//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
569
[1055]570inline void G4VEmModel::SetParticleChange(G4VParticleChange* p,
571 G4VEmFluctuationModel* f = 0)
572{
573 if(p && pParticleChange != p) pParticleChange = p;
574 fluc = f;
575}
[819]576
577//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
578
[1055]579inline void G4VEmModel::SetCurrentCouple(const G4MaterialCutsCouple* p)
580{
581 currentCouple = p;
582}
583
584//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
585
586inline const G4MaterialCutsCouple* G4VEmModel::CurrentCouple() const
587{
588 return currentCouple;
589}
590
591//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
592
593inline void G4VEmModel::SetCurrentElement(const G4Element* elm)
594{
595 currentElement = elm;
596}
597
598//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
599
600inline const G4Element* G4VEmModel::GetCurrentElement() const
601{
602 return currentElement;
603}
604
605//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
606
[819]607#endif
608
Note: See TracBrowser for help on using the repository browser.