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

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

update CVS release candidate geant4.9.3.01

File size: 9.6 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.10 2009/06/11 15:47:08 mantero Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
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   G4cout << G4endl;
88   G4cout << "*******************************************************************************" << G4endl;
89   G4cout << "*******************************************************************************" << G4endl;
90   G4cout << "   The class G4FinalStateElastiChampion is NOT SUPPORTED ANYMORE. " << G4endl;
91   G4cout << "   It will be REMOVED with the next major release of Geant4. " << G4endl;
92   G4cout << "   Please consult: https://twiki.cern.ch/twiki/bin/view/Geant4/LoweProcesses" << G4endl;
93   G4cout << "*******************************************************************************" << G4endl;
94   G4cout << "*******************************************************************************" << G4endl;
95   G4cout << G4endl;
96}
97
98//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99
100G4FinalStateElasticChampion::~G4FinalStateElasticChampion()
101{ 
102  eVecm.clear();
103}
104
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
106
107const G4FinalStateProduct& G4FinalStateElasticChampion::GenerateFinalState(const G4Track& track, const G4Step& )
108{
109  product.Clear();
110
111  G4double k = track.GetDynamicParticle()->GetKineticEnergy();
112
113  if (k>= lowEnergyLimit && k < highEnergyLimit)
114  { 
115    G4double cosTheta = RandomizeCosTheta(k);
116 
117    G4double phi = 2. * pi * G4UniformRand();
118
119    G4ThreeVector zVers = track.GetDynamicParticle()->GetMomentumDirection();
120    G4ThreeVector xVers = zVers.orthogonal();
121    G4ThreeVector yVers = zVers.cross(xVers);
122
123    G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
124    G4double yDir = xDir;
125    xDir *= std::cos(phi);
126    yDir *= std::sin(phi);
127
128    G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
129
130    product.ModifyPrimaryParticle(zPrimeVers,k);
131  }
132
133  if (k<lowEnergyLimit)
134  {
135    product.KillPrimaryParticle();
136  }
137 
138  return product;
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
142
143G4double G4FinalStateElasticChampion::DifferentialCrossSection
144  (G4ParticleDefinition * particleDefinition, G4double k, G4double theta)                                                         
145{
146  G4double sigma = 0.;
147  G4double valueT1 = 0;
148  G4double valueT2 = 0;
149  G4double valueE21 = 0;
150  G4double valueE22 = 0;
151  G4double valueE12 = 0;
152  G4double valueE11 = 0;
153  G4double xs11 = 0;   
154  G4double xs12 = 0; 
155  G4double xs21 = 0; 
156  G4double xs22 = 0; 
157
158  //SI : ensure the correct computation of cross section at the 180*deg limit
159  if (theta==180.) theta=theta-1e-9;
160
161  if (particleDefinition == G4Electron::ElectronDefinition()) 
162  {
163    std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k);
164    std::vector<double>::iterator t1 = t2-1;
165 
166    std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),eVecm[(*t1)].end(), theta);
167    std::vector<double>::iterator e11 = e12-1;
168         
169    std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),eVecm[(*t2)].end(), theta);
170    std::vector<double>::iterator e21 = e22-1;
171               
172    valueT1  =*t1;
173    valueT2  =*t2;
174    valueE21 =*e21;
175    valueE22 =*e22;
176    valueE12 =*e12;
177    valueE11 =*e11;
178
179    xs11 = eDiffCrossSectionData[0][valueT1][valueE11];
180    xs12 = eDiffCrossSectionData[0][valueT1][valueE12];
181    xs21 = eDiffCrossSectionData[0][valueT2][valueE21];
182    xs22 = eDiffCrossSectionData[0][valueT2][valueE22];
183  }
184     
185  G4double xsProduct = xs11 * xs12 * xs21 * xs22;
186 
187  if (xs11==0 || xs12==0 ||xs21==0 ||xs22==0) return (0.);
188     
189  if (xsProduct != 0.)
190  {
191    sigma = QuadInterpolator(  valueE11, valueE12, 
192                               valueE21, valueE22, 
193                               xs11, xs12, 
194                               xs21, xs22, 
195                               valueT1, valueT2, 
196                               k, theta );
197  }
198
199  return sigma;
200}
201
202//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
203
204G4double G4FinalStateElasticChampion::LinLogInterpolate(G4double e1, 
205                                                        G4double e2, 
206                                                        G4double e, 
207                                                        G4double xs1, 
208                                                        G4double xs2)
209{
210  G4double d1 = std::log(xs1);
211  G4double d2 = std::log(xs2);
212  G4double value = std::exp(d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
213  return value;
214}
215
216//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
217
218G4double G4FinalStateElasticChampion::LogLogInterpolate(G4double e1, 
219                                                        G4double e2, 
220                                                        G4double e, 
221                                                        G4double xs1, 
222                                                        G4double xs2)
223{
224  G4double a = (std::log10(xs2)-std::log10(xs1)) / (std::log10(e2)-std::log10(e1));
225  G4double b = std::log10(xs2) - a*std::log10(e2);
226  G4double sigma = a*std::log10(e) + b;
227  G4double value = (std::pow(10.,sigma));
228  return value;
229}
230
231//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
232
233G4double G4FinalStateElasticChampion::QuadInterpolator(G4double e11, G4double e12, 
234                                                       G4double e21, G4double e22, 
235                                                       G4double xs11, G4double xs12, 
236                                                       G4double xs21, G4double xs22, 
237                                                       G4double t1, G4double t2, 
238                                                       G4double t, G4double e)
239{
240// Log-Log
241/*
242  G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
243  G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
244  G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
245*/
246
247// Lin-Log
248  G4double interpolatedvalue1 = LinLogInterpolate(e11, e12, e, xs11, xs12);
249  G4double interpolatedvalue2 = LinLogInterpolate(e21, e22, e, xs21, xs22);
250  G4double value = LinLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
251  return value;
252}
253
254//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
255
256G4double G4FinalStateElasticChampion::RandomizeCosTheta(G4double k) 
257{
258
259 G4double cosTheta = -1;
260 G4double cumul = 0;
261 G4double value = 0;
262 
263 // Number of integration steps in the 0-180 deg range
264 G4int iMax=180;
265 
266 G4double random = G4UniformRand();
267 
268 // Cumulate differential cross section
269 for (G4int i=0; i<iMax; i++) 
270 {
271   cumul = cumul + DifferentialCrossSection(G4Electron::ElectronDefinition(),k/eV,G4double(i)*180./(iMax-1));
272 }
273 
274 // Select theta angle
275 for (G4int i=0; i<iMax; i++) 
276 {
277   value = value + DifferentialCrossSection(G4Electron::ElectronDefinition(),k/eV,G4double(i)*180./(iMax-1)) / cumul;
278   if (random < value) 
279   {
280     cosTheta=std::cos( deg*G4double(i)*180./(iMax-1) ); 
281     break;
282   } 
283 } 
284
285 return cosTheta;
286}
Note: See TracBrowser for help on using the repository browser.