Geant4  10.01.p03
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 87137 2014-11-25 09:12:48Z 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  // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
46  fpWaterDensity = 0;
47 
48  lowEnergyLimit = 2 * eV;
49  highEnergyLimit = 100 * eV;
52  nLevels = 9;
53 
54  verboseLevel = 0;
55  // Verbosity scale:
56  // 0 = nothing
57  // 1 = warning for energy non-conservation
58  // 2 = details of energy budget
59  // 3 = calculation of cross sections, file openings, sampling of atoms
60  // 4 = entering in methods
61 
62  if (verboseLevel > 0)
63  {
64  G4cout << "Sanche Excitation model is constructed " << G4endl<< "Energy range: "
65  << lowEnergyLimit / eV << " eV - "
66  << highEnergyLimit / eV << " eV"
67  << G4endl;
68  }
70  fpWaterDensity = 0;
71 }
72 
73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74 
76 {
77 }
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80 
82  const G4DataVector& /*cuts*/)
83 {
84 
85  if (verboseLevel > 3)
86  G4cout << "Calling G4DNASancheExcitationModel::Initialise()" << G4endl;
87 
88  // Energy limits
89 
91  {
92  G4cout << "G4DNASancheExcitationModel: low energy limit increased from " <<
93  LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
95  }
96 
98  {
99  G4cout << "G4DNASancheExcitationModel: high energy limit decreased from " <<
100  HighEnergyLimit()/eV << " eV to " << highEnergyLimit/eV << " eV" << G4endl;
102  }
103 
104  //
105 
106  if (verboseLevel > 0)
107  {
108  G4cout << "Sanche Excitation model is initialized " << G4endl
109  << "Energy range: "
110  << LowEnergyLimit() / eV << " eV - "
111  << HighEnergyLimit() / eV << " eV"
112  << G4endl;
113  }
114 
115  // Initialize water density pointer
117  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
118 
119  if (isInitialised) {return;} // RETURNS HERE
120 
122  isInitialised = true;
123 
124  char *path = getenv("G4LEDATA");
125  std::ostringstream eFullFileName;
126  eFullFileName << path << "/dna/sigma_excitationvib_e_sanche.dat";
127  std::ifstream input(eFullFileName.str().c_str());
128 
129  if (!input)
130  {
131  G4Exception("G4DNASancheExcitationModel::Initialise","em0003",
132  FatalException,"Missing data file:/dna/sigma_excitationvib_e_sanche.dat");
133  }
134 
135  // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
136  // Added clear for MT
137  tdummyVec.clear();
138  //
139 
140  while(!input.eof())
141  {
142  double t;
143  input>>t;
144  tdummyVec.push_back(t);
145  input>>map1[t][0]>>map1[t][1]>>map1[t][2]>>map1[t][3]>>map1[t][4]>>map1[t][5]>>map1[t][6]>>map1[t][7]>>map1[t][8];
146  //G4cout<<t<<" "<<map1[t][0]<<map1[t][1]<<map1[t][2]<<map1[t][3]<<map1[t][4]<<map1[t][5]<<map1[t][6]<<map1[t][7]<<map1[t][8]<<G4endl;
147  }
148 }
149 
150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
151 
153  const G4ParticleDefinition* particleDefinition,
154  G4double ekin,
155  G4double,
156  G4double)
157 {
158  if (verboseLevel > 3)
159  G4cout << "Calling CrossSectionPerVolume() of G4DNASancheExcitationModel"
160  << G4endl;
161 
162  // Calculate total cross section for model
163 
164  G4double sigma=0;
165 
166  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
167 
168  if(waterDensity!= 0.0)
169  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
170  {
171 
172  if (particleDefinition == G4Electron::ElectronDefinition())
173  {
174  if (ekin >= lowEnergyLimit && ekin < highEnergyLimit)
175  {
176  sigma = Sum(ekin);
177  }
178  }
179 
180  if (verboseLevel > 2)
181  {
182  G4cout << "__________________________________" << G4endl;
183  G4cout << "°°° G4DNASancheExcitationModel - XS INFO START" << G4endl;
184  G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
185  G4cout << "°°° Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
186  G4cout << "°°° Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
187  // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
188  G4cout << "°°° G4DNASancheExcitationModel - XS INFO END" << G4endl;
189  }
190 
191  } // if water
192 
193  // return sigma*2*material->GetAtomicNumDensityVector()[1];
194  return sigma*2*waterDensity;
195  // see papers for factor 2 description
196 
197 }
198 
199 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
200 
203  const G4MaterialCutsCouple*,
204  const G4DynamicParticle* aDynamicElectron,
205  G4double,
206  G4double)
207 {
208 
209  if (verboseLevel > 3)
210  G4cout << "Calling SampleSecondaries() of G4DNASancheExcitationModel"
211  << G4endl;
212 
213  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
214  G4int level = RandomSelect(electronEnergy0);
215  G4double excitationEnergy = VibrationEnergy(level); // levels go from 0 to 8
216  G4double newEnergy = electronEnergy0 - excitationEnergy;
217 
218  /*
219  if (electronEnergy0 < highEnergyLimit)
220  {
221  if (newEnergy >= lowEnergyLimit)
222  {
223  fParticleChangeForGamma->ProposeMomentumDirection(aDynamicElectron->GetMomentumDirection());
224  fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
225  fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy);
226  }
227 
228  else
229  {
230  fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
231  fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
232  }
233  }
234  */
235 
236  if (electronEnergy0 < highEnergyLimit && newEnergy>0.)
237  {
241  }
242 
243  //
244 }
245 
246 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
247 
249  G4int level)
250 {
251  std::vector<double>::iterator t2 = std::upper_bound(tdummyVec.begin(),
252  tdummyVec.end(), t / eV);
253  std::vector<double>::iterator t1 = t2 - 1;
254 
255  double sigma = LinInterpolate((*t1), (*t2), t / eV, map1[*t1][level],
256  map1[*t2][level]);
257  sigma *= 1e-16 * cm * cm;
258  if (sigma == 0.) sigma = 1e-30;
259  return (sigma);
260 }
261 
262 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
263 
265 {
266  G4double energies[9] = { 0.01, 0.024, 0.061, 0.092, 0.204, 0.417, 0.460,
267  0.500, 0.835 };
268  return (energies[level] * eV);
269 }
270 
271 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
272 
274 {
275 
276  // Level Selection Counting can be done here !
277 
278  G4int i = nLevels;
279  G4double value = 0.;
280  std::deque<double> values;
281 
282  while (i > 0)
283  {
284  i--;
285  G4double partial = PartialCrossSection(k, i);
286  values.push_front(partial);
287  value += partial;
288  }
289 
290  value *= G4UniformRand();
291 
292  i = nLevels;
293 
294  while (i > 0)
295  {
296  i--;
297  if (values[i] > value)
298  {
299  //outcount<<i<<" "<<VibrationEnergy(i)<<G4endl;
300  return i;
301  }
302  value -= values[i];
303  }
304 
305  //outcount<<0<<" "<<VibrationEnergy(0)<<G4endl;
306 
307  return 0;
308 }
309 
310 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
311 
313 {
314  G4double totalCrossSection = 0.;
315 
316  for (G4int i = 0; i < nLevels; i++)
317  {
318  totalCrossSection += PartialCrossSection(k, i);
319  }
320  return totalCrossSection;
321 }
322 
323 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
324 
326  G4double e2,
327  G4double e,
328  G4double xs1,
329  G4double xs2)
330 {
331  G4double a = (xs2 - xs1) / (e2 - e1);
332  G4double b = xs2 - a * e2;
333  G4double value = a * e + b;
334  // G4cout<<"interP >> "<<e1<<" "<<e2<<" "<<e<<" "<<xs1<<" "<<xs2<<" "<<a<<" "<<b<<" "<<value<<G4endl;
335 
336  return value;
337 }
338 
static const double cm
Definition: G4SIunits.hh:106
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
G4DNASancheExcitationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNASancheExcitationModel")
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:623
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4double GetKineticEnergy() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:616
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:603
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:707
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:93
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:194
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:714
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134