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

#include <G4DNAELSEPAElasticModel.hh>

+ Inheritance diagram for G4DNAELSEPAElasticModel:

Public Member Functions

 G4DNAELSEPAElasticModel (const G4ParticleDefinition *particle=0, const G4String &nam="DNAELSEPAElasticModel")
 
virtual ~G4DNAELSEPAElasticModel ()
 
virtual void Initialise (const G4ParticleDefinition *particle, const G4DataVector &)
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *particle, G4double ekin, G4double emin, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SetMaximumEnergy (G4double input)
 
void SetKillBelowThreshold (G4double threshold)
 
G4double GetKillBelowThreshold ()
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
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 CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=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 LPMFlag () 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 SetLPMFlag (G4bool val)
 
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)
 
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
 
size_t currentCoupleIndex = 0
 
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 49 of file G4DNAELSEPAElasticModel.hh.

Constructor & Destructor Documentation

◆ G4DNAELSEPAElasticModel()

G4DNAELSEPAElasticModel::G4DNAELSEPAElasticModel ( const G4ParticleDefinition particle = 0,
const G4String nam = "DNAELSEPAElasticModel" 
)

Definition at line 48 of file G4DNAELSEPAElasticModel.cc.

49 :
50G4VEmModel(nam), isInitialised(false)
51{
52 verboseLevel = 0;
53
54 G4ProductionCutsTable* theCoupleTable =
56 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
57
58 for(G4int i=0; i<numOfCouples; ++i)
59 {
60 const G4MaterialCutsCouple* couple =
61 theCoupleTable->GetMaterialCutsCouple(i);
62 const G4Material* material = couple->GetMaterial();
63 G4int nelm = (G4int)material->GetNumberOfElements();
64 const G4ElementVector* theElementVector = material->GetElementVector();
65
66 if(nelm==1)
67 {// Protection: only for single element
68 G4int Z = 79;
69 Z = G4lrint((*theElementVector)[0]->GetZ());
70 // Protection: only for GOLD
71 if (Z==79){
72 fkillBelowEnergy_Au = 10. * eV; // Kills e- tracking
73 flowEnergyLimit = 0 * eV; // Must stay at zero for killing
74 fhighEnergyLimit = 1 * GeV; // Default
75 SetLowEnergyLimit (flowEnergyLimit);
76 SetHighEnergyLimit(fhighEnergyLimit);
77 }else{
78 //continue;
79 }
80 }else{// Protection: H2O only is available
81 if(material->GetName()=="G4_WATER"){
82 flowEnergyLimit = 10. * eV;
83 fhighEnergyLimit = 1 * MeV;
84 SetLowEnergyLimit (flowEnergyLimit);
85 SetHighEnergyLimit(fhighEnergyLimit);
86 }else{
87 //continue;
88 }
89 }
90
91 if (verboseLevel > 0)
92 {
93 G4cout << "ELSEPA Elastic model is constructed for "
94 << material->GetName() << G4endl
95 << "Energy range: "
96 << flowEnergyLimit / eV << " eV - "
97 << fhighEnergyLimit / MeV << " MeV"
98 << G4endl;
99 }
100 }
101
102
104 fpMolDensity = 0;
105
106 fpData_Au=nullptr;
107 fpData_H2O=nullptr;
108}
std::vector< const G4Element * > G4ElementVector
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4String & GetName() const
Definition: G4Material.hh:172
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
int G4lrint(double ad)
Definition: templates.hh:134

◆ ~G4DNAELSEPAElasticModel()

G4DNAELSEPAElasticModel::~G4DNAELSEPAElasticModel ( )
virtual

Definition at line 112 of file G4DNAELSEPAElasticModel.cc.

113{
114 //std::map<G4int,G4DNACrossSectionDataSet*,
115 // std::less<G4String>>::iterator posZ;
116 //for (posZ = tableZData.begin(); posZ != tableZData.end(); ++posZ)
117 //{
118 // G4DNACrossSectionDataSet* table = posZ->second;
119 // delete table;
120 //}
121 //for (posZ = tableZData_Au.begin(); posZ != tableZData_Au.end(); ++posZ)
122 //{
123 // G4DNACrossSectionDataSet* table = posZ->second;
124 // delete table;
125 //}
126 //for (posZ = tableZData_H2O.begin(); posZ != tableZData_H2O.end(); ++posZ)
127 //{
128 // G4DNACrossSectionDataSet* table = posZ->second;
129 // delete table;
130 //}
131
132 if(fpData_Au) delete fpData_Au;
133 if(fpData_H2O) delete fpData_H2O;
134
135 //eEdummyVecZ.clear();
136 //eCumZ.clear();
137 //fAngleDataZ.clear();
138
139 eEdummyVec_Au.clear();
140 eEdummyVec_H2O.clear();
141 eCum_Au.clear();
142 eCum_H2O.clear();
143 fAngleData_Au.clear();
144 fAngleData_H2O.clear();
145}

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4DNAELSEPAElasticModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition particle,
G4double  ekin,
G4double  emin,
G4double  emax 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 388 of file G4DNAELSEPAElasticModel.cc.

