source: trunk/source/processes/electromagnetic/lowenergy/src/G4DNAChampionElasticModel.cc @ 991

Last change on this file since 991 was 968, checked in by garnier, 15 years ago

fichier ajoutes

File size: 14.8 KB
Line 
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,v 1.3 2009/02/16 11:00:11 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28//
29
30#include "G4DNAChampionElasticModel.hh"
31
32//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
33
34using namespace std;
35
36//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
37
38G4DNAChampionElasticModel::G4DNAChampionElasticModel(const G4ParticleDefinition*,
39                                             const G4String& nam)
40:G4VEmModel(nam),isInitialised(false)
41{
42
43  killBelowEnergy = 8.23*eV; // Minimum e- energy for energy loss by excitation
44  lowEnergyLimit = 0 * eV; 
45  lowEnergyLimitOfModel = 7 * eV; // The model lower energy is 7 eV
46  highEnergyLimit = 10 * keV;
47  SetLowEnergyLimit(lowEnergyLimit);
48  SetHighEnergyLimit(highEnergyLimit);
49
50  verboseLevel= 0;
51  // Verbosity scale:
52  // 0 = nothing
53  // 1 = warning for energy non-conservation
54  // 2 = details of energy budget
55  // 3 = calculation of cross sections, file openings, sampling of atoms
56  // 4 = entering in methods
57 
58  G4cout << "Champion Elastic model is constructed " << G4endl
59         << "Energy range: "
60         << lowEnergyLimit / eV << " eV - "
61         << highEnergyLimit / keV << " keV"
62         << G4endl;
63 
64}
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67
68G4DNAChampionElasticModel::~G4DNAChampionElasticModel()
69{ 
70  // For total cross section
71 
72  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
73  for (pos = tableData.begin(); pos != tableData.end(); ++pos)
74  {
75    G4DNACrossSectionDataSet* table = pos->second;
76    delete table;
77  }
78
79   // For final state
80   
81   eVecm.clear();
82
83}
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86
87void G4DNAChampionElasticModel::Initialise(const G4ParticleDefinition* /*particle*/,
88                                       const G4DataVector& /*cuts*/)
89{
90
91  if (verboseLevel > 3)
92    G4cout << "Calling G4DNAChampionElasticModel::Initialise()" << G4endl;
93
94  // Energy limits
95 
96  if (LowEnergyLimit() < lowEnergyLimit)
97  {
98    G4cout << "G4DNAChampionElasticModel: low energy limit increased from " << 
99        LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
100    SetLowEnergyLimit(lowEnergyLimit);
101    }
102
103  if (HighEnergyLimit() > highEnergyLimit)
104  {
105    G4cout << "G4DNAChampionElasticModel: high energy limit decreased from " << 
106        HighEnergyLimit()/keV << " keV to " << highEnergyLimit/keV << " keV" << G4endl;
107    SetHighEnergyLimit(highEnergyLimit);
108  }
109
110  // Reading of data files
111 
112  G4double scaleFactor = 1e-16*cm*cm;
113
114  G4String fileElectron("dna/sigma_elastic_e_champion");
115
116  G4ParticleDefinition* electronDef = G4Electron::ElectronDefinition();
117  G4String electron;
118 
119  if (electronDef != 0)
120  {
121    // For total cross section
122   
123    electron = electronDef->GetParticleName();
124
125    tableFile[electron] = fileElectron;
126
127    G4DNACrossSectionDataSet* tableE = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
128    tableE->LoadData(fileElectron);
129    tableData[electron] = tableE;
130   
131    // For final state
132   
133    char *path = getenv("G4LEDATA");
134 
135    if (!path)
136    G4Exception("G4FinalStateElasticChampion::Initialise: G4LEDATA environment variable not set");
137
138    std::ostringstream eFullFileName;
139    eFullFileName << path << "/dna/sigmadiff_elastic_e_champion.dat";
140    std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
141     
142    if (!eDiffCrossSection) G4Exception("G4DNAChampionElasticModel::Initialise: error opening electron DATA FILE");
143     
144    eTdummyVec.push_back(0.);
145
146    while(!eDiffCrossSection.eof())
147    {
148        double tDummy;
149        double eDummy;
150        eDiffCrossSection>>tDummy>>eDummy;
151
152        // SI : mandatory eVecm initialization
153        if (tDummy != eTdummyVec.back()) 
154        { 
155          eTdummyVec.push_back(tDummy); 
156          eVecm[tDummy].push_back(0.);
157        }
158         
159        eDiffCrossSection>>eDiffCrossSectionData[tDummy][eDummy];
160
161        // SI : only if not end of file reached !
162        if (!eDiffCrossSection.eof()) eDiffCrossSectionData[tDummy][eDummy]*=scaleFactor;
163         
164        if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy);
165         
166    }
167
168    // End final state
169 
170  }
171  else G4Exception("G4DNAChampionElasticModel::Initialise: electron is not defined");
172 
173  if (verboseLevel > 2) 
174    G4cout << "Loaded cross section files for Champion Elastic model" << G4endl;
175
176  G4cout << "Champion Elastic model is initialized " << G4endl
177         << "Energy range: "
178         << LowEnergyLimit() / eV << " eV - "
179         << HighEnergyLimit() / keV << " keV"
180         << G4endl;
181
182  if(!isInitialised) 
183  {
184    isInitialised = true;
185 
186    if(pParticleChange)
187      fParticleChangeForGamma = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
188    else
189      fParticleChangeForGamma = new G4ParticleChangeForGamma();
190  }   
191
192  // InitialiseElementSelectors(particle,cuts);
193
194  // Test if water material
195
196  flagMaterialIsWater= false;
197  densityWater = 0;
198
199  const G4ProductionCutsTable* theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
200
201  if(theCoupleTable) 
202  {
203    G4int numOfCouples = theCoupleTable->GetTableSize();
204 
205    if(numOfCouples>0) 
206    {
207          for (G4int i=0; i<numOfCouples; i++) 
208          {
209            const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
210            const G4Material* material = couple->GetMaterial();
211
212            size_t j = material->GetNumberOfElements();
213            while (j>0)
214            {
215               j--;
216               const G4Element* element(material->GetElement(j));
217               if (element->GetZ() == 8.)
218               {
219                  G4double density = material->GetAtomicNumDensityVector()[j];
220                  if (density > 0.) 
221                  { 
222                    flagMaterialIsWater = true; 
223                    densityWater = density; 
224                   
225                    if (verboseLevel > 3) 
226                    G4cout << "Water material is found with density(cm^-3)=" << density/(cm*cm*cm) << G4endl;
227                  }
228               }
229            }
230 
231          }
232   } // if(numOfCouples>0)
233
234  } // if (theCoupleTable)
235
236}
237
238//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
239
240G4double G4DNAChampionElasticModel::CrossSectionPerVolume(const G4Material*,
241                                           const G4ParticleDefinition* p,
242                                           G4double ekin,
243                                           G4double,
244                                           G4double)
245{
246  if (verboseLevel > 3)
247    G4cout << "Calling CrossSectionPerVolume() of G4DNAChampionElasticModel" << G4endl;
248
249 // Calculate total cross section for model
250
251 G4double sigma=0;
252 
253 if (flagMaterialIsWater)
254 {
255  const G4String& particleName = p->GetParticleName();
256
257  if (ekin >= lowEnergyLimitOfModel && ekin < highEnergyLimit)
258  {
259        std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
260        pos = tableData.find(particleName);
261       
262        if (pos != tableData.end())
263        {
264          G4DNACrossSectionDataSet* table = pos->second;
265          if (table != 0)
266          {
267            sigma = table->FindValue(ekin);
268          }
269        }
270        else
271        {
272            G4Exception("G4DNAChampionElasticModel::ComputeCrossSectionPerVolume: attempting to calculate cross section for wrong particle");
273        }
274  }
275
276  if (verboseLevel > 3)
277  {
278    G4cout << "---> Kinetic energy(eV)=" << ekin/eV << G4endl;
279    G4cout << " - Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
280    G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*densityWater/(1./cm) << G4endl;
281  } 
282
283 } // if (flagMaterialIsWater)
284         
285 return sigma*densityWater;               
286}
287
288//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
289
290void G4DNAChampionElasticModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
291                                              const G4MaterialCutsCouple* /*couple*/,
292                                              const G4DynamicParticle* aDynamicElectron,
293                                              G4double,
294                                              G4double)
295{
296
297  if (verboseLevel > 3)
298    G4cout << "Calling SampleSecondaries() of G4DNAChampionElasticModel" << G4endl;
299
300  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
301 
302  if (electronEnergy0 < killBelowEnergy)
303  {
304    fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
305    fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
306    return ;
307  }
308
309  if (electronEnergy0>= killBelowEnergy && electronEnergy0 < highEnergyLimit)
310  { 
311    G4double cosTheta = RandomizeCosTheta(electronEnergy0);
312 
313    G4double phi = 2. * pi * G4UniformRand();
314
315    G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
316    G4ThreeVector xVers = zVers.orthogonal();
317    G4ThreeVector yVers = zVers.cross(xVers);
318
319    G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
320    G4double yDir = xDir;
321    xDir *= std::cos(phi);
322    yDir *= std::sin(phi);
323
324    G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
325
326    fParticleChangeForGamma->ProposeMomentumDirection(zPrimeVers.unit()) ;
327
328    fParticleChangeForGamma->SetProposedKineticEnergy(electronEnergy0);
329  }
330
331}
332
333//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
334
335G4double G4DNAChampionElasticModel::DifferentialCrossSection
336  (G4ParticleDefinition * particleDefinition, G4double k, G4double theta)                                                         
337{
338
339  G4double sigma = 0.;
340  G4double valueT1 = 0;
341  G4double valueT2 = 0;
342  G4double valueE21 = 0;
343  G4double valueE22 = 0;
344  G4double valueE12 = 0;
345  G4double valueE11 = 0;
346  G4double xs11 = 0;   
347  G4double xs12 = 0; 
348  G4double xs21 = 0; 
349  G4double xs22 = 0; 
350
351  //SI : ensure the correct computation of cross section at the 180*deg limit
352  if (theta==180.) theta=theta-1e-9;
353
354  if (particleDefinition == G4Electron::ElectronDefinition()) 
355  {
356    std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k);
357    std::vector<double>::iterator t1 = t2-1;
358 
359    std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),eVecm[(*t1)].end(), theta);
360    std::vector<double>::iterator e11 = e12-1;
361         
362    std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),eVecm[(*t2)].end(), theta);
363    std::vector<double>::iterator e21 = e22-1;
364               
365    valueT1  =*t1;
366    valueT2  =*t2;
367    valueE21 =*e21;
368    valueE22 =*e22;
369    valueE12 =*e12;
370    valueE11 =*e11;
371
372    xs11 = eDiffCrossSectionData[valueT1][valueE11];
373    xs12 = eDiffCrossSectionData[valueT1][valueE12];
374    xs21 = eDiffCrossSectionData[valueT2][valueE21];
375    xs22 = eDiffCrossSectionData[valueT2][valueE22];
376
377}
378     
379  G4double xsProduct = xs11 * xs12 * xs21 * xs22;
380 
381  if (xs11==0 || xs12==0 ||xs21==0 ||xs22==0) return (0.);
382     
383  if (xsProduct != 0.)
384  {
385    sigma = QuadInterpolator(  valueE11, valueE12, 
386                               valueE21, valueE22, 
387                               xs11, xs12, 
388                               xs21, xs22, 
389                               valueT1, valueT2, 
390                               k, theta );
391  }
392
393  return sigma;
394}
395
396//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
397
398G4double G4DNAChampionElasticModel::LinLogInterpolate(G4double e1, 
399                                                        G4double e2, 
400                                                        G4double e, 
401                                                        G4double xs1, 
402                                                        G4double xs2)
403{
404  G4double d1 = std::log(xs1);
405  G4double d2 = std::log(xs2);
406  G4double value = std::exp(d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
407  return value;
408}
409
410//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
411
412G4double G4DNAChampionElasticModel::LogLogInterpolate(G4double e1, 
413                                                        G4double e2, 
414                                                        G4double e, 
415                                                        G4double xs1, 
416                                                        G4double xs2)
417{
418  G4double a = (std::log10(xs2)-std::log10(xs1)) / (std::log10(e2)-std::log10(e1));
419  G4double b = std::log10(xs2) - a*std::log10(e2);
420  G4double sigma = a*std::log10(e) + b;
421  G4double value = (std::pow(10.,sigma));
422  return value;
423}
424
425//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
426
427G4double G4DNAChampionElasticModel::QuadInterpolator(G4double e11, G4double e12, 
428                                                       G4double e21, G4double e22, 
429                                                       G4double xs11, G4double xs12, 
430                                                       G4double xs21, G4double xs22, 
431                                                       G4double t1, G4double t2, 
432                                                       G4double t, G4double e)
433{
434// Log-Log
435/*
436  G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
437  G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
438  G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
439*/
440
441// Lin-Log
442  G4double interpolatedvalue1 = LinLogInterpolate(e11, e12, e, xs11, xs12);
443  G4double interpolatedvalue2 = LinLogInterpolate(e21, e22, e, xs21, xs22);
444  G4double value = LinLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
445  return value;
446}
447
448//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
449
450G4double G4DNAChampionElasticModel::RandomizeCosTheta(G4double k) 
451{
452 // ***** Similar method as for screened Rutherford scattering
453
454 G4int iMax=180;
455 G4double max=0;
456 G4double tmp=0;
457 
458 // Look for maximum :
459 for (G4int i=0; i<iMax; i++) 
460 {
461   tmp = DifferentialCrossSection(G4Electron::ElectronDefinition(),k/eV,G4double(i)*180./(iMax-1));
462   if (tmp>max) max = tmp;
463 }
464
465 G4double oneOverMax=0;
466 if (max!=0) oneOverMax = 1./max;
467 
468 G4double cosTheta = 0.;
469 G4double fCosTheta = 0.;
470 
471 do
472 {
473    cosTheta = 2. * G4UniformRand() - 1.;
474    fCosTheta = oneOverMax * DifferentialCrossSection(G4Electron::ElectronDefinition(),k/eV,std::acos(cosTheta)*180./pi);
475 }
476 while (fCosTheta < G4UniformRand());
477
478 if (verboseLevel > 3)
479 {
480   G4cout << "---> Cos(theta)=" << cosTheta << G4endl;
481 } 
482
483 return cosTheta; 
484}
Note: See TracBrowser for help on using the repository browser.