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

#include <G4DNABornIonisationModel.hh>

Inheritance diagram for G4DNABornIonisationModel:

Public Member Functions

 G4DNABornIonisationModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="DNABornIonisationModel")
 ~G4DNABornIonisationModel () override
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
G4double ProbabilityDensityFunction (G4double ekin) override
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void StartTracking (G4Track *) override
void SelectFasterComputation (G4bool input)
void SelectStationary (G4bool input)
void SelectSPScaling (G4bool input)
G4DNABornIonisationModeloperator= (const G4DNABornIonisationModel &right)=delete
 G4DNABornIonisationModel (const G4DNABornIonisationModel &)=delete
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 CorrectionsAlongStep (const G4Material *, const G4ParticleDefinition *, const G4double kinEnergy, const G4double cutEnergy, 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
Public Member Functions inherited from G4VSIntegration
 G4VSIntegration ()=default
virtual ~G4VSIntegration ()=default
virtual const G4StringModelName () const
void InitialiseIntegrator (G4double accuracy, G4double fact1, G4double fact2, G4double de, G4double dmin, G4double dmax)
G4double ComputeIntegral (const G4double emin, const G4double emax)
G4double SampleValue ()
 G4VSIntegration (const G4VSIntegration &)=delete
G4VSIntegrationoperator= (const G4VSIntegration &)=delete
G4bool operator== (const G4VSIntegration &right) const =delete
G4bool operator!= (const G4VSIntegration &right) const =delete
void SetVerbose (G4int verb)

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
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 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 45 of file G4DNABornIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4DNABornIonisationModel() [1/2]

G4DNABornIonisationModel::G4DNABornIonisationModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNABornIonisationModel" )

Definition at line 63 of file G4DNABornIonisationModel.cc.

64 :
65 G4VEmModel(nam)
66{
68
69 // Define default angular generator
70 SetAngularDistribution(new G4DNABornAngle());
71
72 fasterCode = G4EmParameters::Instance()->DNAFast();
73
74 if (nullptr == xsdata_p) {
75 isFirst = true;
76 LoadData();
77 }
78}
static G4EmParameters * Instance()
G4bool DNAFast() const
void SetDeexcitationFlag(G4bool val)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
void SetAngularDistribution(G4VEmAngularDistribution *)

Referenced by G4DNABornIonisationModel(), and operator=().

◆ ~G4DNABornIonisationModel()

G4DNABornIonisationModel::~G4DNABornIonisationModel ( )
override

Definition at line 82 of file G4DNABornIonisationModel.cc.

83{
84 if (isFirst) {
85 delete xsdata_e;
86 xsdata_e = nullptr;
87 delete xsdata_p;
88 xsdata_p = nullptr;
89 delete sampling_e;
90 sampling_e = nullptr;
91 delete sampling_p;
92 sampling_p = nullptr;
93 }
94}

◆ G4DNABornIonisationModel() [2/2]

G4DNABornIonisationModel::G4DNABornIonisationModel ( const G4DNABornIonisationModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 196 of file G4DNABornIonisationModel.cc.

199{
200 // check if model is applicable for given material
201 G4double density = (material->GetIndex() < fpWaterDensity->size())
202 ? (*fpWaterDensity)[material->GetIndex()] : 0.0;
203 if (0.0 == density) { return 0.0; }
204
205 // check on kinetic energy (not scaled energy) to stop low-energy ion
206 const G4double xSecMax = 1.e+10*CLHEP::barn;
207 if (ekin < fAbsorptionEnergy) { return xSecMax; }
208
209 G4double e = std::min(ekin, fHighEnergy);
210 G4double sigma = (e > fLowEnergy) ? xsdata->FindValue(e)
211 : xsdata->FindValue(fLowEnergy) * e / fLowEnergy;
212
213 sigma *= density;
214
215 // ICRU49 electronic SP scaling - ZF, SI
216 if (!isElectron && spScaling && e < fpLimitEnergy) {
217 const G4double A = 1.39241700556072800000e-9;
218 const G4double B = -8.52610412942622630000e-2;
219 sigma *= G4Exp(A*(ekin/CLHEP::eV) + B);
220 }
221 if (verbose > 1) {
222 G4cout << "G4DNABornIonisationModel for " << fParticle->GetParticleName()
223 << " Ekin(keV)=" << ekin/CLHEP::keV
224 << " sigma(cm^2)=" << sigma/CLHEP::cm2 << G4endl;
225 }
226 return sigma;
227}
G4double B(G4double temperature)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:132
double G4double
Definition G4Types.hh:83
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
std::size_t GetIndex() const

◆ Initialise()

void G4DNABornIonisationModel::Initialise ( const G4ParticleDefinition * p,
const G4DataVector &  )
overridevirtual

Implements G4VEmModel.

Definition at line 133 of file G4DNABornIonisationModel.cc.

135{
136 if (isInitialised) { return; }
138 isInitialised = true;
139
140 if (p == G4Electron::Electron()) {
141 fParticle = p;
142 xsdata = xsdata_e;
143 sampling = sampling_e;
144 fLowEnergy = 8*CLHEP::eV;
145 fHighEnergy = 1*CLHEP::MeV;
146 feLimitEnergy = 19*CLHEP::eV;
147 fAbsorptionEnergy = 6*CLHEP::eV;
148 fMass = CLHEP::electron_mass_c2;
149 isElectron = true;
150 } else if (p == G4Proton::Proton()) {
151 fParticle = p;
152 xsdata = xsdata_p;
153 sampling = sampling_p;
154 fLowEnergy = 100*CLHEP::keV;
155 fHighEnergy = 100*CLHEP::MeV;
156 fpLimitEnergy = 70*CLHEP::MeV;
157 fAbsorptionEnergy = 50*CLHEP::eV;
158 fMass = CLHEP::proton_mass_c2;
159 isElectron = false;
160 } else {
162 ed << "Born ionisation model is used for " << p->GetParticleName();
163 G4Exception("G4DNABornIonisationModel::Initialise","em0003",
164 FatalException, ed, " it is not available.");
165 }
167
168 // defined stationary mode
170
171 // initialise atomic de-excitation
172 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
173
174 // chemistry
176 if (chem->IsChemistryActivated()) {
177 fChemistry = chem;
178 }
179
180 InitialiseIntegrator(0.1, 0.25, 1.05, 1*CLHEP::eV, 0.2*CLHEP::eV, 10*CLHEP::keV);
181 if (verbose > 1) {
182 G4cout << "Born ionisation model is initialized for "
183 << fParticle->GetParticleName() << G4endl;
184 }
185}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4ParticleChangeForGamma * fParticleChangeForGamma
static G4DNAChemistryManager * Instance()
static G4Electron * Electron()
Definition G4Electron.cc:91
G4bool DNAStationary() const
G4int WorkerVerbose() const
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4String & GetParticleName() const
static G4Proton * Proton()
Definition G4Proton.cc:90
G4ParticleChangeForGamma * GetParticleChangeForGamma()
void InitialiseIntegrator(G4double accuracy, G4double fact1, G4double fact2, G4double de, G4double dmin, G4double dmax)

◆ operator=()

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

◆ ProbabilityDensityFunction()

G4double G4DNABornIonisationModel::ProbabilityDensityFunction ( G4double ekin)
overridevirtual

Implements G4VSIntegration.

Definition at line 380 of file G4DNABornIonisationModel.cc.

381{
382 return sampling->GetValue(fPrimaryEnergy, ekin, fSelectedShell);
383}

◆ SampleSecondaries()

void G4DNABornIonisationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * fvect,
const G4MaterialCutsCouple * couple,
const G4DynamicParticle * dynParticle,
G4double tmin,
G4double maxEnergy )
overridevirtual

Implements G4VEmModel.

Definition at line 231 of file G4DNABornIonisationModel.cc.

235{
236 fPrimaryEnergy = dynParticle->GetKineticEnergy();
237 // proton shoud be stopped - check on kinetic energy
238 // electrons never have such low energy
239 if (fPrimaryEnergy <= fAbsorptionEnergy) {
240 fParticleChangeForGamma->SetProposedKineticEnergy(0.);
241 fParticleChangeForGamma->ProposeTrackStatus(fStopButAlive);
242 fParticleChangeForGamma->ProposeLocalEnergyDeposit(fPrimaryEnergy);
243 return;
244 }
245
246 fSelectedShell = SelectShell();
247 G4double bindingEnergy = waterStructure.IonisationEnergy(fSelectedShell);
248
249 //SI: additional protection if tcs interpolation method is modified
250 if (fPrimaryEnergy < bindingEnergy) { return; }
251
252 // compute max energy
253 if (isElectron) {
254 fMaxEnergy = 0.5*(fPrimaryEnergy - bindingEnergy);
255 } else {
256 G4double tau = fPrimaryEnergy/fMass;
257 fMaxEnergy = 2.0*CLHEP::electron_mass_c2*tau*(tau + 2.0);
258 }
259 // SI: The following protection is necessary to avoid infinite loops :
260 // e- ionisation cross section has non zero partial xs at 18 eV for shell 2.
261 // e- has zero cumulated partial xs at 18 eV for shell 2.
262 // This is due to the fact that the max allowed transfered energy is
263 // (18+10.79)/2=17.025 eV and only transfered energies strictly above this
264 // value have non zero partial cross section starting at transition energy 17.12 eV.
265 if (fasterCode && isElectron && 2 == fSelectedShell && fPrimaryEnergy < feLimitEnergy) {
266 do {
267 fSelectedShell = SelectShell();
268 } while (2 == fSelectedShell);
269 }
270
271 G4double esec = fasterCode ? SampleCumulative() : SampleDifferential();
272 G4double esum = 0.0;
273
274 // sample deexcitation
275 // here we assume that H2O electronic levels are the same as Oxygen.
276 // this can be considered true with a rough 10% error in energy on K-shell,
277 G4int Z = 8;
278 G4ThreeVector deltaDir =
279 GetAngularDistribution()->SampleDirectionForShell(dynParticle, esec, Z,
280 fSelectedShell,
281 couple->GetMaterial());
282
283 // SI: only atomic deexcitation from K shell is considered
284 if (fAtomDeexcitation != nullptr && fSelectedShell == 4) {
285 auto as = G4AtomicShellEnumerator(0);
286 auto ashell = fAtomDeexcitation->GetAtomicShell(Z, as);
287 fAtomDeexcitation->GenerateParticles(fvect, ashell, Z, 0, 0);
288
289 // compute energy sum from de-excitation
290 for (auto const & ptr : *fvect) {
291 esum += ptr->GetKineticEnergy();
292 }
293 }
294 // check energy balance
295 // remaining excitation energy of water molecule
296 G4double exc = std::max(bindingEnergy - esum, 0.0);
297
298 // remaining projectile energy
299 G4double scatteredEnergy = fPrimaryEnergy - bindingEnergy - esec;
300 if (scatteredEnergy < -tolerance || exc < -tolerance) {
301 G4cout << "G4DNABornIonisationModel::SampleSecondaries: "
302 << "final E(keV)=" << scatteredEnergy/CLHEP::keV << " Ein(keV)="
303 << fPrimaryEnergy/CLHEP::keV << " " << fParticle->GetParticleName()
304 << " Edelta(keV)=" << esec/CLHEP::keV << " MeV, Exc(keV)=" << exc/CLHEP::keV
305 << G4endl;
306 }
307 scatteredEnergy = std::max(scatteredEnergy, 0.0);
308
309 // projectile
310 if (!statCode) {
311 fParticleChangeForGamma->SetProposedKineticEnergy(scatteredEnergy);
312 fParticleChangeForGamma->ProposeLocalEnergyDeposit(exc);
313 } else {
314 fParticleChangeForGamma->SetProposedKineticEnergy(fPrimaryEnergy);
315 fParticleChangeForGamma->ProposeLocalEnergyDeposit(fPrimaryEnergy - scatteredEnergy);
316 }
317
318 // delta-electron
319 auto dp = new G4DynamicParticle(G4Electron::Electron(), deltaDir, esec);
320 fvect->push_back(dp);
321
322 // create radical
323 if (nullptr != fChemistry) {
324 fChemistry->CreateWaterMolecule(eIonizedMolecule, fSelectedShell, fTrack);
325 }
326}
@ eIonizedMolecule
CLHEP::Hep3Vector G4ThreeVector
@ fStopButAlive
int G4int
Definition G4Types.hh:85
G4double GetKineticEnergy() const
const G4Material * GetMaterial() const
virtual G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
G4VEmAngularDistribution * GetAngularDistribution()
G4double bindingEnergy(G4int A, G4int Z)

◆ SelectFasterComputation()

void G4DNABornIonisationModel::SelectFasterComputation ( G4bool input)
inline

Definition at line 71 of file G4DNABornIonisationModel.hh.

71{ fasterCode = input; };

◆ SelectSPScaling()

void G4DNABornIonisationModel::SelectSPScaling ( G4bool input)
inline

Definition at line 75 of file G4DNABornIonisationModel.hh.

75{ spScaling = input; };

◆ SelectStationary()

void G4DNABornIonisationModel::SelectStationary ( G4bool input)
inline

Definition at line 73 of file G4DNABornIonisationModel.hh.

73{ statCode = input; };

◆ StartTracking()

void G4DNABornIonisationModel::StartTracking ( G4Track * track)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 189 of file G4DNABornIonisationModel.cc.

190{
191 fTrack = track;
192}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNABornIonisationModel::fParticleChangeForGamma
protected

Definition at line 92 of file G4DNABornIonisationModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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