Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAPTBElasticModel Class Reference

The G4DNAPTBElasticModel class This class implements the elastic model for the DNA materials and precursors. More...

#include <G4DNAPTBElasticModel.hh>

+ Inheritance diagram for G4DNAPTBElasticModel:

Public Member Functions

 G4DNAPTBElasticModel (const G4String &applyToMaterial="all", const G4ParticleDefinition *p=0, const G4String &nam="DNAPTBElasticModel")
 G4DNAPTBElasticModel Constructor.
 
virtual ~G4DNAPTBElasticModel ()
 ~G4DNAPTBElasticModel Destructor
 
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 used have to be added here through the AddCrossSectionData method. Then the LoadCrossSectionData method must be called to trigger the load process. Scale factors to be applied to the cross section can be defined here.
 
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 section value for the current material, particle and energy values. The number of molecule per volume is not used here but in the G4DNAModelInterface class.
 
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 selected (according to the sampling on the calculated path length). Here, the characteristics of the incident and created (if any) particle(s) are set (energy, momentum ...).
 
- Public Member Functions inherited from G4VDNAModel
 G4VDNAModel (const G4String &nam, const G4String &applyToMaterial)
 G4VDNAModel Constructeur of the G4VDNAModel class.
 
virtual ~G4VDNAModel ()
 ~G4VDNAModel
 
virtual void Initialise (const G4ParticleDefinition *particle, const G4DataVector &cuts, G4ParticleChangeForGamma *fpChangeForGamme=nullptr)=0
 Initialise Each model must implement an Initialize method.
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4String &materialName, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)=0
 CrossSectionPerVolume Every model must implement its own CrossSectionPerVolume method. It is used by the process to determine the step path and must return a cross section times a number of molecules per volume unit.
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4String &materialName, const G4DynamicParticle *, G4ParticleChangeForGamma *particleChangeForGamma, G4double tmin=0, G4double tmax=DBL_MAX)=0
 SampleSecondaries Each model must implement SampleSecondaries to decide if a particle will be created after the ModelInterface or if any charateristic of the incident particle will change.
 
G4bool IsMaterialDefine (const G4String &materialName)
 IsMaterialDefine Check if the given material is defined in the simulation.
 
G4bool IsMaterialExistingInModel (const G4String &materialName)
 IsMaterialExistingInModel Check if the given material is defined in the current model class.
 
G4bool IsParticleExistingInModelForMaterial (const G4String &particleName, const G4String &materialName)
 IsParticleExistingInModelForMaterial To check two things: 1- is the material existing in model ? 2- if yes, is the particle defined for that material ?
 
G4String GetName ()
 GetName.
 
G4double GetHighELimit (const G4String &material, const G4String &particle)
 GetHighEnergyLimit.
 
G4double GetLowELimit (const G4String &material, const G4String &particle)
 GetLowEnergyLimit.
 
void SetHighELimit (const G4String &material, const G4String &particle, G4double lim)
 SetHighEnergyLimit.
 
void SetLowELimit (const G4String &material, const G4String &particle, G4double lim)
 SetLowEnergyLimit.
 

Additional Inherited Members

- Protected Types inherited from G4VDNAModel
typedef std::map< G4String, std::map< G4String, G4DNACrossSectionDataSet *, std::less< G4String > > > TableMapData
 
typedef std::map< G4String, std::map< G4String, G4double > > RatioMapData
 
typedef std::map< G4String, G4double >::const_iterator ItCompoMapData
 
- Protected Member Functions inherited from G4VDNAModel
TableMapDataGetTableData ()
 GetTableData.
 
std::vector< G4StringBuildApplyToMatVect (const G4String &materials)
 BuildApplyToMatVect Build the material name vector which is used to know the materials the user want to include in the model.
 
void ReadAndSaveCSFile (const G4String &materialName, const G4String &particleName, const G4String &file, G4double scaleFactor)
 ReadAndSaveCSFile Read and save a "simple" cross section file : use of G4DNACrossSectionDataSet->loadData()
 
G4int RandomSelectShell (G4double k, const G4String &particle, const G4String &materialName)
 RandomSelectShell Method to randomely select a shell from the data table uploaded. The size of the table (number of columns) is used to determine the total number of possible shells.
 
void AddCrossSectionData (G4String materialName, G4String particleName, G4String fileCS, G4String fileDiffCS, G4double scaleFactor)
 AddCrossSectionData Method used during the initialization of the model class to add a new material. It adds a material to the model and fills vectors with informations.
 
void AddCrossSectionData (G4String materialName, G4String particleName, G4String fileCS, G4double scaleFactor)
 AddCrossSectionData Method used during the initialization of the model class to add a new material. It adds a material to the model and fills vectors with informations. Not every model needs differential cross sections.
 
