Geant4 11.2.2
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 Types

using TriDimensionMap
 
using VecMap
 

Public Member Functions

 G4DNAPTBElasticModel (const G4String &applyToMaterial="all", const G4ParticleDefinition *p=nullptr, const G4String &nam="DNAPTBElasticModel")
 G4DNAPTBElasticModel Constructor.
 
 ~G4DNAPTBElasticModel () override=default
 ~G4DNAPTBElasticModel Destructor
 
 G4DNAPTBElasticModel (G4DNAPTBElasticModel &)=delete
 
G4DNAPTBElasticModeloperator= (const G4DNAPTBElasticModel &right)=delete
 
void Initialise (const G4ParticleDefinition *particle, const G4DataVector &) override
 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.
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 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.
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax) override
 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.
 
 ~G4VDNAModel () override
 ~G4VDNAModel
 
G4bool IsMaterialDefine (const size_t &materialID)
 IsMaterialDefine Check if the given material is defined in the simulation.
 
G4bool IsParticleExistingInModelForMaterial (const G4ParticleDefinition *particleName, const size_t &materialID)
 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 size_t &materialID, const G4ParticleDefinition *particle)
 GetHighEnergyLimit.
 
G4double GetLowELimit (const size_t &materialID, const G4ParticleDefinition *particle)
 GetLowEnergyLimit.
 
void SetHighELimit (const size_t &materialID, const G4ParticleDefinition *particle, G4double lim)
 SetHighEnergyLimit.
 
void SetLowELimit (const size_t &materialID, const G4ParticleDefinition *particle, G4double lim)
 SetLowEnergyLimit.
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma = nullptr
 
- Protected Attributes inherited from G4VDNAModel
G4bool isInitialised = false
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Additional Inherited Members

- Protected Types inherited from G4VDNAModel
using MaterialParticleMapData
 
- Protected Member Functions inherited from G4VDNAModel
MaterialParticleMapDataGetData ()
 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 size_t &materialID, const G4ParticleDefinition *p, const G4String &file, const G4double &scaleFactor)
 ReadAndSaveCSFile Read and save a "simple" cross section file : use of G4DNACrossSectionDataSet->loadData()
 
