Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:59:20

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 // This class is used to support PTB models that come from
0028 // M. Bug et al, Rad. Phys and Chem. 130, 459-479 (2017)
0029 //
0030 
0031 #ifndef G4VDNAModel_HH
0032 #define G4VDNAModel_HH
0033 
0034 #ifdef _MSC_VER
0035 #  pragma warning(disable : 4503)
0036 #endif
0037 
0038 #include "G4DNACrossSectionDataSet.hh"
0039 #include "G4DNAMolecularMaterial.hh"
0040 #include "G4LogLogInterpolation.hh"
0041 #include "G4VEmModel.hh"
0042 #include "G4DNAMaterialManager.hh"
0043 /*! \class G4VDNAModel
0044  * \brief The G4VDNAModel class
0045  *
0046  * All the models using the DNA material management should inherit from that class.
0047  * The goal is to allow the use of the material management system with little code interferences
0048  * within the model classes.
0049  */
0050 class G4VDNAModel : public G4VEmModel
0051 {
0052  public:
0053   /*!
0054    * \brief G4VDNAModel
0055    * Constructeur of the  G4VDNAModel class.
0056    * \param nam
0057    * \param applyToMaterial
0058    */
0059   G4VDNAModel(const G4String& nam, const G4String& applyToMaterial = "");
0060 
0061   /*!
0062    * \brief ~G4VDNAModel
0063    */
0064   ~G4VDNAModel() override;
0065 
0066   /*!
0067    * \brief Initialise
0068    * Each model must implement an Initialize method.
0069    * \param particle
0070    * \param cuts
0071    */
0072   void Initialise(const G4ParticleDefinition* particle, const G4DataVector& cuts) override = 0;
0073 
0074   /*!
0075    * \brief CrossSectionPerVolume
0076    * Every model must implement its own CrossSectionPerVolume method.
0077    * It is used by the process to determine the step path and must return a cross section times a
0078    * number of molecules per volume unit. \param material \param materialName \param p \param ekin
0079    * \param emin
0080    * \param emax
0081    * \return crossSection*numberOfMoleculesPerVolumeUnit
0082    */
0083   G4double CrossSectionPerVolume(const G4Material* material, const G4ParticleDefinition* p,
0084     G4double ekin, G4double emin, G4double emax) override = 0;
0085 
0086   /*!
0087    * \brief SampleSecondaries
0088    * Each model must implement SampleSecondaries to decide if a particle will be created after the
0089    * ModelInterface or if any charateristic of the incident particle will change. \param
0090    * materialName \param particleChangeForGamma \param tmin \param tmax
0091    */
0092   void SampleSecondaries(std::vector<G4DynamicParticle*>*, const G4MaterialCutsCouple*,
0093     const G4DynamicParticle*, G4double tmin = 0, G4double tmax = DBL_MAX) override = 0;
0094 
0095   /*!
0096    * \brief IsMaterialDefine
0097    * Check if the given material is defined in the simulation
0098    * \param materialName
0099    * \return true if the material is defined in the simulation
0100    */
0101   G4bool IsMaterialDefine(const size_t& materialID);
0102 
0103   /*!
0104    * \brief IsParticleExistingInModelForMaterial
0105    * To check two things:
0106    * 1- is the material existing in model ?
0107    * 2- if yes, is the particle defined for that material ?
0108    * \param particleName
0109    * \param materialName
0110    * \return true if the particle/material couple is defined in the model
0111    */
0112   G4bool IsParticleExistingInModelForMaterial(
0113     const G4ParticleDefinition* particleName, const size_t& materialID);
0114 
0115   /*!
0116    * \brief GetName
0117    * \return the name of the model
0118    */
0119   G4String GetName()
0120   {
0121     return fName;
0122   }
0123 
0124   /*!
0125    * \brief GetHighEnergyLimit
0126    * \param material
0127    * \param particle
0128    * \return fHighEnergyLimits[material][particle]
0129    */
0130   G4double GetHighELimit(const size_t& materialID, const G4ParticleDefinition* particle)
0131   {
0132     return fHighEnergyLimits[materialID][particle];
0133   }
0134 
0135   /*!
0136    * \brief GetLowEnergyLimit
0137    * \param material
0138    * \param particle
0139    * \return fLowEnergyLimits[material][particle]
0140    */
0141   G4double GetLowELimit(const size_t& materialID, const G4ParticleDefinition* particle)
0142   {
0143     return fLowEnergyLimits[materialID][particle];
0144   }
0145 
0146   /*!
0147    * \brief SetHighEnergyLimit
0148    * \param material
0149    * \param particle
0150    * \param lim
0151    */
0152   void SetHighELimit(const size_t& materialID, const G4ParticleDefinition* particle, G4double lim)
0153   {
0154     fHighEnergyLimits[materialID][particle] = lim;
0155   }
0156 
0157   /*!
0158    * \brief SetLowEnergyLimit
0159    * \param material
0160    * \param particle
0161    * \param lim
0162    */
0163   void SetLowELimit(const size_t& materialID, const G4ParticleDefinition* particle, G4double lim)
0164   {
0165     fLowEnergyLimits[materialID][particle] = lim;
0166   }
0167 
0168  protected:
0169   G4bool isInitialised = false;
0170 
0171   // typedef used to ease the data container reading
0172   //
0173   using MaterialParticleMapData = std::map<size_t /*MaterialsID*/,
0174     std::map<const G4ParticleDefinition*, std::unique_ptr<G4DNACrossSectionDataSet>>>;
0175 
0176   // Getters
0177   //
0178   /*!
0179    * \brief GetTableData
0180    * \return a pointer to a map with the following structure:
0181    * [materialName][particleName]=G4DNACrossSectionDataSet*
0182    */
0183   MaterialParticleMapData* GetData()
0184   {
0185     return &fData;
0186   }
0187 
0188   // Setters
0189   // ... no setters
0190 
0191   /*!
0192    * \brief BuildApplyToMatVect
0193    * Build the material name vector which is used to know the materials the user want to include in
0194    * the model. \param materials \return a vector with all the material names
0195    */
0196   std::vector<G4String> BuildApplyToMatVect(const G4String& materials);
0197 
0198   /*!
0199    * \brief ReadAndSaveCSFile
0200    * Read and save a "simple" cross section file : use of G4DNACrossSectionDataSet->loadData()
0201    * \param materialName
0202    * \param particleName
0203    * \param file
0204    * \param scaleFactor
0205    */
0206   void ReadAndSaveCSFile(const size_t& materialID, const G4ParticleDefinition* p,
0207     const G4String& file, const G4double& scaleFactor);
0208 
0209   /*!
0210    * \brief RandomSelectShell
0211    * Method to randomely select a shell from the data table uploaded.
0212    * The size of the table (number of columns) is used to determine the total number of possible
0213    * shells. \param k \param particle \param materialName \return the selected shell
0214    */
0215   G4int RandomSelectShell(
0216     const G4double& k, const G4ParticleDefinition* particle, const size_t& materialName);
0217 
0218   /*!
0219    * \brief AddCrossSectionData
0220    * Method used during the initialization of the model class to add a new material. It adds a
0221    * material to the model and fills vectors with informations. \param materialName \param
0222    * particleName \param fileCS \param fileDiffCS \param scaleFactor
0223    */
0224   void AddCrossSectionData(const size_t& materialName, const G4ParticleDefinition* particleName,
0225     const G4String& fileCS, const G4String& fileDiffCS, const G4double& scaleFactor);
0226 
0227   /*!
0228    * \brief AddCrossSectionData
0229    * Method used during the initialization of the model class to add a new material. It adds a
0230    * material to the model and fills vectors with informations. Not every model needs differential
0231    * cross sections. \param materialName \param particleName \param fileCS \param scaleFactor
0232    */
0233   void AddCrossSectionData(const size_t& materialName, const G4ParticleDefinition* particleName,
0234     const G4String& fileCS, const G4double& scaleFactor);
0235 
0236   /*!
0237    * \brief LoadCrossSectionData
0238    * Method to loop on all the registered materials in the model and load the corresponding data.
0239    */
0240   void LoadCrossSectionData(const G4ParticleDefinition* particleName);
0241 
0242   /*!
0243    * \brief ReadDiffCSFile
0244    * Virtual method that need to be implemented if one wish to use the differential cross sections.
0245    * The read method for that kind of information is not standardized yet.
0246    * \param materialName
0247    * \param particleName
0248    * \param path
0249    * \param scaleFactor
0250    */
0251   virtual void ReadDiffCSFile(const size_t& materialName, const G4ParticleDefinition* particleName,
0252     const G4String& path, const G4double& scaleFactor);
0253 
0254   /*!
0255    * \brief EnableMaterialAndParticle
0256    * \param materialName
0257    * \param particleName
0258    * Meant to fill fTableData with 0 for the specified material and particle, therefore allowing the
0259    * ModelInterface class to proceed with the material and particle even if no data are registered
0260    * here. The data should obviously be registered somewhere in the child class. This method is here
0261    * to allow an easy use of the no-ModelInterface dna models within the ModelInterface system.
0262    */
0263   void EnableForMaterialAndParticle(const size_t& materialID, const G4ParticleDefinition* p);
0264 
0265  private:
0266 
0267    /*!
0268    * \brief IsMaterialExistingInModel
0269    * Check if the given material is defined in the current model class
0270    * \param materialName
0271    * \return true if the material is defined in the model
0272     */
0273    G4bool IsMaterialExistingInModel(const size_t& materialID);
0274 
0275   G4String fName;  ///< model name
0276   /*!
0277    * \brief fStringOfMaterials
0278    * The user can decide to specify by hand which are the materials the be activated among those
0279    * implemented in the model. If the user does then only the specified materials contained in this
0280    * string variable will be activated. The string is like: mat1/mat2/mat3/mat4
0281    */
0282   const G4String fStringOfMaterials;
0283 
0284   /*!
0285    * \brief fTableData
0286    * It contains the cross section data and can be used like:
0287    * dataTable=fTableData[material][particle]
0288    */
0289   MaterialParticleMapData fData;
0290 
0291   std::vector<size_t> fModelMaterials;  ///< List the materials that can be activated (and will be
0292                                         ///< by default) within the model.
0293   std::vector<const G4ParticleDefinition*>
0294     fModelParticles;  ///< List the particles that can be activated within the model
0295   std::vector<G4String> fModelCSFiles;  ///< List the cross section data files
0296   std::vector<G4String> fModelDiffCSFiles;  ///< List the differential corss section data files
0297   std::vector<G4double>
0298     fModelScaleFactors;  ///< List the model scale factors (they could change with material)
0299 
0300   std::map<size_t, std::map<const G4ParticleDefinition*, G4double>>
0301     fLowEnergyLimits;  ///< List the low energy limits
0302   std::map<size_t, std::map<const G4ParticleDefinition*, G4double>>
0303     fHighEnergyLimits;  ///< List the high energy limits
0304 };
0305 
0306 #endif  // G4VDNAModel_HH