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

The G4DNAPTBExcitationModel class This class implements the PTB excitation model. More...

#include <G4DNAPTBExcitationModel.hh>

+ Inheritance diagram for G4DNAPTBExcitationModel:

Public Member Functions

 G4DNAPTBExcitationModel (const G4String &applyToMaterial="all", const G4ParticleDefinition *p=0, const G4String &nam="DNAPTBExcitationModel")
 G4DNAPTBExcitationModel Constructor.
 
virtual ~G4DNAPTBExcitationModel ()
 ~G4DNAPTBExcitationModel Destructor
 
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.
 
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.
 
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 will be called. The method sets the incident particle characteristics after the ModelInterface.
 
- 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 G4DNAPTBExcitationModel class This class implements the PTB excitation model.

Definition at line 50 of file G4DNAPTBExcitationModel.hh.

Constructor & Destructor Documentation

◆ G4DNAPTBExcitationModel()

G4DNAPTBExcitationModel::G4DNAPTBExcitationModel ( const G4String applyToMaterial = "all",
const G4ParticleDefinition p = 0,
const G4String nam = "DNAPTBExcitationModel" 
)

G4DNAPTBExcitationModel Constructor.

Parameters
applyToMaterial
p
nam

Definition at line 36 of file G4DNAPTBExcitationModel.cc.

38 : G4VDNAModel(nam, applyToMaterial)
39{
40 verboseLevel= 0;
41 // Verbosity scale:
42 // 0 = nothing
43 // 1 = warning for energy non-conservation
44 // 2 = details of energy budget
45 // 3 = calculation of cross sections, file openings, sampling of atoms
46 // 4 = entering in methods
47
48 // initialisation of mean energy loss for each material
49 tableMeanEnergyPTB["THF"] = 8.01*eV;
50 tableMeanEnergyPTB["PY"] = 7.61*eV;
51 tableMeanEnergyPTB["PU"] = 7.61*eV;
52 tableMeanEnergyPTB["TMP"] = 8.01*eV;
53
54 if( verboseLevel>0 )
55 {
56 G4cout << "PTB excitation model is constructed " << G4endl;
57 }
58}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
The G4VDNAModel class.
Definition: G4VDNAModel.hh:50

◆ ~G4DNAPTBExcitationModel()

G4DNAPTBExcitationModel::~G4DNAPTBExcitationModel ( )
virtual

~G4DNAPTBExcitationModel Destructor

Definition at line 62 of file G4DNAPTBExcitationModel.cc.

63{
64
65}

Member Function Documentation

◆ CrossSectionPerVolume()

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

CrossSectionPerVolume Retrieve the cross section corresponding to the current material, particle and energy.

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

Implements G4VDNAModel.

Definition at line 186 of file G4DNAPTBExcitationModel.cc.

192{
193 if (verboseLevel > 3)
194 G4cout << "Calling CrossSectionPerVolume() of G4DNAPTBExcitationModel" << G4endl;
195
196 // Get the name of the current particle
197 G4String particleName = particleDefinition->GetParticleName();
198
199 // initialise variables
200 G4double lowLim = 0;
201 G4double highLim = 0;
202 G4double sigma=0;
203
204 // Get the low energy limit for the current particle
205 lowLim = GetLowELimit(materialName, particleName);
206
207 // Get the high energy limit for the current particle
208 highLim = GetHighELimit(materialName, particleName);
209
210 // Check that we are in the correct energy range
211 if (ekin >= lowLim && ekin < highLim)
212 {
213 // Get the map with all the data tables
214 TableMapData* tableData = GetTableData();
215
216 // Retrieve the cross section value
217 sigma = (*tableData)[materialName][particleName]->FindValue(ekin);
218
219 if (verboseLevel > 2)
220 {
221 G4cout << "__________________________________" << G4endl;
222 G4cout << "°°° G4DNAPTBExcitationModel - XS INFO START" << G4endl;
223 G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
224 G4cout << "°°° Cross section per "<< materialName <<" molecule (cm^2)=" << sigma/cm/cm << G4endl;
225 G4cout << "°°° G4DNAPTBExcitationModel - XS INFO END" << G4endl;
226 }
227
228 }
229
230 // Return the cross section value
231 return sigma;
232}
double G4double
Definition: G4Types.hh:83
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

◆ Initialise()

void G4DNAPTBExcitationModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector = *(new G4DataVector()),
G4ParticleChangeForGamma fpChangeForGamme = nullptr 
)
virtual

Initialise Set the materials for which the model can be used and defined the energy limits.

Implements G4VDNAModel.

Definition at line 69 of file G4DNAPTBExcitationModel.cc.

