Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4ComponentGGHadronNucleusXsc Class Reference

#include <G4ComponentGGHadronNucleusXsc.hh>

Inheritance diagram for G4ComponentGGHadronNucleusXsc:
G4VComponentCrossSection

Public Member Functions

 G4ComponentGGHadronNucleusXsc ()
 
virtual ~G4ComponentGGHadronNucleusXsc ()
 
virtual G4double GetTotalIsotopeCrossSection (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
 
virtual G4double GetTotalElementCrossSection (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
 
virtual G4double GetInelasticIsotopeCrossSection (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
 
virtual G4double GetProductionIsotopeCrossSection (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
 
virtual G4double GetInelasticElementCrossSection (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
 
virtual G4double GetProductionElementCrossSection (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
 
virtual G4double GetElasticElementCrossSection (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
 
virtual G4double GetElasticIsotopeCrossSection (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
 
virtual G4double ComputeQuasiElasticRatio (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
 
G4bool IsIsoApplicable (const G4DynamicParticle *aDP, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
 
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
G4double GetRatioSD (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double GetRatioQE (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double GetHadronNucleonXsc (const G4DynamicParticle *, const G4Element *)
 
G4double GetHadronNucleonXsc (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double GetHadronNucleonXscPDG (const G4DynamicParticle *, const G4Element *)
 
G4double GetHadronNucleonXscPDG (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double GetHadronNucleonXscNS (const G4DynamicParticle *, const G4Element *)
 
G4double GetHadronNucleonXscNS (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double GetKaonNucleonXscVector (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double GetHNinelasticXsc (const G4DynamicParticle *, const G4Element *)
 
G4double GetHNinelasticXsc (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double GetHNinelasticXscVU (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double CalculateEcmValue (const G4double, const G4double, const G4double)
 
G4double CalcMandelstamS (const G4double, const G4double, const G4double)
 
G4double GetNucleusRadius (const G4DynamicParticle *, const G4Element *)
 
G4double GetNucleusRadius (G4int At)
 
virtual void CrossSectionDescription (std::ostream &) const
 
G4double GetElasticGlauberGribov (const G4DynamicParticle *, G4int Z, G4int A)
 
G4double GetInelasticGlauberGribov (const G4DynamicParticle *, G4int Z, G4int A)
 
G4double GetTotalGlauberGribovXsc ()
 
G4double GetElasticGlauberGribovXsc ()
 
G4double GetInelasticGlauberGribovXsc ()
 
G4double GetProductionGlauberGribovXsc ()
 
G4double GetDiffractionGlauberGribovXsc ()
 
G4double GetRadiusConst ()
 
G4double GetParticleBarCorTot (const G4ParticleDefinition *theParticle, G4int Z)
 
G4double GetParticleBarCorIn (const G4ParticleDefinition *theParticle, G4int Z)
 
void SetEnergyLowerLimit (G4double E)
 
- Public Member Functions inherited from G4VComponentCrossSection
 G4VComponentCrossSection (const G4String &nam="")
 
virtual ~G4VComponentCrossSection ()
 
G4double GetTotalElementCrossSection (const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)
 
G4double GetInelasticElementCrossSection (const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)
 
G4double GetElasticElementCrossSection (const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void Description () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
const G4StringGetName () const
 

Detailed Description

Definition at line 51 of file G4ComponentGGHadronNucleusXsc.hh.

Constructor & Destructor Documentation

G4ComponentGGHadronNucleusXsc::G4ComponentGGHadronNucleusXsc ( )

Definition at line 44 of file G4ComponentGGHadronNucleusXsc.cc.

References G4Alpha::Alpha(), G4AntiLambda::AntiLambda(), G4AntiNeutron::AntiNeutron(), G4AntiOmegaMinus::AntiOmegaMinus(), G4AntiProton::AntiProton(), G4AntiSigmaMinus::AntiSigmaMinus(), G4AntiSigmaPlus::AntiSigmaPlus(), G4AntiSigmaZero::AntiSigmaZero(), G4AntiXiMinus::AntiXiMinus(), G4AntiXiZero::AntiXiZero(), G4Deuteron::Deuteron(), G4Gamma::Gamma(), G4He3::He3(), G4KaonMinus::KaonMinus(), G4KaonPlus::KaonPlus(), G4KaonZeroLong::KaonZeroLong(), G4KaonZeroShort::KaonZeroShort(), G4Lambda::Lambda(), G4Neutron::Neutron(), G4OmegaMinus::OmegaMinus(), G4PionMinus::PionMinus(), G4PionPlus::PionPlus(), G4PionZero::PionZero(), G4Proton::Proton(), G4SigmaMinus::SigmaMinus(), G4SigmaPlus::SigmaPlus(), G4SigmaZero::SigmaZero(), G4Triton::Triton(), G4XiMinus::XiMinus(), and G4XiZero::XiZero().

45  : G4VComponentCrossSection("Glauber-Gribov"),
46 // fUpperLimit(100000*GeV),
47  fLowerLimit(10.*MeV),// fLowerLimit(3*GeV),
48  fRadiusConst(1.08*fermi), // 1.1, 1.3 ?
49  fTotalXsc(0.0), fElasticXsc(0.0), fInelasticXsc(0.0), fProductionXsc(0.0),
50  fDiffractionXsc(0.0)
51 // , fHadronNucleonXsc(0.0)
52 {
53  theGamma = G4Gamma::Gamma();
54  theProton = G4Proton::Proton();
55  theNeutron = G4Neutron::Neutron();
56  theAProton = G4AntiProton::AntiProton();
57  theANeutron = G4AntiNeutron::AntiNeutron();
58  thePiPlus = G4PionPlus::PionPlus();
59  thePiMinus = G4PionMinus::PionMinus();
60  thePiZero = G4PionZero::PionZero();
61  theKPlus = G4KaonPlus::KaonPlus();
62  theKMinus = G4KaonMinus::KaonMinus();
65  theL = G4Lambda::Lambda();
66  theAntiL = G4AntiLambda::AntiLambda();
67  theSPlus = G4SigmaPlus::SigmaPlus();
68  theASPlus = G4AntiSigmaPlus::AntiSigmaPlus();
69  theSMinus = G4SigmaMinus::SigmaMinus();
70  theASMinus = G4AntiSigmaMinus::AntiSigmaMinus();
71  theS0 = G4SigmaZero::SigmaZero();
73  theXiMinus = G4XiMinus::XiMinus();
74  theXi0 = G4XiZero::XiZero();
75  theAXiMinus = G4AntiXiMinus::AntiXiMinus();
76  theAXi0 = G4AntiXiZero::AntiXiZero();
77  theOmega = G4OmegaMinus::OmegaMinus();
79  theD = G4Deuteron::Deuteron();
80  theT = G4Triton::Triton();
81  theA = G4Alpha::Alpha();
82  theHe3 = G4He3::He3();
83 
84  hnXsc = new G4HadronNucleonXsc();
85 }
static G4AntiOmegaMinus * AntiOmegaMinus()
static G4OmegaMinus * OmegaMinus()
static G4KaonZeroLong * KaonZeroLong()
static G4AntiSigmaPlus * AntiSigmaPlus()
static G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:99
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
static G4AntiSigmaMinus * AntiSigmaMinus()
static G4XiZero * XiZero()
Definition: G4XiZero.cc:106
static G4KaonZeroShort * KaonZeroShort()
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4XiMinus * XiMinus()
Definition: G4XiMinus.cc:106
static G4AntiXiMinus * AntiXiMinus()
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4PionZero * PionZero()
Definition: G4PionZero.cc:104
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4SigmaMinus * SigmaMinus()
static G4AntiLambda * AntiLambda()
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4AntiSigmaZero * AntiSigmaZero()
G4VComponentCrossSection(const G4String &nam="")
static G4AntiXiZero * AntiXiZero()
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
static G4Lambda * Lambda()
Definition: G4Lambda.cc:108
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
static G4He3 * He3()
Definition: G4He3.cc:94
static G4AntiNeutron * AntiNeutron()
G4ComponentGGHadronNucleusXsc::~G4ComponentGGHadronNucleusXsc ( )
virtual

Definition at line 91 of file G4ComponentGGHadronNucleusXsc.cc.

92 {
93  if (hnXsc) delete hnXsc;
94 }

Member Function Documentation

G4double G4ComponentGGHadronNucleusXsc::CalcMandelstamS ( const G4double  mp,
const G4double  mt,
const G4double  Plab 
)

Definition at line 1469 of file G4ComponentGGHadronNucleusXsc.cc.

Referenced by GetHadronNucleonXsc(), GetHadronNucleonXscNS(), and GetHadronNucleonXscPDG().

1472 {
1473  G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
1474  G4double sMand = mp*mp + mt*mt + 2*Elab*mt ;
1475 
1476  return sMand;
1477 }
double G4double
Definition: G4Types.hh:76
G4double G4ComponentGGHadronNucleusXsc::CalculateEcmValue ( const G4double  mp,
const G4double  mt,
const G4double  Plab 
)

Definition at line 1453 of file G4ComponentGGHadronNucleusXsc.cc.

1456 {
1457  G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
1458  G4double Ecm = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
1459  // G4double Pcm = Plab * mt / Ecm;
1460  // G4double KEcm = std::sqrt ( Pcm * Pcm + mp * mp ) - mp;
1461 
1462  return Ecm ; // KEcm;
1463 }
double G4double
Definition: G4Types.hh:76
G4double G4ComponentGGHadronNucleusXsc::ComputeQuasiElasticRatio ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4int  A 
)
virtual

Reimplemented from G4VComponentCrossSection.

Definition at line 210 of file G4ComponentGGHadronNucleusXsc.cc.

References GetIsoCrossSection().

213 {
214  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
215  kinEnergy);
216  fTotalXsc = GetIsoCrossSection(aDP, Z, A);
217  delete aDP;
218  G4double ratio = 0.;
219 
220  if(fInelasticXsc > 0.)
221  {
222  ratio = (fInelasticXsc - fProductionXsc)/fInelasticXsc;
223  if(ratio < 0.) ratio = 0.;
224  }
225  return ratio;
226 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
double G4double
Definition: G4Types.hh:76
G4ThreeVector G4ParticleMomentum
void G4ComponentGGHadronNucleusXsc::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Definition at line 1483 of file G4ComponentGGHadronNucleusXsc.cc.

1484 {
1485  outFile << "G4ComponentGGHadronNucleusXsc calculates total, inelastic and\n"
1486  << "elastic cross sections for hadron-nucleus cross sections using\n"
1487  << "the Glauber model with Gribov corrections. It is valid for all\n"
1488  << "targets except hydrogen, and for incident p, pbar, n, sigma-,\n"
1489  << "pi+, pi-, K+, K- and gammas with energies above 3 GeV. This is\n"
1490  << "a cross section component which is to be used to build a cross\n"
1491  << "data set.\n";
1492 }
std::ofstream outFile
Definition: GammaRayTel.cc:68
G4double G4ComponentGGHadronNucleusXsc::GetDiffractionGlauberGribovXsc ( )
inline

Definition at line 150 of file G4ComponentGGHadronNucleusXsc.hh.

150 { return fDiffractionXsc; };
G4double G4ComponentGGHadronNucleusXsc::GetElasticElementCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4double  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 182 of file G4ComponentGGHadronNucleusXsc.cc.

References GetIsoCrossSection().

185 {
186  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
187  kinEnergy);
188  fTotalXsc = GetIsoCrossSection(aDP, Z, G4int(A));
189  delete aDP;
190 
191  return fElasticXsc;
192 }
int G4int
Definition: G4Types.hh:78
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4ThreeVector G4ParticleMomentum
G4double G4ComponentGGHadronNucleusXsc::GetElasticGlauberGribov ( const G4DynamicParticle dp,
G4int  Z,
G4int  A 
)
inline

Definition at line 220 of file G4ComponentGGHadronNucleusXsc.hh.

References GetIsoCrossSection().

222 {
223  GetIsoCrossSection(dp, Z, A);
224  return fElasticXsc;
225 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4double G4ComponentGGHadronNucleusXsc::GetElasticGlauberGribovXsc ( )
inline

Definition at line 147 of file G4ComponentGGHadronNucleusXsc.hh.

147 { return fElasticXsc; };
G4double G4ComponentGGHadronNucleusXsc::GetElasticIsotopeCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4int  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 196 of file G4ComponentGGHadronNucleusXsc.cc.

References GetIsoCrossSection().

199 {
200  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
201  kinEnergy);
202  fTotalXsc = GetIsoCrossSection(aDP, Z, A);
203  delete aDP;
204 
205  return fElasticXsc;
206 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4ThreeVector G4ParticleMomentum
G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXsc ( const G4DynamicParticle aParticle,
const G4Element anElement 
)

Definition at line 484 of file G4ComponentGGHadronNucleusXsc.cc.

References G4lrint(), G4Element::GetN(), and G4Element::GetZ().

486 {
487  G4int At = G4lrint(anElement->GetN()); // number of nucleons
488  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
489 
490  return GetHadronNucleonXsc(aParticle, At, Zt);
491 }
G4double GetN() const
Definition: G4Element.hh:134
G4double GetZ() const
Definition: G4Element.hh:131
int G4int
Definition: G4Types.hh:78
G4double GetHadronNucleonXsc(const G4DynamicParticle *, const G4Element *)
int G4lrint(double ad)
Definition: templates.hh:163
G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXsc ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 501 of file G4ComponentGGHadronNucleusXsc.cc.

References CalcMandelstamS(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetMass(), G4DynamicParticle::GetMomentum(), python.hepunit::GeV, CLHEP::Hep3Vector::mag(), and python.hepunit::millibarn.

503 {
504  G4double xsection;
505 
506  //G4double targ_mass = G4NucleiProperties::GetNuclearMass(At, Zt);
507 
508  G4double targ_mass = 0.939*GeV; // ~mean neutron and proton ???
509 
510  G4double proj_mass = aParticle->GetMass();
511  G4double proj_momentum = aParticle->GetMomentum().mag();
512  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
513 
514  sMand /= GeV*GeV; // in GeV for parametrisation
515  proj_momentum /= GeV;
516 
517  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
518 
519  G4double aa = At;
520 
521  if(theParticle == theGamma)
522  {
523  xsection = aa*(0.0677*std::pow(sMand,0.0808) + 0.129*std::pow(sMand,-0.4525));
524  }
525  else if(theParticle == theNeutron) // as proton ???
526  {
527  xsection = aa*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
528  }
529  else if(theParticle == theProton)
530  {
531  xsection = aa*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
532  // xsection = At*( 49.51*std::pow(sMand,-0.097) + 0.314*std::log(sMand)*std::log(sMand) );
533  // xsection = At*( 38.4 + 0.85*std::abs(std::pow(log(sMand),1.47)) );
534  }
535  else if(theParticle == theAProton)
536  {
537  xsection = aa*( 21.70*std::pow(sMand,0.0808) + 98.39*std::pow(sMand,-0.4525));
538  }
539  else if(theParticle == thePiPlus)
540  {
541  xsection = aa*(13.63*std::pow(sMand,0.0808) + 27.56*std::pow(sMand,-0.4525));
542  }
543  else if(theParticle == thePiMinus)
544  {
545  // xsection = At*( 55.2*std::pow(sMand,-0.255) + 0.346*std::log(sMand)*std::log(sMand) );
546  xsection = aa*(13.63*std::pow(sMand,0.0808) + 36.02*std::pow(sMand,-0.4525));
547  }
548  else if(theParticle == theKPlus)
549  {
550  xsection = aa*(11.82*std::pow(sMand,0.0808) + 8.15*std::pow(sMand,-0.4525));
551  }
552  else if(theParticle == theKMinus)
553  {
554  xsection = aa*(11.82*std::pow(sMand,0.0808) + 26.36*std::pow(sMand,-0.4525));
555  }
556  else // as proton ???
557  {
558  xsection = aa*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
559  }
560  xsection *= millibarn;
561  return xsection;
562 }
G4double CalcMandelstamS(const G4double, const G4double, const G4double)
G4ParticleDefinition * GetDefinition() const
int millibarn
Definition: hepunit.py:40
G4double GetMass() const
double G4double
Definition: G4Types.hh:76
double mag() const
G4ThreeVector GetMomentum() const
G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscNS ( const G4DynamicParticle aParticle,
const G4Element anElement 
)

Definition at line 705 of file G4ComponentGGHadronNucleusXsc.cc.

References G4lrint(), G4Element::GetN(), and G4Element::GetZ().

Referenced by GetHNinelasticXsc(), GetIsoCrossSection(), GetRatioQE(), and GetRatioSD().

707 {
708  G4int At = G4lrint(anElement->GetN()); // number of nucleons
709  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
710 
711  return GetHadronNucleonXscNS(aParticle, At, Zt);
712 }
G4double GetN() const
Definition: G4Element.hh:134
G4double GetZ() const
Definition: G4Element.hh:131
int G4int
Definition: G4Types.hh:78
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
int G4lrint(double ad)
Definition: templates.hh:163
G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscNS ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 723 of file G4ComponentGGHadronNucleusXsc.cc.

References CalcMandelstamS(), G4DynamicParticle::GetDefinition(), GetHadronNucleonXscPDG(), G4DynamicParticle::GetMass(), G4DynamicParticle::GetMomentum(), G4ParticleTable::GetParticleTable(), G4DynamicParticle::GetTotalEnergy(), python.hepunit::GeV, CLHEP::Hep3Vector::mag(), python.hepunit::millibarn, G4InuclParticleNames::nn, and G4InuclParticleNames::s0.

725 {
726  G4double xsection(0);
727  // G4double Delta; DHW 19 May 2011: variable set but not used
728  G4double A0, B0;
729  G4double hpXscv(0);
730  G4double hnXscv(0);
731 
732  G4int Nt = At-Zt; // number of neutrons
733  if (Nt < 0) Nt = 0;
734 
735  G4double aa = At;
736  G4double zz = Zt;
737  G4double nn = Nt;
738 
740  GetIonTable()->GetIonMass(Zt, At);
741 
742  targ_mass = 0.939*GeV; // ~mean neutron and proton ???
743 
744  G4double proj_mass = aParticle->GetMass();
745  G4double proj_energy = aParticle->GetTotalEnergy();
746  G4double proj_momentum = aParticle->GetMomentum().mag();
747 
748  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
749 
750  sMand /= GeV*GeV; // in GeV for parametrisation
751  proj_momentum /= GeV;
752  proj_energy /= GeV;
753  proj_mass /= GeV;
754 
755  // General PDG fit constants
756 
757  G4double s0 = 5.38*5.38; // in Gev^2
758  G4double eta1 = 0.458;
759  G4double eta2 = 0.458;
760  G4double B = 0.308;
761 
762 
763  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
764 
765 
766  if(theParticle == theNeutron)
767  {
768  if( proj_momentum >= 373.)
769  {
770  return GetHadronNucleonXscPDG(aParticle,At,Zt);
771  }
772  else if( proj_momentum >= 10.)
773  // if( proj_momentum >= 2.)
774  {
775  // Delta = 1.; // DHW 19 May 2011: variable set but not used
776  // if( proj_energy < 40. ) Delta = 0.916+0.0021*proj_energy;
777 
778  if(proj_momentum >= 10.)
779  {
780  B0 = 7.5;
781  A0 = 100. - B0*std::log(3.0e7);
782 
783  xsection = A0 + B0*std::log(proj_energy) - 11
784  + 103*std::pow(2*0.93827*proj_energy + proj_mass*proj_mass+
785  0.93827*0.93827,-0.165); // mb
786  }
787  xsection *= zz + nn;
788  }
789  else
790  {
791  // nn to be pp
792 
793  if( proj_momentum < 0.73 )
794  {
795  hnXscv = 23 + 50*( std::pow( std::log(0.73/proj_momentum), 3.5 ) );
796  }
797  else if( proj_momentum < 1.05 )
798  {
799  hnXscv = 23 + 40*(std::log(proj_momentum/0.73))*
800  (std::log(proj_momentum/0.73));
801  }
802  else // if( proj_momentum < 10. )
803  {
804  hnXscv = 39.0+
805  75*(proj_momentum - 1.2)/(std::pow(proj_momentum,3.0) + 0.15);
806  }
807  // pn to be np
808 
809  if( proj_momentum < 0.8 )
810  {
811  hpXscv = 33+30*std::pow(std::log(proj_momentum/1.3),4.0);
812  }
813  else if( proj_momentum < 1.4 )
814  {
815  hpXscv = 33+30*std::pow(std::log(proj_momentum/0.95),2.0);
816  }
817  else // if( proj_momentum < 10. )
818  {
819  hpXscv = 33.3+
820  20.8*(std::pow(proj_momentum,2.0)-1.35)/
821  (std::pow(proj_momentum,2.50)+0.95);
822  }
823  xsection = hpXscv*zz + hnXscv*nn;
824  }
825  }
826  else if(theParticle == theProton)
827  {
828  if( proj_momentum >= 373.)
829  {
830  return GetHadronNucleonXscPDG(aParticle,At,Zt);
831  }
832  else if( proj_momentum >= 10.)
833  // if( proj_momentum >= 2.)
834  {
835  // Delta = 1.; DHW 19 May 2011: variable set but not used
836  // if( proj_energy < 40. ) Delta = 0.916+0.0021*proj_energy;
837 
838  if(proj_momentum >= 10.)
839  {
840  B0 = 7.5;
841  A0 = 100. - B0*std::log(3.0e7);
842 
843  xsection = A0 + B0*std::log(proj_energy) - 11
844  + 103*std::pow(2*0.93827*proj_energy + proj_mass*proj_mass+
845  0.93827*0.93827,-0.165); // mb
846  }
847  xsection *= zz + nn;
848  }
849  else
850  {
851  // pp
852 
853  if( proj_momentum < 0.73 )
854  {
855  hpXscv = 23 + 50*( std::pow( std::log(0.73/proj_momentum), 3.5 ) );
856  }
857  else if( proj_momentum < 1.05 )
858  {
859  hpXscv = 23 + 40*(std::log(proj_momentum/0.73))*
860  (std::log(proj_momentum/0.73));
861  }
862  else // if( proj_momentum < 10. )
863  {
864  hpXscv = 39.0+
865  75*(proj_momentum - 1.2)/(std::pow(proj_momentum,3.0) + 0.15);
866  }
867  // pn to be np
868 
869  if( proj_momentum < 0.8 )
870  {
871  hnXscv = 33+30*std::pow(std::log(proj_momentum/1.3),4.0);
872  }
873  else if( proj_momentum < 1.4 )
874  {
875  hnXscv = 33+30*std::pow(std::log(proj_momentum/0.95),2.0);
876  }
877  else // if( proj_momentum < 10. )
878  {
879  hnXscv = 33.3+
880  20.8*(std::pow(proj_momentum,2.0)-1.35)/
881  (std::pow(proj_momentum,2.50)+0.95);
882  }
883  xsection = hpXscv*zz + hnXscv*nn;
884  // xsection = hpXscv*(Zt + Nt);
885  // xsection = hnXscv*(Zt + Nt);
886  }
887  // xsection *= 0.95;
888  }
889  else if( theParticle == theAProton )
890  {
891  // xsection = Zt*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
892  // + 42.53*std::pow(sMand,-eta1) + 33.34*std::pow(sMand,-eta2));
893 
894  // xsection += Nt*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
895  // + 40.15*std::pow(sMand,-eta1) + 30.*std::pow(sMand,-eta2));
896 
897  G4double logP = std::log(proj_momentum);
898 
899  if( proj_momentum <= 1.0 )
900  {
901  xsection = zz*(65.55 + 53.84/(proj_momentum+1.e-6) );
902  }
903  else
904  {
905  xsection = zz*( 41.1 + 77.2*std::pow( proj_momentum, -0.68)
906  + 0.293*logP*logP - 1.82*logP );
907  }
908  if ( nn > 0.)
909  {
910  xsection += nn*( 41.9 + 96.2*std::pow( proj_momentum, -0.99) - 0.154*logP);
911  }
912  else // H
913  {
914  fInelasticXsc = 38.0 + 38.0*std::pow( proj_momentum, -0.96)
915  - 0.169*logP*logP;
916  fInelasticXsc *= millibarn;
917  }
918  }
919  else if( theParticle == thePiPlus )
920  {
921  if(proj_momentum < 0.4)
922  {
923  G4double Ex3 = 180*std::exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.085/0.085);
924  hpXscv = Ex3+20.0;
925  }
926  else if( proj_momentum < 1.15 )
927  {
928  G4double Ex4 = 88*(std::log(proj_momentum/0.75))*(std::log(proj_momentum/0.75));
929  hpXscv = Ex4+14.0;
930  }
931  else if(proj_momentum < 3.5)
932  {
933  G4double Ex1 = 3.2*std::exp(-(proj_momentum-2.55)*(proj_momentum-2.55)/0.55/0.55);
934  G4double Ex2 = 12*std::exp(-(proj_momentum-1.47)*(proj_momentum-1.47)/0.225/0.225);
935  hpXscv = Ex1+Ex2+27.5;
936  }
937  else // if(proj_momentum > 3.5) // mb
938  {
939  hpXscv = 10.6+2.*std::log(proj_energy)+25*std::pow(proj_energy,-0.43);
940  }
941  // pi+n = pi-p??
942 
943  if(proj_momentum < 0.37)
944  {
945  hnXscv = 28.0 + 40*std::exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.07/0.07);
946  }
947  else if(proj_momentum<0.65)
948  {
949  hnXscv = 26+110*(std::log(proj_momentum/0.48))*(std::log(proj_momentum/0.48));
950  }
951  else if(proj_momentum<1.3)
952  {
953  hnXscv = 36.1+
954  10*std::exp(-(proj_momentum-0.72)*(proj_momentum-0.72)/0.06/0.06)+
955  24*std::exp(-(proj_momentum-1.015)*(proj_momentum-1.015)/0.075/0.075);
956  }
957  else if(proj_momentum<3.0)
958  {
959  hnXscv = 36.1+0.079-4.313*std::log(proj_momentum)+
960  3*std::exp(-(proj_momentum-2.1)*(proj_momentum-2.1)/0.4/0.4)+
961  1.5*std::exp(-(proj_momentum-1.4)*(proj_momentum-1.4)/0.12/0.12);
962  }
963  else // mb
964  {
965  hnXscv = 10.6+2*std::log(proj_energy)+30*std::pow(proj_energy,-0.43);
966  }
967  xsection = hpXscv*zz + hnXscv*nn;
968  }
969  else if(theParticle == thePiMinus)
970  {
971  // pi-n = pi+p??
972 
973  if(proj_momentum < 0.4)
974  {
975  G4double Ex3 = 180*std::exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.085/0.085);
976  hnXscv = Ex3+20.0;
977  }
978  else if(proj_momentum < 1.15)
979  {
980  G4double Ex4 = 88*(std::log(proj_momentum/0.75))*(std::log(proj_momentum/0.75));
981  hnXscv = Ex4+14.0;
982  }
983  else if(proj_momentum < 3.5)
984  {
985  G4double Ex1 = 3.2*std::exp(-(proj_momentum-2.55)*(proj_momentum-2.55)/0.55/0.55);
986  G4double Ex2 = 12*std::exp(-(proj_momentum-1.47)*(proj_momentum-1.47)/0.225/0.225);
987  hnXscv = Ex1+Ex2+27.5;
988  }
989  else // if(proj_momentum > 3.5) // mb
990  {
991  hnXscv = 10.6+2.*std::log(proj_energy)+25*std::pow(proj_energy,-0.43);
992  }
993  // pi-p
994 
995  if(proj_momentum < 0.37)
996  {
997  hpXscv = 28.0 + 40*std::exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.07/0.07);
998  }
999  else if(proj_momentum<0.65)
1000  {
1001  hpXscv = 26+110*(std::log(proj_momentum/0.48))*(std::log(proj_momentum/0.48));
1002  }
1003  else if(proj_momentum<1.3)
1004  {
1005  hpXscv = 36.1+
1006  10*std::exp(-(proj_momentum-0.72)*(proj_momentum-0.72)/0.06/0.06)+
1007  24*std::exp(-(proj_momentum-1.015)*(proj_momentum-1.015)/0.075/0.075);
1008  }
1009  else if(proj_momentum<3.0)
1010  {
1011  hpXscv = 36.1+0.079-4.313*std::log(proj_momentum)+
1012  3*std::exp(-(proj_momentum-2.1)*(proj_momentum-2.1)/0.4/0.4)+
1013  1.5*std::exp(-(proj_momentum-1.4)*(proj_momentum-1.4)/0.12/0.12);
1014  }
1015  else // mb
1016  {
1017  hpXscv = 10.6+2*std::log(proj_energy)+30*std::pow(proj_energy,-0.43);
1018  }
1019  xsection = hpXscv*zz + hnXscv*nn;
1020  }
1021  else if(theParticle == theKPlus)
1022  {
1023  xsection = zz*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
1024  + 7.14*std::pow(sMand,-eta1) - 13.45*std::pow(sMand,-eta2));
1025 
1026  xsection += nn*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
1027  + 5.17*std::pow(sMand,-eta1) - 7.23*std::pow(sMand,-eta2));
1028  }
1029  else if(theParticle == theKMinus)
1030  {
1031  xsection = zz*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
1032  + 7.14*std::pow(sMand,-eta1) + 13.45*std::pow(sMand,-eta2));
1033 
1034  xsection += nn*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
1035  + 5.17*std::pow(sMand,-eta1) + 7.23*std::pow(sMand,-eta2));
1036  }
1037  else if(theParticle == theSMinus)
1038  {
1039  xsection = aa*( 35.20 + B*std::pow(std::log(sMand/s0),2.)
1040  - 199.*std::pow(sMand,-eta1) + 264.*std::pow(sMand,-eta2));
1041  }
1042  else if(theParticle == theGamma) // modify later on
1043  {
1044  xsection = aa*( 0.0 + B*std::pow(std::log(sMand/s0),2.)
1045  + 0.032*std::pow(sMand,-eta1) - 0.0*std::pow(sMand,-eta2));
1046 
1047  }
1048  else // as proton ???
1049  {
1050  xsection = zz*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
1051  + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
1052 
1053  xsection += nn*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
1054  + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
1055  }
1056  xsection *= millibarn; // parametrised in mb
1057  return xsection;
1058 }
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4Element *)
G4double GetTotalEnergy() const
G4double CalcMandelstamS(const G4double, const G4double, const G4double)
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
int millibarn
Definition: hepunit.py:40
G4double GetMass() const
static G4ParticleTable * GetParticleTable()
double G4double
Definition: G4Types.hh:76
double mag() const
G4ThreeVector GetMomentum() const
G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscPDG ( const G4DynamicParticle aParticle,
const G4Element anElement 
)

Definition at line 571 of file G4ComponentGGHadronNucleusXsc.cc.

References G4lrint(), G4Element::GetN(), and G4Element::GetZ().

Referenced by GetHadronNucleonXscNS(), and GetKaonNucleonXscVector().

573 {
574  G4int At = G4lrint(anElement->GetN()); // number of nucleons
575  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
576 
577  return GetHadronNucleonXscPDG(aParticle, At, Zt);
578 }
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4Element *)
G4double GetN() const
Definition: G4Element.hh:134
G4double GetZ() const
Definition: G4Element.hh:131
int G4int
Definition: G4Types.hh:78
int G4lrint(double ad)
Definition: templates.hh:163
G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscPDG ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 590 of file G4ComponentGGHadronNucleusXsc.cc.

References CalcMandelstamS(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetMass(), G4DynamicParticle::GetMomentum(), G4ParticleTable::GetParticleTable(), python.hepunit::GeV, CLHEP::Hep3Vector::mag(), python.hepunit::millibarn, G4InuclParticleNames::nn, and G4InuclParticleNames::s0.

592 {
593  G4double xsection;
594 
595  G4int Nt = At-Zt; // number of neutrons
596  if (Nt < 0) Nt = 0;
597 
598  G4double zz = Zt;
599  G4double aa = At;
600  G4double nn = Nt;
601 
603  GetIonTable()->GetIonMass(Zt, At);
604 
605  targ_mass = 0.939*GeV; // ~mean neutron and proton ???
606 
607  G4double proj_mass = aParticle->GetMass();
608  G4double proj_momentum = aParticle->GetMomentum().mag();
609 
610  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
611 
612  sMand /= GeV*GeV; // in GeV for parametrisation
613 
614  // General PDG fit constants
615 
616  G4double s0 = 5.38*5.38; // in Gev^2
617  G4double eta1 = 0.458;
618  G4double eta2 = 0.458;
619  G4double B = 0.308;
620 
621 
622  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
623 
624 
625  if(theParticle == theNeutron) // proton-neutron fit
626  {
627  xsection = zz*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
628  + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
629  xsection += nn*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
630  + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2)); // pp for nn
631  }
632  else if(theParticle == theProton)
633  {
634 
635  xsection = zz*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
636  + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
637 
638  xsection += nn*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
639  + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
640  }
641  else if(theParticle == theAProton)
642  {
643  xsection = zz*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
644  + 42.53*std::pow(sMand,-eta1) + 33.34*std::pow(sMand,-eta2));
645 
646  xsection += nn*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
647  + 40.15*std::pow(sMand,-eta1) + 30.*std::pow(sMand,-eta2));
648  }
649  else if(theParticle == thePiPlus)
650  {
651  xsection = aa*( 20.86 + B*std::pow(std::log(sMand/s0),2.)
652  + 19.24*std::pow(sMand,-eta1) - 6.03*std::pow(sMand,-eta2));
653  }
654  else if(theParticle == thePiMinus)
655  {
656  xsection = aa*( 20.86 + B*std::pow(std::log(sMand/s0),2.)
657  + 19.24*std::pow(sMand,-eta1) + 6.03*std::pow(sMand,-eta2));
658  }
659  else if(theParticle == theKPlus || theParticle == theK0L )
660  {
661  xsection = zz*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
662  + 7.14*std::pow(sMand,-eta1) - 13.45*std::pow(sMand,-eta2));
663 
664  xsection += nn*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
665  + 5.17*std::pow(sMand,-eta1) - 7.23*std::pow(sMand,-eta2));
666  }
667  else if(theParticle == theKMinus || theParticle == theK0S )
668  {
669  xsection = zz*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
670  + 7.14*std::pow(sMand,-eta1) + 13.45*std::pow(sMand,-eta2));
671 
672  xsection += nn*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
673  + 5.17*std::pow(sMand,-eta1) + 7.23*std::pow(sMand,-eta2));
674  }
675  else if(theParticle == theSMinus)
676  {
677  xsection = aa*( 35.20 + B*std::pow(std::log(sMand/s0),2.)
678  - 199.*std::pow(sMand,-eta1) + 264.*std::pow(sMand,-eta2));
679  }
680  else if(theParticle == theGamma) // modify later on
681  {
682  xsection = aa*( 0.0 + B*std::pow(std::log(sMand/s0),2.)
683  + 0.032*std::pow(sMand,-eta1) - 0.0*std::pow(sMand,-eta2));
684 
685  }
686  else // as proton ???
687  {
688  xsection = zz*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
689  + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
690 
691  xsection += nn*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
692  + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
693  }
694  xsection *= millibarn; // parametrised in mb
695  return xsection;
696 }
G4double CalcMandelstamS(const G4double, const G4double, const G4double)
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
int millibarn
Definition: hepunit.py:40
G4double GetMass() const
static G4ParticleTable * GetParticleTable()
double G4double
Definition: G4Types.hh:76
double mag() const
G4ThreeVector GetMomentum() const
G4double G4ComponentGGHadronNucleusXsc::GetHNinelasticXsc ( const G4DynamicParticle aParticle,
const G4Element anElement 
)

Definition at line 1099 of file G4ComponentGGHadronNucleusXsc.cc.

References G4lrint(), G4Element::GetN(), and G4Element::GetZ().

Referenced by GetIsoCrossSection(), and GetRatioQE().

1101 {
1102  G4int At = G4lrint(anElement->GetN()); // number of nucleons
1103  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
1104 
1105  return GetHNinelasticXsc(aParticle, At, Zt);
1106 }
G4double GetN() const
Definition: G4Element.hh:134
G4double GetZ() const
Definition: G4Element.hh:131
int G4int
Definition: G4Types.hh:78
G4double GetHNinelasticXsc(const G4DynamicParticle *, const G4Element *)
int G4lrint(double ad)
Definition: templates.hh:163
G4double G4ComponentGGHadronNucleusXsc::GetHNinelasticXsc ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 1113 of file G4ComponentGGHadronNucleusXsc.cc.

References G4DynamicParticle::GetDefinition(), GetHadronNucleonXscNS(), and GetHNinelasticXscVU().

1115 {
1116  G4ParticleDefinition* hadron = aParticle->GetDefinition();
1117  G4double sumInelastic;
1118  G4int Nt = At - Zt;
1119  if(Nt < 0) Nt = 0;
1120 
1121  if( hadron == theKPlus )
1122  {
1123  sumInelastic = GetHNinelasticXscVU(aParticle, At, Zt);
1124  }
1125  else
1126  {
1127  //sumInelastic = Zt*GetHadronNucleonXscMK(aParticle, theProton);
1128  // sumInelastic += Nt*GetHadronNucleonXscMK(aParticle, theNeutron);
1129  sumInelastic = G4double(Zt)*GetHadronNucleonXscNS(aParticle, 1, 1);
1130  sumInelastic += G4double(Nt)*GetHadronNucleonXscNS(aParticle, 1, 0);
1131  }
1132  return sumInelastic;
1133 }
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
G4double GetHNinelasticXscVU(const G4DynamicParticle *, G4int At, G4int Zt)
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
double G4double
Definition: G4Types.hh:76
G4double G4ComponentGGHadronNucleusXsc::GetHNinelasticXscVU ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 1141 of file G4ComponentGGHadronNucleusXsc.cc.

References G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetMomentum(), G4ParticleDefinition::GetPDGEncoding(), G4DynamicParticle::GetTotalEnergy(), python.hepunit::GeV, CLHEP::Hep3Vector::mag(), and python.hepunit::millibarn.

Referenced by GetHNinelasticXsc().

1143 {
1144  G4int PDGcode = aParticle->GetDefinition()->GetPDGEncoding();
1145  G4int absPDGcode = std::abs(PDGcode);
1146 
1147  G4double Elab = aParticle->GetTotalEnergy();
1148  // (s - 2*0.88*GeV*GeV)/(2*0.939*GeV)/GeV;
1149  G4double Plab = aParticle->GetMomentum().mag();
1150  // std::sqrt(Elab * Elab - 0.88);
1151 
1152  Elab /= GeV;
1153  Plab /= GeV;
1154 
1155  G4double LogPlab = std::log( Plab );
1156  G4double sqrLogPlab = LogPlab * LogPlab;
1157 
1158  //G4cout<<"Plab = "<<Plab<<G4endl;
1159 
1160  G4double NumberOfTargetProtons = G4double(Zt);
1161  G4double NumberOfTargetNucleons = G4double(At);
1162  G4double NumberOfTargetNeutrons = NumberOfTargetNucleons - NumberOfTargetProtons;
1163 
1164  if(NumberOfTargetNeutrons < 0.0) NumberOfTargetNeutrons = 0.0;
1165 
1166  G4double Xtotal, Xelastic, Xinelastic;
1167 
1168  if( absPDGcode > 1000 ) //------Projectile is baryon --------
1169  {
1170  G4double XtotPP = 48.0 + 0. *std::pow(Plab, 0. ) +
1171  0.522*sqrLogPlab - 4.51*LogPlab;
1172 
1173  G4double XtotPN = 47.3 + 0. *std::pow(Plab, 0. ) +
1174  0.513*sqrLogPlab - 4.27*LogPlab;
1175 
1176  G4double XelPP = 11.9 + 26.9*std::pow(Plab,-1.21) +
1177  0.169*sqrLogPlab - 1.85*LogPlab;
1178 
1179  G4double XelPN = 11.9 + 26.9*std::pow(Plab,-1.21) +
1180  0.169*sqrLogPlab - 1.85*LogPlab;
1181 
1182  Xtotal = (NumberOfTargetProtons * XtotPP +
1183  NumberOfTargetNeutrons * XtotPN);
1184 
1185  Xelastic = (NumberOfTargetProtons * XelPP +
1186  NumberOfTargetNeutrons * XelPN);
1187  }
1188  else if( PDGcode == 211 ) //------Projectile is PionPlus -------
1189  {
1190  G4double XtotPiP = 16.4 + 19.3 *std::pow(Plab,-0.42) +
1191  0.19 *sqrLogPlab - 0.0 *LogPlab;
1192 
1193  G4double XtotPiN = 33.0 + 14.0 *std::pow(Plab,-1.36) +
1194  0.456*sqrLogPlab - 4.03*LogPlab;
1195 
1196  G4double XelPiP = 0.0 + 11.4*std::pow(Plab,-0.40) +
1197  0.079*sqrLogPlab - 0.0 *LogPlab;
1198 
1199  G4double XelPiN = 1.76 + 11.2*std::pow(Plab,-0.64) +
1200  0.043*sqrLogPlab - 0.0 *LogPlab;
1201 
1202  Xtotal = ( NumberOfTargetProtons * XtotPiP +
1203  NumberOfTargetNeutrons * XtotPiN );
1204 
1205  Xelastic = ( NumberOfTargetProtons * XelPiP +
1206  NumberOfTargetNeutrons * XelPiN );
1207  }
1208  else if( PDGcode == -211 ) //------Projectile is PionMinus -------
1209  {
1210  G4double XtotPiP = 33.0 + 14.0 *std::pow(Plab,-1.36) +
1211  0.456*sqrLogPlab - 4.03*LogPlab;
1212 
1213  G4double XtotPiN = 16.4 + 19.3 *std::pow(Plab,-0.42) +
1214  0.19 *sqrLogPlab - 0.0 *LogPlab;
1215 
1216  G4double XelPiP = 1.76 + 11.2*std::pow(Plab,-0.64) +
1217  0.043*sqrLogPlab - 0.0 *LogPlab;
1218 
1219  G4double XelPiN = 0.0 + 11.4*std::pow(Plab,-0.40) +
1220  0.079*sqrLogPlab - 0.0 *LogPlab;
1221 
1222  Xtotal = ( NumberOfTargetProtons * XtotPiP +
1223  NumberOfTargetNeutrons * XtotPiN );
1224 
1225  Xelastic = ( NumberOfTargetProtons * XelPiP +
1226  NumberOfTargetNeutrons * XelPiN );
1227  }
1228  else if( PDGcode == 111 ) //------Projectile is PionZero -------
1229  {
1230  G4double XtotPiP =(16.4 + 19.3 *std::pow(Plab,-0.42) +
1231  0.19 *sqrLogPlab - 0.0 *LogPlab + //Pi+
1232  33.0 + 14.0 *std::pow(Plab,-1.36) +
1233  0.456*sqrLogPlab - 4.03*LogPlab)/2; //Pi-
1234 
1235  G4double XtotPiN =(33.0 + 14.0 *std::pow(Plab,-1.36) +
1236  0.456*sqrLogPlab - 4.03*LogPlab + //Pi+
1237  16.4 + 19.3 *std::pow(Plab,-0.42) +
1238  0.19 *sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1239 
1240  G4double XelPiP =( 0.0 + 11.4*std::pow(Plab,-0.40) +
1241  0.079*sqrLogPlab - 0.0 *LogPlab + //Pi+
1242  1.76 + 11.2*std::pow(Plab,-0.64) +
1243  0.043*sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1244 
1245  G4double XelPiN =( 1.76 + 11.2*std::pow(Plab,-0.64) +
1246  0.043*sqrLogPlab - 0.0 *LogPlab + //Pi+
1247  0.0 + 11.4*std::pow(Plab,-0.40) +
1248  0.079*sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1249 
1250  Xtotal = ( NumberOfTargetProtons * XtotPiP +
1251  NumberOfTargetNeutrons * XtotPiN );
1252 
1253  Xelastic = ( NumberOfTargetProtons * XelPiP +
1254  NumberOfTargetNeutrons * XelPiN );
1255  }
1256  else if( PDGcode == 321 ) //------Projectile is KaonPlus -------
1257  {
1258  G4double XtotKP = 18.1 + 0. *std::pow(Plab, 0. ) +
1259  0.26 *sqrLogPlab - 1.0 *LogPlab;
1260  G4double XtotKN = 18.7 + 0. *std::pow(Plab, 0. ) +
1261  0.21 *sqrLogPlab - 0.89*LogPlab;
1262 
1263  G4double XelKP = 5.0 + 8.1*std::pow(Plab,-1.8 ) +
1264  0.16 *sqrLogPlab - 1.3 *LogPlab;
1265 
1266  G4double XelKN = 7.3 + 0. *std::pow(Plab,-0. ) +
1267  0.29 *sqrLogPlab - 2.4 *LogPlab;
1268 
1269  Xtotal = ( NumberOfTargetProtons * XtotKP +
1270  NumberOfTargetNeutrons * XtotKN );
1271 
1272  Xelastic = ( NumberOfTargetProtons * XelKP +
1273  NumberOfTargetNeutrons * XelKN );
1274  }
1275  else if( PDGcode ==-321 ) //------Projectile is KaonMinus ------
1276  {
1277  G4double XtotKP = 32.1 + 0. *std::pow(Plab, 0. ) +
1278  0.66 *sqrLogPlab - 5.6 *LogPlab;
1279  G4double XtotKN = 25.2 + 0. *std::pow(Plab, 0. ) +
1280  0.38 *sqrLogPlab - 2.9 *LogPlab;
1281 
1282  G4double XelKP = 7.3 + 0. *std::pow(Plab,-0. ) +
1283  0.29 *sqrLogPlab - 2.4 *LogPlab;
1284 
1285  G4double XelKN = 5.0 + 8.1*std::pow(Plab,-1.8 ) +
1286  0.16 *sqrLogPlab - 1.3 *LogPlab;
1287 
1288  Xtotal = ( NumberOfTargetProtons * XtotKP +
1289  NumberOfTargetNeutrons * XtotKN );
1290 
1291  Xelastic = ( NumberOfTargetProtons * XelKP +
1292  NumberOfTargetNeutrons * XelKN );
1293  }
1294  else if( PDGcode == 311 ) //------Projectile is KaonZero ------
1295  {
1296  G4double XtotKP = ( 18.1 + 0. *std::pow(Plab, 0. ) +
1297  0.26 *sqrLogPlab - 1.0 *LogPlab + //K+
1298  32.1 + 0. *std::pow(Plab, 0. ) +
1299  0.66 *sqrLogPlab - 5.6 *LogPlab)/2; //K-
1300 
1301  G4double XtotKN = ( 18.7 + 0. *std::pow(Plab, 0. ) +
1302  0.21 *sqrLogPlab - 0.89*LogPlab + //K+
1303  25.2 + 0. *std::pow(Plab, 0. ) +
1304  0.38 *sqrLogPlab - 2.9 *LogPlab)/2; //K-
1305 
1306  G4double XelKP = ( 5.0 + 8.1*std::pow(Plab,-1.8 )
1307  + 0.16 *sqrLogPlab - 1.3 *LogPlab + //K+
1308  7.3 + 0. *std::pow(Plab,-0. ) +
1309  0.29 *sqrLogPlab - 2.4 *LogPlab)/2; //K-
1310 
1311  G4double XelKN = ( 7.3 + 0. *std::pow(Plab,-0. ) +
1312  0.29 *sqrLogPlab - 2.4 *LogPlab + //K+
1313  5.0 + 8.1*std::pow(Plab,-1.8 ) +
1314  0.16 *sqrLogPlab - 1.3 *LogPlab)/2; //K-
1315 
1316  Xtotal = ( NumberOfTargetProtons * XtotKP +
1317  NumberOfTargetNeutrons * XtotKN );
1318 
1319  Xelastic = ( NumberOfTargetProtons * XelKP +
1320  NumberOfTargetNeutrons * XelKN );
1321  }
1322  else //------Projectile is undefined, Nucleon assumed
1323  {
1324  G4double XtotPP = 48.0 + 0. *std::pow(Plab, 0. ) +
1325  0.522*sqrLogPlab - 4.51*LogPlab;
1326 
1327  G4double XtotPN = 47.3 + 0. *std::pow(Plab, 0. ) +
1328  0.513*sqrLogPlab - 4.27*LogPlab;
1329 
1330  G4double XelPP = 11.9 + 26.9*std::pow(Plab,-1.21) +
1331  0.169*sqrLogPlab - 1.85*LogPlab;
1332  G4double XelPN = 11.9 + 26.9*std::pow(Plab,-1.21) +
1333  0.169*sqrLogPlab - 1.85*LogPlab;
1334 
1335  Xtotal = ( NumberOfTargetProtons * XtotPP +
1336  NumberOfTargetNeutrons * XtotPN );
1337 
1338  Xelastic = ( NumberOfTargetProtons * XelPP +
1339  NumberOfTargetNeutrons * XelPN );
1340  }
1341  Xinelastic = Xtotal - Xelastic;
1342 
1343  if( Xinelastic < 0.) Xinelastic = 0.;
1344 
1345  return Xinelastic*= millibarn;
1346 }
G4double GetTotalEnergy() const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
int millibarn
Definition: hepunit.py:40
double G4double
Definition: G4Types.hh:76
double mag() const
G4ThreeVector GetMomentum() const
G4double G4ComponentGGHadronNucleusXsc::GetInelasticElementCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4double  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 154 of file G4ComponentGGHadronNucleusXsc.cc.

References GetIsoCrossSection().

157 {
158  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
159  kinEnergy);
160  fTotalXsc = GetIsoCrossSection(aDP, Z, G4int(A));
161  delete aDP;
162 
163  return fInelasticXsc;
164 }
int G4int
Definition: G4Types.hh:78
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4ThreeVector G4ParticleMomentum
G4double G4ComponentGGHadronNucleusXsc::GetInelasticGlauberGribov ( const G4DynamicParticle dp,
G4int  Z,
G4int  A 
)
inline

Definition at line 231 of file G4ComponentGGHadronNucleusXsc.hh.

References GetIsoCrossSection().

233 {
234  GetIsoCrossSection(dp, Z, A);
235  return fInelasticXsc;
236 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4double G4ComponentGGHadronNucleusXsc::GetInelasticGlauberGribovXsc ( )
inline

Definition at line 148 of file G4ComponentGGHadronNucleusXsc.hh.

148 { return fInelasticXsc; };
G4double G4ComponentGGHadronNucleusXsc::GetInelasticIsotopeCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4int  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 126 of file G4ComponentGGHadronNucleusXsc.cc.

References GetIsoCrossSection().

129 {
130  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
131  kinEnergy);
132  fTotalXsc = GetIsoCrossSection(aDP, Z, A);
133  delete aDP;
134 
135  return fInelasticXsc;
136 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4ThreeVector G4ParticleMomentum
G4double G4ComponentGGHadronNucleusXsc::GetIsoCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
G4int  A,
const G4Isotope iso = 0,
const G4Element elm = 0,
const G4Material mat = 0 
)

Definition at line 270 of file G4ComponentGGHadronNucleusXsc.cc.

References python.hepunit::fermi, G4DynamicParticle::GetDefinition(), G4HadronNucleonXsc::GetHadronNucleonXscNS(), GetHadronNucleonXscNS(), GetHNinelasticXsc(), G4HadronNucleonXsc::GetInelasticHadronNucleonXsc(), G4HadronNucleonXsc::GetKaonNucleonXscGG(), GetNucleusRadius(), GetParticleBarCorIn(), GetParticleBarCorTot(), N, and python.hepunit::pi.

Referenced by ComputeQuasiElasticRatio(), GetElasticElementCrossSection(), GetElasticGlauberGribov(), GetElasticIsotopeCrossSection(), GetInelasticElementCrossSection(), GetInelasticGlauberGribov(), GetInelasticIsotopeCrossSection(), GetProductionElementCrossSection(), GetProductionIsotopeCrossSection(), GetTotalElementCrossSection(), and GetTotalIsotopeCrossSection().

275 {
276  G4double xsection, sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
277  G4double hpInXsc(0.), hnInXsc(0.);
278  G4double R = GetNucleusRadius(A);
279 
280  G4int N = A - Z; // number of neutrons
281  if (N < 0) N = 0;
282 
283  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
284 
285  if( theParticle == theProton ||
286  theParticle == theNeutron ||
287  theParticle == thePiPlus ||
288  theParticle == thePiMinus )
289  {
290  // sigma = GetHadronNucleonXscNS(aParticle, A, Z);
291 
292  sigma = Z*hnXsc->GetHadronNucleonXscNS(aParticle, theProton);
293 
294  hpInXsc = hnXsc->GetInelasticHadronNucleonXsc();
295 
296  sigma += N*hnXsc->GetHadronNucleonXscNS(aParticle, theNeutron);
297 
298  hnInXsc = hnXsc->GetInelasticHadronNucleonXsc();
299 
300  cofInelastic = 2.4;
301  cofTotal = 2.0;
302  }
303  else if( theParticle == theKPlus ||
304  theParticle == theKMinus ||
305  theParticle == theK0S ||
306  theParticle == theK0L )
307  {
308  // sigma = GetKaonNucleonXscVector(aParticle, A, Z);
309 
310  sigma = Z*hnXsc->GetKaonNucleonXscGG(aParticle, theProton);
311 
312  hpInXsc = hnXsc->GetInelasticHadronNucleonXsc();
313 
314  sigma += N*hnXsc->GetKaonNucleonXscGG(aParticle, theNeutron);
315 
316  hnInXsc = hnXsc->GetInelasticHadronNucleonXsc();
317 
318  cofInelastic = 2.2;
319  cofTotal = 2.0;
320  R = 1.3*fermi;
321  R *= std::pow(G4double(A), 0.3333);
322  }
323  else
324  {
325  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
326  cofInelastic = 2.2;
327  cofTotal = 2.0;
328  }
329  // cofInelastic = 2.0;
330 
331  if( A > 1 )
332  {
333  nucleusSquare = cofTotal*pi*R*R; // basically 2piRR
334  ratio = sigma/nucleusSquare;
335 
336  xsection = nucleusSquare*std::log( 1. + ratio );
337 
338  xsection *= GetParticleBarCorTot(theParticle, Z);
339 
340  fTotalXsc = xsection;
341 
342 
343 
344  fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
345 
346  fInelasticXsc *= GetParticleBarCorIn(theParticle, Z);
347 
348  fElasticXsc = fTotalXsc - fInelasticXsc;
349 
350  if(fElasticXsc < 0.) fElasticXsc = 0.;
351 
352  G4double difratio = ratio/(1.+ratio);
353 
354  fDiffractionXsc = 0.5*nucleusSquare*( difratio - std::log( 1. + difratio ) );
355 
356 
357  // sigma = GetHNinelasticXsc(aParticle, A, Z);
358 
359  sigma = Z*hpInXsc + N*hnInXsc;
360 
361  ratio = sigma/nucleusSquare;
362 
363  fProductionXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
364 
365  fProductionXsc *= GetParticleBarCorIn(theParticle, Z);
366 
367  if (fElasticXsc < 0.) fElasticXsc = 0.;
368  }
369  else // H
370  {
371  fTotalXsc = sigma;
372  xsection = sigma;
373 
374  if ( theParticle != theAProton )
375  {
376  sigma = GetHNinelasticXsc(aParticle, A, Z);
377  fInelasticXsc = sigma;
378  fElasticXsc = fTotalXsc - fInelasticXsc;
379  }
380  else
381  {
382  fElasticXsc = fTotalXsc - fInelasticXsc;
383  }
384  if (fElasticXsc < 0.) fElasticXsc = 0.;
385 
386  }
387  return xsection;
388 }
G4double GetNucleusRadius(const G4DynamicParticle *, const G4Element *)
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
G4double GetKaonNucleonXscGG(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double GetHNinelasticXsc(const G4DynamicParticle *, const G4Element *)
G4double GetParticleBarCorIn(const G4ParticleDefinition *theParticle, G4int Z)
**D E S C R I P T I O N
double G4double
Definition: G4Types.hh:76
G4double GetParticleBarCorTot(const G4ParticleDefinition *theParticle, G4int Z)
G4double GetInelasticHadronNucleonXsc()
G4double G4ComponentGGHadronNucleusXsc::GetKaonNucleonXscVector ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 1061 of file G4ComponentGGHadronNucleusXsc.cc.

References G4DynamicParticle::GetDefinition(), GetHadronNucleonXscPDG(), G4DynamicParticle::GetKineticEnergy(), G4HadronNucleonXsc::GetKmNeutronTotXscVector(), G4HadronNucleonXsc::GetKmProtonTotXscVector(), G4HadronNucleonXsc::GetKpNeutronTotXscVector(), G4HadronNucleonXsc::GetKpProtonTotXscVector(), and python.hepunit::GeV.

1063 {
1064  G4double Tkin, logTkin, xsc, xscP, xscN;
1065  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
1066 
1067  G4int Nt = At-Zt; // number of neutrons
1068  if (Nt < 0) Nt = 0;
1069 
1070  Tkin = aParticle->GetKineticEnergy(); // Tkin in MeV
1071 
1072  if( Tkin > 70*GeV ) return GetHadronNucleonXscPDG(aParticle,At,Zt);
1073 
1074  logTkin = std::log(Tkin); // Tkin in MeV!!!
1075 
1076  if( theParticle == theKPlus )
1077  {
1078  xscP = hnXsc->GetKpProtonTotXscVector(logTkin);
1079  xscN = hnXsc->GetKpNeutronTotXscVector(logTkin);
1080  }
1081  else if( theParticle == theKMinus )
1082  {
1083  xscP = hnXsc->GetKmProtonTotXscVector(logTkin);
1084  xscN = hnXsc->GetKmNeutronTotXscVector(logTkin);
1085  }
1086  else // K-zero as half of K+ and K-
1087  {
1088  xscP = (hnXsc->GetKpProtonTotXscVector(logTkin)+hnXsc->GetKmProtonTotXscVector(logTkin))*0.5;
1089  xscN = (hnXsc->GetKpNeutronTotXscVector(logTkin)+hnXsc->GetKmNeutronTotXscVector(logTkin))*0.5;
1090  }
1091  xsc = xscP*Zt + xscN*Nt;
1092  return xsc;
1093 }
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4Element *)
G4double GetKmNeutronTotXscVector(G4double logEnergy)
G4double GetKpProtonTotXscVector(G4double logEnergy)
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
G4double GetKmProtonTotXscVector(G4double logEnergy)
double G4double
Definition: G4Types.hh:76
G4double GetKpNeutronTotXscVector(G4double logEnergy)
G4double G4ComponentGGHadronNucleusXsc::GetNucleusRadius ( const G4DynamicParticle ,
const G4Element anElement 
)

Definition at line 1353 of file G4ComponentGGHadronNucleusXsc.cc.

References G4lrint(), G4Element::GetN(), and G4INCL::Math::oneThird.

Referenced by GetIsoCrossSection(), GetRatioQE(), and GetRatioSD().

1355 {
1356  G4int At = G4lrint(anElement->GetN());
1357  G4double oneThird = 1.0/3.0;
1358  G4double cubicrAt = std::pow(G4double(At), oneThird);
1359 
1360  G4double R; // = fRadiusConst*cubicrAt;
1361  /*
1362  G4double tmp = std::pow( cubicrAt-1., 3.);
1363  tmp += At;
1364  tmp *= 0.5;
1365 
1366  if (At > 20.) // 20.
1367  {
1368  R = fRadiusConst*std::pow (tmp, oneThird);
1369  }
1370  else
1371  {
1372  R = fRadiusConst*cubicrAt;
1373  }
1374  */
1375 
1376  R = fRadiusConst*cubicrAt;
1377 
1378  G4double meanA = 21.;
1379 
1380  G4double tauA1 = 40.;
1381  G4double tauA2 = 10.;
1382  G4double tauA3 = 5.;
1383 
1384  G4double a1 = 0.85;
1385  G4double b1 = 1. - a1;
1386 
1387  G4double b2 = 0.3;
1388  G4double b3 = 4.;
1389 
1390  if (At > 20) // 20.
1391  {
1392  R *= ( a1 + b1*std::exp( -(At - meanA)/tauA1) );
1393  }
1394  else if (At > 3)
1395  {
1396  R *= ( 1.0 + b2*( 1. - std::exp( (At - meanA)/tauA2) ) );
1397  }
1398  else
1399  {
1400  R *= ( 1.0 + b3*( 1. - std::exp( (At - meanA)/tauA3) ) );
1401  }
1402  return R;
1403 
1404 }
G4double GetN() const
Definition: G4Element.hh:134
int G4int
Definition: G4Types.hh:78
int G4lrint(double ad)
Definition: templates.hh:163
double G4double
Definition: G4Types.hh:76
const G4double oneThird
G4double G4ComponentGGHadronNucleusXsc::GetNucleusRadius ( G4int  At)

Definition at line 1410 of file G4ComponentGGHadronNucleusXsc.cc.

References G4INCL::Math::oneThird.

1411 {
1412  G4double oneThird = 1.0/3.0;
1413  G4double cubicrAt = std::pow(G4double(At), oneThird);
1414 
1415  G4double R; // = fRadiusConst*cubicrAt;
1416 
1417  /*
1418  G4double tmp = std::pow( cubicrAt-1., 3.);
1419  tmp += At;
1420  tmp *= 0.5;
1421 
1422  if (At > 20.)
1423  {
1424  R = fRadiusConst*std::pow (tmp, oneThird);
1425  }
1426  else
1427  {
1428  R = fRadiusConst*cubicrAt;
1429  }
1430  */
1431 
1432  R = fRadiusConst*cubicrAt;
1433 
1434  G4double meanA = 20.;
1435  G4double tauA = 20.;
1436 
1437  if (At > 20) // 20.
1438  {
1439  R *= ( 0.8 + 0.2*std::exp( -(G4double(At) - meanA)/tauA) );
1440  }
1441  else
1442  {
1443  R *= ( 1.0 + 0.1*( 1. - std::exp( (G4double(At) - meanA)/tauA) ) );
1444  }
1445 
1446  return R;
1447 }
double G4double
Definition: G4Types.hh:76
const G4double oneThird
G4double G4ComponentGGHadronNucleusXsc::GetParticleBarCorIn ( const G4ParticleDefinition theParticle,
G4int  Z 
)
inline

Definition at line 264 of file G4ComponentGGHadronNucleusXsc.hh.

Referenced by GetIsoCrossSection().

266 {
267  if(Z >= 2 && Z <= 92)
268  {
269  if( theParticle == theProton ) return fProtonBarCorrectionIn[Z];
270  else if( theParticle == theNeutron) return fNeutronBarCorrectionIn[Z];
271  else if( theParticle == thePiPlus ) return fPionPlusBarCorrectionIn[Z];
272  else if( theParticle == thePiMinus) return fPionMinusBarCorrectionIn[Z];
273  else return 1.0;
274  }
275  else return 1.0;
276 }
G4double G4ComponentGGHadronNucleusXsc::GetParticleBarCorTot ( const G4ParticleDefinition theParticle,
G4int  Z 
)
inline

Definition at line 244 of file G4ComponentGGHadronNucleusXsc.hh.

Referenced by GetIsoCrossSection().

246 {
247  if(Z >= 2 && Z <= 92)
248  {
249  if( theParticle == theProton ) return fProtonBarCorrectionTot[Z];
250  else if( theParticle == theNeutron) return fNeutronBarCorrectionTot[Z];
251  else if( theParticle == thePiPlus ) return fPionPlusBarCorrectionTot[Z];
252  else if( theParticle == thePiMinus) return fPionMinusBarCorrectionTot[Z];
253  else return 1.0;
254  }
255  else return 1.0;
256 }
G4double G4ComponentGGHadronNucleusXsc::GetProductionElementCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4double  A 
)
virtual

Definition at line 168 of file G4ComponentGGHadronNucleusXsc.cc.

References GetIsoCrossSection().

171 {
172  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
173  kinEnergy);
174  fTotalXsc = GetIsoCrossSection(aDP, Z, G4int(A));
175  delete aDP;
176 
177  return fProductionXsc;
178 }
int G4int
Definition: G4Types.hh:78
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4ThreeVector G4ParticleMomentum
G4double G4ComponentGGHadronNucleusXsc::GetProductionGlauberGribovXsc ( )
inline

Definition at line 149 of file G4ComponentGGHadronNucleusXsc.hh.

149 { return fProductionXsc; };
G4double G4ComponentGGHadronNucleusXsc::GetProductionIsotopeCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4int  A 
)
virtual

Definition at line 140 of file G4ComponentGGHadronNucleusXsc.cc.

References GetIsoCrossSection().

143 {
144  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
145  kinEnergy);
146  fTotalXsc = GetIsoCrossSection(aDP, Z, A);
147  delete aDP;
148 
149  return fProductionXsc;
150 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4ThreeVector G4ParticleMomentum
G4double G4ComponentGGHadronNucleusXsc::GetRadiusConst ( )
inline

Definition at line 151 of file G4ComponentGGHadronNucleusXsc.hh.

151 { return fRadiusConst; };
G4double G4ComponentGGHadronNucleusXsc::GetRatioQE ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 437 of file G4ComponentGGHadronNucleusXsc.cc.

References G4DynamicParticle::GetDefinition(), GetHadronNucleonXscNS(), GetHNinelasticXsc(), GetNucleusRadius(), and python.hepunit::pi.

438 {
439  G4double sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
440  G4double R = GetNucleusRadius(A);
441 
442  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
443 
444  if( theParticle == theProton ||
445  theParticle == theNeutron ||
446  theParticle == thePiPlus ||
447  theParticle == thePiMinus )
448  {
449  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
450  cofInelastic = 2.4;
451  cofTotal = 2.0;
452  }
453  else
454  {
455  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
456  cofInelastic = 2.2;
457  cofTotal = 2.0;
458  }
459  nucleusSquare = cofTotal*pi*R*R; // basically 2piRR
460  ratio = sigma/nucleusSquare;
461 
462  fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
463 
464  sigma = GetHNinelasticXsc(aParticle, A, Z);
465  ratio = sigma/nucleusSquare;
466 
467  fProductionXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
468 
469  if (fInelasticXsc > fProductionXsc) ratio = (fInelasticXsc-fProductionXsc)/fInelasticXsc;
470  else ratio = 0.;
471  if ( ratio < 0. ) ratio = 0.;
472 
473  return ratio;
474 }
G4double GetNucleusRadius(const G4DynamicParticle *, const G4Element *)
G4ParticleDefinition * GetDefinition() const
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
G4double GetHNinelasticXsc(const G4DynamicParticle *, const G4Element *)
double G4double
Definition: G4Types.hh:76
G4double G4ComponentGGHadronNucleusXsc::GetRatioSD ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 395 of file G4ComponentGGHadronNucleusXsc.cc.

References G4DynamicParticle::GetDefinition(), GetHadronNucleonXscNS(), GetNucleusRadius(), and python.hepunit::pi.

396 {
397  G4double sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
398  G4double R = GetNucleusRadius(A);
399 
400  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
401 
402  if( theParticle == theProton ||
403  theParticle == theNeutron ||
404  theParticle == thePiPlus ||
405  theParticle == thePiMinus )
406  {
407  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
408  cofInelastic = 2.4;
409  cofTotal = 2.0;
410  }
411  else
412  {
413  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
414  cofInelastic = 2.2;
415  cofTotal = 2.0;
416  }
417  nucleusSquare = cofTotal*pi*R*R; // basically 2piRR
418  ratio = sigma/nucleusSquare;
419 
420  fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
421 
422  G4double difratio = ratio/(1.+ratio);
423 
424  fDiffractionXsc = 0.5*nucleusSquare*( difratio - std::log( 1. + difratio ) );
425 
426  if (fInelasticXsc > 0.) ratio = fDiffractionXsc/fInelasticXsc;
427  else ratio = 0.;
428 
429  return ratio;
430 }
G4double GetNucleusRadius(const G4DynamicParticle *, const G4Element *)
G4ParticleDefinition * GetDefinition() const
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
double G4double
Definition: G4Types.hh:76
G4double G4ComponentGGHadronNucleusXsc::GetTotalElementCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4double  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 112 of file G4ComponentGGHadronNucleusXsc.cc.

References GetIsoCrossSection().

115 {
116  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
117  kinEnergy);
118  fTotalXsc = GetIsoCrossSection(aDP, Z, G4int(A));
119  delete aDP;
120 
121  return fTotalXsc;
122 }
int G4int
Definition: G4Types.hh:78
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4ThreeVector G4ParticleMomentum
G4double G4ComponentGGHadronNucleusXsc::GetTotalGlauberGribovXsc ( )
inline

