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

#include <G4DNAIonElasticModel.hh>

+ Inheritance diagram for G4DNAIonElasticModel:

Public Member Functions

 G4DNAIonElasticModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="DNAIonElasticModel")
 
 ~G4DNAIonElasticModel () override
 
G4DNAIonElasticModeloperator= (const G4DNAIonElasticModel &right)=delete
 
 G4DNAIonElasticModel (const G4DNAIonElasticModel &)=delete
 
void Initialise (const G4ParticleDefinition *particuleDefinition, const G4DataVector &) 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 SetKillBelowThreshold (G4double threshold)
 
G4double GetKillBelowThreshold ()
 
void SelectStationary (G4bool input)
 
- 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
 
- 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 46 of file G4DNAIonElasticModel.hh.

Constructor & Destructor Documentation

◆ G4DNAIonElasticModel() [1/2]

G4DNAIonElasticModel::G4DNAIonElasticModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNAIonElasticModel" )

Definition at line 46 of file G4DNAIonElasticModel.cc.

47 :
48 G4VEmModel(nam)
49{
50 killBelowEnergy = 100 * eV;
51 lowEnergyLimit = 0 * eV;
52 highEnergyLimit = 1 * MeV;
53 SetLowEnergyLimit(lowEnergyLimit);
54 SetHighEnergyLimit(highEnergyLimit);
55
56 verboseLevel = 0;
57 // Verbosity scale:
58 // 0 = nothing
59 // 1 = warning for energy non-conservation
60 // 2 = details of energy budget
61 // 3 = calculation of cross sections, file openings, sampling of atoms
62 // 4 = entering in methods
63
64 if(verboseLevel > 0)
65 {
66 G4cout << "Ion elastic model is constructed " << G4endl<< "Energy range: "
67 << lowEnergyLimit / eV << " eV - "
68 << highEnergyLimit / MeV << " MeV"
69 << G4endl;
70 }
71
73 fpMolWaterDensity = nullptr;
74 fpTableData = nullptr;
75 fParticle_Mass = -1;
76
77 // Selection of stationary mode
78
79 statCode = false;
80}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetHighEnergyLimit(G4double)
void SetLowEnergyLimit(G4double)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4DNAIonElasticModel()

G4DNAIonElasticModel::~G4DNAIonElasticModel ( )
override

Definition at line 84 of file G4DNAIonElasticModel.cc.

85{
86 // For total cross section
87 delete fpTableData;
88}

◆ G4DNAIonElasticModel() [2/2]

