Geant4-11
G4DNAPTBElasticModel.hh
Go to the documentation of this file.
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
32#ifndef G4DNAPTBElasticModel_h
33#define G4DNAPTBElasticModel_h 1
34
35#include <map>
37#include "G4VDNAModel.hh"
38#include "G4Electron.hh"
42#include "G4NistManager.hh"
43
49{
50
51public:
52
60 G4DNAPTBElasticModel(const G4String &applyToMaterial = "all", const G4ParticleDefinition* p = 0,
61 const G4String& nam = "DNAPTBElasticModel");
62
67 virtual ~G4DNAPTBElasticModel();
68
76 virtual void Initialise(const G4ParticleDefinition* particle, const G4DataVector&, G4ParticleChangeForGamma* fpChangeForGamme=nullptr);
77
92 const G4String& materialName,
93 const G4ParticleDefinition* p,
94 G4double ekin,
95 G4double emin,
97
107 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
109 const G4String& materialName,
110 const G4DynamicParticle*,
111 G4ParticleChangeForGamma *particleChangeForGamma,
112 G4double tmin,
113 G4double tmax);
114
115protected:
116
117
118
119private:
120
122 std::map<G4String, double > killBelowEnergyTable;
124
125 typedef std::map<G4String, std::map<G4String, std::map<double, std::map<double, double> > > > TriDimensionMap;
127
128 typedef std::map<G4String, std::map<G4String, std::map<double, std::vector<double> > > > VecMap;
130 std::map<G4String, std::map<G4String, std::vector<double> > > tValuesVec;
131
139 void ReadDiffCSFile(const G4String &materialName, const G4String &particleName, const G4String &file, const G4double);
140
151 G4double Theta(G4ParticleDefinition * fParticleDefinition, G4double k, G4double integrDiff, const G4String &materialName);
152
163
174
185
203 G4double e12,
204 G4double e21,
205 G4double e22,
206 G4double x11,
207 G4double x12,
208 G4double x21,
209 G4double x22,
210 G4double t1,
211 G4double t2,
212 G4double t,
213 G4double e);
214
221 G4double RandomizeCosTheta(G4double k, const G4String &materialName);
222
223 // copy constructor and hide assignment operator
224 G4DNAPTBElasticModel(G4DNAPTBElasticModel &); // prevent copy-construction
225 G4DNAPTBElasticModel & operator=(const G4DNAPTBElasticModel &right); // prevent assignement
226};
227
228#endif
static const G4double e1[44]
static const G4double e2[44]
static const G4double emax
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
The G4DNAPTBElasticModel class This class implements the elastic model for the DNA materials and prec...
TriDimensionMap diffCrossSectionData
A map: [materialName][particleName]=DiffCrossSectionTable.
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4String &materialName, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
CrossSectionPerVolume This method is mandatory for any model class. It finds and return the cross sec...
G4DNAPTBElasticModel & operator=(const G4DNAPTBElasticModel &right)
G4DNAPTBElasticModel(G4DNAPTBElasticModel &)
G4double RandomizeCosTheta(G4double k, const G4String &materialName)
RandomizeCosTheta.
G4DNAPTBElasticModel(const G4String &applyToMaterial="all", const G4ParticleDefinition *p=0, const G4String &nam="DNAPTBElasticModel")
G4DNAPTBElasticModel Constructor.
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4String &materialName, const G4DynamicParticle *, G4ParticleChangeForGamma *particleChangeForGamma, G4double tmin, G4double tmax)
SampleSecondaries Method called after CrossSectionPerVolume if the process is the one which is select...
virtual ~G4DNAPTBElasticModel()
~G4DNAPTBElasticModel Destructor
G4double Theta(G4ParticleDefinition *fParticleDefinition, G4double k, G4double integrDiff, const G4String &materialName)
Theta To return an angular theta value from the differential file. This method uses interpolations to...
std::map< G4String, std::map< G4String, std::vector< double > > > tValuesVec
map with vectors containing all the incident (T) energy of the differential file
std::map< G4String, std::map< G4String, std::map< double, std::map< double, double > > > > TriDimensionMap
G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
LogLogInterpolate.
void ReadDiffCSFile(const G4String &materialName, const G4String &particleName, const G4String &file, const G4double)
ReadDiffCSFile Method to read the differential cross section files. This method is not standard yet s...
G4double fKillBelowEnergy
energy kill limit
G4double LinLinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
LinLinInterpolate.
G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double x11, G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, G4double t, G4double e)
QuadInterpolator.
std::map< G4String, double > killBelowEnergyTable
map to save the different energy kill limits for the materials
virtual void Initialise(const G4ParticleDefinition *particle, const G4DataVector &, G4ParticleChangeForGamma *fpChangeForGamme=nullptr)
Initialise Mandatory method for every model class. The material/particle for which the model can be u...
G4double LinLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
LinLogInterpolate.
std::map< G4String, std::map< G4String, std::map< double, std::vector< double > > > > VecMap
G4int verboseLevel
verbose level
The G4VDNAModel class.
Definition: G4VDNAModel.hh:50
string material
Definition: eplot.py:19