Definition at line 146 of file G4ComponentGGHadronNucleusXsc.hh.

146 { return fTotalXsc; };
G4double G4ComponentGGHadronNucleusXsc::GetTotalIsotopeCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4int  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 98 of file G4ComponentGGHadronNucleusXsc.cc.

References GetIsoCrossSection().

101 {
102  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
103  kinEnergy);
104  fTotalXsc = GetIsoCrossSection(aDP, Z, A);
105  delete aDP;
106 
107  return fTotalXsc;
108 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4ThreeVector G4ParticleMomentum
G4bool G4ComponentGGHadronNucleusXsc::IsIsoApplicable ( const G4DynamicParticle aDP,
G4int  Z,
G4int  A,
const G4Element elm = 0,
const G4Material mat = 0 
)

Definition at line 234 of file G4ComponentGGHadronNucleusXsc.cc.

References G4DynamicParticle::GetDefinition(), and G4DynamicParticle::GetKineticEnergy().

238 {
239  G4bool applicable = false;
240  // G4int baryonNumber = aDP->GetDefinition()->GetBaryonNumber();
241  G4double kineticEnergy = aDP->GetKineticEnergy();
242 
243  const G4ParticleDefinition* theParticle = aDP->GetDefinition();
244 
245  if ( ( kineticEnergy >= fLowerLimit &&
246  Z > 1 && // >= He
247  ( theParticle == theAProton ||
248  theParticle == theGamma ||
249  theParticle == theKPlus ||
250  theParticle == theKMinus ||
251  theParticle == theK0L ||
252  theParticle == theK0S ||
253  theParticle == theSMinus ||
254  theParticle == theProton ||
255  theParticle == theNeutron ||
256  theParticle == thePiPlus ||
257  theParticle == thePiMinus ) ) ) applicable = true;
258 
259  return applicable;
260 }
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
bool G4bool
Definition: G4Types.hh:79
double G4double
Definition: G4Types.hh:76
void G4ComponentGGHadronNucleusXsc::SetEnergyLowerLimit ( G4double  E)
inline

Definition at line 156 of file G4ComponentGGHadronNucleusXsc.hh.

156 {fLowerLimit=E;};

The documentation for this class was generated from the following files: