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

#include <G4DNACPA100ElasticModel.hh>

+ Inheritance diagram for G4DNACPA100ElasticModel:

Public Member Functions

 G4DNACPA100ElasticModel (const G4ParticleDefinition *p=0, const G4String &nam="DNACPA100ElasticModel")
 
virtual ~G4DNACPA100ElasticModel ()
 
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 SelectStationary (G4bool input)
 
- 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 54 of file G4DNACPA100ElasticModel.hh.

Constructor & Destructor Documentation

◆ G4DNACPA100ElasticModel()

G4DNACPA100ElasticModel::G4DNACPA100ElasticModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNACPA100ElasticModel" 
)

Definition at line 54 of file G4DNACPA100ElasticModel.cc.

56:G4VEmModel(nam),isInitialised(false)
57{
58
59 SetLowEnergyLimit(11*eV);
60 SetHighEnergyLimit(255955*eV);
61
62 verboseLevel= 0;
63 // Verbosity scale:
64 // 0 = nothing
65 // 1 = warning for energy non-conservation
66 // 2 = details of energy budget
67 // 3 = calculation of cross sections, file openings, sampling of atoms
68 // 4 = entering in methods
69
70#ifdef UEHARA_VERBOSE
71 if( verboseLevel>0 )
72 {
73 G4cout << "CPA100 Elastic model is constructed " << G4endl
74 << "Energy range: "
75 << LowEnergyLimit()/eV << " eV - "
76 << HighEnergyLimit()/ keV << " keV"
77 << G4endl;
78 }
79#endif
80
82 fpMolWaterDensity = 0;
83
84 // Selection of stationary mode
85
86 statCode = false;
87}
#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

◆ ~G4DNACPA100ElasticModel()

G4DNACPA100ElasticModel::~G4DNACPA100ElasticModel ( )
virtual

Definition at line 91 of file G4DNACPA100ElasticModel.cc.

92{
93 // For total cross section
94
95 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
96 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
97 {
98 G4DNACrossSectionDataSet* table = pos->second;
99 delete table;
100 }
101
102 // For final state
103 eVecm.clear();
104
105}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 257 of file G4DNACPA100ElasticModel.cc.

263{
264
265#ifdef UEHARA_VERBOSE
266 if (verboseLevel > 3)
267 G4cout <<
268 "Calling CrossSectionPerVolume() of G4DNACPA100ElasticModel" << G4endl;
269#endif
270
271 // Calculate total cross section for model
272
273 G4double sigma=0;
274
275 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
276
277 const G4String& particleName = p->GetParticleName();
278
279 if (ekin <= HighEnergyLimit() && ekin >= LowEnergyLimit())
280 {
281 //SI : XS must not be zero otherwise sampling of secondaries
282 // method ignored
283
284 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
285 pos = tableData.find(particleName);
286
287 if (pos != tableData.end())
288 {
289 G4DNACrossSectionDataSet* table = pos->second;
290 if (table != 0)
291 {
292 sigma = table->FindValue(ekin);
293 //
294 //Dump in non-MT mode
295 //
296 /*
297 G4double minEnergy = 10.481 * eV;
298 G4double maxEnergy = 255955. * eV;
299 G4int nEnergySteps = 1000;
300 G4double energy(minEnergy);
301 G4double
302 stpEnergy(std::pow(maxEnergy/energy,
303 1./static_cast<G4double>(nEnergySteps-1)));
304 G4int step(nEnergySteps);
305 system ("rm -rf elastic-cpa100.out");
306 FILE* myFile=fopen("elastic-cpa100.out","a");
307 while (step>0)
308 {
309 step--;
310 fprintf (myFile,"%16.9le %16.9le\n",
311 energy/eV,
312 table->FindValue(energy)/(1e-20*m*m));
313 energy*=stpEnergy;
314 }
315 fclose (myFile);
316 abort();
317 */
318 //
319 // end of dump
320 //
321 }
322 }
323 else
324 {
325 G4Exception("G4DNACPA100ElasticModel::ComputeCrossSectionPerVolume",
326 "em0002",
327 FatalException,"Model not applicable to particle type.");
328 }
329 }
330
331#ifdef UEHARA_VERBOSE
332 if (verboseLevel > 2)
333 {
334 G4cout << "__________________________________" << G4endl;
335 G4cout << "G4DNACPA100ElasticModel - XS INFO START" << G4endl;
336 G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
337 G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
338 G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
339 // G4cout << " - Cross section per water molecule (cm^-1)="
340 // << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
341 G4cout << "G4DNACPA100ElasticModel - XS INFO END" << G4endl;
342 }
343#endif
344
345 return sigma*waterDensity;
346}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
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

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 109 of file G4DNACPA100ElasticModel.cc.

