Geant4  10.02
G4BoldyshevTripletModel.hh
Go to the documentation of this file.
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: G.Depaola & F.Longo
27 
28 #ifndef G4BoldyshevTripletModel_h
29 #define G4BoldyshevTripletModel_h 1
30 
31 #include "G4VEmModel.hh"
32 #include "G4Electron.hh"
33 #include "G4Positron.hh"
35 #include "G4LPhysicsFreeVector.hh"
36 #include "G4ProductionCutsTable.hh"
37 
39 {
40 
41 public:
42 
43  G4BoldyshevTripletModel(const G4ParticleDefinition* p = 0, const G4String& nam = "BoldyshevTripletConversion");
44 
45  virtual ~G4BoldyshevTripletModel();
46 
47  virtual void Initialise(const G4ParticleDefinition*,
48  const G4DataVector&);
49 
50  //MT
51  virtual void InitialiseLocal(const G4ParticleDefinition*,
52  G4VEmModel* masterModel);
53 
54  virtual void InitialiseForElement(const G4ParticleDefinition*, G4int Z);
55  //END MT
56 
58  const G4ParticleDefinition*,
59  G4double kinEnergy,
60  G4double Z,
61  G4double A=0,
62  G4double cut=0,
64 
65  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
66  const G4MaterialCutsCouple*,
67  const G4DynamicParticle*,
68  G4double tmin,
69  G4double maxEnergy);
70 
71  virtual G4double MinPrimaryEnergy(const G4Material*,
72  const G4ParticleDefinition*,
73  G4double);
74 
75 private:
76 
77  void ReadData(size_t Z, const char* path = 0);
78 
79  G4double ScreenFunction1(G4double screenVariable);
80  G4double ScreenFunction2(G4double screenVariable);
81 
84 
87 
90 
91  //MT
92  static G4int maxZ;
93  static G4LPhysicsFreeVector* data[100]; // 100 because Z range is 1-99
94  // in LivermoreRayleighModel, 101
95  // because Z range is 1-100
96  //END MT
97 
99 
100  G4double asinh (G4double value);
101 
102 
103 };
104 
105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
106 
108 {
109  G4double out;
110 
111  if (value>0)
112  out = std::log(value+std::sqrt(value*value+1));
113  else
114  out = -std::log(-value+std::sqrt(value*value+1));
115 
116  return out;
117 }
118 
119 
120 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
121 
122 #endif
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4BoldyshevTripletModel(const G4ParticleDefinition *p=0, const G4String &nam="BoldyshevTripletConversion")
G4double ScreenFunction1(G4double screenVariable)
G4double asinh(G4double value)
void ReadData(size_t Z, const char *path=0)
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
int G4int
Definition: G4Types.hh:78
G4ParticleChangeForGamma * fParticleChange
G4BoldyshevTripletModel & operator=(const G4BoldyshevTripletModel &right)
double A(double temperature)
bool G4bool
Definition: G4Types.hh:79
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
static G4LPhysicsFreeVector * data[100]
static const G4double emax
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double)
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
G4double ScreenFunction2(G4double screenVariable)