G4DNAIonElasticModel::G4DNAIonElasticModel ( const G4DNAIonElasticModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 246 of file G4DNAIonElasticModel.cc.

249{
250 if(verboseLevel > 3)
251 {
252 G4cout << "Calling CrossSectionPerVolume() of G4DNAIonElasticModel"
253 << G4endl;
254 }
255
256 // Calculate total cross section for model
257
258 G4double sigma=0;
259
260 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
261
262 const G4String& particleName = p->GetParticleName();
263
264 if (ekin <= highEnergyLimit)
265 {
266 //SI : XS must not be zero otherwise sampling of secondaries method ignored
267 if (ekin < killBelowEnergy) return DBL_MAX;
268 //
269
270 if (fpTableData != nullptr)
271 {
272 sigma = fpTableData->FindValue(ekin);
273 }
274 else
275 {
276 G4Exception("G4DNAIonElasticModel::ComputeCrossSectionPerVolume","em0002",
277 FatalException,"Model not applicable to particle type.");
278 }
279 }
280
281 if (verboseLevel > 2)
282 {
283 G4cout << "__________________________________" << G4endl;
284 G4cout << "G4DNAIonElasticModel - XS INFO START" << G4endl;
285 G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
286 G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
287 G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
288 G4cout << "G4DNAIonElasticModel - XS INFO END" << G4endl;
289 }
290
291 return sigma*waterDensity;
292}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
G4double FindValue(G4double e, G4int componentId=0) const override
std::size_t GetIndex() const
const G4String & GetParticleName() const
#define DBL_MAX
Definition templates.hh:62

◆ GetKillBelowThreshold()

G4double G4DNAIonElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 79 of file G4DNAIonElasticModel.hh.

80 {
81 return killBelowEnergy;
82 }

◆ Initialise()

void G4DNAIonElasticModel::Initialise ( const G4ParticleDefinition * particuleDefinition,
const G4DataVector &  )
overridevirtual

Implements G4VEmModel.

Definition at line 93 of file G4DNAIonElasticModel.cc.

96{
97
98 if(verboseLevel > 3)
99 {
100 G4cout << "Calling G4DNAIonElasticModel::Initialise()" << G4endl;
101 }
102
103 // Energy limits
104
105 if (LowEnergyLimit() < lowEnergyLimit)
106 {
107 G4cout << "G4DNAIonElasticModel: low energy limit increased from " <<
108 LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
109 SetLowEnergyLimit(lowEnergyLimit);
110 }
111
112 if (HighEnergyLimit() > highEnergyLimit)
113 {
114 G4cout << "G4DNAIonElasticModel: high energy limit decreased from " <<
115 HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl;
116 SetHighEnergyLimit(highEnergyLimit);
117 }
118
119 // Reading of data files
120
121 G4double scaleFactor = 1e-16*cm*cm;
122
123 const char *path = G4FindDataDir("G4LEDATA");
124
125 if (path == nullptr)
126 {
127 G4Exception("G4IonElasticModel::Initialise","em0006",
128 FatalException,"G4LEDATA environment variable not set.");
129 return;
130 }
131
132 G4String totalXSFile;
133 std::ostringstream fullFileName;
134
135 G4DNAGenericIonsManager *instance;
137 G4ParticleDefinition* protonDef =
139 G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
140 G4ParticleDefinition* heliumDef = instance->GetIon("helium");
141 G4ParticleDefinition* alphaplusDef = instance->GetIon("alpha+");
142 G4ParticleDefinition* alphaplusplusDef = instance->GetIon("alpha++");
143 G4String proton, hydrogen, helium, alphaplus, alphaplusplus;
144
145 if (
146 (particleDefinition == protonDef && protonDef != nullptr)
147 ||
148 (particleDefinition == hydrogenDef && hydrogenDef != nullptr)
149 )
150 {
151 // For total cross section of p,h
152 fParticle_Mass = 1.;
153 totalXSFile = "dna/sigma_elastic_proton_HTran";
154
155 // For final state
156 fullFileName << path << "/dna/sigmadiff_cumulated_elastic_proton_HTran.dat";
157 }
158
159 if (
160 (particleDefinition == instance->GetIon("helium") && (heliumDef != nullptr))
161 ||
162 (particleDefinition == instance->GetIon("alpha+") && (alphaplusDef != nullptr))
163 ||
164 (particleDefinition == instance->GetIon("alpha++") && (alphaplusplusDef != nullptr))
165 )
166 {
167 // For total cross section of he,he+,he++
168 fParticle_Mass = 4.;
169 totalXSFile = "dna/sigma_elastic_alpha_HTran";
170
171 // For final state
172 fullFileName << path << "/dna/sigmadiff_cumulated_elastic_alpha_HTran.dat";
173 }
174
175 fpTableData = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
176 fpTableData->LoadData(totalXSFile);
177 std::ifstream diffCrossSection(fullFileName.str().c_str());
178
179 if (!diffCrossSection)
180 {
181 G4ExceptionDescription description;
182 description << "Missing data file:"
183 <<fullFileName.str().c_str()<< G4endl;
184 G4Exception("G4IonElasticModel::Initialise","em0003",
186 description);
187 }
188
189 // Added clear for MT
190
191 eTdummyVec.clear();
192 eVecm.clear();
193 fDiffCrossSectionData.clear();
194
195 //
196
197 eTdummyVec.push_back(0.);
198
199 while(!diffCrossSection.eof())
200 {
201 G4double tDummy;
202 G4double eDummy;
203 diffCrossSection>>tDummy>>eDummy;
204
205 // SI : mandatory eVecm initialization
206
207 if (tDummy != eTdummyVec.back())
208 {
209 eTdummyVec.push_back(tDummy);
210 eVecm[tDummy].push_back(0.);
211 }
212
213 diffCrossSection>>fDiffCrossSectionData[tDummy][eDummy];
214
215 if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy);
216
217 }
218
219 // End final state
220 if( verboseLevel>0 )
221 {
222 if (verboseLevel > 2)
223 {
224 G4cout << "Loaded cross section files for ion elastic model" << G4endl;
225 }
226 G4cout << "Ion elastic model is initialized " << G4endl
227 << "Energy range: "
228 << LowEnergyLimit() / eV << " eV - "
229 << HighEnergyLimit() / MeV << " MeV"
230 << G4endl;
231 }
232
233 // Initialize water density pointer
235 fpMolWaterDensity = G4DNAMolecularMaterial::Instance()->
236 GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
237
238 if (isInitialised) return;
240 isInitialised = true;
241}
const char * G4FindDataDir(const char *)
std::ostringstream G4ExceptionDescription
G4bool LoadData(const G4String &argFileName) override
static G4DNAGenericIonsManager * Instance()
G4ParticleDefinition * GetIon(const G4String &name)
static G4DNAMolecularMaterial * Instance()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const

