Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/dna/models/include/G4DNAPTBElasticModel.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 // Authors: S. Meylan and C. Villagrasa (IRSN, France)
 27 // Models come from
 28 // M. Bug et al, Rad. Phys and Chem. 130, 459-479 (2017)
 29 //
 30 
 31 #ifndef G4DNAPTBElasticModel_h
 32 #define G4DNAPTBElasticModel_h 1
 33 
 34 #include "G4DNACrossSectionDataSet.hh"
 35 #include "G4Electron.hh"
 36 #include "G4LogLogInterpolation.hh"
 37 #include "G4NistManager.hh"
 38 #include "G4ParticleChangeForGamma.hh"
 39 #include "G4ProductionCutsTable.hh"
 40 #include "G4VDNAModel.hh"
 41 
 42 #include <map>
 43 
 44 /*!
 45  * \brief The G4DNAPTBElasticModel class
 46  * This class implements the elastic model for the DNA materials and precursors.
 47  */
 48 class G4DNAPTBElasticModel : public G4VDNAModel
 49 {
 50  public:
 51   using TriDimensionMap = std::map<std::size_t,
 52     std::map<const G4ParticleDefinition*, std::map<G4double, std::map<G4double, G4double>>>>;
 53   using VecMap = std::map<std::size_t,
 54     std::map<const G4ParticleDefinition*, std::map<G4double, std::vector<G4double>>>>;
 55   /*!
 56    * \brief G4DNAPTBElasticModel
 57    * Constructor
 58    * \param applyToMaterial
 59    * \param p
 60    * \param nam
 61    */
 62   G4DNAPTBElasticModel(const G4String& applyToMaterial = "all",
 63     const G4ParticleDefinition* p = nullptr, const G4String& nam = "DNAPTBElasticModel");
 64 
 65   /*!
 66    * \brief ~G4DNAPTBElasticModel
 67    * Destructor
 68    */
 69   ~G4DNAPTBElasticModel() override = default;
 70 
 71   // copy constructor and hide assignment operator
 72   G4DNAPTBElasticModel(G4DNAPTBElasticModel&) = delete;  // prevent copy-construction
 73   G4DNAPTBElasticModel& operator=(
 74     const G4DNAPTBElasticModel& right) = delete;  // prevent assignement
 75 
 76   /*!
 77    * \brief Initialise
 78    * Mandatory method for every model class. The material/particle for which the model
 79    * can be used have to be added here through the AddCrossSectionData method.
 80    * Then the LoadCrossSectionData method must be called to trigger the load process.
 81    * Scale factors to be applied to the cross section can be defined here.
 82    */
 83   void Initialise(const G4ParticleDefinition* particle, const G4DataVector&) override;
 84 
 85   /*!
 86    * \brief CrossSectionPerVolume
 87    * This method is mandatory for any model class. It finds and return the cross section value
 88    * for the current material, particle and energy values.
 89    * The number of molecule per volume is not used here but in the G4DNAModelInterface class.
 90    * \param material
 91    * \param materialName
 92    * \param p
 93    * \param ekin
 94    * \param emin
 95    * \param emax
 96    * \return the cross section value
 97    */
 98   G4double CrossSectionPerVolume(const G4Material* material, const G4ParticleDefinition* p,
 99     G4double ekin, G4double emin, G4double emax) override;
100 
101   /*!
102    * \brief SampleSecondaries
103    * Method called after CrossSectionPerVolume if the process is the one which is selected
104    * (according to the sampling on the calculated path length). Here, the characteristics of the
105    * incident and created (if any) particle(s) are set (energy, momentum ...). \param materialName
106    * \param particleChangeForGamma
107    * \param tmin
108    * \param tmax
109    */
110   void SampleSecondaries(std::vector<G4DynamicParticle*>*, const G4MaterialCutsCouple*,
111     const G4DynamicParticle*, G4double tmin, G4double tmax) override;
112 
113  protected:
114   G4ParticleChangeForGamma* fParticleChangeForGamma = nullptr;
115 
116  private:
117   G4int verboseLevel = 0;  ///< verbose level
118   // Verbosity scale:
119   // 0 = nothing
120   // 1 = warning for energy non-conservation
121   // 2 = details of energy budget
122   // 3 = calculation of cross sections, file openings, sampling of atoms
123   // 4 = entering in methods
124   G4double fKillBelowEnergy = 0.;
125   ///< energy kill limit
126   TriDimensionMap diffCrossSectionData;
127   ///< A map: [materialName][particleName]=DiffCrossSectionTable
128   VecMap eValuesVect;
129   /*!< map with vectors containing all the output energy (E) of the diff. file */
130   std::map<std::size_t, std::map<const G4ParticleDefinition*, std::vector<G4double>>> tValuesVec;
131   ///< map with vectors containing all the incident (T) energy of the dif. file
132 
133   G4Material* fpGuanine_PU = nullptr;
134   G4Material* fpTHF = nullptr;
135   G4Material* fpPY = nullptr;
136   G4Material* fpPU = nullptr;
137   G4Material* fpTMP = nullptr;
138   G4Material* fpG4_WATER = nullptr;
139   G4Material* fpBackbone_THF = nullptr;
140   G4Material* fpCytosine_PY = nullptr;
141   G4Material* fpThymine_PY = nullptr;
142   G4Material* fpAdenine_PU = nullptr;
143   G4Material* fpBackbone_TMP = nullptr;
144   G4Material* fpN2 = nullptr;
145   G4DNAPTBElasticModel* fpModelData = nullptr;
146 
147   /*!
148    * \brief ReadDiffCSFile
149    * Method to read the differential cross section files. This method is not standard yet so every
150    * model must implement its own. \param materialName \param particleName \param file
151    */
152   void ReadDiffCSFile(const std::size_t& materialID, const G4ParticleDefinition* particleName,
153     const G4String& file, const G4double&) override;
154 
155   /*!
156    * \brief Theta
157    * To return an angular theta value from the differential file. This method uses interpolations to
158    * calculate the theta value. \param fParticleDefinition \param k \param integrDiff \param
159    * materialName \return a theta value
160    */
161   G4double Theta(
162     const G4ParticleDefinition* p, G4double k, G4double integrDiff, const std::size_t& materialID);
163 
164   /*!
165    * \brief LinLinInterpolate
166    * \param e1
167    * \param e2
168    * \param e
169    * \param xs1
170    * \param xs2
171    * \return
172    */
173   G4double LinLinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
174 
175   /*!
176    * \brief LinLogInterpolate
177    * \param e1
178    * \param e2
179    * \param e
180    * \param xs1
181    * \param xs2
182    * \return
183    */
184   G4double LinLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
185 
186   /*!
187    * \brief LogLogInterpolate
188    * \param e1
189    * \param e2
190    * \param e
191    * \param xs1
192    * \param xs2
193    * \return
194    */
195   G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
196 
197   /*!
198    * \brief QuadInterpolator
199    * \param e11
200    * \param e12
201    * \param e21
202    * \param e22
203    * \param x11
204    * \param x12
205    * \param x21
206    * \param x22
207    * \param t1
208    * \param t2
209    * \param t
210    * \param e
211    * \return
212    */
213   G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double x11,
214     G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, G4double t, G4double e);
215 
216   /*!
217    * \brief RandomizeCosTheta
218    * \param k
219    * \param materialName
220    * \return
221    */
222   G4double RandomizeCosTheta(const G4double& k, const std::size_t& materialName);
223 
224 };
225 
226 #endif
227