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

#include <G4DNATripleIonisationModel.hh>

+ Inheritance diagram for G4DNATripleIonisationModel:

Public Member Functions

 G4DNATripleIonisationModel (const G4ParticleDefinition *p=nullptr, const G4String &model_name="G4DNATripleIonisationModel")
 
 ~G4DNATripleIonisationModel () override=default
 
void Initialise (const G4ParticleDefinition *particle, const G4DataVector &) override
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *pdef, G4double ekin, G4double, G4double) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *vsec, const G4MaterialCutsCouple *couple, const G4DynamicParticle *particle, G4double, G4double) override
 
- Public Member Functions inherited from G4DNADoubleIonisationModel
 G4DNADoubleIonisationModel (const G4ParticleDefinition *p=nullptr, const G4String &model_name="G4DNADoubleIonisationModel")
 
 ~G4DNADoubleIonisationModel () override
 
G4DNADoubleIonisationModeloperator= (const G4DNADoubleIonisationModel &)=delete
 
 G4DNADoubleIonisationModel (const G4DNADoubleIonisationModel &)=delete
 
void SelectStationary (G4bool in)
 
void SelectVerboseLevel (G4int in)
 
void UseChampionAlphaParameter (G4bool in)
 
void SetMultipleIonisationEnergy (G4double in)
 
- 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)
 
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)
 
void SetMasterThread (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4DNADoubleIonisationModel
G4double RandomizeEjectedElectronEnergy (G4ParticleDefinition *pdef, G4double ekin, G4int shell)
 
G4int RandomSelect (G4double energy, G4double scale_param, const G4String &pname)
 
G4double GenerateSecondaries (std::vector< G4DynamicParticle * > *vsec, const G4MaterialCutsCouple *couple, const G4DynamicParticle *particle, G4int ioni_shell, G4double &theta, G4double &phi, G4double &shell_energy)
 
G4double GetLowEnergyLimit (const G4String &pname)
 
G4double GetUppEnergyLimit (const G4String &pname)
 
- 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 G4DNADoubleIonisationModel
G4ParticleChangeForGammaparticle_change_ = nullptr
 
G4bool stat_code_
 
G4VAtomDeexcitationatom_deex_
 
EnergyLimitTable elow_tab_
 
EnergyLimitTable eupp_tab_
 
CrossSectionDataTable xs_tab_
 
G4ParticleDefinitionproton_def_ {nullptr}
 
G4ParticleDefinitionalpha_def_ {nullptr}
 
G4ParticleDefinitioncarbon_def_ {nullptr}
 
const std::vector< G4double > * water_density_
 
G4bool is_initialized_
 
G4int verbose_level_
 
std::map< G4double, G4doublemodel_elow_tab_
 
G4DNAMultipleIonisationManagermioni_manager_ {nullptr}
 
G4bool use_champion_param_
 
G4double energy_threshold_
 
- 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 40 of file G4DNATripleIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4DNATripleIonisationModel()

G4DNATripleIonisationModel::G4DNATripleIonisationModel ( const G4ParticleDefinition * p = nullptr,
const G4String & model_name = "G4DNATripleIonisationModel" )

Definition at line 58 of file G4DNATripleIonisationModel.cc.

60 : G4DNADoubleIonisationModel(p, model_name)
61{
62 // Triple-ionisation energy
63 energy_threshold_ = 65.0 * eV;
64}
G4DNADoubleIonisationModel(const G4ParticleDefinition *p=nullptr, const G4String &model_name="G4DNADoubleIonisationModel")

◆ ~G4DNATripleIonisationModel()

G4DNATripleIonisationModel::~G4DNATripleIonisationModel ( )
overridedefault

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4DNATripleIonisationModel::CrossSectionPerVolume ( const G4Material * material,
const G4ParticleDefinition * pdef,
G4double ekin,
G4double ,
G4double  )
overridevirtual

Reimplemented from G4DNADoubleIonisationModel.

Definition at line 186 of file G4DNATripleIonisationModel.cc.

189{
190
191 if (verbose_level_ > 3) {
192 G4cout << "Calling G4DNATripleIonisationModel::CrossSectionPerVolume()"
193 << G4endl;
194 }
195
196 // Calculate total cross section for model
197
198 if (pdef != proton_def_ && pdef != alpha_def_ && pdef != carbon_def_) {
199 return 0.0;
200 }
201
202 static G4double water_dens = (*water_density_)[material->GetIndex()];
203
204 const auto& pname = pdef->GetParticleName();
205
206 const auto low_energy_lim = GetLowEnergyLimit(pname);
207 const auto upp_energy_lim = GetUppEnergyLimit(pname);
208
209 G4double sigma{0.0};
210 if (ekin <= upp_energy_lim) {
211
212 if (ekin < low_energy_lim) { ekin = low_energy_lim; }
213
214 CrossSectionDataTable::iterator pos = xs_tab_.find(pname);
215 if (pos == xs_tab_.end()) {
216 G4Exception("G4DNATripleIonisationModel::CrossSectionPerVolume",
217 "em0002", FatalException,
218 "Model not applicable to particle type.");
219 }
220
221 G4DNACrossSectionDataSet* table = pos->second;
222 if (table != nullptr) {
223 const auto a = mioni_manager_->GetAlphaParam(ekin);
224 sigma = table->FindValue(ekin) * a * a;
225 }
226
227 }
228
229 if (verbose_level_ > 2) {
230
231 std::stringstream msg;
232
233 msg << "----------------------------------------------------------------\n";
234 msg << " G4DNATripleIonisationModel - XS INFO START\n";
235 msg << " - Kinetic energy(eV): " << ekin/eV << ", Particle : "
236 << pdef->GetParticleName() << "\n";
237 msg << " - Cross section per water molecule (cm^2): "
238 << sigma / cm / cm << "\n";
239 msg << " - Cross section per water molecule (cm^-1): "
240 << sigma * water_dens / (1.0 / cm) << "\n";
241 msg << " G4DNATripleIonisationModel - XS INFO END\n";
242 msg << "----------------------------------------------------------------\n";
243
244 G4cout << msg.str() << G4endl;
245
246 }
247
248 return (sigma * water_dens);
249
250}
@ 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
G4double FindValue(G4double e, G4int componentId=0) const override
G4DNAMultipleIonisationManager * mioni_manager_
G4double GetUppEnergyLimit(const G4String &pname)
G4double GetLowEnergyLimit(const G4String &pname)
std::size_t GetIndex() const
const G4String & GetParticleName() const

◆ Initialise()

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

Reimplemented from G4DNADoubleIonisationModel.

Definition at line 67 of file G4DNATripleIonisationModel.cc.

69{
70 if (verbose_level_ > 3) {
71 G4cout << "Calling G4DNATripleIonisationModel::Initialise()" << G4endl;
72 }
73
77
78 constexpr G4double kScaleFactor = 1.0 * m * m;
79
80 mioni_manager_ = new G4DNAMultipleIonisationManager();
81
82 G4double Z{0.0}, A{0.0};
83 G4String alpha_param_file{"dna/multipleionisation_alphaparam_champion.dat"};
84
85 if (particle == proton_def_) {
86
87 // *************************************************************************
88 // for protons
89 const auto& proton = proton_def_->GetParticleName();
91 eupp_tab_[proton] = 3.0 * MeV;
92
93 // load cross-section data for single ionization process
94 auto xs_proton = new G4DNACrossSectionDataSet(
95 new G4LogLogInterpolation, eV, kScaleFactor);
96 xs_proton->LoadData("dna/sigma_ionisation_p_rudd");
97 xs_tab_[proton] = xs_proton;
98
99 // set energy limits
102
103 if (!use_champion_param_) {
104 alpha_param_file = "dna/multipleionisation_alphaparam_p.dat";
105 }
106
107 Z = static_cast<G4double>(proton_def_->GetAtomicNumber());
108 A = static_cast<G4double>(proton_def_->GetAtomicMass());
109
110 } else if (particle == alpha_def_) {
111
112 //**************************************************************************
113 // for alpha particles
114 const auto& alpha = alpha_def_->GetParticleName();
116 eupp_tab_[alpha] = 23.0 * MeV;
117
118 // load cross-section data for single ionization process
119 auto xs_alpha = new G4DNACrossSectionDataSet(
120 new G4LogLogInterpolation, eV, kScaleFactor);
121 xs_alpha->LoadData("dna/sigma_ionisation_alphaplusplus_rudd");
122 xs_tab_[alpha] = xs_alpha;
123
124 // set energy limits
127
128 if (!use_champion_param_) {
129 alpha_param_file = "dna/multipleionisation_alphaparam_alphaplusplus.dat";
130 }
131
132 Z = static_cast<G4double>(alpha_def_->GetAtomicNumber());
133 A = static_cast<G4double>(alpha_def_->GetAtomicMass());
134
135 } else if (particle == G4GenericIon::GenericIonDefinition()) {
136
137 // *************************************************************************
138 // for carbon ions
139 const auto& carbon = carbon_def_->GetParticleName();
140 elow_tab_[carbon] = model_elow_tab_[5] * carbon_def_->GetAtomicMass();
141 eupp_tab_[carbon] = 120.0 * MeV;
142
143 // load cross-section data for single ionization process
144 auto xs_carbon = new G4DNACrossSectionDataSet(
145 new G4LogLogInterpolation, eV, kScaleFactor);
146 xs_carbon->LoadData("dna/sigma_ionisation_c_rudd");
147 xs_tab_[carbon] = xs_carbon;
148
149 // set energy limits
152
153 if (!use_champion_param_) {
154 alpha_param_file = "dna/multipleionisation_alphaparam_c.dat";
155 }
156
157 Z = static_cast<G4double>(carbon_def_->GetAtomicNumber());
158 A = static_cast<G4double>(carbon_def_->GetAtomicMass());
159
160 }
161
162 // load alpha parameter
163 mioni_manager_->LoadAlphaParam(alpha_param_file, Z, A);
164
165 if (verbose_level_ > 0) {
166 G4cout << "G4DNATripleIonisationModel is initialized " << G4endl
167 << "Energy range: "
168 << LowEnergyLimit() / eV << " eV - "
169 << HighEnergyLimit() / keV << " keV for "
170 << particle->GetParticleName()
171 << G4endl;
172 }
173
175 G4Material::GetMaterial("G4_WATER"));
176
178
179 if (is_initialized_) { return; }
180
182 is_initialized_ = true;
183}
const G4double A[17]
std::map< G4double, G4double > model_elow_tab_
const std::vector< G4double > * water_density_
G4ParticleChangeForGamma * particle_change_
static G4DNAGenericIonsManager * Instance()
G4ParticleDefinition * GetIon(const G4String &name)
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()
static G4GenericIon * GenericIonDefinition()
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
static G4IonTable * GetIonTable()
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
static G4Proton * ProtonDefinition()
Definition G4Proton.cc:85
void SetHighEnergyLimit(G4double)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
void SetLowEnergyLimit(G4double)

