Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/dna/models/include/G4DNARPWBAIonisationModel.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 // Reference:
 27 //    A.D. Dominguez-Munoz, M.I. Gallardo, M.C. Bordage,
 28 //    Z. Francis, S. Incerti, M.A. Cortes-Giraldo,
 29 //    Radiat. Phys. Chem. 199 (2022) 110363.
 30 //
 31 // Class authors:
 32 //    A.D. Dominguez-Munoz
 33 //    M.A. Cortes-Giraldo (miancortes -at- us.es)
 34 //
 35 // Class creation: 2022-03-03
 36 //
 37 //
 38 
 39 #ifndef G4DNARPWBAIonisationModel_h
 40 #define G4DNARPWBAIonisationModel_h 1
 41 #include "G4VEmModel.hh"
 42 #include "G4ParticleChangeForGamma.hh"
 43 #include "G4ProductionCutsTable.hh"
 44 #include "G4DNACrossSectionDataSet.hh"
 45 #include "G4Electron.hh"
 46 #include "G4Proton.hh"
 47 #include "G4DNAGenericIonsManager.hh"
 48 #include "G4LogLogInterpolation.hh"
 49 #include "G4DNAWaterIonisationStructure.hh"
 50 #include "G4VAtomDeexcitation.hh"
 51 #include "G4NistManager.hh"
 52 
 53 class G4DNARPWBAIonisationModel : public G4VEmModel
 54 {
 55  public:
 56   G4DNARPWBAIonisationModel(const G4ParticleDefinition* p = nullptr,
 57                             const G4String& nam = "DNARPWBAIonisationModel");
 58   ~G4DNARPWBAIonisationModel() override;
 59   G4DNARPWBAIonisationModel& operator=(const G4DNARPWBAIonisationModel& right) =
 60     delete;
 61   G4DNARPWBAIonisationModel(const G4DNARPWBAIonisationModel&) = delete;
 62   void Initialise(const G4ParticleDefinition*,
 63                   const G4DataVector& = *(new G4DataVector())) override;
 64 
 65   G4double CrossSectionPerVolume(const G4Material* material,
 66                                  const G4ParticleDefinition* p, G4double ekin,
 67                                  G4double emin, G4double emax) override;
 68   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
 69                          const G4MaterialCutsCouple*, const G4DynamicParticle*,
 70                          G4double tmin, G4double maxEnergy) override;
 71   G4double GetPartialCrossSection(const G4Material*, G4int /*level*/,
 72                                   const G4ParticleDefinition*,
 73                                   G4double /*kineticEnergy*/) override;
 74   G4double DifferentialCrossSection(const G4double& k, const G4double& energyTransfer,
 75                                     const G4int& shell);
 76   G4double TransferedEnergy(G4double incomingParticleEnergy, G4int shell,
 77                             const G4double& random);
 78   inline void SelectFasterComputation(const G4bool& input);
 79   inline void SelectStationary(const G4bool& input);
 80   inline void SelectSPScaling(const G4bool& input);
 81 
 82  protected:
 83   G4ParticleChangeForGamma* fParticleChangeForGamma = nullptr;
 84 
 85  private:
 86   using MapData = std::map<G4String, std::unique_ptr<G4DNACrossSectionDataSet>,
 87                            std::less<G4String>>;
 88   using TriDimensionMap = std::map<G4double, std::map<G4double, G4double>>;
 89   using VecMap          = std::map<G4double, std::vector<G4double>>;
 90   // methods
 91   G4double RandomizeEjectedElectronEnergy(const G4double&, const G4int& );
 92   G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(const G4double& incomingParticleEnergy,
 93     const G4int& shell);
 94   G4double Interpolate(const G4double& e1, const G4double& e2, const G4double& e, const G4double& xs1,
 95                        const G4double& xs2);
 96   G4double QuadInterpolator(const G4double& e11, const G4double& e12, const G4double& e21,
 97                             const G4double& e22, const G4double& x11, const G4double& x12,
 98                             const G4double& x21, const G4double& x22, const G4double& t1,
 99                             const G4double& t2, const G4double& t, const G4double& e);
100   //shoud change to array ?
101   G4int RandomSelect(G4double energy);
102   void InitialiseForProton(const G4ParticleDefinition*);
103   G4bool InEnergyLimit(const G4double&);
104 
105   //members
106   G4bool fasterCode = false;
107   G4bool statCode   = false;
108   G4bool spScaling  = true;
109   // Water density table
110   const std::vector<G4double>* fpMolWaterDensity = nullptr;
111   // Deexcitation manager to produce fluo photons and e-
112   G4VAtomDeexcitation* fAtomDeexcitation = nullptr;
113   G4double lowEnergyLimit = 0;
114   G4double highEnergyLimit = 0;
115   G4bool isInitialised = false;
116   G4int verboseLevel   = 0;
117   // Cross section
118 
119   std::unique_ptr<G4DNACrossSectionDataSet> fpTotalCrossSection;
120   // Final state
121   G4DNAWaterIonisationStructure waterStructure;
122   TriDimensionMap eDiffCrossSectionData[6];
123   TriDimensionMap eNrjTransfData[6];  // for cumulated dcs
124   TriDimensionMap pDiffCrossSectionData[6];
125   TriDimensionMap pNrjTransfData[6];  // for cumulated dcs
126   std::vector<G4double> eTdummyVec;
127   std::vector<G4double> pTdummyVec;
128   VecMap eVecm;
129   VecMap pVecm;
130   VecMap eProbaShellMap[6];  // for cumulated dcs
131   VecMap pProbaShellMap[6];  // for cumulated dcs
132   const G4ParticleDefinition* fProtonDef = G4Proton::ProtonDefinition();
133 };
134 
135 inline void G4DNARPWBAIonisationModel::SelectFasterComputation(
136   const G4bool& input)
137 {
138   fasterCode = input;
139 }
140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
141 
142 inline void G4DNARPWBAIonisationModel::SelectStationary(const G4bool& input)
143 {
144   statCode = input;
145 }
146 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
147 
148 inline void G4DNARPWBAIonisationModel::SelectSPScaling(const G4bool& input)
149 {
150   spScaling = input;
151 }
152 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
153 
154 #endif
155