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

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 18.8 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//
[1340]26// $Id: G4VEmModel.hh,v 1.77 2010/10/14 16:27:35 vnivanch Exp $
27// GEANT4 tag $Name: emutils-V09-03-23 $
[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)
[1340]69// 13-10-10 Added G4VEmAngularDistribution (VI)
[819]70//
71// Class Description:
72//
73// Abstract interface to energy loss models
74
75// -------------------------------------------------------------------
76//
77
78#ifndef G4VEmModel_h
79#define G4VEmModel_h 1
80
81#include "globals.hh"
82#include "G4DynamicParticle.hh"
83#include "G4ParticleDefinition.hh"
84#include "G4MaterialCutsCouple.hh"
85#include "G4Material.hh"
86#include "G4Element.hh"
87#include "G4ElementVector.hh"
88#include "G4DataVector.hh"
89#include "G4VEmFluctuationModel.hh"
[1340]90#include "G4VEmAngularDistribution.hh"
[961]91#include "G4EmElementSelector.hh"
[819]92#include "Randomize.hh"
[961]93#include <vector>
[819]94
95class G4PhysicsTable;
96class G4Region;
97class G4VParticleChange;
[1055]98class G4ParticleChangeForLoss;
99class G4ParticleChangeForGamma;
[819]100class G4Track;
101
102class G4VEmModel
103{
104
105public:
106
107  G4VEmModel(const G4String& nam);
108
109  virtual ~G4VEmModel();
110
111  //------------------------------------------------------------------------
[961]112  // Virtual methods to be implemented for any concrete model
[819]113  //------------------------------------------------------------------------
114
[961]115  virtual void Initialise(const G4ParticleDefinition*, 
116                          const G4DataVector&) = 0;
[819]117
118  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
119                                 const G4MaterialCutsCouple*,
120                                 const G4DynamicParticle*,
121                                 G4double tmin = 0.0,
122                                 G4double tmax = DBL_MAX) = 0;
123
124  //------------------------------------------------------------------------
125  // Methods with standard implementation; may be overwritten if needed
126  //------------------------------------------------------------------------
127
[961]128  // main method to compute dEdx
[819]129  virtual G4double ComputeDEDXPerVolume(const G4Material*,
130                                        const G4ParticleDefinition*,
131                                        G4double kineticEnergy,
132                                        G4double cutEnergy = DBL_MAX);
133
[961]134  // main method to compute cross section per Volume
135  virtual G4double CrossSectionPerVolume(const G4Material*,
136                                         const G4ParticleDefinition*,
137                                         G4double kineticEnergy,
138                                         G4double cutEnergy = 0.0,
139                                         G4double maxEnergy = DBL_MAX);
[819]140
[1196]141  // main method to compute cross section per atom
[819]142  virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
143                                              G4double kinEnergy, 
144                                              G4double Z, 
[961]145                                              G4double A = 0., /* amu */
[819]146                                              G4double cutEnergy = 0.0,
147                                              G4double maxEnergy = DBL_MAX);
[961]148                                                                     
149  // min cut in kinetic energy allowed by the model
[1340]150  // obsolete method will be removed
[961]151  virtual G4double MinEnergyCut(const G4ParticleDefinition*,
152                                const G4MaterialCutsCouple*);
[819]153
[961]154  // Compute effective ion charge square
[1315]155  virtual G4double ChargeSquareRatio(const G4Track&);
156
157  // Compute effective ion charge square
[961]158  virtual G4double GetChargeSquareRatio(const G4ParticleDefinition*,
159                                        const G4Material*,
160                                        G4double kineticEnergy);
161
162  // Compute ion charge
163  virtual G4double GetParticleCharge(const G4ParticleDefinition*,
164                                     const G4Material*,
165                                     G4double kineticEnergy);
166
[1196]167  // add correction to energy loss and compute non-ionizing energy loss
[961]168  virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*,
169                                    const G4DynamicParticle*,
170                                    G4double& eloss,
171                                    G4double& niel,
172                                    G4double length);
173
[1055]174  // sample PIXE deexcitation
175  virtual void SampleDeexcitationAlongStep(const G4Material*,
176                                           const G4Track&,
177                                           G4double& eloss);
178
[1196]179  // initilisation at run time for a given material
[1055]180  virtual void SetupForMaterial(const G4ParticleDefinition*,
181                                const G4Material*,
182                                G4double kineticEnergy);
183
[1196]184  // add a region for the model
185  virtual void DefineForRegion(const G4Region*);
186
[819]187protected:
188
[1055]189  // initialisation of the ParticleChange for the model
190  G4ParticleChangeForLoss* GetParticleChangeForLoss();
191
192  // initialisation of the ParticleChange for the model
193  G4ParticleChangeForGamma* GetParticleChangeForGamma();
194
[961]195  // kinematically allowed max kinetic energy of a secondary
[819]196  virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
197                                      G4double kineticEnergy); 
198
199public:
200
201  //------------------------------------------------------------------------
202  // Generic methods common to all models
203  //------------------------------------------------------------------------
204
[961]205  // should be called at initialisation to build element selectors
206  void InitialiseElementSelectors(const G4ParticleDefinition*, 
207                                  const G4DataVector&);
[819]208
[1055]209  // dEdx per unit length
210  inline G4double ComputeDEDX(const G4MaterialCutsCouple*,
211                              const G4ParticleDefinition*,
212                              G4double kineticEnergy,
213                              G4double cutEnergy = DBL_MAX);
214
215  // cross section per volume
216  inline G4double CrossSection(const G4MaterialCutsCouple*,
[1196]217                               const G4ParticleDefinition*,
218                               G4double kineticEnergy,
219                               G4double cutEnergy = 0.0,
220                               G4double maxEnergy = DBL_MAX);
[1055]221
[961]222  // compute mean free path via cross section per volume
[1055]223  inline G4double ComputeMeanFreePath(const G4ParticleDefinition*,
224                                      G4double kineticEnergy,
225                                      const G4Material*,   
226                                      G4double cutEnergy = 0.0,
227                                      G4double maxEnergy = DBL_MAX);
[819]228
[961]229  // generic cross section per element
230  inline G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
231                                             const G4Element*,
232                                             G4double kinEnergy, 
233                                             G4double cutEnergy = 0.0,
234                                             G4double maxEnergy = DBL_MAX);
[819]235
[1315]236  // select isotope in order to have precise mass of the nucleus
237  inline G4int SelectIsotopeNumber(const G4Element*);
238
[961]239  // atom can be selected effitiantly if element selectors are initialised
240  inline const G4Element* SelectRandomAtom(const G4MaterialCutsCouple*,
241                                           const G4ParticleDefinition*,
242                                           G4double kineticEnergy,
243                                           G4double cutEnergy = 0.0,
244                                           G4double maxEnergy = DBL_MAX);
[819]245
[1055]246  // to select atom cross section per volume is recomputed for each element
[1315]247  const G4Element* SelectRandomAtom(const G4Material*,
248                                    const G4ParticleDefinition*,
249                                    G4double kineticEnergy,
250                                    G4double cutEnergy = 0.0,
251                                    G4double maxEnergy = DBL_MAX);
[819]252
[961]253  //------------------------------------------------------------------------
254  // Get/Set methods
255  //------------------------------------------------------------------------
[819]256
[1340]257  void SetParticleChange(G4VParticleChange*, G4VEmFluctuationModel* f=0);
258
[961]259  inline G4VEmFluctuationModel* GetModelOfFluctuations();
260
[1340]261  inline G4VEmAngularDistribution* GetAngularDistribution();
262
263  inline void SetAngularDistribution(G4VEmAngularDistribution*);
264
[961]265  inline G4double HighEnergyLimit() const;
266
267  inline G4double LowEnergyLimit() const;
268
269  inline G4double PolarAngleLimit() const;
270
271  inline G4double SecondaryThreshold() const;
272
273  inline G4bool LPMFlag() const;
274
[1055]275  inline G4bool DeexcitationFlag() const;
276
[961]277  inline void SetHighEnergyLimit(G4double);
278
279  inline void SetLowEnergyLimit(G4double);
280
[1196]281  inline void SetActivationHighEnergyLimit(G4double);
282
283  inline void SetActivationLowEnergyLimit(G4double);
284
285  inline G4bool IsActive(G4double kinEnergy);
286
[961]287  inline void SetPolarAngleLimit(G4double);
288
289  inline void SetSecondaryThreshold(G4double);
290
291  inline void SetLPMFlag(G4bool val);
292
[1055]293  inline void SetDeexcitationFlag(G4bool val);
294
[961]295  inline G4double MaxSecondaryKinEnergy(const G4DynamicParticle* dynParticle);
296
297  inline const G4String& GetName() const;
298
[1055]299  inline void SetCurrentCouple(const G4MaterialCutsCouple*);
300
[1315]301  inline const G4Element* GetCurrentElement() const;
302
[819]303protected:
304
[1055]305  inline const G4MaterialCutsCouple* CurrentCouple() const;
[819]306
[961]307  inline void SetCurrentElement(const G4Element*);
[819]308
309private:
310
311  //  hide assignment operator
312  G4VEmModel & operator=(const  G4VEmModel &right);
313  G4VEmModel(const  G4VEmModel&);
314
[961]315  // ======== Parameters of the class fixed at construction =========
316
[1340]317  G4VEmFluctuationModel* flucModel;
318  G4VEmAngularDistribution* anglModel;
[961]319  const G4String   name;
320
321  // ======== Parameters of the class fixed at initialisation =======
322
[819]323  G4double        lowLimit;
324  G4double        highLimit;
[1196]325  G4double        eMinActive;
326  G4double        eMaxActive;
[961]327  G4double        polarAngleLimit;
328  G4double        secondaryThreshold;
329  G4bool          theLPMflag;
[819]330
[961]331  G4int           nSelectors;
332  std::vector<G4EmElementSelector*> elmSelectors;
[819]333
334protected:
335
336  G4VParticleChange*  pParticleChange;
[1340]337  //  G4bool              nuclearStopping;
[961]338
339  // ======== Cashed values - may be state dependent ================
340
341private:
342
[1055]343  const G4MaterialCutsCouple* currentCouple;
344  const G4Element*            currentElement;
345
[961]346  G4int                  nsec;
[1055]347  G4bool                 flagDeexcitation;
[961]348  std::vector<G4double>  xsec;
349
[819]350};
351
[1315]352// ======== Run time inline methods ================
353
354inline void G4VEmModel::SetCurrentCouple(const G4MaterialCutsCouple* p)
355{
356  currentCouple = p;
357}
358
[819]359//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[1315]360
361inline const G4MaterialCutsCouple* G4VEmModel::CurrentCouple() const
362{
363  return currentCouple;
364}
365
[819]366//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
367
[1315]368inline void G4VEmModel::SetCurrentElement(const G4Element* elm)
369{
370  currentElement = elm;
371}
372
373//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
374
375inline const G4Element* G4VEmModel::GetCurrentElement() const
376{
377  return currentElement;
378}
379
380//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
381
382inline 
383G4double G4VEmModel::MaxSecondaryKinEnergy(const G4DynamicParticle* dynPart)
384{
[1340]385  return MaxSecondaryEnergy(dynPart->GetParticleDefinition(),
[1315]386                            dynPart->GetKineticEnergy());
387}
388
389//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
390
[1055]391inline G4double G4VEmModel::ComputeDEDX(const G4MaterialCutsCouple* c,
392                                        const G4ParticleDefinition* p,
393                                        G4double kinEnergy,
394                                        G4double cutEnergy)
[819]395{
[1055]396  currentCouple = c;
397  return ComputeDEDXPerVolume(c->GetMaterial(),p,kinEnergy,cutEnergy);
[819]398}
399
400//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
401
[1055]402inline G4double G4VEmModel::CrossSection(const G4MaterialCutsCouple* c,
403                                         const G4ParticleDefinition* p,
404                                         G4double kinEnergy,
405                                         G4double cutEnergy,
406                                         G4double maxEnergy)
[819]407{
[1055]408  currentCouple = c;
409  return CrossSectionPerVolume(c->GetMaterial(),p,kinEnergy,cutEnergy,maxEnergy);
[819]410}
411
412//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
413
[1055]414inline G4double G4VEmModel::ComputeMeanFreePath(const G4ParticleDefinition* p,
415                                                G4double ekin,
416                                                const G4Material* material,     
417                                                G4double emin,
418                                                G4double emax)
[819]419{
[1055]420  G4double mfp = DBL_MAX;
421  G4double cross = CrossSectionPerVolume(material,p,ekin,emin,emax);
[1315]422  if (cross > DBL_MIN) { mfp = 1./cross; }
[1055]423  return mfp;
[819]424}
425
426//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
427
[961]428inline G4double G4VEmModel::ComputeCrossSectionPerAtom(
429                const G4ParticleDefinition* part,
430                const G4Element* elm,
431                G4double kinEnergy, 
432                G4double cutEnergy,
433                G4double maxEnergy)
[819]434{
[961]435  currentElement = elm;
436  return ComputeCrossSectionPerAtom(part,kinEnergy,elm->GetZ(),elm->GetN(),
437                                    cutEnergy,maxEnergy);
[819]438}
439
440//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
441
[961]442inline 
443const G4Element* G4VEmModel::SelectRandomAtom(const G4MaterialCutsCouple* couple,
444                                              const G4ParticleDefinition* p,
445                                              G4double kinEnergy,
446                                              G4double cutEnergy,
447                                              G4double maxEnergy)
[819]448{
[1055]449  currentCouple = couple;
[961]450  if(nSelectors > 0) {
451    currentElement = 
452      elmSelectors[couple->GetIndex()]->SelectRandomAtom(kinEnergy);
453  } else {
454    currentElement = SelectRandomAtom(couple->GetMaterial(),p,kinEnergy,
455                                      cutEnergy,maxEnergy);
456  }
457  return currentElement;
[819]458}
459
460//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
461
[961]462inline G4int G4VEmModel::SelectIsotopeNumber(const G4Element* elm)
463{
[1196]464  currentElement = elm;
[961]465  G4int N = G4int(elm->GetN() + 0.5);
466  G4int ni = elm->GetNumberOfIsotopes();
467  if(ni > 0) {
468    G4int idx = 0;
469    if(ni > 1) {
[1196]470      G4double* ab = elm->GetRelativeAbundanceVector();
[961]471      G4double x = G4UniformRand();
[1315]472      for(; idx<ni; ++idx) {
[961]473        x -= ab[idx];
[1315]474        if (x <= 0.0) { break; }
[961]475      }
[1315]476      if(idx >= ni) { idx = ni - 1; }
[961]477    }
478    N = elm->GetIsotope(idx)->GetN();
479  }
480  return N;
481}
482
[1315]483// ======== Get/Set inline methods used at initialisation ================
[961]484
[1055]485inline G4VEmFluctuationModel* G4VEmModel::GetModelOfFluctuations()
[819]486{
[1340]487  return flucModel;
[819]488}
489
490//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
491
[1340]492inline G4VEmAngularDistribution* G4VEmModel::GetAngularDistribution()
493{
494  return anglModel;
495}
496
497//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
498
499inline void G4VEmModel::SetAngularDistribution(G4VEmAngularDistribution* p)
500{
501  anglModel = p;
502}
503
504//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
505
[1055]506inline G4double G4VEmModel::HighEnergyLimit() const
[819]507{
[1055]508  return highLimit;
[819]509}
510
511//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
512
[1055]513inline G4double G4VEmModel::LowEnergyLimit() const
[819]514{
[1055]515  return lowLimit;
[819]516}
517
518//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
519
[1055]520inline G4double G4VEmModel::PolarAngleLimit() const
[819]521{
[1055]522  return polarAngleLimit;
[819]523}
524
525//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
526
[1055]527inline G4double G4VEmModel::SecondaryThreshold() const
[819]528{
[1055]529  return secondaryThreshold;
[819]530}
531
[1055]532//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
533
534inline G4bool G4VEmModel::LPMFlag() const 
535{
536  return theLPMflag;
537}
538
539//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
540
541inline G4bool G4VEmModel::DeexcitationFlag() const 
542{
543  return flagDeexcitation;
544}
545
[819]546//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[1055]547
548inline void G4VEmModel::SetHighEnergyLimit(G4double val)
549{
550  highLimit = val;
551}
552
[819]553//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
554
[1055]555inline void G4VEmModel::SetLowEnergyLimit(G4double val)
556{
557  lowLimit = val;
558}
[819]559
560//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
561
[1196]562inline void G4VEmModel::SetActivationHighEnergyLimit(G4double val)
563{
564  eMaxActive = val;
565}
566
567//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
568
569inline void G4VEmModel::SetActivationLowEnergyLimit(G4double val)
570{
571  eMinActive = val;
572}
573
574//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
575
576inline G4bool G4VEmModel::IsActive(G4double kinEnergy)
577{
578  return (kinEnergy >= eMinActive && kinEnergy <= eMaxActive);
579}
580
581//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
582
[1055]583inline void G4VEmModel::SetPolarAngleLimit(G4double val)
[819]584{
[1055]585  polarAngleLimit = val;
[819]586}
587
588//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
589
[1055]590inline void G4VEmModel::SetSecondaryThreshold(G4double val) 
[819]591{
[1055]592  secondaryThreshold = val;
[819]593}
594
[1055]595//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
596
597inline void G4VEmModel::SetLPMFlag(G4bool val) 
598{
599  theLPMflag = val;
600}
601
602//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
603
604inline void G4VEmModel::SetDeexcitationFlag(G4bool val) 
605{
606  flagDeexcitation = val;
607}
608
[819]609//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
610
[1055]611inline const G4String& G4VEmModel::GetName() const 
612{
613  return name;
614}
[819]615
616//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
617
618#endif
619
Note: See TracBrowser for help on using the repository browser.