394{
395
396 if (verboseLevel > 3)
397 {
398 G4cout <<
399 "Calling CrossSectionPerVolume() of G4DNAELSEPAElasticModel"
400 << G4endl;
401 }
402
403 G4double atomicNDensity=0.0;
404 G4double sigma=0;
405
406 const G4ElementVector* theElementVector = material->GetElementVector();
407 std::size_t nelm = material->GetNumberOfElements();
408 if (nelm==1) // Protection: only for single element
409 {
410 // Protection: only for GOLD
411 if (material->GetZ()!=79) return 0.0;
412
413 G4int Z = G4lrint((*theElementVector)[0]->GetZ());
414
415 const G4String& particleName = particle->GetParticleName();
416 atomicNDensity = material->GetAtomicNumDensityVector()[0];
417 if(atomicNDensity!= 0.0)
418 {
419 if (ekin < fhighEnergyLimit)
420 {
421 if (ekin < fkillBelowEnergy_Au) return DBL_MAX;
422
423 //std::map< G4int,G4DNACrossSectionDataSet*,
424 // std::less<G4String> >::iterator pos;
425 ////pos = tableZData_Au.find(0);
426 //pos = tableZData.find(Z);
427 //
428 ////if (pos != tableZData_Au.end())
429 //if (pos != tableZData.end())
430 //{
431 // G4DNACrossSectionDataSet* table = pos->second;
432 // if (table != 0)
433 // {
434 // // XS takes its 10 eV value below 10 eV for GOLD
435 // if (ekin < 10*eV) sigma = table->FindValue(10*eV);
436 // else sigma = table->FindValue(ekin);
437 // }
438 //}
439 //else
440 //{
441 // G4Exception("G4DNAELSEPAElasticModel::ComputeCrossSectionPerVolume",
442 // "em0006",FatalException,"Model not applicable to particle type.");
443 //}
444
445 if (ekin < 10*eV) sigma = fpData_Au->FindValue(10*eV);
446 else sigma = fpData_Au->FindValue(ekin);
447 }
448 }
449 if (verboseLevel > 2)
450 {
451 G4cout << "__________________________________" << G4endl;
452 G4cout << "=== G4DNAELSEPAElasticModel - XS INFO START" << G4endl;
453 G4cout << "=== Material is made of one element with Z =" << Z << G4endl;
454 G4cout << "=== Kinetic energy(eV)=" << ekin/eV << " particle : "
455 << particleName << G4endl;
456 G4cout << "=== Cross section per atom for Z="<<Z<<" is (cm^2)"
457 << sigma/cm/cm << G4endl;
458 G4cout << "=== Cross section per atom for Z="<<Z<<" is (cm^-1)="
459 << sigma*atomicNDensity/(1./cm) << G4endl;
460 G4cout << "=== G4DNAELSEPAElasticModel - XS INFO END" << G4endl;
461 }
462 }
463 else
464 {
465 fpMolDensity =
467 GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
468 atomicNDensity = (*fpMolDensity)[material->GetIndex()];
469 if(atomicNDensity!= 0.0)
470 {
471 if (ekin < HighEnergyLimit() && ekin >= LowEnergyLimit())
472 {
473 //std::map< G4int,G4DNACrossSectionDataSet*,
474 //std::less<G4String> >::iterator pos;
475 ////pos = tableZData_H2O.find(0); // the data is stored as Z=0
476 //pos = tableZData.find(0); // the data is stored as Z=0
477 ////SI : XS must not be zero
478 //// otherwise sampling of secondaries method ignored
479 ////if (pos != tableZData_H2O.end())
480 //if (pos != tableZData.end())
481 //{
482 // G4DNACrossSectionDataSet* table = pos->second;
483 // if (table != 0)
484 // {
485 // sigma = table->FindValue(ekin);
486 // }
487 //}
488
489 sigma = fpData_H2O->FindValue(ekin);
490 }
491 }
492 if (verboseLevel > 2)
493 {
494 G4cout << "__________________________________" << G4endl;
495 G4cout << "=== G4DNAELSEPAElasticModel - XS INFO START" << G4endl;
496 G4cout << "=== Kinetic energy(eV)=" << ekin/eV
497 << " particle : " << particle->GetParticleName() << G4endl;
498 G4cout << "=== Cross section per water molecule (cm^2)="
499 << sigma/cm/cm << G4endl;
500 G4cout << "=== Cross section per water molecule (cm^-1)="
501 << sigma*atomicNDensity/(1./cm) << G4endl;
502 G4cout << "=== G4DNAELSEPAElasticModel - XS INFO END" << G4endl;
503 }
504 }
505
506 return sigma*atomicNDensity;
507}
double G4double
Definition: G4Types.hh:83
virtual G4double FindValue(G4double e, G4int componentId=0) const
static G4DNAMolecularMaterial * Instance()
G4double GetZ() const
Definition: G4Material.cc:745
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:211
size_t GetIndex() const
Definition: G4Material.hh:255
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:691
const G4String & GetParticleName() const
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
#define DBL_MAX
Definition: templates.hh:62

