Geant4  10.02.p02
G4DNAChampionElasticModel.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 // $Id: G4DNAChampionElasticModel.cc 87375 2014-12-02 08:17:28Z gcosmo $
27 //
28 
30 #include "G4PhysicalConstants.hh"
31 #include "G4SystemOfUnits.hh"
33 
34 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
35 
36 using namespace std;
37 
38 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
39 
41  const G4String& nam) :
42  G4VEmModel(nam), isInitialised(false)
43 {
44 // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
45 
46  killBelowEnergy = 7.4 * eV;
47  lowEnergyLimit = 0 * eV;
48  highEnergyLimit = 1. * MeV;
51 
52  verboseLevel = 0;
53  // Verbosity scale:
54  // 0 = nothing
55  // 1 = warning for energy non-conservation
56  // 2 = details of energy budget
57  // 3 = calculation of cross sections, file openings, sampling of atoms
58  // 4 = entering in methods
59 
60  if (verboseLevel > 0)
61  {
62  G4cout << "Champion Elastic model is constructed " << G4endl<< "Energy range: "
63  << lowEnergyLimit / eV << " eV - "
64  << highEnergyLimit / MeV << " MeV"
65  << G4endl;
66  }
69 }
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
72 
74 {
75  // For total cross section
76 
77  std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
78  for (pos = tableData.begin(); pos != tableData.end(); ++pos)
79  {
80  G4DNACrossSectionDataSet* table = pos->second;
81  delete table;
82  }
83 
84  // For final state
85 
86  eVecm.clear();
87 
88 }
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91 
93  const G4DataVector& /*cuts*/)
94 {
95 
96  if (verboseLevel > 3)
97  G4cout << "Calling G4DNAChampionElasticModel::Initialise()" << G4endl;
98 
99  // Energy limits
100 
102  {
103  G4cout << "G4DNAChampionElasticModel: low energy limit increased from " <<
104  LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
106  }
107 
109  {
110  G4cout << "G4DNAChampionElasticModel: high energy limit decreased from " <<
111  HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl;
113  }
114 
115  // Reading of data files
116 
117  G4double scaleFactor = 1e-16*cm*cm;
118 
119  G4String fileElectron("dna/sigma_elastic_e_champion");
120 
123 
124  // *** ELECTRON
125 
126  // For total cross section
127 
128  electron = electronDef->GetParticleName();
129 
130  tableFile[electron] = fileElectron;
131 
132  G4DNACrossSectionDataSet* tableE =
133  new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
134  tableE->LoadData(fileElectron);
135  tableData[electron] = tableE;
136 
137  // For final state
138 
139  char *path = getenv("G4LEDATA");
140 
141  if (!path)
142  {
143  G4Exception("G4ChampionElasticModel::Initialise","em0006",
144  FatalException,"G4LEDATA environment variable not set.");
145  return;
146  }
147 
148  std::ostringstream eFullFileName;
149  eFullFileName << path << "/dna/sigmadiff_cumulated_elastic_e_champion_hp.dat";
150  std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
151 
152  if (!eDiffCrossSection)
153  G4Exception("G4DNAChampionElasticModel::Initialise","em0003",
155  "Missing data file:/dna/sigmadiff_cumulated_elastic_e_champion_hp.dat; please use G4EMLOW6.36 and above.");
156 
157  // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
158  // Added clear for MT
159 
160  eTdummyVec.clear();
161  eVecm.clear();
162  eDiffCrossSectionData.clear();
163 
164  //
165 
166  eTdummyVec.push_back(0.);
167 
168  while(!eDiffCrossSection.eof())
169  {
170  double tDummy;
171  double eDummy;
172  eDiffCrossSection>>tDummy>>eDummy;
173 
174  // SI : mandatory eVecm initialization
175 
176  if (tDummy != eTdummyVec.back())
177  {
178  eTdummyVec.push_back(tDummy);
179  eVecm[tDummy].push_back(0.);
180  }
181 
182  eDiffCrossSection>>eDiffCrossSectionData[tDummy][eDummy];
183 
184  if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy);
185 
186  }
187 
188  // End final state
189 
190  if (verboseLevel > 2)
191  G4cout << "Loaded cross section files for Champion Elastic model" << G4endl;
192 
193  if( verboseLevel>0 )
194  {
195  G4cout << "Champion Elastic model is initialized " << G4endl
196  << "Energy range: "
197  << LowEnergyLimit() / eV << " eV - "
198  << HighEnergyLimit() / MeV << " MeV"
199  << G4endl;
200  }
201 
202  // Initialize water density pointer
205 
206  if (isInitialised)
207  { return;}
209  isInitialised = true;
210 
211 }
212 
213 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
214 
216  const G4ParticleDefinition* p,
217  G4double ekin,
218  G4double,
219  G4double)
220 {
221  if (verboseLevel > 3)
222  G4cout << "Calling CrossSectionPerVolume() of G4DNAChampionElasticModel"
223  << G4endl;
224 
225  // Calculate total cross section for model
226 
227  G4double sigma=0;
228 
229  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
230 
231  if(waterDensity!= 0.0)
232 // if (material == nistwater || material->GetBaseMaterial() == nistwater)
233  {
234  const G4String& particleName = p->GetParticleName();
235 
236  if (ekin < highEnergyLimit)
237  {
238  //SI : XS must not be zero otherwise sampling of secondaries method ignored
239  if (ekin < killBelowEnergy) return DBL_MAX;
240  //
241 
242  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
243  pos = tableData.find(particleName);
244 
245  if (pos != tableData.end())
246  {
247  G4DNACrossSectionDataSet* table = pos->second;
248  if (table != 0)
249  {
250  sigma = table->FindValue(ekin);
251  }
252  }
253  else
254  {
255  G4Exception("G4DNAChampionElasticModel::ComputeCrossSectionPerVolume","em0002",
256  FatalException,"Model not applicable to particle type.");
257  }
258  }
259 
260  if (verboseLevel > 2)
261  {
262  G4cout << "__________________________________" << G4endl;
263  G4cout << "=== G4DNAChampionElasticModel - XS INFO START" << G4endl;
264  G4cout << "=== Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
265  G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
266  G4cout << "=== Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
267  // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
268  G4cout << "=== G4DNAChampionElasticModel - XS INFO END" << G4endl;
269  }
270 
271  }
272 
273  return sigma*waterDensity;
274 // return sigma*material->GetAtomicNumDensityVector()[1];
275 }
276 
277 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
278 
279 void G4DNAChampionElasticModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
280  const G4MaterialCutsCouple* /*couple*/,
281  const G4DynamicParticle* aDynamicElectron,
282  G4double,
283  G4double)
284 {
285 
286  if (verboseLevel > 3)
287  G4cout << "Calling SampleSecondaries() of G4DNAChampionElasticModel" << G4endl;
288 
289  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
290 
291  if (electronEnergy0 < killBelowEnergy)
292  {
296  return;
297  }
298 
299  if (electronEnergy0>= killBelowEnergy && electronEnergy0 < highEnergyLimit)
300  {
301 
302  G4double cosTheta = RandomizeCosTheta(electronEnergy0);
303 
304  G4double phi = 2. * pi * G4UniformRand();
305 
306  G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
307  G4ThreeVector xVers = zVers.orthogonal();
308  G4ThreeVector yVers = zVers.cross(xVers);
309 
310  G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
311  G4double yDir = xDir;
312  xDir *= std::cos(phi);
313  yDir *= std::sin(phi);
314 
315  G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
316 
318 
320  }
321 
322 }
323 
324 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
325 
327  G4double k,
328  G4double integrDiff)
329 {
330  G4double theta = 0.;
331  G4double valueT1 = 0;
332  G4double valueT2 = 0;
333  G4double valueE21 = 0;
334  G4double valueE22 = 0;
335  G4double valueE12 = 0;
336  G4double valueE11 = 0;
337  G4double xs11 = 0;
338  G4double xs12 = 0;
339  G4double xs21 = 0;
340  G4double xs22 = 0;
341 
342  if (particleDefinition == G4Electron::ElectronDefinition())
343  {
344  std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),
345  eTdummyVec.end(), k);
346  std::vector<double>::iterator t1 = t2 - 1;
347 
348  std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),
349  eVecm[(*t1)].end(),
350  integrDiff);
351  std::vector<double>::iterator e11 = e12 - 1;
352 
353  std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),
354  eVecm[(*t2)].end(),
355  integrDiff);
356  std::vector<double>::iterator e21 = e22 - 1;
357 
358  valueT1 = *t1;
359  valueT2 = *t2;
360  valueE21 = *e21;
361  valueE22 = *e22;
362  valueE12 = *e12;
363  valueE11 = *e11;
364 
365  xs11 = eDiffCrossSectionData[valueT1][valueE11];
366  xs12 = eDiffCrossSectionData[valueT1][valueE12];
367  xs21 = eDiffCrossSectionData[valueT2][valueE21];
368  xs22 = eDiffCrossSectionData[valueT2][valueE22];
369  }
370 
371  if (xs11 == 0 && xs12 == 0 && xs21 == 0 && xs22 == 0) return (0.);
372 
373  theta = QuadInterpolator(valueE11, valueE12, valueE21, valueE22, xs11, xs12,
374  xs21, xs22, valueT1, valueT2, k, integrDiff);
375 
376  return theta;
377 }
378 
379 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
380 
382  G4double e2,
383  G4double e,
384  G4double xs1,
385  G4double xs2)
386 {
387  G4double d1 = std::log(xs1);
388  G4double d2 = std::log(xs2);
389  G4double value = std::exp(d1 + (d2 - d1) * (e - e1) / (e2 - e1));
390  return value;
391 }
392 
393 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
394 
396  G4double e2,
397  G4double e,
398  G4double xs1,
399  G4double xs2)
400 {
401  G4double d1 = xs1;
402  G4double d2 = xs2;
403  G4double value = (d1 + (d2 - d1) * (e - e1) / (e2 - e1));
404  return value;
405 }
406 
407 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
408 
410  G4double e2,
411  G4double e,
412  G4double xs1,
413  G4double xs2)
414 {
415  G4double a = (std::log10(xs2) - std::log10(xs1))
416  / (std::log10(e2) - std::log10(e1));
417  G4double b = std::log10(xs2) - a * std::log10(e2);
418  G4double sigma = a * std::log10(e) + b;
419  G4double value = (std::pow(10., sigma));
420  return value;
421 }
422 
423 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
424 
426  G4double e12,
427  G4double e21,
428  G4double e22,
429  G4double xs11,
430  G4double xs12,
431  G4double xs21,
432  G4double xs22,
433  G4double t1,
434  G4double t2,
435  G4double t,
436  G4double e)
437 {
438  // Log-Log
439  /*
440  G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
441  G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
442  G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
443 
444 
445  // Lin-Log
446  G4double interpolatedvalue1 = LinLogInterpolate(e11, e12, e, xs11, xs12);
447  G4double interpolatedvalue2 = LinLogInterpolate(e21, e22, e, xs21, xs22);
448  G4double value = LinLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
449  */
450 
451  // Lin-Lin
452  G4double interpolatedvalue1 = LinLinInterpolate(e11, e12, e, xs11, xs12);
453  G4double interpolatedvalue2 = LinLinInterpolate(e21, e22, e, xs21, xs22);
454  G4double value = LinLinInterpolate(t1, t2, t, interpolatedvalue1,
455  interpolatedvalue2);
456 
457  return value;
458 }
459 
460 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
461 
463 {
464 
465  G4double integrdiff = 0;
466  G4double uniformRand = G4UniformRand();
467  integrdiff = uniformRand;
468 
469  G4double theta = 0.;
470  G4double cosTheta = 0.;
471  theta = Theta(G4Electron::ElectronDefinition(), k / eV, integrdiff);
472 
473  cosTheta = std::cos(theta * pi / 180);
474 
475  return cosTheta;
476 }
477 
478 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
479 
481 {
482  killBelowEnergy = threshold;
483 
484  if (threshold < 7.4 * eV)
485  {
486  G4cout<< "*** WARNING : the G4DNAChampionElasticModel class is not "
487  "activated below 7.4 eV !" << G4endl;
488  G4cout<< "*** NOTE: if you are using G4EmDNAChemistry, do not worry, this is"
489  " the expected behavior" << G4endl;
490  }
491 }
static const double cm
Definition: G4SIunits.hh:118
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
static const double MeV
Definition: G4SIunits.hh:211
static const G4double d1
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
G4DNAChampionElasticModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAChampionElasticModel")
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
size_t GetIndex() const
Definition: G4Material.hh:262
void SetKillBelowThreshold(G4double threshold)
static const G4double e2
virtual G4bool LoadData(const G4String &argFileName)
G4double a
Definition: TRTMaterials.hh:39
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4double LinLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4String & GetParticleName() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
const G4ThreeVector & GetMomentumDirection() const
virtual G4double FindValue(G4double e, G4int componentId=0) const
G4double Theta(G4ParticleDefinition *aParticleDefinition, G4double k, G4double integrDiff)
static const G4double e1
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
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 const double eV
Definition: G4SIunits.hh:212
static G4DNAMolecularMaterial * Instance()
static const double pi
Definition: G4SIunits.hh:74
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
static const G4double d2
G4double LinLinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
#define DBL_MAX
Definition: templates.hh:83
G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
static const G4double pos
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
const std::vector< G4double > * fpMolWaterDensity