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.61 2008/08/05 07:37:39 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-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 | |
---|
81 | G4HadronElastic::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 | nnans = 0; |
---|
104 | npos = 0; |
---|
105 | nneg = 0; |
---|
106 | neneg = 0; |
---|
107 | } |
---|
108 | |
---|
109 | G4HadronElastic::~G4HadronElastic() |
---|
110 | { |
---|
111 | delete hElastic; |
---|
112 | if( (nnans + npos + nneg + neneg) > 0 ) { |
---|
113 | G4cout << "### G4HadronElastic destructor Warnings: "; |
---|
114 | if(nnans > 0) G4cout << "### N(nans) = " << nnans; |
---|
115 | if(npos > 0) G4cout << "### N(cost > 1)= " << npos; |
---|
116 | if(nneg > 0) G4cout << "### N(cost <-1)= " << nneg; |
---|
117 | if(neneg > 0) G4cout << "### N(E < 0)= " << neneg; |
---|
118 | G4cout << "###" << G4endl; |
---|
119 | } |
---|
120 | } |
---|
121 | |
---|
122 | G4VQCrossSection* G4HadronElastic::GetCS() |
---|
123 | { |
---|
124 | return qCManager; |
---|
125 | } |
---|
126 | |
---|
127 | G4ElasticHadrNucleusHE* G4HadronElastic::GetHElastic() |
---|
128 | { |
---|
129 | return hElastic; |
---|
130 | } |
---|
131 | |
---|
132 | G4HadFinalState* 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= " |
---|
213 | << N << " pdg= " << projPDG |
---|
214 | << " mom(GeV)= " << plab/GeV << " " << qCManager << G4endl; |
---|
215 | } |
---|
216 | if(Z == 1 && N == 2) N = 1; |
---|
217 | else if(Z == 2 && N == 1) N = 2; |
---|
218 | G4double cs = qCManager->GetCrossSection(false,plab,Z,N,projPDG); |
---|
219 | |
---|
220 | // check if cross section is reasonable |
---|
221 | if(cs > 0.0) t = qCManager->GetExchangeT(Z,N,projPDG); |
---|
222 | else if(plab > plabLowLimit) gtype = fLElastic; |
---|
223 | else gtype = fSWave; |
---|
224 | } |
---|
225 | |
---|
226 | if(gtype == fLElastic) { |
---|
227 | G4double g2 = GeV*GeV; |
---|
228 | t = g2*SampleT(tmax/g2,m1,m2,aTarget); |
---|
229 | } |
---|
230 | |
---|
231 | // use mean atomic number |
---|
232 | if(gtype == fHElastic) { |
---|
233 | t = hElastic->SampleT(theParticle,plab,Z,A); |
---|
234 | } |
---|
235 | |
---|
236 | // NaN finder |
---|
237 | if(!(t < 0.0 || t >= 0.0)) { |
---|
238 | if (verboseLevel > 0) { |
---|
239 | G4cout << "G4HadronElastic:WARNING: Z= " << Z << " N= " |
---|
240 | << N << " pdg= " << projPDG |
---|
241 | << " mom(GeV)= " << plab/GeV |
---|
242 | << " the model type " << gtype; |
---|
243 | if(gtype == fQElastic) G4cout << " CHIPS "; |
---|
244 | else if(gtype == fLElastic) G4cout << " LElastic "; |
---|
245 | else if(gtype == fHElastic) G4cout << " HElastic "; |
---|
246 | G4cout << " S-wave will be sampled" |
---|
247 | << G4endl; |
---|
248 | } |
---|
249 | t = 0.0; |
---|
250 | nnans++; |
---|
251 | } |
---|
252 | |
---|
253 | if(gtype == fSWave) t = G4UniformRand()*tmax; |
---|
254 | |
---|
255 | if(verboseLevel>1) { |
---|
256 | G4cout <<"type= " << gtype <<" t= " << t << " tmax= " << tmax |
---|
257 | << " ptot= " << ptot << G4endl; |
---|
258 | } |
---|
259 | // Sampling in CM system |
---|
260 | G4double phi = G4UniformRand()*twopi; |
---|
261 | G4double cost = 1. - 2.0*t/tmax; |
---|
262 | G4double sint; |
---|
263 | |
---|
264 | // problem in sampling |
---|
265 | if(cost >= 1.0) { |
---|
266 | cost = 1.0; |
---|
267 | sint = 0.0; |
---|
268 | npos++; |
---|
269 | } else if(cost < -1 ) { |
---|
270 | /* |
---|
271 | G4cout << "G4HadronElastic:WARNING: Z= " << Z << " N= " |
---|
272 | << N << " " << aParticle->GetDefinition()->GetParticleName() |
---|
273 | << " mom(GeV)= " << plab/GeV |
---|
274 | << " the model type " << gtype; |
---|
275 | if(gtype == fQElastic) G4cout << " CHIPS "; |
---|
276 | else if(gtype == fLElastic) G4cout << " LElastic "; |
---|
277 | else if(gtype == fHElastic) G4cout << " HElastic "; |
---|
278 | G4cout << " cost= " << cost |
---|
279 | << G4endl; |
---|
280 | */ |
---|
281 | cost = 1.0; |
---|
282 | sint = 0.0; |
---|
283 | nneg++; |
---|
284 | |
---|
285 | // normal situation |
---|
286 | } else { |
---|
287 | sint = std::sqrt((1.0-cost)*(1.0+cost)); |
---|
288 | } |
---|
289 | if (verboseLevel>1) { |
---|
290 | G4cout << "cos(t)=" << cost << " std::sin(t)=" << sint << G4endl; |
---|
291 | } |
---|
292 | G4ThreeVector v1(sint*std::cos(phi),sint*std::sin(phi),cost); |
---|
293 | v1 *= ptot; |
---|
294 | G4LorentzVector nlv1(v1.x(),v1.y(),v1.z(),std::sqrt(ptot*ptot + m1*m1)); |
---|
295 | |
---|
296 | nlv1.boost(bst); |
---|
297 | |
---|
298 | G4double eFinal = nlv1.e() - m1; |
---|
299 | if (verboseLevel > 1) { |
---|
300 | G4cout << "Scattered: " |
---|
301 | << nlv1<<" m= " << m1 << " ekin(MeV)= " << eFinal |
---|
302 | << " Proj: 4-mom " << lv1 |
---|
303 | <<G4endl; |
---|
304 | } |
---|
305 | if(eFinal <= lowestEnergyLimit) { |
---|
306 | if(eFinal < 0.0 && verboseLevel > 0) { |
---|
307 | neneg++; |
---|
308 | G4cout << "G4HadronElastic WARNING ekin= " << eFinal |
---|
309 | << " after scattering of " |
---|
310 | << aParticle->GetDefinition()->GetParticleName() |
---|
311 | << " p(GeV/c)= " << plab |
---|
312 | << " on " << theDef->GetParticleName() |
---|
313 | << G4endl; |
---|
314 | } |
---|
315 | theParticleChange.SetEnergyChange(0.0); |
---|
316 | nlv1 = G4LorentzVector(0.0,0.0,0.0,m1); |
---|
317 | |
---|
318 | } else { |
---|
319 | theParticleChange.SetMomentumChange(nlv1.vect().unit()); |
---|
320 | theParticleChange.SetEnergyChange(eFinal); |
---|
321 | } |
---|
322 | |
---|
323 | G4LorentzVector nlv0 = lv - nlv1; |
---|
324 | G4double erec = nlv0.e() - m2; |
---|
325 | if (verboseLevel > 1) { |
---|
326 | G4cout << "Recoil: " |
---|
327 | << nlv0<<" m= " << m2 << " ekin(MeV)= " << erec |
---|
328 | <<G4endl; |
---|
329 | } |
---|
330 | if(erec > lowEnergyRecoilLimit) { |
---|
331 | G4DynamicParticle * aSec = new G4DynamicParticle(theDef, nlv0); |
---|
332 | theParticleChange.AddSecondary(aSec); |
---|
333 | } else { |
---|
334 | if(erec < 0.0) erec = 0.0; |
---|
335 | theParticleChange.SetLocalEnergyDeposit(erec); |
---|
336 | } |
---|
337 | |
---|
338 | return &theParticleChange; |
---|
339 | } |
---|
340 | |
---|
341 | G4double |
---|
342 | G4HadronElastic::SampleT(G4double tmax, G4double, G4double, G4double atno2) |
---|
343 | { |
---|
344 | // G4cout << "Entering elastic scattering 2"<<G4endl; |
---|
345 | // Compute the direction of elastic scattering. |
---|
346 | // It is planned to replace this code with a method based on |
---|
347 | // parameterized functions and a Monte Carlo method to invert the CDF. |
---|
348 | |
---|
349 | // G4double ran = G4UniformRand(); |
---|
350 | G4double aa, bb, cc, dd, rr; |
---|
351 | if (atno2 <= 62.) { |
---|
352 | aa = std::pow(atno2, 1.63); |
---|
353 | bb = 14.5*std::pow(atno2, 0.66); |
---|
354 | cc = 1.4*std::pow(atno2, 0.33); |
---|
355 | dd = 10.; |
---|
356 | } else { |
---|
357 | aa = std::pow(atno2, 1.33); |
---|
358 | bb = 60.*std::pow(atno2, 0.33); |
---|
359 | cc = 0.4*std::pow(atno2, 0.40); |
---|
360 | dd = 10.; |
---|
361 | } |
---|
362 | aa = aa/bb; |
---|
363 | cc = cc/dd; |
---|
364 | G4double ran, t1, t2; |
---|
365 | do { |
---|
366 | ran = G4UniformRand(); |
---|
367 | t1 = -std::log(ran)/bb; |
---|
368 | t2 = -std::log(ran)/dd; |
---|
369 | } while(t1 > tmax || t2 > tmax); |
---|
370 | |
---|
371 | rr = (aa + cc)*ran; |
---|
372 | |
---|
373 | if (verboseLevel > 1) { |
---|
374 | G4cout << "DoIt: aa,bb,cc,dd,rr" << G4endl; |
---|
375 | G4cout << aa << " " << bb << " " << cc << " " << dd << " " << rr << G4endl; |
---|
376 | G4cout << "t1,Fctcos " << t1 << " " << Fctcos(t1, aa, bb, cc, dd, rr) << G4endl; |
---|
377 | G4cout << "t2,Fctcos " << t2 << " " << Fctcos(t2, aa, bb, cc, dd, rr) << G4endl; |
---|
378 | } |
---|
379 | G4double eps = 0.001; |
---|
380 | G4int ind1 = 10; |
---|
381 | G4double t = 0.0; |
---|
382 | G4int ier1; |
---|
383 | ier1 = Rtmi(&t, t1, t2, eps, ind1, |
---|
384 | aa, bb, cc, dd, rr); |
---|
385 | if (verboseLevel > 1) { |
---|
386 | G4cout << "From Rtmi, ier1=" << ier1 << " t= " << t << G4endl; |
---|
387 | G4cout << "t, Fctcos " << t << " " << Fctcos(t, aa, bb, cc, dd, rr) << G4endl; |
---|
388 | } |
---|
389 | if (ier1 != 0) t = 0.25*(3.*t1 + t2); |
---|
390 | if (verboseLevel > 1) { |
---|
391 | G4cout << "t, Fctcos " << t << " " << Fctcos(t, aa, bb, cc, dd, rr) << |
---|
392 | G4endl; |
---|
393 | } |
---|
394 | return t; |
---|
395 | } |
---|
396 | |
---|
397 | // The following is a "translation" of a root-finding routine |
---|
398 | // from GEANT3.21/GHEISHA. Some of the labelled block structure has |
---|
399 | // been retained for clarity. This routine will not be needed after |
---|
400 | // the planned revisions to DoIt(). |
---|
401 | |
---|
402 | G4int |
---|
403 | G4HadronElastic::Rtmi(G4double* x, G4double xli, G4double xri, G4double eps, |
---|
404 | G4int iend, |
---|
405 | G4double aa, G4double bb, G4double cc, G4double dd, |
---|
406 | G4double rr) |
---|
407 | { |
---|
408 | G4int ier = 0; |
---|
409 | G4double xl = xli; |
---|
410 | G4double xr = xri; |
---|
411 | *x = xl; |
---|
412 | G4double tol = *x; |
---|
413 | G4double f = Fctcos(tol, aa, bb, cc, dd, rr); |
---|
414 | if (f == 0.) return ier; |
---|
415 | G4double fl, fr; |
---|
416 | fl = f; |
---|
417 | *x = xr; |
---|
418 | tol = *x; |
---|
419 | f = Fctcos(tol, aa, bb, cc, dd, rr); |
---|
420 | if (f == 0.) return ier; |
---|
421 | fr = f; |
---|
422 | |
---|
423 | // Error return in case of wrong input data |
---|
424 | if (fl*fr >= 0.) { |
---|
425 | ier = 2; |
---|
426 | return ier; |
---|
427 | } |
---|
428 | |
---|
429 | // Basic assumption fl*fr less than 0 is satisfied. |
---|
430 | // Generate tolerance for function values. |
---|
431 | G4int i = 0; |
---|
432 | G4double tolf = 100.*eps; |
---|
433 | |
---|
434 | // Start iteration loop |
---|
435 | label4: |
---|
436 | i++; |
---|
437 | |
---|
438 | // Start bisection loop |
---|
439 | for (G4int k = 1; k <= iend; k++) { |
---|
440 | *x = 0.5*(xl + xr); |
---|
441 | tol = *x; |
---|
442 | f = Fctcos(tol, aa, bb, cc, dd, rr); |
---|
443 | if (f == 0.) return 0; |
---|
444 | if (f*fr < 0.) { // Interchange xl and xr in order to get the |
---|
445 | tol = xl; // same Sign in f and fr |
---|
446 | xl = xr; |
---|
447 | xr = tol; |
---|
448 | tol = fl; |
---|
449 | fl = fr; |
---|
450 | fr = tol; |
---|
451 | } |
---|
452 | tol = f - fl; |
---|
453 | G4double a = f*tol; |
---|
454 | a = a + a; |
---|
455 | if (a < fr*(fr - fl) && i <= iend) goto label17; |
---|
456 | xr = *x; |
---|
457 | fr = f; |
---|
458 | |
---|
459 | // Test on satisfactory accuracy in bisection loop |
---|
460 | tol = eps; |
---|
461 | a = std::abs(xr); |
---|
462 | if (a > 1.) tol = tol*a; |
---|
463 | if (std::abs(xr - xl) <= tol && std::abs(fr - fl) <= tolf) goto label14; |
---|
464 | } |
---|
465 | // End of bisection loop |
---|
466 | |
---|
467 | // No convergence after iend iteration steps followed by iend |
---|
468 | // successive steps of bisection or steadily increasing function |
---|
469 | // values at right bounds. Error return. |
---|
470 | ier = 1; |
---|
471 | |
---|
472 | label14: |
---|
473 | if (std::abs(fr) > std::abs(fl)) { |
---|
474 | *x = xl; |
---|
475 | f = fl; |
---|
476 | } |
---|
477 | return ier; |
---|
478 | |
---|
479 | // Computation of iterated x-value by inverse parabolic interp |
---|
480 | label17: |
---|
481 | G4double a = fr - f; |
---|
482 | G4double dx = (*x - xl)*fl*(1. + f*(a - tol)/(a*(fr - fl)))/tol; |
---|
483 | G4double xm = *x; |
---|
484 | G4double fm = f; |
---|
485 | *x = xl - dx; |
---|
486 | tol = *x; |
---|
487 | f = Fctcos(tol, aa, bb, cc, dd, rr); |
---|
488 | if (f == 0.) return ier; |
---|
489 | |
---|
490 | // Test on satisfactory accuracy in iteration loop |
---|
491 | tol = eps; |
---|
492 | a = std::abs(*x); |
---|
493 | if (a > 1) tol = tol*a; |
---|
494 | if (std::abs(dx) <= tol && std::abs(f) <= tolf) return ier; |
---|
495 | |
---|
496 | // Preparation of next bisection loop |
---|
497 | if (f*fl < 0.) { |
---|
498 | xr = *x; |
---|
499 | fr = f; |
---|
500 | } |
---|
501 | else { |
---|
502 | xl = *x; |
---|
503 | fl = f; |
---|
504 | xr = xm; |
---|
505 | fr = fm; |
---|
506 | } |
---|
507 | goto label4; |
---|
508 | } |
---|
509 | |
---|
510 | // Test function for root-finder |
---|
511 | G4double |
---|
512 | G4HadronElastic::Fctcos(G4double t, |
---|
513 | G4double aa, G4double bb, G4double cc, G4double dd, |
---|
514 | G4double rr) |
---|
515 | { |
---|
516 | const G4double expxl = -82.; |
---|
517 | const G4double expxu = 82.; |
---|
518 | |
---|
519 | G4double test1 = -bb*t; |
---|
520 | if (test1 > expxu) test1 = expxu; |
---|
521 | if (test1 < expxl) test1 = expxl; |
---|
522 | |
---|
523 | G4double test2 = -dd*t; |
---|
524 | if (test2 > expxu) test2 = expxu; |
---|
525 | if (test2 < expxl) test2 = expxl; |
---|
526 | |
---|
527 | return aa*std::exp(test1) + cc*std::exp(test2) - rr; |
---|
528 | } |
---|
529 | |
---|
530 | |
---|