◆ SampleSecondaries()

void G4DNATripleIonisationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * vsec,
const G4MaterialCutsCouple * couple,
const G4DynamicParticle * particle,
G4double ,
G4double  )
overridevirtual

Reimplemented from G4DNADoubleIonisationModel.

Definition at line 253 of file G4DNATripleIonisationModel.cc.

256{
257
258 if (verbose_level_ > 3) {
259 G4cout << "Calling SampleSecondaries() of G4DNATripleIonisationModel"
260 << G4endl;
261 }
262
263 // get the definition for this parent particle
264 auto pdef = particle->GetDefinition();
265
266 // get kinetic energy
267 auto ekin = particle->GetKineticEnergy();
268
269 // get particle name
270 const auto& pname = pdef->GetParticleName();
271
272 // get energy limits
273 const auto low_energy_lim = GetLowEnergyLimit(pname);
274
275 // ***************************************************************************
276 // stop the transportation process of this parent particle
277 // if its kinetic energy is below the lower limit
278 if (ekin < low_energy_lim) {
279 particle_change_->SetProposedKineticEnergy(0.0);
280 particle_change_->ProposeTrackStatus(fStopAndKill);
281 particle_change_->ProposeLocalEnergyDeposit(ekin);
282 return;
283 }
284 // ***************************************************************************
285
286 constexpr G4int kNumSecondaries = 3;
287 constexpr G4double kDeltaTheta = pi * 0.666666667;
288
289 G4int ioni_shell[kNumSecondaries] = {0, 0, 0};
290 G4double shell_energy[kNumSecondaries];
291
292 auto scale_param = mioni_manager_->GetAlphaParam(ekin);
293 scale_param *= scale_param;
294
295 G4bool is_continue{true};
296 while (1) {
297 ioni_shell[0] = RandomSelect(ekin, scale_param, pname);
298 ioni_shell[1] = RandomSelect(ekin, scale_param, pname);
299 ioni_shell[2] = RandomSelect(ekin, scale_param, pname);
300 is_continue = (ioni_shell[0] == ioni_shell[1] &&
301 ioni_shell[1] == ioni_shell[2]);
302 if (!is_continue) { break; }
303 }
304
305 G4double tot_ioni_energy{0.0};
306 for (int i = 0; i < kNumSecondaries; i++) {
307 shell_energy[i] = ::water_structure.IonisationEnergy(ioni_shell[i]);
308 tot_ioni_energy += shell_energy[i];
309 }
310
311 if (ekin < tot_ioni_energy || tot_ioni_energy < energy_threshold_) {
312 return;
313 }
314
315 // generate secondary electrons
316 G4double theta{0.0}, phi{0.0}, tot_ekin2{0.0};
317 for (int i = 0; i < kNumSecondaries; i++) {
318 tot_ekin2 += GenerateSecondaries(vsec, couple, particle, ioni_shell[i],
319 theta, phi, shell_energy[i]);
320 theta += kDeltaTheta;
321 }
322
323 // This should never happen
324 if (mioni_manager_->CheckShellEnergy(eTripleIonisedMolecule, shell_energy)) {
325 G4Exception("G4DNATripleIonisatioModel::SampleSecondaries()",
326 "em2050", FatalException, "Negative local energy deposit");
327 }
328
329 // ***************************************************************************
330 // update kinematics for this parent particle
331 const auto primary_dir = particle->GetMomentumDirection();
332 particle_change_->ProposeMomentumDirection(primary_dir);
333
334 const auto scattered_energy = ekin - tot_ioni_energy - tot_ekin2;
335
336 // update total amount of shell energy
337 tot_ioni_energy = shell_energy[0] + shell_energy[1] + shell_energy[2];
338
339 if (stat_code_) {
340 particle_change_->SetProposedKineticEnergy(ekin);
341 particle_change_->ProposeLocalEnergyDeposit(
342 ekin - scattered_energy);
343 } else {
344 particle_change_->SetProposedKineticEnergy(scattered_energy);
345 particle_change_->ProposeLocalEnergyDeposit(tot_ioni_energy);
346 }
347
348 // ***************************************************************************
349 // generate triple-ionized water molecules (H2O^3+)
350 const auto the_track = particle_change_->GetCurrentTrack();
351 mioni_manager_->CreateMultipleIonisedWaterMolecule(
352 eTripleIonisedMolecule, ioni_shell, the_track);
353 // ***************************************************************************
354
355}
@ fStopAndKill
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4int RandomSelect(G4double energy, G4double scale_param, const G4String &pname)
G4double GenerateSecondaries(std::vector< G4DynamicParticle * > *vsec, const G4MaterialCutsCouple *couple, const G4DynamicParticle *particle, G4int ioni_shell, G4double &theta, G4double &phi, G4double &shell_energy)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4double pi

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