112{
113
114#ifdef UEHARA_VERBOSE
115 if (verboseLevel > 3)
116 G4cout << "Calling G4DNACPA100ElasticModel::Initialise()" << G4endl;
117#endif
118
119 if(particle->GetParticleName() != "e-")
120 {
121 G4Exception("*** WARNING: the G4DNACPA100ElasticModel is "
122 "not intented to be used with another particle than the electron",
123 "",FatalException,"") ;
124 }
125
126 // Energy limits
127
128 if (LowEnergyLimit() < 11.*eV)
129 {
130 G4cout << "G4DNACPA100ElasticModel: low energy limit increased from " <<
131 LowEnergyLimit()/eV << " eV to " << 11 << " eV" << G4endl;
132 SetLowEnergyLimit(11.*eV);
133 }
134
135 if (HighEnergyLimit() > 255955.*eV)
136 {
137 G4cout << "G4DNACPA100ElasticModel: high energy limit decreased from " <<
138 HighEnergyLimit()/keV << " keV to " << 255.955 << " keV"
139 << G4endl;
140 SetHighEnergyLimit(255955.*eV);
141 }
142
143 // Reading of data files
144
145 G4double scaleFactor = 1e-20*m*m;
146
147 G4String fileElectron("dna/sigma_elastic_e_cpa100");
148
151
152 // *** ELECTRON
153
154 // For total cross section
155
156 electron = electronDef->GetParticleName();
157
158 tableFile[electron] = fileElectron;
159
160 G4DNACrossSectionDataSet* tableE =
162 eV,scaleFactor );
163
164 /*
165 G4DNACrossSectionDataSet* tableE =
166 new G4DNACrossSectionDataSet(new G4DNACPA100LogLogInterpolation,
167 eV,scaleFactor );
168 */
169
170 tableE->LoadData(fileElectron);
171
172 tableData[electron] = tableE;
173
174 // For final state
175
176 char *path = getenv("G4LEDATA");
177
178 if (!path)
179 {
180 G4Exception("G4DNACPA100ElasticModel::Initialise","em0006",
181 FatalException,"G4LEDATA environment variable not set.");
182 return;
183 }
184
185 std::ostringstream eFullFileName;
186
187 eFullFileName << path
188 << "/dna/sigmadiff_cumulated_elastic_e_cpa100.dat";
189
190 std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
191
192 if (!eDiffCrossSection)
193 G4Exception("G4DNACPA100ElasticModel::Initialise","em0003",
195 "Missing data file:/dna/sigmadiff_cumulated_elastic_e_cpa100.dat");
196
197 // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
198 // Added clear for MT
199
200 eTdummyVec.clear();
201 eVecm.clear();
202 eDiffCrossSectionData.clear();
203
204 //
205
206 eTdummyVec.push_back(0.);
207
208 while(!eDiffCrossSection.eof())
209 {
210 G4double tDummy;
211 G4double eDummy;
212 eDiffCrossSection>>tDummy>>eDummy;
213
214 // SI : mandatory eVecm initialization
215
216 if (tDummy != eTdummyVec.back())
217 {
218 eTdummyVec.push_back(tDummy);
219 eVecm[tDummy].push_back(0.);
220 }
221
222 eDiffCrossSection>>eDiffCrossSectionData[tDummy][eDummy];
223
224 if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy);
225 }
226
227 // End final state
228
229#ifdef UEHARA_VERBOSE
230 if (verboseLevel > 2)
231 G4cout << "Loaded cross section files for CPA100 Elastic model" << G4endl;
232#endif
233
234#ifdef UEHARA_VERBOSE
235 if( verboseLevel>0 )
236 {
237 G4cout << "CPA100 Elastic model is initialized " << G4endl
238 << "Energy range: "
239 << LowEnergyLimit() / eV << " eV - "
240 << HighEnergyLimit() / keV << " keV"
241 << G4endl;
242 }
243#endif
244
245 // Initialize water density pointer
246 fpMolWaterDensity = G4DNAMolecularMaterial::Instance()
248
249 if (isInitialised) { return; }
251 isInitialised = true;
252
253}
virtual G4bool LoadData(const G4String &argFileName)
const std::vector< G4double > * GetNumMolPerVolTableFor(const G4Material *) const
Retrieve a table of molecular densities (number of molecules per unit volume) in the G4 unit system f...
static G4DNAMolecularMaterial * Instance()
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:651
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 350 of file G4DNACPA100ElasticModel.cc.

