Geant4  10.02.p02
G4DNASancheExcitationModel.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: G4DNASancheExcitationModel.cc 93616 2015-10-27 08:59:17Z gcosmo $
27 //
28 
29 // Created by Z. Francis
30 
32 #include "G4SystemOfUnits.hh"
34 
35 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
36 
37 using namespace std;
38 
39 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
40 
42  const G4String& nam) :
43  G4VEmModel(nam), isInitialised(false)
44 {
45  fpWaterDensity = 0;
46 
47  lowEnergyLimit = 2 * eV;
48  highEnergyLimit = 100 * eV;
51  nLevels = 9;
52 
53  verboseLevel = 0;
54  // Verbosity scale:
55  // 0 = nothing
56  // 1 = warning for energy non-conservation
57  // 2 = details of energy budget
58  // 3 = calculation of cross sections, file openings, sampling of atoms
59  // 4 = entering in methods
60 
61  if (verboseLevel > 0)
62  {
63  G4cout << "Sanche Excitation model is constructed " << G4endl<< "Energy range: "
64  << lowEnergyLimit / eV << " eV - "
65  << highEnergyLimit / eV << " eV"
66  << G4endl;
67  }
69  fpWaterDensity = 0;
70 }
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73 
75 {
76 }
77 
78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
79 
81  const G4DataVector& /*cuts*/)
82 {
83 
84  if (verboseLevel > 3)
85  G4cout << "Calling G4DNASancheExcitationModel::Initialise()" << G4endl;
86 
87  // Energy limits
88 
90  {
91  G4cout << "G4DNASancheExcitationModel: low energy limit increased from " <<
92  LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
94  }
95 
97  {
98  G4cout << "G4DNASancheExcitationModel: high energy limit decreased from " <<
99  HighEnergyLimit()/eV << " eV to " << highEnergyLimit/eV << " eV" << G4endl;
101  }
102 
103  //
104 
105  if (verboseLevel > 0)
106  {
107  G4cout << "Sanche Excitation model is initialized " << G4endl
108  << "Energy range: "
109  << LowEnergyLimit() / eV << " eV - "
110  << HighEnergyLimit() / eV << " eV"
111  << G4endl;
112  }
113 
114  // Initialize water density pointer
116  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
117 
118  if (isInitialised) {return;} // RETURNS HERE
119 
121  isInitialised = true;
122 
123  char *path = getenv("G4LEDATA");
124  std::ostringstream eFullFileName;
125  eFullFileName << path << "/dna/sigma_excitationvib_e_sanche.dat";
126  std::ifstream input(eFullFileName.str().c_str());
127 
128  if (!input)
129  {
130  G4Exception("G4DNASancheExcitationModel::Initialise","em0003",
131  FatalException,"Missing data file:/dna/sigma_excitationvib_e_sanche.dat");
132  }
133 
134  // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
135  // Added clear for MT
136  tdummyVec.clear();
137  //
138 
139  double t;
140  double xs;
141 
142  while(!input.eof())
143  {
144  input>>t;
145  tdummyVec.push_back(t);
146 
147  fEnergyLevelXS.push_back(std::vector<double>());
148  fEnergyTotalXS.push_back(0);
149  std::vector<double>& levelXS = fEnergyLevelXS.back();
150  levelXS.reserve(9);
151 
152 // G4cout<<t;
153 
154  for(size_t i = 0 ; i < 9 ;++i)
155  {
156  input>>xs;
157  levelXS.push_back(xs);
158  fEnergyTotalXS.back() += xs;
159 // G4cout <<" " << levelXS[i];
160  }
161 
162 // G4cout << G4endl;
163  }
164 }
165 
166 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
167 
169  const G4ParticleDefinition* particleDefinition,
170  G4double ekin,
171  G4double,
172  G4double)
173 {
174  if (verboseLevel > 3)
175  G4cout << "Calling CrossSectionPerVolume() of G4DNASancheExcitationModel"
176  << G4endl;
177 
178  // Calculate total cross section for model
179 
180  G4double sigma=0;
181 
182  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
183 
184  if(waterDensity!= 0.0)
185  {
186 
187 // if (particleDefinition == G4Electron::ElectronDefinition()) // pas besoin
188  {
189  if (ekin >= lowEnergyLimit && ekin < highEnergyLimit)
190  // for now this is necessary
191  {
192  //sigma = Sum(ekin);
193  sigma = TotalCrossSection(ekin);
194  }
195  }
196 
197  if (verboseLevel > 2)
198  {
199  G4cout << "__________________________________" << G4endl;
200  G4cout << "=== G4DNASancheExcitationModel - XS INFO START" << G4endl;
201  G4cout << "=== Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
202  G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
203  G4cout << "=== Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
204  G4cout << "=== G4DNASancheExcitationModel - XS INFO END" << G4endl;
205  }
206 
207  } // if water
208 
209  return sigma*2*waterDensity;
210  // see papers for factor 2 description
211 
212 }
213 
214 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
215 
218  const G4MaterialCutsCouple*,
219  const G4DynamicParticle* aDynamicElectron,
220  G4double,
221  G4double)
222 {
223 
224  if (verboseLevel > 3)
225  G4cout << "Calling SampleSecondaries() of G4DNASancheExcitationModel"
226  << G4endl;
227 
228  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
229  G4int level = RandomSelect(electronEnergy0);
230  G4double excitationEnergy = VibrationEnergy(level); // levels go from 0 to 8
231  G4double newEnergy = electronEnergy0 - excitationEnergy;
232 
233  /*
234  if (electronEnergy0 < highEnergyLimit)
235  {
236  if (newEnergy >= lowEnergyLimit)
237  {
238  fParticleChangeForGamma->ProposeMomentumDirection(aDynamicElectron->GetMomentumDirection());
239  fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
240  fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy);
241  }
242 
243  else
244  {
245  fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
246  fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
247  }
248  }
249  */
250 
251  if (electronEnergy0 < highEnergyLimit && newEnergy>0.)
252  {
256  }
257 
258  //
259 }
260 
261 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
262 
264  G4int level)
265 {
266  std::vector<double>::iterator t2 = std::upper_bound(tdummyVec.begin(),
267  tdummyVec.end(), t / eV);
268  std::vector<double>::iterator t1 = t2 - 1;
269 
270  size_t i1 = t1 - tdummyVec.begin();
271  size_t i2 = t2 - tdummyVec.begin();
272 
273  double sigma = LinInterpolate((*t1), (*t2),
274  t / eV,
275  fEnergyLevelXS[i1][level],
276  fEnergyLevelXS[i2][level]);
277 
278  static const double conv_factor = 1e-16 * cm * cm;
279 
280  sigma *= conv_factor;
281  if (sigma == 0.) sigma = 1e-30;
282  return (sigma);
283 }
284 
285 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
286 
288 {
289  std::vector<double>::iterator t2 = std::upper_bound(tdummyVec.begin(),
290  tdummyVec.end(), t / eV);
291  std::vector<double>::iterator t1 = t2 - 1;
292 
293  size_t i1 = t1 - tdummyVec.begin();
294  size_t i2 = t2 - tdummyVec.begin();
295 
296  double sigma = LinInterpolate((*t1), (*t2),
297  t / eV,
298  fEnergyTotalXS[i1],
299  fEnergyTotalXS[i2]);
300 
301  static const double conv_factor = 1e-16 * cm * cm;
302 
303  sigma *= conv_factor;
304  if (sigma == 0.) sigma = 1e-30;
305  return (sigma);
306 }
307 
308 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
309 
311 {
312  static G4double energies[9] = { 0.01, 0.024, 0.061, 0.092, 0.204, 0.417, 0.460,
313  0.500, 0.835 };
314  return (energies[level] * eV);
315 }
316 
317 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
318 
320 {
321 
322  // Level Selection Counting can be done here !
323 
324  G4int i = nLevels;
325  G4double value = 0.;
326  std::deque<double> values;
327 
328  while (i > 0)
329  {
330  i--;
331  G4double partial = PartialCrossSection(k, i);
332  values.push_front(partial);
333  value += partial;
334  }
335 
336  value *= G4UniformRand();
337 
338  i = nLevels;
339 
340  while (i > 0)
341  {
342  i--;
343  if (values[i] > value)
344  {
345  //outcount<<i<<" "<<VibrationEnergy(i)<<G4endl;
346  return i;
347  }
348  value -= values[i];
349  }
350 
351  //outcount<<0<<" "<<VibrationEnergy(0)<<G4endl;
352 
353  return 0;
354 }
355 
356 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
357 
359 {
360  G4double totalCrossSection = 0.;
361 
362  for (G4int i = 0; i < nLevels; i++)
363  {
364  totalCrossSection += PartialCrossSection(k, i);
365  }
366 
367  return totalCrossSection;
368 }
369 
370 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
371 
373  G4double e2,
374  G4double e,
375  G4double xs1,
376  G4double xs2)
377 {
378  G4double a = (xs2 - xs1) / (e2 - e1);
379  G4double b = xs2 - a * e2;
380  G4double value = a * e + b;
381  // G4cout<<"interP >> "<<e1<<" "<<e2<<" "<<e<<" "<<xs1<<" "<<xs2<<" "<<a<<" "<<b<<" "<<value<<G4endl;
382 
383  return value;
384 }
385 
static const double cm
Definition: G4SIunits.hh:118
G4DNASancheExcitationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNASancheExcitationModel")
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4double GetKineticEnergy() const
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
G4double PartialCrossSection(G4double energy, G4int level)
static const G4double e2
G4double a
Definition: TRTMaterials.hh:39
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4String & GetParticleName() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
G4double LinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const std::vector< G4double > * fpWaterDensity
const G4ThreeVector & GetMomentumDirection() const
G4ParticleChangeForGamma * fParticleChangeForGamma
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
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
static G4DNAMolecularMaterial * Instance()
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
std::vector< std::vector< double > > fEnergyLevelXS