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

Last change on this file since 961 was 961, checked in by garnier, 15 years ago

update processes

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