Geant4 11.3.0
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=nullptr, const G4String &nam="DNAELSEPAElasticModel")
 
 ~G4DNAELSEPAElasticModel () override
 
G4DNAELSEPAElasticModeloperator= (const G4DNAELSEPAElasticModel &right)=delete
 
 G4DNAELSEPAElasticModel (const G4DNAELSEPAElasticModel &)=delete
 
void Initialise (const G4ParticleDefinition *particle, const G4DataVector &) override
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *particle, G4double ekin, G4double emin, G4double emax) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SetMaximumEnergy (G4double input)
 
void SetKillBelowThreshold (G4double threshold)
 
G4double GetKillBelowThreshold ()
 
- 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
 

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 49 of file G4DNAELSEPAElasticModel.hh.

Constructor & Destructor Documentation

◆ G4DNAELSEPAElasticModel() [1/2]

G4DNAELSEPAElasticModel::G4DNAELSEPAElasticModel ( const G4ParticleDefinition * particle = nullptr,
const G4String & nam = "DNAELSEPAElasticModel" )

Definition at line 49 of file G4DNAELSEPAElasticModel.cc.

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

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

◆ ~G4DNAELSEPAElasticModel()

G4DNAELSEPAElasticModel::~G4DNAELSEPAElasticModel ( )
override

Definition at line 117 of file G4DNAELSEPAElasticModel.cc.

118{
119 delete fpData_Au;
120 delete fpData_H2O;
121
122 eEdummyVec_Au.clear();
123 eEdummyVec_H2O.clear();
124 eCum_Au.clear();
125 eCum_H2O.clear();
126 fAngleData_Au.clear();
127 fAngleData_H2O.clear();
128}

◆ G4DNAELSEPAElasticModel() [2/2]

G4DNAELSEPAElasticModel::G4DNAELSEPAElasticModel ( const G4DNAELSEPAElasticModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 336 of file G4DNAELSEPAElasticModel.cc.

342{
343
344 if (verboseLevel > 3)
345 {
346 G4cout <<
347 "Calling CrossSectionPerVolume() of G4DNAELSEPAElasticModel"
348 << G4endl;
349 }
350
351 G4double atomicNDensity=0.0;
352 G4double sigma=0;
353
354 std::size_t nelm = material->GetNumberOfElements();
355 if (nelm==1) // Protection: only for single element
356 {
357 // Protection: only for GOLD
358 if (material->GetZ()!=79) return 0.0;
359
360 const G4ElementVector* theElementVector = material->GetElementVector();
361 G4int Z = G4lrint((*theElementVector)[0]->GetZ());
362
363 const G4String& particleName = particle->GetParticleName();
364 atomicNDensity = material->GetAtomicNumDensityVector()[0];
365 if(atomicNDensity!= 0.0)
366 {
367 if (ekin < fhighEnergyLimit)
368 {
369 if (ekin < fkillBelowEnergy_Au) return DBL_MAX;
370
371 if (ekin < 10*eV) sigma = fpData_Au->FindValue(10*eV);
372 else sigma = fpData_Au->FindValue(ekin);
373 }
374 }
375 if (verboseLevel > 2)
376 {
377 G4cout << "__________________________________" << G4endl;
378 G4cout << "=== G4DNAELSEPAElasticModel - XS INFO START" << G4endl;
379 G4cout << "=== Material is made of one element with Z =" << Z << G4endl;
380 G4cout << "=== Kinetic energy(eV)=" << ekin/eV << " particle : "
381 << particleName << G4endl;
382 G4cout << "=== Cross section per atom for Z="<<Z<<" is (cm^2)"
383 << sigma/cm/cm << G4endl;
384 G4cout << "=== Cross section per atom for Z="<<Z<<" is (cm^-1)="
385 << sigma*atomicNDensity/(1./cm) << G4endl;
386 G4cout << "=== G4DNAELSEPAElasticModel - XS INFO END" << G4endl;
387 }
388 }
389 else
390 {
391 atomicNDensity = (*fpMolDensity)[material->GetIndex()];
392 if(atomicNDensity!= 0.0)
393 {
394 if (ekin < HighEnergyLimit() && ekin >= LowEnergyLimit())
395 {
396 sigma = fpData_H2O->FindValue(ekin);
397 }
398 }
399 if (verboseLevel > 2)
400 {
401 G4cout << "__________________________________" << G4endl;
402 G4cout << "=== G4DNAELSEPAElasticModel - XS INFO START" << G4endl;
403 G4cout << "=== Kinetic energy(eV)=" << ekin/eV
404 << " particle : " << particle->GetParticleName() << G4endl;
405 G4cout << "=== Cross section per water molecule (cm^2)="
406 << sigma/cm/cm << G4endl;
407 G4cout << "=== Cross section per water molecule (cm^-1)="
408 << sigma*atomicNDensity/(1./cm) << G4endl;
409 G4cout << "=== G4DNAELSEPAElasticModel - XS INFO END" << G4endl;
410 }
411 }
412
413 return sigma*atomicNDensity;
414}
double G4double
Definition G4Types.hh:83
G4double GetZ() const
const G4double * GetAtomicNumDensityVector() const
std::size_t GetIndex() const
const G4String & GetParticleName() const
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
#define DBL_MAX
Definition templates.hh:62

◆ GetKillBelowThreshold()

G4double G4DNAELSEPAElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 82 of file G4DNAELSEPAElasticModel.hh.

82{return fkillBelowEnergy_Au;}

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 132 of file G4DNAELSEPAElasticModel.cc.

134{
135 if (verboseLevel > 3)
136 G4cout << "Calling G4DNAELSEPAElasticModel::Initialise()" << G4endl;
137
138 if (isInitialised) {return;}
139
140 if(particle->GetParticleName() != "e-")
141 {
142 G4Exception("G4DNAELSEPAElasticModel::Initialise","em0001",
143 FatalException,"Model not applicable to particle type.");
144 return;
145 }
146
147 G4ProductionCutsTable* theCoupleTable =
149 auto numOfCouples = (G4int)theCoupleTable->GetTableSize();
150
151 // UNIT OF TCS
152 G4double scaleFactor = 1.*cm*cm;
153
154 fpData_Au=nullptr;
155 fpData_H2O=nullptr;
156 fpBaseWater = G4Material::GetMaterial("G4_WATER");
157
158 for(G4int i=0; i<numOfCouples; ++i)
159 {
160 const G4MaterialCutsCouple* couple =
161 theCoupleTable->GetMaterialCutsCouple(i);
162 const G4Material* material = couple->GetMaterial()->GetBaseMaterial();
163 if(!material) material = couple->GetMaterial();
164
165 auto nelm = (G4int)material->GetNumberOfElements();
166 if (nelm==1){// Protection: only for single element
167 const G4ElementVector* theElementVector = material->GetElementVector();
168 G4int Z = G4lrint((*theElementVector)[0]->GetZ());
169 if (Z!=79)// Protection: only for GOLD
170 {
171 continue;
172 }
173
174 if (Z>0)
175 {
176 G4String fileZElectron("dna/sigma_elastic_e_elsepa_Z");
177 std::ostringstream oss;
178 oss.str("");
179 oss.clear(stringstream::goodbit);
180 oss << Z;
181 fileZElectron += oss.str()+"_muffintin";
182
183 fpData_Au = new G4DNACrossSectionDataSet(new G4LogLogInterpolation,
184 eV,
185 scaleFactor );
186 fpData_Au->LoadData(fileZElectron);
187
188 std::ostringstream eFullFileNameZ;
189 const char *path = G4EmParameters::Instance()->GetDirLEDATA();
190
191 if (path == nullptr)
192 {
193 G4Exception("G4DNAELSEPAElasticModel::Initialise","em0002",
194 FatalException,"G4LEDATA environment variable not set.");
195 return;
196 }
197
198 eFullFileNameZ.str("");
199 eFullFileNameZ.clear(stringstream::goodbit);
200
201 eFullFileNameZ
202 << path
203 << "/dna/sigmadiff_cumulated_elastic_e_elsepa_Z"
204 << Z << "_muffintin.dat";
205
206 std::ifstream eDiffCrossSectionZ(eFullFileNameZ.str().c_str());
207
208 if (!eDiffCrossSectionZ)
209 {
210 G4Exception("G4DNAELSEPAElasticModel::Initialise","em0003",
211 FatalException,"Missing data file for cumulated DCS");
212 return;
213 }
214
215 eEdummyVec_Au.clear();
216 eCum_Au.clear();
217 fAngleData_Au.clear();
218
219 eEdummyVec_Au.push_back(0.);
220 do
221 {
222 G4double eDummy;
223 G4double cumDummy;
224 eDiffCrossSectionZ>>eDummy>>cumDummy;
225 if (eDummy != eEdummyVec_Au.back())
226 {
227 eEdummyVec_Au.push_back(eDummy);
228 eCum_Au[eDummy].push_back(0.);
229 }
230 eDiffCrossSectionZ>>fAngleData_Au[eDummy][cumDummy];
231 if (cumDummy != eCum_Au[eDummy].back())
232 {
233 eCum_Au[eDummy].push_back(cumDummy);
234 }
235 }while(!eDiffCrossSectionZ.eof());
236 }
237
238 }else{// Protection: H2O only is available
239 if(material == fpBaseWater && !fpData_H2O){
240 if (LowEnergyLimit() < 10*eV)
241 {
242 G4cout<<"G4DNAELSEPAElasticModel: low energy limit increased from "
243 << LowEnergyLimit()/eV << " eV to " << 10 << " eV"
244 << G4endl;
245 SetLowEnergyLimit(10.*eV);
246 }
247
248 if (HighEnergyLimit() > 1.*MeV)
249 {
250 G4cout<<"G4DNAELSEPAElasticModel: high energy limit decreased from "
251 << HighEnergyLimit()/MeV << " MeV to " << 1. << " MeV"
252 << G4endl;
253 SetHighEnergyLimit(1.*MeV);
254 }
255
256 G4String fileZElectron("dna/sigma_elastic_e_elsepa_muffin");
257
258 fpData_H2O = new G4DNACrossSectionDataSet(new G4LogLogInterpolation,
259 eV,
260 scaleFactor );
261 fpData_H2O->LoadData(fileZElectron);
262
263 std::ostringstream eFullFileNameZ;
264
265 const char *path = G4EmParameters::Instance()->GetDirLEDATA();
266 if (path == nullptr)
267 {
268 G4Exception("G4DNAELSEPAElasticModel::Initialise","em0004",
269 FatalException,"G4LEDATA environment variable not set.");
270 return;
271 }
272
273 eFullFileNameZ.str("");
274 eFullFileNameZ.clear(stringstream::goodbit);
275
276 eFullFileNameZ
277 << path
278 << "/dna/sigmadiff_cumulated_elastic_e_elsepa_muffin.dat";
279
280 std::ifstream eDiffCrossSectionZ(eFullFileNameZ.str().c_str());
281
282 if (!eDiffCrossSectionZ)
283 G4Exception("G4DNAELSEPAElasticModel::Initialise","em0005",
285 "Missing data file for cumulated DCS");
286
287 eEdummyVec_H2O.clear();
288 eCum_H2O.clear();
289 fAngleData_H2O.clear();
290
291 eEdummyVec_H2O.push_back(0.);
292
293 do
294 {
295 G4double eDummy;
296 G4double cumDummy;
297 eDiffCrossSectionZ>>eDummy>>cumDummy;
298 if (eDummy != eEdummyVec_H2O.back())
299 {
300 eEdummyVec_H2O.push_back(eDummy);
301 eCum_H2O[eDummy].push_back(0.);
302 }
303 eDiffCrossSectionZ>>fAngleData_H2O[eDummy][cumDummy];
304 if (cumDummy != eCum_H2O[eDummy].back()){
305 eCum_H2O[eDummy].push_back(cumDummy);
306 }
307 }while(!eDiffCrossSectionZ.eof());
308 }
309 }
310 if (verboseLevel > 2)
311 G4cout << "Loaded cross section files of ELSEPA Elastic model for"
312 << material->GetName() << G4endl;
313
314 if( verboseLevel>0 )
315 {
316 G4cout << "ELSEPA elastic model is initialized " << G4endl
317 << "Energy range: "
318 << LowEnergyLimit() / eV << " eV - "
319 << HighEnergyLimit()/ MeV << " MeV"
320 << G4endl;
321 }
322 } // Loop on couples
323
324
326
327 fpMolDensity =
329 GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
330
331 isInitialised = true;
332}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
static G4DNAMolecularMaterial * Instance()
static G4EmParameters * Instance()
const G4String & GetDirLEDATA() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 418 of file G4DNAELSEPAElasticModel.cc.

424{
425
426 if (verboseLevel > 3){
427 G4cout <<
428 "Calling SampleSecondaries() of G4DNAELSEPAElasticModel"
429 << G4endl;
430 }
431
432 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
433
434 const G4Material* material = couple->GetMaterial()->GetBaseMaterial();
435 if(!material) material = couple->GetMaterial();
436
437 std::size_t nelm = material->GetNumberOfElements();
438 if (nelm==1) // Protection: only for single element
439 {
440 const G4ElementVector* theElementVector = material->GetElementVector();
441 G4int Z = G4lrint((*theElementVector)[0]->GetZ());
442 if (Z!=79) return;
443 if (electronEnergy0 < fkillBelowEnergy_Au)
444 {
445 fParticleChangeForGamma->SetProposedKineticEnergy(0.);
446 fParticleChangeForGamma->ProposeMomentumDirection(G4ThreeVector(0,0,0));
447 fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
448 fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
449 return;
450 }
451
452 if(electronEnergy0>= fkillBelowEnergy_Au && electronEnergy0 < fhighEnergyLimit)
453 {
454 G4double cosTheta = 0;
455 if (electronEnergy0>=10*eV)
456 {
457 cosTheta = RandomizeCosTheta(Z,electronEnergy0);
458 }
459 else
460 {
461 cosTheta = RandomizeCosTheta(Z,10*eV);
462 }
463
464 G4double phi = 2. * CLHEP::pi * G4UniformRand();
465
466 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
467 G4ThreeVector xVers = zVers.orthogonal();
468 G4ThreeVector yVers = zVers.cross(xVers);
469
470 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
471 G4double yDir = xDir;
472 xDir *= std::cos(phi);
473 yDir *= std::sin(phi);
474
475 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
476 fParticleChangeForGamma->ProposeMomentumDirection(zPrimeVers.unit());
477 fParticleChangeForGamma->SetProposedKineticEnergy(electronEnergy0);
478
479 }
480 }
481 else
482 {
483 if(material == fpBaseWater)
484 {
485 //The data for water is stored as Z=0
486 G4double cosTheta = RandomizeCosTheta(0,electronEnergy0);
487
488 G4double phi = 2. * pi * G4UniformRand();
489
490 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
491 G4ThreeVector xVers = zVers.orthogonal();
492 G4ThreeVector yVers = zVers.cross(xVers);
493
494 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
495 G4double yDir = xDir;
496 xDir *= std::cos(phi);
497 yDir *= std::sin(phi);
498
499 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
500 fParticleChangeForGamma->ProposeMomentumDirection(zPrimeVers.unit());
501 fParticleChangeForGamma->SetProposedKineticEnergy(electronEnergy0);
502 }
503 }
504}
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
const G4double pi

◆ SetKillBelowThreshold()

void G4DNAELSEPAElasticModel::SetKillBelowThreshold ( G4double threshold)

Definition at line 718 of file G4DNAELSEPAElasticModel.cc.

719{
720 fkillBelowEnergy_Au = threshold;
721
722 if (threshold < 10 * eV)
723 {
724 G4cout<< "*** WARNING : the G4DNAELSEPAElasticModel model is not "
725 "defined below 10 eV !" << G4endl;
726 }
727}

◆ SetMaximumEnergy()

void G4DNAELSEPAElasticModel::SetMaximumEnergy ( G4double input)
inline

Definition at line 77 of file G4DNAELSEPAElasticModel.hh.

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

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAELSEPAElasticModel::fParticleChangeForGamma
protected

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