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

#include <G4DNAPTBExcitationModel.hh>

+ Inheritance diagram for G4DNAPTBExcitationModel:

Public Types

using MapMeanEnergy = std::map<std::size_t, G4double>
 

Public Member Functions

 G4DNAPTBExcitationModel (const G4String &applyToMaterial="all", const G4ParticleDefinition *p=nullptr, const G4String &nam="DNAPTBExcitationModel")
 G4DNAPTBExcitationModel Constructor.
 
 ~G4DNAPTBExcitationModel () override=default
 ~G4DNAPTBExcitationModel Destructor
 
 G4DNAPTBExcitationModel (const G4DNAPTBExcitationModel &)=delete
 
G4DNAPTBExcitationModeloperator= (const G4DNAPTBExcitationModel &right)=delete
 
void Initialise (const G4ParticleDefinition *particle, const G4DataVector &) override
 Initialise Set the materials for which the model can be used and defined the energy limits.
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 CrossSectionPerVolume Retrieve the cross section corresponding to the current material, particle and energy.
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax) override
 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.
 
 ~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
 

Public Attributes

G4ParticleChangeForGammafParticleChangeForGamma = nullptr
 

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 *)
 
- 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
 

Detailed Description

Definition at line 50 of file G4DNAPTBExcitationModel.hh.

Member Typedef Documentation

◆ MapMeanEnergy

using G4DNAPTBExcitationModel::MapMeanEnergy = std::map<std::size_t, G4double>

Definition at line 53 of file G4DNAPTBExcitationModel.hh.

Constructor & Destructor Documentation

◆ G4DNAPTBExcitationModel() [1/2]

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

G4DNAPTBExcitationModel Constructor.

Parameters
applyToMaterial
p
nam

Definition at line 38 of file G4DNAPTBExcitationModel.cc.

40 : G4VDNAModel(nam, applyToMaterial)
41{
42 fpTHF = G4Material::GetMaterial("THF", false);
43 fpPY = G4Material::GetMaterial("PY", false);
44 fpPU = G4Material::GetMaterial("PU", false);
45 fpTMP = G4Material::GetMaterial("TMP", false);
46 fpG4_WATER = G4Material::GetMaterial("G4_WATER", false);
47 fpBackbone_THF = G4Material::GetMaterial("backbone_THF", false);
48 fpCytosine_PY = G4Material::GetMaterial("cytosine_PY", false);
49 fpThymine_PY = G4Material::GetMaterial("thymine_PY", false);
50 fpAdenine_PU = G4Material::GetMaterial("adenine_PU", false);
51 fpBackbone_TMP = G4Material::GetMaterial("backbone_TMP", false);
52 fpGuanine_PU = G4Material::GetMaterial("guanine_PU", false);
53 fpN2 = G4Material::GetMaterial("N2", false);
54 // initialisation of mean energy loss for each material
55
56 if (fpTHF != nullptr) {
57 fTableMeanEnergyPTB[fpTHF->GetIndex()] = 8.01 * eV;
58 }
59
60 if (fpPY != nullptr) {
61 fTableMeanEnergyPTB[fpPY->GetIndex()] = 7.61 * eV;
62 }
63
64 if (fpPU != nullptr) {
65 fTableMeanEnergyPTB[fpPU->GetIndex()] = 7.61 * eV;
66 }
67 if (fpTMP != nullptr) {
68 fTableMeanEnergyPTB[fpTMP->GetIndex()] = 8.01 * eV;
69 }
70
71 if (verboseLevel > 0) {
72 G4cout << "PTB excitation model is constructed " << G4endl;
73 }
74}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
std::size_t GetIndex() const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
G4VDNAModel(const G4String &nam, const G4String &applyToMaterial="")
G4VDNAModel Constructeur of the G4VDNAModel class.

◆ ~G4DNAPTBExcitationModel()

G4DNAPTBExcitationModel::~G4DNAPTBExcitationModel ( )
overridedefault

~G4DNAPTBExcitationModel Destructor

◆ G4DNAPTBExcitationModel() [2/2]

