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

#include <G4DNAQuadrupleIonisationModel.hh>

+ Inheritance diagram for G4DNAQuadrupleIonisationModel:

Public Member Functions

 G4DNAQuadrupleIonisationModel (const G4ParticleDefinition *p=nullptr, const G4String &model_name="G4DNAQuadrupleIonisationModel")
 
 ~G4DNAQuadrupleIonisationModel () 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 G4DNAQuadrupleIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4DNAQuadrupleIonisationModel()

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

Definition at line 60 of file G4DNAQuadrupleIonisationModel.cc.

62 : G4DNADoubleIonisationModel(p, model_name)
63{
64 // Quadruple-ionisation energy
65 energy_threshold_ = 88.0 * eV;
66}
G4DNADoubleIonisationModel(const G4ParticleDefinition *p=nullptr, const G4String &model_name="G4DNADoubleIonisationModel")

◆ ~G4DNAQuadrupleIonisationModel()

G4DNAQuadrupleIonisationModel::~G4DNAQuadrupleIonisationModel ( )
overridedefault

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4DNADoubleIonisationModel.

Definition at line 188 of file G4DNAQuadrupleIonisationModel.cc.

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

Reimplemented from G4DNADoubleIonisationModel.

Definition at line 69 of file G4DNAQuadrupleIonisationModel.cc.

71{
72 if (verbose_level_ > 3) {
73 G4cout << "Calling G4DNAQuadrupleIonisationModel::Initialise()" << G4endl;
74 }
75
79
80 constexpr G4double kScaleFactor = 1.0 * m * m;
81
82 mioni_manager_ = new G4DNAMultipleIonisationManager();
83
84 G4double Z{0.0}, A{0.0};
85 G4String alpha_param_file{"dna/multipleionisation_alphaparam_champion.dat"};
86
87 if (particle == proton_def_) {
88
89 // *************************************************************************
90 // for protons
91 const auto& proton = proton_def_->GetParticleName();
93 eupp_tab_[proton] = 3.0 * MeV;
94
95 // load cross-section data for single ionization process
96 auto xs_proton = new G4DNACrossSectionDataSet(
97 new G4LogLogInterpolation, eV, kScaleFactor);
98 xs_proton->LoadData("dna/sigma_ionisation_p_rudd");
99 xs_tab_[proton] = xs_proton;
100
101 // set energy limits
104
105 if (!use_champion_param_) {
106 alpha_param_file = "dna/multipleionisation_alphaparam_p.dat";
107 }
108
109 Z = static_cast<G4double>(proton_def_->GetAtomicNumber());
110 A = static_cast<G4double>(proton_def_->GetAtomicMass());
111
112 } else if (particle == alpha_def_) {
113
114 //**************************************************************************
115 // for alpha particles
116 const auto& alpha = alpha_def_->GetParticleName();
118 eupp_tab_[alpha] = 23.0 * MeV;
119
120 // load cross-section data for single ionization process
121 auto xs_alpha = new G4DNACrossSectionDataSet(
122 new G4LogLogInterpolation, eV, kScaleFactor);
123 xs_alpha->LoadData("dna/sigma_ionisation_alphaplusplus_rudd");
124 xs_tab_[alpha] = xs_alpha;
125
126 // set energy limits
129
130 if (!use_champion_param_) {
131 alpha_param_file = "dna/multipleionisation_alphaparam_alphaplusplus.dat";
132 }
133
134 Z = static_cast<G4double>(alpha_def_->GetAtomicNumber());
135 A = static_cast<G4double>(alpha_def_->GetAtomicMass());
136
137 } else if (particle == G4GenericIon::GenericIonDefinition()) {
138
139 // *************************************************************************
140 // for carbon ions
141 const auto& carbon = carbon_def_->GetParticleName();
142 elow_tab_[carbon] = model_elow_tab_[5] * carbon_def_->GetAtomicMass();
143 eupp_tab_[carbon] = 120.0 * MeV;
144
145 // load cross-section data for single ionization process
146 auto xs_carbon = new G4DNACrossSectionDataSet(
147 new G4LogLogInterpolation, eV, kScaleFactor);
148 xs_carbon->LoadData("dna/sigma_ionisation_c_rudd");
149 xs_tab_[carbon] = xs_carbon;
150
151 // set energy limits
154
155 if (!use_champion_param_) {
156 alpha_param_file = "dna/multipleionisation_alphaparam_c.dat";
157 }
158
159 Z = static_cast<G4double>(carbon_def_->GetAtomicNumber());
160 A = static_cast<G4double>(carbon_def_->GetAtomicMass());
161
162 }
163
164 // load alpha parameter
165 mioni_manager_->LoadAlphaParam(alpha_param_file, Z, A);
166
167 if (verbose_level_ > 0) {
168 G4cout << "G4DNAQuadrupleIonisationModel is initialized " << G4endl
169 << "Energy range: "
170 << LowEnergyLimit() / eV << " eV - "
171 << HighEnergyLimit() / keV << " keV for "
172 << particle->GetParticleName()
173 << G4endl;
174 }
175
177 G4Material::GetMaterial("G4_WATER"));
178
180
181 if (is_initialized_) { return; }
182
184 is_initialized_ = true;
185}
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 G4DNAQuadrupleIonisationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * vsec,
const G4MaterialCutsCouple * couple,
const G4DynamicParticle * particle,
G4double ,
G4double  )
overridevirtual

Reimplemented from G4DNADoubleIonisationModel.

Definition at line 256 of file G4DNAQuadrupleIonisationModel.cc.

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