00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030 #include "G4ComponentGGHadronNucleusXsc.hh"
00031
00032 #include "G4PhysicalConstants.hh"
00033 #include "G4SystemOfUnits.hh"
00034 #include "G4ParticleTable.hh"
00035 #include "G4IonTable.hh"
00036 #include "G4ParticleDefinition.hh"
00037 #include "G4DynamicParticle.hh"
00038 #include "G4HadronNucleonXsc.hh"
00039
00040
00042
00043
00044 G4ComponentGGHadronNucleusXsc::G4ComponentGGHadronNucleusXsc()
00045 : G4VComponentCrossSection("Glauber-Gribov"),
00046 fUpperLimit(100000*GeV), fLowerLimit(10.*MeV),
00047 fRadiusConst(1.08*fermi),
00048 fTotalXsc(0.0), fElasticXsc(0.0), fInelasticXsc(0.0), fProductionXsc(0.0),
00049 fDiffractionXsc(0.0), fHadronNucleonXsc(0.0)
00050 {
00051 theGamma = G4Gamma::Gamma();
00052 theProton = G4Proton::Proton();
00053 theNeutron = G4Neutron::Neutron();
00054 theAProton = G4AntiProton::AntiProton();
00055 theANeutron = G4AntiNeutron::AntiNeutron();
00056 thePiPlus = G4PionPlus::PionPlus();
00057 thePiMinus = G4PionMinus::PionMinus();
00058 thePiZero = G4PionZero::PionZero();
00059 theKPlus = G4KaonPlus::KaonPlus();
00060 theKMinus = G4KaonMinus::KaonMinus();
00061 theK0S = G4KaonZeroShort::KaonZeroShort();
00062 theK0L = G4KaonZeroLong::KaonZeroLong();
00063 theL = G4Lambda::Lambda();
00064 theAntiL = G4AntiLambda::AntiLambda();
00065 theSPlus = G4SigmaPlus::SigmaPlus();
00066 theASPlus = G4AntiSigmaPlus::AntiSigmaPlus();
00067 theSMinus = G4SigmaMinus::SigmaMinus();
00068 theASMinus = G4AntiSigmaMinus::AntiSigmaMinus();
00069 theS0 = G4SigmaZero::SigmaZero();
00070 theAS0 = G4AntiSigmaZero::AntiSigmaZero();
00071 theXiMinus = G4XiMinus::XiMinus();
00072 theXi0 = G4XiZero::XiZero();
00073 theAXiMinus = G4AntiXiMinus::AntiXiMinus();
00074 theAXi0 = G4AntiXiZero::AntiXiZero();
00075 theOmega = G4OmegaMinus::OmegaMinus();
00076 theAOmega = G4AntiOmegaMinus::AntiOmegaMinus();
00077 theD = G4Deuteron::Deuteron();
00078 theT = G4Triton::Triton();
00079 theA = G4Alpha::Alpha();
00080 theHe3 = G4He3::He3();
00081
00082 hnXsc = new G4HadronNucleonXsc();
00083 }
00084
00086
00087
00088
00089 G4ComponentGGHadronNucleusXsc::~G4ComponentGGHadronNucleusXsc()
00090 {
00091 if (hnXsc) delete hnXsc;
00092 }
00093
00095
00096 G4double G4ComponentGGHadronNucleusXsc::GetTotalIsotopeCrossSection(const G4ParticleDefinition* aParticle,
00097 G4double kinEnergy,
00098 G4int Z, G4int A)
00099 {
00100 G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
00101 kinEnergy);
00102 fTotalXsc = GetIsoCrossSection(aDP, Z, A);
00103 delete aDP;
00104
00105 return fTotalXsc;
00106 }
00107
00109
00110 G4double G4ComponentGGHadronNucleusXsc::GetTotalElementCrossSection(const G4ParticleDefinition* aParticle,
00111 G4double kinEnergy,
00112 G4int Z, G4double A)
00113 {
00114 G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
00115 kinEnergy);
00116 fTotalXsc = GetIsoCrossSection(aDP, Z, G4int(A));
00117 delete aDP;
00118
00119 return fTotalXsc;
00120 }
00121
00123
00124 G4double G4ComponentGGHadronNucleusXsc::GetInelasticIsotopeCrossSection(const G4ParticleDefinition* aParticle,
00125 G4double kinEnergy,
00126 G4int Z, G4int A)
00127 {
00128 G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
00129 kinEnergy);
00130 fTotalXsc = GetIsoCrossSection(aDP, Z, A);
00131 delete aDP;
00132
00133 return fInelasticXsc;
00134 }
00135
00137
00138 G4double G4ComponentGGHadronNucleusXsc::GetInelasticElementCrossSection(const G4ParticleDefinition* aParticle,
00139 G4double kinEnergy,
00140 G4int Z, G4double A)
00141 {
00142 G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
00143 kinEnergy);
00144 fTotalXsc = GetIsoCrossSection(aDP, Z, G4int(A));
00145 delete aDP;
00146
00147 return fInelasticXsc;
00148 }
00149
00151
00152 G4double G4ComponentGGHadronNucleusXsc::GetElasticElementCrossSection(const G4ParticleDefinition* aParticle,
00153 G4double kinEnergy,
00154 G4int Z, G4double A)
00155 {
00156 G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
00157 kinEnergy);
00158 fTotalXsc = GetIsoCrossSection(aDP, Z, G4int(A));
00159 delete aDP;
00160
00161 return fElasticXsc;
00162 }
00163
00165
00166 G4double G4ComponentGGHadronNucleusXsc::GetElasticIsotopeCrossSection(const G4ParticleDefinition* aParticle,
00167 G4double kinEnergy,
00168 G4int Z, G4int A)
00169 {
00170 G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
00171 kinEnergy);
00172 fTotalXsc = GetIsoCrossSection(aDP, Z, A);
00173 delete aDP;
00174
00175 return fElasticXsc;
00176 }
00177
00179
00180 G4double G4ComponentGGHadronNucleusXsc::ComputeQuasiElasticRatio(const G4ParticleDefinition* aParticle,
00181 G4double kinEnergy,
00182 G4int Z, G4int A)
00183 {
00184 G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
00185 kinEnergy);
00186 fTotalXsc = GetIsoCrossSection(aDP, Z, A);
00187 delete aDP;
00188 G4double ratio = 0.;
00189
00190 if(fInelasticXsc > 0.)
00191 {
00192 ratio = (fInelasticXsc - fProductionXsc)/fInelasticXsc;
00193 if(ratio < 0.) ratio = 0.;
00194 }
00195 return ratio;
00196 }
00197
00198
00199
00200
00202
00203 G4bool
00204 G4ComponentGGHadronNucleusXsc::IsIsoApplicable(const G4DynamicParticle* aDP,
00205 G4int Z, G4int ,
00206 const G4Element*,
00207 const G4Material*)
00208 {
00209 G4bool applicable = false;
00210
00211 G4double kineticEnergy = aDP->GetKineticEnergy();
00212
00213 const G4ParticleDefinition* theParticle = aDP->GetDefinition();
00214
00215 if ( ( kineticEnergy >= fLowerLimit &&
00216 Z > 1 &&
00217 ( theParticle == theAProton ||
00218 theParticle == theGamma ||
00219 theParticle == theKPlus ||
00220 theParticle == theKMinus ||
00221 theParticle == theK0L ||
00222 theParticle == theK0S ||
00223 theParticle == theSMinus ||
00224 theParticle == theProton ||
00225 theParticle == theNeutron ||
00226 theParticle == thePiPlus ||
00227 theParticle == thePiMinus ) ) ) applicable = true;
00228
00229 return applicable;
00230 }
00231
00233
00234
00235
00236
00237
00238
00239 G4double
00240 G4ComponentGGHadronNucleusXsc::GetIsoCrossSection(const G4DynamicParticle* aParticle,
00241 G4int Z, G4int A,
00242 const G4Isotope*,
00243 const G4Element*,
00244 const G4Material*)
00245 {
00246 G4double xsection, sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
00247 G4double hpInXsc(0.), hnInXsc(0.);
00248 G4double R = GetNucleusRadius(A);
00249
00250 G4int N = A - Z;
00251 if (N < 0) N = 0;
00252
00253 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
00254
00255 if( theParticle == theProton ||
00256 theParticle == theNeutron ||
00257 theParticle == thePiPlus ||
00258 theParticle == thePiMinus )
00259 {
00260
00261
00262 sigma = Z*hnXsc->GetHadronNucleonXscNS(aParticle, theProton);
00263
00264 hpInXsc = hnXsc->GetInelasticHadronNucleonXsc();
00265
00266 sigma += N*hnXsc->GetHadronNucleonXscNS(aParticle, theNeutron);
00267
00268 hnInXsc = hnXsc->GetInelasticHadronNucleonXsc();
00269
00270 cofInelastic = 2.4;
00271 cofTotal = 2.0;
00272 }
00273 else if( theParticle == theKPlus ||
00274 theParticle == theKMinus ||
00275 theParticle == theK0S ||
00276 theParticle == theK0L )
00277 {
00278 sigma = GetKaonNucleonXscVector(aParticle, A, Z);
00279 cofInelastic = 2.2;
00280 cofTotal = 2.0;
00281 R = 1.3*fermi;
00282 R *= std::pow(G4double(A), 0.3333);
00283 }
00284 else
00285 {
00286 sigma = GetHadronNucleonXscNS(aParticle, A, Z);
00287 cofInelastic = 2.2;
00288 cofTotal = 2.0;
00289 }
00290
00291
00292 if( A > 1 )
00293 {
00294 nucleusSquare = cofTotal*pi*R*R;
00295 ratio = sigma/nucleusSquare;
00296
00297 xsection = nucleusSquare*std::log( 1. + ratio );
00298
00299 xsection *= GetParticleBarCorTot(theParticle, Z);
00300
00301 fTotalXsc = xsection;
00302
00303
00304
00305 fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
00306
00307 fInelasticXsc *= GetParticleBarCorIn(theParticle, Z);
00308
00309 fElasticXsc = fTotalXsc - fInelasticXsc;
00310
00311 if(fElasticXsc < 0.) fElasticXsc = 0.;
00312
00313 G4double difratio = ratio/(1.+ratio);
00314
00315 fDiffractionXsc = 0.5*nucleusSquare*( difratio - std::log( 1. + difratio ) );
00316
00317
00318
00319
00320 sigma = Z*hpInXsc + N*hnInXsc;
00321
00322 ratio = sigma/nucleusSquare;
00323
00324 fProductionXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
00325
00326 if (fElasticXsc < 0.) fElasticXsc = 0.;
00327 }
00328 else
00329 {
00330 fTotalXsc = sigma;
00331 xsection = sigma;
00332
00333 if ( theParticle != theAProton )
00334 {
00335 sigma = GetHNinelasticXsc(aParticle, A, Z);
00336 fInelasticXsc = sigma;
00337 fElasticXsc = fTotalXsc - fInelasticXsc;
00338 }
00339 else
00340 {
00341 fElasticXsc = fTotalXsc - fInelasticXsc;
00342 }
00343 if (fElasticXsc < 0.) fElasticXsc = 0.;
00344
00345 }
00346 return xsection;
00347 }
00348
00350
00351
00352
00353 G4double G4ComponentGGHadronNucleusXsc::
00354 GetRatioSD(const G4DynamicParticle* aParticle, G4int A, G4int Z)
00355 {
00356 G4double sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
00357 G4double R = GetNucleusRadius(A);
00358
00359 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
00360
00361 if( theParticle == theProton ||
00362 theParticle == theNeutron ||
00363 theParticle == thePiPlus ||
00364 theParticle == thePiMinus )
00365 {
00366 sigma = GetHadronNucleonXscNS(aParticle, A, Z);
00367 cofInelastic = 2.4;
00368 cofTotal = 2.0;
00369 }
00370 else
00371 {
00372 sigma = GetHadronNucleonXscNS(aParticle, A, Z);
00373 cofInelastic = 2.2;
00374 cofTotal = 2.0;
00375 }
00376 nucleusSquare = cofTotal*pi*R*R;
00377 ratio = sigma/nucleusSquare;
00378
00379 fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
00380
00381 G4double difratio = ratio/(1.+ratio);
00382
00383 fDiffractionXsc = 0.5*nucleusSquare*( difratio - std::log( 1. + difratio ) );
00384
00385 if (fInelasticXsc > 0.) ratio = fDiffractionXsc/fInelasticXsc;
00386 else ratio = 0.;
00387
00388 return ratio;
00389 }
00390
00392
00393
00394
00395 G4double G4ComponentGGHadronNucleusXsc::
00396 GetRatioQE(const G4DynamicParticle* aParticle, G4int A, G4int Z)
00397 {
00398 G4double sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
00399 G4double R = GetNucleusRadius(A);
00400
00401 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
00402
00403 if( theParticle == theProton ||
00404 theParticle == theNeutron ||
00405 theParticle == thePiPlus ||
00406 theParticle == thePiMinus )
00407 {
00408 sigma = GetHadronNucleonXscNS(aParticle, A, Z);
00409 cofInelastic = 2.4;
00410 cofTotal = 2.0;
00411 }
00412 else
00413 {
00414 sigma = GetHadronNucleonXscNS(aParticle, A, Z);
00415 cofInelastic = 2.2;
00416 cofTotal = 2.0;
00417 }
00418 nucleusSquare = cofTotal*pi*R*R;
00419 ratio = sigma/nucleusSquare;
00420
00421 fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
00422
00423 sigma = GetHNinelasticXsc(aParticle, A, Z);
00424 ratio = sigma/nucleusSquare;
00425
00426 fProductionXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
00427
00428 if (fInelasticXsc > fProductionXsc) ratio = (fInelasticXsc-fProductionXsc)/fInelasticXsc;
00429 else ratio = 0.;
00430 if ( ratio < 0. ) ratio = 0.;
00431
00432 return ratio;
00433 }
00434
00436
00437
00438
00439
00440
00441
00442 G4double
00443 G4ComponentGGHadronNucleusXsc::GetHadronNucleonXsc(const G4DynamicParticle* aParticle,
00444 const G4Element* anElement)
00445 {
00446 G4int At = G4lrint(anElement->GetN());
00447 G4int Zt = G4lrint(anElement->GetZ());
00448
00449 return GetHadronNucleonXsc(aParticle, At, Zt);
00450 }
00451
00453
00454
00455
00456
00457
00458
00459 G4double
00460 G4ComponentGGHadronNucleusXsc::GetHadronNucleonXsc(const G4DynamicParticle* aParticle,
00461 G4int At, G4int )
00462 {
00463 G4double xsection;
00464
00465
00466
00467 G4double targ_mass = 0.939*GeV;
00468
00469 G4double proj_mass = aParticle->GetMass();
00470 G4double proj_momentum = aParticle->GetMomentum().mag();
00471 G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
00472
00473 sMand /= GeV*GeV;
00474 proj_momentum /= GeV;
00475
00476 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
00477
00478 G4double aa = At;
00479
00480 if(theParticle == theGamma)
00481 {
00482 xsection = aa*(0.0677*std::pow(sMand,0.0808) + 0.129*std::pow(sMand,-0.4525));
00483 }
00484 else if(theParticle == theNeutron)
00485 {
00486 xsection = aa*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
00487 }
00488 else if(theParticle == theProton)
00489 {
00490 xsection = aa*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
00491
00492
00493 }
00494 else if(theParticle == theAProton)
00495 {
00496 xsection = aa*( 21.70*std::pow(sMand,0.0808) + 98.39*std::pow(sMand,-0.4525));
00497 }
00498 else if(theParticle == thePiPlus)
00499 {
00500 xsection = aa*(13.63*std::pow(sMand,0.0808) + 27.56*std::pow(sMand,-0.4525));
00501 }
00502 else if(theParticle == thePiMinus)
00503 {
00504
00505 xsection = aa*(13.63*std::pow(sMand,0.0808) + 36.02*std::pow(sMand,-0.4525));
00506 }
00507 else if(theParticle == theKPlus)
00508 {
00509 xsection = aa*(11.82*std::pow(sMand,0.0808) + 8.15*std::pow(sMand,-0.4525));
00510 }
00511 else if(theParticle == theKMinus)
00512 {
00513 xsection = aa*(11.82*std::pow(sMand,0.0808) + 26.36*std::pow(sMand,-0.4525));
00514 }
00515 else
00516 {
00517 xsection = aa*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
00518 }
00519 xsection *= millibarn;
00520 return xsection;
00521 }
00522
00523
00525
00526
00527
00528
00529 G4double
00530 G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscPDG(const G4DynamicParticle* aParticle,
00531 const G4Element* anElement)
00532 {
00533 G4int At = G4lrint(anElement->GetN());
00534 G4int Zt = G4lrint(anElement->GetZ());
00535
00536 return GetHadronNucleonXscPDG(aParticle, At, Zt);
00537 }
00538
00539
00540
00541
00543
00544
00545
00546
00547
00548 G4double
00549 G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscPDG(const G4DynamicParticle* aParticle,
00550 G4int At, G4int Zt)
00551 {
00552 G4double xsection;
00553
00554 G4int Nt = At-Zt;
00555 if (Nt < 0) Nt = 0;
00556
00557 G4double zz = Zt;
00558 G4double aa = At;
00559 G4double nn = Nt;
00560
00561 G4double targ_mass = G4ParticleTable::GetParticleTable()->
00562 GetIonTable()->GetIonMass(Zt, At);
00563
00564 targ_mass = 0.939*GeV;
00565
00566 G4double proj_mass = aParticle->GetMass();
00567 G4double proj_momentum = aParticle->GetMomentum().mag();
00568
00569 G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
00570
00571 sMand /= GeV*GeV;
00572
00573
00574
00575 G4double s0 = 5.38*5.38;
00576 G4double eta1 = 0.458;
00577 G4double eta2 = 0.458;
00578 G4double B = 0.308;
00579
00580
00581 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
00582
00583
00584 if(theParticle == theNeutron)
00585 {
00586 xsection = zz*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
00587 + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
00588 xsection += nn*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
00589 + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
00590 }
00591 else if(theParticle == theProton)
00592 {
00593
00594 xsection = zz*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
00595 + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
00596
00597 xsection += nn*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
00598 + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
00599 }
00600 else if(theParticle == theAProton)
00601 {
00602 xsection = zz*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
00603 + 42.53*std::pow(sMand,-eta1) + 33.34*std::pow(sMand,-eta2));
00604
00605 xsection += nn*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
00606 + 40.15*std::pow(sMand,-eta1) + 30.*std::pow(sMand,-eta2));
00607 }
00608 else if(theParticle == thePiPlus)
00609 {
00610 xsection = aa*( 20.86 + B*std::pow(std::log(sMand/s0),2.)
00611 + 19.24*std::pow(sMand,-eta1) - 6.03*std::pow(sMand,-eta2));
00612 }
00613 else if(theParticle == thePiMinus)
00614 {
00615 xsection = aa*( 20.86 + B*std::pow(std::log(sMand/s0),2.)
00616 + 19.24*std::pow(sMand,-eta1) + 6.03*std::pow(sMand,-eta2));
00617 }
00618 else if(theParticle == theKPlus || theParticle == theK0L )
00619 {
00620 xsection = zz*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
00621 + 7.14*std::pow(sMand,-eta1) - 13.45*std::pow(sMand,-eta2));
00622
00623 xsection += nn*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
00624 + 5.17*std::pow(sMand,-eta1) - 7.23*std::pow(sMand,-eta2));
00625 }
00626 else if(theParticle == theKMinus || theParticle == theK0S )
00627 {
00628 xsection = zz*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
00629 + 7.14*std::pow(sMand,-eta1) + 13.45*std::pow(sMand,-eta2));
00630
00631 xsection += nn*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
00632 + 5.17*std::pow(sMand,-eta1) + 7.23*std::pow(sMand,-eta2));
00633 }
00634 else if(theParticle == theSMinus)
00635 {
00636 xsection = aa*( 35.20 + B*std::pow(std::log(sMand/s0),2.)
00637 - 199.*std::pow(sMand,-eta1) + 264.*std::pow(sMand,-eta2));
00638 }
00639 else if(theParticle == theGamma)
00640 {
00641 xsection = aa*( 0.0 + B*std::pow(std::log(sMand/s0),2.)
00642 + 0.032*std::pow(sMand,-eta1) - 0.0*std::pow(sMand,-eta2));
00643
00644 }
00645 else
00646 {
00647 xsection = zz*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
00648 + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
00649
00650 xsection += nn*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
00651 + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
00652 }
00653 xsection *= millibarn;
00654 return xsection;
00655 }
00656
00657
00659
00660
00661
00662
00663 G4double
00664 G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscNS(const G4DynamicParticle* aParticle,
00665 const G4Element* anElement)
00666 {
00667 G4int At = G4lrint(anElement->GetN());
00668 G4int Zt = G4lrint(anElement->GetZ());
00669
00670 return GetHadronNucleonXscNS(aParticle, At, Zt);
00671 }
00672
00673
00674
00675
00677
00678
00679
00680
00681 G4double
00682 G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscNS(const G4DynamicParticle* aParticle,
00683 G4int At, G4int Zt)
00684 {
00685 G4double xsection(0);
00686
00687 G4double A0, B0;
00688 G4double hpXscv(0);
00689 G4double hnXscv(0);
00690
00691 G4int Nt = At-Zt;
00692 if (Nt < 0) Nt = 0;
00693
00694 G4double aa = At;
00695 G4double zz = Zt;
00696 G4double nn = Nt;
00697
00698 G4double targ_mass = G4ParticleTable::GetParticleTable()->
00699 GetIonTable()->GetIonMass(Zt, At);
00700
00701 targ_mass = 0.939*GeV;
00702
00703 G4double proj_mass = aParticle->GetMass();
00704 G4double proj_energy = aParticle->GetTotalEnergy();
00705 G4double proj_momentum = aParticle->GetMomentum().mag();
00706
00707 G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
00708
00709 sMand /= GeV*GeV;
00710 proj_momentum /= GeV;
00711 proj_energy /= GeV;
00712 proj_mass /= GeV;
00713
00714
00715
00716 G4double s0 = 5.38*5.38;
00717 G4double eta1 = 0.458;
00718 G4double eta2 = 0.458;
00719 G4double B = 0.308;
00720
00721
00722 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
00723
00724
00725 if(theParticle == theNeutron)
00726 {
00727 if( proj_momentum >= 373.)
00728 {
00729 return GetHadronNucleonXscPDG(aParticle,At,Zt);
00730 }
00731 else if( proj_momentum >= 10.)
00732
00733 {
00734
00735
00736
00737 if(proj_momentum >= 10.)
00738 {
00739 B0 = 7.5;
00740 A0 = 100. - B0*std::log(3.0e7);
00741
00742 xsection = A0 + B0*std::log(proj_energy) - 11
00743 + 103*std::pow(2*0.93827*proj_energy + proj_mass*proj_mass+
00744 0.93827*0.93827,-0.165);
00745 }
00746 xsection *= zz + nn;
00747 }
00748 else
00749 {
00750
00751
00752 if( proj_momentum < 0.73 )
00753 {
00754 hnXscv = 23 + 50*( std::pow( std::log(0.73/proj_momentum), 3.5 ) );
00755 }
00756 else if( proj_momentum < 1.05 )
00757 {
00758 hnXscv = 23 + 40*(std::log(proj_momentum/0.73))*
00759 (std::log(proj_momentum/0.73));
00760 }
00761 else
00762 {
00763 hnXscv = 39.0+
00764 75*(proj_momentum - 1.2)/(std::pow(proj_momentum,3.0) + 0.15);
00765 }
00766
00767
00768 if( proj_momentum < 0.8 )
00769 {
00770 hpXscv = 33+30*std::pow(std::log(proj_momentum/1.3),4.0);
00771 }
00772 else if( proj_momentum < 1.4 )
00773 {
00774 hpXscv = 33+30*std::pow(std::log(proj_momentum/0.95),2.0);
00775 }
00776 else
00777 {
00778 hpXscv = 33.3+
00779 20.8*(std::pow(proj_momentum,2.0)-1.35)/
00780 (std::pow(proj_momentum,2.50)+0.95);
00781 }
00782 xsection = hpXscv*zz + hnXscv*nn;
00783 }
00784 }
00785 else if(theParticle == theProton)
00786 {
00787 if( proj_momentum >= 373.)
00788 {
00789 return GetHadronNucleonXscPDG(aParticle,At,Zt);
00790 }
00791 else if( proj_momentum >= 10.)
00792
00793 {
00794
00795
00796
00797 if(proj_momentum >= 10.)
00798 {
00799 B0 = 7.5;
00800 A0 = 100. - B0*std::log(3.0e7);
00801
00802 xsection = A0 + B0*std::log(proj_energy) - 11
00803 + 103*std::pow(2*0.93827*proj_energy + proj_mass*proj_mass+
00804 0.93827*0.93827,-0.165);
00805 }
00806 xsection *= zz + nn;
00807 }
00808 else
00809 {
00810
00811
00812 if( proj_momentum < 0.73 )
00813 {
00814 hpXscv = 23 + 50*( std::pow( std::log(0.73/proj_momentum), 3.5 ) );
00815 }
00816 else if( proj_momentum < 1.05 )
00817 {
00818 hpXscv = 23 + 40*(std::log(proj_momentum/0.73))*
00819 (std::log(proj_momentum/0.73));
00820 }
00821 else
00822 {
00823 hpXscv = 39.0+
00824 75*(proj_momentum - 1.2)/(std::pow(proj_momentum,3.0) + 0.15);
00825 }
00826
00827
00828 if( proj_momentum < 0.8 )
00829 {
00830 hnXscv = 33+30*std::pow(std::log(proj_momentum/1.3),4.0);
00831 }
00832 else if( proj_momentum < 1.4 )
00833 {
00834 hnXscv = 33+30*std::pow(std::log(proj_momentum/0.95),2.0);
00835 }
00836 else
00837 {
00838 hnXscv = 33.3+
00839 20.8*(std::pow(proj_momentum,2.0)-1.35)/
00840 (std::pow(proj_momentum,2.50)+0.95);
00841 }
00842 xsection = hpXscv*zz + hnXscv*nn;
00843
00844
00845 }
00846
00847 }
00848 else if( theParticle == theAProton )
00849 {
00850
00851
00852
00853
00854
00855
00856 G4double logP = std::log(proj_momentum);
00857
00858 if( proj_momentum <= 1.0 )
00859 {
00860 xsection = zz*(65.55 + 53.84/(proj_momentum+1.e-6) );
00861 }
00862 else
00863 {
00864 xsection = zz*( 41.1 + 77.2*std::pow( proj_momentum, -0.68)
00865 + 0.293*logP*logP - 1.82*logP );
00866 }
00867 if ( nn > 0.)
00868 {
00869 xsection += nn*( 41.9 + 96.2*std::pow( proj_momentum, -0.99) - 0.154*logP);
00870 }
00871 else
00872 {
00873 fInelasticXsc = 38.0 + 38.0*std::pow( proj_momentum, -0.96)
00874 - 0.169*logP*logP;
00875 fInelasticXsc *= millibarn;
00876 }
00877 }
00878 else if( theParticle == thePiPlus )
00879 {
00880 if(proj_momentum < 0.4)
00881 {
00882 G4double Ex3 = 180*std::exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.085/0.085);
00883 hpXscv = Ex3+20.0;
00884 }
00885 else if( proj_momentum < 1.15 )
00886 {
00887 G4double Ex4 = 88*(std::log(proj_momentum/0.75))*(std::log(proj_momentum/0.75));
00888 hpXscv = Ex4+14.0;
00889 }
00890 else if(proj_momentum < 3.5)
00891 {
00892 G4double Ex1 = 3.2*std::exp(-(proj_momentum-2.55)*(proj_momentum-2.55)/0.55/0.55);
00893 G4double Ex2 = 12*std::exp(-(proj_momentum-1.47)*(proj_momentum-1.47)/0.225/0.225);
00894 hpXscv = Ex1+Ex2+27.5;
00895 }
00896 else
00897 {
00898 hpXscv = 10.6+2.*std::log(proj_energy)+25*std::pow(proj_energy,-0.43);
00899 }
00900
00901
00902 if(proj_momentum < 0.37)
00903 {
00904 hnXscv = 28.0 + 40*std::exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.07/0.07);
00905 }
00906 else if(proj_momentum<0.65)
00907 {
00908 hnXscv = 26+110*(std::log(proj_momentum/0.48))*(std::log(proj_momentum/0.48));
00909 }
00910 else if(proj_momentum<1.3)
00911 {
00912 hnXscv = 36.1+
00913 10*std::exp(-(proj_momentum-0.72)*(proj_momentum-0.72)/0.06/0.06)+
00914 24*std::exp(-(proj_momentum-1.015)*(proj_momentum-1.015)/0.075/0.075);
00915 }
00916 else if(proj_momentum<3.0)
00917 {
00918 hnXscv = 36.1+0.079-4.313*std::log(proj_momentum)+
00919 3*std::exp(-(proj_momentum-2.1)*(proj_momentum-2.1)/0.4/0.4)+
00920 1.5*std::exp(-(proj_momentum-1.4)*(proj_momentum-1.4)/0.12/0.12);
00921 }
00922 else
00923 {
00924 hnXscv = 10.6+2*std::log(proj_energy)+30*std::pow(proj_energy,-0.43);
00925 }
00926 xsection = hpXscv*zz + hnXscv*nn;
00927 }
00928 else if(theParticle == thePiMinus)
00929 {
00930
00931
00932 if(proj_momentum < 0.4)
00933 {
00934 G4double Ex3 = 180*std::exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.085/0.085);
00935 hnXscv = Ex3+20.0;
00936 }
00937 else if(proj_momentum < 1.15)
00938 {
00939 G4double Ex4 = 88*(std::log(proj_momentum/0.75))*(std::log(proj_momentum/0.75));
00940 hnXscv = Ex4+14.0;
00941 }
00942 else if(proj_momentum < 3.5)
00943 {
00944 G4double Ex1 = 3.2*std::exp(-(proj_momentum-2.55)*(proj_momentum-2.55)/0.55/0.55);
00945 G4double Ex2 = 12*std::exp(-(proj_momentum-1.47)*(proj_momentum-1.47)/0.225/0.225);
00946 hnXscv = Ex1+Ex2+27.5;
00947 }
00948 else
00949 {
00950 hnXscv = 10.6+2.*std::log(proj_energy)+25*std::pow(proj_energy,-0.43);
00951 }
00952
00953
00954 if(proj_momentum < 0.37)
00955 {
00956 hpXscv = 28.0 + 40*std::exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.07/0.07);
00957 }
00958 else if(proj_momentum<0.65)
00959 {
00960 hpXscv = 26+110*(std::log(proj_momentum/0.48))*(std::log(proj_momentum/0.48));
00961 }
00962 else if(proj_momentum<1.3)
00963 {
00964 hpXscv = 36.1+
00965 10*std::exp(-(proj_momentum-0.72)*(proj_momentum-0.72)/0.06/0.06)+
00966 24*std::exp(-(proj_momentum-1.015)*(proj_momentum-1.015)/0.075/0.075);
00967 }
00968 else if(proj_momentum<3.0)
00969 {
00970 hpXscv = 36.1+0.079-4.313*std::log(proj_momentum)+
00971 3*std::exp(-(proj_momentum-2.1)*(proj_momentum-2.1)/0.4/0.4)+
00972 1.5*std::exp(-(proj_momentum-1.4)*(proj_momentum-1.4)/0.12/0.12);
00973 }
00974 else
00975 {
00976 hpXscv = 10.6+2*std::log(proj_energy)+30*std::pow(proj_energy,-0.43);
00977 }
00978 xsection = hpXscv*zz + hnXscv*nn;
00979 }
00980 else if(theParticle == theKPlus)
00981 {
00982 xsection = zz*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
00983 + 7.14*std::pow(sMand,-eta1) - 13.45*std::pow(sMand,-eta2));
00984
00985 xsection += nn*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
00986 + 5.17*std::pow(sMand,-eta1) - 7.23*std::pow(sMand,-eta2));
00987 }
00988 else if(theParticle == theKMinus)
00989 {
00990 xsection = zz*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
00991 + 7.14*std::pow(sMand,-eta1) + 13.45*std::pow(sMand,-eta2));
00992
00993 xsection += nn*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
00994 + 5.17*std::pow(sMand,-eta1) + 7.23*std::pow(sMand,-eta2));
00995 }
00996 else if(theParticle == theSMinus)
00997 {
00998 xsection = aa*( 35.20 + B*std::pow(std::log(sMand/s0),2.)
00999 - 199.*std::pow(sMand,-eta1) + 264.*std::pow(sMand,-eta2));
01000 }
01001 else if(theParticle == theGamma)
01002 {
01003 xsection = aa*( 0.0 + B*std::pow(std::log(sMand/s0),2.)
01004 + 0.032*std::pow(sMand,-eta1) - 0.0*std::pow(sMand,-eta2));
01005
01006 }
01007 else
01008 {
01009 xsection = zz*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
01010 + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
01011
01012 xsection += nn*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
01013 + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
01014 }
01015 xsection *= millibarn;
01016 return xsection;
01017 }
01018
01019 G4double
01020 G4ComponentGGHadronNucleusXsc::GetKaonNucleonXscVector(const G4DynamicParticle* aParticle,
01021 G4int At, G4int Zt)
01022 {
01023 G4double Tkin, logTkin, xsc, xscP, xscN;
01024 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
01025
01026 G4int Nt = At-Zt;
01027 if (Nt < 0) Nt = 0;
01028
01029 Tkin = aParticle->GetKineticEnergy();
01030
01031 if( Tkin > 70*GeV ) return GetHadronNucleonXscPDG(aParticle,At,Zt);
01032
01033 logTkin = std::log(Tkin);
01034
01035 if( theParticle == theKPlus )
01036 {
01037 xscP = hnXsc->GetKpProtonTotXscVector(logTkin);
01038 xscN = hnXsc->GetKpNeutronTotXscVector(logTkin);
01039 }
01040 else if( theParticle == theKMinus )
01041 {
01042 xscP = hnXsc->GetKmProtonTotXscVector(logTkin);
01043 xscN = hnXsc->GetKmNeutronTotXscVector(logTkin);
01044 }
01045 else
01046 {
01047 xscP = (hnXsc->GetKpProtonTotXscVector(logTkin)+hnXsc->GetKmProtonTotXscVector(logTkin))*0.5;
01048 xscN = (hnXsc->GetKpNeutronTotXscVector(logTkin)+hnXsc->GetKmNeutronTotXscVector(logTkin))*0.5;
01049 }
01050 xsc = xscP*Zt + xscN*Nt;
01051 return xsc;
01052 }
01054
01055
01056
01057 G4double
01058 G4ComponentGGHadronNucleusXsc::GetHNinelasticXsc(const G4DynamicParticle* aParticle,
01059 const G4Element* anElement)
01060 {
01061 G4int At = G4lrint(anElement->GetN());
01062 G4int Zt = G4lrint(anElement->GetZ());
01063
01064 return GetHNinelasticXsc(aParticle, At, Zt);
01065 }
01066
01068
01069
01070
01071 G4double
01072 G4ComponentGGHadronNucleusXsc::GetHNinelasticXsc(const G4DynamicParticle* aParticle,
01073 G4int At, G4int Zt)
01074 {
01075 G4ParticleDefinition* hadron = aParticle->GetDefinition();
01076 G4double sumInelastic;
01077 G4int Nt = At - Zt;
01078 if(Nt < 0) Nt = 0;
01079
01080 if( hadron == theKPlus )
01081 {
01082 sumInelastic = GetHNinelasticXscVU(aParticle, At, Zt);
01083 }
01084 else
01085 {
01086
01087
01088 sumInelastic = G4double(Zt)*GetHadronNucleonXscNS(aParticle, 1, 1);
01089 sumInelastic += G4double(Nt)*GetHadronNucleonXscNS(aParticle, 1, 0);
01090 }
01091 return sumInelastic;
01092 }
01093
01094
01096
01097
01098
01099 G4double
01100 G4ComponentGGHadronNucleusXsc::GetHNinelasticXscVU(const G4DynamicParticle* aParticle,
01101 G4int At, G4int Zt)
01102 {
01103 G4int PDGcode = aParticle->GetDefinition()->GetPDGEncoding();
01104 G4int absPDGcode = std::abs(PDGcode);
01105
01106 G4double Elab = aParticle->GetTotalEnergy();
01107
01108 G4double Plab = aParticle->GetMomentum().mag();
01109
01110
01111 Elab /= GeV;
01112 Plab /= GeV;
01113
01114 G4double LogPlab = std::log( Plab );
01115 G4double sqrLogPlab = LogPlab * LogPlab;
01116
01117
01118
01119 G4double NumberOfTargetProtons = G4double(Zt);
01120 G4double NumberOfTargetNucleons = G4double(At);
01121 G4double NumberOfTargetNeutrons = NumberOfTargetNucleons - NumberOfTargetProtons;
01122
01123 if(NumberOfTargetNeutrons < 0.0) NumberOfTargetNeutrons = 0.0;
01124
01125 G4double Xtotal, Xelastic, Xinelastic;
01126
01127 if( absPDGcode > 1000 )
01128 {
01129 G4double XtotPP = 48.0 + 0. *std::pow(Plab, 0. ) +
01130 0.522*sqrLogPlab - 4.51*LogPlab;
01131
01132 G4double XtotPN = 47.3 + 0. *std::pow(Plab, 0. ) +
01133 0.513*sqrLogPlab - 4.27*LogPlab;
01134
01135 G4double XelPP = 11.9 + 26.9*std::pow(Plab,-1.21) +
01136 0.169*sqrLogPlab - 1.85*LogPlab;
01137
01138 G4double XelPN = 11.9 + 26.9*std::pow(Plab,-1.21) +
01139 0.169*sqrLogPlab - 1.85*LogPlab;
01140
01141 Xtotal = (NumberOfTargetProtons * XtotPP +
01142 NumberOfTargetNeutrons * XtotPN);
01143
01144 Xelastic = (NumberOfTargetProtons * XelPP +
01145 NumberOfTargetNeutrons * XelPN);
01146 }
01147 else if( PDGcode == 211 )
01148 {
01149 G4double XtotPiP = 16.4 + 19.3 *std::pow(Plab,-0.42) +
01150 0.19 *sqrLogPlab - 0.0 *LogPlab;
01151
01152 G4double XtotPiN = 33.0 + 14.0 *std::pow(Plab,-1.36) +
01153 0.456*sqrLogPlab - 4.03*LogPlab;
01154
01155 G4double XelPiP = 0.0 + 11.4*std::pow(Plab,-0.40) +
01156 0.079*sqrLogPlab - 0.0 *LogPlab;
01157
01158 G4double XelPiN = 1.76 + 11.2*std::pow(Plab,-0.64) +
01159 0.043*sqrLogPlab - 0.0 *LogPlab;
01160
01161 Xtotal = ( NumberOfTargetProtons * XtotPiP +
01162 NumberOfTargetNeutrons * XtotPiN );
01163
01164 Xelastic = ( NumberOfTargetProtons * XelPiP +
01165 NumberOfTargetNeutrons * XelPiN );
01166 }
01167 else if( PDGcode == -211 )
01168 {
01169 G4double XtotPiP = 33.0 + 14.0 *std::pow(Plab,-1.36) +
01170 0.456*sqrLogPlab - 4.03*LogPlab;
01171
01172 G4double XtotPiN = 16.4 + 19.3 *std::pow(Plab,-0.42) +
01173 0.19 *sqrLogPlab - 0.0 *LogPlab;
01174
01175 G4double XelPiP = 1.76 + 11.2*std::pow(Plab,-0.64) +
01176 0.043*sqrLogPlab - 0.0 *LogPlab;
01177
01178 G4double XelPiN = 0.0 + 11.4*std::pow(Plab,-0.40) +
01179 0.079*sqrLogPlab - 0.0 *LogPlab;
01180
01181 Xtotal = ( NumberOfTargetProtons * XtotPiP +
01182 NumberOfTargetNeutrons * XtotPiN );
01183
01184 Xelastic = ( NumberOfTargetProtons * XelPiP +
01185 NumberOfTargetNeutrons * XelPiN );
01186 }
01187 else if( PDGcode == 111 )
01188 {
01189 G4double XtotPiP =(16.4 + 19.3 *std::pow(Plab,-0.42) +
01190 0.19 *sqrLogPlab - 0.0 *LogPlab +
01191 33.0 + 14.0 *std::pow(Plab,-1.36) +
01192 0.456*sqrLogPlab - 4.03*LogPlab)/2;
01193
01194 G4double XtotPiN =(33.0 + 14.0 *std::pow(Plab,-1.36) +
01195 0.456*sqrLogPlab - 4.03*LogPlab +
01196 16.4 + 19.3 *std::pow(Plab,-0.42) +
01197 0.19 *sqrLogPlab - 0.0 *LogPlab)/2;
01198
01199 G4double XelPiP =( 0.0 + 11.4*std::pow(Plab,-0.40) +
01200 0.079*sqrLogPlab - 0.0 *LogPlab +
01201 1.76 + 11.2*std::pow(Plab,-0.64) +
01202 0.043*sqrLogPlab - 0.0 *LogPlab)/2;
01203
01204 G4double XelPiN =( 1.76 + 11.2*std::pow(Plab,-0.64) +
01205 0.043*sqrLogPlab - 0.0 *LogPlab +
01206 0.0 + 11.4*std::pow(Plab,-0.40) +
01207 0.079*sqrLogPlab - 0.0 *LogPlab)/2;
01208
01209 Xtotal = ( NumberOfTargetProtons * XtotPiP +
01210 NumberOfTargetNeutrons * XtotPiN );
01211
01212 Xelastic = ( NumberOfTargetProtons * XelPiP +
01213 NumberOfTargetNeutrons * XelPiN );
01214 }
01215 else if( PDGcode == 321 )
01216 {
01217 G4double XtotKP = 18.1 + 0. *std::pow(Plab, 0. ) +
01218 0.26 *sqrLogPlab - 1.0 *LogPlab;
01219 G4double XtotKN = 18.7 + 0. *std::pow(Plab, 0. ) +
01220 0.21 *sqrLogPlab - 0.89*LogPlab;
01221
01222 G4double XelKP = 5.0 + 8.1*std::pow(Plab,-1.8 ) +
01223 0.16 *sqrLogPlab - 1.3 *LogPlab;
01224
01225 G4double XelKN = 7.3 + 0. *std::pow(Plab,-0. ) +
01226 0.29 *sqrLogPlab - 2.4 *LogPlab;
01227
01228 Xtotal = ( NumberOfTargetProtons * XtotKP +
01229 NumberOfTargetNeutrons * XtotKN );
01230
01231 Xelastic = ( NumberOfTargetProtons * XelKP +
01232 NumberOfTargetNeutrons * XelKN );
01233 }
01234 else if( PDGcode ==-321 )
01235 {
01236 G4double XtotKP = 32.1 + 0. *std::pow(Plab, 0. ) +
01237 0.66 *sqrLogPlab - 5.6 *LogPlab;
01238 G4double XtotKN = 25.2 + 0. *std::pow(Plab, 0. ) +
01239 0.38 *sqrLogPlab - 2.9 *LogPlab;
01240
01241 G4double XelKP = 7.3 + 0. *std::pow(Plab,-0. ) +
01242 0.29 *sqrLogPlab - 2.4 *LogPlab;
01243
01244 G4double XelKN = 5.0 + 8.1*std::pow(Plab,-1.8 ) +
01245 0.16 *sqrLogPlab - 1.3 *LogPlab;
01246
01247 Xtotal = ( NumberOfTargetProtons * XtotKP +
01248 NumberOfTargetNeutrons * XtotKN );
01249
01250 Xelastic = ( NumberOfTargetProtons * XelKP +
01251 NumberOfTargetNeutrons * XelKN );
01252 }
01253 else if( PDGcode == 311 )
01254 {
01255 G4double XtotKP = ( 18.1 + 0. *std::pow(Plab, 0. ) +
01256 0.26 *sqrLogPlab - 1.0 *LogPlab +
01257 32.1 + 0. *std::pow(Plab, 0. ) +
01258 0.66 *sqrLogPlab - 5.6 *LogPlab)/2;
01259
01260 G4double XtotKN = ( 18.7 + 0. *std::pow(Plab, 0. ) +
01261 0.21 *sqrLogPlab - 0.89*LogPlab +
01262 25.2 + 0. *std::pow(Plab, 0. ) +
01263 0.38 *sqrLogPlab - 2.9 *LogPlab)/2;
01264
01265 G4double XelKP = ( 5.0 + 8.1*std::pow(Plab,-1.8 )
01266 + 0.16 *sqrLogPlab - 1.3 *LogPlab +
01267 7.3 + 0. *std::pow(Plab,-0. ) +
01268 0.29 *sqrLogPlab - 2.4 *LogPlab)/2;
01269
01270 G4double XelKN = ( 7.3 + 0. *std::pow(Plab,-0. ) +
01271 0.29 *sqrLogPlab - 2.4 *LogPlab +
01272 5.0 + 8.1*std::pow(Plab,-1.8 ) +
01273 0.16 *sqrLogPlab - 1.3 *LogPlab)/2;
01274
01275 Xtotal = ( NumberOfTargetProtons * XtotKP +
01276 NumberOfTargetNeutrons * XtotKN );
01277
01278 Xelastic = ( NumberOfTargetProtons * XelKP +
01279 NumberOfTargetNeutrons * XelKN );
01280 }
01281 else
01282 {
01283 G4double XtotPP = 48.0 + 0. *std::pow(Plab, 0. ) +
01284 0.522*sqrLogPlab - 4.51*LogPlab;
01285
01286 G4double XtotPN = 47.3 + 0. *std::pow(Plab, 0. ) +
01287 0.513*sqrLogPlab - 4.27*LogPlab;
01288
01289 G4double XelPP = 11.9 + 26.9*std::pow(Plab,-1.21) +
01290 0.169*sqrLogPlab - 1.85*LogPlab;
01291 G4double XelPN = 11.9 + 26.9*std::pow(Plab,-1.21) +
01292 0.169*sqrLogPlab - 1.85*LogPlab;
01293
01294 Xtotal = ( NumberOfTargetProtons * XtotPP +
01295 NumberOfTargetNeutrons * XtotPN );
01296
01297 Xelastic = ( NumberOfTargetProtons * XelPP +
01298 NumberOfTargetNeutrons * XelPN );
01299 }
01300 Xinelastic = Xtotal - Xelastic;
01301
01302 if( Xinelastic < 0.) Xinelastic = 0.;
01303
01304 return Xinelastic*= millibarn;
01305 }
01306
01308
01309
01310
01311 G4double
01312 G4ComponentGGHadronNucleusXsc::GetNucleusRadius(const G4DynamicParticle* ,
01313 const G4Element* anElement)
01314 {
01315 G4int At = G4lrint(anElement->GetN());
01316 G4double oneThird = 1.0/3.0;
01317 G4double cubicrAt = std::pow(G4double(At), oneThird);
01318
01319 G4double R;
01320
01321
01322
01323
01324
01325
01326
01327
01328
01329
01330
01331
01332
01333
01334
01335 R = fRadiusConst*cubicrAt;
01336
01337 G4double meanA = 21.;
01338
01339 G4double tauA1 = 40.;
01340 G4double tauA2 = 10.;
01341 G4double tauA3 = 5.;
01342
01343 G4double a1 = 0.85;
01344 G4double b1 = 1. - a1;
01345
01346 G4double b2 = 0.3;
01347 G4double b3 = 4.;
01348
01349 if (At > 20)
01350 {
01351 R *= ( a1 + b1*std::exp( -(At - meanA)/tauA1) );
01352 }
01353 else if (At > 3)
01354 {
01355 R *= ( 1.0 + b2*( 1. - std::exp( (At - meanA)/tauA2) ) );
01356 }
01357 else
01358 {
01359 R *= ( 1.0 + b3*( 1. - std::exp( (At - meanA)/tauA3) ) );
01360 }
01361 return R;
01362
01363 }
01365
01366
01367
01368 G4double
01369 G4ComponentGGHadronNucleusXsc::GetNucleusRadius(G4int At)
01370 {
01371 G4double oneThird = 1.0/3.0;
01372 G4double cubicrAt = std::pow(G4double(At), oneThird);
01373
01374 G4double R;
01375
01376
01377
01378
01379
01380
01381
01382
01383
01384
01385
01386
01387
01388
01389
01390
01391 R = fRadiusConst*cubicrAt;
01392
01393 G4double meanA = 20.;
01394 G4double tauA = 20.;
01395
01396 if (At > 20)
01397 {
01398 R *= ( 0.8 + 0.2*std::exp( -(G4double(At) - meanA)/tauA) );
01399 }
01400 else
01401 {
01402 R *= ( 1.0 + 0.1*( 1. - std::exp( (G4double(At) - meanA)/tauA) ) );
01403 }
01404
01405 return R;
01406 }
01407
01409
01410
01411
01412 G4double G4ComponentGGHadronNucleusXsc::CalculateEcmValue( const G4double mp ,
01413 const G4double mt ,
01414 const G4double Plab )
01415 {
01416 G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
01417 G4double Ecm = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
01418
01419
01420
01421 return Ecm ;
01422 }
01423
01425
01426
01427
01428 G4double G4ComponentGGHadronNucleusXsc::CalcMandelstamS( const G4double mp ,
01429 const G4double mt ,
01430 const G4double Plab )
01431 {
01432 G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
01433 G4double sMand = mp*mp + mt*mt + 2*Elab*mt ;
01434
01435 return sMand;
01436 }
01437
01439
01440
01441
01442 void G4ComponentGGHadronNucleusXsc::CrossSectionDescription(std::ostream& outFile) const
01443 {
01444 outFile << "G4ComponentGGHadronNucleusXsc calculates total, inelastic and\n"
01445 << "elastic cross sections for hadron-nucleus cross sections using\n"
01446 << "the Glauber model with Gribov corrections. It is valid for all\n"
01447 << "targets except hydrogen, and for incident p, pbar, n, sigma-,\n"
01448 << "pi+, pi-, K+, K- and gammas with energies above 3 GeV. This is\n"
01449 << "a cross section component which is to be used to build a cross\n"
01450 << "data set.\n";
01451 }
01452
01453
01455
01456
01457
01458 const G4double G4ComponentGGHadronNucleusXsc::fNeutronBarCorrectionTot[93] = {
01459
01460 1.0, 1.0, 1.118517e+00, 1.082002e+00, 1.116171e+00, 1.078747e+00, 1.061315e+00,
01461 1.058205e+00, 1.082663e+00, 1.068500e+00, 1.076912e+00, 1.083475e+00, 1.079117e+00,
01462 1.071856e+00, 1.071990e+00, 1.073774e+00, 1.079356e+00, 1.081314e+00, 1.082056e+00,
01463 1.090772e+00, 1.096776e+00, 1.095828e+00, 1.097678e+00, 1.099157e+00, 1.103677e+00,
01464 1.105132e+00, 1.109806e+00, 1.110816e+00, 1.117378e+00, 1.115165e+00, 1.115710e+00,
01465 1.111855e+00, 1.110482e+00, 1.110112e+00, 1.106676e+00, 1.108706e+00, 1.105549e+00,
01466 1.106318e+00, 1.106242e+00, 1.107672e+00, 1.107342e+00, 1.108119e+00, 1.106655e+00,
01467 1.102588e+00, 1.096657e+00, 1.092920e+00, 1.086629e+00, 1.083592e+00, 1.076030e+00,
01468 1.083777e+00, 1.089460e+00, 1.086545e+00, 1.079924e+00, 1.082218e+00, 1.077798e+00,
01469 1.077062e+00, 1.072825e+00, 1.072241e+00, 1.072104e+00, 1.072490e+00, 1.069829e+00,
01470 1.070398e+00, 1.065458e+00, 1.064968e+00, 1.060524e+00, 1.060048e+00, 1.057620e+00,
01471 1.056428e+00, 1.055366e+00, 1.055017e+00, 1.052304e+00, 1.051767e+00, 1.049728e+00,
01472 1.048745e+00, 1.047399e+00, 1.045876e+00, 1.042972e+00, 1.041824e+00, 1.039993e+00,
01473 1.039021e+00, 1.036627e+00, 1.034176e+00, 1.032526e+00, 1.033633e+00, 1.036107e+00,
01474 1.037803e+00, 1.031266e+00, 1.032991e+00, 1.033284e+00, 1.035015e+00, 1.033945e+00,
01475 1.037075e+00, 1.034721e+00
01476
01477 };
01478
01479 const G4double G4ComponentGGHadronNucleusXsc::fNeutronBarCorrectionIn[93] = {
01480
01481 1.0, 1.0, 1.167421e+00, 1.156250e+00, 1.205364e+00, 1.154225e+00, 1.120391e+00,
01482 1.124632e+00, 1.129460e+00, 1.107863e+00, 1.102152e+00, 1.104593e+00, 1.100285e+00,
01483 1.098450e+00, 1.092677e+00, 1.101124e+00, 1.106461e+00, 1.115049e+00, 1.123903e+00,
01484 1.126661e+00, 1.131259e+00, 1.133949e+00, 1.134185e+00, 1.133767e+00, 1.132813e+00,
01485 1.131515e+00, 1.130338e+00, 1.134171e+00, 1.139206e+00, 1.141474e+00, 1.142189e+00,
01486 1.140725e+00, 1.140100e+00, 1.139848e+00, 1.137674e+00, 1.138645e+00, 1.136339e+00,
01487 1.136439e+00, 1.135946e+00, 1.136431e+00, 1.135702e+00, 1.135703e+00, 1.134113e+00,
01488 1.131935e+00, 1.128381e+00, 1.126373e+00, 1.122453e+00, 1.120908e+00, 1.115953e+00,
01489 1.115947e+00, 1.114426e+00, 1.111749e+00, 1.106207e+00, 1.107494e+00, 1.103622e+00,
01490 1.102576e+00, 1.098816e+00, 1.097889e+00, 1.097306e+00, 1.097130e+00, 1.094578e+00,
01491 1.094552e+00, 1.090222e+00, 1.089358e+00, 1.085409e+00, 1.084560e+00, 1.082182e+00,
01492 1.080773e+00, 1.079464e+00, 1.078724e+00, 1.076121e+00, 1.075235e+00, 1.073159e+00,
01493 1.071920e+00, 1.070395e+00, 1.069503e+00, 1.067525e+00, 1.066919e+00, 1.065779e+00,
01494 1.065319e+00, 1.063730e+00, 1.062092e+00, 1.061085e+00, 1.059908e+00, 1.059815e+00,
01495 1.059109e+00, 1.051920e+00, 1.051258e+00, 1.049473e+00, 1.048823e+00, 1.045984e+00,
01496 1.046435e+00, 1.042614e+00
01497
01498 };
01499
01500 const G4double G4ComponentGGHadronNucleusXsc::fProtonBarCorrectionTot[93] = {
01501
01502 1.0, 1.0,
01503 1.118515e+00, 1.082000e+00, 1.116169e+00, 1.078745e+00, 1.061313e+00, 1.058203e+00,
01504 1.082661e+00, 1.068498e+00, 1.076910e+00, 1.083474e+00, 1.079115e+00, 1.071854e+00,
01505 1.071988e+00, 1.073772e+00, 1.079355e+00, 1.081312e+00, 1.082054e+00, 1.090770e+00,
01506 1.096774e+00, 1.095827e+00, 1.097677e+00, 1.099156e+00, 1.103676e+00, 1.105130e+00,
01507 1.109805e+00, 1.110814e+00, 1.117377e+00, 1.115163e+00, 1.115708e+00, 1.111853e+00,
01508 1.110480e+00, 1.110111e+00, 1.106674e+00, 1.108705e+00, 1.105548e+00, 1.106317e+00,
01509 1.106241e+00, 1.107671e+00, 1.107341e+00, 1.108118e+00, 1.106654e+00, 1.102586e+00,
01510 1.096655e+00, 1.092918e+00, 1.086628e+00, 1.083590e+00, 1.076028e+00, 1.083776e+00,
01511 1.089458e+00, 1.086543e+00, 1.079923e+00, 1.082216e+00, 1.077797e+00, 1.077061e+00,
01512 1.072824e+00, 1.072239e+00, 1.072103e+00, 1.072488e+00, 1.069828e+00, 1.070396e+00,
01513 1.065456e+00, 1.064966e+00, 1.060523e+00, 1.060047e+00, 1.057618e+00, 1.056427e+00,
01514 1.055365e+00, 1.055016e+00, 1.052303e+00, 1.051766e+00, 1.049727e+00, 1.048743e+00,
01515 1.047397e+00, 1.045875e+00, 1.042971e+00, 1.041823e+00, 1.039992e+00, 1.039019e+00,
01516 1.036626e+00, 1.034175e+00, 1.032525e+00, 1.033632e+00, 1.036106e+00, 1.037802e+00,
01517 1.031265e+00, 1.032990e+00, 1.033283e+00, 1.035014e+00, 1.033944e+00, 1.037074e+00,
01518 1.034720e+00
01519
01520 };
01521
01522 const G4double G4ComponentGGHadronNucleusXsc::fProtonBarCorrectionIn[93] = {
01523
01524 1.0, 1.0,
01525 1.167419e+00, 1.156248e+00, 1.205362e+00, 1.154224e+00, 1.120390e+00, 1.124630e+00,
01526 1.129459e+00, 1.107861e+00, 1.102151e+00, 1.104591e+00, 1.100284e+00, 1.098449e+00,
01527 1.092675e+00, 1.101122e+00, 1.106460e+00, 1.115048e+00, 1.123902e+00, 1.126659e+00,
01528 1.131258e+00, 1.133948e+00, 1.134183e+00, 1.133766e+00, 1.132812e+00, 1.131514e+00,
01529 1.130337e+00, 1.134170e+00, 1.139205e+00, 1.141472e+00, 1.142188e+00, 1.140724e+00,
01530 1.140099e+00, 1.139847e+00, 1.137672e+00, 1.138644e+00, 1.136338e+00, 1.136438e+00,
01531 1.135945e+00, 1.136429e+00, 1.135701e+00, 1.135702e+00, 1.134112e+00, 1.131934e+00,
01532 1.128380e+00, 1.126371e+00, 1.122452e+00, 1.120907e+00, 1.115952e+00, 1.115946e+00,
01533 1.114425e+00, 1.111748e+00, 1.106205e+00, 1.107493e+00, 1.103621e+00, 1.102575e+00,
01534 1.098815e+00, 1.097888e+00, 1.097305e+00, 1.097129e+00, 1.094577e+00, 1.094551e+00,
01535 1.090221e+00, 1.089357e+00, 1.085408e+00, 1.084559e+00, 1.082181e+00, 1.080772e+00,
01536 1.079463e+00, 1.078723e+00, 1.076120e+00, 1.075234e+00, 1.073158e+00, 1.071919e+00,
01537 1.070394e+00, 1.069502e+00, 1.067524e+00, 1.066918e+00, 1.065778e+00, 1.065318e+00,
01538 1.063729e+00, 1.062091e+00, 1.061084e+00, 1.059907e+00, 1.059814e+00, 1.059108e+00,
01539 1.051919e+00, 1.051257e+00, 1.049472e+00, 1.048822e+00, 1.045983e+00, 1.046434e+00,
01540 1.042613e+00
01541
01542 };
01543
01544
01545 const G4double G4ComponentGGHadronNucleusXsc::fPionPlusBarCorrectionTot[93] = {
01546
01547 1.0, 1.0,
01548 1.075927e+00, 1.074407e+00, 1.126098e+00, 1.100127e+00, 1.089742e+00, 1.083536e+00,
01549 1.089988e+00, 1.103566e+00, 1.096922e+00, 1.126573e+00, 1.132734e+00, 1.136512e+00,
01550 1.136629e+00, 1.133086e+00, 1.132428e+00, 1.129299e+00, 1.125622e+00, 1.126992e+00,
01551 1.127840e+00, 1.162670e+00, 1.160392e+00, 1.157864e+00, 1.157227e+00, 1.154627e+00,
01552 1.192555e+00, 1.197243e+00, 1.197911e+00, 1.200326e+00, 1.220053e+00, 1.215019e+00,
01553 1.211703e+00, 1.209080e+00, 1.204248e+00, 1.203328e+00, 1.198671e+00, 1.196840e+00,
01554 1.194392e+00, 1.193037e+00, 1.190408e+00, 1.188583e+00, 1.206127e+00, 1.210028e+00,
01555 1.206434e+00, 1.204456e+00, 1.200547e+00, 1.199058e+00, 1.200174e+00, 1.200276e+00,
01556 1.198912e+00, 1.213048e+00, 1.207160e+00, 1.208020e+00, 1.203814e+00, 1.202380e+00,
01557 1.198306e+00, 1.197002e+00, 1.196027e+00, 1.195449e+00, 1.192563e+00, 1.192135e+00,
01558 1.187556e+00, 1.186308e+00, 1.182124e+00, 1.180900e+00, 1.178224e+00, 1.176471e+00,
01559 1.174811e+00, 1.173702e+00, 1.170827e+00, 1.169581e+00, 1.167205e+00, 1.165626e+00,
01560 1.180244e+00, 1.177626e+00, 1.175121e+00, 1.173903e+00, 1.172192e+00, 1.171128e+00,
01561 1.168997e+00, 1.166826e+00, 1.164130e+00, 1.165412e+00, 1.165504e+00, 1.165020e+00,
01562 1.158462e+00, 1.158014e+00, 1.156519e+00, 1.156081e+00, 1.153602e+00, 1.154190e+00,
01563 1.152974e+00
01564
01565 };
01566
01567 const G4double G4ComponentGGHadronNucleusXsc::fPionPlusBarCorrectionIn[93] = {
01568
01569 1.0, 1.0,
01570 1.140246e+00, 1.097872e+00, 1.104301e+00, 1.068722e+00, 1.044495e+00, 1.062622e+00,
01571 1.047987e+00, 1.037032e+00, 1.035686e+00, 1.042870e+00, 1.052222e+00, 1.065100e+00,
01572 1.070480e+00, 1.078286e+00, 1.081488e+00, 1.089713e+00, 1.099105e+00, 1.098003e+00,
01573 1.102175e+00, 1.117707e+00, 1.121734e+00, 1.125229e+00, 1.126457e+00, 1.128905e+00,
01574 1.137312e+00, 1.126263e+00, 1.126459e+00, 1.115191e+00, 1.116986e+00, 1.117184e+00,
01575 1.117037e+00, 1.116777e+00, 1.115858e+00, 1.115745e+00, 1.114489e+00, 1.113993e+00,
01576 1.113226e+00, 1.112818e+00, 1.111890e+00, 1.111238e+00, 1.111209e+00, 1.111775e+00,
01577 1.110256e+00, 1.109414e+00, 1.107647e+00, 1.106980e+00, 1.106096e+00, 1.107331e+00,
01578 1.107849e+00, 1.106407e+00, 1.103426e+00, 1.103896e+00, 1.101756e+00, 1.101031e+00,
01579 1.098915e+00, 1.098260e+00, 1.097768e+00, 1.097487e+00, 1.095964e+00, 1.095773e+00,
01580 1.093348e+00, 1.092687e+00, 1.090465e+00, 1.089821e+00, 1.088394e+00, 1.087462e+00,
01581 1.086571e+00, 1.085997e+00, 1.084451e+00, 1.083798e+00, 1.082513e+00, 1.081670e+00,
01582 1.080735e+00, 1.075659e+00, 1.074341e+00, 1.073689e+00, 1.072787e+00, 1.072237e+00,
01583 1.071107e+00, 1.069955e+00, 1.064856e+00, 1.065873e+00, 1.065938e+00, 1.065694e+00,
01584 1.062192e+00, 1.061967e+00, 1.061180e+00, 1.060960e+00, 1.059646e+00, 1.059975e+00,
01585 1.059658e+00
01586
01587 };
01588
01589
01590 const G4double G4ComponentGGHadronNucleusXsc::fPionMinusBarCorrectionTot[93] = {
01591
01592 1.0, 1.0,
01593 1.075927e+00, 1.077959e+00, 1.129145e+00, 1.102088e+00, 1.089765e+00, 1.083542e+00,
01594 1.089995e+00, 1.104895e+00, 1.097154e+00, 1.127663e+00, 1.133063e+00, 1.137425e+00,
01595 1.136724e+00, 1.133859e+00, 1.132498e+00, 1.130276e+00, 1.127896e+00, 1.127656e+00,
01596 1.127905e+00, 1.164210e+00, 1.162259e+00, 1.160075e+00, 1.158978e+00, 1.156649e+00,
01597 1.194157e+00, 1.199177e+00, 1.198983e+00, 1.202325e+00, 1.221967e+00, 1.217548e+00,
01598 1.214389e+00, 1.211760e+00, 1.207335e+00, 1.206081e+00, 1.201766e+00, 1.199779e+00,
01599 1.197283e+00, 1.195706e+00, 1.193071e+00, 1.191115e+00, 1.208838e+00, 1.212681e+00,
01600 1.209235e+00, 1.207163e+00, 1.203451e+00, 1.201807e+00, 1.203283e+00, 1.203388e+00,
01601 1.202244e+00, 1.216509e+00, 1.211066e+00, 1.211504e+00, 1.207539e+00, 1.205991e+00,
01602 1.202143e+00, 1.200724e+00, 1.199595e+00, 1.198815e+00, 1.196025e+00, 1.195390e+00,
01603 1.191137e+00, 1.189791e+00, 1.185888e+00, 1.184575e+00, 1.181996e+00, 1.180229e+00,
01604 1.178545e+00, 1.177355e+00, 1.174616e+00, 1.173312e+00, 1.171016e+00, 1.169424e+00,
01605 1.184120e+00, 1.181478e+00, 1.179085e+00, 1.177817e+00, 1.176124e+00, 1.175003e+00,
01606 1.172947e+00, 1.170858e+00, 1.168170e+00, 1.169397e+00, 1.169304e+00, 1.168706e+00,
01607 1.162774e+00, 1.162217e+00, 1.160740e+00, 1.160196e+00, 1.157857e+00, 1.158220e+00,
01608 1.157267e+00
01609 };
01610
01611
01612 const G4double G4ComponentGGHadronNucleusXsc::fPionMinusBarCorrectionIn[93] = {
01613
01614 1.0, 1.0,
01615 1.140246e+00, 1.100898e+00, 1.106773e+00, 1.070289e+00, 1.044514e+00, 1.062628e+00,
01616 1.047992e+00, 1.038041e+00, 1.035862e+00, 1.043679e+00, 1.052466e+00, 1.065780e+00,
01617 1.070551e+00, 1.078869e+00, 1.081541e+00, 1.090455e+00, 1.100847e+00, 1.098511e+00,
01618 1.102226e+00, 1.118865e+00, 1.123143e+00, 1.126904e+00, 1.127785e+00, 1.130444e+00,
01619 1.138502e+00, 1.127678e+00, 1.127244e+00, 1.116634e+00, 1.118347e+00, 1.118988e+00,
01620 1.118957e+00, 1.118696e+00, 1.118074e+00, 1.117722e+00, 1.116717e+00, 1.116111e+00,
01621 1.115311e+00, 1.114745e+00, 1.113814e+00, 1.113069e+00, 1.113141e+00, 1.113660e+00,
01622 1.112249e+00, 1.111343e+00, 1.109718e+00, 1.108942e+00, 1.108310e+00, 1.109549e+00,
01623 1.110227e+00, 1.108846e+00, 1.106183e+00, 1.106354e+00, 1.104388e+00, 1.103583e+00,
01624 1.101632e+00, 1.100896e+00, 1.100296e+00, 1.099873e+00, 1.098420e+00, 1.098082e+00,
01625 1.095892e+00, 1.095162e+00, 1.093144e+00, 1.092438e+00, 1.091083e+00, 1.090142e+00,
01626 1.089236e+00, 1.088604e+00, 1.087159e+00, 1.086465e+00, 1.085239e+00, 1.084388e+00,
01627 1.083473e+00, 1.078373e+00, 1.077136e+00, 1.076450e+00, 1.075561e+00, 1.074973e+00,
01628 1.073898e+00, 1.072806e+00, 1.067706e+00, 1.068684e+00, 1.068618e+00, 1.068294e+00,
01629 1.065241e+00, 1.064939e+00, 1.064166e+00, 1.063872e+00, 1.062659e+00, 1.062828e+00,
01630 1.062699e+00
01631
01632 };
01633
01634
01635
01636