Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/dna/models/include/G4DNADiracRMatrixExcitationModel.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 // Created on 2016/04/08
 27 //
 28 // Authors: D. Sakata, S. Incerti
 29 //
 30 // This class perform electric excitation for electron transportation,
 31 // based on Dirac B-Spline R-Matrix Model and scaled experimental data.
 32 // See following reference paper 
 33 // Phys.Rev.A77,062711(2008) and Phys.Rev.A78,042713(2008)  
 34 
 35 #ifndef G4DNADiracRMatrixExcitationModel_h
 36 #define G4DNADiracRMatrixExcitationModel_h 1
 37 
 38 #include "G4VEmModel.hh"
 39 #include "G4ParticleChangeForGamma.hh"
 40 #include "G4ProductionCutsTable.hh"
 41 #include "G4VAtomDeexcitation.hh"
 42 
 43 #include "G4LogLogInterpolation.hh"
 44 #include "G4Electron.hh"
 45 #include "G4Proton.hh"
 46 #include "G4NistManager.hh"
 47 
 48 #include "G4DNACrossSectionDataSet.hh"
 49 
 50 class G4DNADiracRMatrixExcitationModel: public G4VEmModel
 51 {
 52 
 53 public:
 54 
 55   G4DNADiracRMatrixExcitationModel(const G4ParticleDefinition* p = nullptr,
 56                         const G4String& nam = "DNADiracRMatrixExcitationModel");
 57 
 58   ~G4DNADiracRMatrixExcitationModel() override;
 59 
 60   G4DNADiracRMatrixExcitationModel & operator
 61                            =(const  G4DNADiracRMatrixExcitationModel &right) = delete;
 62   G4DNADiracRMatrixExcitationModel(const  G4DNADiracRMatrixExcitationModel&) = delete;
 63 
 64 
 65   void Initialise(const G4ParticleDefinition*,
 66                           const G4DataVector& = *(new G4DataVector())) override;
 67 
 68   G4double CrossSectionPerVolume(const G4Material* material,
 69                                          const G4ParticleDefinition* p,
 70                                          G4double ekin,
 71                                          G4double emin,
 72                                          G4double emax) override;
 73 
 74   virtual G4double GetExtendedTotalCrossSection  (const G4Material* material,
 75                                           const G4ParticleDefinition*,
 76                                           G4double kineticEnergy);
 77   
 78   virtual G4double GetExtendedPartialCrossSection(const G4Material* material,
 79                                           G4int level,
 80                                           const G4ParticleDefinition*,
 81                                           G4double kineticEnergy);
 82 
 83   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
 84                                  const G4MaterialCutsCouple*,
 85                                  const G4DynamicParticle*,
 86                                  G4double tmin,
 87                                  G4double maxEnergy) override;
 88 
 89   inline  void SelectStationary(G4bool input);
 90 
 91 protected:
 92 
 93   G4ParticleChangeForGamma* fParticleChangeForGamma;
 94 
 95 private:
 96 
 97   const G4double paramFuncTCS_5dto6s1[3]={-3e-50      , 9.46358e-16,  1.4237  }; // y = [0]+[1]/pow(x-[2],2)
 98   const G4double paramFuncTCS_5dto6s2[3]={-3e-50      , 4.24498e-15, -0.674543}; // y = [0]+[1]/pow(x-[2],2)
 99   const G4double paramFuncTCS_6sto6p1[3]={ 1.50018e-26, 2.459e-15  ,-40.8088  }; // y = [0]+[1]*log(x-[2])/(x-[2])
100   const G4double paramFuncTCS_6sto6p2[3]={ 1.26684e-25, 3.97221e-15,-55.6954  }; // y = [0]+[1]*log(x-[2])/(x-[2])
101   const G4double ExcitationEnergyAu[4]={ 2.66 , 1.14 , 4.63 ,  5.11}; 
102   // [eV] 5dto6s1,6sto6p1,6sto6p2
103 
104   G4double fLowEnergyLimit=0.;
105   G4double fExperimentalEnergyLimit=0.;
106   G4double fHighEnergyLimit=0.;
107 
108   G4bool   isInitialised=false;
109   G4bool   statCode=false;
110   G4int    verboseLevel=0;
111 
112   G4String                     fTableFile="";
113   G4DNACrossSectionDataSet*    fTableData=nullptr;
114   const std::vector<G4double>* fpMaterialDensity=nullptr;
115   const  G4ParticleDefinition* fParticleDefinition=nullptr;
116 
117   G4int RandomSelect(const G4Material* material,
118                      const G4ParticleDefinition*,
119                      G4double kineticEnergy);
120   
121    
122 };
123 
124 inline void G4DNADiracRMatrixExcitationModel::SelectStationary(G4bool input)
125 {
126   statCode = input;
127 }
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
130 
131 #endif
132