G4DNAPTBExcitationModel::G4DNAPTBExcitationModel ( const G4DNAPTBExcitationModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

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

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 204 of file G4DNAPTBExcitationModel.cc.

208{
209 // Get the name of the current particle
210 G4String particleName = p->GetParticleName();
211 const std::size_t& MatID = material->GetIndex();
212 // initialise variables
213 G4double lowLim;
214 G4double highLim;
215 G4double sigma = 0;
216
217 // Get the low energy limit for the current particle
218 lowLim = fpModelData->GetLowELimit(MatID, p);
219
220 // Get the high energy limit for the current particle
221 highLim = fpModelData->GetHighELimit(MatID, p);
222
223 // Check that we are in the correct energy range
224 if (ekin >= lowLim && ekin < highLim) {
225 // Get the map with all the data tables
226 auto Data = fpModelData->GetData();
227
228 if ((*Data)[MatID][p] == nullptr) {
229 G4Exception("G4DNAPTBExcitationModel::CrossSectionPerVolume", "em00236", FatalException,
230 "No model is registered");
231 }
232 // Retrieve the cross section value
233 sigma = (*Data)[MatID][p]->FindValue(ekin);
234
235 if (verboseLevel > 2) {
236 G4cout << "__________________________________" << G4endl;
237 G4cout << "°°° G4DNAPTBExcitationModel - XS INFO START" << G4endl;
238 G4cout << "°°° Kinetic energy(eV)=" << ekin / eV << " particle : " << particleName << G4endl;
239 G4cout << "°°° Cross section per " << MatID << " ID molecule (cm^2)=" << sigma / cm / cm
240 << G4endl;
241 G4cout << "°°° G4DNAPTBExcitationModel - XS INFO END" << G4endl;
242 }
243 }
244
245 // Return the cross section value
246 auto MolDensity =
248 return sigma * MolDensity;
249}
@ 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.

◆ Initialise()

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

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

Implements G4VDNAModel.

Definition at line 78 of file G4DNAPTBExcitationModel.cc.

80{
81 if (isInitialised) {
82 return;
83 }
84 if (verboseLevel > 3)
85 {
86 G4cout << "Calling G4DNAPTBExcitationModel::Initialise()" << G4endl;
87 }
88
89 if (particle != G4Electron::ElectronDefinition()) {
90 std::ostringstream oss;
91 oss << " Model is not applied for this particle " << particle->GetParticleName();
92 G4Exception("G4DNAPTBExcitationModel::Initialise", "PTB001", FatalException, oss.str().c_str());
93 }
94
95 G4double scaleFactor = 1e-16 * cm * cm;
96 G4double scaleFactorBorn = (1.e-22 / 3.343) * m * m;
97
98 //*******************************************************
99 // Cross section data
100 //*******************************************************
101 std::size_t index;
102 if (fpTHF != nullptr) {
103 index = fpTHF->GetIndex();
104 AddCrossSectionData(index, particle, "dna/sigma_excitation_e-_PTB_THF", scaleFactor);
105 SetLowELimit(index, particle, 9. * eV);
106 SetHighELimit(index, particle, 1. * keV);
107 }
108 if (fpPY != nullptr) {
109 index = fpPY->GetIndex();
110 AddCrossSectionData(index, particle, "dna/sigma_excitation_e-_PTB_PY", scaleFactor);
111 SetLowELimit(index, particle, 9. * eV);
112 SetHighELimit(index, particle, 1. * keV);
113 }
114
115 if (fpPU != nullptr) {
116 index = fpPU->GetIndex();
117 AddCrossSectionData(index, particle, "dna/sigma_excitation_e-_PTB_PU", scaleFactor);
118 SetLowELimit(index, particle, 9. * eV);
119 SetHighELimit(index, particle, 1. * keV);
120 }
121
122 if (fpTMP != nullptr) {
123 index = fpTMP->GetIndex();
124 AddCrossSectionData(index, particle, "dna/sigma_excitation_e-_PTB_TMP", scaleFactor);
125 SetLowELimit(index, particle, 9. * eV);
126 SetHighELimit(index, particle, 1. * keV);
127 }
128 if (fpG4_WATER != nullptr) {
129 index = fpG4_WATER->GetIndex();
130 AddCrossSectionData(index, particle, "dna/sigma_excitation_e_born", scaleFactorBorn);
131 SetLowELimit(index, particle, 9. * eV);
132 SetHighELimit(index, particle, 1. * keV);
133 }
134 // DNA materials
135 //
136 if (fpBackbone_THF != nullptr) {
137 index = fpBackbone_THF->GetIndex();
138 AddCrossSectionData(index, particle, "dna/sigma_excitation_e-_PTB_THF", scaleFactor * 33. / 30);
139 SetLowELimit(index, particle, 9. * eV);
140 SetHighELimit(index, particle, 1. * keV);
141 }
142 if (fpCytosine_PY != nullptr) {
143 index = fpCytosine_PY->GetIndex();
144 AddCrossSectionData(index, particle, "dna/sigma_excitation_e-_PTB_PY", scaleFactor * 42. / 30);
145 SetLowELimit(index, particle, 9. * eV);
146 SetHighELimit(index, particle, 1. * keV);
147 }
148 if (fpThymine_PY != nullptr) {
149 index = fpThymine_PY->GetIndex();
150 AddCrossSectionData(index, particle, "dna/sigma_excitation_e-_PTB_PY", scaleFactor * 48. / 30);
151 SetLowELimit(index, particle, 9. * eV);
152 SetHighELimit(index, particle, 1. * keV);
153 }
154 if (fpAdenine_PU != nullptr) {
155 index = fpAdenine_PU->GetIndex();
156 AddCrossSectionData(index, particle, "dna/sigma_excitation_e-_PTB_PU", scaleFactor * 50. / 44);
157 SetLowELimit(index, particle, 9. * eV);
158 SetHighELimit(index, particle, 1. * keV);
159 }
160 if (fpGuanine_PU != nullptr) {
161 index = fpGuanine_PU->GetIndex();
162 AddCrossSectionData(index, particle, "dna/sigma_excitation_e-_PTB_PU", scaleFactor * 56. / 44);
163 SetLowELimit(index, particle, 9. * eV);
164 SetHighELimit(index, particle, 1. * keV);
165 }
166 if (fpBackbone_TMP != nullptr) {
167 index = fpBackbone_TMP->GetIndex();
168 AddCrossSectionData(index, particle, "dna/sigma_excitation_e-_PTB_TMP", scaleFactor * 33. / 50);
169 SetLowELimit(index, particle, 9. * eV);
170 SetHighELimit(index, particle, 1. * keV);
171 }
172 // MPietrzak, adding paths for N2
173 if (fpN2 != nullptr) {
174 index = fpN2->GetIndex();
175 AddCrossSectionData(index, particle, "dna/sigma_excitation_e-_PTB_N2", scaleFactor);
176 SetLowELimit(index, particle, 13. * eV);
177 SetHighELimit(index, particle, 1.02 * MeV);
178 }
180 // Load data
181 LoadCrossSectionData(particle);
183 fpModelData = this;
184 }
185 else {
186 auto dataModel = dynamic_cast<G4DNAPTBExcitationModel*>(
188 if (dataModel == nullptr) {
189 G4cout << "G4DNAPTBExcitationModel::Initialise:: not good modelData" << G4endl;
190 G4Exception("G4DNAPTBExcitationModel::Initialise", "PTB0006", FatalException,
191 "not good modelData");
192 }
193 else {
194 fpModelData = dataModel;
195 }
196 }
197
199 isInitialised = true;
200}
static G4DNAMaterialManager * Instance()
void SetMasterDataModel(const DNAModelType &t, G4VEmModel *m)
G4VEmModel * GetModel(const DNAModelType &t)
G4ParticleChangeForGamma * fParticleChangeForGamma
static G4Electron * ElectronDefinition()
Definition G4Electron.cc:86
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=()

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

◆ SampleSecondaries()

void G4DNAPTBExcitationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * fvect,
const G4MaterialCutsCouple * couple,
const G4DynamicParticle * aDynamicParticle,
G4double tmin,
G4double tmax )
overridevirtual

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

