Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4LivermoreGammaConversionModel.hh
이 파일의 문서화 페이지로 가기
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 // Author: Sebastien Incerti
27 // 22 January 2012
28 // on base of G4LivermoreGammaConversionModel (original version)
29 // and G4LivermoreRayleighModel (MT version)
30 
31 #ifndef G4LivermoreGammaConversionModel_h
32 #define G4LivermoreGammaConversionModel_h 1
33 
34 #include "G4VEmModel.hh"
35 #include "G4Log.hh"
36 
39 class G4PhysicsLogVector;
40 
42 {
43 
44 public:
45 
47  const G4ParticleDefinition* p = nullptr,
48  const G4String& nam = "LivermoreConversion");
49 
51 
52  virtual void Initialise(const G4ParticleDefinition*,
53  const G4DataVector&);
54 
55  virtual void InitialiseLocal(const G4ParticleDefinition*,
56  G4VEmModel* masterModel);
57 
58  virtual void InitialiseForElement(const G4ParticleDefinition*, G4int Z);
59 
61  const G4ParticleDefinition*,
62  G4double kinEnergy,
63  G4double Z,
64  G4double A=0.0,
65  G4double cut=0.0,
67 
68  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
69  const G4MaterialCutsCouple*,
70  const G4DynamicParticle*,
71  G4double tmin,
72  G4double maxEnergy);
73 
74  virtual G4double MinPrimaryEnergy(const G4Material*,
75  const G4ParticleDefinition*,
76  G4double);
77 
78 private:
79 
80  void ReadData(size_t Z, const char* path = nullptr);
82 
83  inline G4double ScreenFunction1(G4double screenVariable);
84  inline G4double ScreenFunction2(G4double screenVariable);
85 
87  (const G4LivermoreGammaConversionModel &right) = delete;
89 
93 
96  static G4int maxZ;
97 
98  static G4LPhysicsFreeVector* data[100]; // 100 because Z range is 1-99
99  static G4PhysicsLogVector* probTriplet[100]; //
100 
102 };
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105 
106 inline G4double
108 {
109  // Compute the value of the screening function 3*phi1 - phi2
110  return (screenVariable > 1.)
111  ? 42.24 - 8.368 * G4Log(screenVariable + 0.952)
112  : 42.392 - screenVariable * (7.796 - 1.961 * screenVariable);
113 }
114 
115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
116 
117 inline G4double
119 {
120  // Compute the value of the screening function 1.5*phi1 - 0.5*phi2
121  return (screenVariable > 1.)
122  ? 42.24 - 8.368 * G4Log(screenVariable + 0.952)
123  : 41.405 - screenVariable * (5.828 - 0.8945 * screenVariable);
124 }
125 
126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
127 
128 #endif
void ReadData(size_t Z, const char *path=nullptr)
G4double ScreenFunction2(G4double screenVariable)
const char * p
Definition: xmltok.h:285
static const G4double emax
G4double G4Log(G4double x)
Definition: G4Log.hh:230
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0.0, G4double cut=0.0, G4double emax=DBL_MAX)
const XML_Char const XML_Char * data
Definition: expat.h:268
Float_t Z
void InitialiseProbability(const G4ParticleDefinition *, G4int Z)
double G4double
Definition: G4Types.hh:76
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double)
double A(double temperature)
int G4int
Definition: G4Types.hh:78
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4LivermoreGammaConversionModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="LivermoreConversion")
G4double ScreenFunction1(G4double screenVariable)
#define DBL_MAX
Definition: templates.hh:83
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
static G4PhysicsLogVector * probTriplet[100]