source: trunk/source/processes/hadronic/models/coherent_elastic/include/G4DiffuseElastic.hh @ 1007

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

update to geant4.9.2

File size: 16.0 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: G4DiffuseElastic.hh,v 1.13 2007/11/06 17:01:20 grichine Exp $
28// GEANT4 tag $Name: geant4-09-02 $
29//
30//
31// G4 Model: optical elastic scattering with 4-momentum balance
32//
33// Class Description
34// Final state production model for hadron nuclear elastic scattering;
35// Class Description - End
36//
37//
38// 24.05.07 V. Grichine first implementation for hadron (no Coulomb) elastic scattering
39// 04.09.07 V. Grichine implementation for Coulomb elastic scattering
40
41
42#ifndef G4DiffuseElastic_h
43#define G4DiffuseElastic_h 1
44 
45#include "globals.hh"
46#include "G4HadronicInteraction.hh"
47#include "G4HadProjectile.hh"
48#include "G4Nucleus.hh"
49
50using namespace std;
51
52class G4ParticleDefinition;
53class G4PhysicsTable;
54class G4PhysicsLogVector;
55
56class G4DiffuseElastic : public G4HadronicInteraction
57{
58public:
59
60  G4DiffuseElastic();
61
62  G4DiffuseElastic(const G4ParticleDefinition* aParticle);
63
64
65
66
67
68  virtual ~G4DiffuseElastic();
69
70  void Initialise();
71
72  void InitialiseOnFly(G4double Z, G4double A);
73
74  void BuildAngleTable();
75
76 
77  G4HadFinalState * ApplyYourself(const G4HadProjectile & aTrack, 
78                                  G4Nucleus & targetNucleus);
79
80
81  void SetPlabLowLimit(G4double value);
82
83  void SetHEModelLowLimit(G4double value);
84
85  void SetQModelLowLimit(G4double value);
86
87  void SetLowestEnergyLimit(G4double value);
88
89  void SetRecoilKinEnergyLimit(G4double value);
90
91  G4double SampleT(const G4ParticleDefinition* aParticle, 
92                         G4double p, G4double A);
93
94  G4double SampleTableT(const G4ParticleDefinition* aParticle, 
95                         G4double p, G4double Z, G4double A);
96
97  G4double SampleThetaCMS(const G4ParticleDefinition* aParticle, G4double p, G4double A);
98
99  G4double SampleTableThetaCMS(const G4ParticleDefinition* aParticle, G4double p, 
100                                     G4double Z, G4double A);
101
102  G4double GetScatteringAngle(G4int iMomentum, G4int iAngle, G4double position);
103
104  G4double SampleThetaLab(const G4HadProjectile* aParticle, 
105                                G4double tmass, G4double A);
106
107  G4double GetDiffuseElasticXsc( const G4ParticleDefinition* particle, 
108                                 G4double theta, 
109                                 G4double momentum, 
110                                 G4double A         );
111
112  G4double GetInvElasticXsc( const G4ParticleDefinition* particle, 
113                                 G4double theta, 
114                                 G4double momentum, 
115                                 G4double A, G4double Z );
116
117  G4double GetDiffuseElasticSumXsc( const G4ParticleDefinition* particle, 
118                                 G4double theta, 
119                                 G4double momentum, 
120                                 G4double A, G4double Z );
121
122  G4double GetInvElasticSumXsc( const G4ParticleDefinition* particle, 
123                                 G4double tMand, 
124                                 G4double momentum, 
125                                 G4double A, G4double Z );
126
127  G4double IntegralElasticProb( const G4ParticleDefinition* particle, 
128                                 G4double theta, 
129                                 G4double momentum, 
130                                 G4double A            );
131 
132
133  G4double GetCoulombElasticXsc( const G4ParticleDefinition* particle, 
134                                 G4double theta, 
135                                 G4double momentum, 
136                                 G4double Z         );
137
138  G4double GetInvCoulombElasticXsc( const G4ParticleDefinition* particle, 
139                                 G4double tMand, 
140                                 G4double momentum, 
141                                 G4double A, G4double Z         );
142
143  G4double GetCoulombTotalXsc( const G4ParticleDefinition* particle, 
144                                 G4double momentum, G4double Z       );
145
146  G4double GetCoulombIntegralXsc( const G4ParticleDefinition* particle, 
147                                 G4double momentum, G4double Z, 
148                                 G4double theta1, G4double theta2         );
149
150
151  G4double CalculateParticleBeta( const G4ParticleDefinition* particle, 
152                                        G4double momentum    );
153
154  G4double CalculateZommerfeld( G4double beta, G4double Z1, G4double Z2 );
155
156  G4double CalculateAm( G4double momentum, G4double n, G4double Z);
157
158  G4double CalculateNuclearRad( G4double A);
159
160  G4double ThetaCMStoThetaLab(const G4DynamicParticle* aParticle, 
161                                G4double tmass, G4double thetaCMS);
162
163  G4double ThetaLabToThetaCMS(const G4DynamicParticle* aParticle, 
164                                G4double tmass, G4double thetaLab);
165
166
167
168
169
170  G4double BesselJzero(G4double z);
171  G4double BesselJone(G4double z);
172  G4double DampFactor(G4double z);
173  G4double BesselOneByArg(G4double z);
174
175  G4double GetDiffElasticProb(G4double theta);
176  G4double GetDiffElasticSumProb(G4double theta);
177  G4double GetIntegrandFunction(G4double theta);
178
179
180  G4double GetNuclearRadius(){return fNuclearRadius;};
181
182private:
183
184
185  G4ParticleDefinition* theProton;
186  G4ParticleDefinition* theNeutron;
187  G4ParticleDefinition* theDeuteron;
188  G4ParticleDefinition* theAlpha;
189
190  const G4ParticleDefinition* thePionPlus;
191  const G4ParticleDefinition* thePionMinus;
192
193  G4double lowEnergyRecoilLimit; 
194  G4double lowEnergyLimitHE; 
195  G4double lowEnergyLimitQ; 
196  G4double lowestEnergyLimit; 
197  G4double plabLowLimit;
198
199  G4int fEnergyBin;
200  G4int fAngleBin;
201
202  G4PhysicsLogVector*           fEnergyVector;
203  G4PhysicsTable*               fAngleTable;
204  std::vector<G4PhysicsTable*>  fAngleBank;
205
206  std::vector<G4double> fElementNumberVector;
207  std::vector<G4String> fElementNameVector;
208
209  const G4ParticleDefinition* fParticle;
210  G4double fWaveVector;
211  G4double fAtomicWeight;
212  G4double fAtomicNumber;
213  G4double fNuclearRadius;
214  G4double fBeta;
215  G4double fZommerfeld;
216  G4double fAm;
217  G4bool fAddCoulomb;
218
219};
220
221
222inline void G4DiffuseElastic::SetRecoilKinEnergyLimit(G4double value)
223{
224  lowEnergyRecoilLimit = value;
225}
226
227inline void G4DiffuseElastic::SetPlabLowLimit(G4double value)
228{
229  plabLowLimit = value;
230}
231
232inline void G4DiffuseElastic::SetHEModelLowLimit(G4double value)
233{
234  lowEnergyLimitHE = value;
235}
236
237inline void G4DiffuseElastic::SetQModelLowLimit(G4double value)
238{
239  lowEnergyLimitQ = value;
240}
241
242inline void G4DiffuseElastic::SetLowestEnergyLimit(G4double value)
243{
244  lowestEnergyLimit = value;
245}
246
247
248/////////////////////////////////////////////////////////////
249//
250// Bessel J0 function based on rational approximation from
251// J.F. Hart, Computer Approximations, New York, Willey 1968, p. 141
252
253inline G4double G4DiffuseElastic::BesselJzero(G4double value)
254{
255  G4double modvalue, value2, fact1, fact2, arg, shift, bessel;
256
257  modvalue = fabs(value);
258
259  if ( value < 8.0 && value > -8.0 )
260  {
261    value2 = value*value;
262
263    fact1  = 57568490574.0 + value2*(-13362590354.0 
264                           + value2*( 651619640.7 
265                           + value2*(-11214424.18 
266                           + value2*( 77392.33017 
267                           + value2*(-184.9052456   ) ) ) ) );
268
269    fact2  = 57568490411.0 + value2*( 1029532985.0 
270                           + value2*( 9494680.718
271                           + value2*(59272.64853
272                           + value2*(267.8532712 
273                           + value2*1.0               ) ) ) );
274
275    bessel = fact1/fact2;
276  } 
277  else 
278  {
279    arg    = 8.0/modvalue;
280
281    value2 = arg*arg;
282
283    shift  = modvalue-0.785398164;
284
285    fact1  = 1.0 + value2*(-0.1098628627e-2 
286                 + value2*(0.2734510407e-4
287                 + value2*(-0.2073370639e-5 
288                 + value2*0.2093887211e-6    ) ) );
289
290    fact2  = -0.1562499995e-1 + value2*(0.1430488765e-3
291                              + value2*(-0.6911147651e-5 
292                              + value2*(0.7621095161e-6
293                              - value2*0.934945152e-7    ) ) );
294
295    bessel = sqrt(0.636619772/modvalue)*(cos(shift)*fact1 - arg*sin(shift)*fact2 );
296  }
297  return bessel;
298}
299
300/////////////////////////////////////////////////////////////
301//
302// Bessel J1 function based on rational approximation from
303// J.F. Hart, Computer Approximations, New York, Willey 1968, p. 141
304
305inline G4double G4DiffuseElastic::BesselJone(G4double value)
306{
307  G4double modvalue, value2, fact1, fact2, arg, shift, bessel;
308
309  modvalue = fabs(value);
310
311  if ( modvalue < 8.0 ) 
312  {
313    value2 = value*value;
314
315    fact1  = value*(72362614232.0 + value2*(-7895059235.0 
316                                  + value2*( 242396853.1
317                                  + value2*(-2972611.439 
318                                  + value2*( 15704.48260 
319                                  + value2*(-30.16036606  ) ) ) ) ) );
320
321    fact2  = 144725228442.0 + value2*(2300535178.0 
322                            + value2*(18583304.74
323                            + value2*(99447.43394 
324                            + value2*(376.9991397 
325                            + value2*1.0             ) ) ) );
326    bessel = fact1/fact2;
327  } 
328  else 
329  {
330    arg    = 8.0/modvalue;
331
332    value2 = arg*arg;
333
334    shift  = modvalue - 2.356194491;
335
336    fact1  = 1.0 + value2*( 0.183105e-2 
337                 + value2*(-0.3516396496e-4
338                 + value2*(0.2457520174e-5 
339                 + value2*(-0.240337019e-6          ) ) ) );
340
341    fact2 = 0.04687499995 + value2*(-0.2002690873e-3
342                          + value2*( 0.8449199096e-5
343                          + value2*(-0.88228987e-6
344                          + value2*0.105787412e-6       ) ) );
345
346    bessel = sqrt( 0.636619772/modvalue)*(cos(shift)*fact1 - arg*sin(shift)*fact2);
347
348    if (value < 0.0) bessel = -bessel;
349  }
350  return bessel;
351}
352
353////////////////////////////////////////////////////////////////////
354//
355// damp factor in diffraction x/sh(x), x was already *pi
356
357inline G4double G4DiffuseElastic::DampFactor(G4double x)
358{
359  G4double df;
360  G4double f2 = 2., f3 = 6., f4 = 24.; // first factorials
361
362  // x *= pi;
363
364  if( std::fabs(x) < 0.01 )
365  { 
366    df = 1./(1. + x/f2 + x*x/f3 + x*x*x/f4);
367  }
368  else
369  {
370    df = x/std::sinh(x); 
371  }
372  return df;
373}
374
375
376////////////////////////////////////////////////////////////////////
377//
378// return J1(x)/x with special case for small x
379
380inline G4double G4DiffuseElastic::BesselOneByArg(G4double x)
381{
382  G4double x2, result;
383 
384  if( std::fabs(x) < 0.01 )
385  { 
386   x     *= 0.5;
387   x2     = x*x;
388   result = 2. - x2 + x2*x2/6.;
389  }
390  else
391  {
392    result = BesselJone(x)/x; 
393  }
394  return result;
395}
396
397////////////////////////////////////////////////////////////////////
398//
399// return particle beta
400
401inline  G4double G4DiffuseElastic::CalculateParticleBeta( const G4ParticleDefinition* particle, 
402                                        G4double momentum    )
403{
404  G4double mass = particle->GetPDGMass();
405  G4double a    = momentum/mass;
406  fBeta         = a/std::sqrt(1+a*a);
407
408  return fBeta; 
409}
410
411////////////////////////////////////////////////////////////////////
412//
413// return Zommerfeld parameter for Coulomb scattering
414
415inline  G4double G4DiffuseElastic::CalculateZommerfeld( G4double beta, G4double Z1, G4double Z2 )
416{
417  fZommerfeld = fine_structure_const*Z1*Z2/beta;
418
419  return fZommerfeld; 
420}
421
422////////////////////////////////////////////////////////////////////
423//
424// return Wentzel correction for Coulomb scattering
425
426inline  G4double G4DiffuseElastic::CalculateAm( G4double momentum, G4double n, G4double Z)
427{
428  G4double k   = momentum/hbarc;
429  G4double ch  = 1.13 + 3.76*n*n;
430  G4double zn  = 1.77*k*std::pow(Z,-1./3.)*Bohr_radius;
431  G4double zn2 = zn*zn;
432  fAm          = ch/zn2;
433
434  return fAm;
435}
436
437////////////////////////////////////////////////////////////////////
438//
439// calculate nuclear radius for different atomic weights using different approximations
440
441inline  G4double G4DiffuseElastic::CalculateNuclearRad( G4double A)
442{
443  G4double r0;
444
445  if(A < 50.)
446  {
447    if(A > 10.) r0  = 1.16*( 1 - std::pow(A, -2./3.) )*fermi;   // 1.08*fermi;
448    else        r0  = 1.1*fermi;
449
450    fNuclearRadius = r0*std::pow(A, 1./3.);
451  }
452  else
453  {
454    r0 = 1.7*fermi;
455
456    fNuclearRadius = r0*std::pow(A, 0.27);
457  }
458  return fNuclearRadius;
459}
460
461////////////////////////////////////////////////////////////////////
462//
463// return Coulomb scattering differential xsc with Wentzel correction 
464
465inline  G4double G4DiffuseElastic::GetCoulombElasticXsc( const G4ParticleDefinition* particle, 
466                                 G4double theta, 
467                                 G4double momentum, 
468                                 G4double Z         )
469{
470  G4double sinHalfTheta  = std::sin(0.5*theta);
471  G4double sinHalfTheta2 = sinHalfTheta*sinHalfTheta;
472  G4double beta          = CalculateParticleBeta( particle, momentum);
473  G4double z             = particle->GetPDGCharge();
474  G4double n             = CalculateZommerfeld( beta, z, Z );
475  G4double am            = CalculateAm( momentum, n, Z);
476  G4double k             = momentum/hbarc;
477  G4double ch            = 0.5*n/k;
478  G4double ch2           = ch*ch;
479  G4double xsc           = ch2/(sinHalfTheta2+am)/(sinHalfTheta2+am);
480
481  return xsc;
482}
483
484
485////////////////////////////////////////////////////////////////////
486//
487// return Coulomb scattering total xsc with Wentzel correction 
488
489inline  G4double G4DiffuseElastic::GetCoulombTotalXsc( const G4ParticleDefinition* particle, 
490                                                             G4double momentum, G4double Z  )
491{
492  G4double beta          = CalculateParticleBeta( particle, momentum);
493  G4cout<<"beta = "<<beta<<G4endl;
494  G4double z             = particle->GetPDGCharge();
495  G4double n             = CalculateZommerfeld( beta, z, Z );
496  G4cout<<"fZomerfeld = "<<n<<G4endl;
497  G4double am            = CalculateAm( momentum, n, Z);
498  G4cout<<"cof Am = "<<am<<G4endl;
499  G4double k             = momentum/hbarc;
500  G4cout<<"k = "<<k*fermi<<" 1/fermi"<<G4endl;
501  G4cout<<"k*Bohr_radius = "<<k*Bohr_radius<<G4endl;
502  G4double ch            = n/k;
503  G4double ch2           = ch*ch;
504  G4double xsc           = ch2*pi/(am +am*am);
505
506  return xsc;
507}
508
509////////////////////////////////////////////////////////////////////
510//
511// return Coulomb scattering xsc with Wentzel correction  integrated between
512// theta1 and < theta2
513
514inline  G4double G4DiffuseElastic::GetCoulombIntegralXsc( const G4ParticleDefinition* particle, 
515                                 G4double momentum, G4double Z, 
516                                 G4double theta1, G4double theta2 )
517{
518  G4double c1 = std::cos(theta1);
519  G4cout<<"c1 = "<<c1<<G4endl;
520  G4double c2 = std::cos(theta2);
521  G4cout<<"c2 = "<<c2<<G4endl;
522  G4double beta          = CalculateParticleBeta( particle, momentum);
523  // G4cout<<"beta = "<<beta<<G4endl;
524  G4double z             = particle->GetPDGCharge();
525  G4double n             = CalculateZommerfeld( beta, z, Z );
526  // G4cout<<"fZomerfeld = "<<n<<G4endl;
527  G4double am            = CalculateAm( momentum, n, Z);
528  // G4cout<<"cof Am = "<<am<<G4endl;
529  G4double k             = momentum/hbarc;
530  // G4cout<<"k = "<<k*fermi<<" 1/fermi"<<G4endl;
531  // G4cout<<"k*Bohr_radius = "<<k*Bohr_radius<<G4endl;
532  G4double ch            = n/k;
533  G4double ch2           = ch*ch;
534  am *= 2.;
535  G4double xsc           = ch2*twopi*(c1-c2);
536           xsc          /= (1 - c1 + am)*(1 - c2 + am);
537
538  return xsc;
539}
540
541#endif
Note: See TracBrowser for help on using the repository browser.