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

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

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

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