Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/dna/models/include/G4DNASancheExcitationModel.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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 29 
 30 // Created by Z. Francis
 31 
 32 #ifndef G4DNASancheExcitationModel_h
 33 #define G4DNASancheExcitationModel_h 1
 34 
 35 #include <deque>
 36 #include <CLHEP/Units/SystemOfUnits.h>
 37 
 38 #include "G4VEmModel.hh"
 39 #include "G4ParticleChangeForGamma.hh"
 40 #include "G4Electron.hh"
 41 #include "G4NistManager.hh"
 42 
 43 class G4DNASancheExcitationModel : public G4VEmModel
 44 {
 45 public:
 46   G4DNASancheExcitationModel(const G4ParticleDefinition* p = nullptr,
 47                              const G4String& nam = "DNASancheExcitationModel");
 48 
 49   ~G4DNASancheExcitationModel() override;
 50 
 51   G4DNASancheExcitationModel & operator=(const G4DNASancheExcitationModel &right) = delete;
 52   G4DNASancheExcitationModel(const G4DNASancheExcitationModel&) = delete;
 53 
 54   void Initialise(const G4ParticleDefinition*,
 55                           const G4DataVector&) override;
 56 
 57   G4double CrossSectionPerVolume(const G4Material* material,
 58                                          const G4ParticleDefinition* p,
 59                                          G4double ekin,
 60                                          G4double emin,
 61                                          G4double emax) override;
 62 
 63   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
 64                                  const G4MaterialCutsCouple*,
 65                                  const G4DynamicParticle*,
 66                                  G4double tmin,
 67                                  G4double maxEnergy) override;
 68 
 69   // Cross section
 70 
 71   G4double PartialCrossSection(G4double energy, G4int level);
 72   G4double TotalCrossSection(G4double t);
 73 
 74   inline void ExtendLowEnergyLimit(G4double /*threshold*/);
 75 
 76   inline void SetVerboseLevel(G4int verbose)
 77   {
 78     verboseLevel = verbose;
 79   }
 80 
 81   inline void SelectStationary(G4bool input); 
 82 
 83 protected:
 84 
 85   G4ParticleChangeForGamma* fParticleChangeForGamma;
 86 
 87 private:
 88 
 89   G4bool statCode;
 90 
 91   // Water density table
 92   const std::vector<G4double>* fpWaterDensity;
 93 
 94   G4bool isInitialised{false};
 95   G4int verboseLevel;
 96 
 97   // Cross section
 98 
 99   G4int RandomSelect(G4double energy);
100   G4int nLevels;
101   G4double VibrationEnergy(G4int level);
102   G4double Sum(G4double k);
103   G4double LinInterpolate(G4double e1,
104                           G4double e2,
105                           G4double e,
106                           G4double xs1,
107                           G4double xs2);
108 
109   //
110 //  typedef std::map<double, std::map<double, double> > TriDimensionMap;
111 //  TriDimensionMap map1;
112   std::vector<G4double> tdummyVec;
113   std::vector<std::vector<G4double>> fEnergyLevelXS;
114   std::vector<G4double> fEnergyTotalXS;
115 
116 };
117 
118 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
119 
120 inline void G4DNASancheExcitationModel::ExtendLowEnergyLimit(G4double threshold)
121 {
122   if(threshold < 2 * CLHEP::eV)
123     G4Exception("*** WARNING : the G4DNASancheExcitationModel class is not "
124                 "validated below 2 eV !",
125                 "", JustWarning, "");
126   
127   SetLowEnergyLimit(threshold);
128 }
129 
130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
131 
132 inline void G4DNASancheExcitationModel::SelectStationary (G4bool input)
133 { 
134     statCode = input; 
135 }    
136 
137 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
138 
139 #endif
140