◆ GetKillBelowThreshold()

G4double G4DNAELSEPAElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 79 of file G4DNAELSEPAElasticModel.hh.

79{return fkillBelowEnergy_Au;}

◆ Initialise()

void G4DNAELSEPAElasticModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 149 of file G4DNAELSEPAElasticModel.cc.

151{
152 if (verboseLevel > 3)
153 G4cout << "Calling G4DNAELSEPAElasticModel::Initialise()" << G4endl;
154
155 if (isInitialised) {return;}
156
157 if(particle->GetParticleName() != "e-")
158 {
159 G4Exception("G4DNAELSEPAElasticModel::Initialise","em0001",
160 FatalException,"Model not applicable to particle type.");
161 return;
162 }
163
164 G4ProductionCutsTable* theCoupleTable =
166 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
167
168 // UNIT OF TCS
169 G4double scaleFactor = 1.*cm*cm;
170
171 //tableZData.clear();
172 //tableZData_Au.clear();
173 //tableZData_H2O.clear();
174
175 fpData_Au=nullptr;
176 fpData_H2O=nullptr;
177
178 for(G4int i=0; i<numOfCouples; ++i)
179 {
180 const G4MaterialCutsCouple* couple =
181 theCoupleTable->GetMaterialCutsCouple(i);
182 const G4Material* material = couple->GetMaterial();
183 const G4ElementVector* theElementVector = material->GetElementVector();
184
185 G4int nelm = (G4int)material->GetNumberOfElements();
186 if (nelm==1){// Protection: only for single element
187 G4int Z = G4lrint((*theElementVector)[0]->GetZ());
188 if (Z!=79)// Protection: only for GOLD
189 {
190 continue;
191 }
192
193 if (Z>0)
194 {
195 G4String fileZElectron("dna/sigma_elastic_e_elsepa_Z");
196 std::ostringstream oss;
197 oss.str("");
198 oss.clear(stringstream::goodbit);
199 oss << Z;
200 fileZElectron += oss.str()+"_muffintin";
201
202 //G4DNACrossSectionDataSet* tableZE =
203 // new G4DNACrossSectionDataSet
204 // (new G4LogLogInterpolation, eV,scaleFactor );
205 //tableZE->LoadData(fileZElectron);
206 ////tableZData_Au[0] = tableZE;
207 //tableZData[Z] = tableZE;
208
210 eV,
211 scaleFactor );
212 fpData_Au->LoadData(fileZElectron);
213
214 std::ostringstream eFullFileNameZ;
215 const char *path = G4FindDataDir("G4LEDATA");
216 if (!path)
217 {
218 G4Exception("G4DNAELSEPAElasticModel::Initialise","em0002",
219 FatalException,"G4LEDATA environment variable not set.");
220 return;
221 }
222
223 eFullFileNameZ.str("");
224 eFullFileNameZ.clear(stringstream::goodbit);
225
226 eFullFileNameZ
227 << path
228 << "/dna/sigmadiff_cumulated_elastic_e_elsepa_Z"
229 << Z << "_muffintin.dat";
230
231 std::ifstream eDiffCrossSectionZ(eFullFileNameZ.str().c_str());
232
233 if (!eDiffCrossSectionZ)
234 {
235 G4Exception("G4DNAELSEPAElasticModel::Initialise","em0003",
236 FatalException,"Missing data file for cumulated DCS");
237 return;
238 }
239
240 //eEdummyVecZ.clear();
241 //eCumZ.clear();
242 //fAngleDataZ.clear();
243
244 eEdummyVec_Au.clear();
245 eCum_Au.clear();
246 fAngleData_Au.clear();
247
248 //eEdummyVecZ[Z].push_back(0.);
249 eEdummyVec_Au.push_back(0.);
250 do
251 {
252 G4double eDummy;
253 G4double cumDummy;
254 eDiffCrossSectionZ>>eDummy>>cumDummy;
255 //if (eDummy != eEdummyVecZ[Z].back())
256 if (eDummy != eEdummyVec_Au.back())
257 {
258
259 //eEdummyVecZ[Z].push_back(eDummy);
260 eEdummyVec_Au.push_back(eDummy);
261 //eCumZ[Z][eDummy].push_back(0.);
262 eCum_Au[eDummy].push_back(0.);
263 }
264 //eDiffCrossSectionZ>>fAngleDataZ[Z][eDummy][cumDummy];
265 eDiffCrossSectionZ>>fAngleData_Au[eDummy][cumDummy];
266 //if (cumDummy != eCumZ[Z][eDummy].back())
267 if (cumDummy != eCum_Au[eDummy].back())
268 {
269 //eCumZ[Z][eDummy].push_back(cumDummy);
270 eCum_Au[eDummy].push_back(cumDummy);
271 }
272 }while(!eDiffCrossSectionZ.eof());
273 }
274
275 }else{// Protection: H2O only is available
276 if(material->GetName()=="G4_WATER"){
277 if (LowEnergyLimit() < 10*eV)
278 {
279 G4cout<<"G4DNAELSEPAElasticModel: low energy limit increased from "
280 << LowEnergyLimit()/eV << " eV to " << 10 << " eV"
281 << G4endl;
282 SetLowEnergyLimit(10.*eV);
283 }
284
285 if (HighEnergyLimit() > 1.*MeV)
286 {
287 G4cout<<"G4DNAELSEPAElasticModel: high energy limit decreased from "
288 << HighEnergyLimit()/MeV << " MeV to " << 1. << " MeV"
289 << G4endl;
290 SetHighEnergyLimit(1.*MeV);
291 }
292
293 G4String fileZElectron("dna/sigma_elastic_e_elsepa_muffin");
294
295 //G4DNACrossSectionDataSet* tableZE =
296 // new G4DNACrossSectionDataSet(
297 // new G4LogLogInterpolation, eV,scaleFactor );
298 //tableZE->LoadData(fileZElectron);
299 ////tableZData_H2O[0] = tableZE;
300 //tableZData[0] = tableZE;
301
303 eV,
304 scaleFactor );
305 fpData_H2O->LoadData(fileZElectron);
306
307 std::ostringstream eFullFileNameZ;
308
309 const char *path = G4FindDataDir("G4LEDATA");
310 if (!path)
311 {
312 G4Exception("G4DNAELSEPAElasticModel::Initialise","em0004",
313 FatalException,"G4LEDATA environment variable not set.");
314 return;
315 }
316
317 eFullFileNameZ.str("");
318 eFullFileNameZ.clear(stringstream::goodbit);
319
320 eFullFileNameZ
321 << path
322 << "/dna/sigmadiff_cumulated_elastic_e_elsepa_muffin.dat";
323
324 std::ifstream eDiffCrossSectionZ(eFullFileNameZ.str().c_str());
325
326 if (!eDiffCrossSectionZ)
327 G4Exception("G4DNAELSEPAElasticModel::Initialise","em0005",
329 "Missing data file for cumulated DCS");
330
331 //eEdummyVecZ.clear();
332 //eCumZ.clear();
333 //fAngleDataZ.clear();
334
335 eEdummyVec_H2O.clear();
336 eCum_H2O.clear();
337 fAngleData_H2O.clear();
338
339 //eEdummyVecZ[0].push_back(0.);
340 eEdummyVec_H2O.push_back(0.);
341
342 do
343 {
344 G4double eDummy;
345 G4double cumDummy;
346 eDiffCrossSectionZ>>eDummy>>cumDummy;
347 //if (eDummy != eEdummyVecZ[0].back())
348 if (eDummy != eEdummyVec_H2O.back())
349 {
350 //eEdummyVecZ[0].push_back(eDummy);
351 eEdummyVec_H2O.push_back(eDummy);
352 //eCumZ[0][eDummy].push_back(0.);
353 eCum_H2O[eDummy].push_back(0.);
354 }
355 //eDiffCrossSectionZ>>fAngleDataZ[0][eDummy][cumDummy];
356 eDiffCrossSectionZ>>fAngleData_H2O[eDummy][cumDummy];
357 //if (cumDummy != eCumZ[0][eDummy].back()){
358 if (cumDummy != eCum_H2O[eDummy].back()){
359 //eCumZ[0][eDummy].push_back(cumDummy);
360 eCum_H2O[eDummy].push_back(cumDummy);
361 }
362 }while(!eDiffCrossSectionZ.eof());
363 }
364 }
365 if (verboseLevel > 2)
366 G4cout << "Loaded cross section files of ELSEPA Elastic model for"
367 << material->GetName() << G4endl;
368
369 if( verboseLevel>0 )
370 {
371 G4cout << "ELSEPA elastic model is initialized " << G4endl
372 << "Energy range: "
373 << LowEnergyLimit() / eV << " eV - "
374 << HighEnergyLimit()/ MeV << " MeV"
375 << G4endl;
376 }
377 } // Loop on couples
378
379
381 fpMolDensity = 0;
382
383 isInitialised = true;
384}
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
virtual G4bool LoadData(const G4String &argFileName)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124