void LoadCrossSectionData (const G4String &particleName)
 LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresponding data.
 
virtual void ReadDiffCSFile (const G4String &materialName, const G4String &particleName, const G4String &path, const G4double scaleFactor)
 ReadDiffCSFile Virtual method that need to be implemented if one wish to use the differential cross sections. The read method for that kind of information is not standardized yet.
 
void EnableForMaterialAndParticle (const G4String &materialName, const G4String &particleName)
 EnableMaterialAndParticle.
 

Detailed Description

The G4DNAPTBElasticModel class This class implements the elastic model for the DNA materials and precursors.

Definition at line 48 of file G4DNAPTBElasticModel.hh.

Constructor & Destructor Documentation

◆ G4DNAPTBElasticModel()

G4DNAPTBElasticModel::G4DNAPTBElasticModel ( const G4String applyToMaterial = "all",
const G4ParticleDefinition p = 0,
const G4String nam = "DNAPTBElasticModel" 
)

G4DNAPTBElasticModel Constructor.

Parameters
applyToMaterial
p
nam

Definition at line 38 of file G4DNAPTBElasticModel.cc.

40 : G4VDNAModel(nam, applyToMaterial)
41{
42 fKillBelowEnergy = 10*eV; // will be override by the limits defined for each material
43
44 verboseLevel= 0;
45 // Verbosity scale:
46 // 0 = nothing
47 // 1 = warning for energy non-conservation
48 // 2 = details of energy budget
49 // 3 = calculation of cross sections, file openings, sampling of atoms
50 // 4 = entering in methods
51
52 if( verboseLevel>0 )
53 {
54 G4cout << "PTB Elastic model is constructed " << G4endl;
55 }
56}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
The G4VDNAModel class.
Definition: G4VDNAModel.hh:50

◆ ~G4DNAPTBElasticModel()

G4DNAPTBElasticModel::~G4DNAPTBElasticModel ( )
virtual

~G4DNAPTBElasticModel Destructor

Definition at line 60 of file G4DNAPTBElasticModel.cc.

61{
62
63}

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4DNAPTBElasticModel::CrossSectionPerVolume ( const G4Material material,
const G4String materialName,
const G4ParticleDefinition p,
G4double  ekin,
G4double  emin,
G4double  emax 
)
virtual

CrossSectionPerVolume This method is mandatory for any model class. It finds and return the cross section value for the current material, particle and energy values. The number of molecule per volume is not used here but in the G4DNAModelInterface class.

Parameters
material
materialName
p
ekin
emin
emax
Returns
the cross section value

Implements G4VDNAModel.

Definition at line 305 of file G4DNAPTBElasticModel.cc.

311{
312 if (verboseLevel > 3)
313 G4cout << "Calling CrossSectionPerVolume() of G4DNAPTBElasticModel" << G4endl;
314
315 // Get the name of the current particle
316 const G4String& particleName = p->GetParticleName();
317
318 // set killBelowEnergy value for current material
319 fKillBelowEnergy = GetLowELimit(materialName, particleName);
320
321 // initialise the return value (cross section) to zero
322 G4double sigma(0);
323
324 // check if we are below the high energy limit
325 if (ekin < GetHighELimit(materialName, particleName) )
326 {
327 // This is used to kill the particle if its kinetic energy is below fKillBelowEnergy.
328 // If the energy is lower then we return a maximum cross section and thus the SampleSecondaries method will be called for sure.
329 // SampleSecondaries will remove the particle from the simulation.
330 //
331 //SI : XS must not be zero otherwise sampling of secondaries method ignored
332 if (ekin < fKillBelowEnergy) return DBL_MAX;
333
334 // Get the tables with the cross section data
335 TableMapData* tableData = GetTableData();
336
337 // Retrieve the cross section value
338 sigma = (*tableData)[materialName][particleName]->FindValue(ekin);
339 }
340
341 if (verboseLevel > 2)
342 {
343 G4cout << "__________________________________" << G4endl;
344 G4cout << "°°° G4DNAPTBElasticModel - XS INFO START" << G4endl;
345 G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
346 G4cout << "°°° Cross section per molecule (cm^2)=" << sigma/cm/cm << G4endl;
347 G4cout << "°°° G4DNAPTBElasticModel - XS INFO END" << G4endl;
348 }
349
350 // Return the cross section
351 return sigma;
352}
double G4double
Definition: G4Types.hh:83
const G4String & GetParticleName() const
TableMapData * GetTableData()
GetTableData.
Definition: G4VDNAModel.hh:193
std::map< G4String, std::map< G4String, G4DNACrossSectionDataSet *, std::less< G4String > > > TableMapData
Definition: G4VDNAModel.hh:183
G4double GetHighELimit(const G4String &material, const G4String &particle)
GetHighEnergyLimit.
Definition: G4VDNAModel.hh:153
G4double GetLowELimit(const G4String &material, const G4String &particle)
GetLowEnergyLimit.
Definition: G4VDNAModel.hh:161
#define DBL_MAX
Definition: templates.hh:62

