source: trunk/source/processes/electromagnetic/lowenergy/src/G4FinalStateElasticChampion.cc @ 989

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

fichier ajoutes

File size: 8.9 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: G4FinalStateElasticChampion.cc,v 1.7 2008/12/10 18:25:28 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28// -------------------------------------------------------------------
29
30
31#include "G4FinalStateElasticChampion.hh"
32
33//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
34
35G4FinalStateElasticChampion::G4FinalStateElasticChampion()
36{
37  G4ParticleDefinition* electronDef = G4Electron::ElectronDefinition();
38 
39  lowEnergyLimit = 8.23* eV; // SI : i/o of 7. * eV;
40  highEnergyLimit = 10. * keV;
41 
42  G4double scaleFactor = 1e-16*cm*cm;
43
44  char *path = getenv("G4LEDATA");
45 
46  if (!path)
47    G4Exception("G4FinalStateElasticChampion: constructor: G4LEDATA environment variable not set");
48
49  if (electronDef != 0)
50  {
51    std::ostringstream eFullFileName;
52    eFullFileName << path << "/dna/sigmadiff_elastic_e_champion.dat";
53    std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
54     
55    if (!eDiffCrossSection) G4Exception("G4FinalStateElasticChampion: constructor: error opening electron DATA FILE");
56     
57    eTdummyVec.push_back(0.);
58
59    while(!eDiffCrossSection.eof())
60    {
61          double tDummy;
62          double eDummy;
63          eDiffCrossSection>>tDummy>>eDummy;
64
65          // SI : mandatory eVecm initialization
66          if (tDummy != eTdummyVec.back()) 
67          { 
68            eTdummyVec.push_back(tDummy); 
69            eVecm[tDummy].push_back(0.);
70          }
71         
72          eDiffCrossSection>>eDiffCrossSectionData[0][tDummy][eDummy];
73
74          // SI : only if not end of file reached !
75          if (!eDiffCrossSection.eof()) eDiffCrossSectionData[0][tDummy][eDummy]*=scaleFactor;
76         
77          if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy);
78         
79    }
80
81  }
82  else
83  {
84      G4Exception("G4FinalStateElastiChampion : constructor: electron is not defined");
85  }
86}
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89
90G4FinalStateElasticChampion::~G4FinalStateElasticChampion()
91{ 
92  eVecm.clear();
93}
94
95//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
96
97const G4FinalStateProduct& G4FinalStateElasticChampion::GenerateFinalState(const G4Track& track, const G4Step& )
98{
99  product.Clear();
100
101  G4double k = track.GetDynamicParticle()->GetKineticEnergy();
102
103  if (k>= lowEnergyLimit && k < highEnergyLimit)
104  { 
105    G4double cosTheta = RandomizeCosTheta(k);
106 
107    G4double phi = 2. * pi * G4UniformRand();
108
109    G4ThreeVector zVers = track.GetDynamicParticle()->GetMomentumDirection();
110    G4ThreeVector xVers = zVers.orthogonal();
111    G4ThreeVector yVers = zVers.cross(xVers);
112
113    G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
114    G4double yDir = xDir;
115    xDir *= std::cos(phi);
116    yDir *= std::sin(phi);
117
118    G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
119
120    product.ModifyPrimaryParticle(zPrimeVers,k);
121  }
122
123  if (k<lowEnergyLimit)
124  {
125    product.KillPrimaryParticle();
126  }
127 
128  return product;
129}
130
131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
132
133G4double G4FinalStateElasticChampion::DifferentialCrossSection
134  (G4ParticleDefinition * particleDefinition, G4double k, G4double theta)                                                         
135{
136  G4double sigma = 0.;
137  G4double valueT1 = 0;
138  G4double valueT2 = 0;
139  G4double valueE21 = 0;
140  G4double valueE22 = 0;
141  G4double valueE12 = 0;
142  G4double valueE11 = 0;
143  G4double xs11 = 0;   
144  G4double xs12 = 0; 
145  G4double xs21 = 0; 
146  G4double xs22 = 0; 
147
148  //SI : ensure the correct computation of cross section at the 180*deg limit
149  if (theta==180.) theta=theta-1e-9;
150
151  if (particleDefinition == G4Electron::ElectronDefinition()) 
152  {
153    std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k);
154    std::vector<double>::iterator t1 = t2-1;
155 
156    std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),eVecm[(*t1)].end(), theta);
157    std::vector<double>::iterator e11 = e12-1;
158         
159    std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),eVecm[(*t2)].end(), theta);
160    std::vector<double>::iterator e21 = e22-1;
161               
162    valueT1  =*t1;
163    valueT2  =*t2;
164    valueE21 =*e21;
165    valueE22 =*e22;
166    valueE12 =*e12;
167    valueE11 =*e11;
168
169    xs11 = eDiffCrossSectionData[0][valueT1][valueE11];
170    xs12 = eDiffCrossSectionData[0][valueT1][valueE12];
171    xs21 = eDiffCrossSectionData[0][valueT2][valueE21];
172    xs22 = eDiffCrossSectionData[0][valueT2][valueE22];
173  }
174     
175  G4double xsProduct = xs11 * xs12 * xs21 * xs22;
176 
177  if (xs11==0 || xs12==0 ||xs21==0 ||xs22==0) return (0.);
178     
179  if (xsProduct != 0.)
180  {
181    sigma = QuadInterpolator(  valueE11, valueE12, 
182                               valueE21, valueE22, 
183                               xs11, xs12, 
184                               xs21, xs22, 
185                               valueT1, valueT2, 
186                               k, theta );
187  }
188
189  return sigma;
190}
191
192//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
193
194G4double G4FinalStateElasticChampion::LinLogInterpolate(G4double e1, 
195                                                        G4double e2, 
196                                                        G4double e, 
197                                                        G4double xs1, 
198                                                        G4double xs2)
199{
200  G4double d1 = std::log(xs1);
201  G4double d2 = std::log(xs2);
202  G4double value = std::exp(d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
203  return value;
204}
205
206//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
207
208G4double G4FinalStateElasticChampion::LogLogInterpolate(G4double e1, 
209                                                        G4double e2, 
210                                                        G4double e, 
211                                                        G4double xs1, 
212                                                        G4double xs2)
213{
214  G4double a = (std::log10(xs2)-std::log10(xs1)) / (std::log10(e2)-std::log10(e1));
215  G4double b = std::log10(xs2) - a*std::log10(e2);
216  G4double sigma = a*std::log10(e) + b;
217  G4double value = (std::pow(10.,sigma));
218  return value;
219}
220
221//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
222
223G4double G4FinalStateElasticChampion::QuadInterpolator(G4double e11, G4double e12, 
224                                                       G4double e21, G4double e22, 
225                                                       G4double xs11, G4double xs12, 
226                                                       G4double xs21, G4double xs22, 
227                                                       G4double t1, G4double t2, 
228                                                       G4double t, G4double e)
229{
230// Log-Log
231/*
232  G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
233  G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
234  G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
235*/
236
237// Lin-Log
238  G4double interpolatedvalue1 = LinLogInterpolate(e11, e12, e, xs11, xs12);
239  G4double interpolatedvalue2 = LinLogInterpolate(e21, e22, e, xs21, xs22);
240  G4double value = LinLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
241  return value;
242}
243
244//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
245
246G4double G4FinalStateElasticChampion::RandomizeCosTheta(G4double k) 
247{
248
249 G4double cosTheta = -1;
250 G4double cumul = 0;
251 G4double value = 0;
252 
253 // Number of integration steps in the 0-180 deg range
254 G4int iMax=180;
255 
256 G4double random = G4UniformRand();
257 
258 // Cumulate differential cross section
259 for (G4int i=0; i<iMax; i++) 
260 {
261   cumul = cumul + DifferentialCrossSection(G4Electron::ElectronDefinition(),k/eV,G4double(i)*180./(iMax-1));
262 }
263 
264 // Select theta angle
265 for (G4int i=0; i<iMax; i++) 
266 {
267   value = value + DifferentialCrossSection(G4Electron::ElectronDefinition(),k/eV,G4double(i)*180./(iMax-1)) / cumul;
268   if (random < value) 
269   {
270     cosTheta=std::cos( deg*G4double(i)*180./(iMax-1) ); 
271     break;
272   } 
273 } 
274
275 return cosTheta;
276}
Note: See TracBrowser for help on using the repository browser.