Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4HadronNucleonXsc.cc
Go to the documentation of this file.
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// 14.03.07 V. Grichine - first implementation
27//
28
29#include "G4HadronNucleonXsc.hh"
30
32#include "G4SystemOfUnits.hh"
33#include "G4ParticleTable.hh"
34#include "G4IonTable.hh"
36#include "G4HadTmpUtil.hh"
37
38
40: fUpperLimit( 10000 * GeV ),
41 fLowerLimit( 0.03 * MeV ),
42 fTotalXsc(0.0), fElasticXsc(0.0), fInelasticXsc(0.0), fHadronNucleonXsc(0.0)
43{
44 theGamma = G4Gamma::Gamma();
45 theProton = G4Proton::Proton();
46 theNeutron = G4Neutron::Neutron();
47 theAProton = G4AntiProton::AntiProton();
48 theANeutron = G4AntiNeutron::AntiNeutron();
49 thePiPlus = G4PionPlus::PionPlus();
50 thePiMinus = G4PionMinus::PionMinus();
51 thePiZero = G4PionZero::PionZero();
52 theKPlus = G4KaonPlus::KaonPlus();
53 theKMinus = G4KaonMinus::KaonMinus();
56 theL = G4Lambda::Lambda();
57 theAntiL = G4AntiLambda::AntiLambda();
58 theSPlus = G4SigmaPlus::SigmaPlus();
60 theSMinus = G4SigmaMinus::SigmaMinus();
62 theS0 = G4SigmaZero::SigmaZero();
64 theXiMinus = G4XiMinus::XiMinus();
65 theXi0 = G4XiZero::XiZero();
66 theAXiMinus = G4AntiXiMinus::AntiXiMinus();
67 theAXi0 = G4AntiXiZero::AntiXiZero();
68 theOmega = G4OmegaMinus::OmegaMinus();
70 theD = G4Deuteron::Deuteron();
71 theT = G4Triton::Triton();
72 theA = G4Alpha::Alpha();
73 theHe3 = G4He3::He3();
74
76}
77
78
80{}
81
82void G4HadronNucleonXsc::CrossSectionDescription(std::ostream& outFile) const
83{
84 outFile << "G4HadronNucleonXsc calculates the total, inelastic and elastic\n"
85 << "hadron-nucleon cross sections using several different\n"
86 << "parameterizations within the Glauber-Gribov approach. It is\n"
87 << "valid for all incident gammas and long-lived hadrons at\n"
88 << "energies above 30 keV. This is a cross section component which\n"
89 << "is to be used to build a cross section data set.\n";
90}
91
92G4bool
94 const G4Element* anElement)
95{
96 G4int Z = G4lrint(anElement->GetZ());
97 G4int A = G4lrint(anElement->GetN());
98 return IsIsoApplicable(aDP, Z, A);
99}
100
101////////////////////////////////////////////////////////////////////////////////////////
102//
103
104G4bool
106 G4int Z, G4int)
107{
108 G4bool applicable = false;
109 // G4int baryonNumber = aDP->GetDefinition()->GetBaryonNumber();
110 G4double kineticEnergy = aDP->GetKineticEnergy();
111
112 const G4ParticleDefinition* theParticle = aDP->GetDefinition();
113
114 if ( ( kineticEnergy >= fLowerLimit &&
115 Z > 1 && // >= He
116 ( theParticle == theAProton ||
117 theParticle == theGamma ||
118 theParticle == theKPlus ||
119 theParticle == theKMinus ||
120 theParticle == theSMinus) ) ||
121
122 ( kineticEnergy >= 0.1*fLowerLimit &&
123 Z > 1 && // >= He
124 ( theParticle == theProton ||
125 theParticle == theNeutron ||
126 theParticle == thePiPlus ||
127 theParticle == thePiMinus ) ) ) applicable = true;
128
129 return applicable;
130}
131
132
133/////////////////////////////////////////////////////////////////////////////////////
134//
135// Returns hadron-nucleon Xsc according to differnt parametrisations:
136// [2] E. Levin, hep-ph/9710546
137// [3] U. Dersch, et al, hep-ex/9910052
138// [4] M.J. Longo, et al, Phys.Rev.Lett. 33 (1974) 725
139
142 const G4ParticleDefinition* nucleon )
143{
144 G4double xsection;
145
146
147 G4double targ_mass = 0.939*GeV; // ~mean neutron and proton ???
148
149 G4double proj_mass = aParticle->GetMass();
150 G4double proj_momentum = aParticle->GetMomentum().mag();
151 G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
152
153 sMand /= GeV*GeV; // in GeV for parametrisation
154 proj_momentum /= GeV;
155
156 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
157
158 G4bool pORn = (nucleon == theProton || nucleon == theNeutron );
159
160
161 if(theParticle == theGamma && pORn )
162 {
163 xsection = (0.0677*std::pow(sMand,0.0808) + 0.129*std::pow(sMand,-0.4525));
164 }
165 else if(theParticle == theNeutron && pORn ) // as proton ???
166 {
167 xsection = (21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
168 }
169 else if(theParticle == theProton && pORn )
170 {
171 xsection = (21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
172
173 // xsection = At*( 49.51*std::pow(sMand,-0.097) + 0.314*std::log(sMand)*std::log(sMand) );
174 // xsection = At*( 38.4 + 0.85*std::abs(std::pow(log(sMand),1.47)) );
175 }
176 else if(theParticle == theAProton && pORn )
177 {
178 xsection = ( 21.70*std::pow(sMand,0.0808) + 98.39*std::pow(sMand,-0.4525));
179 }
180 else if(theParticle == thePiPlus && pORn )
181 {
182 xsection = (13.63*std::pow(sMand,0.0808) + 27.56*std::pow(sMand,-0.4525));
183 }
184 else if(theParticle == thePiMinus && pORn )
185 {
186 // xsection = At*( 55.2*std::pow(sMand,-0.255) + 0.346*std::log(sMand)*std::log(sMand) );
187 xsection = (13.63*std::pow(sMand,0.0808) + 36.02*std::pow(sMand,-0.4525));
188 }
189 else if(theParticle == theKPlus && pORn )
190 {
191 xsection = (11.82*std::pow(sMand,0.0808) + 8.15*std::pow(sMand,-0.4525));
192 }
193 else if(theParticle == theKMinus && pORn )
194 {
195 xsection = (11.82*std::pow(sMand,0.0808) + 26.36*std::pow(sMand,-0.4525));
196 }
197 else // as proton ???
198 {
199 xsection = (21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
200 }
201 xsection *= millibarn;
202
203 fTotalXsc = xsection;
204 fInelasticXsc = 0.83*xsection;
205 fElasticXsc = fTotalXsc - fInelasticXsc;
206 if (fElasticXsc < 0.)fElasticXsc = 0.;
207
208 return xsection;
209}
210
211
212/////////////////////////////////////////////////////////////////////////////////////
213//
214// Returns hadron-nucleon Xsc according to PDG parametrisation (2005):
215// http://pdg.lbl.gov/2006/reviews/hadronicrpp.pdf
216// At = number of nucleons, Zt = number of protons
217
220 const G4ParticleDefinition* nucleon )
221{
222 G4double xsection(0);
223 G4int Zt=1, Nt=1, At=1;
224
225 G4double targ_mass = 0.939*GeV; // ~mean neutron and proton ???
226
227 G4double proj_mass = aParticle->GetMass();
228 G4double proj_momentum = aParticle->GetMomentum().mag();
229
230 G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
231
232 sMand /= GeV*GeV; // in GeV for parametrisation
233
234 // General PDG fit constants
235
236 G4double s0 = 5.38*5.38; // in Gev^2
237 G4double eta1 = 0.458;
238 G4double eta2 = 0.458;
239 G4double B = 0.308;
240
241 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
242
243 G4bool pORn = (nucleon == theProton || nucleon == theNeutron );
244 G4bool proton = (nucleon == theProton);
245 G4bool neutron = (nucleon == theNeutron);
246
247 if(theParticle == theNeutron) // proton-neutron fit
248 {
249 if ( proton )
250 {
251 xsection = Zt*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
252 + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));// on p
253 }
254 if ( neutron )
255 {
256 xsection = Nt*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
257 + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2)); // on n pp for nn
258 }
259 }
260 else if(theParticle == theProton)
261 {
262 if ( proton )
263 {
264 xsection = Zt*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
265 + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
266 }
267 if ( neutron )
268 {
269 xsection = Nt*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
270 + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
271 }
272 }
273 else if(theParticle == theAProton)
274 {
275 if ( proton )
276 {
277 xsection = Zt*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
278 + 42.53*std::pow(sMand,-eta1) + 33.34*std::pow(sMand,-eta2));
279 }
280 if ( neutron )
281 {
282 xsection = Nt*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
283 + 40.15*std::pow(sMand,-eta1) + 30.*std::pow(sMand,-eta2));
284 }
285 }
286 else if(theParticle == thePiPlus && pORn )
287 {
288 xsection = At*( 20.86 + B*std::pow(std::log(sMand/s0),2.)
289 + 19.24*std::pow(sMand,-eta1) - 6.03*std::pow(sMand,-eta2));
290 }
291 else if(theParticle == thePiMinus && pORn )
292 {
293 xsection = At*( 20.86 + B*std::pow(std::log(sMand/s0),2.)
294 + 19.24*std::pow(sMand,-eta1) + 6.03*std::pow(sMand,-eta2));
295 }
296 else if(theParticle == theKPlus)
297 {
298 if ( proton )
299 {
300 xsection = Zt*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
301 + 7.14*std::pow(sMand,-eta1) - 13.45*std::pow(sMand,-eta2));
302 }
303 if ( neutron )
304 {
305 xsection = Nt*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
306 + 5.17*std::pow(sMand,-eta1) - 7.23*std::pow(sMand,-eta2));
307 }
308 }
309 else if(theParticle == theKMinus)
310 {
311 if ( proton )
312 {
313 xsection = Zt*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
314 + 7.14*std::pow(sMand,-eta1) + 13.45*std::pow(sMand,-eta2));
315 }
316 if ( neutron )
317 {
318 xsection = Nt*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
319 + 5.17*std::pow(sMand,-eta1) + 7.23*std::pow(sMand,-eta2) );
320 }
321 }
322 else if(theParticle == theSMinus && pORn )
323 {
324 xsection = At*( 35.20 + B*std::pow(std::log(sMand/s0),2.)
325 - 199.*std::pow(sMand,-eta1) + 264.*std::pow(sMand,-eta2) );
326 }
327 else if(theParticle == theGamma && pORn ) // modify later on
328 {
329 xsection = At*( 0.0 + B*std::pow(std::log(sMand/s0),2.)
330 + 0.032*std::pow(sMand,-eta1) - 0.0*std::pow(sMand,-eta2) );
331
332 }
333 else // as proton ???
334 {
335 if ( proton )
336 {
337 xsection = Zt*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
338 + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2) );
339 }
340 if ( neutron )
341 {
342 xsection = Nt*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
343 + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
344 }
345 }
346 xsection *= millibarn; // parametrised in mb
347
348 fTotalXsc = xsection;
349 fInelasticXsc = 0.75*xsection;
350 fElasticXsc = fTotalXsc - fInelasticXsc;
351 if (fElasticXsc < 0.) fElasticXsc = 0.;
352
353 return xsection;
354}
355
356
357/////////////////////////////////////////////////////////////////////////////////////
358//
359// Returns hadron-nucleon cross-section based on N. Starkov parametrisation of
360// data from mainly http://wwwppds.ihep.su:8001/c5-6A.html database
361
364 const G4ParticleDefinition* nucleon )
365{
366 G4double xsection(0);
367
368 G4double A0, B0;
369 G4double hpXsc(0);
370 G4double hnXsc(0);
371
372
373 G4double tM = 0.939*GeV; // ~mean neutron and proton ???
374
375 G4double pM = aParticle->GetMass();
376 G4double pE = aParticle->GetTotalEnergy();
377 G4double pLab = aParticle->GetMomentum().mag();
378
379 G4double sMand = CalcMandelstamS ( pM , tM , pLab );
380
381 sMand /= GeV*GeV; // in GeV for parametrisation
382 pLab /= GeV;
383 pE /= GeV;
384 pM /= GeV;
385
386 G4double logP = std::log(pLab);
387
388
389 // General PDG fit constants
390
391 G4double s0 = 5.38*5.38; // in Gev^2
392 G4double eta1 = 0.458;
393 G4double eta2 = 0.458;
394 G4double B = 0.308;
395
396
397 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
398
399 G4bool pORn = (nucleon == theProton || nucleon == theNeutron );
400 G4bool proton = (nucleon == theProton);
401 G4bool neutron = (nucleon == theNeutron);
402
403 if( theParticle == theNeutron && pORn )
404 {
405 if( pLab >= 373.)
406 {
407 xsection = GetHadronNucleonXscPDG(aParticle, nucleon)/millibarn;
408
409 fElasticXsc = 6.5 + 0.308*std::pow(std::log(sMand/400.),1.65) + 9.19*std::pow(sMand,-0.458);
410
411 fTotalXsc = xsection;
412
413 }
414 else if( pLab >= 100.)
415 {
416 B0 = 7.5;
417 A0 = 100. - B0*std::log(3.0e7);
418
419 xsection = A0 + B0*std::log(pE) - 11
420 // + 103*std::pow(2*0.93827*pE + pM*pM+0.93827*0.93827,-0.165); // mb
421 + 103*std::pow(sMand,-0.165); // mb
422
423 fElasticXsc = 5.53 + 0.308*std::pow(std::log(sMand/28.9),1.1) + 9.19*std::pow(sMand,-0.458);
424
425 fTotalXsc = xsection;
426 }
427 else if( pLab >= 10.)
428 {
429 B0 = 7.5;
430 A0 = 100. - B0*std::log(3.0e7);
431
432 xsection = A0 + B0*std::log(pE) - 11
433 + 103*std::pow(2*0.93827*pE + pM*pM+
434 0.93827*0.93827,-0.165); // mb
435 fTotalXsc = xsection;
436 fElasticXsc = 6 + 20/( (logP-0.182)*(logP-0.182) + 1.0 );
437 }
438 else // pLab < 10 GeV/c
439 {
440 if( neutron ) // nn to be pp
441 {
442 if( pLab < 0.4 )
443 {
444 hnXsc = 23 + 50*( std::pow( std::log(0.73/pLab), 3.5 ) );
445 fElasticXsc = hnXsc;
446 }
447 else if( pLab < 0.73 )
448 {
449 hnXsc = 23 + 50*( std::pow( std::log(0.73/pLab), 3.5 ) );
450 fElasticXsc = hnXsc;
451 }
452 else if( pLab < 1.05 )
453 {
454 hnXsc = 23 + 40*(std::log(pLab/0.73))*
455 (std::log(pLab/0.73));
456 fElasticXsc = 23 + 20*(std::log(pLab/0.73))*
457 (std::log(pLab/0.73));
458 }
459 else // 1.05 - 10 GeV/c
460 {
461 hnXsc = 39.0+75*(pLab - 1.2)/(std::pow(pLab,3.0) + 0.15);
462
463 fElasticXsc = 6 + 20/( (logP-0.182)*(logP-0.182) + 1.0 );
464 }
465 fTotalXsc = hnXsc;
466 }
467 if( proton ) // pn to be np
468 {
469 if( pLab < 0.02 )
470 {
471 hpXsc = 4100+30*std::pow(std::log(1.3/pLab),3.6); // was as pLab < 0.8
472 fElasticXsc = hpXsc;
473 }
474 else if( pLab < 0.8 )
475 {
476 hpXsc = 33+30*std::pow(std::log(pLab/1.3),4.0);
477 fElasticXsc = hpXsc;
478 }
479 else if( pLab < 1.05 )
480 {
481 hpXsc = 33+30*std::pow(std::log(pLab/0.95),2.0);
482 fElasticXsc = 6 + 52/( std::log(0.511/pLab)*std::log(0.511/pLab) + 1.6 );
483 }
484 else if( pLab < 1.4 )
485 {
486 hpXsc = 33+30*std::pow(std::log(pLab/0.95),2.0);
487 fElasticXsc = 6 + 52/( std::log(0.511/pLab)*std::log(0.511/pLab) + 1.6 );
488 }
489 else // 1.4 < pLab < 10. )
490 {
491 hpXsc = 33.3 + 20.8*(std::pow(pLab,2.0) - 1.35)/(std::pow(pLab,2.50) + 0.95);
492
493 fElasticXsc = 6 + 20/( (logP-0.182)*(logP-0.182) + 1.0 );
494 }
495 fTotalXsc = hpXsc;
496 }
497 }
498 }
499 else if( theParticle == theProton && pORn ) ////// proton //////////////////////////////////////////////
500 {
501 if( pLab >= 373.) // pdg due to TOTEM data
502 {
503 xsection = GetHadronNucleonXscPDG(aParticle, nucleon)/millibarn;
504
505 fElasticXsc = 6.5 + 0.308*std::pow(std::log(sMand/400.),1.65) + 9.19*std::pow(sMand,-0.458);
506
507 fTotalXsc = xsection;
508 }
509 else if( pLab >= 100.)
510 {
511 B0 = 7.5;
512 A0 = 100. - B0*std::log(3.0e7);
513
514 xsection = A0 + B0*std::log(pE) - 11 + 103*std::pow(sMand,-0.165); // mb
515
516 fElasticXsc = 5.53 + 0.308*std::pow(std::log(sMand/28.9),1.1) + 9.19*std::pow(sMand,-0.458);
517
518 fTotalXsc = xsection;
519 }
520 else if( pLab >= 10.)
521 {
522 B0 = 7.5;
523 A0 = 100. - B0*std::log(3.0e7);
524
525 xsection = A0 + B0*std::log(pE) - 11 + 103*std::pow(sMand,-0.165); // mb
526
527 fElasticXsc = 6 + 20/( (logP-0.182)*(logP-0.182) + 1.0 );
528
529 fTotalXsc = xsection;
530 }
531 else
532 {
533 // pp
534
535 if( proton )
536 {
537 if( pLab < 0.4 )
538 {
539 hpXsc = 23 + 50*( std::pow( std::log(0.73/pLab), 3.5 ) );
540 fElasticXsc = hpXsc;
541 }
542 else if( pLab < 0.73 )
543 {
544 hpXsc = 23 + 50*( std::pow( std::log(0.73/pLab), 3.5 ) );
545 fElasticXsc = hpXsc;
546 }
547 else if( pLab < 1.05 )
548 {
549 hpXsc = 23 + 40*(std::log(pLab/0.73))*
550 (std::log(pLab/0.73));
551 fElasticXsc = 23 + 20*(std::log(pLab/0.73))*
552 (std::log(pLab/0.73));
553 }
554 else // 1.05 - 10 GeV/c
555 {
556 hpXsc = 39.0+75*(pLab - 1.2)/(std::pow(pLab,3.0) + 0.15);
557
558 fElasticXsc = 6 + 20/( (logP-0.182)*(logP-0.182) + 1.0 );
559 }
560 fTotalXsc = hpXsc;
561 }
562 if( neutron ) // pn to be np
563 {
564 if( pLab < 0.02 )
565 {
566 hnXsc = 4100+30*std::pow(std::log(1.3/pLab),3.6); // was as pLab < 0.8
567 fElasticXsc = hnXsc;
568 }
569 else if( pLab < 0.8 )
570 {
571 hnXsc = 33+30*std::pow(std::log(pLab/1.3),4.0);
572 fElasticXsc = hnXsc;
573 }
574 else if( pLab < 1.05 )
575 {
576 hnXsc = 33+30*std::pow(std::log(pLab/0.95),2.0);
577 fElasticXsc = 6 + 52/( std::log(0.511/pLab)*std::log(0.511/pLab) + 1.6 );
578 }
579 else if( pLab < 1.4 )
580 {
581 hnXsc = 33+30*std::pow(std::log(pLab/0.95),2.0);
582 fElasticXsc = 6 + 52/( std::log(0.511/pLab)*std::log(0.511/pLab) + 1.6 );
583 }
584 else // 1.4 < pLab < 10. )
585 {
586 hnXsc = 33.3 + 20.8*(std::pow(pLab,2.0) - 1.35)/(std::pow(pLab,2.50) + 0.95);
587
588 fElasticXsc = 6 + 20/( (logP-0.182)*(logP-0.182) + 1.0 );
589 }
590 fTotalXsc = hnXsc;
591 }
592 }
593 }
594 else if( theParticle == theAProton && pORn ) /////////////////// p_bar ///////////////////////////
595 {
596 if( proton )
597 {
598 xsection = 35.45 + B*std::pow(std::log(sMand/s0),2.)
599 + 42.53*std::pow(sMand,-eta1) + 33.34*std::pow(sMand,-eta2);
600 }
601 if( neutron ) // ???
602 {
603 xsection = 35.80 + B*std::pow(std::log(sMand/s0),2.)
604 + 40.15*std::pow(sMand,-eta1) + 30.*std::pow(sMand,-eta2);
605 }
606 fTotalXsc = xsection;
607 }
608 else if( theParticle == thePiPlus && pORn ) // pi+ /////////////////////////////////////////////
609 {
610 if( proton ) // pi+ p
611 {
612 if( pLab < 0.28 )
613 {
614 hpXsc = 10./((logP + 1.273)*(logP + 1.273) + 0.05);
615 fElasticXsc = hpXsc;
616 }
617 else if( pLab < 0.4 )
618 {
619 hpXsc = 14./( (logP + 1.273)*(logP + 1.273) + 0.07);
620 fElasticXsc = hpXsc;
621 }
622 else if( pLab < 0.68 )
623 {
624 hpXsc = 14./( (logP + 1.273)*(logP + 1.273) + 0.07);
625 fElasticXsc = hpXsc;
626 }
627 else if( pLab < 0.85 )
628 {
629 G4double Ex4 = 88*(std::log(pLab/0.77))*(std::log(pLab/0.77));
630 hpXsc = Ex4 + 14.9;
631 fElasticXsc = hpXsc*std::exp(-3.*(pLab - 0.68));
632 }
633 else if( pLab < 1.15 )
634 {
635 G4double Ex4 = 88*(std::log(pLab/0.77))*(std::log(pLab/0.77));
636 hpXsc = Ex4 + 14.9;
637
638 fElasticXsc = 6.0 + 1.4/(( pLab - 1.4)*( pLab - 1.4) + 0.1);
639 }
640 else if( pLab < 1.4) // ns original
641 {
642 G4double Ex1 = 3.2*std::exp(-(pLab-2.55)*(pLab-2.55)/0.55/0.55);
643 G4double Ex2 = 12*std::exp(-(pLab-1.47)*(pLab-1.47)/0.225/0.225);
644 hpXsc = Ex1 + Ex2 + 27.5;
645 fElasticXsc = 6.0 + 1.4/(( pLab - 1.4)*( pLab - 1.4) + 0.1);
646 }
647 else if( pLab < 2.0 ) // ns original
648 {
649 G4double Ex1 = 3.2*std::exp(-(pLab-2.55)*(pLab-2.55)/0.55/0.55);
650 G4double Ex2 = 12*std::exp(-(pLab-1.47)*(pLab-1.47)/0.225/0.225);
651 hpXsc = Ex1 + Ex2 + 27.5;
652 fElasticXsc = 3.0 + 1.36/( (logP - 0.336)*(logP - 0.336) + 0.08);
653 }
654 else if( pLab < 3.5 ) // ns original
655 {
656 G4double Ex1 = 3.2*std::exp(-(pLab-2.55)*(pLab-2.55)/0.55/0.55);
657 G4double Ex2 = 12*std::exp(-(pLab-1.47)*(pLab-1.47)/0.225/0.225);
658 hpXsc = Ex1 + Ex2 + 27.5;
659 fElasticXsc = 3.0 + 6.20/( (logP - 0.336)*(logP - 0.336) + 0.8);
660 }
661 else if( pLab < 200. ) // my
662 {
663 hpXsc = 10.6 + 2.*std::log(pE) + 25*std::pow(pE, -0.43 ); // ns original
664 // hpXsc = GetHadronNucleonXscPDG(aParticle, nucleon )/millibarn;
665 fElasticXsc = 3.0 + 6.20/( (logP - 0.336)*(logP - 0.336) + 0.8);
666 }
667 else // pLab > 100 // my
668 {
669 hpXsc = GetHadronNucleonXscPDG(aParticle, nucleon )/millibarn;
670 fElasticXsc = 3.0 + 6.20/( (logP - 0.336)*(logP - 0.336) + 0.8);
671 }
672 fTotalXsc = hpXsc;
673 }
674 if( neutron ) // pi+ n = pi- p??
675 {
676 if( pLab < 0.28 )
677 {
678 hnXsc = 0.288/((pLab - 0.28)*(pLab - 0.28) + 0.004);
679 fElasticXsc = 1.8/((logP + 1.273)*(logP + 1.273) + 0.07);
680 }
681 else if( pLab < 0.395676 ) // first peak
682 {
683 hnXsc = 0.648/((pLab - 0.28)*(pLab - 0.28) + 0.009);
684 fElasticXsc = 0.257/((pLab - 0.28)*(pLab - 0.28) + 0.01);
685 }
686 else if( pLab < 0.5 )
687 {
688 hnXsc = 26 + 110*(std::log(pLab/0.48))*(std::log(pLab/0.48));
689 fElasticXsc = 0.37*hnXsc;
690 }
691 else if( pLab < 0.65 )
692 {
693 hnXsc = 26 + 110*(std::log(pLab/0.48))*(std::log(pLab/0.48));
694 fElasticXsc = 0.95/((pLab - 0.72)*(pLab - 0.72) + 0.049);
695 }
696 else if( pLab < 0.72 )
697 {
698 hnXsc = 36.1 + 10*std::exp(-(pLab-0.72)*(pLab-0.72)/0.06/0.06)+
699 24*std::exp(-(pLab-1.015)*(pLab-1.015)/0.075/0.075);
700 fElasticXsc = 0.95/((pLab - 0.72)*(pLab - 0.72) + 0.049);
701 }
702 else if( pLab < 0.88 )
703 {
704 hnXsc = 36.1 + 10*std::exp(-(pLab-0.72)*(pLab-0.72)/0.06/0.06)+
705 24*std::exp(-(pLab-1.015)*(pLab-1.015)/0.075/0.075);
706 fElasticXsc = 0.95/((pLab - 0.72)*(pLab - 0.72) + 0.049);
707 }
708 else if( pLab < 1.03 )
709 {
710 hnXsc = 36.1 + 10*std::exp(-(pLab-0.72)*(pLab-0.72)/0.06/0.06)+
711 24*std::exp(-(pLab-1.015)*(pLab-1.015)/0.075/0.075);
712 fElasticXsc = 2.0 + 0.4/((pLab - 1.03)*(pLab - 1.03) + 0.016);
713 }
714 else if( pLab < 1.15 )
715 {
716 hnXsc = 36.1 + 10*std::exp(-(pLab-0.72)*(pLab-0.72)/0.06/0.06)+
717 24*std::exp(-(pLab-1.015)*(pLab-1.015)/0.075/0.075);
718 fElasticXsc = 2.0 + 0.4/((pLab - 1.03)*(pLab - 1.03) + 0.016);
719 }
720 else if( pLab < 1.3 )
721 {
722 hnXsc = 36.1 + 10*std::exp(-(pLab-0.72)*(pLab-0.72)/0.06/0.06)+
723 24*std::exp(-(pLab-1.015)*(pLab-1.015)/0.075/0.075);
724 fElasticXsc = 3. + 13./pLab;
725 }
726 else if( pLab < 2.6 ) // < 3.0) // ns original
727 {
728 hnXsc = 36.1 + 0.079-4.313*std::log(pLab)+
729 3*std::exp(-(pLab-2.1)*(pLab-2.1)/0.4/0.4)+
730 1.5*std::exp(-(pLab-1.4)*(pLab-1.4)/0.12/0.12);
731 fElasticXsc = 3. + 13./pLab;
732 }
733 else if( pLab < 20. ) // < 3.0) // ns original
734 {
735 hnXsc = 36.1 + 0.079 - 4.313*std::log(pLab)+
736 3*std::exp(-(pLab-2.1)*(pLab-2.1)/0.4/0.4)+
737 1.5*std::exp(-(pLab-1.4)*(pLab-1.4)/0.12/0.12);
738 fElasticXsc = 3. + 13./pLab;
739 }
740 else // mb
741 {
742 hnXsc = GetHadronNucleonXscPDG(aParticle, nucleon )/millibarn;
743 fElasticXsc = 3. + 13./pLab;
744 }
745 fTotalXsc = hnXsc;
746 }
747 }
748 else if( theParticle == thePiMinus && pORn ) /// pi- ////////////////////////////////////////////
749 {
750 if( neutron ) // pi- n = pi+ p??
751 {
752 if( pLab < 0.28 )
753 {
754 hnXsc = 10./((logP + 1.273)*(logP + 1.273) + 0.05);
755 fElasticXsc = hnXsc;
756 }
757 else if( pLab < 0.4 )
758 {
759 hnXsc = 14./( (logP + 1.273)*(logP + 1.273) + 0.07);
760 fElasticXsc = hnXsc;
761 }
762 else if( pLab < 0.68 )
763 {
764 hnXsc = 14./( (logP + 1.273)*(logP + 1.273) + 0.07);
765 fElasticXsc = hnXsc;
766 }
767 else if( pLab < 0.85 )
768 {
769 G4double Ex4 = 88*(std::log(pLab/0.77))*(std::log(pLab/0.77));
770 hnXsc = Ex4 + 14.9;
771 fElasticXsc = hnXsc*std::exp(-3.*(pLab - 0.68));
772 }
773 else if( pLab < 1.15 )
774 {
775 G4double Ex4 = 88*(std::log(pLab/0.77))*(std::log(pLab/0.77));
776 hnXsc = Ex4 + 14.9;
777
778 fElasticXsc = 6.0 + 1.4/(( pLab - 1.4)*( pLab - 1.4) + 0.1);
779 }
780 else if( pLab < 1.4) // ns original
781 {
782 G4double Ex1 = 3.2*std::exp(-(pLab-2.55)*(pLab-2.55)/0.55/0.55);
783 G4double Ex2 = 12*std::exp(-(pLab-1.47)*(pLab-1.47)/0.225/0.225);
784 hnXsc = Ex1 + Ex2 + 27.5;
785 fElasticXsc = 6.0 + 1.4/(( pLab - 1.4)*( pLab - 1.4) + 0.1);
786 }
787 else if( pLab < 2.0 ) // ns original
788 {
789 G4double Ex1 = 3.2*std::exp(-(pLab-2.55)*(pLab-2.55)/0.55/0.55);
790 G4double Ex2 = 12*std::exp(-(pLab-1.47)*(pLab-1.47)/0.225/0.225);
791 hnXsc = Ex1 + Ex2 + 27.5;
792 fElasticXsc = 3.0 + 1.36/( (logP - 0.336)*(logP - 0.336) + 0.08);
793 }
794 else if( pLab < 3.5 ) // ns original
795 {
796 G4double Ex1 = 3.2*std::exp(-(pLab-2.55)*(pLab-2.55)/0.55/0.55);
797 G4double Ex2 = 12*std::exp(-(pLab-1.47)*(pLab-1.47)/0.225/0.225);
798 hnXsc = Ex1 + Ex2 + 27.5;
799 fElasticXsc = 3.0 + 6.20/( (logP - 0.336)*(logP - 0.336) + 0.8);
800 }
801 else if( pLab < 200. ) // my
802 {
803 hnXsc = 10.6 + 2.*std::log(pE) + 25*std::pow(pE, -0.43 ); // ns original
804 fElasticXsc = 3.0 + 6.20/( (logP - 0.336)*(logP - 0.336) + 0.8);
805 }
806 else // pLab > 100 // my
807 {
808 hnXsc = GetHadronNucleonXscPDG(aParticle, nucleon )/millibarn;
809 fElasticXsc = 3.0 + 6.20/( (logP - 0.336)*(logP - 0.336) + 0.8);
810 }
811 fTotalXsc = hnXsc;
812 }
813 if( proton ) // pi- p
814 {
815 if( pLab < 0.28 )
816 {
817 hpXsc = 0.288/((pLab - 0.28)*(pLab - 0.28) + 0.004);
818 fElasticXsc = 1.8/((logP + 1.273)*(logP + 1.273) + 0.07);
819 }
820 else if( pLab < 0.395676 ) // first peak
821 {
822 hpXsc = 0.648/((pLab - 0.28)*(pLab - 0.28) + 0.009);
823 fElasticXsc = 0.257/((pLab - 0.28)*(pLab - 0.28) + 0.01);
824 }
825 else if( pLab < 0.5 )
826 {
827 hpXsc = 26 + 110*(std::log(pLab/0.48))*(std::log(pLab/0.48));
828 fElasticXsc = 0.37*hpXsc;
829 }
830 else if( pLab < 0.65 )
831 {
832 hpXsc = 26 + 110*(std::log(pLab/0.48))*(std::log(pLab/0.48));
833 fElasticXsc = 0.95/((pLab - 0.72)*(pLab - 0.72) + 0.049);
834 }
835 else if( pLab < 0.72 )
836 {
837 hpXsc = 36.1+
838 10*std::exp(-(pLab-0.72)*(pLab-0.72)/0.06/0.06)+
839 24*std::exp(-(pLab-1.015)*(pLab-1.015)/0.075/0.075);
840 fElasticXsc = 0.95/((pLab - 0.72)*(pLab - 0.72) + 0.049);
841 }
842 else if( pLab < 0.88 )
843 {
844 hpXsc = 36.1+
845 10*std::exp(-(pLab-0.72)*(pLab-0.72)/0.06/0.06)+
846 24*std::exp(-(pLab-1.015)*(pLab-1.015)/0.075/0.075);
847 fElasticXsc = 0.95/((pLab - 0.72)*(pLab - 0.72) + 0.049);
848 }
849 else if( pLab < 1.03 )
850 {
851 hpXsc = 36.1+
852 10*std::exp(-(pLab-0.72)*(pLab-0.72)/0.06/0.06)+
853 24*std::exp(-(pLab-1.015)*(pLab-1.015)/0.075/0.075);
854 fElasticXsc = 2.0 + 0.4/((pLab - 1.03)*(pLab - 1.03) + 0.016);
855 }
856 else if( pLab < 1.15 )
857 {
858 hpXsc = 36.1+
859 10*std::exp(-(pLab-0.72)*(pLab-0.72)/0.06/0.06)+
860 24*std::exp(-(pLab-1.015)*(pLab-1.015)/0.075/0.075);
861 fElasticXsc = 2.0 + 0.4/((pLab - 1.03)*(pLab - 1.03) + 0.016);
862 }
863 else if( pLab < 1.3 )
864 {
865 hpXsc = 36.1+
866 10*std::exp(-(pLab-0.72)*(pLab-0.72)/0.06/0.06)+
867 24*std::exp(-(pLab-1.015)*(pLab-1.015)/0.075/0.075);
868 fElasticXsc = 3. + 13./pLab;
869 }
870 else if( pLab < 2.6 ) // < 3.0) // ns original
871 {
872 hpXsc = 36.1+0.079-4.313*std::log(pLab)+
873 3*std::exp(-(pLab-2.1)*(pLab-2.1)/0.4/0.4)+
874 1.5*std::exp(-(pLab-1.4)*(pLab-1.4)/0.12/0.12);
875 fElasticXsc = 3. +13./pLab; // *std::log(pLab*6.79);
876 }
877 else // mb
878 {
879 hpXsc = GetHadronNucleonXscPDG(aParticle, nucleon )/millibarn;
880 fElasticXsc = 3. + 13./pLab;
881 }
882 fTotalXsc = hpXsc;
883 }
884 }
885 else if( theParticle == theKPlus && pORn )
886 {
887 if( proton )
888 {
889 xsection = 17.91 + B*std::pow(std::log(sMand/s0),2.)
890 + 7.14*std::pow(sMand,-eta1) - 13.45*std::pow(sMand,-eta2);
891 }
892 if( neutron )
893 {
894 xsection = 17.87 + B*std::pow(std::log(sMand/s0),2.)
895 + 5.17*std::pow(sMand,-eta1) - 7.23*std::pow(sMand,-eta2);
896 }
897 fTotalXsc = xsection;
898 }
899 else if( theParticle == theKMinus && pORn )
900 {
901 if( proton )
902 {
903 xsection = 17.91 + B*std::pow(std::log(sMand/s0),2.)
904 + 7.14*std::pow(sMand,-eta1) + 13.45*std::pow(sMand,-eta2);
905 }
906 if( neutron )
907 {
908 xsection = 17.87 + B*std::pow(std::log(sMand/s0),2.)
909 + 5.17*std::pow(sMand,-eta1) + 7.23*std::pow(sMand,-eta2);
910 }
911 fTotalXsc = xsection;
912 }
913 else if( theParticle == theSMinus && pORn )
914 {
915 xsection = 35.20 + B*std::pow(std::log(sMand/s0),2.)
916 - 199.*std::pow(sMand,-eta1) + 264.*std::pow(sMand,-eta2);
917 }
918 else if( theParticle == theGamma && pORn ) // modify later on
919 {
920 xsection = 0.0 + B*std::pow(std::log(sMand/s0),2.)
921 + 0.032*std::pow(sMand,-eta1) - 0.0*std::pow(sMand,-eta2);
922 fTotalXsc = xsection;
923 }
924 else // other then p,n,pi+,pi-,K+,K- as proton ???
925 {
926 if( proton )
927 {
928 xsection = 35.45 + B*std::pow(std::log(sMand/s0),2.)
929 + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2);
930 }
931 if( neutron )
932 {
933 xsection += 35.80 + B*std::pow(std::log(sMand/s0),2.)
934 + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2);
935 }
936 fTotalXsc = xsection;
937 }
938 fTotalXsc *= millibarn; // parametrised in mb
939 fElasticXsc *= millibarn; // parametrised in mb
940
941 if( proton && aParticle->GetDefinition()->GetPDGCharge() > 0. )
942 {
943 G4double cB = GetCoulombBarrier(aParticle, nucleon);
944 fTotalXsc *= cB;
945 fElasticXsc *= cB;
946 }
947 fInelasticXsc = fTotalXsc - fElasticXsc;
948 if( fInelasticXsc < 0. ) fInelasticXsc = 0.;
949
950 // G4cout<<fTotalXsc/millibarn<<"; "<<fElasticXsc/millibarn<<"; "<<fInelasticXsc/millibarn<<G4endl;
951
952 return fTotalXsc;
953}
954
955/////////////////////////////////////////////////////////////////////////////////////
956//
957// Returns hadron-nucleon cross-section based on V. Uzjinsky parametrisation of
958// data from G4FTFCrossSection class
959
962 const G4ParticleDefinition* nucleon )
963{
964 G4int PDGcode = aParticle->GetDefinition()->GetPDGEncoding();
965 G4int absPDGcode = std::abs(PDGcode);
966 G4double Elab = aParticle->GetTotalEnergy();
967 // (s - 2*0.88*GeV*GeV)/(2*0.939*GeV)/GeV;
968 G4double Plab = aParticle->GetMomentum().mag();
969 // std::sqrt(Elab * Elab - 0.88);
970
971 Elab /= GeV;
972 Plab /= GeV;
973
974 G4double LogPlab = std::log( Plab );
975 G4double sqrLogPlab = LogPlab * LogPlab;
976
977 G4bool pORn = (nucleon == theProton || nucleon == theNeutron );
978 G4bool proton = (nucleon == theProton);
979 G4bool neutron = (nucleon == theNeutron);
980
981
982 if( absPDGcode > 1000 && pORn ) //------Projectile is baryon -
983 {
984 if(proton)
985 {
986 fTotalXsc = 48.0 + 0. *std::pow(Plab, 0. ) + 0.522*sqrLogPlab - 4.51*LogPlab;
987 fElasticXsc = 11.9 + 26.9*std::pow(Plab,-1.21) + 0.169*sqrLogPlab - 1.85*LogPlab;
988 }
989 if(neutron)
990 {
991 fTotalXsc = 47.3 + 0. *std::pow(Plab, 0. ) + 0.513*sqrLogPlab - 4.27*LogPlab;
992 fElasticXsc = 11.9 + 26.9*std::pow(Plab,-1.21) + 0.169*sqrLogPlab - 1.85*LogPlab;
993 }
994 }
995 else if( PDGcode == 211 && pORn ) //------Projectile is PionPlus ----
996 {
997 if(proton)
998 {
999 fTotalXsc = 16.4 + 19.3 *std::pow(Plab,-0.42) + 0.19 *sqrLogPlab - 0.0 *LogPlab;
1000 fElasticXsc = 0.0 + 11.4*std::pow(Plab,-0.40) + 0.079*sqrLogPlab - 0.0 *LogPlab;
1001 }
1002 if(neutron)
1003 {
1004 fTotalXsc = 33.0 + 14.0 *std::pow(Plab,-1.36) + 0.456*sqrLogPlab - 4.03*LogPlab;
1005 fElasticXsc = 1.76 + 11.2*std::pow(Plab,-0.64) + 0.043*sqrLogPlab - 0.0 *LogPlab;
1006 }
1007 }
1008 else if( PDGcode == -211 && pORn ) //------Projectile is PionMinus ----
1009 {
1010 if(proton)
1011 {
1012 fTotalXsc = 33.0 + 14.0 *std::pow(Plab,-1.36) + 0.456*sqrLogPlab - 4.03*LogPlab;
1013 fElasticXsc = 1.76 + 11.2*std::pow(Plab,-0.64) + 0.043*sqrLogPlab - 0.0 *LogPlab;
1014 }
1015 if(neutron)
1016 {
1017 fTotalXsc = 16.4 + 19.3 *std::pow(Plab,-0.42) + 0.19 *sqrLogPlab - 0.0 *LogPlab;
1018 fElasticXsc = 0.0 + 11.4*std::pow(Plab,-0.40) + 0.079*sqrLogPlab - 0.0 *LogPlab;
1019 }
1020 }
1021 else if( PDGcode == 111 && pORn ) //------Projectile is PionZero --
1022 {
1023 if(proton)
1024 {
1025 fTotalXsc = (16.4 + 19.3 *std::pow(Plab,-0.42) + 0.19 *sqrLogPlab - 0.0 *LogPlab + //Pi+
1026 33.0 + 14.0 *std::pow(Plab,-1.36) + 0.456*sqrLogPlab - 4.03*LogPlab)/2; //Pi-
1027
1028 fElasticXsc = ( 0.0 + 11.4*std::pow(Plab,-0.40) + 0.079*sqrLogPlab - 0.0 *LogPlab + //Pi+
1029 1.76 + 11.2*std::pow(Plab,-0.64) + 0.043*sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1030
1031 }
1032 if(neutron)
1033 {
1034 fTotalXsc = (33.0 + 14.0 *std::pow(Plab,-1.36) + 0.456*sqrLogPlab - 4.03*LogPlab + //Pi+
1035 16.4 + 19.3 *std::pow(Plab,-0.42) + 0.19 *sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1036 fElasticXsc = ( 1.76 + 11.2*std::pow(Plab,-0.64) + 0.043*sqrLogPlab - 0.0 *LogPlab + //Pi+
1037 0.0 + 11.4*std::pow(Plab,-0.40) + 0.079*sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1038 }
1039 }
1040 else if( PDGcode == 321 && pORn ) //------Projectile is KaonPlus --
1041 {
1042 if(proton)
1043 {
1044 fTotalXsc = 18.1 + 0. *std::pow(Plab, 0. ) + 0.26 *sqrLogPlab - 1.0 *LogPlab;
1045 fElasticXsc = 5.0 + 8.1*std::pow(Plab,-1.8 ) + 0.16 *sqrLogPlab - 1.3 *LogPlab;
1046 }
1047 if(neutron)
1048 {
1049 fTotalXsc = 18.7 + 0. *std::pow(Plab, 0. ) + 0.21 *sqrLogPlab - 0.89*LogPlab;
1050 fElasticXsc = 7.3 + 0. *std::pow(Plab,-0. ) + 0.29 *sqrLogPlab - 2.4 *LogPlab;
1051 }
1052 }
1053 else if( PDGcode ==-321 && pORn ) //------Projectile is KaonMinus ----
1054 {
1055 if(proton)
1056 {
1057 fTotalXsc = 32.1 + 0. *std::pow(Plab, 0. ) + 0.66*sqrLogPlab - 5.6*LogPlab;
1058 fElasticXsc = 7.3 + 0. *std::pow(Plab,-0. ) + 0.29*sqrLogPlab - 2.4*LogPlab;
1059 }
1060 if(neutron)
1061 {
1062 fTotalXsc = 25.2 + 0. *std::pow(Plab, 0. ) + 0.38*sqrLogPlab - 2.9*LogPlab;
1063 fElasticXsc = 5.0 + 8.1*std::pow(Plab,-1.8 ) + 0.16*sqrLogPlab - 1.3*LogPlab;
1064 }
1065 }
1066 else if( PDGcode == 311 && pORn ) //------Projectile is KaonZero -----
1067 {
1068 if(proton)
1069 {
1070 fTotalXsc = ( 18.1 + 0. *std::pow(Plab, 0. ) + 0.26 *sqrLogPlab - 1.0 *LogPlab + //K+
1071 32.1 + 0. *std::pow(Plab, 0. ) + 0.66 *sqrLogPlab - 5.6 *LogPlab)/2; //K-
1072 fElasticXsc = ( 5.0 + 8.1*std::pow(Plab,-1.8 ) + 0.16 *sqrLogPlab - 1.3 *LogPlab + //K+
1073 7.3 + 0. *std::pow(Plab,-0. ) + 0.29 *sqrLogPlab - 2.4 *LogPlab)/2; //K-
1074 }
1075 if(neutron)
1076 {
1077 fTotalXsc = ( 18.7 + 0. *std::pow(Plab, 0. ) + 0.21 *sqrLogPlab - 0.89*LogPlab + //K+
1078 25.2 + 0. *std::pow(Plab, 0. ) + 0.38 *sqrLogPlab - 2.9 *LogPlab)/2; //K-
1079 fElasticXsc = ( 7.3 + 0. *std::pow(Plab,-0. ) + 0.29 *sqrLogPlab - 2.4 *LogPlab + //K+
1080 5.0 + 8.1*std::pow(Plab,-1.8 ) + 0.16 *sqrLogPlab - 1.3 *LogPlab)/2; //K-
1081 }
1082 }
1083 else //------Projectile is undefined, Nucleon assumed
1084 {
1085 if(proton)
1086 {
1087 fTotalXsc = 48.0 + 0. *std::pow(Plab, 0. ) + 0.522*sqrLogPlab - 4.51*LogPlab;
1088 fElasticXsc = 11.9 + 26.9*std::pow(Plab,-1.21) + 0.169*sqrLogPlab - 1.85*LogPlab;
1089 }
1090 if(neutron)
1091 {
1092 fTotalXsc = 47.3 + 0. *std::pow(Plab, 0. ) + 0.513*sqrLogPlab - 4.27*LogPlab;
1093 fElasticXsc = 11.9 + 26.9*std::pow(Plab,-1.21) + 0.169*sqrLogPlab - 1.85*LogPlab;
1094 }
1095 }
1096 fTotalXsc *= millibarn;
1097 fElasticXsc *= millibarn;
1098 fInelasticXsc = fTotalXsc - fElasticXsc;
1099 if (fInelasticXsc < 0.) fInelasticXsc = 0.;
1100
1101 return fTotalXsc;
1102}
1103
1104////////////////////////////////////////////////////////////////////////////////////
1105//
1106//
1107
1109 const G4double mt ,
1110 const G4double Plab )
1111{
1112 G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
1113 G4double Ecm = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
1114 // G4double Pcm = Plab * mt / Ecm;
1115 // G4double KEcm = std::sqrt ( Pcm * Pcm + mp * mp ) - mp;
1116
1117 return Ecm ; // KEcm;
1118}
1119
1120
1121////////////////////////////////////////////////////////////////////////////////////
1122//
1123//
1124
1126 const G4double mt ,
1127 const G4double Plab )
1128{
1129 G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
1130 G4double sMand = mp*mp + mt*mt + 2*Elab*mt ;
1131
1132 return sMand;
1133}
1134
1135
1136///////////////////////////////////////////////////////////////////////////////
1137//
1138//
1139
1141 const G4ParticleDefinition* nucleon )
1142{
1143 G4double ratio;
1144
1145 G4double tR = 0.895*fermi, pR;
1146
1147 if ( aParticle->GetDefinition() == theProton ) pR = 0.895*fermi;
1148 else if( aParticle->GetDefinition() == thePiPlus ) pR = 0.663*fermi;
1149 else if( aParticle->GetDefinition() == theKPlus ) pR = 0.340*fermi;
1150 else pR = 0.500*fermi;
1151
1152 G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
1153 G4double tZ = nucleon->GetPDGCharge();
1154
1155 G4double pTkin = aParticle->GetKineticEnergy();
1156
1157 G4double pM = aParticle->GetDefinition()->GetPDGMass();
1158 G4double tM = nucleon->GetPDGMass();
1159
1160 G4double pElab = pTkin + pM;
1161
1162 G4double totEcm = std::sqrt(pM*pM + tM*tM + 2.*pElab*tM);
1163
1164 G4double totTcm = totEcm - pM -tM;
1165
1166 G4double bC = fine_structure_const*hbarc*pZ*tZ;
1167 bC /= pR + tR;
1168 bC /= 2.; // 4., 2. parametrisation cof ??? vmg
1169
1170 // G4cout<<"pTkin = "<<pTkin/GeV<<"; pPlab = "
1171 // <<pPlab/GeV<<"; bC = "<<bC/GeV<<"; pTcm = "<<pTcm/GeV<<G4endl;
1172
1173 if( totTcm <= bC ) ratio = 0.;
1174 else ratio = 1. - bC/totTcm;
1175
1176 // if(ratio < DBL_MIN) ratio = DBL_MIN;
1177 if( ratio < 0.) ratio = 0.;
1178
1179 // G4cout <<"ratio = "<<ratio<<G4endl;
1180 return ratio;
1181}
1182
1183
1184
1185
1186
1187////////////////////////////////////////////////////////////////////////////////////
1188//
1189// Initialaise K(p,m)-(p,n) total cross section vectors
1190
1191
1193{
1194 G4int i = 0, iMax;
1195 G4double tmpxsc[106];
1196
1197 // Kp-proton tot xsc
1198
1199 iMax = 66;
1200 fKpProtonTotXscVector = G4LPhysicsFreeVector(iMax, fKpProtonTotTkin[0], fKpProtonTotTkin[iMax-1]);
1201 fKpProtonTotXscVector.SetSpline(true);
1202
1203 for( i = 0; i < iMax; i++)
1204 {
1205 tmpxsc[i] = 0.;
1206
1207 if( i == 0 || i == iMax-1 ) tmpxsc[i] = fKpProtonTotXsc[i];
1208 else tmpxsc[i] = 0.5*(fKpProtonTotXsc[i-1]+fKpProtonTotXsc[i+1]);
1209
1210 fKpProtonTotXscVector.PutValues(size_t(i), fKpProtonTotTkin[i], tmpxsc[i]*millibarn);
1211 }
1212
1213 // Kp-neutron tot xsc
1214
1215 iMax = 75;
1216 fKpNeutronTotXscVector = G4LPhysicsFreeVector(iMax, fKpNeutronTotTkin[0], fKpNeutronTotTkin[iMax-1]);
1217 fKpNeutronTotXscVector.SetSpline(true);
1218
1219 for( i = 0; i < iMax; i++)
1220 {
1221 tmpxsc[i] = 0.;
1222 if( i == 0 || i == iMax-1 ) tmpxsc[i] = fKpNeutronTotXsc[i];
1223 else tmpxsc[i] = 0.5*(fKpNeutronTotXsc[i-1]+fKpNeutronTotXsc[i+1]);
1224
1225 fKpNeutronTotXscVector.PutValues(size_t(i), fKpNeutronTotTkin[i], tmpxsc[i]*millibarn);
1226 }
1227
1228 // Km-proton tot xsc
1229
1230 iMax = 106;
1231 fKmProtonTotXscVector = G4LPhysicsFreeVector(iMax, fKmProtonTotTkin[0], fKmProtonTotTkin[iMax-1]);
1232 fKmProtonTotXscVector.SetSpline(true);
1233
1234 for( i = 0; i < iMax; i++)
1235 {
1236 tmpxsc[i] = 0.;
1237
1238 if( i == 0 || i == iMax-1 ) tmpxsc[i] = fKmProtonTotXsc[i];
1239 else tmpxsc[i] = 0.5*(fKmProtonTotXsc[i-1]+fKmProtonTotXsc[i+1]);
1240
1241 fKmProtonTotXscVector.PutValues(size_t(i), fKmProtonTotTkin[i], tmpxsc[i]*millibarn);
1242 }
1243
1244 // Km-neutron tot xsc
1245
1246 iMax = 68;
1247 fKmNeutronTotXscVector = G4LPhysicsFreeVector(iMax, fKmNeutronTotTkin[0], fKmNeutronTotTkin[iMax-1]);
1248 fKmNeutronTotXscVector.SetSpline(true);
1249
1250 for( i = 0; i < iMax; i++)
1251 {
1252 tmpxsc[i] = 0.;
1253 if( i == 0 || i == iMax-1 ) tmpxsc[i] = fKmNeutronTotXsc[i];
1254 else tmpxsc[i] = 0.5*(fKmNeutronTotXsc[i-1]+fKmNeutronTotXsc[i+1]);
1255
1256 fKmNeutronTotXscVector.PutValues(size_t(i), fKmNeutronTotTkin[i], tmpxsc[i]*millibarn);
1257 }
1258}
1259
1260///////////////////////////////////////////////////////
1261//
1262// K-nucleon tot xsc (mb) fit data, std::log(Tkin(MeV))
1263
1264const G4double G4HadronNucleonXsc::fKpProtonTotXsc[66] = {
12650.000000e+00, 1.592400e-01, 3.184700e-01, 7.961800e-01, 1.433120e+00, 2.070060e+00,
12662.866240e+00, 3.582800e+00, 4.378980e+00, 5.015920e+00, 5.573250e+00, 6.449040e+00,
12677.404460e+00, 8.200640e+00, 8.837580e+00, 9.713380e+00, 1.027070e+01, 1.090764e+01,
12681.130573e+01, 1.170382e+01, 1.242038e+01, 1.281847e+01, 1.321656e+01, 1.337580e+01,
12691.345541e+01, 1.329618e+01, 1.265924e+01, 1.242038e+01, 1.250000e+01, 1.305732e+01,
12701.369427e+01, 1.425159e+01, 1.544586e+01, 1.648089e+01, 1.751592e+01, 1.791401e+01,
12711.791401e+01, 1.775478e+01, 1.751592e+01, 1.735669e+01, 1.719745e+01, 1.711783e+01,
12721.703822e+01, 1.695860e+01, 1.695860e+01, 1.695860e+01, 1.695860e+01, 1.687898e+01,
12731.687898e+01, 1.703822e+01, 1.719745e+01, 1.735669e+01, 1.751592e+01, 1.767516e+01,
12741.783439e+01, 1.799363e+01, 1.815287e+01, 1.839172e+01, 1.855095e+01, 1.871019e+01,
12751.886943e+01, 1.918790e+01, 1.942675e+01, 1.966561e+01, 2.006369e+01, 2.054140e+01
1276}; // 66
1277
1278
1279const G4double G4HadronNucleonXsc::fKpProtonTotTkin[66] = {
12802.100000e+00, 2.180770e+00, 2.261540e+00, 2.396150e+00, 2.476920e+00, 2.557690e+00,
12812.557690e+00, 2.584620e+00, 2.638460e+00, 2.665380e+00, 2.719230e+00, 2.746150e+00,
12822.800000e+00, 2.853850e+00, 2.934620e+00, 3.042310e+00, 3.150000e+00, 3.311540e+00,
12833.446150e+00, 3.607690e+00, 3.930770e+00, 4.226920e+00, 4.361540e+00, 4.846150e+00,
12844.980770e+00, 5.088460e+00, 5.465380e+00, 5.653850e+00, 5.950000e+00, 6.084620e+00,
12856.246150e+00, 6.300000e+00, 6.380770e+00, 6.515380e+00, 6.730770e+00, 6.838460e+00,
12867.000000e+00, 7.161540e+00, 7.323080e+00, 7.457690e+00, 7.619230e+00, 7.780770e+00,
12877.915380e+00, 8.130770e+00, 8.265380e+00, 8.453850e+00, 8.642310e+00, 8.803850e+00,
12889.019230e+00, 9.234620e+00, 9.530770e+00, 9.773080e+00, 1.001538e+01, 1.017692e+01,
12891.033846e+01, 1.058077e+01, 1.082308e+01, 1.098462e+01, 1.114615e+01, 1.138846e+01,
12901.160385e+01, 1.173846e+01, 1.192692e+01, 1.216923e+01, 1.238461e+01, 1.257308e+01
1291}; // 66
1292
1293const G4double G4HadronNucleonXsc::fKpNeutronTotXsc[75] = {
12943.980900e-01, 3.184700e-01, 3.184700e-01, 3.980900e-01, 3.980900e-01, 3.980900e-01,
12953.980900e-01, 3.980900e-01, 3.980900e-01, 4.777100e-01, 3.980900e-01, 3.980900e-01,
12964.777100e-01, 5.573200e-01, 1.035030e+00, 1.512740e+00, 2.149680e+00, 2.786620e+00,
12973.503180e+00, 4.219750e+00, 5.015920e+00, 5.652870e+00, 6.289810e+00, 7.245220e+00,
12988.121020e+00, 8.837580e+00, 9.633760e+00, 1.042994e+01, 1.114650e+01, 1.194268e+01,
12991.265924e+01, 1.329618e+01, 1.393312e+01, 1.449045e+01, 1.496815e+01, 1.552548e+01,
13001.592357e+01, 1.664013e+01, 1.727707e+01, 1.783439e+01, 1.831210e+01, 1.902866e+01,
13011.902866e+01, 1.878981e+01, 1.847134e+01, 1.831210e+01, 1.807325e+01, 1.791401e+01,
13021.783439e+01, 1.767516e+01, 1.759554e+01, 1.743631e+01, 1.743631e+01, 1.751592e+01,
13031.743631e+01, 1.735669e+01, 1.751592e+01, 1.759554e+01, 1.767516e+01, 1.783439e+01,
13041.783439e+01, 1.791401e+01, 1.815287e+01, 1.823248e+01, 1.847134e+01, 1.878981e+01,
13051.894905e+01, 1.902866e+01, 1.934713e+01, 1.966561e+01, 1.990446e+01, 2.014331e+01,
13062.030255e+01, 2.046178e+01, 2.085987e+01
1307}; // 75
1308
1309const G4double G4HadronNucleonXsc::fKpNeutronTotTkin[75] = {
13102.692000e-02, 1.615400e-01, 2.961500e-01, 4.576900e-01, 6.461500e-01, 7.538500e-01,
13118.884600e-01, 1.103850e+00, 1.211540e+00, 1.400000e+00, 1.561540e+00, 1.776920e+00,
13121.992310e+00, 2.126920e+00, 2.342310e+00, 2.423080e+00, 2.557690e+00, 2.692310e+00,
13132.800000e+00, 2.988460e+00, 3.203850e+00, 3.365380e+00, 3.500000e+00, 3.688460e+00,
13143.850000e+00, 4.011540e+00, 4.173080e+00, 4.415380e+00, 4.630770e+00, 4.873080e+00,
13155.061540e+00, 5.276920e+00, 5.492310e+00, 5.707690e+00, 5.896150e+00, 6.030770e+00,
13166.138460e+00, 6.219230e+00, 6.273080e+00, 6.326920e+00, 6.407690e+00, 6.650000e+00,
13176.784620e+00, 7.026920e+00, 7.242310e+00, 7.350000e+00, 7.484620e+00, 7.619230e+00,
13187.807690e+00, 7.915380e+00, 8.050000e+00, 8.211540e+00, 8.453850e+00, 8.588460e+00,
13198.830770e+00, 9.073080e+00, 9.288460e+00, 9.476920e+00, 9.665380e+00, 9.826920e+00,
13201.004231e+01, 1.031154e+01, 1.052692e+01, 1.071538e+01, 1.095769e+01, 1.120000e+01,
13211.138846e+01, 1.155000e+01, 1.176538e+01, 1.190000e+01, 1.214231e+01, 1.222308e+01,
13221.238461e+01, 1.246538e+01, 1.265385e+01
1323}; // 75
1324
1325const G4double G4HadronNucleonXsc::fKmProtonTotXsc[106] = {
13261.136585e+02, 9.749129e+01, 9.275262e+01, 8.885017e+01, 8.334146e+01, 7.955401e+01,
13277.504530e+01, 7.153658e+01, 6.858537e+01, 6.674913e+01, 6.525784e+01, 6.448781e+01,
13286.360279e+01, 6.255401e+01, 6.127526e+01, 6.032404e+01, 5.997910e+01, 5.443554e+01,
13295.376307e+01, 5.236934e+01, 5.113937e+01, 5.090941e+01, 4.967944e+01, 4.844948e+01,
13304.705575e+01, 4.638327e+01, 4.571080e+01, 4.475958e+01, 4.397213e+01, 4.257840e+01,
13314.102090e+01, 4.090592e+01, 3.906969e+01, 3.839721e+01, 3.756097e+01, 3.644599e+01,
13323.560976e+01, 3.533101e+01, 3.533101e+01, 3.644599e+01, 3.811847e+01, 3.839721e+01,
13333.979094e+01, 4.090592e+01, 4.257840e+01, 4.341463e+01, 4.425087e+01, 4.564460e+01,
13344.759582e+01, 4.703833e+01, 4.843206e+01, 4.787457e+01, 4.452962e+01, 4.202090e+01,
13354.034843e+01, 3.839721e+01, 3.616725e+01, 3.365854e+01, 3.170732e+01, 3.087108e+01,
13363.170732e+01, 3.254355e+01, 3.310104e+01, 3.254355e+01, 3.142857e+01, 3.059233e+01,
13372.947735e+01, 2.891986e+01, 2.836237e+01, 2.752613e+01, 2.696864e+01, 2.641115e+01,
13382.501742e+01, 2.473868e+01, 2.418118e+01, 2.362369e+01, 2.334495e+01, 2.278746e+01,
13392.250871e+01, 2.222997e+01, 2.167247e+01, 2.139373e+01, 2.139373e+01, 2.139373e+01,
13402.111498e+01, 2.083624e+01, 2.055749e+01, 2.083624e+01, 2.055749e+01, 2.083624e+01,
13412.083624e+01, 2.055749e+01, 2.055749e+01, 2.055749e+01, 2.027875e+01, 2.000000e+01,
13422.055749e+01, 2.027875e+01, 2.083624e+01, 2.083624e+01, 2.055749e+01, 2.083624e+01,
13432.083624e+01, 2.083624e+01, 2.139373e+01, 2.139373e+01
1344}; // 106
1345
1346const G4double G4HadronNucleonXsc::fKmProtonTotTkin[106] = {
13474.017980e+00, 4.125840e+00, 4.179780e+00, 4.251690e+00, 4.287640e+00, 4.341570e+00,
13484.395510e+00, 4.467420e+00, 4.503370e+00, 4.575280e+00, 4.683150e+00, 4.737080e+00,
13494.773030e+00, 4.826970e+00, 4.880900e+00, 4.916850e+00, 4.952810e+00, 4.988760e+00,
13504.988760e+00, 5.006740e+00, 5.006740e+00, 5.042700e+00, 5.078650e+00, 5.114610e+00,
13515.132580e+00, 5.150560e+00, 5.186520e+00, 5.204490e+00, 5.276400e+00, 5.348310e+00,
13525.366290e+00, 5.384270e+00, 5.456180e+00, 5.564040e+00, 5.600000e+00, 5.671910e+00,
13535.743820e+00, 5.833710e+00, 5.905620e+00, 5.977530e+00, 6.085390e+00, 6.085390e+00,
13546.157300e+00, 6.175280e+00, 6.211240e+00, 6.229210e+00, 6.247190e+00, 6.337080e+00,
13556.391010e+00, 6.516850e+00, 6.462920e+00, 6.498880e+00, 6.570790e+00, 6.606740e+00,
13566.660670e+00, 6.678650e+00, 6.696630e+00, 6.732580e+00, 6.804490e+00, 6.876400e+00,
13576.948310e+00, 7.020220e+00, 7.074160e+00, 7.182020e+00, 7.235960e+00, 7.289890e+00,
13587.397750e+00, 7.523600e+00, 7.631460e+00, 7.757300e+00, 7.901120e+00, 8.062920e+00,
13598.260670e+00, 8.386520e+00, 8.530340e+00, 8.674160e+00, 8.817980e+00, 8.943820e+00,
13609.087640e+00, 9.267420e+00, 9.429210e+00, 9.573030e+00, 9.698880e+00, 9.896630e+00,
13611.002247e+01, 1.016629e+01, 1.031011e+01, 1.048989e+01, 1.063371e+01, 1.077753e+01,
13621.095730e+01, 1.108315e+01, 1.120899e+01, 1.135281e+01, 1.149663e+01, 1.162247e+01,
13631.174831e+01, 1.187416e+01, 1.200000e+01, 1.212584e+01, 1.221573e+01, 1.234157e+01,
13641.239551e+01, 1.250337e+01, 1.261124e+01, 1.273708e+01
1365}; // 106
1366
1367const G4double G4HadronNucleonXsc::fKmNeutronTotXsc[68] = {
13682.621810e+01, 2.741123e+01, 2.868413e+01, 2.963889e+01, 3.067343e+01, 3.178759e+01,
13693.282148e+01, 3.417466e+01, 3.536778e+01, 3.552620e+01, 3.544576e+01, 3.496756e+01,
13703.433030e+01, 3.401166e+01, 3.313537e+01, 3.257772e+01, 3.178105e+01, 3.138264e+01,
13713.074553e+01, 2.970952e+01, 2.891301e+01, 2.827542e+01, 2.787700e+01, 2.715978e+01,
13722.660181e+01, 2.612394e+01, 2.564574e+01, 2.516721e+01, 2.421098e+01, 2.365235e+01,
13732.317366e+01, 2.261437e+01, 2.237389e+01, 2.205427e+01, 2.181395e+01, 2.165357e+01,
13742.149335e+01, 2.133297e+01, 2.109232e+01, 2.093128e+01, 2.069030e+01, 2.052992e+01,
13752.028927e+01, 2.012824e+01, 1.996737e+01, 1.996590e+01, 1.988530e+01, 1.964432e+01,
13761.948361e+01, 1.940236e+01, 1.940040e+01, 1.931882e+01, 1.947593e+01, 1.947429e+01,
13771.939320e+01, 1.939157e+01, 1.946922e+01, 1.962715e+01, 1.970481e+01, 1.970301e+01,
13781.993958e+01, 2.009669e+01, 2.025380e+01, 2.033178e+01, 2.049003e+01, 2.064747e+01,
13792.080540e+01, 2.096333e+01
1380}; // 68
1381
1382const G4double G4HadronNucleonXsc::fKmNeutronTotTkin[68] = {
13835.708500e+00, 5.809560e+00, 5.896270e+00, 5.954120e+00, 5.997630e+00, 6.041160e+00,
13846.142160e+00, 6.171410e+00, 6.272470e+00, 6.344390e+00, 6.416230e+00, 6.459180e+00,
13856.487690e+00, 6.501940e+00, 6.544740e+00, 6.573280e+00, 6.616110e+00, 6.644710e+00,
13866.658840e+00, 6.744700e+00, 6.773150e+00, 6.830410e+00, 6.859010e+00, 6.916240e+00,
13876.973530e+00, 6.987730e+00, 7.030670e+00, 7.102360e+00, 7.173880e+00, 7.288660e+00,
13887.374720e+00, 7.547000e+00, 7.690650e+00, 7.791150e+00, 7.920420e+00, 8.020980e+00,
13898.107160e+00, 8.207720e+00, 8.365740e+00, 8.523790e+00, 8.710560e+00, 8.811110e+00,
13908.969140e+00, 9.127190e+00, 9.270860e+00, 9.400230e+00, 9.486440e+00, 9.673210e+00,
13919.802510e+00, 9.946220e+00, 1.011870e+01, 1.029116e+01, 1.047808e+01, 1.062181e+01,
13921.075114e+01, 1.089488e+01, 1.106739e+01, 1.118244e+01, 1.135496e+01, 1.151307e+01,
13931.171439e+01, 1.190130e+01, 1.208822e+01, 1.223199e+01, 1.231829e+01, 1.247646e+01,
13941.259150e+01, 1.270655e+01
1395}; // 68
1396
1397
1398
1399
1400
1401//
1402//
1403///////////////////////////////////////////////////////////////////////////////////////
@ neutron
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
double mag() const
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4AntiLambda * AntiLambda()
static G4AntiNeutron * AntiNeutron()
static G4AntiOmegaMinus * AntiOmegaMinus()
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4AntiSigmaMinus * AntiSigmaMinus()
static G4AntiSigmaPlus * AntiSigmaPlus()
static G4AntiSigmaZero * AntiSigmaZero()
static G4AntiXiMinus * AntiXiMinus()
static G4AntiXiZero * AntiXiZero()
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
G4double GetMass() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4ThreeVector GetMomentum() const
G4double GetZ() const
Definition: G4Element.hh:131
G4double GetN() const
Definition: G4Element.hh:134
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
virtual G4bool IsApplicable(const G4DynamicParticle *aDP, const G4Element *)
virtual G4bool IsIsoApplicable(const G4DynamicParticle *aDP, G4int Z, G4int A)
G4double CalcMandelstamS(const G4double, const G4double, const G4double)
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4ParticleDefinition *)
void CrossSectionDescription(std::ostream &) const
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double GetHadronNucleonXscVU(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double CalculateEcmValue(const G4double, const G4double, const G4double)
G4double GetHadronNucleonXscEL(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double GetCoulombBarrier(const G4DynamicParticle *aParticle, const G4ParticleDefinition *nucleon)
static G4He3 * He3()
Definition: G4He3.cc:94
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
static G4KaonZeroLong * KaonZeroLong()
static G4KaonZeroShort * KaonZeroShort()
void PutValues(size_t binNumber, G4double binValue, G4double dataValue)
static G4Lambda * Lambda()
Definition: G4Lambda.cc:108
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4OmegaMinus * OmegaMinus()
G4double GetPDGCharge() const
void SetSpline(G4bool)
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4PionZero * PionZero()
Definition: G4PionZero.cc:104
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4SigmaMinus * SigmaMinus()
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
static G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:99
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4XiMinus * XiMinus()
Definition: G4XiMinus.cc:106
static G4XiZero * XiZero()
Definition: G4XiZero.cc:106
int G4lrint(double ad)
Definition: templates.hh:163