G4int RandomSelectShell (const G4double &k, const G4ParticleDefinition *particle, const size_t &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 (const size_t &materialName, const G4ParticleDefinition *particleName, const G4String &fileCS, const G4String &fileDiffCS, const 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 (const size_t &materialName, const G4ParticleDefinition *particleName, const G4String &fileCS, const 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 G4ParticleDefinition *particleName)
 LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresponding data.
 
virtual void ReadDiffCSFile (const size_t &materialName, const G4ParticleDefinition *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 size_t &materialID, const G4ParticleDefinition *p)
 EnableMaterialAndParticle.
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

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.

Member Typedef Documentation

◆ TriDimensionMap

Initial value:
std::map<std::size_t,
std::map<const G4ParticleDefinition*, std::map<G4double, std::map<G4double, G4double>>>>

Definition at line 51 of file G4DNAPTBElasticModel.hh.

◆ VecMap

Initial value:
std::map<std::size_t,
std::map<const G4ParticleDefinition*, std::map<G4double, std::vector<G4double>>>>

Definition at line 53 of file G4DNAPTBElasticModel.hh.

Constructor & Destructor Documentation

◆ G4DNAPTBElasticModel() [1/2]

G4DNAPTBElasticModel::G4DNAPTBElasticModel ( const G4String & applyToMaterial = "all",
const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNAPTBElasticModel" )

G4DNAPTBElasticModel Constructor.

Parameters
applyToMaterial
p
nam

Definition at line 39 of file G4DNAPTBElasticModel.cc.

41 : G4VDNAModel(nam, applyToMaterial)
42{
43 if (verboseLevel > 0) {
44 G4cout << "PTB Elastic model is constructed : " << G4endl;
45 }
46 fpTHF = G4Material::GetMaterial("THF", false);
47 fpPY = G4Material::GetMaterial("PY", false);
48 fpPU = G4Material::GetMaterial("PU", false);
49 fpTMP = G4Material::GetMaterial("TMP", false);
50 fpG4_WATER = G4Material::GetMaterial("G4_WATER", false);
51 fpBackbone_THF = G4Material::GetMaterial("backbone_THF", false);
52 fpCytosine_PY = G4Material::GetMaterial("cytosine_PY", false);
53 fpThymine_PY = G4Material::GetMaterial("thymine_PY", false);
54 fpAdenine_PU = G4Material::GetMaterial("adenine_PU", false);
55 fpBackbone_TMP = G4Material::GetMaterial("backbone_TMP", false);
56 fpGuanine_PU = G4Material::GetMaterial("guanine_PU", false);
57 fpN2 = G4Material::GetMaterial("N2", false);
58}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
G4VDNAModel(const G4String &nam, const G4String &applyToMaterial="")
G4VDNAModel Constructeur of the G4VDNAModel class.

◆ ~G4DNAPTBElasticModel()

G4DNAPTBElasticModel::~G4DNAPTBElasticModel ( )
overridedefault

~G4DNAPTBElasticModel Destructor

◆ G4DNAPTBElasticModel() [2/2]

G4DNAPTBElasticModel::G4DNAPTBElasticModel ( G4DNAPTBElasticModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4DNAPTBElasticModel::CrossSectionPerVolume ( const G4Material * material,
const G4ParticleDefinition * p,
G4double ekin,
G4double emin,
G4double emax )
overridevirtual

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 306 of file G4DNAPTBElasticModel.cc.

309{
310 if (verboseLevel > 3){
311 G4cout << "Calling CrossSectionPerVolume() of G4DNAPTBElasticModel" << G4endl;
312 }
313
314 // Get the name of the current particle
315 const G4String& particleName = p->GetParticleName();
316 const std::size_t& materialID = pMaterial->GetIndex();
317
318 // set killBelowEnergy value for current material
319 fKillBelowEnergy = fpModelData->GetLowELimit(materialID, p);
320 // initialise the return value (cross section) to zero
321 G4double sigma = 0.;
322
323 // check if we are below the high energy limit
324 if (ekin < fpModelData->GetHighELimit(materialID, p)) {
325 // This is used to kill the particle if its kinetic energy is below fKillBelowEnergy.
326 // If the energy is lower then we return a maximum cross section and thus the SampleSecondaries
327 // method will be called for sure. SampleSecondaries will remove the particle from the
328 // simulation.
329 //
330 // SI : XS must not be zero otherwise sampling of secondaries method ignored
331 if (ekin < fKillBelowEnergy) {
332 return DBL_MAX;
333 }
334
335 // Get the tables with the cross section data
336 auto tableData = fpModelData->GetData();
337 if ((*tableData)[materialID][p] == nullptr) {
338 G4Exception("G4DNAPTBElasticModel::CrossSectionPerVolume", "em00236", FatalException,
339 "No model is registered");
340 }
341 // Retrieve the cross section value
342 sigma = (*tableData)[materialID][p]->FindValue(ekin);
343 }
344
345 if (verboseLevel > 2) {
346 G4cout << "__________________________________" << G4endl;
347 G4cout << "°°° G4DNAPTBElasticModel - XS INFO START" << G4endl;
348 G4cout << "°°° Kinetic energy(eV)=" << ekin / eV << " particle : " << particleName << G4endl;
349 G4cout << "°°° Cross section per molecule (cm^2)=" << sigma / cm / cm << G4endl;
350 G4cout << "°°° G4DNAPTBElasticModel - XS INFO END" << G4endl;
351 }
352
353 // Return the cross section
354 auto MolDensity =
356 return sigma * MolDensity;
357}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
const std::vector< G4double > * GetNumMolPerVolTableFor(const G4Material *) const
Retrieve a table of molecular densities (number of molecules per unit volume) in the G4 unit system f...
static G4DNAMolecularMaterial * Instance()
const G4String & GetParticleName() const
G4double GetHighELimit(const size_t &materialID, const G4ParticleDefinition *particle)
GetHighEnergyLimit.
G4double GetLowELimit(const size_t &materialID, const G4ParticleDefinition *particle)
GetLowEnergyLimit.
MaterialParticleMapData * GetData()
GetTableData.
#define DBL_MAX
Definition templates.hh:62

◆ Initialise()

void G4DNAPTBElasticModel::Initialise ( const G4ParticleDefinition * particle,
const G4DataVector &  )
overridevirtual

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 62 of file G4DNAPTBElasticModel.cc.

64{
65 if (isInitialised) {
66 return;
67 }
68 if (verboseLevel > 3)
69 {
70 G4cout << "Calling G4DNAPTBElasticModel::Initialise()" << G4endl;
71 }
72 if (particle != G4Electron::ElectronDefinition()) {
73 std::ostringstream oss;
74 oss << " Model is not applied for this particle " << particle->GetParticleName();
75 G4Exception("G4DNAPTBElasticModel::G4DNAPTBElasticModel", "PTB001", FatalException,
76 oss.str().c_str());
77 }
78 G4double scaleFactor = 1e-16 * cm * cm;
79 //*******************************************************
80 // Cross section data
81 //*******************************************************
82
83 std::size_t index;
84 // MPietrzak, adding paths for N2
85 if (fpN2 != nullptr) {
86 index = fpN2->GetIndex();
87 AddCrossSectionData(index, particle, "dna/sigma_elastic_e-_PTB_N2",
88 "dna/sigmadiff_cumulated_elastic_e-_PTB_N2", scaleFactor);
89 SetLowELimit(index, particle, 10 * eV);
90 SetHighELimit(index, particle, 1.02 * MeV);
91 }
92 // MPietrzak
93
94 if (fpTHF != nullptr) {
95 index = fpTHF->GetIndex();
96 AddCrossSectionData(index, particle, "dna/sigma_elastic_e-_PTB_THF",
97 "dna/sigmadiff_cumulated_elastic_e-_PTB_THF", scaleFactor);
98 SetLowELimit(index, particle, 10 * eV);
99 SetHighELimit(index, particle, 1 * keV);
100 }
101
102 if (fpPY != nullptr) {
103 index = fpPY->GetIndex();
104 AddCrossSectionData(index, particle, "dna/sigma_elastic_e-_PTB_PY",
105 "dna/sigmadiff_cumulated_elastic_e-_PTB_PY", scaleFactor);
106 SetLowELimit(index, particle, 10 * eV);
107 SetHighELimit(index, particle, 1 * keV);
108 }
109
110 if (fpPU != nullptr) {
111 index = fpPU->GetIndex();
112 AddCrossSectionData(index, particle, "dna/sigma_elastic_e-_PTB_PU",
113 "dna/sigmadiff_cumulated_elastic_e-_PTB_PU", scaleFactor);
114 SetLowELimit(index, particle, 10 * eV);
115 SetHighELimit(index, particle, 1 * keV);
116 }
117
118 if (fpTMP != nullptr) {
119 index = fpTMP->GetIndex();
120 AddCrossSectionData(index, particle, "dna/sigma_elastic_e-_PTB_TMP",
121 "dna/sigmadiff_cumulated_elastic_e-_PTB_TMP", scaleFactor);
122 SetLowELimit(index, particle, 10 * eV);
123 SetHighELimit(index, particle, 1 * keV);
124 }
125 //????
126 if (fpG4_WATER != nullptr) {
127 index = fpG4_WATER->GetIndex();
128 AddCrossSectionData(index, particle, "dna/sigma_elastic_e_champion",
129 "dna/sigmadiff_cumulated_elastic_e_champion", scaleFactor);
130 SetLowELimit(index, particle, 10 * eV);
131 SetHighELimit(index, particle, 1 * keV);
132 }
133 // DNA materials
134 //
135 if (fpBackbone_THF != nullptr) {
136 index = fpBackbone_THF->GetIndex();
137 AddCrossSectionData(index, particle, "dna/sigma_elastic_e-_PTB_THF",
138 "dna/sigmadiff_cumulated_elastic_e-_PTB_THF", scaleFactor * 33. / 30);
139 SetLowELimit(index, particle, 10 * eV);
140 SetHighELimit(index, particle, 1 * keV);
141 }
142
143 if (fpCytosine_PY != nullptr) {
144 index = fpCytosine_PY->GetIndex();
145 AddCrossSectionData(index, particle, "dna/sigma_elastic_e-_PTB_PY",
146 "dna/sigmadiff_cumulated_elastic_e-_PTB_PY", scaleFactor * 42. / 30);
147 SetLowELimit(index, particle, 10 * eV);
148 SetHighELimit(index, particle, 1 * keV);
149 }
150
151 if (fpThymine_PY != nullptr) {
152 index = fpThymine_PY->GetIndex();
153 AddCrossSectionData(index, particle, "dna/sigma_elastic_e-_PTB_PY",
154 "dna/sigmadiff_cumulated_elastic_e-_PTB_PY", scaleFactor * 48. / 30);
155 SetLowELimit(index, particle, 10 * eV);
156 SetHighELimit(index, particle, 1 * keV);
157 }
158
159 if (fpAdenine_PU != nullptr) {
160 index = fpAdenine_PU->GetIndex();
161 AddCrossSectionData(index, particle, "dna/sigma_elastic_e-_PTB_PU",
162 "dna/sigmadiff_cumulated_elastic_e-_PTB_PU", scaleFactor * 50. / 44);
163 SetLowELimit(index, particle, 10 * eV);
164 SetHighELimit(index, particle, 1 * keV);
165 }
166 if (fpGuanine_PU != nullptr) {
167 index = fpGuanine_PU->GetIndex();
168 AddCrossSectionData(index, particle, "dna/sigma_elastic_e-_PTB_PU",
169 "dna/sigmadiff_cumulated_elastic_e-_PTB_PU", scaleFactor * 56. / 44);
170 SetLowELimit(index, particle, 10 * eV);
171 SetHighELimit(index, particle, 1 * keV);
172 }
173
174 if (fpBackbone_TMP != nullptr) {
175 index = fpBackbone_TMP->GetIndex();
176 AddCrossSectionData(index, particle, "dna/sigma_elastic_e-_PTB_TMP",
177 "dna/sigmadiff_cumulated_elastic_e-_PTB_TMP", scaleFactor * 33. / 50);
178 SetLowELimit(index, particle, 10 * eV);
179 SetHighELimit(index, particle, 1 * keV);
180 }
181
183 // Load the data
184 LoadCrossSectionData(particle);
186 fpModelData = this;
187 }
188 else {
189 auto dataModel = dynamic_cast<G4DNAPTBElasticModel*>(
191 if (dataModel == nullptr) {
192 G4cout << "G4DNAPTBElasticModel::Initialise:: not good modelData" << G4endl;
193 G4Exception("G4DNAPTBElasticModel::Initialise", "PTB0006", FatalException,
194 "not good modelData");
195 }
196 else {
197 fpModelData = dataModel;
198 }
199 }
200
201 if (verboseLevel > 2) {
202 G4cout << "Loaded cross section files for PTB Elastic model" << G4endl;
203 }
204
206 isInitialised = true;
207}
static G4DNAMaterialManager * Instance()
void SetMasterDataModel(const DNAModelType &t, G4VEmModel *m)
G4VEmModel * GetModel(const DNAModelType &t)
The G4DNAPTBElasticModel class This class implements the elastic model for the DNA materials and prec...
G4ParticleChangeForGamma * fParticleChangeForGamma
static G4Electron * ElectronDefinition()
Definition G4Electron.cc:86
std::size_t GetIndex() const
void LoadCrossSectionData(const G4ParticleDefinition *particleName)
LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresp...
void SetLowELimit(const size_t &materialID, const G4ParticleDefinition *particle, G4double lim)
SetLowEnergyLimit.
void SetHighELimit(const size_t &materialID, const G4ParticleDefinition *particle, G4double lim)
SetHighEnergyLimit.
void AddCrossSectionData(const size_t &materialName, const G4ParticleDefinition *particleName, const G4String &fileCS, const G4String &fileDiffCS, const G4double &scaleFactor)
AddCrossSectionData Method used during the initialization of the model class to add a new material....
G4bool isInitialised
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4bool IsLocked() const

◆ operator=()

G4DNAPTBElasticModel & G4DNAPTBElasticModel::operator= ( const G4DNAPTBElasticModel & right)
delete

◆ SampleSecondaries()

void G4DNAPTBElasticModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * ,
const G4MaterialCutsCouple * couple,
const G4DynamicParticle * aDynamicElectron,
G4double tmin,
G4double tmax )
overridevirtual

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 361 of file G4DNAPTBElasticModel.cc.

365{
366 if (verboseLevel > 3) {
367 G4cout << "Calling SampleSecondaries() of G4DNAPTBElasticModel" << G4endl;
368 }
369
370 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
371 const std::size_t& materialID = couple->GetIndex();
372 auto p = aDynamicElectron->GetParticleDefinition();
373
374 // set killBelowEnergy value for material
375 fKillBelowEnergy = fpModelData->GetLowELimit(materialID, p);
376
377 // If the particle (electron here) energy is below the kill limit then we remove it from the
378 // simulation
379 if (electronEnergy0 < fKillBelowEnergy) {
383 }
384 // If we are above the kill limite and below the high limit then we proceed
385 else if (electronEnergy0 >= fKillBelowEnergy && electronEnergy0 < GetHighELimit(materialID, p)) {
386 // Random sampling of the cosTheta
387 G4double cosTheta = fpModelData->RandomizeCosTheta(electronEnergy0, materialID);
388
389 // Random sampling of phi
390 G4double phi = 2. * CLHEP::pi * G4UniformRand();
391
392 auto zVers = aDynamicElectron->GetMomentumDirection();
393 auto xVers = zVers.orthogonal();
394 auto yVers = zVers.cross(xVers);
395
396 G4double xDir = std::sqrt(1. - cosTheta * cosTheta);
397 G4double yDir = xDir;
398 xDir *= std::cos(phi);
399 yDir *= std::sin(phi);
400
401 // Particle direction after ModelInterface
402 G4ThreeVector zPrikeVers((xDir * xVers + yDir * yVers + cosTheta * zVers));
403
404 // Give the new direction
406
407 // Update the energy which does not change here
409 }
410}
@ 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)

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAPTBElasticModel::fParticleChangeForGamma = nullptr
protected

Definition at line 114 of file G4DNAPTBElasticModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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