Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/coherent_elastic/src/G4AntiNuclElastic.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 //
 27 // Geant4 Header : G4AntiNuclElastic
 28 //
 29 // 
 30 
 31 #include "G4AntiNuclElastic.hh"
 32 
 33 #include "G4PhysicalConstants.hh"
 34 #include "G4SystemOfUnits.hh"
 35 #include "G4ParticleTable.hh"
 36 #include "G4ParticleDefinition.hh"
 37 #include "G4IonTable.hh"
 38 #include "Randomize.hh"
 39 #include "G4AntiProton.hh"
 40 #include "G4AntiNeutron.hh"
 41 #include "G4AntiDeuteron.hh"
 42 #include "G4AntiAlpha.hh"
 43 #include "G4AntiTriton.hh"
 44 #include "G4AntiHe3.hh"
 45 #include "G4Proton.hh"
 46 #include "G4Neutron.hh"
 47 #include "G4Deuteron.hh"
 48 #include "G4Alpha.hh"
 49 #include "G4Pow.hh"
 50 #include "G4Exp.hh"
 51 #include "G4Log.hh"
 52 
 53 #include "G4NucleiProperties.hh"
 54 #include "G4CrossSectionDataSetRegistry.hh"
 55 
 56 G4AntiNuclElastic::G4AntiNuclElastic() 
 57   : G4HadronElastic("AntiAElastic")
 58 {
 59   //V.Ivanchenko commented out 
 60   //SetMinEnergy( 0.1*GeV );
 61   //SetMaxEnergy( 10.*TeV );
 62 
 63   theAProton       = G4AntiProton::AntiProton();
 64   theANeutron      = G4AntiNeutron::AntiNeutron();
 65   theADeuteron     = G4AntiDeuteron::AntiDeuteron();
 66   theATriton       = G4AntiTriton::AntiTriton();
 67   theAAlpha        = G4AntiAlpha::AntiAlpha();
 68   theAHe3          = G4AntiHe3::AntiHe3();
 69 
 70   theProton   = G4Proton::Proton();
 71   theNeutron  = G4Neutron::Neutron();
 72   theDeuteron = G4Deuteron::Deuteron();
 73   theAlpha    = G4Alpha::Alpha();
 74 
 75   G4CrossSectionDataSetRegistry* reg = G4CrossSectionDataSetRegistry::Instance();
 76   cs = static_cast<G4ComponentAntiNuclNuclearXS*>(reg->GetComponentCrossSection("AntiAGlauber"));
 77   if(!cs) { cs = new G4ComponentAntiNuclNuclearXS(); }
 78 
 79   fParticle = 0;
 80   fWaveVector = 0.;
 81   fBeta = 0.;
 82   fZommerfeld = 0.;
 83   fAm = 0.;
 84   fTetaCMS = 0.;    
 85   fRa = 0.;
 86   fRef = 0.;
 87   fceff = 0.;
 88   fptot = 0.;
 89   fTmax = 0.;
 90   fThetaLab = 0.;
 91 }
 92 
 93 /////////////////////////////////////////////////////////////////////////
 94 G4AntiNuclElastic::~G4AntiNuclElastic()
 95 {}
 96 
 97 ////////////////////////////////////////////////////////////////////////
 98 // sample momentum transfer in the CMS system 
 99 G4double G4AntiNuclElastic::SampleInvariantT(const G4ParticleDefinition* particle, 
100                G4double Plab,  G4int Z, G4int A)
101 {
102   G4double T;
103   G4double Mproj = particle->GetPDGMass();
104   G4LorentzVector Pproj(0.,0.,Plab,std::sqrt(Plab*Plab+Mproj*Mproj));
105   G4double ctet1 = GetcosTeta1(Plab, A); 
106 
107   G4double energy=Pproj.e()-Mproj;   
108   
109   const G4ParticleDefinition* theParticle = particle;
110 
111   G4ParticleDefinition * theTargetDef = 0;
112 
113   if      (Z == 1 && A == 1) theTargetDef = theProton;
114   else if (Z == 1 && A == 2) theTargetDef = theDeuteron;
115   else if (Z == 1 && A == 3) theTargetDef = G4Triton::Triton();
116   else if (Z == 2 && A == 3) theTargetDef = G4He3::He3();
117   else if (Z == 2 && A == 4) theTargetDef = theAlpha;
118 
119 
120   G4double TargMass =G4NucleiProperties::GetNuclearMass(A,Z); 
121 
122   //transform to CMS
123 
124   G4LorentzVector lv(0.0,0.0,0.0,TargMass);   
125   lv += Pproj;
126   G4double S = lv.mag2()/(GeV*GeV);
127 
128   G4ThreeVector bst = lv.boostVector();
129   Pproj.boost(-bst);
130 
131   G4ThreeVector p1 = Pproj.vect();
132   G4double ptot    = p1.mag();
133 
134   fbst = bst;
135   fptot= ptot;
136   fTmax = 4.0*ptot*ptot;  // In (MeV/c)^2
137 
138   if(Plab < (std::abs(particle->GetBaryonNumber())*100)*MeV)
139   {return fTmax*G4UniformRand();}
140   
141   // Calculation of NN collision properties
142   G4double PlabPerN = Plab/std::abs(theParticle->GetBaryonNumber());
143   G4double NucleonMass = 0.5*( theProton->GetPDGMass() + theNeutron->GetPDGMass() );
144   G4double PrNucleonMass(0.);  // Projectile average nucleon mass
145   if( std::abs(theParticle->GetBaryonNumber()) == 1 ) { PrNucleonMass = theParticle->GetPDGMass(); }
146   else                                                { PrNucleonMass = NucleonMass;               }
147   G4double energyPerN = std::sqrt( sqr(PlabPerN) + sqr(PrNucleonMass));
148   energyPerN -= PrNucleonMass;
149   //---
150 
151   G4double  Z1 = particle->GetPDGCharge();
152   G4double  Z2 = Z;  
153   
154   G4double beta = CalculateParticleBeta(particle, ptot); 
155   G4double n  = CalculateZommerfeld( beta,  Z1,  Z2 );
156   G4double Am = CalculateAm( ptot,  n,  Z2 );
157   fWaveVector = ptot;     //   /hbarc; 
158       
159   G4LorentzVector Fproj(0.,0.,0.,0.);  
160   const G4double mevToBarn = 0.38938e+6;
161   G4double XsCoulomb = mevToBarn*sqr(n/fWaveVector)*pi*(1+ctet1)/(1.+Am)/(1.+2.*Am-ctet1);
162 
163   G4double XsElastHadronic =cs->GetElasticElementCrossSection(particle, energy, Z, (G4double)A);
164   G4double XsTotalHadronic =cs->GetTotalElementCrossSection(particle, energy, Z, (G4double)A);
165 
166   XsElastHadronic/=millibarn; XsTotalHadronic/=millibarn;
167 
168   G4double CoulombProb =  XsCoulomb/(XsCoulomb+XsElastHadronic);
169 
170   if(G4UniformRand() < CoulombProb)
171   {  // Simulation of Coulomb scattering
172 
173     G4double phi = twopi * G4UniformRand();
174     G4double Ksi =  G4UniformRand();
175 
176     G4double par1 = 2.*(1.+Am)/(1.+ctet1);
177 
178     // ////sample ThetaCMS in Coulomb part
179 
180     G4double cosThetaCMS = (par1*ctet1- Ksi*(1.+2.*Am))/(par1-Ksi);
181 
182     G4double PtZ=ptot*cosThetaCMS;
183     Fproj.setPz(PtZ);
184     G4double PtProjCMS = ptot*std::sqrt(1.0 - cosThetaCMS*cosThetaCMS);
185     G4double PtX= PtProjCMS * std::cos(phi);
186     G4double PtY= PtProjCMS * std::sin(phi);
187     Fproj.setPx(PtX);     
188     Fproj.setPy(PtY);
189     Fproj.setE(std::sqrt(PtX*PtX+PtY*PtY+PtZ*PtZ+Mproj*Mproj));    
190     T =  -(Pproj-Fproj).mag2();      
191   } 
192   else 
193   {  
194     // Simulation of strong interaction scattering
195 
196     G4double Qmax = 2.*ptot/197.33;   // in fm^-1
197 
198     G4double Amag      = 1.0;  // A1 in Majorant funct:A1*exp(-q*A2)
199     G4double SlopeMag  = 0.5;  // A2 in Majorant funct:A1*exp(-q*A2)
200     
201     G4double sig_pbarp = cs->GetAntiHadronNucleonTotCrSc(theAProton,energyPerN);  //mb
202 
203     fRa = 1.113*G4Pow::GetInstance()->Z13(A) -                     
204           0.227/G4Pow::GetInstance()->Z13(A);                      
205     if(A == 3) fRa=1.81;
206     if(A == 4) fRa=1.37;
207  
208     if((A>=12.) && (A<27) ) fRa=fRa*0.85;
209     if((A>=27.) && (A<48) ) fRa=fRa*0.90;
210     if((A>=48.) && (A<65) ) fRa=fRa*0.95;
211 
212     G4double Ref2 = XsTotalHadronic/10./2./pi;  // in fm^2
213     G4double ceff2 = 0.0;
214     G4double rho = 0.0;
215 
216     if  ((theParticle == theAProton) || (theParticle == theANeutron))  
217     {
218       if(theTargetDef == theProton)
219       { 
220         // Determination of the real part of Pbar+N amplitude
221         if(Plab < 610.) 
222         { rho = 1.3347-10.342*Plab/1000.+22.277*Plab/1000.*Plab/1000.-
223                 13.634*Plab/1000.*Plab/1000.*Plab/1000. ;}
224         if((Plab < 5500.)&&(Plab >= 610.) )
225         { rho = 0.22; }
226         if((Plab >= 5500.)&&(Plab < 12300.) )
227         { rho = -0.32; }
228         if( Plab >= 12300.)
229         { rho = 0.135-2.26/(std::sqrt(S)) ;}
230         Ref2  = 0.35 + 0.9/std::sqrt(std::sqrt(S-4.*0.88))+0.04*G4Log(S) ;
231         ceff2 = 0.375 - 2./S + 0.44/(sqr(S-4.)+1.5) ;
232         Ref2  =Ref2*Ref2;
233         ceff2 = ceff2*ceff2;
234       }
235 
236       if( (Z==1)&&(A==2) )
237       {
238         Ref2 = fRa*fRa - 0.28 + 0.019 * sig_pbarp + 2.06e-6 * sig_pbarp*sig_pbarp;  
239         ceff2 = 0.297 + 7.853e-04*sig_pbarp + 0.2899*G4Exp(-0.03*sig_pbarp);
240       }
241       if( (Z==1)&&(A==3) )
242       { 
243         Ref2 = fRa*fRa - 1.36 + 0.025 * sig_pbarp - 3.69e-7 * sig_pbarp*sig_pbarp;
244         ceff2 = 0.149 + 7.091e-04*sig_pbarp + 0.3743*G4Exp(-0.03*sig_pbarp);
245       }
246       if( (Z==2)&&(A==3) )
247       { 
248         Ref2 = fRa*fRa - 1.36 + 0.025 * sig_pbarp - 3.69e-7 * sig_pbarp*sig_pbarp;
249         ceff2 = 0.149 + 7.091e-04*sig_pbarp + 0.3743*G4Exp(-0.03*sig_pbarp);
250       }  
251       if( (Z==2)&&(A==4) )
252       { 
253         Ref2 = fRa*fRa -0.46 +0.03*sig_pbarp - 2.98e-6*sig_pbarp*sig_pbarp;
254         ceff2= 0.078 + 6.657e-4*sig_pbarp + 0.3359*G4Exp(-0.03*sig_pbarp);
255       }
256       if(Z>2) 
257       { 
258         Ref2 = fRa*fRa +2.48*0.01*sig_pbarp*fRa - 2.23e-6*sig_pbarp*sig_pbarp*fRa*fRa; 
259         ceff2 = 0.16+3.3e-4*sig_pbarp+0.35*G4Exp(-0.03*sig_pbarp);
260       }
261     }  // End of if ((theParticle == theAProton) || (theParticle == theANeutron))
262 
263     if (theParticle == theADeuteron)
264     {
265       if(theTargetDef == theProton)
266       { 
267         ceff2 = 0.297 + 7.853e-04*sig_pbarp + 0.2899*G4Exp(-0.03*sig_pbarp);
268       }
269       if(theTargetDef == theDeuteron)
270       { 
271         ceff2 = 0.65 + 3.0e-4*sig_pbarp + 0.55 * G4Exp(-0.03*sig_pbarp);
272       }
273       if( (theTargetDef == G4Triton::Triton()) || (theTargetDef == G4He3::He3() ) )
274       {
275         ceff2 = 0.57 + 2.5e-4*sig_pbarp + 0.65 * G4Exp(-0.02*sig_pbarp);
276       }
277       if(theTargetDef == theAlpha)
278       {
279         ceff2 = 0.40 + 3.5e-4 *sig_pbarp + 0.45 * G4Exp(-0.02*sig_pbarp);
280       }
281       if(Z>2)
282       {
283         ceff2 = 0.38 + 2.0e-4 *sig_pbarp + 0.5 * G4Exp(-0.03*sig_pbarp);
284       }
285     }
286 
287     if( (theParticle ==theAHe3) || (theParticle ==theATriton) )
288     {
289       if(theTargetDef == theProton)
290       {
291         ceff2 = 0.149 + 7.091e-04*sig_pbarp + 0.3743*G4Exp(-0.03*sig_pbarp);         
292       }
293       if(theTargetDef == theDeuteron)
294       {
295         ceff2 = 0.57 + 2.5e-4*sig_pbarp + 0.65 * G4Exp(-0.02*sig_pbarp);
296       }
297       if( (theTargetDef == G4Triton::Triton()) || (theTargetDef == G4He3::He3() ) )
298       {
299         ceff2 = 0.39 + 2.7e-4*sig_pbarp + 0.7 * G4Exp(-0.02*sig_pbarp);
300       }
301       if(theTargetDef == theAlpha)
302       {
303         ceff2 = 0.24 + 3.5e-4*sig_pbarp + 0.75 * G4Exp(-0.03*sig_pbarp);
304       }
305       if(Z>2)
306       {
307         ceff2 = 0.26 + 2.2e-4*sig_pbarp + 0.33*G4Exp(-0.03*sig_pbarp);
308       }  
309     }
310 
311     if ( (theParticle == theAAlpha) || (ceff2 == 0.0) )
312     {
313       if(theTargetDef == theProton)
314       {
315         ceff2= 0.078 + 6.657e-4*sig_pbarp + 0.3359*G4Exp(-0.03*sig_pbarp);   
316       }
317       if(theTargetDef == theDeuteron)  
318       {
319         ceff2 = 0.40 + 3.5e-4 *sig_pbarp + 0.45 * G4Exp(-0.02*sig_pbarp);
320       }   
321       if( (theTargetDef == G4Triton::Triton()) || (theTargetDef == G4He3::He3() ) )
322       {
323         ceff2 = 0.24 + 3.5e-4*sig_pbarp + 0.75 * G4Exp(-0.03*sig_pbarp);   
324       }
325       if(theTargetDef == theAlpha)   
326       {
327         ceff2 = 0.17 + 3.5e-4*sig_pbarp + 0.45 * G4Exp(-0.03*sig_pbarp);
328       }
329       if(Z>2)
330       {
331         ceff2 = 0.22 + 2.0e-4*sig_pbarp + 0.2 * G4Exp(-0.03*sig_pbarp);
332       }
333     }
334 
335     fRef=std::sqrt(Ref2);
336     fceff = std::sqrt(ceff2);     
337 
338     G4double Q = 0.0 ;
339     G4double BracFunct;
340 
341     const G4int maxNumberOfLoops = 10000;
342     G4int loopCounter = 0;
343     do 
344     {
345       Q = -G4Log(1.-(1.- G4Exp(-SlopeMag * Qmax))* G4UniformRand() )/SlopeMag;
346       G4double x = fRef * Q;
347       BracFunct = ( ( sqr(BesselOneByArg(x))+sqr(rho/2. * BesselJzero(x)) )
348         *   sqr(DampFactor(pi*fceff*Q))) /(Amag*G4Exp(-SlopeMag*Q));
349       BracFunct = BracFunct * Q;
350     } 
351     while ( (G4UniformRand()>BracFunct) && 
352             ++loopCounter < maxNumberOfLoops );  /* Loop checking, 10.08.2015, A.Ribon */
353     if ( loopCounter >= maxNumberOfLoops ) {
354       fTetaCMS = 0.0;
355       return 0.0;
356     }
357 
358     T= sqr(Q); 
359     T*=3.893913e+4;  // fm^(-2) -> MeV^2
360 
361   }  // End of simulation of strong interaction scattering
362 
363   return T;
364 }
365 
366 /////////////////////////////////////////////////////////////////////
367 //  Sample of Theta in CMS
368  G4double G4AntiNuclElastic::SampleThetaCMS(const G4ParticleDefinition* p, G4double plab,
369                                                                          G4int Z, G4int A)
370 { 
371   G4double T;
372   T =  SampleInvariantT( p, plab,  Z,  A);
373 
374    // NaN finder
375   if(!(T < 0.0 || T >= 0.0))
376   {
377     if (verboseLevel > 0)
378     {
379       G4cout << "G4DiffuseElastic:WARNING: A = " << A
380              << " mom(GeV)= " << plab/GeV
381              << " S-wave will be sampled"
382              << G4endl;
383     }
384     T = G4UniformRand()*fTmax;
385  
386   }
387 
388   if(fptot > 0.)
389   {
390    G4double cosTet=1.0-T/(2.*fptot*fptot);
391    if(cosTet >  1.0 ) cosTet= 1.;
392    if(cosTet < -1.0 ) cosTet=-1.;
393    fTetaCMS=std::acos(cosTet); 
394    return fTetaCMS;
395   } else
396   {
397    return 2.*G4UniformRand()-1.;
398   }
399 }  
400 
401 
402 /////////////////////////////////////////////////////////////////////
403 //  Sample of Theta in Lab System
404  G4double G4AntiNuclElastic::SampleThetaLab(const G4ParticleDefinition* p, G4double plab,
405                                                                          G4int Z, G4int A)
406 { 
407   G4double T; 
408   T = SampleInvariantT( p, plab,  Z,  A);
409 
410  // NaN finder
411   if(!(T < 0.0 || T >= 0.0))
412   {
413     if (verboseLevel > 0)               
414     {
415       G4cout << "G4DiffuseElastic:WARNING: A = " << A
416              << " mom(GeV)= " << plab/GeV
417              << " S-wave will be sampled"
418              << G4endl;
419     }
420     T = G4UniformRand()*fTmax;
421   }
422 
423   G4double phi  = G4UniformRand()*twopi;
424 
425   G4double cost(1.);
426   if(fTmax > 0.) {cost = 1. - 2.0*T/fTmax;}
427 
428   G4double sint;
429   if( cost >= 1.0 )
430   {
431     cost = 1.0;
432     sint = 0.0;
433   }
434   else if( cost <= -1.0)
435   {
436     cost = -1.0;
437     sint =  0.0;
438   }
439   else
440   {
441     sint = std::sqrt((1.0-cost)*(1.0+cost));
442   }
443 
444   G4double m1 = p->GetPDGMass();
445   G4ThreeVector v(sint*std::cos(phi),sint*std::sin(phi),cost);
446   v *= fptot;
447   G4LorentzVector nlv(v.x(),v.y(),v.z(),std::sqrt(fptot*fptot + m1*m1));
448    
449   nlv.boost(fbst);
450    
451   G4ThreeVector np = nlv.vect();
452   G4double theta = np.theta();
453   fThetaLab = theta; 
454 
455   return theta;
456 }
457 
458 ////////////////////////////////////////////////////////////////////
459 //   Calculation of Damp factor
460  G4double G4AntiNuclElastic::DampFactor(G4double x)
461 {
462    G4double df;
463    G4double  f3 = 6.; // first factorials
464 
465  if( std::fabs(x) < 0.01 )
466   { 
467       df=1./(1.+x*x/f3); 
468  }
469   else
470   {
471     df = x/std::sinh(x); 
472   }
473   return df;
474 }
475 
476 
477 /////////////////////////////////////////////////////////////////////////////////
478 //  Calculation of particle velocity Beta
479 
480  G4double G4AntiNuclElastic::CalculateParticleBeta( const G4ParticleDefinition* particle, 
481                                   G4double momentum    )
482 {
483   G4double mass = particle->GetPDGMass();
484   G4double a    = momentum/mass;
485   fBeta         = a/std::sqrt(1+a*a);
486 
487   return fBeta; 
488 }
489 
490 
491 ///////////////////////////////////////////////////////////////////////////////////
492 //   Calculation of parameter Zommerfeld
493 
494  G4double G4AntiNuclElastic::CalculateZommerfeld( G4double beta, G4double Z1, G4double Z2 )
495 {
496   fZommerfeld = fine_structure_const*Z1*Z2/beta;
497 
498   return fZommerfeld; 
499 }
500 
501 ////////////////////////////////////////////////////////////////////////////////////
502 //  
503 G4double G4AntiNuclElastic::CalculateAm( G4double momentum, G4double n, G4double Z)
504 {
505   G4double k   = momentum/hbarc;
506   G4double ch  = 1.13 + 3.76*n*n;
507   G4double zn  = 1.77*k/G4Pow::GetInstance()->A13(Z)*Bohr_radius; 
508   G4double zn2 = zn*zn;
509   fAm          = ch/zn2;
510 
511   return fAm;
512 }
513 
514 /////////////////////////////////////////////////////////////
515 //
516 // Bessel J0 function based on rational approximation from
517 // J.F. Hart, Computer Approximations, New York, Willey 1968, p. 141
518    
519 G4double G4AntiNuclElastic::BesselJzero(G4double value)
520 {  
521   G4double modvalue, value2, fact1, fact2, arg, shift, bessel;
522 
523   modvalue = std::fabs(value);
524 
525   if ( value < 8.0 && value > -8.0 )
526   {
527     value2 = value*value;
528                  
529     fact1  = 57568490574.0 + value2*(-13362590354.0
530                            + value2*( 651619640.7  
531                            + value2*(-11214424.18
532                            + value2*( 77392.33017
533                            + value2*(-184.9052456   ) ) ) ) );
534                               
535     fact2  = 57568490411.0 + value2*( 1029532985.0
536                            + value2*( 9494680.718
537                            + value2*(59272.64853
538                            + value2*(267.8532712
539                            + value2*1.0               ) ) ) );
540  
541     bessel = fact1/fact2;
542   }
543   else
544   {
545     arg    = 8.0/modvalue;
546 
547     value2 = arg*arg;
548  
549     shift  = modvalue-0.785398164;
550 
551     fact1  = 1.0 + value2*(-0.1098628627e-2
552                  + value2*(0.2734510407e-4
553                  + value2*(-0.2073370639e-5
554                  + value2*0.2093887211e-6    ) ) );
555   fact2  = -0.1562499995e-1 + value2*(0.1430488765e-3
556                               + value2*(-0.6911147651e-5
557                               + value2*(0.7621095161e-6
558                               - value2*0.934945152e-7    ) ) );
559 
560     bessel = std::sqrt(0.636619772/modvalue)*(std::cos(shift)*fact1 - arg*std::sin(shift)*fact2);
561   }
562   return bessel;
563 }
564 
565 
566 //////////////////////////////////////////////////////////////////////////////
567 // Bessel J1 function based on rational approximation from
568 // J.F. Hart, Computer Approximations, New York, Willey 1968, p. 141
569     
570  G4double G4AntiNuclElastic::BesselJone(G4double value)
571 {
572   G4double modvalue, value2, fact1, fact2, arg, shift, bessel;
573                           
574   modvalue = std::fabs(value);
575                  
576   if ( modvalue < 8.0 )
577   {
578     value2 = value*value;
579     fact1  = value*(72362614232.0 + value2*(-7895059235.0
580                                   + value2*( 242396853.1
581                                   + value2*(-2972611.439
582                                   + value2*( 15704.48260
583                                   + value2*(-30.16036606  ) ) ) ) ) );
584     
585     fact2  = 144725228442.0 + value2*(2300535178.0
586                             + value2*(18583304.74
587                             + value2*(99447.43394
588                             + value2*(376.9991397
589                             + value2*1.0             ) ) ) );
590     bessel = fact1/fact2;
591   }
592   else
593   {
594     arg    = 8.0/modvalue;  
595   value2 = arg*arg;
596 
597     shift  = modvalue - 2.356194491;
598 
599     fact1  = 1.0 + value2*( 0.183105e-2
600                  + value2*(-0.3516396496e-4
601                  + value2*(0.2457520174e-5
602                  + value2*(-0.240337019e-6          ) ) ) );
603 
604     fact2 = 0.04687499995 + value2*(-0.2002690873e-3
605                           + value2*( 0.8449199096e-5
606                           + value2*(-0.88228987e-6
607                           + value2*0.105787412e-6       ) ) );
608 
609     bessel = std::sqrt( 0.636619772/modvalue)*(std::cos(shift)*fact1 - arg*std::sin(shift)*fact2);
610     if (value < 0.0) bessel = -bessel;
611   }
612   return bessel;
613 }
614 
615 ////////////////////////////////////////////////////////////////////////////////
616 // return J1(x)/x with special case for small x 
617 G4double G4AntiNuclElastic::BesselOneByArg(G4double x)
618 {
619   G4double x2, result;
620  
621   if( std::fabs(x) < 0.01 )
622   {
623    x  *= 0.5;
624    x2 = x*x;
625    result = (2.- x2 + x2*x2/6.)/4.;
626   }
627   else
628   {
629     result = BesselJone(x)/x;
630   }
631   return result;
632 } 
633 
634 /////////////////////////////////////////////////////////////////////////////////
635 // return  angle from which Coulomb scattering is calculated 
636 G4double G4AntiNuclElastic::GetcosTeta1(G4double plab, G4int A)
637 {
638 
639 // G4double p0 =G4LossTableManager::Instance()->FactorForAngleLimit()*CLHEP::hbarc/CLHEP::fermi;
640   G4double p0 = 1.*hbarc/fermi;
641 //G4double cteta1 = 1.0 - p0*p0/2.0 * pow(A,2./3.)/(plab*plab);
642   G4double cteta1 = 1.0 - p0*p0/2.0 * G4Pow::GetInstance()->Z23(A)/(plab*plab);
643 //////////////////
644   if(cteta1 < -1.) cteta1 = -1.0;
645   return cteta1;
646 }
647 
648 
649 
650 
651 
652 
653 
654