Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNAUeharaScreenedRutherfordElasticModel.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 #ifndef G4DNAUeharaScreenedRutherfordElasticModel_h
27 #define G4DNAUeharaScreenedRutherfordElasticModel_h 1
28 
30 
31 #include "G4VEmModel.hh"
33 #include "G4ProductionCutsTable.hh"
34 #include "G4NistManager.hh"
35 
37 {
38 public:
40  const G4String& nam = "DNAUeharaScreenedRutherfordElasticModel");
41 
43 
44  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
45 
46  virtual G4double CrossSectionPerVolume(const G4Material* material,
47  const G4ParticleDefinition* p,
48  G4double ekin,
49  G4double emin,
50  G4double emax);
51 
52  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
53  const G4MaterialCutsCouple*,
54  const G4DynamicParticle*,
55  G4double tmin,
56  G4double maxEnergy);
57 
58  inline void SelectFasterComputation(G4bool input);
59 
60  //---
61  // kept for backward compatibility
62  inline void SetKillBelowThreshold (G4double threshold);
64  inline void SelectHighEnergyLimit(G4double threshold);
65  //---
66 
67 private:
68  G4double intermediateEnergyLimit;
69 
70  // -- Brenner & Zaider
71  std::vector<G4double> betaCoeff;
72  std::vector<G4double> deltaCoeff;
73  std::vector<G4double> gamma035_10Coeff;
74  std::vector<G4double> gamma10_100Coeff;
75  std::vector<G4double> gamma100_200Coeff;
76 
77  // -- Water density table
78  const std::vector<G4double>* fpWaterDensity;
79 
80 protected:
82 
83 private:
84  G4int verboseLevel;
85  G4bool fasterCode;
86  G4bool isInitialised;
87 
88  // -- Cross section
89  G4double RutherfordCrossSection(G4double energy, G4double z);
90  G4double ScreeningFactor(G4double energy, G4double z);
91 
92  // -- Final state according to Brenner & Zaider
93  G4double BrennerZaiderRandomizeCosTheta(G4double k);
94  G4double CalculatePolynomial(G4double k,
95  std::vector<G4double>& vec);
96 
97  // -- Final state according to Screened Rutherford
98  G4double ScreenedRutherfordRandomizeCosTheta(G4double k,
99  G4double z);
100 
101  //
103  operator=(const G4DNAUeharaScreenedRutherfordElasticModel &right);
106 };
107 
108 
111 {
112  fasterCode = input;
113 }
114 
115 //---
116 // kept for backward compatibility
117 
118 inline void
120  G4double threshold)
121 {
122  if(threshold > 10. * CLHEP::keV)
123  {
124  G4Exception (
125  "*** WARNING : the G4DNAUeharaScreenedRutherfordElasticModel class is "
126  "used above 10 keV !",
127  "", JustWarning, "");
128  }
129 
130  SetHighEnergyLimit(threshold);
131 }
132 
133 inline void
135 {
136  G4ExceptionDescription errMsg;
137  errMsg << "*** WARNING : "
138  << "G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold"
139  << "is deprecated, the kill threshold won't be taken into account";
140 
141  G4Exception (
142  "G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold",
143  "DEPRECATED", JustWarning, errMsg);
144 }
145 
146 inline G4double
148 {
149  G4ExceptionDescription errMsg;
150  errMsg << "*** WARNING : "
151  << "G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold"
152  << "is deprecated, the returned value is nonsense";
153 
154  G4Exception (
155  "G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold",
156  "DEPRECATED", JustWarning, errMsg);
157 
158  return -1;
159 }
160 //---
161 
162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
163 
164 #endif
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static constexpr double keV
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
const char * p
Definition: xmltok.h:285
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
int G4int
Definition: G4Types.hh:78
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
G4DNAUeharaScreenedRutherfordElasticModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAUeharaScreenedRutherfordElasticModel")
bool G4bool
Definition: G4Types.hh:79
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const G4double emax
G4double energy(const ThreeVector &p, const G4double m)
double G4double
Definition: G4Types.hh:76