◆ operator=()

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

◆ SampleSecondaries()

void G4DNAIonElasticModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * ,
const G4MaterialCutsCouple * ,
const G4DynamicParticle * aDynamicParticle,
G4double tmin,
G4double maxEnergy )
overridevirtual

Implements G4VEmModel.

Definition at line 297 of file G4DNAIonElasticModel.cc.

301{
302
303 if(verboseLevel > 3)
304 {
305 G4cout << "Calling SampleSecondaries() of G4DNAIonElasticModel" << G4endl;
306 }
307
308 G4double particleEnergy0 = aDynamicParticle->GetKineticEnergy();
309
310 if (particleEnergy0 < killBelowEnergy)
311 {
315 return;
316 }
317
318 if (particleEnergy0>= killBelowEnergy && particleEnergy0 <= highEnergyLimit)
319 {
320 G4double water_mass = 18.;
321
322 G4double thetaCM = RandomizeThetaCM(particleEnergy0, aDynamicParticle->GetDefinition());
323
324 //HT:convert to laboratory system
325
326 G4double theta = std::atan(std::sin(thetaCM*CLHEP::pi/180)
327 /(fParticle_Mass/water_mass+std::cos(thetaCM*CLHEP::pi/180)));
328
329 G4double cosTheta= std::cos(theta);
330
331 //
332
333 G4double phi = 2. * CLHEP::pi * G4UniformRand();
334
335 G4ThreeVector zVers = aDynamicParticle->GetMomentumDirection();
336 G4ThreeVector xVers = zVers.orthogonal();
337 G4ThreeVector yVers = zVers.cross(xVers);
338
339 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
340 G4double yDir = xDir;
341 xDir *= std::cos(phi);
342 yDir *= std::sin(phi);
343
344 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
345
347
348 G4double depositEnergyCM = 0;
349
350 //HT: deposited energy
351 depositEnergyCM = 4. * particleEnergy0 * fParticle_Mass * water_mass *
352 (1-std::cos(thetaCM*CLHEP::pi/180))
353 / (2 * std::pow((fParticle_Mass+water_mass),2));
354
355 //SI: added protection particleEnergy0 >= depositEnergyCM
356 if (!statCode && (particleEnergy0 >= depositEnergyCM) )
357
358 fParticleChangeForGamma->SetProposedKineticEnergy(particleEnergy0 - depositEnergyCM);
359
361
363 }
364}
@ fStopAndKill
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SelectStationary()

void G4DNAIonElasticModel::SelectStationary ( G4bool input)
inline

Definition at line 146 of file G4DNAIonElasticModel.hh.

147{
148 statCode = input;
149}

◆ SetKillBelowThreshold()

void G4DNAIonElasticModel::SetKillBelowThreshold ( G4double threshold)

Definition at line 506 of file G4DNAIonElasticModel.cc.

507{
508 killBelowEnergy = threshold;
509
510 if(killBelowEnergy < 100 * eV)
511 {
512 G4cout << "*** WARNING : the G4DNAIonElasticModel class is not "
513 "activated below 100 eV !"
514 << G4endl;
515 }
516}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAIonElasticModel::fParticleChangeForGamma
protected

Definition at line 88 of file G4DNAIonElasticModel.hh.

Referenced by G4DNAIonElasticModel(), Initialise(), and SampleSecondaries().


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