Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:58:09

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 // Authors: S. Meylan and C. Villagrasa (IRSN, France)
0027 // Models come from
0028 // M. Bug et al, Rad. Phys and Chem. 130, 459-479 (2017)
0029 //
0030 
0031 #ifndef G4DNAPTBElasticModel_h
0032 #define G4DNAPTBElasticModel_h 1
0033 
0034 #include "G4DNACrossSectionDataSet.hh"
0035 #include "G4Electron.hh"
0036 #include "G4LogLogInterpolation.hh"
0037 #include "G4NistManager.hh"
0038 #include "G4ParticleChangeForGamma.hh"
0039 #include "G4ProductionCutsTable.hh"
0040 #include "G4VDNAModel.hh"
0041 
0042 #include <map>
0043 
0044 /*!
0045  * \brief The G4DNAPTBElasticModel class
0046  * This class implements the elastic model for the DNA materials and precursors.
0047  */
0048 class G4DNAPTBElasticModel : public G4VDNAModel
0049 {
0050  public:
0051   using TriDimensionMap = std::map<std::size_t,
0052     std::map<const G4ParticleDefinition*, std::map<G4double, std::map<G4double, G4double>>>>;
0053   using VecMap = std::map<std::size_t,
0054     std::map<const G4ParticleDefinition*, std::map<G4double, std::vector<G4double>>>>;
0055   /*!
0056    * \brief G4DNAPTBElasticModel
0057    * Constructor
0058    * \param applyToMaterial
0059    * \param p
0060    * \param nam
0061    */
0062   G4DNAPTBElasticModel(const G4String& applyToMaterial = "all",
0063     const G4ParticleDefinition* p = nullptr, const G4String& nam = "DNAPTBElasticModel");
0064 
0065   /*!
0066    * \brief ~G4DNAPTBElasticModel
0067    * Destructor
0068    */
0069   ~G4DNAPTBElasticModel() override = default;
0070 
0071   // copy constructor and hide assignment operator
0072   G4DNAPTBElasticModel(G4DNAPTBElasticModel&) = delete;  // prevent copy-construction
0073   G4DNAPTBElasticModel& operator=(
0074     const G4DNAPTBElasticModel& right) = delete;  // prevent assignement
0075 
0076   /*!
0077    * \brief Initialise
0078    * Mandatory method for every model class. The material/particle for which the model
0079    * can be used have to be added here through the AddCrossSectionData method.
0080    * Then the LoadCrossSectionData method must be called to trigger the load process.
0081    * Scale factors to be applied to the cross section can be defined here.
0082    */
0083   void Initialise(const G4ParticleDefinition* particle, const G4DataVector&) override;
0084 
0085   /*!
0086    * \brief CrossSectionPerVolume
0087    * This method is mandatory for any model class. It finds and return the cross section value
0088    * for the current material, particle and energy values.
0089    * The number of molecule per volume is not used here but in the G4DNAModelInterface class.
0090    * \param material
0091    * \param materialName
0092    * \param p
0093    * \param ekin
0094    * \param emin
0095    * \param emax
0096    * \return the cross section value
0097    */
0098   G4double CrossSectionPerVolume(const G4Material* material, const G4ParticleDefinition* p,
0099     G4double ekin, G4double emin, G4double emax) override;
0100 
0101   /*!
0102    * \brief SampleSecondaries
0103    * Method called after CrossSectionPerVolume if the process is the one which is selected
0104    * (according to the sampling on the calculated path length). Here, the characteristics of the
0105    * incident and created (if any) particle(s) are set (energy, momentum ...). \param materialName
0106    * \param particleChangeForGamma
0107    * \param tmin
0108    * \param tmax
0109    */
0110   void SampleSecondaries(std::vector<G4DynamicParticle*>*, const G4MaterialCutsCouple*,
0111     const G4DynamicParticle*, G4double tmin, G4double tmax) override;
0112 
0113  protected:
0114   G4ParticleChangeForGamma* fParticleChangeForGamma = nullptr;
0115 
0116  private:
0117   G4int verboseLevel = 0;  ///< verbose level
0118   // Verbosity scale:
0119   // 0 = nothing
0120   // 1 = warning for energy non-conservation
0121   // 2 = details of energy budget
0122   // 3 = calculation of cross sections, file openings, sampling of atoms
0123   // 4 = entering in methods
0124   G4double fKillBelowEnergy = 0.;
0125   ///< energy kill limit
0126   TriDimensionMap diffCrossSectionData;
0127   ///< A map: [materialName][particleName]=DiffCrossSectionTable
0128   VecMap eValuesVect;
0129   /*!< map with vectors containing all the output energy (E) of the diff. file */
0130   std::map<std::size_t, std::map<const G4ParticleDefinition*, std::vector<G4double>>> tValuesVec;
0131   ///< map with vectors containing all the incident (T) energy of the dif. file
0132 
0133   G4Material* fpGuanine_PU = nullptr;
0134   G4Material* fpTHF = nullptr;
0135   G4Material* fpPY = nullptr;
0136   G4Material* fpPU = nullptr;
0137   G4Material* fpTMP = nullptr;
0138   G4Material* fpG4_WATER = nullptr;
0139   G4Material* fpBackbone_THF = nullptr;
0140   G4Material* fpCytosine_PY = nullptr;
0141   G4Material* fpThymine_PY = nullptr;
0142   G4Material* fpAdenine_PU = nullptr;
0143   G4Material* fpBackbone_TMP = nullptr;
0144   G4Material* fpN2 = nullptr;
0145   G4DNAPTBElasticModel* fpModelData = nullptr;
0146 
0147   /*!
0148    * \brief ReadDiffCSFile
0149    * Method to read the differential cross section files. This method is not standard yet so every
0150    * model must implement its own. \param materialName \param particleName \param file
0151    */
0152   void ReadDiffCSFile(const std::size_t& materialID, const G4ParticleDefinition* particleName,
0153     const G4String& file, const G4double&) override;
0154 
0155   /*!
0156    * \brief Theta
0157    * To return an angular theta value from the differential file. This method uses interpolations to
0158    * calculate the theta value. \param fParticleDefinition \param k \param integrDiff \param
0159    * materialName \return a theta value
0160    */
0161   G4double Theta(
0162     const G4ParticleDefinition* p, G4double k, G4double integrDiff, const std::size_t& materialID);
0163 
0164   /*!
0165    * \brief LinLinInterpolate
0166    * \param e1
0167    * \param e2
0168    * \param e
0169    * \param xs1
0170    * \param xs2
0171    * \return
0172    */
0173   G4double LinLinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
0174 
0175   /*!
0176    * \brief LinLogInterpolate
0177    * \param e1
0178    * \param e2
0179    * \param e
0180    * \param xs1
0181    * \param xs2
0182    * \return
0183    */
0184   G4double LinLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
0185 
0186   /*!
0187    * \brief LogLogInterpolate
0188    * \param e1
0189    * \param e2
0190    * \param e
0191    * \param xs1
0192    * \param xs2
0193    * \return
0194    */
0195   G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
0196 
0197   /*!
0198    * \brief QuadInterpolator
0199    * \param e11
0200    * \param e12
0201    * \param e21
0202    * \param e22
0203    * \param x11
0204    * \param x12
0205    * \param x21
0206    * \param x22
0207    * \param t1
0208    * \param t2
0209    * \param t
0210    * \param e
0211    * \return
0212    */
0213   G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double x11,
0214     G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, G4double t, G4double e);
0215 
0216   /*!
0217    * \brief RandomizeCosTheta
0218    * \param k
0219    * \param materialName
0220    * \return
0221    */
0222   G4double RandomizeCosTheta(const G4double& k, const std::size_t& materialName);
0223 
0224 };
0225 
0226 #endif