355{
356#ifdef UEHARA_VERBOSE
357 if (verboseLevel > 3)
358 G4cout << "Calling SampleSecondaries() of G4DNACPA100ElasticModel" << G4endl;
359#endif
360
361 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
362
363 G4double cosTheta = RandomizeCosTheta(electronEnergy0);
364 G4double phi = 2. * pi * G4UniformRand();
365
366 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
367
368 //G4ThreeVector xVers = zVers.orthogonal();
369 //G4ThreeVector yVers = zVers.cross(xVers);
370 //G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
371 //G4double yDir = xDir;
372 //xDir *= std::cos(phi);
373 //yDir *= std::sin(phi);
374
375 // Computation of scattering angles (from Subroutine DIRAN in CPA100)
376
377 G4double CT1, ST1, CF1, SF1, CT2, ST2, CF2, SF2;
378 G4double sinTheta = std::sqrt (1-cosTheta*cosTheta);
379
380 CT1=0;
381 ST1=0;
382 CF1=0;
383 SF1=0;
384 CT2=0;
385 ST2=0;
386 CF2=0;
387 SF2=0;
388
389 CT1 = zVers.z();
390 ST1=std::sqrt(1.-CT1*CT1);
391
392 if (ST1!=0) CF1 = zVers.x()/ST1; else CF1 = std::cos(2. * pi * G4UniformRand());
393 if (ST1!=0) SF1 = zVers.y()/ST1; else SF1 = std::sqrt(1.-CF1*CF1);
394
395 G4double A3, A4, A5, A2, A1;
396 A3=0;
397 A4=0;
398 A5=0;
399 A2=0;
400 A1=0;
401
402 A3 = sinTheta*std::cos(phi);
403 A4 = A3*CT1 + ST1*cosTheta;
404 A5 = sinTheta * std::sin(phi);
405 A2 = A4 * SF1 + A5 * CF1;
406 A1 = A4 * CF1 - A5 * SF1;
407
408 CT2 = CT1*cosTheta - ST1*A3;
409 ST2 = std::sqrt(1.-CT2*CT2);
410
411 if (ST2==0) ST2=1E-6;
412 CF2 = A1/ST2;
413 SF2 = A2/ST2;
414
415 /*
416 G4cout << "CT1=" << CT1 << G4endl;
417 G4cout << "ST1=" << ST1 << G4endl;
418 G4cout << "CF1=" << CF1 << G4endl;
419 G4cout << "SF1=" << SF1 << G4endl;
420 G4cout << "cosTheta=" << cosTheta << G4endl;
421 G4cout << "sinTheta=" << sinTheta << G4endl;
422 G4cout << "cosPhi=" << std::cos(phi) << G4endl;
423 G4cout << "sinPhi=" << std::sin(phi) << G4endl;
424 G4cout << "CT2=" << CT2 << G4endl;
425 G4cout << "ST2=" << ST2 << G4endl;
426 G4cout << "CF2=" << CF2 << G4endl;
427 G4cout << "SF2=" << SF2 << G4endl;
428 */
429
430 G4ThreeVector zPrimeVers(ST2*CF2,ST2*SF2,CT2);
431
432 //
433
435
436 if (!statCode)
437
439 (electronEnergy0-1.214E-4*(1.-cosTheta)*electronEnergy0);
440
442
443 //
444
445 fParticleChangeForGamma->ProposeLocalEnergyDeposit(1.214E-4*(1.-cosTheta)*electronEnergy0);
446
447}
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
double x() const
double y() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
const G4double pi

◆ SelectStationary()

void G4DNACPA100ElasticModel::SelectStationary ( G4bool  input)
inline

Definition at line 146 of file G4DNACPA100ElasticModel.hh.

147{
148 statCode = input;
149}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNACPA100ElasticModel::fParticleChangeForGamma
protected

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