source: trunk/source/processes/hadronic/models/coherent_elastic/src/G4HadronElastic.cc @ 847

Last change on this file since 847 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 15.2 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: G4HadronElastic.cc,v 1.56.2.1 2008/04/23 14:14:55 gcosmo Exp $
27// GEANT4 tag $Name: geant4-09-01-patch-02 $
28//
29//
30// Physics model class G4HadronElastic (derived from G4LElastic)
31//
32//
33// G4 Model: Low-energy Elastic scattering with 4-momentum balance
34// F.W. Jones, TRIUMF, 04-JUN-96
35// Uses  G4ElasticHadrNucleusHE and G4VQCrossSection
36//
37//
38// 25-JUN-98 FWJ: replaced missing Initialize for ParticleChange.
39// 09-Set-05 V.Ivanchenko HARP version of the model: fix scattering
40//           on hydrogen, use relativistic Lorentz transformation
41// 24-Nov-05 V.Ivanchenko sample cost in center of mass reference system
42// 03-Dec-05 V.Ivanchenko add protection to initial momentum 20 MeV/c in
43//           center of mass system (before it was in lab system)
44//           below model is not valid
45// 14-Dec-05 V.Ivanchenko change protection to cos(theta) < -1 and
46//           rename the class
47// 13-Apr-06 V.Ivanchenko move to coherent_elastic subdirectory; remove
48//           charge exchange; remove limitation on incident momentum;
49//           add s-wave regim below some momentum       
50// 24-Apr-06 V.Ivanchenko add neutron scattering on hydrogen from CHIPS
51// 07-Jun-06 V.Ivanchenko fix problem of rotation
52// 25-Jul-06 V.Ivanchenko add 19 MeV low energy, below which S-wave is sampled
53// 02-Aug-06 V.Ivanchenko introduce energy cut on the aria of S-wave for pions
54// 24-Aug-06 V.Ivanchenko switch on G4ElasticHadrNucleusHE
55// 31-Aug-06 V.Ivanchenko do not sample sacttering for particles with kinetic
56//                        energy below 10 keV
57// 16-Nov-06 V.Ivanchenko Simplify logic of choosing of the model for sampling
58// 30-Mar-07 V.Ivanchenko lowEnergyLimitQ=0, lowEnergyLimitHE = 1.0*GeV,
59//                        lowestEnergyLimit= 0
60// 04-May-07 V.Ivanchenko do not use HE model for hydrogen target to avoid NaN;
61//                        use QElastic for p, n incident for any energy for
62//                        p and He targets only 
63// 11-May-07 V.Ivanchenko remove unused method Defs1
64//
65
66#include "G4HadronElastic.hh"
67#include "G4ParticleTable.hh"
68#include "G4ParticleDefinition.hh"
69#include "G4IonTable.hh"
70#include "G4QElasticCrossSection.hh"
71#include "G4VQCrossSection.hh"
72#include "G4ElasticHadrNucleusHE.hh"
73#include "Randomize.hh"
74#include "G4Proton.hh"
75#include "G4Neutron.hh"
76#include "G4Deuteron.hh"
77#include "G4Alpha.hh"
78#include "G4PionPlus.hh"
79#include "G4PionMinus.hh"
80
81G4HadronElastic::G4HadronElastic(G4ElasticHadrNucleusHE* HModel) 
82  : G4HadronicInteraction("G4HadronElastic"), hElastic(HModel)
83{
84  SetMinEnergy( 0.0*GeV );
85  SetMaxEnergy( 100.*TeV );
86  verboseLevel= 0;
87  lowEnergyRecoilLimit = 100.*keV; 
88  lowEnergyLimitQ  = 0.0*GeV; 
89  lowEnergyLimitHE = 1.0*GeV; 
90  lowestEnergyLimit= 1.e-6*eV; 
91  plabLowLimit     = 20.0*MeV;
92
93  qCManager   = G4QElasticCrossSection::GetPointer();
94  if(!hElastic) hElastic = new G4ElasticHadrNucleusHE();
95
96  theProton   = G4Proton::Proton();
97  theNeutron  = G4Neutron::Neutron();
98  theDeuteron = G4Deuteron::Deuteron();
99  theAlpha    = G4Alpha::Alpha();
100  thePionPlus = G4PionPlus::PionPlus();
101  thePionMinus= G4PionMinus::PionMinus();
102}
103
104G4HadronElastic::~G4HadronElastic()
105{
106  delete hElastic;
107}
108
109G4VQCrossSection* G4HadronElastic::GetCS()
110{
111  return qCManager;
112}
113
114G4ElasticHadrNucleusHE* G4HadronElastic::GetHElastic()
115{
116  return hElastic;
117}
118
119G4HadFinalState* G4HadronElastic::ApplyYourself(
120                 const G4HadProjectile& aTrack, G4Nucleus& targetNucleus)
121{
122  theParticleChange.Clear();
123
124  const G4HadProjectile* aParticle = &aTrack;
125  G4double ekin = aParticle->GetKineticEnergy();
126  if(ekin <= lowestEnergyLimit) {
127    theParticleChange.SetEnergyChange(ekin);
128    theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
129    return &theParticleChange;
130  }
131
132  G4double aTarget = targetNucleus.GetN();
133  G4double zTarget = targetNucleus.GetZ();
134
135  G4double plab = aParticle->GetTotalMomentum();
136  if (verboseLevel >1) {
137    G4cout << "G4HadronElastic::DoIt: Incident particle plab=" 
138           << plab/GeV << " GeV/c " 
139           << " ekin(MeV) = " << ekin/MeV << "  " 
140           << aParticle->GetDefinition()->GetParticleName() << G4endl;
141  }
142  // Scattered particle referred to axis of incident particle
143  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
144  G4double m1 = theParticle->GetPDGMass();
145
146  G4int Z = static_cast<G4int>(zTarget+0.5);
147  G4int A = static_cast<G4int>(aTarget+0.5);
148  G4int N = A - Z;
149  G4int projPDG = theParticle->GetPDGEncoding();
150  if (verboseLevel>1) 
151    G4cout << "G4HadronElastic for " << theParticle->GetParticleName()
152           << " PDGcode= " << projPDG << " on nucleus Z= " << Z
153           << " A= " << A << " N= " << N
154           << G4endl;
155
156  G4ParticleDefinition * theDef = 0;
157
158  if(Z == 1 && A == 1)       theDef = theProton;
159  else if (Z == 1 && A == 2) theDef = theDeuteron;
160  else if (Z == 1 && A == 3) theDef = G4Triton::Triton();
161  else if (Z == 2 && A == 3) theDef = G4He3::He3();
162  else if (Z == 2 && A == 4) theDef = theAlpha;
163  else theDef = G4ParticleTable::GetParticleTable()->FindIon(Z,A,0,Z);
164 
165  G4double m2 = theDef->GetPDGMass();
166  G4LorentzVector lv1 = aParticle->Get4Momentum();
167  G4LorentzVector lv(0.0,0.0,0.0,m2);   
168  lv += lv1;
169
170  G4ThreeVector bst = lv.boostVector();
171  lv1.boost(-bst);
172
173  G4ThreeVector p1 = lv1.vect();
174  G4double ptot = p1.mag();
175  G4double tmax = 4.0*ptot*ptot;
176  G4double t = 0.0;
177
178  // Choose generator
179  G4ElasticGenerator gtype = fLElastic;
180
181  // Q-elastic for p,n scattering on H and He
182  if (theParticle == theProton || theParticle == theNeutron)
183    //     && Z <= 2 && ekin >= lowEnergyLimitQ) 
184    gtype = fQElastic;
185
186  else {
187    // S-wave for very low energy
188    if(plab < plabLowLimit) gtype = fSWave;
189    // HE-elastic for energetic projectile mesons
190    else if(ekin >= lowEnergyLimitHE && theParticle->GetBaryonNumber() == 0) 
191      gtype = fHElastic;
192  }
193
194  //
195  // Sample t
196  //
197  if(gtype == fQElastic) {
198    if (verboseLevel >1) {
199      G4cout << "G4HadronElastic: Z= " << Z << " N= " 
200             << N << " pdg= " <<  projPDG
201             << " mom(GeV)= " << plab/GeV << "  " << qCManager << G4endl; 
202    }
203    if(Z == 1 && N == 2) N = 1;
204    else if(Z == 2 && N == 1) N = 2;
205    G4double cs = qCManager->GetCrossSection(false,plab,Z,N,projPDG);
206
207    // check if cross section is reasonable
208    if(cs > 0.0) t = qCManager->GetExchangeT(Z,N,projPDG);
209    else if(plab > plabLowLimit) gtype = fLElastic;
210    else gtype = fSWave;
211  }
212
213  if(gtype == fLElastic) {
214    t = GeV*GeV*SampleT(ptot,m1,m2,aTarget);
215  }
216
217  // use mean atomic number
218  if(gtype == fHElastic) {
219    t = hElastic->SampleT(theParticle,plab,Z,A);
220  }
221
222  // NaN finder
223  if(!(t < 0.0 || t >= 0.0)) {
224    if (verboseLevel > 0) {
225      G4cout << "G4HadronElastic:WARNING: Z= " << Z << " N= " 
226             << N << " pdg= " <<  projPDG
227             << " mom(GeV)= " << plab/GeV
228             << " the model type " << gtype;
229      if(gtype ==  fQElastic) G4cout << " CHIPS ";
230      else if(gtype ==  fLElastic) G4cout << " LElastic ";
231      else if(gtype ==  fHElastic) G4cout << " HElastic ";
232      G4cout << " S-wave will be sampled" 
233             << G4endl; 
234    }
235    t = 0.0;
236  }
237
238  if(gtype == fSWave) t = G4UniformRand()*tmax;
239
240  if(verboseLevel>1)
241    G4cout <<"type= " << gtype <<" t= " << t << " tmax= " << tmax
242           << " ptot= " << ptot << G4endl;
243
244  // Sampling in CM system
245  G4double phi  = G4UniformRand()*twopi;
246  G4double cost = 1. - 2.0*t/tmax;
247  G4double sint;
248
249  if( cost >= 1.0 || cost < -1 ) 
250  {
251    cost = 1.0;
252    sint = 0.0;
253  }
254  else 
255  {
256    sint = std::sqrt((1.0-cost)*(1.0+cost));
257  }   
258  if (verboseLevel>1) 
259    G4cout << "cos(t)=" << cost << " std::sin(t)=" << sint << G4endl;
260
261  G4ThreeVector v1(sint*std::cos(phi),sint*std::sin(phi),cost);
262  v1 *= ptot;
263  G4LorentzVector nlv1(v1.x(),v1.y(),v1.z(),std::sqrt(ptot*ptot + m1*m1));
264
265  nlv1.boost(bst); 
266
267  G4double eFinal = nlv1.e() - m1;
268  if (verboseLevel > 1) {
269    G4cout << "Scattered: "
270           << nlv1<<" m= " << m1 << " ekin(MeV)= " << eFinal
271           << " Proj: 4-mom " << lv1
272           <<G4endl;
273  }
274  if(eFinal <= lowestEnergyLimit) {
275    if(eFinal < 0.0) {
276      G4cout << "G4HadronElastic WARNING ekin= " << eFinal
277             << " after scattering of " 
278             << aParticle->GetDefinition()->GetParticleName()
279             << " p(GeV/c)= " << plab
280             << " on " << theDef->GetParticleName()
281             << G4endl;
282    }
283    theParticleChange.SetEnergyChange(0.0);
284    nlv1 = G4LorentzVector(0.0,0.0,0.0,m1);
285
286  } else {
287    theParticleChange.SetMomentumChange(nlv1.vect().unit());
288    theParticleChange.SetEnergyChange(eFinal);
289  } 
290
291  G4LorentzVector nlv0 = lv - nlv1;
292  G4double erec =  nlv0.e() - m2;
293  if (verboseLevel > 1) {
294    G4cout << "Recoil: "
295           << nlv0<<" m= " << m2 << " ekin(MeV)= " << erec
296           <<G4endl;
297  }
298  if(erec > lowEnergyRecoilLimit) {
299    G4DynamicParticle * aSec = new G4DynamicParticle(theDef, nlv0);
300    theParticleChange.AddSecondary(aSec);
301  } else {
302    if(erec < 0.0) erec = 0.0;
303    theParticleChange.SetLocalEnergyDeposit(erec);
304  }
305
306  return &theParticleChange;
307}
308
309G4double
310G4HadronElastic::SampleT(G4double, G4double, G4double, G4double atno2)
311{
312  // G4cout << "Entering elastic scattering 2"<<G4endl;
313  // Compute the direction of elastic scattering.
314  // It is planned to replace this code with a method based on
315  // parameterized functions and a Monte Carlo method to invert the CDF.
316
317  G4double ran = G4UniformRand();
318  G4double aa, bb, cc, dd, rr;
319  if (atno2 <= 62.) {
320    aa = std::pow(atno2, 1.63);
321    bb = 14.5*std::pow(atno2, 0.66);
322    cc = 1.4*std::pow(atno2, 0.33);
323    dd = 10.;
324  } else {
325    aa = std::pow(atno2, 1.33);
326    bb = 60.*std::pow(atno2, 0.33);
327    cc = 0.4*std::pow(atno2, 0.40);
328    dd = 10.;
329  }
330  aa = aa/bb;
331  cc = cc/dd;
332  rr = (aa + cc)*ran;
333  if (verboseLevel > 1) {
334    G4cout << "DoIt: aa,bb,cc,dd,rr" << G4endl;
335    G4cout << aa << " " << bb << " " << cc << " " << dd << " " << rr << G4endl;
336  }
337  G4double t1 = -std::log(ran)/bb;
338  G4double t2 = -std::log(ran)/dd;
339  if (verboseLevel > 1) {
340    G4cout << "t1,Fctcos " << t1 << " " << Fctcos(t1, aa, bb, cc, dd, rr) << G4endl;
341    G4cout << "t2,Fctcos " << t2 << " " << Fctcos(t2, aa, bb, cc, dd, rr) << G4endl;
342  }
343  G4double eps = 0.001;
344  G4int ind1 = 10;
345  G4double t = 0.0;
346  G4int ier1;
347  ier1 = Rtmi(&t, t1, t2, eps, ind1,
348              aa, bb, cc, dd, rr);
349  if (verboseLevel > 1) {
350    G4cout << "From Rtmi, ier1=" << ier1 << G4endl;
351    G4cout << "t, Fctcos " << t << " " << Fctcos(t, aa, bb, cc, dd, rr) << G4endl;
352  }
353  if (ier1 != 0) t = 0.25*(3.*t1 + t2);
354  if (verboseLevel > 1) {
355      G4cout << "t, Fctcos " << t << " " << Fctcos(t, aa, bb, cc, dd, rr) << 
356              G4endl;
357  }
358  return t;
359}
360
361// The following is a "translation" of a root-finding routine
362// from GEANT3.21/GHEISHA.  Some of the labelled block structure has
363// been retained for clarity.  This routine will not be needed after
364// the planned revisions to DoIt().
365
366G4int
367G4HadronElastic::Rtmi(G4double* x, G4double xli, G4double xri, G4double eps, 
368                      G4int iend, 
369                      G4double aa, G4double bb, G4double cc, G4double dd, 
370                      G4double rr)
371{
372   G4int ier = 0;
373   G4double xl = xli;
374   G4double xr = xri;
375   *x = xl;
376   G4double tol = *x;
377   G4double f = Fctcos(tol, aa, bb, cc, dd, rr);
378   if (f == 0.) return ier;
379   G4double fl, fr;
380   fl = f;
381   *x = xr;
382   tol = *x;
383   f = Fctcos(tol, aa, bb, cc, dd, rr);
384   if (f == 0.) return ier;
385   fr = f;
386
387// Error return in case of wrong input data
388   if (fl*fr >= 0.) {
389      ier = 2;
390      return ier;
391   }
392
393// Basic assumption fl*fr less than 0 is satisfied.
394// Generate tolerance for function values.
395   G4int i = 0;
396   G4double tolf = 100.*eps;
397
398// Start iteration loop
399label4:
400   i++;
401
402// Start bisection loop
403   for (G4int k = 1; k <= iend; k++) {
404      *x = 0.5*(xl + xr);
405      tol = *x;
406      f = Fctcos(tol, aa, bb, cc, dd, rr);
407      if (f == 0.) return 0;
408      if (f*fr < 0.) {      // Interchange xl and xr in order to get the
409         tol = xl;          // same Sign in f and fr
410         xl = xr;
411         xr = tol;
412         tol = fl;
413         fl = fr;
414         fr = tol;
415      }
416      tol = f - fl;
417      G4double a = f*tol;
418      a = a + a;
419      if (a < fr*(fr - fl) && i <= iend) goto label17;
420      xr = *x;
421      fr = f;
422
423// Test on satisfactory accuracy in bisection loop
424      tol = eps;
425      a = std::abs(xr);
426      if (a > 1.) tol = tol*a;
427      if (std::abs(xr - xl) <= tol && std::abs(fr - fl) <= tolf) goto label14;
428   }
429// End of bisection loop
430
431// No convergence after iend iteration steps followed by iend
432// successive steps of bisection or steadily increasing function
433// values at right bounds.  Error return.
434   ier = 1;
435
436label14:
437   if (std::abs(fr) > std::abs(fl)) {
438      *x = xl;
439      f = fl;
440   }
441   return ier;
442
443// Computation of iterated x-value by inverse parabolic interp
444label17:
445   G4double a = fr - f;
446   G4double dx = (*x - xl)*fl*(1. + f*(a - tol)/(a*(fr - fl)))/tol;
447   G4double xm = *x;
448   G4double fm = f;
449   *x = xl - dx;
450   tol = *x;
451   f = Fctcos(tol, aa, bb, cc, dd, rr);
452   if (f == 0.) return ier;
453
454// Test on satisfactory accuracy in iteration loop
455   tol = eps;
456   a = std::abs(*x);
457   if (a > 1) tol = tol*a;
458   if (std::abs(dx) <= tol && std::abs(f) <= tolf) return ier;
459
460// Preparation of next bisection loop
461   if (f*fl < 0.) {
462      xr = *x;
463      fr = f;
464   }
465   else {
466      xl = *x;
467      fl = f;
468      xr = xm;
469      fr = fm;
470   }
471   goto label4;
472}
473
474// Test function for root-finder
475G4double
476G4HadronElastic::Fctcos(G4double t, 
477                        G4double aa, G4double bb, G4double cc, G4double dd, 
478                        G4double rr)
479{
480   const G4double expxl = -82.;
481   const G4double expxu = 82.;
482
483   G4double test1 = -bb*t;
484   if (test1 > expxu) test1 = expxu;
485   if (test1 < expxl) test1 = expxl;
486
487   G4double test2 = -dd*t;
488   if (test2 > expxu) test2 = expxu;
489   if (test2 < expxl) test2 = expxl;
490
491   return aa*std::exp(test1) + cc*std::exp(test2) - rr;
492}
493
494
Note: See TracBrowser for help on using the repository browser.