Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4DNAPTBIonisationModel.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 G4DNAPTBIONISATIONMODEL_h
32 #define G4DNAPTBIONISATIONMODEL_h 1
33 
34 #include "G4VDNAModel.hh"
36 #include "G4ProductionCutsTable.hh"
37 
39 #include "G4Electron.hh"
40 #include "G4Proton.hh"
42 
43 #include "G4LogLogInterpolation.hh"
44 
46 #include "G4DNAPTBAugerModel.hh"
47 #include "G4NistManager.hh"
48 
54 {
55 
56 public:
65  G4DNAPTBIonisationModel(const G4String &applyToMaterial = "all",
66  const G4ParticleDefinition* p = 0,
67  const G4String &nam = "DNAPTBIonisationModel",
68  const G4bool isAuger=true);
69 
74  virtual ~G4DNAPTBIonisationModel();
75 
81  virtual void Initialise(const G4ParticleDefinition* particle, const G4DataVector& = *(new G4DataVector()), G4ParticleChangeForGamma* fpChangeForGamme=nullptr);
82 
95  virtual G4double CrossSectionPerVolume(const G4Material* material,
96  const G4String& materialName,
97  const G4ParticleDefinition* p,
98  G4double ekin,
99  G4double emin,
100  G4double emax);
101 
112  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
113  const G4MaterialCutsCouple*,
114  const G4String& materialName,
115  const G4DynamicParticle*,
116  G4ParticleChangeForGamma *particleChangeForGamma,
117  G4double tmin,
118  G4double tmax);
119 
120 protected:
121 
122 private:
123 
125 
127 
130  typedef std::map<G4String, std::map<G4String, std::map<double, std::map<double, std::map<double, double> > > > > TriDimensionMap;
133  std::map<G4String, std::map<G4String, std::vector<double> > > fTMapWithVec;
134  typedef std::map<G4String, std::map<G4String, std::map<double, std::vector<double> > > > VecMap;
136  typedef std::map<G4String, std::map<G4String, std::map<double, std::map<double, std::vector<double> > > > > VecMapWithShell;
138 
139  G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell, const G4String& materialName);
140  double DifferentialCrossSection(G4ParticleDefinition * aParticleDefinition, G4double k, G4double energyTransfer, G4int shell, const G4String &materialName);
141 
151  G4double RandomizeEjectedElectronEnergyFromCumulated(G4ParticleDefinition *particleDefinition, G4double k, G4int shell, const G4String& materialName);
152 
162  void RandomizeEjectedElectronDirection(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4double
163  outgoingParticleEnergy, G4double & cosTheta, G4double & phi );
172  void ReadDiffCSFile(const G4String &materialName, const G4String &particleName, const G4String &file, const G4double scaleFactor);
173 
201 
202  // copy constructor and hide assignment operator
203  G4DNAPTBIonisationModel(const G4DNAPTBIonisationModel&); // prevent copy-construction
204  G4DNAPTBIonisationModel & operator=(const G4DNAPTBIonisationModel &right); // prevent assignement
205 };
206 
207 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
208 
209 #endif
std::map< G4String, std::map< G4String, std::map< double, std::vector< double > > > > VecMap
virtual void Initialise(const G4ParticleDefinition *particle, const G4DataVector &=*(new G4DataVector()), G4ParticleChangeForGamma *fpChangeForGamme=nullptr)
Initialise Method called once at the beginning of the simulation. It is used to setup the list of the...
The G4DNAPTBIonisationModel class Implements the PTB ionisation model.
The G4DNAPTBAugerModel class Implement the PTB Auger model.
TTree * t1
Definition: plottest35.C:26
virtual ~G4DNAPTBIonisationModel()
~G4DNAPTBIonisationModel Destructor
const char * p
Definition: xmltok.h:285
std::map< G4String, std::map< G4String, std::map< double, std::map< double, std::vector< double > > > > > VecMapWithShell
G4DNAPTBIonisationStructure ptbStructure
static const G4double emax
std::map< G4String, std::map< G4String, std::vector< double > > > fTMapWithVec
G4DNAPTBIonisationModel & operator=(const G4DNAPTBIonisationModel &right)
std::map< G4String, std::map< G4String, std::map< double, std::map< double, std::map< double, double > > > > > TriDimensionMap
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4double RandomizeEjectedElectronEnergyFromCumulated(G4ParticleDefinition *particleDefinition, G4double k, G4int shell, const G4String &materialName)
RandomizeEjectedElectronEnergyFromCumulated Uses the cumulated tables to find the energy of the eject...
G4DNAPTBIonisationModel(const G4String &applyToMaterial="all", const G4ParticleDefinition *p=0, const G4String &nam="DNAPTBIonisationModel", const G4bool isAuger=true)
G4DNAPTBIonisationModel Constructor.
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 SampleSecondaries will be call...
double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell, const G4String &materialName)
TTree * t2
Definition: plottest35.C:36
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4String &materialName, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
CrossSectionPerVolume Mandatory for every model the CrossSectionPerVolume method is in charge of retu...
void ReadDiffCSFile(const G4String &materialName, const G4String &particleName, const G4String &file, const G4double scaleFactor)
ReadDiffCSFile Method to read the differential cross section files.
G4DNAPTBAugerModel * fDNAPTBAugerModel
PTB Auger model instanciated in the constructor and deleted in the destructor of the class...
int G4int
Definition: G4Types.hh:78
G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
LogLogInterpolate.
void RandomizeEjectedElectronDirection(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4double outgoingParticleEnergy, G4double &cosTheta, G4double &phi)
RandomizeEjectedElectronDirection Method to calculate the ejected electron direction.
G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double xs11, G4double xs12, G4double xs21, G4double xs22, G4double t1, G4double t2, G4double t, G4double e)
QuadInterpolator.
TFile * file
The G4VDNAModel class.
Definition: G4VDNAModel.hh:49
G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, const G4String &materialName)