Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/dna/models/include/G4DNAMillerGreenExcitationModel.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 //
 27 
 28 #ifndef G4DNAMillerGreenExcitationModel_h
 29 #define G4DNAMillerGreenExcitationModel_h 1
 30 
 31 #include "G4VEmModel.hh"
 32 #include "G4ParticleChangeForGamma.hh"
 33 #include "G4ProductionCutsTable.hh"
 34 
 35 #include "G4Proton.hh"
 36 #include "G4DNAGenericIonsManager.hh"
 37 #include "G4DNAWaterExcitationStructure.hh"
 38 #include "Randomize.hh"
 39 #include "G4NistManager.hh"
 40 
 41 #include <deque>
 42 
 43 class G4DNAMillerGreenExcitationModel : public G4VEmModel
 44 {
 45 
 46 public:
 47 
 48   G4DNAMillerGreenExcitationModel(const G4ParticleDefinition* p = nullptr, 
 49               const G4String& nam = "DNAMillerGreenExcitationModel");
 50 
 51   ~G4DNAMillerGreenExcitationModel() override;
 52    
 53   G4DNAMillerGreenExcitationModel & operator=(const  G4DNAMillerGreenExcitationModel &right) = delete;
 54   G4DNAMillerGreenExcitationModel(const  G4DNAMillerGreenExcitationModel&) = delete;
 55 
 56 
 57   void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
 58 
 59   G4double CrossSectionPerVolume(const G4Material* material,
 60              const G4ParticleDefinition* p,
 61              G4double ekin,
 62              G4double emin,
 63              G4double emax) override;
 64 
 65   G4double GetPartialCrossSection(const G4Material*,
 66                                           G4int /*level*/,
 67                                           const G4ParticleDefinition*,
 68                                           G4double /*kineticEnergy*/) override;
 69 
 70   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
 71          const G4MaterialCutsCouple*,
 72          const G4DynamicParticle*,
 73          G4double tmin,
 74          G4double maxEnergy) override;
 75 
 76   inline void SelectStationary(G4bool input); 
 77 
 78 protected:
 79 
 80   G4ParticleChangeForGamma* fParticleChangeForGamma;
 81 
 82 private:
 83 
 84   G4bool statCode;
 85 
 86   // Water density table
 87   const std::vector<G4double>* fpMolWaterDensity;
 88 
 89   std::map<G4String,G4double,std::less<G4String> > lowEnergyLimit;
 90   std::map<G4String,G4double,std::less<G4String> > highEnergyLimit;
 91 
 92   G4bool isInitialised{false};
 93   G4int verboseLevel;
 94   
 95   // Cross section
 96 
 97   // Partial cross section
 98   
 99   G4double PartialCrossSection(G4double energy,G4int level, const G4ParticleDefinition* particle);
100 
101   G4double Sum(G4double energy, const G4ParticleDefinition* particle);
102 
103   G4int RandomSelect(G4double energy, const G4ParticleDefinition* particle);
104 
105   G4int nLevels;
106 
107   G4DNAWaterExcitationStructure waterExcitation;
108   
109   G4double S_1s(G4double t, 
110     G4double energyTransferred, 
111     G4double slaterEffectiveCharge,
112     G4double shellNumber);
113 
114   G4double S_2s(G4double t, 
115     G4double energyTransferred,  
116     G4double slaterEffectiveCharge, 
117     G4double shellNumber);
118 
119   G4double S_2p(G4double t, 
120     G4double energyTransferred, 
121     G4double slaterEffectiveCharge, 
122     G4double shellNumber);
123 
124   G4double R(G4double t, 
125        G4double energyTransferred, 
126        G4double slaterEffectiveCharge, 
127        G4double shellNumber);
128 
129   G4double kineticEnergyCorrection[4]; // 4 is the particle type index
130   G4double slaterEffectiveCharge[3][4];
131   G4double sCoefficient[3][4];
132 
133   //
134   // Reusable particle definitions
135   G4ParticleDefinition* protonDef = nullptr;
136   G4ParticleDefinition* hydrogenDef = nullptr;
137   G4ParticleDefinition* alphaPlusPlusDef = nullptr;
138   G4ParticleDefinition* alphaPlusDef = nullptr;
139   G4ParticleDefinition* heliumDef = nullptr;
140   
141 };
142 
143 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
144 
145 inline void G4DNAMillerGreenExcitationModel::SelectStationary (G4bool input)
146 { 
147     statCode = input; 
148 }    
149 
150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
151 
152 #endif
153