◆ Initialise()

void G4DNAPTBElasticModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector ,
G4ParticleChangeForGamma fpChangeForGamme = nullptr 
)
virtual

Initialise Mandatory method for every model class. The material/particle for which the model can be used have to be added here through the AddCrossSectionData method. Then the LoadCrossSectionData method must be called to trigger the load process. Scale factors to be applied to the cross section can be defined here.

Implements G4VDNAModel.

Definition at line 67 of file G4DNAPTBElasticModel.cc.

69{
70 if (verboseLevel > 3)
71 G4cout << "Calling G4DNAPTBElasticModel::Initialise()" << G4endl;
72
73 G4double scaleFactor = 1e-16*cm*cm;
74
76
77 //*******************************************************
78 // Cross section data
79 //*******************************************************
80
81 if(particle == electronDef)
82 {
83 G4String particleName = particle->GetParticleName();
84
85 // MPietrzak, adding paths for N2
87 particleName,
88 "dna/sigma_elastic_e-_PTB_N2",
89 "dna/sigmadiff_cumulated_elastic_e-_PTB_N2",
90 scaleFactor);
91 SetLowELimit("N2", particleName, 10*eV);
92 SetHighELimit("N2", particleName, 1.02*MeV);
93 // MPietrzak
94
96 particleName,
97 "dna/sigma_elastic_e-_PTB_THF",
98 "dna/sigmadiff_cumulated_elastic_e-_PTB_THF",
99 scaleFactor);
100 SetLowELimit("THF", particleName, 10*eV);
101 SetHighELimit("THF", particleName, 1*keV);
102
104 particleName,
105 "dna/sigma_elastic_e-_PTB_PY",
106 "dna/sigmadiff_cumulated_elastic_e-_PTB_PY",
107 scaleFactor);
108 SetLowELimit("PY", particleName, 10*eV);
109 SetHighELimit("PY", particleName, 1*keV);
110
112 particleName,
113 "dna/sigma_elastic_e-_PTB_PU",
114 "dna/sigmadiff_cumulated_elastic_e-_PTB_PU",
115 scaleFactor);
116 SetLowELimit("PU", particleName, 10*eV);
117 SetHighELimit("PU", particleName, 1*keV);
118
120 particleName,
121 "dna/sigma_elastic_e-_PTB_TMP",
122 "dna/sigmadiff_cumulated_elastic_e-_PTB_TMP",
123 scaleFactor);
124 SetLowELimit("TMP", particleName, 10*eV);
125 SetHighELimit("TMP", particleName, 1*keV);
126
127 AddCrossSectionData("G4_WATER",
128 particleName,
129 "dna/sigma_elastic_e_champion",
130 "dna/sigmadiff_cumulated_elastic_e_champion",
131 scaleFactor);
132 SetLowELimit("G4_WATER", particleName, 10*eV);
133 SetHighELimit("G4_WATER", particleName, 1*keV);
134
135 // DNA materials
136 //
137 AddCrossSectionData("backbone_THF",
138 particleName,
139 "dna/sigma_elastic_e-_PTB_THF",
140 "dna/sigmadiff_cumulated_elastic_e-_PTB_THF",
141 scaleFactor*33./30);
142 SetLowELimit("backbone_THF", particleName, 10*eV);
143 SetHighELimit("backbone_THF", particleName, 1*keV);
144
145 AddCrossSectionData("cytosine_PY",
146 particleName,
147 "dna/sigma_elastic_e-_PTB_PY",
148 "dna/sigmadiff_cumulated_elastic_e-_PTB_PY",
149 scaleFactor*42./30);
150 SetLowELimit("cytosine_PY", particleName, 10*eV);
151 SetHighELimit("cytosine_PY", particleName, 1*keV);
152
153 AddCrossSectionData("thymine_PY",
154 particleName,
155 "dna/sigma_elastic_e-_PTB_PY",
156 "dna/sigmadiff_cumulated_elastic_e-_PTB_PY",
157 scaleFactor*48./30);
158 SetLowELimit("thymine_PY", particleName, 10*eV);
159 SetHighELimit("thymine_PY", particleName, 1*keV);
160
161 AddCrossSectionData("adenine_PU",
162 particleName,
163 "dna/sigma_elastic_e-_PTB_PU",
164 "dna/sigmadiff_cumulated_elastic_e-_PTB_PU",
165 scaleFactor*50./44);
166 SetLowELimit("adenine_PU", particleName, 10*eV);
167 SetHighELimit("adenine_PU", particleName, 1*keV);
168
169 AddCrossSectionData("guanine_PU",
170 particleName,
171 "dna/sigma_elastic_e-_PTB_PU",
172 "dna/sigmadiff_cumulated_elastic_e-_PTB_PU",
173 scaleFactor*56./44);
174 SetLowELimit("guanine_PU", particleName, 10*eV);
175 SetHighELimit("guanine_PU", particleName, 1*keV);
176
177 AddCrossSectionData("backbone_TMP",
178 particleName,
179 "dna/sigma_elastic_e-_PTB_TMP",
180 "dna/sigmadiff_cumulated_elastic_e-_PTB_TMP",
181 scaleFactor*33./50);
182 SetLowELimit("backbone_TMP", particleName, 10*eV);
183 SetHighELimit("backbone_TMP", particleName, 1*keV);
184 }
185
186 //*******************************************************
187 // Load the data
188 //*******************************************************
189
191
192 //*******************************************************
193 // Verbose output
194 //*******************************************************
195
196 if (verboseLevel > 2)
197 G4cout << "Loaded cross section files for PTB Elastic model" << G4endl;
198
199 if( verboseLevel>0 )
200 {
201 G4cout << "PTB Elastic model is initialized " << G4endl;
202 }
203}
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
void SetHighELimit(const G4String &material, const G4String &particle, G4double lim)
SetHighEnergyLimit.
Definition: G4VDNAModel.hh:169
void AddCrossSectionData(G4String materialName, G4String particleName, G4String fileCS, G4String fileDiffCS, G4double scaleFactor)
AddCrossSectionData Method used during the initialization of the model class to add a new material....
Definition: G4VDNAModel.cc:58
void SetLowELimit(const G4String &material, const G4String &particle, G4double lim)
SetLowEnergyLimit.
Definition: G4VDNAModel.hh:177
void LoadCrossSectionData(const G4String &particleName)
LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresp...
Definition: G4VDNAModel.cc:75