71{
72 if (verboseLevel > 3)
73 G4cout << "Calling G4DNAPTBExcitationModel::Initialise()" << G4endl;
74
75 G4double scaleFactor = 1e-16*cm*cm;
76 G4double scaleFactorBorn = (1.e-22 / 3.343) * m*m;
77
79
80 //*******************************************************
81 // Cross section data
82 //*******************************************************
83
84 if(particle == electronDef)
85 {
86 G4String particleName = particle->GetParticleName();
87
89 particleName,
90 "dna/sigma_excitation_e-_PTB_THF",
91 scaleFactor);
92 SetLowELimit("THF", particleName, 9.*eV);
93 SetHighELimit("THF", particleName, 1.*keV);
94
96 particleName,
97 "dna/sigma_excitation_e-_PTB_PY",
98 scaleFactor);
99 SetLowELimit("PY", particleName, 9.*eV);
100 SetHighELimit("PY", particleName, 1.*keV);
101
103 particleName,
104 "dna/sigma_excitation_e-_PTB_PU",
105 scaleFactor);
106 SetLowELimit("PU", particleName, 9.*eV);
107 SetHighELimit("PU", particleName, 1.*keV);
108
110 particleName,
111 "dna/sigma_excitation_e-_PTB_TMP",
112 scaleFactor);
113 SetLowELimit("TMP", particleName, 9.*eV);
114 SetHighELimit("TMP", particleName, 1.*keV);
115
116 AddCrossSectionData("G4_WATER",
117 particleName,
118 "dna/sigma_excitation_e_born",
119 scaleFactorBorn);
120 SetLowELimit("G4_WATER", particleName, 9.*eV);
121 SetHighELimit("G4_WATER", particleName, 1.*keV);
122
123 // DNA materials
124 //
125 AddCrossSectionData("backbone_THF",
126 particleName,
127 "dna/sigma_excitation_e-_PTB_THF",
128 scaleFactor*33./30);
129 SetLowELimit("backbone_THF", particleName, 9.*eV);
130 SetHighELimit("backbone_THF", particleName, 1.*keV);
131
132 AddCrossSectionData("cytosine_PY",
133 particleName,
134 "dna/sigma_excitation_e-_PTB_PY",
135 scaleFactor*42./30);
136 SetLowELimit("cytosine_PY", particleName, 9.*eV);
137 SetHighELimit("cytosine_PY", particleName, 1.*keV);
138
139 AddCrossSectionData("thymine_PY",
140 particleName,
141 "dna/sigma_excitation_e-_PTB_PY",
142 scaleFactor*48./30);
143 SetLowELimit("thymine_PY", particleName, 9.*eV);
144 SetHighELimit("thymine_PY", particleName, 1.*keV);
145
146 AddCrossSectionData("adenine_PU",
147 particleName,
148 "dna/sigma_excitation_e-_PTB_PU",
149 scaleFactor*50./44);
150 SetLowELimit("adenine_PU", particleName, 9.*eV);
151 SetHighELimit("adenine_PU", particleName, 1.*keV);
152
153 AddCrossSectionData("guanine_PU",
154 particleName,
155 "dna/sigma_excitation_e-_PTB_PU",
156 scaleFactor*56./44);
157 SetLowELimit("guanine_PU", particleName, 9.*eV);
158 SetHighELimit("guanine_PU", particleName, 1.*keV);
159
160 AddCrossSectionData("backbone_TMP",
161 particleName,
162 "dna/sigma_excitation_e-_PTB_TMP",
163 scaleFactor*33./50);
164 SetLowELimit("backbone_TMP", particleName, 9.*eV);
165 SetHighELimit("backbone_TMP", particleName, 1.*keV);
166 }
167
168 //*******************************************************
169 // Load data
170 //*******************************************************
171
173
174 //*******************************************************
175 // Verbose
176 //*******************************************************
177
178 if( verboseLevel>0 )
179 {
180 G4cout << "PTB excitation model is initialized " << G4endl;
181 }
182}
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
const G4String & GetParticleName() const
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 G4DNAPTBExcitationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4String materialName,
const G4DynamicParticle aDynamicParticle,
G4ParticleChangeForGamma particleChangeForGamma,
G4double  tmin,
G4double  tmax 
)
virtual

SampleSecondaries If the model is selected for the ModelInterface then the SampleSecondaries method will be called. The method sets the incident particle characteristics after the ModelInterface.

Parameters
materialName
particleChangeForGamma
tmin
tmax

Implements G4VDNAModel.

Definition at line 236 of file G4DNAPTBExcitationModel.cc.

243{
244 if (verboseLevel > 3)
245 G4cout << "Calling SampleSecondaries() of G4DNAPTBExcitationModel" << G4endl;
246
247 // Get the incident particle kinetic energy
248 G4double k = aDynamicParticle->GetKineticEnergy();
249
250 if(materialName!="G4_WATER")
251 {
252 // Retrieve the excitation energy for the current material
253 G4double excitationEnergy = tableMeanEnergyPTB[materialName];
254
255 // Calculate the new energy of the particle
256 G4double newEnergy = k - excitationEnergy;
257
258 // Check that the new energy is above zero before applying it the particle.
259 // Otherwise, do nothing.
260 if (newEnergy > 0)
261 {
262 particleChangeForGamma->ProposeMomentumDirection(aDynamicParticle->GetMomentumDirection());
263 particleChangeForGamma->SetProposedKineticEnergy(newEnergy);
264 particleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy);
265 }
266 }
267 else
268 {
269 const G4String& particleName = aDynamicParticle->GetDefinition()->GetParticleName();
270
271 G4int level = RandomSelectShell(k,particleName, materialName);
272 G4double excitationEnergy = waterStructure.ExcitationEnergy(level);
273 G4double newEnergy = k - excitationEnergy;
274
275 if (newEnergy > 0)
276 {
277 particleChangeForGamma->ProposeMomentumDirection(aDynamicParticle->GetMomentumDirection());
278 particleChangeForGamma->SetProposedKineticEnergy(newEnergy);
279 particleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy);
280 }
281
282 const G4Track * theIncomingTrack = particleChangeForGamma->GetCurrentTrack();
284 level,
285 theIncomingTrack);
286 }
287}
@ eExcitedMolecule
int G4int
Definition: G4Types.hh:85
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4int RandomSelectShell(G4double k, const G4String &particle, const G4String &materialName)
RandomSelectShell Method to randomely select a shell from the data table uploaded....
Definition: G4VDNAModel.cc:182
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

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