◆ SampleSecondaries()

void G4DNAELSEPAElasticModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicElectron,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 511 of file G4DNAELSEPAElasticModel.cc.

517{
518
519 if (verboseLevel > 3){
520 G4cout <<
521 "Calling SampleSecondaries() of G4DNAELSEPAElasticModel"
522 << G4endl;
523 }
524
525 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
526
527 const G4Material* material = couple->GetMaterial();
528 const G4ElementVector* theElementVector = material->GetElementVector();
529 std::size_t nelm = material->GetNumberOfElements();
530 if (nelm==1) // Protection: only for single element
531 {
532 G4int Z = G4lrint((*theElementVector)[0]->GetZ());
533 if (Z!=79) return;
534 if (electronEnergy0 < fkillBelowEnergy_Au)
535 {
540 return;
541 }
542
543 if(electronEnergy0>= fkillBelowEnergy_Au && electronEnergy0 < fhighEnergyLimit)
544 {
545 G4double cosTheta = 0;
546 if (electronEnergy0>=10*eV)
547 {
548 cosTheta = RandomizeCosTheta(Z,electronEnergy0);
549 }
550 else
551 {
552 cosTheta = RandomizeCosTheta(Z,10*eV);
553 }
554
555 G4double phi = 2. * CLHEP::pi * G4UniformRand();
556
557 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
558 G4ThreeVector xVers = zVers.orthogonal();
559 G4ThreeVector yVers = zVers.cross(xVers);
560
561 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
562 G4double yDir = xDir;
563 xDir *= std::cos(phi);
564 yDir *= std::sin(phi);
565
566 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
569
570 }
571 }
572 else
573 {
574 if(material->GetName()=="G4_WATER")
575 {
576 //The data for water is stored as Z=0
577 G4double cosTheta = RandomizeCosTheta(0,electronEnergy0);
578
579 G4double phi = 2. * pi * G4UniformRand();
580
581 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
582 G4ThreeVector xVers = zVers.orthogonal();
583 G4ThreeVector yVers = zVers.cross(xVers);
584
585 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
586 G4double yDir = xDir;
587 xDir *= std::cos(phi);
588 yDir *= std::sin(phi);
589
590 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
593 }
594 }
595}
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
const G4double pi

◆ SetKillBelowThreshold()

void G4DNAELSEPAElasticModel::SetKillBelowThreshold ( G4double  threshold)

Definition at line 819 of file G4DNAELSEPAElasticModel.cc.

820{
821 fkillBelowEnergy_Au = threshold;
822
823 if (threshold < 10 * eV)
824 {
825 G4cout<< "*** WARNING : the G4DNAELSEPAElasticModel model is not "
826 "defined below 10 eV !" << G4endl;
827 }
828}

◆ SetMaximumEnergy()

void G4DNAELSEPAElasticModel::SetMaximumEnergy ( G4double  input)
inline

Definition at line 74 of file G4DNAELSEPAElasticModel.hh.

75 {fhighEnergyLimit = input; SetHighEnergyLimit(input);};

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAELSEPAElasticModel::fParticleChangeForGamma
protected

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