◆ SampleSecondaries()

void G4DNAPTBElasticModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4String materialName,
const G4DynamicParticle aDynamicElectron,
G4ParticleChangeForGamma particleChangeForGamma,
G4double  tmin,
G4double  tmax 
)
virtual

SampleSecondaries Method called after CrossSectionPerVolume if the process is the one which is selected (according to the sampling on the calculated path length). Here, the characteristics of the incident and created (if any) particle(s) are set (energy, momentum ...).

Parameters
materialName
particleChangeForGamma
tmin
tmax

Implements G4VDNAModel.

Definition at line 356 of file G4DNAPTBElasticModel.cc.

363{
364 if (verboseLevel > 3)
365 G4cout << "Calling SampleSecondaries() of G4DNAPTBElasticModel" << G4endl;
366
367 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
368
369 const G4String& particleName = aDynamicElectron->GetParticleDefinition()->GetParticleName();
370
371 // set killBelowEnergy value for material
372 fKillBelowEnergy = GetLowELimit(materialName, particleName);
373
374 // If the particle (electron here) energy is below the kill limit then we remove it from the simulation
375 if (electronEnergy0 < fKillBelowEnergy)
376 {
377 particleChangeForGamma->SetProposedKineticEnergy(0.);
378 particleChangeForGamma->ProposeTrackStatus(fStopAndKill);
379 particleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
380 }
381 // If we are above the kill limite and below the high limit then we proceed
382 else if (electronEnergy0>= fKillBelowEnergy && electronEnergy0 < GetHighELimit(materialName, particleName) )
383 {
384 // Random sampling of the cosTheta
385 G4double cosTheta = RandomizeCosTheta(electronEnergy0, materialName);
386
387 // Random sampling of phi
388 G4double phi = 2. * pi * G4UniformRand();
389
390 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
391 G4ThreeVector xVers = zVers.orthogonal();
392 G4ThreeVector yVers = zVers.cross(xVers);
393
394 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
395 G4double yDir = xDir;
396 xDir *= std::cos(phi);
397 yDir *= std::sin(phi);
398
399 // Particle direction after ModelInterface
400 G4ThreeVector zPrikeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
401
402 // Give the new direction
403 particleChangeForGamma->ProposeMomentumDirection(zPrikeVers.unit()) ;
404
405 // Update the energy which does not change here
406 particleChangeForGamma->SetProposedKineticEnergy(electronEnergy0);
407 }
408}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
const G4ThreeVector & GetMomentumDirection() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
const G4double pi

The documentation for this class was generated from the following files: