Geant4 10.7.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 *p=0, const G4String &nam="DNAELSEPAElasticModel")
 
virtual ~G4DNAELSEPAElasticModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
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 *, G4double &eloss, G4double &niel, G4double length)
 
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 ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

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

Constructor & Destructor Documentation

◆ G4DNAELSEPAElasticModel()

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

Definition at line 43 of file G4DNAELSEPAElasticModel.cc.

44 :
45 G4VEmModel(nam), isInitialised(false)
46{
47 SetLowEnergyLimit(10. * eV);
48 SetHighEnergyLimit(1. * MeV);
49
50 verboseLevel = 0;
51 // Verbosity scale:
52 // 0 = nothing
53 // 1 = warning for energy non-conservation
54 // 2 = details of energy budget
55 // 3 = calculation of cross sections, file openings, sampling of atoms
56 // 4 = entering in methods
57
58#ifdef ELSEPA_VERBOSE
59 if (verboseLevel > 0)
60 {
61 G4cout << "ELSEPA Elastic model is constructed "
62 << G4endl
63 << "Energy range: "
64 << LowEnergyLimit() / eV << " eV - "
65 << HighEnergyLimit() / MeV << " MeV"
66 << G4endl;
67 }
68#endif
69
71 fpMolWaterDensity = 0;
72 fpData = 0;
73}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764

◆ ~G4DNAELSEPAElasticModel()

G4DNAELSEPAElasticModel::~G4DNAELSEPAElasticModel ( )
virtual

Definition at line 77 of file G4DNAELSEPAElasticModel.cc.

78{
79 // For total cross section
80 if(fpData) delete fpData;
81
82 // For final state
83 eVecm.clear();
84}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 235 of file G4DNAELSEPAElasticModel.cc.

245{
246#ifdef ELSEPA_VERBOSE
247 if (verboseLevel > 3)
248 {
249 G4cout << "Calling CrossSectionPerVolume() of G4DNAELSEPAElasticModel"
250 << G4endl;
251 }
252#endif
253
254 // Calculate total cross section for model
255
256 G4double sigma = 0.;
257 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
258
259 if(waterDensity!= 0.0)
260 {
261 if (ekin < HighEnergyLimit() && ekin >= LowEnergyLimit())
262 {
263 //SI : XS must not be zero otherwise sampling of secondaries method ignored
264 //
265 sigma = fpData->FindValue(ekin);
266 }
267
268#ifdef ELSEPA_VERBOSE
269 if (verboseLevel > 2)
270 {
271 G4cout << "__________________________________" << G4endl;
272 G4cout << "=== G4DNAELSEPAElasticModel - XS INFO START" << G4endl;
273 G4cout << "=== Kinetic energy(eV)=" << ekin/eV << " particle : " << p->GetParticleName() << G4endl;
274 G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
275 G4cout << "=== Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
276 G4cout << "=== G4DNAELSEPAElasticModel - XS INFO END" << G4endl;
277 }
278#endif
279 }
280
281 return sigma*waterDensity;
282}
double G4double
Definition: G4Types.hh:83
virtual G4double FindValue(G4double e, G4int componentId=0) const
size_t GetIndex() const
Definition: G4Material.hh:258
const G4String & GetParticleName() const

◆ GetKillBelowThreshold()

G4double G4DNAELSEPAElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 68 of file G4DNAELSEPAElasticModel.hh.

69 {
71 errMsg << "The method G4DNAELSEPAElasticModel::"
72 "GetKillBelowThreshold is deprecated";
73
74 G4Exception("G4DNAELSEPAElasticModel::GetKillBelowThreshold",
75 "deprecated",
77 errMsg);
78 return 0.;
79 }
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 88 of file G4DNAELSEPAElasticModel.cc.