param particleChangeForGamma

Parameters
tmin

param tmax

Implements G4VDNAModel.

Definition at line 253 of file G4DNAPTBExcitationModel.cc.

257{
258 const std::size_t& materialID = (std::size_t)couple->GetIndex();
259
260 // Get the incident particle kinetic energy
261 G4double k = aDynamicParticle->GetKineticEnergy();
262 // Get the particle name
263 const auto& particle = aDynamicParticle->GetDefinition();
264 // Get the energy limits
265 G4double lowLim = fpModelData->GetLowELimit(materialID, particle);
266 G4double highLim = fpModelData->GetHighELimit(materialID, particle);
267
268 // Check if we are in the correct energy range
269 if (k >= lowLim && k < highLim) {
270 if (fpN2 != nullptr && materialID == fpN2->GetIndex()) {
271 // Retrieve the excitation energy for the current material
272 G4int level = fpModelData->RandomSelectShell(k, particle, materialID);
273 G4double excitationEnergy = ptbExcitationStructure.ExcitationEnergy(level, fpN2->GetIndex());
274
275 // Calculate the new energy of the particle
276 G4double newEnergy = k - excitationEnergy;
277
278 // Check that the new energy is above zero before applying it the particle.
279 // Otherwise, do nothing.
280 if (newEnergy > 0) {
284 G4double ioniThres = ptbIonisationStructure.IonisationEnergy(0, fpN2->GetIndex());
285 // if excitation energy greater than ionisation threshold, then autoionisaiton
286 if ((excitationEnergy > ioniThres) && (G4UniformRand() < 0.5)) {
288 // energy of ejected electron
289 G4double secondaryKinetic = excitationEnergy - ioniThres;
290 // random direction
291 G4double cosTheta = 2 * G4UniformRand() - 1., phi = CLHEP::twopi * G4UniformRand();
292 G4double sinTheta = std::sqrt(1. - cosTheta * cosTheta);
293 G4double ux = sinTheta * std::cos(phi), uy = sinTheta * std::sin(phi), uz = cosTheta;
294 G4ThreeVector deltaDirection(ux, uy, uz);
295 // Create the new particle with its characteristics
296 auto dp = new G4DynamicParticle(G4Electron::Electron(), deltaDirection, secondaryKinetic);
297 fvect->push_back(dp);
298 }
299 }
300 else {
301 G4ExceptionDescription description;
302 description << "Kinetic energy <= 0 at " << fpN2->GetName() << " material !!!";
303 G4Exception("G4DNAPTBExcitationModel::SampleSecondaries", "", FatalException, description);
304 }
305 }
306 else if (fpG4_WATER == nullptr || materialID != fpG4_WATER->GetIndex()) {
307 // Retrieve the excitation energy for the current material
308 G4double excitationEnergy = fTableMeanEnergyPTB[materialID];
309 // Calculate the new energy of the particle
310 G4double newEnergy = k - excitationEnergy;
311 // Check that the new energy is above zero before applying it the particle.
312 // Otherwise, do nothing.
313 if (newEnergy > 0) {
317 }
318 else {
319 G4ExceptionDescription description;
320 description << "Kinetic energy <= 0 at " << materialID << " index material !!!";
321 G4Exception("G4DNAPTBExcitationModel::SampleSecondaries", "", FatalException, description);
322 }
323 }
324 else {
325 G4int level = RandomSelectShell(k, particle, materialID);
326 G4double excitationEnergy = waterStructure.ExcitationEnergy(level);
327 G4double newEnergy = k - excitationEnergy;
328
329 if (newEnergy > 0) {
333 const G4Track* theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
335 theIncomingTrack);
336 }
337 else {
338 G4ExceptionDescription description;
339 description << "Kinetic energy <= 0 at " << materialID << " ID material !!!";
340 G4Exception("G4DNAPTBExcitationModel::SampleSecondaries", "", FatalException, description);
341 }
342 }
343 }
344}
@ eExcitedMolecule
std::ostringstream G4ExceptionDescription
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
G4double ExcitationEnergy(const G4int &ExcLevel, const size_t &materialID)
G4double IonisationEnergy(G4int level, const size_t &materialName)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition G4Electron.cc:91
const G4String & GetName() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
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)

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAPTBExcitationModel::fParticleChangeForGamma = nullptr

Definition at line 104 of file G4DNAPTBExcitationModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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