Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4DNAPTBExcitationModel.hh
이 파일의 문서화 페이지로 가기
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 // Authors: S. Meylan and C. Villagrasa (IRSN, France)
27 // Models come from
28 // M. Bug et al, Rad. Phys and Chem. 130, 459-479 (2017)
29 //
30 
31 #ifndef G4DNAPTBExcitationModel_h
32 #define G4DNAPTBExcitationModel_h 1
33 
34 #include "G4VDNAModel.hh"
36 #include "G4ProductionCutsTable.hh"
37 
39 #include "G4LogLogInterpolation.hh"
40 #include "G4Electron.hh"
41 #include "G4Proton.hh"
42 #include "G4NistManager.hh"
43 
45 
51 {
52 
53 public:
54 
62  G4DNAPTBExcitationModel(const G4String &applyToMaterial = "all", const G4ParticleDefinition* p = 0,
63  const G4String& nam = "DNAPTBExcitationModel");
64 
69  virtual ~G4DNAPTBExcitationModel();
70 
75  virtual void Initialise(const G4ParticleDefinition* particle, const G4DataVector& = *(new G4DataVector()), G4ParticleChangeForGamma* fpChangeForGamme=nullptr);
76 
88  virtual G4double CrossSectionPerVolume(const G4Material* material,
89  const G4String& materialName,
90  const G4ParticleDefinition* p,
91  G4double ekin,
92  G4double emin,
93  G4double emax);
94 
104  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
105  const G4MaterialCutsCouple*,
106  const G4String& materialName,
107  const G4DynamicParticle*,
108  G4ParticleChangeForGamma *particleChangeForGamma,
109  G4double tmin,
110  G4double tmax);
111 
112 protected:
113 
114 private:
115 
117 
119 
120  typedef std::map<G4String,G4double,std::less<G4String> > MapMeanEnergy;
122 
123  // copy constructor and hide assignment operator
124  G4DNAPTBExcitationModel(const G4DNAPTBExcitationModel&); // prevent copy-construction
125  G4DNAPTBExcitationModel & operator=(const G4DNAPTBExcitationModel &right); // prevent assignement
126 };
127 
128 #endif
const char * p
Definition: xmltok.h:285
virtual ~G4DNAPTBExcitationModel()
~G4DNAPTBExcitationModel Destructor
static const G4double emax
G4DNAPTBExcitationModel & operator=(const G4DNAPTBExcitationModel &right)
double G4double
Definition: G4Types.hh:76
G4DNAWaterExcitationStructure waterStructure
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4String &materialName, const G4DynamicParticle *, G4ParticleChangeForGamma *particleChangeForGamma, G4double tmin, G4double tmax)
SampleSecondaries If the model is selected for the ModelInterface then the SampleSecondaries method w...
int G4int
Definition: G4Types.hh:78
G4DNAPTBExcitationModel(const G4String &applyToMaterial="all", const G4ParticleDefinition *p=0, const G4String &nam="DNAPTBExcitationModel")
G4DNAPTBExcitationModel Constructor.
std::map< G4String, G4double, std::less< G4String > > MapMeanEnergy
The G4DNAPTBExcitationModel class This class implements the PTB excitation model. ...
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4String &materialName, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
CrossSectionPerVolume Retrieve the cross section corresponding to the current material, particle and energy.
MapMeanEnergy tableMeanEnergyPTB
map: [materialName]=energyValue
The G4VDNAModel class.
Definition: G4VDNAModel.hh:49
virtual void Initialise(const G4ParticleDefinition *particle, const G4DataVector &=*(new G4DataVector()), G4ParticleChangeForGamma *fpChangeForGamme=nullptr)
Initialise Set the materials for which the model can be used and defined the energy limits...