source: trunk/source/processes/hadronic/models/coherent_elastic/include/G4NuclNuclDiffuseElastic.hh @ 1340

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 42.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//
27// $Id: G4NuclNuclDiffuseElastic.hh,v 1.13 2010/09/28 16:28:58 gcosmo Exp $
28// GEANT4 tag $Name: geant4-09-03-ref-09 $
29//
30//
31// G4 Model: optical elastic scattering with 4-momentum balance
32//
33// Class Description
34// Final state production model for nucleus-nucleus elastic scattering;
35// Coulomb amplitude is not considered as correction
36// (as in G4DiffuseElastic)
37// Class Description - End
38//
39//
40// 17.03.09 V. Grichine implementation for Coulomb elastic scattering
41
42
43#ifndef G4NuclNuclDiffuseElastic_h
44#define G4NuclNuclDiffuseElastic_h 1
45 
46#include "globals.hh"
47#include <complex>
48#include "G4Integrator.hh"
49
50#include "G4HadronicInteraction.hh"
51#include "G4HadProjectile.hh"
52#include "G4Nucleus.hh"
53
54using namespace std;
55
56class G4ParticleDefinition;
57class G4PhysicsTable;
58class G4PhysicsLogVector;
59
60class G4NuclNuclDiffuseElastic : public G4HadronicInteraction
61{
62public:
63
64  G4NuclNuclDiffuseElastic();
65
66  G4NuclNuclDiffuseElastic(const G4ParticleDefinition* aParticle);
67
68
69
70
71
72  virtual ~G4NuclNuclDiffuseElastic();
73
74  void Initialise();
75
76  void InitialiseOnFly(G4double Z, G4double A);
77
78  void BuildAngleTable();
79
80 
81  G4HadFinalState * ApplyYourself(const G4HadProjectile & aTrack, 
82                                  G4Nucleus & targetNucleus);
83
84
85  void SetPlabLowLimit(G4double value);
86
87  void SetHEModelLowLimit(G4double value);
88
89  void SetQModelLowLimit(G4double value);
90
91  void SetLowestEnergyLimit(G4double value);
92
93  void SetRecoilKinEnergyLimit(G4double value);
94
95  G4double SampleT(const G4ParticleDefinition* aParticle, 
96                         G4double p, G4double A);
97
98  G4double SampleTableT(const G4ParticleDefinition* aParticle, 
99                         G4double p, G4double Z, G4double A);
100
101  G4double SampleThetaCMS(const G4ParticleDefinition* aParticle, G4double p, G4double A);
102
103  G4double SampleTableThetaCMS(const G4ParticleDefinition* aParticle, G4double p, 
104                                     G4double Z, G4double A);
105
106  G4double GetScatteringAngle(G4int iMomentum, G4int iAngle, G4double position);
107
108  G4double SampleThetaLab(const G4HadProjectile* aParticle, 
109                                G4double tmass, G4double A);
110
111  G4double GetDiffuseElasticXsc( const G4ParticleDefinition* particle, 
112                                 G4double theta, 
113                                 G4double momentum, 
114                                 G4double A         );
115
116  G4double GetInvElasticXsc( const G4ParticleDefinition* particle, 
117                                 G4double theta, 
118                                 G4double momentum, 
119                                 G4double A, G4double Z );
120
121  G4double GetDiffuseElasticSumXsc( const G4ParticleDefinition* particle, 
122                                 G4double theta, 
123                                 G4double momentum, 
124                                 G4double A, G4double Z );
125
126  G4double GetInvElasticSumXsc( const G4ParticleDefinition* particle, 
127                                 G4double tMand, 
128                                 G4double momentum, 
129                                 G4double A, G4double Z );
130
131  G4double IntegralElasticProb( const G4ParticleDefinition* particle, 
132                                 G4double theta, 
133                                 G4double momentum, 
134                                 G4double A            );
135 
136
137  G4double GetCoulombElasticXsc( const G4ParticleDefinition* particle, 
138                                 G4double theta, 
139                                 G4double momentum, 
140                                 G4double Z         );
141
142  G4double GetInvCoulombElasticXsc( const G4ParticleDefinition* particle, 
143                                 G4double tMand, 
144                                 G4double momentum, 
145                                 G4double A, G4double Z         );
146
147  G4double GetCoulombTotalXsc( const G4ParticleDefinition* particle, 
148                                 G4double momentum, G4double Z       );
149
150  G4double GetCoulombIntegralXsc( const G4ParticleDefinition* particle, 
151                                 G4double momentum, G4double Z, 
152                                 G4double theta1, G4double theta2         );
153
154
155  G4double CalculateParticleBeta( const G4ParticleDefinition* particle, 
156                                        G4double momentum    );
157
158  G4double CalculateZommerfeld( G4double beta, G4double Z1, G4double Z2 );
159
160  G4double CalculateAm( G4double momentum, G4double n, G4double Z);
161
162  G4double CalculateNuclearRad( G4double A);
163
164  G4double ThetaCMStoThetaLab(const G4DynamicParticle* aParticle, 
165                                G4double tmass, G4double thetaCMS);
166
167  G4double ThetaLabToThetaCMS(const G4DynamicParticle* aParticle, 
168                                G4double tmass, G4double thetaLab);
169
170  void TestAngleTable(const G4ParticleDefinition* theParticle, G4double partMom,
171                      G4double Z, G4double A);
172
173
174
175  G4double BesselJzero(G4double z);
176  G4double BesselJone(G4double z);
177  G4double DampFactor(G4double z);
178  G4double BesselOneByArg(G4double z);
179
180  G4double GetDiffElasticProb(G4double theta);
181  G4double GetDiffElasticSumProb(G4double theta);
182  G4double GetDiffElasticSumProbA(G4double alpha);
183  G4double GetIntegrandFunction(G4double theta);
184
185  G4double GetNuclearRadius(){return fNuclearRadius;};
186
187
188  // Technical math functions for strong Coulomb contribution
189
190  G4complex GammaLogarithm(G4complex xx);
191  G4complex GammaLogB2n(G4complex xx);
192
193  G4double  GetErf(G4double x);
194
195  G4complex GetErfcComp(G4complex z, G4int nMax);
196  G4complex GetErfcSer(G4complex z, G4int nMax);
197  G4complex GetErfcInt(G4complex z); // , G4int nMax);
198
199  G4complex GetErfComp(G4complex z, G4int nMax);  // AandS algorithm != Ser, Int
200  G4complex GetErfSer(G4complex z, G4int nMax);
201
202  G4double GetExpCos(G4double x);
203  G4double GetExpSin(G4double x);
204  G4complex GetErfInt(G4complex z); // , G4int nMax);
205
206  G4double GetLegendrePol(G4int n, G4double x);
207
208  G4complex TestErfcComp(G4complex z, G4int nMax);
209  G4complex TestErfcSer(G4complex z, G4int nMax);
210  G4complex TestErfcInt(G4complex z); // , G4int nMax);
211
212  G4complex CoulombAmplitude(G4double theta);
213  void CalculateCoulombPhaseZero();
214  G4double CalculateCoulombPhase(G4int n);
215  void CalculateRutherfordAnglePar();
216
217  G4double ProfileNear(G4double theta);
218  G4double ProfileFar(G4double theta);
219
220  G4complex PhaseNear(G4double theta);
221  G4complex PhaseFar(G4double theta);
222
223  G4complex GammaLess(G4double theta);
224  G4complex GammaMore(G4double theta);
225
226  G4complex AmplitudeNear(G4double theta);
227  G4complex AmplitudeFar(G4double theta);
228  G4complex Amplitude(G4double theta);
229  G4double  AmplitudeMod2(G4double theta);
230
231  G4complex AmplitudeGla(G4double theta);
232  G4double  AmplitudeGlaMod2(G4double theta);
233
234  G4complex AmplitudeGG(G4double theta);
235  G4double  AmplitudeGGMod2(G4double theta);
236
237  void      InitParameters(const G4ParticleDefinition* theParticle, 
238                              G4double partMom, G4double Z, G4double A); 
239
240  void      InitParametersGla(const G4DynamicParticle* aParticle, 
241                              G4double partMom, G4double Z, G4double A);
242
243  G4double GetHadronNucleonXscNS( G4ParticleDefinition* pParticle, 
244                                                 G4double pTkin, 
245                                  G4ParticleDefinition* tParticle);
246
247  G4double CalcMandelstamS( const G4double mp , 
248                                                       const G4double mt , 
249                                                    const G4double Plab );
250
251  G4double GetProfileLambda(){return fProfileLambda;};
252
253  void SetProfileLambda(G4double pl) {fProfileLambda = pl;};
254  void SetProfileDelta(G4double pd) {fProfileDelta = pd;};
255  void SetProfileAlpha(G4double pa){fProfileAlpha = pa;};
256  void SetCofLambda(G4double pa){fCofLambda = pa;};
257  void SetCofAlpha(G4double pa){fCofAlpha = pa;};
258  void SetCofDelta(G4double pa){fCofDelta = pa;};
259  void SetCofPhase(G4double pa){fCofPhase = pa;};
260  void SetCofFar(G4double pa){fCofFar = pa;};
261  void SetEtaRatio(G4double pa){fEtaRatio = pa;};
262  void SetMaxL(G4int l){fMaxL = l;};
263
264private:
265
266
267  G4ParticleDefinition* theProton;
268  G4ParticleDefinition* theNeutron;
269  G4ParticleDefinition* theDeuteron;
270  G4ParticleDefinition* theAlpha;
271
272  const G4ParticleDefinition* thePionPlus;
273  const G4ParticleDefinition* thePionMinus;
274
275  G4double lowEnergyRecoilLimit; 
276  G4double lowEnergyLimitHE; 
277  G4double lowEnergyLimitQ; 
278  G4double lowestEnergyLimit; 
279  G4double plabLowLimit;
280
281  G4int fEnergyBin;
282  G4int fAngleBin;
283
284  G4PhysicsLogVector*           fEnergyVector;
285  G4PhysicsTable*               fAngleTable;
286  std::vector<G4PhysicsTable*>  fAngleBank;
287
288  std::vector<G4double> fElementNumberVector;
289  std::vector<G4String> fElementNameVector;
290
291  const G4ParticleDefinition* fParticle;
292  G4double fWaveVector;
293  G4double fAtomicWeight;
294  G4double fAtomicNumber;
295
296  G4double fNuclearRadius1;
297  G4double fNuclearRadius2;
298  G4double fNuclearRadius;
299  G4double fNuclearRadiusSquare;
300
301  G4double fBeta;
302  G4double fZommerfeld;
303  G4double fAm;
304  G4bool   fAddCoulomb;
305
306  G4double fCoulombPhase0;
307  G4double fHalfRutThetaTg;
308  G4double fHalfRutThetaTg2;
309  G4double fRutherfordTheta;
310
311  G4double fProfileLambda;
312  G4double fProfileDelta;
313  G4double fProfileAlpha;
314
315  G4double fCofLambda;
316  G4double fCofAlpha;
317  G4double fCofDelta;
318  G4double fCofPhase;
319  G4double fCofFar;
320
321  G4int    fMaxL;
322  G4double fSumSigma;
323  G4double fEtaRatio;
324
325  G4double fReZ;
326
327};
328
329
330inline void G4NuclNuclDiffuseElastic::SetRecoilKinEnergyLimit(G4double value)
331{
332  lowEnergyRecoilLimit = value;
333}
334
335inline void G4NuclNuclDiffuseElastic::SetPlabLowLimit(G4double value)
336{
337  plabLowLimit = value;
338}
339
340inline void G4NuclNuclDiffuseElastic::SetHEModelLowLimit(G4double value)
341{
342  lowEnergyLimitHE = value;
343}
344
345inline void G4NuclNuclDiffuseElastic::SetQModelLowLimit(G4double value)
346{
347  lowEnergyLimitQ = value;
348}
349
350inline void G4NuclNuclDiffuseElastic::SetLowestEnergyLimit(G4double value)
351{
352  lowestEnergyLimit = value;
353}
354
355
356/////////////////////////////////////////////////////////////
357//
358// Bessel J0 function based on rational approximation from
359// J.F. Hart, Computer Approximations, New York, Willey 1968, p. 141
360
361inline G4double G4NuclNuclDiffuseElastic::BesselJzero(G4double value)
362{
363  G4double modvalue, value2, fact1, fact2, arg, shift, bessel;
364
365  modvalue = fabs(value);
366
367  if ( value < 8.0 && value > -8.0 )
368  {
369    value2 = value*value;
370
371    fact1  = 57568490574.0 + value2*(-13362590354.0 
372                           + value2*( 651619640.7 
373                           + value2*(-11214424.18 
374                           + value2*( 77392.33017 
375                           + value2*(-184.9052456   ) ) ) ) );
376
377    fact2  = 57568490411.0 + value2*( 1029532985.0 
378                           + value2*( 9494680.718
379                           + value2*(59272.64853
380                           + value2*(267.8532712 
381                           + value2*1.0               ) ) ) );
382
383    bessel = fact1/fact2;
384  } 
385  else 
386  {
387    arg    = 8.0/modvalue;
388
389    value2 = arg*arg;
390
391    shift  = modvalue-0.785398164;
392
393    fact1  = 1.0 + value2*(-0.1098628627e-2 
394                 + value2*(0.2734510407e-4
395                 + value2*(-0.2073370639e-5 
396                 + value2*0.2093887211e-6    ) ) );
397
398    fact2  = -0.1562499995e-1 + value2*(0.1430488765e-3
399                              + value2*(-0.6911147651e-5 
400                              + value2*(0.7621095161e-6
401                              - value2*0.934945152e-7    ) ) );
402
403    bessel = sqrt(0.636619772/modvalue)*(cos(shift)*fact1 - arg*sin(shift)*fact2 );
404  }
405  return bessel;
406}
407
408/////////////////////////////////////////////////////////////
409//
410// Bessel J1 function based on rational approximation from
411// J.F. Hart, Computer Approximations, New York, Willey 1968, p. 141
412
413inline G4double G4NuclNuclDiffuseElastic::BesselJone(G4double value)
414{
415  G4double modvalue, value2, fact1, fact2, arg, shift, bessel;
416
417  modvalue = fabs(value);
418
419  if ( modvalue < 8.0 ) 
420  {
421    value2 = value*value;
422
423    fact1  = value*(72362614232.0 + value2*(-7895059235.0 
424                                  + value2*( 242396853.1
425                                  + value2*(-2972611.439 
426                                  + value2*( 15704.48260 
427                                  + value2*(-30.16036606  ) ) ) ) ) );
428
429    fact2  = 144725228442.0 + value2*(2300535178.0 
430                            + value2*(18583304.74
431                            + value2*(99447.43394 
432                            + value2*(376.9991397 
433                            + value2*1.0             ) ) ) );
434    bessel = fact1/fact2;
435  } 
436  else 
437  {
438    arg    = 8.0/modvalue;
439
440    value2 = arg*arg;
441
442    shift  = modvalue - 2.356194491;
443
444    fact1  = 1.0 + value2*( 0.183105e-2 
445                 + value2*(-0.3516396496e-4
446                 + value2*(0.2457520174e-5 
447                 + value2*(-0.240337019e-6          ) ) ) );
448
449    fact2 = 0.04687499995 + value2*(-0.2002690873e-3
450                          + value2*( 0.8449199096e-5
451                          + value2*(-0.88228987e-6
452                          + value2*0.105787412e-6       ) ) );
453
454    bessel = sqrt( 0.636619772/modvalue)*(cos(shift)*fact1 - arg*sin(shift)*fact2);
455
456    if (value < 0.0) bessel = -bessel;
457  }
458  return bessel;
459}
460
461////////////////////////////////////////////////////////////////////
462//
463// damp factor in diffraction x/sh(x), x was already *pi
464
465inline G4double G4NuclNuclDiffuseElastic::DampFactor(G4double x)
466{
467  G4double df;
468  G4double f2 = 2., f3 = 6., f4 = 24.; // first factorials
469
470  // x *= pi;
471
472  if( std::fabs(x) < 0.01 )
473  { 
474    df = 1./(1. + x/f2 + x*x/f3 + x*x*x/f4);
475  }
476  else
477  {
478    df = x/std::sinh(x); 
479  }
480  return df;
481}
482
483
484////////////////////////////////////////////////////////////////////
485//
486// return J1(x)/x with special case for small x
487
488inline G4double G4NuclNuclDiffuseElastic::BesselOneByArg(G4double x)
489{
490  G4double x2, result;
491 
492  if( std::fabs(x) < 0.01 )
493  { 
494   x     *= 0.5;
495   x2     = x*x;
496   result = 2. - x2 + x2*x2/6.;
497  }
498  else
499  {
500    result = BesselJone(x)/x; 
501  }
502  return result;
503}
504
505////////////////////////////////////////////////////////////////////
506//
507// return particle beta
508
509inline  G4double G4NuclNuclDiffuseElastic::CalculateParticleBeta( const G4ParticleDefinition* particle, 
510                                        G4double momentum    )
511{
512  G4double mass = particle->GetPDGMass();
513  G4double a    = momentum/mass;
514  fBeta         = a/std::sqrt(1+a*a);
515
516  return fBeta; 
517}
518
519////////////////////////////////////////////////////////////////////
520//
521// return Zommerfeld parameter for Coulomb scattering
522
523inline  G4double G4NuclNuclDiffuseElastic::CalculateZommerfeld( G4double beta, G4double Z1, G4double Z2 )
524{
525  fZommerfeld = fine_structure_const*Z1*Z2/beta;
526
527  return fZommerfeld; 
528}
529
530////////////////////////////////////////////////////////////////////
531//
532// return Wentzel correction for Coulomb scattering
533
534inline  G4double G4NuclNuclDiffuseElastic::CalculateAm( G4double momentum, G4double n, G4double Z)
535{
536  G4double k   = momentum/hbarc;
537  G4double ch  = 1.13 + 3.76*n*n;
538  G4double zn  = 1.77*k*std::pow(Z,-1./3.)*Bohr_radius;
539  G4double zn2 = zn*zn;
540  fAm          = ch/zn2;
541
542  return fAm;
543}
544
545////////////////////////////////////////////////////////////////////
546//
547// calculate nuclear radius for different atomic weights using different approximations
548
549inline  G4double G4NuclNuclDiffuseElastic::CalculateNuclearRad( G4double A)
550{
551  G4double r0, radius;
552
553  if( A < 50. )
554  {
555    if( A > 10. ) r0  = 1.16*( 1 - std::pow(A, -2./3.) )*fermi;   // 1.08*fermi;
556    else          r0  = 1.1*fermi;
557
558    radius = r0*std::pow(A, 1./3.);
559  }
560  else
561  {
562    r0 = 1.7*fermi;   // 1.7*fermi;
563
564    radius = r0*std::pow(A, 0.27); // 0.27);
565  }
566  return radius;
567}
568
569////////////////////////////////////////////////////////////////////
570//
571// return Coulomb scattering differential xsc with Wentzel correction 
572
573inline  G4double G4NuclNuclDiffuseElastic::GetCoulombElasticXsc( const G4ParticleDefinition* particle, 
574                                 G4double theta, 
575                                 G4double momentum, 
576                                 G4double Z         )
577{
578  G4double sinHalfTheta  = std::sin(0.5*theta);
579  G4double sinHalfTheta2 = sinHalfTheta*sinHalfTheta;
580  G4double beta          = CalculateParticleBeta( particle, momentum);
581  G4double z             = particle->GetPDGCharge();
582  G4double n             = CalculateZommerfeld( beta, z, Z );
583  G4double am            = CalculateAm( momentum, n, Z);
584  G4double k             = momentum/hbarc;
585  G4double ch            = 0.5*n/k;
586  G4double ch2           = ch*ch;
587  G4double xsc           = ch2/(sinHalfTheta2+am)/(sinHalfTheta2+am);
588
589  return xsc;
590}
591
592
593////////////////////////////////////////////////////////////////////
594//
595// return Coulomb scattering total xsc with Wentzel correction 
596
597inline  G4double G4NuclNuclDiffuseElastic::GetCoulombTotalXsc( const G4ParticleDefinition* particle, 
598                                                             G4double momentum, G4double Z  )
599{
600  G4double beta          = CalculateParticleBeta( particle, momentum);
601  G4cout<<"beta = "<<beta<<G4endl;
602  G4double z             = particle->GetPDGCharge();
603  G4double n             = CalculateZommerfeld( beta, z, Z );
604  G4cout<<"fZomerfeld = "<<n<<G4endl;
605  G4double am            = CalculateAm( momentum, n, Z);
606  G4cout<<"cof Am = "<<am<<G4endl;
607  G4double k             = momentum/hbarc;
608  G4cout<<"k = "<<k*fermi<<" 1/fermi"<<G4endl;
609  G4cout<<"k*Bohr_radius = "<<k*Bohr_radius<<G4endl;
610  G4double ch            = n/k;
611  G4double ch2           = ch*ch;
612  G4double xsc           = ch2*pi/(am +am*am);
613
614  return xsc;
615}
616
617////////////////////////////////////////////////////////////////////
618//
619// return Coulomb scattering xsc with Wentzel correction  integrated between
620// theta1 and < theta2
621
622inline  G4double G4NuclNuclDiffuseElastic::GetCoulombIntegralXsc( const G4ParticleDefinition* particle, 
623                                 G4double momentum, G4double Z, 
624                                 G4double theta1, G4double theta2 )
625{
626  G4double c1 = std::cos(theta1);
627  G4cout<<"c1 = "<<c1<<G4endl;
628  G4double c2 = std::cos(theta2);
629  G4cout<<"c2 = "<<c2<<G4endl;
630  G4double beta          = CalculateParticleBeta( particle, momentum);
631  // G4cout<<"beta = "<<beta<<G4endl;
632  G4double z             = particle->GetPDGCharge();
633  G4double n             = CalculateZommerfeld( beta, z, Z );
634  // G4cout<<"fZomerfeld = "<<n<<G4endl;
635  G4double am            = CalculateAm( momentum, n, Z);
636  // G4cout<<"cof Am = "<<am<<G4endl;
637  G4double k             = momentum/hbarc;
638  // G4cout<<"k = "<<k*fermi<<" 1/fermi"<<G4endl;
639  // G4cout<<"k*Bohr_radius = "<<k*Bohr_radius<<G4endl;
640  G4double ch            = n/k;
641  G4double ch2           = ch*ch;
642  am *= 2.;
643  G4double xsc           = ch2*twopi*(c1-c2);
644           xsc          /= (1 - c1 + am)*(1 - c2 + am);
645
646  return xsc;
647}
648
649///////////////////////////////////////////////////////////////////
650//
651// For the calculation of arg Gamma(z) one needs complex extension
652// of ln(Gamma(z))
653
654inline G4complex G4NuclNuclDiffuseElastic::GammaLogarithm(G4complex zz)
655{
656  static G4double cof[6] = { 76.18009172947146,     -86.50532032941677,
657                             24.01409824083091,      -1.231739572450155,
658                              0.1208650973866179e-2, -0.5395239384953e-5  } ;
659  register G4int j;
660  G4complex z = zz - 1.0;
661  G4complex tmp = z + 5.5;
662  tmp -= (z + 0.5) * std::log(tmp);
663  G4complex ser = G4complex(1.000000000190015,0.);
664
665  for ( j = 0; j <= 5; j++ )
666  {
667    z += 1.0;
668    ser += cof[j]/z;
669  }
670  return -tmp + std::log(2.5066282746310005*ser);
671}
672
673///////////////////////////////////////////////////////////////////
674//
675// For the calculation of arg Gamma(z) one needs complex extension
676// of ln(Gamma(z)) here is approximate algorithm
677
678inline G4complex G4NuclNuclDiffuseElastic::GammaLogB2n(G4complex z)
679{
680  G4complex z1 = 12.*z;
681  G4complex z2 = z*z;
682  G4complex z3 = z2*z;
683  G4complex z5 = z2*z3;
684  G4complex z7 = z2*z5;
685
686  z3 *= 360.;
687  z5 *= 1260.;
688  z7 *= 1680.;
689
690  G4complex result  = (z-0.5)*std::log(z)-z+0.5*std::log(twopi);
691            result += 1./z1 - 1./z3 +1./z5 -1./z7;
692  return result;
693}
694
695/////////////////////////////////////////////////////////////////
696//
697//
698
699inline G4double  G4NuclNuclDiffuseElastic::GetErf(G4double x)
700{
701  G4double t, z, tmp, result;
702
703  z   = std::fabs(x);
704  t   = 1.0/(1.0+0.5*z);
705
706  tmp = t*exp(-z*z-1.26551223+t*(1.00002368+t*(0.37409196+t*(0.09678418+
707                t*(-0.18628806+t*(0.27886807+t*(-1.13520398+t*(1.48851587+
708                t*(-0.82215223+t*0.17087277)))))))));
709
710  if( x >= 0.) result = 1. - tmp;
711  else         result = 1. + tmp; 
712   
713  return result;
714}
715
716/////////////////////////////////////////////////////////////////
717//
718//
719
720inline G4complex G4NuclNuclDiffuseElastic::GetErfcComp(G4complex z, G4int nMax)
721{
722  G4complex erfcz = 1. - GetErfComp( z, nMax);
723  return erfcz;
724}
725
726/////////////////////////////////////////////////////////////////
727//
728//
729
730inline G4complex G4NuclNuclDiffuseElastic::GetErfcSer(G4complex z, G4int nMax)
731{
732  G4complex erfcz = 1. - GetErfSer( z, nMax);
733  return erfcz;
734}
735
736/////////////////////////////////////////////////////////////////
737//
738//
739
740inline G4complex G4NuclNuclDiffuseElastic::GetErfcInt(G4complex z) // , G4int nMax)
741{
742  G4complex erfcz = 1. - GetErfInt( z); // , nMax);
743  return erfcz;
744}
745
746inline  G4double G4NuclNuclDiffuseElastic::GetLegendrePol(G4int n, G4double theta)
747{
748  G4double legPol, epsilon = 1.e-6;
749  G4double x = std::cos(theta);
750
751  if     ( n  < 0 ) legPol = 0.;
752  else if( n == 0 ) legPol = 1.;
753  else if( n == 1 ) legPol = x;
754  else if( n == 2 ) legPol = (3.*x*x-1.)/2.;
755  else if( n == 3 ) legPol = (5.*x*x*x-3.*x)/2.;
756  else if( n == 4 ) legPol = (35.*x*x*x*x-30.*x*x+3.)/8.;
757  else if( n == 5 ) legPol = (63.*x*x*x*x*x-70.*x*x*x+15.*x)/8.;
758  else if( n == 6 ) legPol = (231.*x*x*x*x*x*x-315.*x*x*x*x+105.*x*x-5.)/16.;
759  else           
760  {
761    // legPol = ( (2*n-1)*x*GetLegendrePol(n-1,x) - (n-1)*GetLegendrePol(n-2,x) )/n;
762
763    legPol = std::sqrt( 2./(n*pi*std::sin(theta+epsilon)) )*std::sin( (n+0.5)*theta+0.25*pi );
764  }
765  return legPol; 
766}
767
768
769
770/////////////////////////////////////////////////////////////////
771//
772//
773
774inline G4complex G4NuclNuclDiffuseElastic::TestErfcComp(G4complex z, G4int nMax)
775{
776  G4complex miz = G4complex( z.imag(), -z.real() ); 
777  G4complex erfcz = 1. - GetErfComp( miz, nMax);
778  G4complex w = std::exp(-z*z)*erfcz;
779  return w;
780}
781
782/////////////////////////////////////////////////////////////////
783//
784//
785
786inline G4complex G4NuclNuclDiffuseElastic::TestErfcSer(G4complex z, G4int nMax)
787{
788  G4complex miz = G4complex( z.imag(), -z.real() ); 
789  G4complex erfcz = 1. - GetErfSer( miz, nMax);
790  G4complex w = std::exp(-z*z)*erfcz;
791  return w;
792}
793
794/////////////////////////////////////////////////////////////////
795//
796//
797
798inline G4complex G4NuclNuclDiffuseElastic::TestErfcInt(G4complex z) // , G4int nMax)
799{
800  G4complex miz = G4complex( z.imag(), -z.real() ); 
801  G4complex erfcz = 1. - GetErfInt( miz); // , nMax);
802  G4complex w = std::exp(-z*z)*erfcz;
803  return w;
804}
805
806/////////////////////////////////////////////////////////////////
807//
808//
809
810inline G4complex G4NuclNuclDiffuseElastic::GetErfComp(G4complex z, G4int nMax)
811{
812  G4int n;
813  G4double n2, cofn, shny, chny, fn, gn;
814
815  G4double x = z.real();
816  G4double y = z.imag();
817
818  G4double outRe = 0., outIm = 0.;
819
820  G4double twox  = 2.*x;
821  G4double twoxy = twox*y;
822  G4double twox2 = twox*twox;
823
824  G4double cof1 = std::exp(-x*x)/pi;
825
826  G4double cos2xy = std::cos(twoxy);
827  G4double sin2xy = std::sin(twoxy);
828
829  G4double twoxcos2xy = twox*cos2xy;
830  G4double twoxsin2xy = twox*sin2xy;
831
832  for( n = 1; n <= nMax; n++)
833  {
834    n2   = n*n;
835
836    cofn = std::exp(-0.5*n2)/(n2+twox2);  // /(n2+0.5*twox2);
837
838    chny = std::cosh(n*y);
839    shny = std::sinh(n*y);
840
841    fn  = twox - twoxcos2xy*chny + n*sin2xy*shny;
842    gn  =        twoxsin2xy*chny + n*cos2xy*shny;
843
844    fn *= cofn;
845    gn *= cofn;
846
847    outRe += fn;
848    outIm += gn;
849  }
850  outRe *= 2*cof1;
851  outIm *= 2*cof1;
852
853  if(std::abs(x) < 0.0001)
854  {
855    outRe += GetErf(x);
856    outIm += cof1*y;
857  }
858  else
859  {
860    outRe += GetErf(x) + cof1*(1-cos2xy)/twox;
861    outIm += cof1*sin2xy/twox;
862  }
863  return G4complex(outRe, outIm);
864}
865
866/////////////////////////////////////////////////////////////////
867//
868//
869
870inline G4complex G4NuclNuclDiffuseElastic::GetErfSer(G4complex z, G4int nMax)
871{
872  G4int n;
873  G4double a =1., b = 1., tmp;
874  G4complex sum = z, d = z;
875
876  for( n = 1; n <= nMax; n++)
877  {
878    a *= 2.;
879    b *= 2.*n +1.;
880    d *= z*z;
881
882    tmp = a/b;
883
884    sum += tmp*d;
885  }
886  sum *= 2.*std::exp(-z*z)/std::sqrt(pi);
887
888  return sum;
889}
890
891/////////////////////////////////////////////////////////////////////
892
893inline  G4double  G4NuclNuclDiffuseElastic::GetExpCos(G4double x)
894{
895  G4double result;
896
897  result = std::exp(x*x-fReZ*fReZ);
898  result *= std::cos(2.*x*fReZ);
899  return result;
900}
901
902/////////////////////////////////////////////////////////////////////
903
904inline  G4double  G4NuclNuclDiffuseElastic::GetExpSin(G4double x)
905{
906  G4double result;
907
908  result = std::exp(x*x-fReZ*fReZ);
909  result *= std::sin(2.*x*fReZ);
910  return result;
911}
912
913
914
915/////////////////////////////////////////////////////////////////
916//
917//
918
919inline G4complex G4NuclNuclDiffuseElastic::GetErfInt(G4complex z) // , G4int nMax)
920{
921  G4double outRe, outIm;
922
923  G4double x = z.real();
924  G4double y = z.imag();
925  fReZ       = x;
926
927  G4Integrator<G4NuclNuclDiffuseElastic,G4double(G4NuclNuclDiffuseElastic::*)(G4double)> integral;
928
929  outRe = integral.Legendre96(this,&G4NuclNuclDiffuseElastic::GetExpSin, 0., y );
930  outIm = integral.Legendre96(this,&G4NuclNuclDiffuseElastic::GetExpCos, 0., y );
931
932  outRe *= 2./sqrt(pi);
933  outIm *= 2./sqrt(pi);
934
935  outRe += GetErf(x);
936
937  return G4complex(outRe, outIm);
938}
939
940
941/////////////////////////////////////////////////////////////////
942//
943//
944
945inline  G4complex G4NuclNuclDiffuseElastic::CoulombAmplitude(G4double theta)
946{
947  G4complex ca;
948 
949  G4double sinHalfTheta  = std::sin(0.5*theta);
950  G4double sinHalfTheta2 = sinHalfTheta*sinHalfTheta; 
951  sinHalfTheta2         += fAm;
952
953  G4double order         = 2.*fCoulombPhase0 - fZommerfeld*std::log(sinHalfTheta2);
954  G4complex z            = G4complex(0., order);
955  ca                     = std::exp(z);
956
957  ca                    *= -fZommerfeld/(2.*fWaveVector*sinHalfTheta2);
958
959  return ca; 
960}
961
962/////////////////////////////////////////////////////////////////
963//
964//
965
966
967inline  void G4NuclNuclDiffuseElastic::CalculateCoulombPhaseZero()
968{
969  G4complex z        = G4complex(1,fZommerfeld); 
970  // G4complex gammalog = GammaLogarithm(z);
971  G4complex gammalog = GammaLogB2n(z);
972  fCoulombPhase0     = gammalog.imag();
973}
974
975/////////////////////////////////////////////////////////////////
976//
977//
978
979
980inline  G4double G4NuclNuclDiffuseElastic::CalculateCoulombPhase(G4int n)
981{
982  G4complex z        = G4complex(1. + n, fZommerfeld); 
983  // G4complex gammalog = GammaLogarithm(z);
984  G4complex gammalog = GammaLogB2n(z);
985  return gammalog.imag();
986}
987
988
989/////////////////////////////////////////////////////////////////
990//
991//
992
993
994inline  void G4NuclNuclDiffuseElastic::CalculateRutherfordAnglePar()
995{
996  fHalfRutThetaTg   = fZommerfeld/fProfileLambda;  // (fWaveVector*fNuclearRadius);
997  fRutherfordTheta  = 2.*std::atan(fHalfRutThetaTg);
998  fHalfRutThetaTg2  = fHalfRutThetaTg*fHalfRutThetaTg;
999  G4cout<<"fRutherfordTheta = "<<fRutherfordTheta/degree<<" degree"<<G4endl;
1000
1001}
1002
1003/////////////////////////////////////////////////////////////////
1004//
1005//
1006
1007inline   G4double G4NuclNuclDiffuseElastic::ProfileNear(G4double theta)
1008{
1009  G4double dTheta = fRutherfordTheta - theta;
1010  G4double result = 0., argument = 0.;
1011
1012  if(std::abs(dTheta) < 0.001) result = fProfileAlpha*fProfileDelta;
1013  else
1014  {
1015    argument = fProfileDelta*dTheta;
1016    result   = pi*argument*std::exp(fProfileAlpha*argument);
1017    result  /= std::sinh(pi*argument);
1018    result  -= 1.;
1019    result  /= dTheta;
1020  }
1021  return result;
1022}
1023
1024/////////////////////////////////////////////////////////////////
1025//
1026//
1027
1028inline   G4double G4NuclNuclDiffuseElastic::ProfileFar(G4double theta)
1029{
1030  G4double dTheta   = fRutherfordTheta + theta;
1031  G4double argument = fProfileDelta*dTheta;
1032
1033  G4double result   = pi*argument*std::exp(fProfileAlpha*argument);
1034  result           /= std::sinh(pi*argument);
1035  result           /= dTheta;
1036
1037  return result;
1038}
1039
1040/////////////////////////////////////////////////////////////////
1041//
1042//
1043
1044inline   G4complex G4NuclNuclDiffuseElastic::PhaseNear(G4double theta)
1045{
1046  G4double twosigma = 2.*fCoulombPhase0; 
1047  twosigma -= fZommerfeld*std::log(fHalfRutThetaTg2/(1.+fHalfRutThetaTg2));
1048  twosigma += fRutherfordTheta*fZommerfeld/fHalfRutThetaTg - halfpi;
1049  twosigma -= fProfileLambda*theta - 0.25*pi;
1050
1051  twosigma *= fCofPhase;
1052
1053  G4complex z = G4complex(0., twosigma);
1054
1055  return std::exp(z);
1056}
1057
1058/////////////////////////////////////////////////////////////////
1059//
1060//
1061
1062inline   G4complex G4NuclNuclDiffuseElastic::PhaseFar(G4double theta)
1063{
1064  G4double twosigma = 2.*fCoulombPhase0; 
1065  twosigma -= fZommerfeld*std::log(fHalfRutThetaTg2/(1.+fHalfRutThetaTg2));
1066  twosigma += fRutherfordTheta*fZommerfeld/fHalfRutThetaTg - halfpi;
1067  twosigma += fProfileLambda*theta - 0.25*pi;
1068
1069  twosigma *= fCofPhase;
1070
1071  G4complex z = G4complex(0., twosigma);
1072
1073  return std::exp(z);
1074}
1075
1076/////////////////////////////////////////////////////////////////
1077//
1078//
1079
1080
1081inline G4complex G4NuclNuclDiffuseElastic::GammaLess(G4double theta)
1082{
1083  G4double sinThetaR      = 2.*fHalfRutThetaTg/(1. + fHalfRutThetaTg2);
1084  G4double cosHalfThetaR2 = 1./(1. + fHalfRutThetaTg2);
1085
1086  G4double u              = std::sqrt(0.5*fProfileLambda/sinThetaR);
1087  G4double kappa          = u/std::sqrt(pi);
1088  G4double dTheta         = theta - fRutherfordTheta;
1089  u                      *= dTheta;
1090  G4double u2             = u*u;
1091  G4double u2m2p3         = u2*2./3.;
1092
1093  G4complex im            = G4complex(0.,1.);
1094  G4complex order         = G4complex(u,u);
1095  order                  /= std::sqrt(2.);
1096
1097  G4complex gamma         = pi*kappa*GetErfcInt(-order)*std::exp(im*(u*u+0.25*pi));
1098  G4complex a0            = 0.5*(1. + 4.*(1.+im*u2)*cosHalfThetaR2/3.)/sinThetaR;
1099  G4complex a1            = 0.5*(1. + 2.*(1.+im*u2m2p3)*cosHalfThetaR2)/sinThetaR;
1100  G4complex out           = gamma*(1. - a1*dTheta) - a0;
1101
1102  return out;
1103}
1104
1105/////////////////////////////////////////////////////////////////
1106//
1107//
1108
1109inline G4complex G4NuclNuclDiffuseElastic::GammaMore(G4double theta)
1110{
1111  G4double sinThetaR      = 2.*fHalfRutThetaTg/(1. + fHalfRutThetaTg2);
1112  G4double cosHalfThetaR2 = 1./(1. + fHalfRutThetaTg2);
1113
1114  G4double u              = std::sqrt(0.5*fProfileLambda/sinThetaR);
1115  G4double kappa          = u/std::sqrt(pi);
1116  G4double dTheta         = theta - fRutherfordTheta;
1117  u                      *= dTheta;
1118  G4double u2             = u*u;
1119  G4double u2m2p3         = u2*2./3.;
1120
1121  G4complex im            = G4complex(0.,1.);
1122  G4complex order         = G4complex(u,u);
1123  order                  /= std::sqrt(2.);
1124  G4complex gamma         = pi*kappa*GetErfcInt(order)*std::exp(im*(u*u+0.25*pi));
1125  G4complex a0            = 0.5*(1. + 4.*(1.+im*u2)*cosHalfThetaR2/3.)/sinThetaR;
1126  G4complex a1            = 0.5*(1. + 2.*(1.+im*u2m2p3)*cosHalfThetaR2)/sinThetaR;
1127  G4complex out           = -gamma*(1. - a1*dTheta) - a0;
1128
1129  return out;
1130}
1131
1132/////////////////////////////////////////////////////////////////
1133//
1134//
1135
1136inline  G4complex G4NuclNuclDiffuseElastic::AmplitudeNear(G4double theta)
1137{
1138  G4double kappa = std::sqrt(0.5*fProfileLambda/std::sin(theta)/pi);
1139  G4complex out = G4complex(kappa/fWaveVector,0.);
1140
1141  out *= PhaseNear(theta);
1142
1143  if( theta <= fRutherfordTheta )
1144  {
1145    out *= GammaLess(theta) + ProfileNear(theta);
1146    // out *= GammaMore(theta) + ProfileNear(theta);
1147    out += CoulombAmplitude(theta);
1148  }
1149  else
1150  {
1151    out *= GammaMore(theta) + ProfileNear(theta);
1152    // out *= GammaLess(theta) + ProfileNear(theta);
1153  }
1154  return out;
1155}
1156
1157/////////////////////////////////////////////////////////////////
1158//
1159//
1160
1161inline  G4complex G4NuclNuclDiffuseElastic::AmplitudeFar(G4double theta)
1162{
1163  G4double kappa = std::sqrt(0.5*fProfileLambda/std::sin(theta)/pi);
1164  G4complex out = G4complex(kappa/fWaveVector,0.);
1165  out *= ProfileFar(theta);
1166  out *= PhaseFar(theta);
1167  return out;
1168}
1169
1170
1171/////////////////////////////////////////////////////////////////
1172//
1173//
1174
1175inline  G4complex G4NuclNuclDiffuseElastic::Amplitude(G4double theta)
1176{
1177 
1178  G4complex out = AmplitudeNear(theta) + fCofFar*AmplitudeFar(theta);
1179  // G4complex out = AmplitudeNear(theta);
1180  // G4complex out = AmplitudeFar(theta);
1181  return    out;
1182}
1183
1184/////////////////////////////////////////////////////////////////
1185//
1186//
1187
1188inline  G4double  G4NuclNuclDiffuseElastic::AmplitudeMod2(G4double theta)
1189{
1190  G4complex out = Amplitude(theta);
1191  G4double mod2 = out.real()*out.real() + out.imag()*out.imag();
1192  return   mod2;
1193}
1194
1195/////////////////////////////////////////////////////////////////
1196//
1197//
1198
1199inline  G4complex G4NuclNuclDiffuseElastic::AmplitudeGla(G4double theta)
1200{
1201  G4int n;
1202  G4double T12b, b, b2; // cosTheta = std::cos(theta);
1203  G4complex out = G4complex(0.,0.), shiftC, shiftN; 
1204  G4complex im  = G4complex(0.,1.);
1205
1206  for( n = 0; n < fMaxL; n++)
1207  {
1208    shiftC = std::exp( im*2.*CalculateCoulombPhase(n) );
1209    // b = ( fZommerfeld + std::sqrt( fZommerfeld*fZommerfeld + n*(n+1) ) )/fWaveVector;
1210    b = ( std::sqrt( G4double(n*(n+1)) ) )/fWaveVector;
1211    b2 = b*b;
1212    T12b = fSumSigma*std::exp(-b2/fNuclearRadiusSquare)/pi/fNuclearRadiusSquare;         
1213    shiftN = std::exp( -0.5*(1.-im*fEtaRatio)*T12b ) - 1.;
1214    out +=  (2.*n+1.)*shiftC*shiftN*GetLegendrePol(n, theta);   
1215  }
1216  out /= 2.*im*fWaveVector;
1217  out += CoulombAmplitude(theta);
1218  return out;
1219}
1220
1221/////////////////////////////////////////////////////////////////
1222//
1223//
1224
1225inline  G4double  G4NuclNuclDiffuseElastic::AmplitudeGlaMod2(G4double theta)
1226{
1227  G4complex out = AmplitudeGla(theta);
1228  G4double mod2 = out.real()*out.real() + out.imag()*out.imag();
1229  return   mod2;
1230}
1231
1232
1233/////////////////////////////////////////////////////////////////
1234//
1235//
1236
1237inline  G4complex G4NuclNuclDiffuseElastic::AmplitudeGG(G4double theta)
1238{
1239  G4int n;
1240  G4double T12b, a, aTemp, b2, sinThetaH = std::sin(0.5*theta);
1241  G4double sinThetaH2 = sinThetaH*sinThetaH;
1242  G4complex out = G4complex(0.,0.); 
1243  G4complex im  = G4complex(0.,1.);
1244
1245  a  = -fSumSigma/twopi/fNuclearRadiusSquare;
1246  b2 = fWaveVector*fWaveVector*fNuclearRadiusSquare*sinThetaH2;
1247
1248  aTemp = a;
1249
1250  for( n = 1; n < fMaxL; n++)
1251  {   
1252    T12b   = aTemp*std::exp(-b2/n)/n;         
1253    aTemp *= a; 
1254    out   += T12b; 
1255    G4cout<<"out = "<<out<<G4endl; 
1256  }
1257  out *= -4.*im*fWaveVector/pi;
1258  out += CoulombAmplitude(theta);
1259  return out;
1260}
1261
1262/////////////////////////////////////////////////////////////////
1263//
1264//
1265
1266inline  G4double  G4NuclNuclDiffuseElastic::AmplitudeGGMod2(G4double theta)
1267{
1268  G4complex out = AmplitudeGG(theta);
1269  G4double mod2 = out.real()*out.real() + out.imag()*out.imag();
1270  return   mod2;
1271}
1272
1273
1274///////////////////////////////////////////////////////////////////////////////
1275//
1276// Test for given particle and element table of momentum, angle probability.
1277// For the partMom in CMS.
1278
1279inline void G4NuclNuclDiffuseElastic::InitParameters(const G4ParticleDefinition* theParticle, 
1280                                          G4double partMom, G4double Z, G4double A) 
1281{
1282  fAtomicNumber  = Z;     // atomic number
1283  fAtomicWeight  = A;     // number of nucleons
1284
1285  fNuclearRadius2 = CalculateNuclearRad(fAtomicWeight);
1286  G4double A1     = G4double( theParticle->GetBaryonNumber() );   
1287  fNuclearRadius1 = CalculateNuclearRad(A1);
1288  // fNuclearRadius = std::sqrt(fNuclearRadius1*fNuclearRadius1+fNuclearRadius2*fNuclearRadius2);
1289  fNuclearRadius = fNuclearRadius1 + fNuclearRadius2;
1290
1291  G4double a = 0.;
1292  G4double z = theParticle->GetPDGCharge();
1293  G4double m1 = theParticle->GetPDGMass();
1294
1295  fWaveVector = partMom/hbarc;
1296
1297  G4double lambda = fCofLambda*fWaveVector*fNuclearRadius;
1298  G4cout<<"kR = "<<lambda<<G4endl;
1299
1300  if( z )
1301  {
1302    a           = partMom/m1; // beta*gamma for m1
1303    fBeta       = a/std::sqrt(1+a*a);
1304    fZommerfeld = CalculateZommerfeld( fBeta, z, fAtomicNumber);
1305    fAm         = CalculateAm( partMom, fZommerfeld, fAtomicNumber);
1306  }
1307  fProfileLambda = lambda*std::sqrt(1.-2*fZommerfeld/lambda);
1308  G4cout<<"fProfileLambda = "<<fProfileLambda<<G4endl;
1309  fProfileDelta  = fCofDelta*fProfileLambda;
1310  fProfileAlpha  = fCofAlpha*fProfileLambda;
1311
1312  CalculateCoulombPhaseZero();
1313  CalculateRutherfordAnglePar();
1314
1315  return;
1316}
1317
1318
1319///////////////////////////////////////////////////////////////////////////////
1320//
1321// Test for given particle and element table of momentum, angle probability.
1322// For the partMom in CMS.
1323
1324inline void G4NuclNuclDiffuseElastic::InitParametersGla(const G4DynamicParticle* aParticle, 
1325                                          G4double partMom, G4double Z, G4double A) 
1326{
1327  fAtomicNumber  = Z;     // target atomic number
1328  fAtomicWeight  = A;     // target number of nucleons
1329
1330  fNuclearRadius2 = CalculateNuclearRad(fAtomicWeight); // target nucleus radius
1331  G4double A1     = G4double( aParticle->GetDefinition()->GetBaryonNumber() );   
1332  fNuclearRadius1 = CalculateNuclearRad(A1); // projectile nucleus radius
1333  fNuclearRadiusSquare = fNuclearRadius1*fNuclearRadius1+fNuclearRadius2*fNuclearRadius2;
1334 
1335
1336  G4double a  = 0., kR12;
1337  G4double z  = aParticle->GetDefinition()->GetPDGCharge();
1338  G4double m1 = aParticle->GetDefinition()->GetPDGMass();
1339
1340  fWaveVector = partMom/hbarc;
1341
1342  G4double pN = A1 - z;
1343  if( pN < 0. ) pN = 0.;
1344
1345  G4double tN = A - Z;
1346  if( tN < 0. ) tN = 0.;
1347
1348  G4double pTkin = aParticle->GetKineticEnergy(); 
1349  pTkin /= A1;
1350
1351
1352  fSumSigma = (Z*z+pN*tN)*GetHadronNucleonXscNS(theProton, pTkin, theProton) +
1353              (z*tN+pN*Z)*GetHadronNucleonXscNS(theProton, pTkin, theNeutron);
1354
1355  G4cout<<"fSumSigma = "<<fSumSigma/millibarn<<" mb"<<G4endl;
1356  G4cout<<"pi*R2 = "<<pi*fNuclearRadiusSquare/millibarn<<" mb"<<G4endl;
1357  kR12 = fWaveVector*std::sqrt(fNuclearRadiusSquare);
1358  G4cout<<"k*sqrt(R2) = "<<kR12<<" "<<G4endl;
1359  fMaxL = (G4int(kR12)+1)*4;
1360  G4cout<<"fMaxL = "<<fMaxL<<" "<<G4endl;
1361
1362  if( z )
1363  {
1364      a           = partMom/m1; // beta*gamma for m1
1365      fBeta       = a/std::sqrt(1+a*a);
1366      fZommerfeld = CalculateZommerfeld( fBeta, z, fAtomicNumber);
1367      fAm         = CalculateAm( partMom, fZommerfeld, fAtomicNumber);
1368  }
1369
1370  CalculateCoulombPhaseZero();
1371 
1372
1373  return;
1374}
1375
1376
1377/////////////////////////////////////////////////////////////////////////////////////
1378//
1379// Returns nucleon-nucleon cross-section based on N. Starkov parametrisation of
1380// data from mainly http://wwwppds.ihep.su:8001/c5-6A.html database
1381// projectile nucleon is pParticle with pTkin shooting target nucleon tParticle
1382
1383inline G4double
1384G4NuclNuclDiffuseElastic::GetHadronNucleonXscNS( G4ParticleDefinition* pParticle, 
1385                                                 G4double pTkin, 
1386                                                 G4ParticleDefinition* tParticle)
1387{
1388  G4double xsection(0), Delta, A0, B0;
1389  G4double hpXsc(0);
1390  G4double hnXsc(0);
1391
1392
1393  G4double targ_mass     = tParticle->GetPDGMass();
1394  G4double proj_mass     = pParticle->GetPDGMass(); 
1395
1396  G4double proj_energy   = proj_mass + pTkin; 
1397  G4double proj_momentum = std::sqrt(pTkin*(pTkin+2*proj_mass));
1398
1399  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
1400
1401  sMand         /= GeV*GeV;  // in GeV for parametrisation
1402  proj_momentum /= GeV;
1403  proj_energy   /= GeV;
1404  proj_mass     /= GeV;
1405  G4double logS = std::log(sMand);
1406
1407  // General PDG fit constants
1408
1409
1410  // fEtaRatio=Re[f(0)]/Im[f(0)]
1411
1412  if( proj_momentum >= 1.2 )
1413  {
1414    fEtaRatio  = 0.13*(logS - 5.8579332)*std::pow(sMand,-0.18);
1415  }       
1416  else if( proj_momentum >= 0.6 )
1417  { 
1418    fEtaRatio = -75.5*(std::pow(proj_momentum,0.25)-0.95)/
1419          (std::pow(3*proj_momentum,2.2)+1);     
1420  }
1421  else 
1422  {
1423    fEtaRatio = 15.5*proj_momentum/(27*proj_momentum*proj_momentum*proj_momentum+2);
1424  }
1425  G4cout<<"fEtaRatio = "<<fEtaRatio<<G4endl;
1426
1427  // xsc
1428 
1429  if( proj_momentum >= 10. ) // high energy: pp = nn = np
1430    // if( proj_momentum >= 2.)
1431  {
1432    Delta = 1.;
1433
1434    if( proj_energy < 40. ) Delta = 0.916+0.0021*proj_energy;
1435
1436    if( proj_momentum >= 10.)
1437    {
1438        B0 = 7.5;
1439        A0 = 100. - B0*std::log(3.0e7);
1440
1441        xsection = A0 + B0*std::log(proj_energy) - 11
1442                  + 103*std::pow(2*0.93827*proj_energy + proj_mass*proj_mass+
1443                     0.93827*0.93827,-0.165);        //  mb
1444    }
1445  }
1446  else // low energy pp = nn != np
1447  {
1448      if(pParticle == tParticle) // pp or nn      // nn to be pp
1449      {
1450        if( proj_momentum < 0.73 )
1451        {
1452          hnXsc = 23 + 50*( std::pow( std::log(0.73/proj_momentum), 3.5 ) );
1453        }
1454        else if( proj_momentum < 1.05  )
1455        {
1456          hnXsc = 23 + 40*(std::log(proj_momentum/0.73))*
1457                         (std::log(proj_momentum/0.73));
1458        }
1459        else  // if( proj_momentum < 10.  )
1460        {
1461          hnXsc = 39.0 + 
1462              75*(proj_momentum - 1.2)/(std::pow(proj_momentum,3.0) + 0.15);
1463        }
1464        xsection = hnXsc;
1465      }
1466      else  // pn to be np
1467      {
1468        if( proj_momentum < 0.8 )
1469        {
1470          hpXsc = 33+30*std::pow(std::log(proj_momentum/1.3),4.0);
1471        }     
1472        else if( proj_momentum < 1.4 )
1473        {
1474          hpXsc = 33+30*std::pow(std::log(proj_momentum/0.95),2.0);
1475        }
1476        else    // if( proj_momentum < 10.  )
1477        {
1478          hpXsc = 33.3+
1479              20.8*(std::pow(proj_momentum,2.0)-1.35)/
1480                 (std::pow(proj_momentum,2.50)+0.95);
1481        }
1482        xsection = hpXsc;
1483      }
1484  }
1485  xsection *= millibarn; // parametrised in mb
1486  G4cout<<"xsection = "<<xsection/millibarn<<" mb"<<G4endl;
1487  return xsection;
1488}
1489
1490////////////////////////////////////////////////////////////////////////////////////
1491//
1492//
1493
1494inline G4double G4NuclNuclDiffuseElastic::CalcMandelstamS( const G4double mp , 
1495                                                       const G4double mt , 
1496                                                       const G4double Plab )
1497{
1498  G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
1499  G4double sMand  = mp*mp + mt*mt + 2*Elab*mt ;
1500
1501  return sMand;
1502}
1503
1504
1505
1506#endif
Note: See TracBrowser for help on using the repository browser.