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

#include <G4DNACPA100IonisationModel.hh>

+ Inheritance diagram for G4DNACPA100IonisationModel:

Public Member Functions

 G4DNACPA100IonisationModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="DNACPA100IonisationModel")
 
 ~G4DNACPA100IonisationModel () override=default
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 Initialise Each model must implement an Initialize method.
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 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.
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 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.
 
G4double DifferentialCrossSection (PartKineticInMat info, const G4double &energyTransfer)
 
void SelectFasterComputation (G4bool input)
 
void SelectUseDcs (G4bool input)
 
void SelectStationary (G4bool input)
 
G4DNACPA100IonisationModeloperator= (const G4DNACPA100IonisationModel &right)=delete
 
 G4DNACPA100IonisationModel (const G4DNACPA100IonisationModel &)=delete
 
void ReadDiffCSFile (const std::size_t &materialID, const G4ParticleDefinition *p, const G4String &file, const G4double &scaleFactor) override
 
- 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

Definition at line 56 of file G4DNACPA100IonisationModel.hh.

Constructor & Destructor Documentation

◆ G4DNACPA100IonisationModel() [1/2]

G4DNACPA100IonisationModel::G4DNACPA100IonisationModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNACPA100IonisationModel" )
explicit

Definition at line 62 of file G4DNACPA100IonisationModel.cc.

64 : G4VDNAModel(nam, "all")
65{
66 fpGuanine = G4Material::GetMaterial("G4_GUANINE", false);
67 fpG4_WATER = G4Material::GetMaterial("G4_WATER", false);
68 fpDeoxyribose = G4Material::GetMaterial("G4_DEOXYRIBOSE", false);
69 fpCytosine = G4Material::GetMaterial("G4_CYTOSINE", false);
70 fpThymine = G4Material::GetMaterial("G4_THYMINE", false);
71 fpAdenine = G4Material::GetMaterial("G4_ADENINE", false);
72 fpPhosphate = G4Material::GetMaterial("G4_PHOSPHORIC_ACID", false);
73 fpParticle = G4Electron::ElectronDefinition();
74}
static G4Electron * ElectronDefinition()
Definition G4Electron.cc:86
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
G4VDNAModel(const G4String &nam, const G4String &applyToMaterial="")
G4VDNAModel Constructeur of the G4VDNAModel class.

◆ ~G4DNACPA100IonisationModel()

G4DNACPA100IonisationModel::~G4DNACPA100IonisationModel ( )
overridedefault

◆ G4DNACPA100IonisationModel() [2/2]

G4DNACPA100IonisationModel::G4DNACPA100IonisationModel ( const G4DNACPA100IonisationModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

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

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.

Parameters
material

param materialName

Parameters
p

param ekin

Parameters
emin
emax
Returns
crossSection*numberOfMoleculesPerVolumeUnit

Implements G4VDNAModel.

Definition at line 205 of file G4DNACPA100IonisationModel.cc.

208{
209 // initialise the cross section value (output value)
210 G4double sigma(0);
211
212 // Get the current particle name
213 const G4String& particleName = p->GetParticleName();
214
215 if (p != fpParticle) {
216 G4Exception("G4DNACPA100IonisationModel::CrossSectionPerVolume", "em00223", FatalException,
217 "No model is registered for this particle");
218 }
219
220 auto matID = material->GetIndex();
221
222 // Set the low and high energy limits
223 G4double lowLim = fpModelData->GetLowELimit(matID, p);
224 G4double highLim = fpModelData->GetHighELimit(matID, p);
225
226 // Check that we are in the correct energy range
227 if (ekin >= lowLim && ekin < highLim) {
228 // Get the map with all the model data tables
229 auto tableData = fpModelData->GetData();
230
231 if ((*tableData)[matID][p] == nullptr) {
232 G4Exception("G4DNACPA100IonisationModel::CrossSectionPerVolume", "em00236", FatalException,
233 "No model is registered");
234 }
235 else {
236 sigma = (*tableData)[matID][p]->FindValue(ekin);
237 }
238
239 if (verboseLevel > 2) {
240 auto MolDensity =
242 G4cout << "__________________________________" << G4endl;
243 G4cout << "°°° G4DNACPA100IonisationModel - XS INFO START" << G4endl;
244 G4cout << "°°° Kinetic energy(eV)=" << ekin / eV << " particle : " << particleName << G4endl;
245 G4cout << "°°° lowLim (eV) = " << lowLim / eV << " highLim (eV) : " << highLim / eV << G4endl;
246 G4cout << "°°° Materials = " << (*G4Material::GetMaterialTable())[matID]->GetName() << G4endl;
247 G4cout << "°°° Cross section per " << matID << " index molecule (cm^2)=" << sigma / cm / cm
248 << G4endl;
249 G4cout << "°°° Cross section per Phosphate molecule (cm^-1)="
250 << sigma * MolDensity / (1. / cm) << G4endl;
251 G4cout << "°°° G4DNACPA100IonisationModel - XS INFO END" << G4endl;
252 }
253 }
254
255 auto MolDensity = (*G4DNAMolecularMaterial::Instance()->GetNumMolPerVolTableFor(material))[matID];
256 return sigma * MolDensity;
257}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
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()
std::size_t GetIndex() const
static G4MaterialTable * GetMaterialTable()
const G4String & GetParticleName() const
G4String GetName()
GetName.
G4double GetHighELimit(const size_t &materialID, const G4ParticleDefinition *particle)
GetHighEnergyLimit.
G4double GetLowELimit(const size_t &materialID, const G4ParticleDefinition *particle)
GetLowEnergyLimit.
MaterialParticleMapData * GetData()
GetTableData.

◆ DifferentialCrossSection()

G4double G4DNACPA100IonisationModel::DifferentialCrossSection ( PartKineticInMat info,
const G4double & energyTransfer )

Definition at line 465 of file G4DNACPA100IonisationModel.cc.

467{
468 auto MatID = std::get<0>(info);
469 auto k = std::get<1>(info) / eV; // in eV unit
470 auto shell = std::get<2>(info);
471 G4double sigma = 0.;
472 G4double shellEnergy = iStructure.IonisationEnergy(shell, MatID);
473 G4double kSE(energyTransfer - shellEnergy);
474
475 if (energyTransfer >= shellEnergy) {
476 G4double valueT1 = 0;
477 G4double valueT2 = 0;
478 G4double valueE21 = 0;
479 G4double valueE22 = 0;
480 G4double valueE12 = 0;
481 G4double valueE11 = 0;
482
483 G4double xs11 = 0;
484 G4double xs12 = 0;
485 G4double xs21 = 0;
486 G4double xs22 = 0;
487
488 auto t2 = std::upper_bound(fTMapWithVec[MatID][fpParticle].begin(),
489 fTMapWithVec[MatID][fpParticle].end(), k);
490 auto t1 = t2 - 1;
491
492 if (kSE <= fEMapWithVector[MatID][fpParticle][(*t1)].back()
493 && kSE <= fEMapWithVector[MatID][fpParticle][(*t2)].back())
494 {
495 auto e12 = std::upper_bound(fEMapWithVector[MatID][fpParticle][(*t1)].begin(),
496 fEMapWithVector[MatID][fpParticle][(*t1)].end(), kSE);
497 auto e11 = e12 - 1;
498
499 auto e22 = std::upper_bound(fEMapWithVector[MatID][fpParticle][(*t2)].begin(),
500 fEMapWithVector[MatID][fpParticle][(*t2)].end(), kSE);
501 auto e21 = e22 - 1;
502
503 valueT1 = *t1;
504 valueT2 = *t2;
505 valueE21 = *e21;
506 valueE22 = *e22;
507 valueE12 = *e12;
508 valueE11 = *e11;
509
510 xs11 = diffCrossSectionData[MatID][fpParticle][shell][valueT1][valueE11];
511 xs12 = diffCrossSectionData[MatID][fpParticle][shell][valueT1][valueE12];
512 xs21 = diffCrossSectionData[MatID][fpParticle][shell][valueT2][valueE21];
513 xs22 = diffCrossSectionData[MatID][fpParticle][shell][valueT2][valueE22];
514 }
515
516 G4double xsProduct = xs11 * xs12 * xs21 * xs22;
517
518 if (xsProduct != 0.) {
519 sigma = QuadInterpolator(valueE11, valueE12, valueE21, valueE22, xs11, xs12, xs21, xs22,
520 valueT1, valueT2, k, kSE);
521 }
522 }
523
524 return sigma;
525}
G4double IonisationEnergy(const std::size_t &level, const std::size_t &MatID)

◆ Initialise()

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

Initialise Each model must implement an Initialize method.

Parameters
particle
cuts

Implements G4VDNAModel.

Definition at line 78 of file G4DNACPA100IonisationModel.cc.

80{
81 if (isInitialised) {
82 return;
83 }
84 if (verboseLevel > 3) {
85 G4cout << "Calling G4DNACPA100IonisationModel::Initialise()" << G4endl;
86 }
87
89 if (p != fpParticle) {
90 std::ostringstream oss;
91 oss << " Model is not applied for this particle " << p->GetParticleName();
92 G4Exception("G4DNACPA100IonisationModel::G4DNACPA100IonisationModel", "CPA001",
93 FatalException, oss.str().c_str());
94 }
95
96 const char* path = G4FindDataDir("G4LEDATA");
97
98 if (path == nullptr) {
99 G4Exception("G4DNACPA100IonisationModel::Initialise", "em0006", FatalException,
100 "G4LEDATA environment variable not set.");
101 return;
102 }
103
104 std::size_t index;
105 if (fpG4_WATER != nullptr) {
106 index = fpG4_WATER->GetIndex();
107 G4String eFullFileName = "";
108 fasterCode ? eFullFileName = "/dna/sigmadiff_cumulated_ionisation_e_cpa100_rel"
109 : eFullFileName = "/dna/sigmadiff_ionisation_e_cpa100_rel";
110 AddCrossSectionData(index, p, "dna/sigma_ionisation_e_cpa100_form_rel", eFullFileName,
111 1.e-20 * m * m);
112 SetLowELimit(index, p, 11 * eV);
113 SetHighELimit(index, p, 255955 * eV);
114 }
115 if (fpGuanine != nullptr) {
116 index = fpGuanine->GetIndex();
117 G4String eFullFileName = "";
118 if(useDcs) {
119 fasterCode ? eFullFileName = "/dna/sigmadiff_cumulated_elastic_e_cpa100_guanine"
120 : eFullFileName = "/dna/sigmadiff_ionisation_e_cpa100_guanine";
121 }
122 AddCrossSectionData(index, p, "dna/sigma_ionisation_e_cpa100_guanine", eFullFileName,
123 1. * cm * cm);
124 SetLowELimit(index, p, 11 * eV);
125 SetHighELimit(index, p, 1 * MeV);
126 }
127 if (fpDeoxyribose != nullptr) {
128 index = fpDeoxyribose->GetIndex();
129 G4String eFullFileName = "";
130 if(useDcs) {
131 eFullFileName = "/dna/sigmadiff_cumulated_ionisation_e_cpa100_deoxyribose";
132 }
133 AddCrossSectionData(index, p, "dna/sigma_ionisation_e_cpa100_deoxyribose", eFullFileName,
134 1. * cm * cm);
135 SetLowELimit(index, p, 11 * eV);
136 SetHighELimit(index, p, 1 * MeV);
137 }
138 if (fpCytosine != nullptr) {
139 index = fpCytosine->GetIndex();
140 G4String eFullFileName = "";
141 if(useDcs) {
142 fasterCode ? eFullFileName = "/dna/sigmadiff_cumulated_ionisation_e_cpa100_cytosine"
143 : eFullFileName = "/dna/sigmadiff_ionisation_e_cpa100_cytosine";
144 }
145 AddCrossSectionData(index, p, "dna/sigma_ionisation_e_cpa100_cytosine", eFullFileName,
146 1. * cm * cm);
147 SetLowELimit(index, p, 11 * eV);
148 SetHighELimit(index, p, 1 * MeV);
149 }
150 if (fpThymine != nullptr) {
151 index = fpThymine->GetIndex();
152 G4String eFullFileName = "";
153 if(useDcs) {
154 fasterCode ? eFullFileName = "/dna/sigmadiff_cumulated_ionisation_e_cpa100_thymine"
155 : eFullFileName = "/dna/sigmadiff_ionisation_e_cpa100_thymine";
156 }
157 AddCrossSectionData(index, p, "dna/sigma_ionisation_e_cpa100_thymine", eFullFileName,
158 1. * cm * cm);
159 SetLowELimit(index, p, 11 * eV);
160 SetHighELimit(index, p, 1 * MeV);
161 }
162 if (fpAdenine != nullptr) {
163 index = fpAdenine->GetIndex();
164 G4String eFullFileName = "";
165 if(useDcs) {
166 fasterCode ? eFullFileName = "/dna/sigmadiff_cumulated_ionisation_e_cpa100_adenine"
167 : eFullFileName = "/dna/sigmadiff_ionisation_e_cpa100_adenine";
168 }
169 AddCrossSectionData(index, p, "dna/sigma_ionisation_e_cpa100_adenine", eFullFileName,
170 1. * cm * cm);
171 SetLowELimit(index, p, 11 * eV);
172 SetHighELimit(index, p, 1 * MeV);
173 }
174 if (fpPhosphate != nullptr) {
175 index = fpPhosphate->GetIndex();
176 G4String eFullFileName = "";
177 if(useDcs) {
178 eFullFileName = "dna/sigmadiff_cumulated_ionisation_e_cpa100_phosphoric_acid";
179 }
180 AddCrossSectionData(index, p, "dna/sigma_ionisation_e_cpa100_phosphoric_acid",eFullFileName,
181 1. * cm * cm);
182 SetLowELimit(index, p, 11 * eV);
183 SetHighELimit(index, p, 1 * MeV);
184 }
187 fpModelData = this;
188 }
189 else {
190 auto dataModel = dynamic_cast<G4DNACPA100IonisationModel*>(
192 if (dataModel == nullptr) {
193 G4cout << "G4DNACPA100IonisationModel::CrossSectionPerVolume:: not good modelData" << G4endl;
194 throw;
195 }
196 fpModelData = dataModel;
197 }
198
199 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
200
202 isInitialised = true;
203}
const char * G4FindDataDir(const char *)
G4ParticleChangeForGamma * fParticleChangeForGamma
static G4DNAMaterialManager * Instance()
void SetMasterDataModel(const DNAModelType &t, G4VEmModel *m)
G4VEmModel * GetModel(const DNAModelType &t)
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
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....
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4bool IsLocked() const

◆ operator=()

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

◆ ReadDiffCSFile()

void G4DNACPA100IonisationModel::ReadDiffCSFile ( const std::size_t & materialID,
const G4ParticleDefinition * p,
const G4String & file,
const G4double & scaleFactor )
override

Definition at line 802 of file G4DNACPA100IonisationModel.cc.

805{
806 const char* path = G4FindDataDir("G4LEDATA");
807 if (path == nullptr) {
808 G4Exception("G4DNACPA100IonisationModel::ReadAllDiffCSFiles", "em0006", FatalException,
809 "G4LEDATA environment variable was not set.");
810 return;
811 }
812
813 std::ostringstream fullFileName;
814 fullFileName << path << "/" << file << ".dat";
815
816 std::ifstream diffCrossSection(fullFileName.str().c_str());
817 std::stringstream endPath;
818 if (!diffCrossSection) {
819 endPath << "Missing data file: " << file;
820 G4Exception("G4DNACPA100IonisationModel::Initialise", "em0003", FatalException,
821 endPath.str().c_str());
822 }
823
824 // load data from the file
825 fTMapWithVec[materialID][p].push_back(0.);
826
827 G4String line;
828
829 while (!diffCrossSection.eof()) {
830 G4double T, E;
831 diffCrossSection >> T >> E;
832
833 if (T != fTMapWithVec[materialID][p].back()) {
834 fTMapWithVec[materialID][p].push_back(T);
835 }
836
837 // T is incident energy, E is the energy transferred
838 if (T != fTMapWithVec[materialID][p].back()) {
839 fTMapWithVec[materialID][p].push_back(T);
840 }
841
842 auto eshell = (G4int)iStructure.NumberOfLevels(materialID);
843 for (G4int shell = 0; shell < eshell; ++shell) {
844 diffCrossSection >> diffCrossSectionData[materialID][p][shell][T][E];
845 if (fasterCode) {
846 fEnergySecondaryData[materialID][p][shell][T]
847 [diffCrossSectionData[materialID][p][shell][T][E]] = E;
848
849 fProbaShellMap[materialID][p][shell][T].push_back(
850 diffCrossSectionData[materialID][p][shell][T][E]);
851 }
852 else {
853 diffCrossSectionData[materialID][p][shell][T][E] *= scaleFactor;
854 fEMapWithVector[materialID][p][T].push_back(E);
855 }
856 }
857 }
858}
int G4int
Definition G4Types.hh:85
std::size_t NumberOfLevels(const std::size_t &MatID)

◆ SampleSecondaries()

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

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.

Parameters
materialName

param particleChangeForGamma

Parameters
tmin

param tmax

Implements G4VDNAModel.

Definition at line 261 of file G4DNACPA100IonisationModel.cc.

265{
266 if (verboseLevel > 3) {
267 G4cout << "Calling SampleSecondaries() of G4DNACPA100IonisationModel" << G4endl;
268 }
269 auto k = particle->GetKineticEnergy();
270
271 const G4Material* material = couple->GetMaterial();
272
273 auto MatID = material->GetIndex();
274
275 auto p = particle->GetDefinition();
276
277 auto lowLim = fpModelData->GetLowELimit(MatID, p);
278 auto highLim = fpModelData->GetHighELimit(MatID, p);
279
280 // Check if we are in the correct energy range
281 if (k >= lowLim && k < highLim) {
282 const auto& primaryDirection = particle->GetMomentumDirection();
283 auto particleMass = particle->GetDefinition()->GetPDGMass();
284 auto totalEnergy = k + particleMass;
285 auto pSquare = k * (totalEnergy + particleMass);
286 auto totalMomentum = std::sqrt(pSquare);
287 G4int shell = -1;
288 G4double bindingEnergy, secondaryKinetic;
289 shell = fpModelData->RandomSelectShell(k, p, MatID);
290 bindingEnergy = iStructure.IonisationEnergy(shell, MatID);
291
292 if (k < bindingEnergy) {
293 return;
294 }
295
296 auto info = std::make_tuple(MatID, k, shell);
297
298 secondaryKinetic = -1000 * eV;
299 if (fpG4_WATER->GetIndex() != MatID) {//for DNA material useDcs = false
300 secondaryKinetic = fpModelData->RandomizeEjectedElectronEnergyFromanalytical(info);
301 }else if(fasterCode){
302 secondaryKinetic = fpModelData->RandomizeEjectedElectronEnergyFromCumulatedDcs(info);
303 }else{
304 secondaryKinetic = fpModelData->RandomizeEjectedElectronEnergy(info);
305 }
306
307 G4double cosTheta = 0.;
308 G4double phi = 0.;
309 RandomizeEjectedElectronDirection(particle->GetDefinition(), k, secondaryKinetic, cosTheta,
310 phi);
311
312 G4double sinTheta = std::sqrt(1. - cosTheta * cosTheta);
313 G4double dirX = sinTheta * std::cos(phi);
314 G4double dirY = sinTheta * std::sin(phi);
315 G4double dirZ = cosTheta;
316 G4ThreeVector deltaDirection(dirX, dirY, dirZ);
317 deltaDirection.rotateUz(primaryDirection);
318
319 // SI - For atom. deexc. tagging - 23/05/2017
320 if (secondaryKinetic > 0) {
321 auto dp = new G4DynamicParticle(G4Electron::Electron(), deltaDirection, secondaryKinetic);
322 fvect->push_back(dp);
323 }
324
325 if (particle->GetDefinition() != fpParticle) {
327 }
328 else {
329 G4double deltaTotalMomentum =
330 std::sqrt(secondaryKinetic * (secondaryKinetic + 2. * electron_mass_c2));
331 G4double finalPx =
332 totalMomentum * primaryDirection.x() - deltaTotalMomentum * deltaDirection.x();
333 G4double finalPy =
334 totalMomentum * primaryDirection.y() - deltaTotalMomentum * deltaDirection.y();
335 G4double finalPz =
336 totalMomentum * primaryDirection.z() - deltaTotalMomentum * deltaDirection.z();
337 G4double finalMomentum = std::sqrt(finalPx * finalPx + finalPy * finalPy + finalPz * finalPz);
338 finalPx /= finalMomentum;
339 finalPy /= finalMomentum;
340 finalPz /= finalMomentum;
341
342 G4ThreeVector direction;
343 direction.set(finalPx, finalPy, finalPz);
344
346 }
347
348 // SI - For atom. deexc. tagging - 23/05/2017
349
350 // AM: sample deexcitation
351 // here we assume that H_{2}O electronic levels are the same of Oxigen.
352 // this can be considered true with a rough 10% error in energy on K-shell,
353
354 G4double scatteredEnergy = k - bindingEnergy - secondaryKinetic;
355
356 // SI: only atomic deexcitation from K shell is considered
357 // Hoang: only for water
358 if (fpG4_WATER != nullptr && material == G4Material::GetMaterial("G4_WATER")) {
359 std::size_t secNumberInit = 0; // need to know at a certain point the energy of secondaries
360 std::size_t secNumberFinal = 0; // So I'll make the diference and then sum the energies
361 if ((fAtomDeexcitation != nullptr) && shell == 4) {
362 G4int Z = 8;
363 auto Kshell = fAtomDeexcitation->GetAtomicShell(Z, G4AtomicShellEnumerator(0));
364 secNumberInit = fvect->size();
365 fAtomDeexcitation->GenerateParticles(fvect, Kshell, Z, 0, 0);
366 secNumberFinal = fvect->size();
367 if (secNumberFinal > secNumberInit) {
368 for (std::size_t i = secNumberInit; i < secNumberFinal; ++i) {
369 // Check if there is enough residual energy
370 if (bindingEnergy >= ((*fvect)[i])->GetKineticEnergy()) {
371 // Ok, this is a valid secondary: keep it
372 bindingEnergy -= ((*fvect)[i])->GetKineticEnergy();
373 }
374 else {
375 // Invalid secondary: not enough energy to create it!
376 // Keep its energy in the local deposit
377 delete (*fvect)[i];
378 (*fvect)[i] = nullptr;
379 }
380 }
381 }
382 }
383 }
384
385 // This should never happen
386 if (bindingEnergy < 0.0) {
387 G4Exception("G4DNACPA100IonisatioModel1::SampleSecondaries()", "em2050", FatalException,
388 "Negative local energy deposit");
389 }
390 if (!statCode) {
393 }
394 else {
397 }
398
399 // only water for chemistry
400 if (fpG4_WATER != nullptr && material == G4Material::GetMaterial("G4_WATER")) {
401 const G4Track* theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
403 theIncomingTrack);
404 }
405 }
406}
@ eIonizedMolecule
Hep3Vector unit() const
void set(double x, double y, double z)
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
static G4Electron * Electron()
Definition G4Electron.cc:91
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
G4int RandomSelectShell(const G4double &k, const G4ParticleDefinition *particle, const size_t &materialName)
RandomSelectShell Method to randomely select a shell from the data table uploaded....
const G4Track * GetCurrentTrack() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)

◆ SelectFasterComputation()

void G4DNACPA100IonisationModel::SelectFasterComputation ( G4bool input)
inline

Definition at line 159 of file G4DNACPA100IonisationModel.hh.

160{
161 fasterCode = input;
162}

◆ SelectStationary()

void G4DNACPA100IonisationModel::SelectStationary ( G4bool input)
inline

Definition at line 175 of file G4DNACPA100IonisationModel.hh.

176{
177 statCode = input;
178}

◆ SelectUseDcs()

void G4DNACPA100IonisationModel::SelectUseDcs ( G4bool input)
inline

Definition at line 166 of file G4DNACPA100IonisationModel.hh.

167{
168 useDcs = input;
169}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNACPA100IonisationModel::fParticleChangeForGamma = nullptr
protected

Definition at line 104 of file G4DNACPA100IonisationModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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