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

Last change on this file since 1200 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

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