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

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

update to geant4.9.2

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