90{
91#ifdef ELSEPA_VERBOSE
92 if (verboseLevel > 3)
93 {
94 G4cout << "Calling G4DNAELSEPAElasticModel::Initialise()" << G4endl;
95 }
96#endif
97
98 if(particle->GetParticleName() != "e-")
99 {
100 G4Exception("G4DNAELSEPAElasticModel::Initialise",
101 "em0002",
103 "Model not applicable to particle type.");
104 }
105
106 // Energy limits
107
108 if (LowEnergyLimit() < 10*eV)
109 {
110 G4cout << "G4DNAELSEPAElasticModel: low energy limit increased from "
111 << LowEnergyLimit()/eV << " eV to " << 10 << " eV"
112 << G4endl;
113 SetLowEnergyLimit(10.*eV);
114 }
115
116 if (HighEnergyLimit() > 1.*MeV)
117 {
118 G4cout << "G4DNAELSEPAElasticModel: high energy limit decreased from "
119 << HighEnergyLimit()/MeV << " MeV to " << 1. << " MeV"
120 << G4endl;
121 SetHighEnergyLimit(1.*MeV);
122 }
123
124 if (isInitialised) { return; }
125
126 // *** ELECTRON
127 // For total cross section
128 // Reading of data files
129
130 G4double scaleFactor = 1*cm*cm;
131
132 G4String fileElectron("dna/sigma_elastic_e_elsepa_muffin");
133
134// Alternative option
135// G4String fileElectron("dna/sigma_elastic_e_elsepa_free");
136
138 eV,
139 scaleFactor );
140 fpData->LoadData(fileElectron);
141 // For final state
142
143 char *path = getenv("G4LEDATA");
144
145 if (!path)
146 {
147 G4Exception("G4ELSEPAElasticModel::Initialise",
148 "em0006",
150 "G4LEDATA environment variable not set.");
151 return;
152 }
153
154 std::ostringstream eFullFileName;
155
156// Alternative option
157// eFullFileName << path << "/dna/sigmadiff_cumulated_elastic_e_elsepa_free.dat";
158
159 eFullFileName << path << "/dna/sigmadiff_cumulated_elastic_e_elsepa_muffin.dat";
160 std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
161
162 if (!eDiffCrossSection)
163 {
165 errMsg << "Missing data file:/dna/sigmadiff_cumulated_elastic_e_elsepa_muffin.dat; "
166 << "please use G4EMLOW7.8 and above.";
167
168 G4Exception("G4DNAELSEPAElasticModel::Initialise",
169 "em0003",
171 errMsg);
172 }
173
174 // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
175 // Added clear for MT
176
177 eTdummyVec.clear();
178 eVecm.clear();
179 eDiffCrossSectionData.clear();
180
181 //
182
183 eTdummyVec.push_back(0.);
184
185 while(!eDiffCrossSection.eof())
186 {
187 double tDummy;
188 double eDummy;
189 eDiffCrossSection >> tDummy >> eDummy;
190
191 // SI : mandatory eVecm initialization
192
193 if (tDummy != eTdummyVec.back())
194 {
195 eTdummyVec.push_back(tDummy);
196 eVecm[tDummy].push_back(0.);
197 }
198
199 eDiffCrossSection >> eDiffCrossSectionData[tDummy][eDummy];
200
201 if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy);
202 }
203
204 // End final state
205#ifdef ELSEPA_VERBOSE
206 if (verboseLevel>0)
207 {
208 if (verboseLevel > 2)
209 {
210 G4cout << "Loaded cross section files for ELSEPA Elastic model" << G4endl;
211 }
212
213 G4cout << "ELSEPA Elastic model is initialized " << G4endl
214 << "Energy range: "
215 << LowEnergyLimit() / eV << " eV - "
216 << HighEnergyLimit() / MeV << " MeV"
217 << G4endl;
218 }
219#endif
220
221 // Initialize water density pointer
223
224 fpMolWaterDensity = G4DNAMolecularMaterial::Instance()->
225 GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
226
228 isInitialised = true;
229
230}
@ FatalException
virtual G4bool LoadData(const G4String &argFileName)
static G4DNAMolecularMaterial * Instance()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:651
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 286 of file G4DNAELSEPAElasticModel.cc.

291{
292
293
294#ifdef ELSEPA_VERBOSE
295 if (verboseLevel > 3)
296 {
297 G4cout << "Calling SampleSecondaries() of G4DNAELSEPAElasticModel" << G4endl;
298 }
299#endif
300
301 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
302
303// if (electronEnergy0 < HighEnergyLimit()) // necessaire ?
304 {
305 G4double cosTheta = RandomizeCosTheta(electronEnergy0);
306
307 G4double phi = 2. * pi * G4UniformRand();
308
309 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
310 G4ThreeVector xVers = zVers.orthogonal();
311 G4ThreeVector yVers = zVers.cross(xVers);
312
313 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
314 G4double yDir = xDir;
315 xDir *= std::cos(phi);
316 yDir *= std::sin(phi);
317
318 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
319
321
323 // necessaire ?
324 }
325}
#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(G4double Px, G4double Py, G4double Pz)
const G4double pi

◆ SetKillBelowThreshold()

void G4DNAELSEPAElasticModel::SetKillBelowThreshold ( G4double  threshold)

Definition at line 484 of file G4DNAELSEPAElasticModel.cc.

485{
487 errMsg << "The method G4DNAELSEPAElasticModel::SetKillBelowThreshold is deprecated";
488
489 G4Exception("G4DNAELSEPAElasticModel::SetKillBelowThreshold",
490 "deprecated",
492 errMsg);
493}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAELSEPAElasticModel::fParticleChangeForGamma
protected

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