Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/dna/models/include/G4DNAUeharaScreenedRutherfordElasticModel.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 
 29 #include <CLHEP/Units/SystemOfUnits.h>
 30 
 31 #include "G4VEmModel.hh"
 32 #include "G4ParticleChangeForGamma.hh"
 33 #include "G4ProductionCutsTable.hh"
 34 #include "G4NistManager.hh"
 35 
 36 class G4DNAUeharaScreenedRutherfordElasticModel : public G4VEmModel
 37 {
 38 public:
 39   G4DNAUeharaScreenedRutherfordElasticModel(const G4ParticleDefinition* p = nullptr, 
 40               const G4String& nam = "DNAUeharaScreenedRutherfordElasticModel");
 41 
 42   ~G4DNAUeharaScreenedRutherfordElasticModel() override = default;
 43 
 44   void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
 45 
 46   G4double CrossSectionPerVolume(const G4Material* material,
 47          const G4ParticleDefinition* p,
 48          G4double ekin,
 49          G4double emin,
 50          G4double emax) override;
 51 
 52   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
 53        const G4MaterialCutsCouple*,
 54        const G4DynamicParticle*,
 55        G4double tmin,
 56        G4double maxEnergy) override;
 57 
 58   inline void SelectFasterComputation(G4bool input);
 59   
 60   //---
 61   // kept for backward compatibility
 62   inline void SetKillBelowThreshold (G4double threshold);
 63   inline G4double GetKillBelowThreshold();
 64   inline void SelectHighEnergyLimit(G4double threshold);
 65   //---
 66 
 67   //
 68   G4DNAUeharaScreenedRutherfordElasticModel&
 69     operator=(const G4DNAUeharaScreenedRutherfordElasticModel &right) = delete;
 70   G4DNAUeharaScreenedRutherfordElasticModel(
 71       const G4DNAUeharaScreenedRutherfordElasticModel&) = delete;
 72 
 73 private:
 74 
 75   // -- Cross section
 76   G4double RutherfordCrossSection(G4double energy, G4double z);
 77   G4double ScreeningFactor(G4double energy, G4double z);
 78   
 79   // -- Final state according to Brenner & Zaider
 80   G4double BrennerZaiderRandomizeCosTheta(G4double k);
 81   G4double CalculatePolynomial(G4double k,
 82                                std::vector<G4double>& vec);
 83    
 84   // -- Final state according to Screened Rutherford
 85   G4double ScreenedRutherfordRandomizeCosTheta(G4double k,
 86                                                G4double z);
 87 
 88 protected:
 89   G4ParticleChangeForGamma* fParticleChangeForGamma = nullptr;
 90 
 91 private:
 92   G4double iLowEnergyLimit;
 93   G4double intermediateEnergyLimit;
 94   G4double iHighEnergyLimit;
 95   
 96   // -- Brenner & Zaider
 97   std::vector<G4double> betaCoeff;
 98   std::vector<G4double> deltaCoeff;
 99   std::vector<G4double> gamma035_10Coeff;
100   std::vector<G4double> gamma10_100Coeff;
101   std::vector<G4double> gamma100_200Coeff;
102 
103   // -- Water density table
104   const std::vector<G4double>* fpWaterDensity = nullptr;
105   
106   G4int verboseLevel;
107   G4bool isInitialised = false;
108   // Selection of computation method
109   // We do not recommend "true" usage with the current cumul. proba. settings
110   G4bool fasterCode = false;
111 };
112  
113 
114 inline void G4DNAUeharaScreenedRutherfordElasticModel::
115 SelectFasterComputation(G4bool input)
116 { 
117   fasterCode = input; 
118 }
119 
120 //---
121 // kept for backward compatibility
122 
123 inline void
124 G4DNAUeharaScreenedRutherfordElasticModel::SelectHighEnergyLimit(
125     G4double threshold)
126 {
127   if(threshold > 10. * CLHEP::keV)
128   {
129     G4Exception (
130         "*** WARNING : the G4DNAUeharaScreenedRutherfordElasticModel class is "
131         "used above 10 keV !",
132         "", JustWarning, "");
133   }
134 
135   SetHighEnergyLimit(threshold);
136 }
137 
138 inline void
139 G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold(G4double)
140 {
141   G4ExceptionDescription errMsg;
142   errMsg << "*** WARNING : "
143       << "G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold"
144       << "is deprecated, the kill threshold won't be taken into account";
145 
146   G4Exception (
147       "G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold",
148       "DEPRECATED", JustWarning, errMsg);
149 }
150 
151 inline G4double
152 G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold()
153 {
154   G4ExceptionDescription errMsg;
155   errMsg << "*** WARNING : "
156       << "G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold"
157       << "is deprecated, the returned value is nonsense";
158 
159   G4Exception (
160       "G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold",
161       "DEPRECATED", JustWarning, errMsg);
162 
163   return -1;
164 }
165 //---
166 
167 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
168 
169 #endif
170