Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/dna/models/include/G4DNADoubleIonisationModel.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 // G4DNADoubleIonisationModel.hh
 27 //
 28 //  Created at 2024/04/03 (Thu.)
 29 //  Author: Shogo OKADA @KEK-CRC (shogo.okada@kek.jp)
 30 //
 31 //  Reference: J.Meesungnoen et. al, DOI: 10.1021/jp058037z
 32 //
 33 
 34 #ifndef G4DNA_DOUBLE_IONISATION_MODEL_HH_
 35 #define G4DNA_DOUBLE_IONISATION_MODEL_HH_
 36 
 37 #include "G4VEmModel.hh"
 38 #include "G4ParticleChangeForGamma.hh"
 39 #include "G4ProductionCutsTable.hh"
 40 
 41 #include "G4DNAGenericIonsManager.hh"
 42 #include "G4DNACrossSectionDataSet.hh"
 43 #include "G4Electron.hh"
 44 #include "G4Proton.hh"
 45 #include "G4LogLogInterpolation.hh"
 46 
 47 #include "G4DNAWaterIonisationStructure.hh"
 48 #include "G4VAtomDeexcitation.hh"
 49 #include "G4NistManager.hh"
 50 #include "G4DNAMultipleIonisationManager.hh"
 51 
 52 using EnergyLimitTable = std::map<G4String, G4double, std::less<G4String>>;
 53 
 54 using CrossSectionDataTable = std::map<G4String, G4DNACrossSectionDataSet*,
 55                                        std::less<G4String>>;
 56 
 57 //==============================================================================
 58 
 59 class G4DNADoubleIonisationModel : public G4VEmModel {
 60 public:
 61 
 62   // constructor
 63   G4DNADoubleIonisationModel(
 64     const G4ParticleDefinition* p = nullptr,
 65     const G4String& model_name = "G4DNADoubleIonisationModel");
 66 
 67   // destructor
 68   ~G4DNADoubleIonisationModel() override;
 69 
 70   G4DNADoubleIonisationModel& operator=(
 71                              const G4DNADoubleIonisationModel&) = delete;
 72   G4DNADoubleIonisationModel(const G4DNADoubleIonisationModel&) = delete;
 73 
 74   void Initialise(
 75     const G4ParticleDefinition* particle, const G4DataVector&) override;
 76 
 77   G4double CrossSectionPerVolume(
 78     const G4Material* material, const G4ParticleDefinition* pdef,
 79     G4double ekin, G4double, G4double) override;
 80 
 81   void SampleSecondaries(
 82     std::vector<G4DynamicParticle*>* vsec, const G4MaterialCutsCouple* couple,
 83     const G4DynamicParticle* particle, G4double, G4double) override;
 84 
 85   void SelectStationary(G4bool in);
 86 
 87   void SelectVerboseLevel(G4int in);
 88 
 89   void UseChampionAlphaParameter(G4bool in);
 90 
 91   void SetMultipleIonisationEnergy(G4double in);
 92 
 93 protected:
 94 
 95   G4double RandomizeEjectedElectronEnergy(
 96     G4ParticleDefinition* pdef, G4double ekin, G4int shell);
 97 
 98   G4ParticleChangeForGamma* particle_change_ = nullptr;
 99 
100   G4int RandomSelect(G4double energy, G4double scale_param,
101                      const G4String& pname);
102 
103   G4double GenerateSecondaries(std::vector<G4DynamicParticle*>* vsec,
104                                const G4MaterialCutsCouple* couple,
105                                const G4DynamicParticle* particle,
106                                G4int ioni_shell,
107                                G4double& theta, G4double& phi,
108                                G4double& shell_energy);
109 
110   G4double GetLowEnergyLimit(const G4String& pname);
111 
112   G4double GetUppEnergyLimit(const G4String& pname);
113 
114   G4bool stat_code_;
115 
116   G4VAtomDeexcitation* atom_deex_;
117 
118   EnergyLimitTable elow_tab_;
119   EnergyLimitTable eupp_tab_;
120 
121   CrossSectionDataTable xs_tab_;
122 
123   G4ParticleDefinition* proton_def_{nullptr};
124   G4ParticleDefinition* alpha_def_{nullptr};
125   G4ParticleDefinition* carbon_def_{nullptr};
126 
127   const std::vector<G4double>* water_density_;
128 
129   G4bool is_initialized_;
130   G4int verbose_level_;
131 
132   std::map<G4double, G4double> model_elow_tab_;
133 
134   G4DNAMultipleIonisationManager* mioni_manager_{nullptr};
135 
136   G4bool use_champion_param_;
137 
138   G4double energy_threshold_;
139 };
140 
141 //==============================================================================
142 inline void G4DNADoubleIonisationModel::SelectStationary(G4bool in)
143 {
144   stat_code_ = in;
145 }
146 
147 //------------------------------------------------------------------------------
148 inline void G4DNADoubleIonisationModel::SelectVerboseLevel(G4int in)
149 {
150   verbose_level_ = in;
151 }
152 
153 //------------------------------------------------------------------------------
154 inline void G4DNADoubleIonisationModel::UseChampionAlphaParameter(G4bool in)
155 {
156   use_champion_param_ = in;
157 }
158 
159 //------------------------------------------------------------------------------
160 inline void G4DNADoubleIonisationModel::SetMultipleIonisationEnergy(G4double in)
161 {
162   energy_threshold_ = in;
163 }
164 
165 #endif // G4DNA_DOUBLE_IONISATION_MODEL_HH_
166