source: trunk/source/processes/electromagnetic/lowenergy/src/G4AnalyticalEcpssrLiCrossSection.cc @ 1347

Last change on this file since 1347 was 1347, checked in by garnier, 13 years ago

geant4 tag 9.4

  • Property svn:executable set to *
File size: 32.3 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: G4AnalyticalEcpssrLiCrossSection.cc,v 1.4 2010/11/22 17:25:45 mantero Exp $
27// GEANT4 tag $Name: geant4-09-04-ref-00 $
28
29#include "globals.hh"
30#include "G4AnalyticalEcpssrLiCrossSection.hh"
31#include "G4AtomicTransitionManager.hh"
32#include "G4NistManager.hh"
33#include "G4Proton.hh"
34#include "G4Alpha.hh"
35#include <math.h>
36#include <iostream>
37#include "G4LinLogInterpolation.hh"
38
39//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
40
41G4AnalyticalEcpssrLiCrossSection::G4AnalyticalEcpssrLiCrossSection()
42{
43
44    // Storing FLi data needed for 0.2 to 3.0  velocities region
45
46    char *path = getenv("G4LEDATA");
47
48    if (!path)
49    G4Exception("G4ecpssrLCrossSection::G4AnalyticalEcpssrLiCrossSection() G4LEDDATA environment variable not set");
50
51    std::ostringstream fileName1;
52    std::ostringstream fileName2;
53   
54    fileName1 << path << "/pixe/uf/FL1.dat";
55    fileName2 << path << "/pixe/uf/FL2.dat";
56
57    // Reading of FL1.dat
58
59    std::ifstream FL1(fileName1.str().c_str());
60    if (!FL1) G4Exception("G4ecpssrLCrossSection::G4AnalyticalEcpssrLiCrossSection() error opening FL1 data file");
61 
62    dummyVec1.push_back(0.);
63
64    while(!FL1.eof())
65    {
66        double x1;
67        double y1;
68
69        FL1>>x1>>y1;
70
71        //  Mandatory vector initialization
72        if (x1 != dummyVec1.back())
73        {
74          dummyVec1.push_back(x1);
75          aVecMap1[x1].push_back(-1.);
76        }
77
78        FL1>>FL1Data[x1][y1];
79
80        if (y1 != aVecMap1[x1].back()) aVecMap1[x1].push_back(y1);
81    }
82
83    // Reading of FL2.dat
84   
85    std::ifstream FL2(fileName2.str().c_str());
86    if (!FL2) G4Exception("G4ecpssrLCrossSection::G4AnalyticalEcpssrLiCrossSection() error opening FL2 data file");
87   
88    dummyVec2.push_back(0.);
89
90    while(!FL2.eof())
91    {
92        double x2;
93        double y2;
94
95        FL2>>x2>>y2;
96
97        //  Mandatory vector initialization
98        if (x2 != dummyVec2.back())
99        {
100          dummyVec2.push_back(x2);
101          aVecMap2[x2].push_back(-1.);
102        }
103
104        FL2>>FL2Data[x2][y2];
105
106        if (y2 != aVecMap2[x2].back()) aVecMap2[x2].push_back(y2);
107    }
108
109    // Verbose level
110    verboseLevel=0;
111
112}
113
114//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
115
116G4AnalyticalEcpssrLiCrossSection::~G4AnalyticalEcpssrLiCrossSection()
117{}
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
120
121G4double G4AnalyticalEcpssrLiCrossSection::ExpIntFunction(G4int n,G4double x)
122
123{
124// this function allows fast evaluation of the n order exponential integral function En(x)
125
126  G4int i;
127  G4int ii;
128  G4int nm1;
129  G4double a;
130  G4double b;
131  G4double c;
132  G4double d;
133  G4double del;
134  G4double fact;
135  G4double h;
136  G4double psi;
137  G4double ans = 0;
138  const G4double euler= 0.5772156649;
139  const G4int maxit= 100;
140  const G4double fpmin = 1.0e-30;
141  const G4double eps = 1.0e-7;
142  nm1=n-1;
143  if (n<0 || x<0.0 || (x==0.0 && (n==0 || n==1)))
144  G4cout << "bad arguments in ExpIntFunction" << G4endl;
145  else {
146       if (n==0) ans=std::exp(-x)/x;
147        else {
148           if (x==0.0) ans=1.0/nm1;
149              else {
150                   if (x > 1.0) {
151                                b=x+n;
152                                c=1.0/fpmin;
153                                d=1.0/b;
154                                h=d;
155                                for (i=1;i<=maxit;i++) {
156                                  a=-i*(nm1+i);
157                                  b +=2.0;
158                                  d=1.0/(a*d+b);
159                                  c=b+a/c;
160                                  del=c*d;
161                                  h *=del;
162                                      if (std::fabs(del-1.0) < eps) {
163                                        ans=h*std::exp(-x);
164                                        return ans;
165                                      }
166                                }
167                   } else {
168                     ans = (nm1!=0 ? 1.0/nm1 : -std::log(x)-euler);
169                     fact=1.0;
170                     for (i=1;i<=maxit;i++) {
171                       fact *=-x/i;
172                       if (i !=nm1) del = -fact/(i-nm1);
173                       else {
174                         psi = -euler;
175                         for (ii=1;ii<=nm1;ii++) psi +=1.0/ii;
176                         del=fact*(-std::log(x)+psi);
177                       }
178                       ans += del;
179                       if (std::fabs(del) < std::fabs(ans)*eps) return ans;
180                     }
181                   }
182              }
183        }
184  }
185return ans;
186}
187
188//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
189
190G4double G4AnalyticalEcpssrLiCrossSection::CalculateL1CrossSection(G4int zTarget,G4double massIncident, G4double energyIncident)
191{
192
193 //this L1-CrossSection calculation method is done according to Werner Brandt and Grzegorz Lapicki, Phys.Rev.A20 N2 (1979),
194 //and using data tables of O. Benka et al. At.Data Nucl.Data Tables Vol.22 No.3 (1978).
195
196  G4NistManager* massManager = G4NistManager::Instance();
197
198  G4AtomicTransitionManager*  transitionManager =  G4AtomicTransitionManager::Instance();
199
200  G4double  zIncident = 0;
201  G4Proton* aProtone = G4Proton::Proton();
202  G4Alpha* aAlpha = G4Alpha::Alpha();
203
204  if (massIncident == aProtone->GetPDGMass() )
205
206   zIncident = (aProtone->GetPDGCharge())/eplus;
207
208  else
209    {
210      if (massIncident == aAlpha->GetPDGMass())
211
212          zIncident  = (aAlpha->GetPDGCharge())/eplus;
213
214      else
215        {
216          G4cout << "*** WARNING in G4AnalyticalEcpssrLiCrossSection::CalculateL1CrossSection : Proton or Alpha incident particles only. " << G4endl;
217          G4cout << massIncident << ", " << aAlpha->GetPDGMass() << " (alpha)" << aProtone->GetPDGMass() << " (proton)" << G4endl;
218          return 0;
219        }
220    }
221
222  G4double l1BindingEnergy = transitionManager->Shell(zTarget,1)->BindingEnergy(); //Observed binding energy of L1-subshell
223
224  G4double massTarget = (massManager->GetAtomicMassAmu(zTarget))*amu_c2;
225
226  G4double systemMass =((massIncident*massTarget)/(massIncident+massTarget))/electron_mass_c2; //Mass of the system (projectile, target)
227
228  const G4double zlshell= 4.15;
229
230  G4double screenedzTarget = zTarget-zlshell; //Effective nuclear charge as seen by electrons in L1-sub shell
231
232  const G4double rydbergMeV= 13.6056923e-6;
233
234  const G4double nl= 2.;
235
236  G4double tetal1 = (l1BindingEnergy*nl*nl)/((screenedzTarget*screenedzTarget)*rydbergMeV); //Screening parameter
237
238    if (verboseLevel>0) G4cout << "  tetal1=" <<  tetal1<< G4endl;
239
240  G4double reducedEnergy = (energyIncident*electron_mass_c2)/(massIncident*rydbergMeV*screenedzTarget*screenedzTarget);
241
242  const G4double bohrPow2Barn=(Bohr_radius*Bohr_radius)/barn ; //Bohr radius of hydrogen
243
244  G4double sigma0 = 8.*pi*(zIncident*zIncident)*bohrPow2Barn*std::pow(screenedzTarget,-4.);
245
246  G4double velocityl1 = CalculateVelocity(1, zTarget, massIncident, energyIncident); // Scaled velocity
247
248    if (verboseLevel>0) G4cout << "  velocityl1=" << velocityl1<< G4endl;
249
250  const G4double l1AnalyticalApproximation= 1.5;
251
252  G4double x1 =(nl*l1AnalyticalApproximation)/velocityl1;
253
254   if (verboseLevel>0) G4cout << "  x1=" << x1<< G4endl;
255
256  G4double electrIonizationEnergyl1=0.;
257
258  if ( x1<=0.035)  electrIonizationEnergyl1= 0.75*pi*(std::log(1./(x1*x1))-1.);
259  else
260    {
261      if ( x1<=3.)
262        electrIonizationEnergyl1 =std::exp(-2.*x1)/(0.031+(0.213*std::pow(x1,0.5))+(0.005*x1)-(0.069*std::pow(x1,3./2.))+(0.324*x1*x1));
263      else
264        {if ( x1<=11.) electrIonizationEnergyl1 =2.*std::exp(-2.*x1)/std::pow(x1,1.6);}
265    }
266
267  G4double hFunctionl1 =(electrIonizationEnergyl1*2.*nl)/(tetal1*std::pow(velocityl1,3)); //takes into account the polarization effect
268
269    if (verboseLevel>0) G4cout << "  hFunctionl1=" << hFunctionl1<< G4endl;
270
271  G4double gFunctionl1 = (1.+(9.*velocityl1)+(31.*velocityl1*velocityl1)+(49.*std::pow(velocityl1,3.))+(162.*std::pow(velocityl1,4.))+(63.*std::pow(velocityl1,5.))+(18.*std::pow(velocityl1,6.))+(1.97*std::pow(velocityl1,7.)))/std::pow(1.+velocityl1,9.);//takes into account the reduced binding effect
272
273    if (verboseLevel>0) G4cout << "  gFunctionl1=" << gFunctionl1<< G4endl;
274
275  G4double sigmaPSS_l1 = 1.+(((2.*zIncident)/(screenedzTarget*tetal1))*(gFunctionl1-hFunctionl1)); //Binding-polarization factor
276
277    if (verboseLevel>0) G4cout << "sigmaPSS_l1 =" << sigmaPSS_l1<< G4endl;
278
279const G4double cNaturalUnit= 137.;
280 
281G4double yl1Formula=0.4*(screenedzTarget/cNaturalUnit)*(screenedzTarget/cNaturalUnit)/(nl*velocityl1/sigmaPSS_l1);
282   
283G4double l1relativityCorrection = std::pow((1.+(1.1*yl1Formula*yl1Formula)),0.5)+yl1Formula; // Relativistic correction parameter
284 
285//G4double reducedVelocity_l1 = velocityl1*std::pow(l1relativityCorrection,0.5); //Reduced velocity parameter
286
287
288G4double L1etaOverTheta2;
289
290G4double  universalFunction_l1 = 0.;
291
292G4double sigmaPSSR_l1;
293
294if ( velocityl1 <5.  )
295  {
296
297L1etaOverTheta2 =(reducedEnergy* l1relativityCorrection)/((tetal1*sigmaPSS_l1)*(tetal1*sigmaPSS_l1));
298   
299 if ( ((tetal1*sigmaPSS_l1) >=0.2) && ((tetal1*sigmaPSS_l1) <=2.6670) && (L1etaOverTheta2>=0.1e-3) && (L1etaOverTheta2<=0.866e2) )
300 
301universalFunction_l1 = FunctionFL1((tetal1*sigmaPSS_l1),L1etaOverTheta2);
302
303 if (verboseLevel>0) G4cout << "at low velocity range, universalFunction_l1  =" << universalFunction_l1 << G4endl;
304
305sigmaPSSR_l1 = (sigma0/(tetal1*sigmaPSS_l1))*universalFunction_l1;// Plane-wave Born -Aproximation L1-subshell ionisation Cross Section
306
307 if (verboseLevel>0) G4cout << "  at low velocity range, sigma PWBA L1 CS  = " << sigmaPSSR_l1<< G4endl;
308
309
310}
311
312else
313
314{
315
316 L1etaOverTheta2 = reducedEnergy/(tetal1*tetal1);
317
318 if ( (tetal1 >=0.2) && (tetal1 <=2.6670) && (L1etaOverTheta2>=0.1e-3) && (L1etaOverTheta2<=0.866e2) ) 
319   
320universalFunction_l1 = FunctionFL1(tetal1,L1etaOverTheta2);
321
322if (verboseLevel>0) G4cout << "at medium and high velocity range, universalFunction_l1  =" << universalFunction_l1 << G4endl;
323
324sigmaPSSR_l1 = (sigma0/tetal1)*universalFunction_l1;// Plane-wave Born -Aproximation L1-subshell ionisation Cross Section
325
326 if (verboseLevel>0) G4cout << "  sigma PWBA L1 CS at medium and high velocity range = " << sigmaPSSR_l1<< G4endl;   
327}
328 
329 G4double pssDeltal1 = (4./(systemMass *sigmaPSS_l1*tetal1))*(sigmaPSS_l1/velocityl1)*(sigmaPSS_l1/velocityl1);
330
331    if (verboseLevel>0) G4cout << "  pssDeltal1=" << pssDeltal1<< G4endl;
332
333  G4double energyLossl1 = std::pow(1-pssDeltal1,0.5);
334
335    if (verboseLevel>0) G4cout << "  energyLossl1=" << energyLossl1<< G4endl;
336
337  G4double coulombDeflectionl1 = 
338(8.*pi*zIncident/systemMass)*std::pow(tetal1*sigmaPSS_l1,-2.)*std::pow(velocityl1/sigmaPSS_l1,-3.)*(zTarget/screenedzTarget);
339
340 G4double cParameterl1 =2.* coulombDeflectionl1/(energyLossl1*(energyLossl1+1.));
341
342  G4double coulombDeflectionFunction_l1 = 9.*ExpIntFunction(10,cParameterl1); //Coulomb-deflection effect correction
343
344 if (verboseLevel>0) G4cout << "  coulombDeflectionFunction_l1 =" << coulombDeflectionFunction_l1 << G4endl;
345
346G4double crossSection_L1 = coulombDeflectionFunction_l1 * sigmaPSSR_l1;
347
348  //ECPSSR L1 -subshell cross section is estimated at perturbed-stationnairy-state(PSS)
349  //and reduced by the energy-loss(E),the Coulomb deflection(C),and the relativity(R) effects
350
351    if (verboseLevel>0) G4cout << "  crossSection_L1 =" << crossSection_L1 << G4endl;
352
353  if (crossSection_L1 >= 0) {
354
355    return crossSection_L1 * barn;
356  }
357 
358  else {return 0;}
359}
360
361//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
362
363G4double G4AnalyticalEcpssrLiCrossSection::CalculateL2CrossSection(G4int zTarget,G4double massIncident, G4double energyIncident)
364
365{
366
367  // this L2-CrossSection calculation method is done according to Werner Brandt and Grzegorz Lapicki, Phys.Rev.A20 N2 (1979),
368  // and using data tables of O. Benka et al. At.Data Nucl.Data Tables Vol.22 No.3 (1978).
369
370  G4NistManager* massManager = G4NistManager::Instance();
371
372  G4AtomicTransitionManager*  transitionManager =  G4AtomicTransitionManager::Instance();
373
374  G4double  zIncident = 0;
375 
376  G4Proton* aProtone = G4Proton::Proton();
377  G4Alpha* aAlpha = G4Alpha::Alpha();
378
379  if (massIncident == aProtone->GetPDGMass() )
380
381   zIncident = (aProtone->GetPDGCharge())/eplus;
382
383  else
384    {
385      if (massIncident == aAlpha->GetPDGMass())
386
387          zIncident  = (aAlpha->GetPDGCharge())/eplus;
388
389      else
390        {
391          G4cout << "*** WARNING in G4AnalyticalEcpssrLiCrossSection::CalculateL2CrossSection : Proton or Alpha incident particles only. " << G4endl;
392          G4cout << massIncident << ", " << aAlpha->GetPDGMass() << " (alpha)" << aProtone->GetPDGMass() << " (proton)" << G4endl;
393          return 0;
394        }
395    }
396
397  G4double l2BindingEnergy = transitionManager->Shell(zTarget,2)->BindingEnergy(); //Observed binding energy of L2-subshell
398
399  G4double massTarget = (massManager->GetAtomicMassAmu(zTarget))*amu_c2;
400
401  G4double systemMass =((massIncident*massTarget)/(massIncident+massTarget))/electron_mass_c2; //Mass of the system (projectile, target)
402
403  const G4double zlshell= 4.15;
404
405  G4double screenedzTarget = zTarget-zlshell; //Effective nuclear charge as seen by electrons in L2-subshell
406
407  const G4double rydbergMeV= 13.6056923e-6;
408
409  const G4double nl= 2.;
410
411  G4double tetal2 = (l2BindingEnergy*nl*nl)/((screenedzTarget*screenedzTarget)*rydbergMeV); //Screening parameter
412
413    if (verboseLevel>0) G4cout << "  tetal2=" <<  tetal2<< G4endl;
414
415  G4double reducedEnergy = (energyIncident*electron_mass_c2)/(massIncident*rydbergMeV*screenedzTarget*screenedzTarget); 
416
417  const G4double bohrPow2Barn=(Bohr_radius*Bohr_radius)/barn ; //Bohr radius of hydrogen
418
419  G4double sigma0 = 8.*pi*(zIncident*zIncident)*bohrPow2Barn*std::pow(screenedzTarget,-4.);
420
421  G4double velocityl2 = CalculateVelocity(2,  zTarget, massIncident, energyIncident); // Scaled velocity
422
423    if (verboseLevel>0) G4cout << "  velocityl2=" << velocityl2<< G4endl;
424
425  const G4double l23AnalyticalApproximation= 1.25;
426
427  G4double x2 = (nl*l23AnalyticalApproximation)/velocityl2;
428   
429    if (verboseLevel>0) G4cout << "  x2=" << x2<< G4endl;
430
431  G4double electrIonizationEnergyl2=0.;
432
433  if ( x2<=0.035)  electrIonizationEnergyl2= 0.75*pi*(std::log(1./(x2*x2))-1.);
434  else
435    {
436      if ( x2<=3.)
437        electrIonizationEnergyl2 =std::exp(-2.*x2)/(0.031+(0.210*std::pow(x2,0.5))+(0.005*x2)-(0.069*std::pow(x2,3./2.))+(0.324*x2*x2));
438      else
439        {if ( x2<=11.) electrIonizationEnergyl2 =2.*std::exp(-2.*x2)/std::pow(x2,1.6);  }
440    }
441
442 G4double hFunctionl2 =(electrIonizationEnergyl2*2.*nl)/(tetal2*std::pow(velocityl2,3)); //takes into account  the polarization effect
443
444   if (verboseLevel>0) G4cout << "  hFunctionl2=" << hFunctionl2<< G4endl;
445
446  G4double gFunctionl2 = (1.+(10.*velocityl2)+(45.*velocityl2*velocityl2)+(102.*std::pow(velocityl2,3.))+(331.*std::pow(velocityl2,4.))+(6.7*std::pow(velocityl2,5.))+(58.*std::pow(velocityl2,6.))+(7.8*std::pow(velocityl2,7.))+ (0.888*std::pow(velocityl2,8.)) )/std::pow(1.+velocityl2,10.);
447  //takes into account the reduced binding effect
448
449    if (verboseLevel>0) G4cout << "  gFunctionl2=" << gFunctionl2<< G4endl;
450
451  G4double sigmaPSS_l2 = 1.+(((2.*zIncident)/(screenedzTarget*tetal2))*(gFunctionl2-hFunctionl2)); //Binding-polarization factor
452
453    if (verboseLevel>0) G4cout << "  sigmaPSS_l2=" << sigmaPSS_l2<< G4endl;
454
455const G4double cNaturalUnit= 137.;
456     
457G4double yl2Formula=0.15*(screenedzTarget/cNaturalUnit)*(screenedzTarget/cNaturalUnit)/(velocityl2/sigmaPSS_l2);
458 
459G4double l2relativityCorrection = std::pow((1.+(1.1*yl2Formula*yl2Formula)),0.5)+yl2Formula; // Relativistic correction parameter
460   
461
462G4double L2etaOverTheta2;
463
464G4double  universalFunction_l2 = 0.;
465
466G4double sigmaPSSR_l2 ;
467 
468if ( velocityl2 < 5. )
469  {
470
471L2etaOverTheta2 = (reducedEnergy*l2relativityCorrection)/((sigmaPSS_l2*tetal2)*(sigmaPSS_l2*tetal2));
472
473if ( (tetal2*sigmaPSS_l2>=0.2) && (tetal2*sigmaPSS_l2<=2.6670) && (L2etaOverTheta2>=0.1e-3) && (L2etaOverTheta2<=0.866e2) )
474
475universalFunction_l2 = FunctionFL2((tetal2*sigmaPSS_l2),L2etaOverTheta2);
476
477sigmaPSSR_l2 = (sigma0/(tetal2*sigmaPSS_l2))*universalFunction_l2;
478
479 if (verboseLevel>0) G4cout << "  sigma PWBA L2 CS at low velocity range = " << sigmaPSSR_l2<< G4endl;
480
481}   
482
483else
484
485{ 
486   
487L2etaOverTheta2 = reducedEnergy /(tetal2*tetal2);
488
489if ( (tetal2>=0.2) && (tetal2<=2.6670) && (L2etaOverTheta2>=0.1e-3) && (L2etaOverTheta2<=0.866e2) )
490
491universalFunction_l2 = FunctionFL2((tetal2),L2etaOverTheta2);
492   
493sigmaPSSR_l2 = (sigma0/tetal2)*universalFunction_l2;
494
495 if (verboseLevel>0) G4cout << "  sigma PWBA L2 CS at medium and high velocity range = " << sigmaPSSR_l2<< G4endl;
496
497}
498 
499G4double pssDeltal2 = (4./(systemMass*sigmaPSS_l2*tetal2))*(sigmaPSS_l2/velocityl2)*(sigmaPSS_l2/velocityl2);
500
501G4double energyLossl2 = std::pow(1-pssDeltal2,0.5);
502 
503  if (verboseLevel>0) G4cout << "  energyLossl2=" << energyLossl2<< G4endl;
504
505  G4double coulombDeflectionl2
506=(8.*pi*zIncident/systemMass)*std::pow(tetal2*sigmaPSS_l2,-2.)*std::pow(velocityl2/sigmaPSS_l2,-3.)*(zTarget/screenedzTarget);
507
508  G4double cParameterl2 = 2.*coulombDeflectionl2/(energyLossl2*(energyLossl2+1.));
509
510  G4double coulombDeflectionFunction_l2 = 11.*ExpIntFunction(12,cParameterl2); //Coulomb-deflection effect correction
511
512    if (verboseLevel>0) G4cout << "  coulombDeflectionFunction_l2 =" << coulombDeflectionFunction_l2 << G4endl;
513
514  G4double crossSection_L2 = coulombDeflectionFunction_l2 * sigmaPSSR_l2;
515  //ECPSSR L2 -subshell cross section is estimated at perturbed-stationnairy-state(PSS)
516  //and reduced by the energy-loss(E),the Coulomb deflection(C),and the relativity(R) effects
517
518    if (verboseLevel>0) G4cout << "  crossSection_L2 =" << crossSection_L2 << G4endl;
519
520  if (crossSection_L2 >= 0) {
521     return crossSection_L2 * barn;
522      }
523  else {return 0;}
524}
525
526//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
527
528
529G4double G4AnalyticalEcpssrLiCrossSection::CalculateL3CrossSection(G4int zTarget,G4double massIncident, G4double energyIncident)
530
531{
532
533  //this L3-CrossSection calculation method is done according to Werner Brandt and Grzegorz Lapicki, Phys.Rev.A20 N2 (1979),
534  //and using data tables of O. Benka et al. At.Data Nucl.Data Tables Vol.22 No.3 (1978).
535
536  G4NistManager* massManager = G4NistManager::Instance();
537
538  G4AtomicTransitionManager*  transitionManager =  G4AtomicTransitionManager::Instance();
539
540  G4double  zIncident = 0;
541
542  G4Proton* aProtone = G4Proton::Proton();
543  G4Alpha* aAlpha = G4Alpha::Alpha();
544
545  if (massIncident == aProtone->GetPDGMass() )
546
547   zIncident = (aProtone->GetPDGCharge())/eplus;
548
549  else
550    {
551      if (massIncident == aAlpha->GetPDGMass())
552
553          zIncident  = (aAlpha->GetPDGCharge())/eplus;
554
555      else
556        {
557          G4cout << "*** WARNING in G4AnalyticalEcpssrLiCrossSection::CalculateL3CrossSection : Proton or Alpha incident particles only. " << G4endl;
558          G4cout << massIncident << ", " << aAlpha->GetPDGMass() << " (alpha)" << aProtone->GetPDGMass() << " (proton)" << G4endl;
559          return 0;
560        }
561    }
562
563  G4double l3BindingEnergy = transitionManager->Shell(zTarget,3)->BindingEnergy();
564
565  G4double massTarget = (massManager->GetAtomicMassAmu(zTarget))*amu_c2;
566
567  G4double systemMass =((massIncident*massTarget)/(massIncident+massTarget))/electron_mass_c2;//Mass of the system (projectile, target)
568
569  const G4double zlshell= 4.15;
570
571  G4double screenedzTarget = zTarget-zlshell;//Effective nuclear charge as seen by electrons in L3-subshell
572
573  const G4double rydbergMeV= 13.6056923e-6;
574
575  const G4double nl= 2.;
576
577  G4double tetal3 = (l3BindingEnergy*nl*nl)/((screenedzTarget*screenedzTarget)*rydbergMeV);//Screening parameter
578
579    if (verboseLevel>0) G4cout << "  tetal3=" <<  tetal3<< G4endl;
580
581  G4double reducedEnergy = (energyIncident*electron_mass_c2)/(massIncident*rydbergMeV*screenedzTarget*screenedzTarget);
582
583  const G4double bohrPow2Barn=(Bohr_radius*Bohr_radius)/barn ;//Bohr radius of hydrogen
584
585  G4double sigma0 = 8.*pi*(zIncident*zIncident)*bohrPow2Barn*std::pow(screenedzTarget,-4.);
586
587  G4double velocityl3 = CalculateVelocity(3, zTarget, massIncident, energyIncident);// Scaled velocity
588
589    if (verboseLevel>0) G4cout << "  velocityl3=" << velocityl3<< G4endl;
590
591  const G4double l23AnalyticalApproximation= 1.25;
592
593  G4double x3 = (nl*l23AnalyticalApproximation)/velocityl3;
594
595    if (verboseLevel>0) G4cout << "  x3=" << x3<< G4endl;
596
597  G4double electrIonizationEnergyl3=0.;
598
599  if ( x3<=0.035)  electrIonizationEnergyl3= 0.75*pi*(std::log(1./(x3*x3))-1.);
600    else
601    {
602      if ( x3<=3.) electrIonizationEnergyl3 =std::exp(-2.*x3)/(0.031+(0.210*std::pow(x3,0.5))+(0.005*x3)-(0.069*std::pow(x3,3./2.))+(0.324*x3*x3));
603      else
604        {
605          if ( x3<=11.) electrIonizationEnergyl3 =2.*std::exp(-2.*x3)/std::pow(x3,1.6);}
606    }
607
608  G4double hFunctionl3 =(electrIonizationEnergyl3*2.*nl)/(tetal3*std::pow(velocityl3,3));//takes into account the polarization effect
609
610    if (verboseLevel>0) G4cout << "  hFunctionl3=" << hFunctionl3<< G4endl;
611
612  G4double gFunctionl3 = (1.+(10.*velocityl3)+(45.*velocityl3*velocityl3)+(102.*std::pow(velocityl3,3.))+(331.*std::pow(velocityl3,4.))+(6.7*std::pow(velocityl3,5.))+(58.*std::pow(velocityl3,6.))+(7.8*std::pow(velocityl3,7.))+ (0.888*std::pow(velocityl3,8.)) )/std::pow(1.+velocityl3,10.);
613  //takes into account the reduced binding effect
614
615    if (verboseLevel>0) G4cout << "  gFunctionl3=" << gFunctionl3<< G4endl;
616
617  G4double sigmaPSS_l3 = 1.+(((2.*zIncident)/(screenedzTarget*tetal3))*(gFunctionl3-hFunctionl3));//Binding-polarization factor
618
619    if (verboseLevel>0) G4cout << "sigmaPSS_l3 =" << sigmaPSS_l3<< G4endl;
620
621const G4double cNaturalUnit= 137.;
622       
623G4double yl3Formula=0.15*(screenedzTarget/cNaturalUnit)*(screenedzTarget/cNaturalUnit)/(velocityl3/sigmaPSS_l3);
624       
625G4double l3relativityCorrection = std::pow((1.+(1.1*yl3Formula*yl3Formula)),0.5)+yl3Formula; // Relativistic correction parameter
626 
627G4double L3etaOverTheta2;
628 
629G4double  universalFunction_l3 = 0.;
630 
631G4double sigmaPSSR_l3;
632
633if ( velocityl3 < 5. )
634  {
635
636L3etaOverTheta2 = (reducedEnergy* l3relativityCorrection)/((sigmaPSS_l3*tetal3)*(sigmaPSS_l3*tetal3));
637
638if ( (tetal3*sigmaPSS_l3>=0.2) && (tetal3*sigmaPSS_l3<=2.6670) && (L3etaOverTheta2>=0.1e-3) && (L3etaOverTheta2<=0.866e2) )
639
640universalFunction_l3 = 2.*FunctionFL2((tetal3*sigmaPSS_l3), L3etaOverTheta2  );
641
642sigmaPSSR_l3 = (sigma0/(tetal3*sigmaPSS_l3))*universalFunction_l3;
643
644 if (verboseLevel>0) G4cout << "  sigma PWBA L3 CS at low velocity range = " << sigmaPSSR_l3<< G4endl;
645
646}
647
648else 
649
650{
651
652 L3etaOverTheta2 = reducedEnergy/(tetal3*tetal3);
653
654if ( (tetal3>=0.2) && (tetal3<=2.6670) && (L3etaOverTheta2>=0.1e-3) && (L3etaOverTheta2<=0.866e2) ) 
655     
656universalFunction_l3 = 2.*FunctionFL2(tetal3, L3etaOverTheta2  );
657
658 sigmaPSSR_l3 = (sigma0/tetal3)*universalFunction_l3;
659
660 if (verboseLevel>0) G4cout << "  sigma PWBA L3 CS at medium and high velocity range = " << sigmaPSSR_l3<< G4endl;
661
662}
663 
664G4double pssDeltal3 = (4./(systemMass*sigmaPSS_l3*tetal3))*(sigmaPSS_l3/velocityl3)*(sigmaPSS_l3/velocityl3);
665
666if (verboseLevel>0) G4cout << "  pssDeltal3=" << pssDeltal3<< G4endl;
667
668G4double energyLossl3 = std::pow(1-pssDeltal3,0.5);
669
670if (verboseLevel>0) G4cout << "  energyLossl3=" << energyLossl3<< G4endl;
671
672  G4double coulombDeflectionl3 = 
673(8.*pi*zIncident/systemMass)*std::pow(tetal3*sigmaPSS_l3,-2.)*std::pow(velocityl3/sigmaPSS_l3,-3.)*(zTarget/screenedzTarget);
674
675  G4double cParameterl3 = 2.*coulombDeflectionl3/(energyLossl3*(energyLossl3+1.));
676
677  G4double coulombDeflectionFunction_l3 = 11.*ExpIntFunction(12,cParameterl3);//Coulomb-deflection effect correction
678
679    if (verboseLevel>0) G4cout << "  coulombDeflectionFunction_l3 =" << coulombDeflectionFunction_l3 << G4endl;
680
681  G4double crossSection_L3 =  coulombDeflectionFunction_l3 * sigmaPSSR_l3;
682  //ECPSSR L3 -subshell cross section is estimated at perturbed-stationnairy-state(PSS)
683  //and reduced by the energy-loss(E),the Coulomb deflection(C),and the relativity(R) effects
684
685    if (verboseLevel>0) G4cout << "  crossSection_L3 =" << crossSection_L3 << G4endl;
686 
687  if (crossSection_L3 >= 0) {
688    return crossSection_L3 * barn;
689  }
690  else {return 0;}
691}
692
693//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
694
695G4double G4AnalyticalEcpssrLiCrossSection::CalculateVelocity(G4int subShell, G4int zTarget, G4double massIncident,  G4double energyIncident)
696
697{
698
699  G4AtomicTransitionManager*  transitionManager =  G4AtomicTransitionManager::Instance();
700
701  G4double liBindingEnergy = transitionManager->Shell(zTarget,subShell)->BindingEnergy();
702
703  G4Proton* aProtone = G4Proton::Proton();
704  G4Alpha* aAlpha = G4Alpha::Alpha();
705
706  if (!((massIncident == aProtone->GetPDGMass()) || (massIncident == aAlpha->GetPDGMass())))
707    {
708      G4cout << "*** WARNING in G4AnalyticalEcpssrLiCrossSection::CalculateVelocity : Proton or Alpha incident particles only. " << G4endl;
709      G4cout << massIncident << ", " << aAlpha->GetPDGMass() << " (alpha)" << aProtone->GetPDGMass() << " (proton)" << G4endl;
710      return 0;
711    }
712
713  const G4double zlshell= 4.15;
714
715  G4double screenedzTarget = zTarget- zlshell;
716
717  const G4double rydbergMeV= 13.6056923e-6;
718
719  const G4double nl= 2.;
720
721  G4double tetali = (liBindingEnergy*nl*nl)/(screenedzTarget*screenedzTarget*rydbergMeV);
722
723  G4double reducedEnergy = (energyIncident*electron_mass_c2)/(massIncident*rydbergMeV*screenedzTarget*screenedzTarget);
724
725  G4double velocity = 2.*nl*std::pow(reducedEnergy,0.5)/tetali;
726
727  return velocity;
728}
729
730//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
731
732G4double G4AnalyticalEcpssrLiCrossSection::FunctionFL1(G4double k, G4double theta)
733{
734
735  G4double sigma = 0.;
736  G4double valueT1 = 0;
737  G4double valueT2 = 0;
738  G4double valueE21 = 0;
739  G4double valueE22 = 0;
740  G4double valueE12 = 0;
741  G4double valueE11 = 0;
742  G4double xs11 = 0;
743  G4double xs12 = 0;
744  G4double xs21 = 0;
745  G4double xs22 = 0;
746
747  // PROTECTION TO ALLOW INTERPOLATION AT MINIMUM AND MAXIMUM Eta/Theta2 values
748
749  if (
750       theta==8.66e-4 ||
751       theta==8.66e-3 ||
752       theta==8.66e-2 ||
753       theta==8.66e-1 ||
754       theta==8.66e+00 ||
755       theta==8.66e+01
756  ) theta=theta-1e-12;
757
758  if (
759       theta==1.e-4 ||
760       theta==1.e-3 ||
761       theta==1.e-2 ||
762       theta==1.e-1 ||
763       theta==1.e+00 ||
764       theta==1.e+01
765  ) theta=theta+1e-12;
766
767  // END PROTECTION
768
769    std::vector<double>::iterator t2 = std::upper_bound(dummyVec1.begin(),dummyVec1.end(), k);
770    std::vector<double>::iterator t1 = t2-1;
771
772    std::vector<double>::iterator e12 = std::upper_bound(aVecMap1[(*t1)].begin(),aVecMap1[(*t1)].end(), theta);
773    std::vector<double>::iterator e11 = e12-1;
774
775    std::vector<double>::iterator e22 = std::upper_bound(aVecMap1[(*t2)].begin(),aVecMap1[(*t2)].end(), theta);
776    std::vector<double>::iterator e21 = e22-1;
777
778    valueT1  =*t1;
779    valueT2  =*t2;
780    valueE21 =*e21;
781    valueE22 =*e22;
782    valueE12 =*e12;
783    valueE11 =*e11;
784
785    xs11 = FL1Data[valueT1][valueE11];
786    xs12 = FL1Data[valueT1][valueE12];
787    xs21 = FL1Data[valueT2][valueE21];
788    xs22 = FL1Data[valueT2][valueE22];
789
790  if (verboseLevel>0) 
791    G4cout
792    << valueT1 << " "
793    << valueT2 << " "
794    << valueE11 << " "
795    << valueE12 << " "
796    << valueE21 << " "
797    << valueE22 << " "
798    << xs11 << " "
799    << xs12 << " "
800    << xs21 << " "
801    << xs22 << " "
802    << G4endl;
803
804  G4double xsProduct = xs11 * xs12 * xs21 * xs22;
805
806  if (xs11==0 || xs12==0 ||xs21==0 ||xs22==0) return (0.);
807
808  if (xsProduct != 0.)
809  {
810    sigma = QuadInterpolator(  valueE11, valueE12,
811                               valueE21, valueE22,
812                               xs11, xs12,
813                               xs21, xs22,
814                               valueT1, valueT2,
815                               k, theta );
816  }
817
818  return sigma;
819}
820
821//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
822
823G4double G4AnalyticalEcpssrLiCrossSection::FunctionFL2(G4double k, G4double theta)
824{
825
826  G4double sigma = 0.;
827  G4double valueT1 = 0;
828  G4double valueT2 = 0;
829  G4double valueE21 = 0;
830  G4double valueE22 = 0;
831  G4double valueE12 = 0;
832  G4double valueE11 = 0;
833  G4double xs11 = 0;
834  G4double xs12 = 0;
835  G4double xs21 = 0;
836  G4double xs22 = 0;
837
838  // PROTECTION TO ALLOW INTERPOLATION AT MINIMUM AND MAXIMUM Eta/Theta2 values
839
840  if (
841       theta==8.66e-4 ||
842       theta==8.66e-3 ||
843       theta==8.66e-2 ||
844       theta==8.66e-1 ||
845       theta==8.66e+00 ||
846       theta==8.66e+01
847  ) theta=theta-1e-12;
848
849  if (
850       theta==1.e-4 ||
851       theta==1.e-3 ||
852       theta==1.e-2 ||
853       theta==1.e-1 ||
854       theta==1.e+00 ||
855       theta==1.e+01
856  ) theta=theta+1e-12;
857
858  // END PROTECTION
859
860    std::vector<double>::iterator t2 = std::upper_bound(dummyVec2.begin(),dummyVec2.end(), k);
861    std::vector<double>::iterator t1 = t2-1;
862
863    std::vector<double>::iterator e12 = std::upper_bound(aVecMap2[(*t1)].begin(),aVecMap2[(*t1)].end(), theta);
864    std::vector<double>::iterator e11 = e12-1;
865
866    std::vector<double>::iterator e22 = std::upper_bound(aVecMap2[(*t2)].begin(),aVecMap2[(*t2)].end(), theta);
867    std::vector<double>::iterator e21 = e22-1;
868
869    valueT1  =*t1;
870    valueT2  =*t2;
871    valueE21 =*e21;
872    valueE22 =*e22;
873    valueE12 =*e12;
874    valueE11 =*e11;
875
876    xs11 = FL2Data[valueT1][valueE11];
877    xs12 = FL2Data[valueT1][valueE12];
878    xs21 = FL2Data[valueT2][valueE21];
879    xs22 = FL2Data[valueT2][valueE22];
880
881  if (verboseLevel>0) 
882    G4cout
883    << valueT1 << " "
884    << valueT2 << " "
885    << valueE11 << " "
886    << valueE12 << " "
887    << valueE21 << " "
888    << valueE22 << " "
889    << xs11 << " "
890    << xs12 << " "
891    << xs21 << " "
892    << xs22 << " "
893    << G4endl;
894
895  G4double xsProduct = xs11 * xs12 * xs21 * xs22;
896
897  if (xs11==0 || xs12==0 ||xs21==0 ||xs22==0) return (0.);
898
899  if (xsProduct != 0.)
900  {
901    sigma = QuadInterpolator(  valueE11, valueE12,
902                               valueE21, valueE22,
903                               xs11, xs12,
904                               xs21, xs22,
905                               valueT1, valueT2,
906                               k, theta );
907  }
908
909  return sigma;
910}
911
912//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
913
914G4double G4AnalyticalEcpssrLiCrossSection::LinLinInterpolate(G4double e1,
915                                                        G4double e2,
916                                                        G4double e,
917                                                        G4double xs1,
918                                                        G4double xs2)
919{
920  G4double value = xs1 + (xs2 - xs1)*(e - e1)/ (e2 - e1);
921  return value;
922}
923
924//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
925
926G4double G4AnalyticalEcpssrLiCrossSection::LinLogInterpolate(G4double e1,
927                                                        G4double e2,
928                                                        G4double e,
929                                                        G4double xs1,
930                                                        G4double xs2)
931{
932  G4double d1 = std::log(xs1);
933  G4double d2 = std::log(xs2);
934  G4double value = std::exp(d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
935  return value;
936}
937
938//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
939
940G4double G4AnalyticalEcpssrLiCrossSection::LogLogInterpolate(G4double e1,
941                                                        G4double e2,
942                                                        G4double e,
943                                                        G4double xs1,
944                                                        G4double xs2)
945{
946  G4double a = (std::log10(xs2)-std::log10(xs1)) / (std::log10(e2)-std::log10(e1));
947  G4double b = std::log10(xs2) - a*std::log10(e2);
948  G4double sigma = a*std::log10(e) + b;
949  G4double value = (std::pow(10.,sigma));
950  return value;
951}
952
953//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
954
955G4double G4AnalyticalEcpssrLiCrossSection::QuadInterpolator(G4double e11, G4double e12,
956                                                       G4double e21, G4double e22,
957                                                       G4double xs11, G4double xs12,
958                                                       G4double xs21, G4double xs22,
959                                                       G4double t1, G4double t2,
960                                                       G4double t, G4double e)
961{
962// Log-Log
963  G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
964  G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
965  G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
966
967/*
968// Lin-Log
969  G4double interpolatedvalue1 = LinLogInterpolate(e11, e12, e, xs11, xs12);
970  G4double interpolatedvalue2 = LinLogInterpolate(e21, e22, e, xs21, xs22);
971  G4double value = LinLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
972*/
973
974/*
975// Lin-Lin
976  G4double interpolatedvalue1 = LinLinInterpolate(e11, e12, e, xs11, xs12);
977  G4double interpolatedvalue2 = LinLinInterpolate(e21, e22, e, xs21, xs22);
978  G4double value = LinLinInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
979*/
980  return value;
981
982}
983
Note: See TracBrowser for help on using the repository browser.