Geant4  10.02
G4DNAIonElasticModel.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // Author: H. N. Tran (Ton Duc Thang University)
27 // p, H, He, He+ and He++ models are assumed identical
28 // NIMB 343, 132-137 (2015)
29 //
30 // The Geant4-DNA web site is available at http://geant4-dna.org
31 //
32 
33 #include "G4DNAIonElasticModel.hh"
34 #include "G4PhysicalConstants.hh"
35 #include "G4SystemOfUnits.hh"
37 #include "G4ParticleTable.hh"
38 
39 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
40 
41 using namespace std;
42 
43 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
44 
46  const G4String& nam) :
47  G4VEmModel(nam), isInitialised(false)
48 {
49  //nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
50 
51  killBelowEnergy = 100 * eV;
52  lowEnergyLimit = 0 * eV;
53  highEnergyLimit = 1 * MeV;
56 
57  verboseLevel = 0;
58  // Verbosity scale:
59  // 0 = nothing
60  // 1 = warning for energy non-conservation
61  // 2 = details of energy budget
62  // 3 = calculation of cross sections, file openings, sampling of atoms
63  // 4 = entering in methods
64 
65  if(verboseLevel > 0)
66  {
67  G4cout << "Ion elastic model is constructed " << G4endl<< "Energy range: "
68  << lowEnergyLimit / eV << " eV - "
69  << highEnergyLimit / MeV << " MeV"
70  << G4endl;
71  }
74  fpTableData = 0;
75  fParticle_Mass = -1;
76 
77  // Selection of stationary mode
78 
79  statCode = false;
80 }
81 
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83 
85 {
86  // For total cross section
87  if(fpTableData) delete fpTableData;
88 }
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91 
92 void
94  const G4ParticleDefinition* particleDefinition,
95  const G4DataVector& /*cuts*/)
96 {
97 
98  if(verboseLevel > 3)
99  {
100  G4cout << "Calling G4DNAIonElasticModel::Initialise()" << G4endl;
101  }
102 
103  // Energy limits
104 
106  {
107  G4cout << "G4DNAIonElasticModel: low energy limit increased from " <<
108  LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
110  }
111 
113  {
114  G4cout << "G4DNAIonElasticModel: high energy limit decreased from " <<
115  HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl;
117  }
118 
119  // Reading of data files
120 
121  G4double scaleFactor = 1e-16*cm*cm;
122 
123  char *path = getenv("G4LEDATA");
124 
125  if (!path)
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 
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 != 0)
147  ||
148  (particleDefinition == hydrogenDef && hydrogenDef != 0)
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)
161  ||
162  (particleDefinition == instance->GetIon("alpha+") && alphaplusDef)
163  ||
164  (particleDefinition == instance->GetIon("alpha++") && alphaplusplusDef)
165  )
166  {
167  fParticle_Mass = 4.;
168  // For total cross section of he,he+,he++
169  totalXSFile = "dna/sigma_elastic_alpha_HTran";
170 
171  // For final state
172  fullFileName << path << "/dna/sigmadiff_cumulated_elastic_alpha_HTran.dat";
173  }
174 
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  double tDummy;
202  double 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
236  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
237 
238  if (isInitialised)
239  { return;}
241  isInitialised = true;
242 }
243 
244 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
245 
246 G4double
248  const G4ParticleDefinition* p,
249  G4double ekin, G4double, G4double)
250 {
251  if(verboseLevel > 3)
252  {
253  G4cout << "Calling CrossSectionPerVolume() of G4DNAIonElasticModel"
254  << G4endl;
255  }
256 
257  // Calculate total cross section for model
258 
259  G4double sigma=0;
260 
261  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
262 
263  if(waterDensity!= 0.0)
264  {
265  const G4String& particleName = p->GetParticleName();
266 
267  if (ekin < highEnergyLimit)
268  {
269  //SI : XS must not be zero otherwise sampling of secondaries method ignored
270  if (ekin < killBelowEnergy) return DBL_MAX;
271  //
272 
273  if (fpTableData != 0) // MK: is this check necessary?
274  {
275  sigma = fpTableData->FindValue(ekin);
276  }
277  else
278  {
279  G4Exception("G4DNAIonElasticModel::ComputeCrossSectionPerVolume","em0002",
280  FatalException,"Model not applicable to particle type.");
281  }
282  }
283 
284  if (verboseLevel > 2)
285  {
286  G4cout << "__________________________________" << G4endl;
287  G4cout << "G4DNAIonElasticModel - XS INFO START" << G4endl;
288  G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
289  G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
290  G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
291  G4cout << "G4DNAIonElasticModel - XS INFO END" << G4endl;
292  }
293 
294  }
295 
296  return sigma*waterDensity;
297 }
298 
299 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
300 
301 void
303  std::vector<G4DynamicParticle*>* /*fvect*/,
304  const G4MaterialCutsCouple* /*couple*/,
305  const G4DynamicParticle* aDynamicParticle, G4double, G4double)
306 {
307 
308  if(verboseLevel > 3)
309  {
310  G4cout << "Calling SampleSecondaries() of G4DNAIonElasticModel" << G4endl;
311  }
312 
313  G4double particleEnergy0 = aDynamicParticle->GetKineticEnergy();
314 
315  if (particleEnergy0 < killBelowEnergy)
316  {
320  return;
321  }
322 
323  if (particleEnergy0>= killBelowEnergy && particleEnergy0 < highEnergyLimit)
324  {
325  G4double water_mass = 18.;
326 
327  G4double thetaCM = RandomizeThetaCM(particleEnergy0, aDynamicParticle->GetDefinition());
328 
329  //HT:convert to laboratory system
330 
331  G4double theta = std::atan(std::sin(thetaCM*CLHEP::pi/180)
332  /(fParticle_Mass/water_mass+std::cos(thetaCM*CLHEP::pi/180)));
333 
334  G4double cosTheta= std::cos(theta);
335 
336  //
337 
338  G4double phi = 2. * CLHEP::pi * G4UniformRand();
339 
340  G4ThreeVector zVers = aDynamicParticle->GetMomentumDirection();
341  G4ThreeVector xVers = zVers.orthogonal();
342  G4ThreeVector yVers = zVers.cross(xVers);
343 
344  G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
345  G4double yDir = xDir;
346  xDir *= std::cos(phi);
347  yDir *= std::sin(phi);
348 
349  G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
350 
352 
353  G4double depositEnergyCM = 0;
354 
355  //HT: deposited energy
356  depositEnergyCM = 4. * particleEnergy0 * fParticle_Mass * water_mass *
357  (1-std::cos(thetaCM*CLHEP::pi/180))
358  / (2 * std::pow((fParticle_Mass+water_mass),2));
359 
360  if (!statCode) fParticleChangeForGamma->SetProposedKineticEnergy(particleEnergy0 - depositEnergyCM);
361 
362  else fParticleChangeForGamma->SetProposedKineticEnergy(particleEnergy0);
363 
365  }
366 }
367 
368 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
369 
370 G4double
372  G4double k, G4double integrDiff)
373 {
374  G4double theta = 0.;
375  G4double valueT1 = 0;
376  G4double valueT2 = 0;
377  G4double valueE21 = 0;
378  G4double valueE22 = 0;
379  G4double valueE12 = 0;
380  G4double valueE11 = 0;
381  G4double xs11 = 0;
382  G4double xs12 = 0;
383  G4double xs21 = 0;
384  G4double xs22 = 0;
385 
386  std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),
387  eTdummyVec.end(), k);
388  std::vector<double>::iterator t1 = t2 - 1;
389 
390  std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),
391  eVecm[(*t1)].end(),
392  integrDiff);
393  std::vector<double>::iterator e11 = e12 - 1;
394 
395  std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),
396  eVecm[(*t2)].end(),
397  integrDiff);
398  std::vector<double>::iterator e21 = e22 - 1;
399 
400  valueT1 = *t1;
401  valueT2 = *t2;
402  valueE21 = *e21;
403  valueE22 = *e22;
404  valueE12 = *e12;
405  valueE11 = *e11;
406 
407  xs11 = fDiffCrossSectionData[valueT1][valueE11];
408  xs12 = fDiffCrossSectionData[valueT1][valueE12];
409  xs21 = fDiffCrossSectionData[valueT2][valueE21];
410  xs22 = fDiffCrossSectionData[valueT2][valueE22];
411 
412  if(xs11 == 0 && xs12 == 0 && xs21 == 0 && xs22 == 0) return (0.);
413 
414  theta = QuadInterpolator(valueE11, valueE12, valueE21, valueE22, xs11, xs12,
415  xs21, xs22, valueT1, valueT2, k, integrDiff);
416 
417  return theta;
418 }
419 
420 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
421 
422 G4double
424  G4double xs1, G4double xs2)
425 {
426  G4double d1 = xs1;
427  G4double d2 = xs2;
428  G4double value = (d1 + (d2 - d1) * (e - e1) / (e2 - e1));
429  return value;
430 }
431 
432 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
433 
434 G4double
436  G4double xs1, G4double xs2)
437 {
438  G4double d1 = std::log(xs1);
439  G4double d2 = std::log(xs2);
440  G4double value = std::exp(d1 + (d2 - d1) * (e - e1) / (e2 - e1));
441  return value;
442 }
443 
444 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
445 
446 G4double
448  G4double xs1, G4double xs2)
449 {
450  G4double a = (std::log10(xs2) - std::log10(xs1))
451  / (std::log10(e2) - std::log10(e1));
452  G4double b = std::log10(xs2) - a * std::log10(e2);
453  G4double sigma = a * std::log10(e) + b;
454  G4double value = (std::pow(10., sigma));
455  return value;
456 }
457 
458 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
459 
460 G4double
462  G4double e21, G4double e22,
463  G4double xs11, G4double xs12,
464  G4double xs21, G4double xs22,
465  G4double t1, G4double t2, G4double t,
466  G4double e)
467 {
468  // Log-Log
469  /*
470  G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
471  G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
472  G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
473  */
474 
475  // Lin-Log
476  /*
477  G4double interpolatedvalue1 = LinLogInterpolate(e11, e12, e, xs11, xs12);
478  G4double interpolatedvalue2 = LinLogInterpolate(e21, e22, e, xs21, xs22);
479  G4double value = LinLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
480  */
481 
482 // Lin-Lin
483  G4double interpolatedvalue1 = LinLinInterpolate(e11, e12, e, xs11, xs12);
484  G4double interpolatedvalue2 = LinLinInterpolate(e21, e22, e, xs21, xs22);
485  G4double value = LinLinInterpolate(t1, t2, t, interpolatedvalue1,
486  interpolatedvalue2);
487 
488  return value;
489 }
490 
491 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
492 
493 G4double
495  G4double k, G4ParticleDefinition * particleDefinition)
496 {
497  G4double integrdiff = G4UniformRand();
498  return Theta(particleDefinition, k / eV, integrdiff);
499 }
500 
501 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
502 
503 void
505 {
506  killBelowEnergy = threshold;
507 
508  if(killBelowEnergy < 100 * eV)
509  {
510  G4cout << "*** WARNING : the G4DNAIonElasticModel class is not "
511  "activated below 100 eV !"
512  << G4endl;
513  }
514 }
515 
static const double cm
Definition: G4SIunits.hh:118
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
static const double MeV
Definition: G4SIunits.hh:211
static const G4double d1
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
size_t GetIndex() const
Definition: G4Material.hh:262
TriDimensionMap fDiffCrossSectionData
std::vector< double > eTdummyVec
static const G4double e2
virtual G4bool LoadData(const G4String &argFileName)
G4double LinLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
G4ParticleDefinition * GetDefinition() const
G4double a
Definition: TRTMaterials.hh:39
G4double RandomizeThetaCM(G4double k, G4ParticleDefinition *aParticleDefinition)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4String & GetParticleName() const
void SetKillBelowThreshold(G4double threshold)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
G4double LinLinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4DNACrossSectionDataSet * fpTableData
G4ParticleChangeForGamma * fParticleChangeForGamma
const G4ThreeVector & GetMomentumDirection() const
virtual void Initialise(const G4ParticleDefinition *particuleDefinition, const G4DataVector &)
const std::vector< G4double > * fpMolWaterDensity
static G4DNAGenericIonsManager * Instance(void)
virtual G4double FindValue(G4double e, G4int componentId=0) const
static const G4double e1
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const double eV
Definition: G4SIunits.hh:212
static G4DNAMolecularMaterial * Instance()
static const double pi
Definition: G4SIunits.hh:74
static G4ParticleTable * GetParticleTable()
G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
#define G4endl
Definition: G4ios.hh:61
G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double x11, G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, G4double t, G4double e)
static MCTruthManager * instance
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4double Theta(G4ParticleDefinition *aParticleDefinition, G4double k, G4double integrDiff)
static const G4double d2
#define DBL_MAX
Definition: templates.hh:83
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
G4DNAIonElasticModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAIonElasticModel")
G4ParticleDefinition * GetIon(const G4String &name)