Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ComponentGGHadronNucleusXsc Class Reference

#include <G4ComponentGGHadronNucleusXsc.hh>

Inheritance diagram for G4ComponentGGHadronNucleusXsc:
Collaboration diagram for G4ComponentGGHadronNucleusXsc:

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 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)
 
G4double GetAxsc2piR2 ()
 
G4double GetModelInLog ()
 
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
 

Static Public Member Functions

static const char * Default_Name ()
 

Detailed Description

Definition at line 51 of file G4ComponentGGHadronNucleusXsc.hh.

Constructor & Destructor Documentation

G4ComponentGGHadronNucleusXsc::G4ComponentGGHadronNucleusXsc ( )

Definition at line 46 of file G4ComponentGGHadronNucleusXsc.cc.

48 // fUpperLimit(100000*GeV),
49  fLowerLimit(10.*MeV),// fLowerLimit(3*GeV),
50  fRadiusConst(1.08*fermi), // 1.1, 1.3 ?
51  fTotalXsc(0.0), fElasticXsc(0.0), fInelasticXsc(0.0), fProductionXsc(0.0),
52  fDiffractionXsc(0.0), fAxsc2piR2(0.0),fModelInLog(0.0)
53 // , fHadronNucleonXsc(0.0)
54 {
55  theGamma = G4Gamma::Gamma();
56  theProton = G4Proton::Proton();
57  theNeutron = G4Neutron::Neutron();
58  theAProton = G4AntiProton::AntiProton();
59  theANeutron = G4AntiNeutron::AntiNeutron();
60  thePiPlus = G4PionPlus::PionPlus();
61  thePiMinus = G4PionMinus::PionMinus();
62  thePiZero = G4PionZero::PionZero();
63  theKPlus = G4KaonPlus::KaonPlus();
64  theKMinus = G4KaonMinus::KaonMinus();
67  theL = G4Lambda::Lambda();
68  theAntiL = G4AntiLambda::AntiLambda();
69  theSPlus = G4SigmaPlus::SigmaPlus();
70  theASPlus = G4AntiSigmaPlus::AntiSigmaPlus();
71  theSMinus = G4SigmaMinus::SigmaMinus();
72  theASMinus = G4AntiSigmaMinus::AntiSigmaMinus();
73  theS0 = G4SigmaZero::SigmaZero();
75  theXiMinus = G4XiMinus::XiMinus();
76  theXi0 = G4XiZero::XiZero();
77  theAXiMinus = G4AntiXiMinus::AntiXiMinus();
78  theAXi0 = G4AntiXiZero::AntiXiZero();
79  theOmega = G4OmegaMinus::OmegaMinus();
81  theD = G4Deuteron::Deuteron();
82  theT = G4Triton::Triton();
83  theA = G4Alpha::Alpha();
84  theHe3 = G4He3::He3();
85 
86  hnXsc = new G4HadronNucleonXsc();
87 }
static G4AntiOmegaMinus * AntiOmegaMinus()
static G4OmegaMinus * OmegaMinus()
static G4KaonZeroLong * KaonZeroLong()
static G4AntiSigmaPlus * AntiSigmaPlus()
static G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:102
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:108
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 constexpr double MeV
Definition: G4SIunits.hh:214
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 constexpr double fermi
Definition: G4SIunits.hh:103
static G4He3 * He3()
Definition: G4He3.cc:94
static G4AntiNeutron * AntiNeutron()

Here is the call graph for this function:

G4ComponentGGHadronNucleusXsc::~G4ComponentGGHadronNucleusXsc ( )
virtual

Definition at line 93 of file G4ComponentGGHadronNucleusXsc.cc.

94 {
95  if (hnXsc) delete hnXsc;
96 }

Member Function Documentation

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

Definition at line 1505 of file G4ComponentGGHadronNucleusXsc.cc.

1508 {
1509  G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
1510  G4double sMand = mp*mp + mt*mt + 2*Elab*mt ;
1511 
1512  return sMand;
1513 }
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

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

Definition at line 1489 of file G4ComponentGGHadronNucleusXsc.cc.

1492 {
1493  G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
1494  G4double Ecm = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
1495  // G4double Pcm = Plab * mt / Ecm;
1496  // G4double KEcm = std::sqrt ( Pcm * Pcm + mp * mp ) - mp;
1497 
1498  return Ecm ; // KEcm;
1499 }
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 212 of file G4ComponentGGHadronNucleusXsc.cc.

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

Here is the call graph for this function:

void G4ComponentGGHadronNucleusXsc::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Definition at line 1519 of file G4ComponentGGHadronNucleusXsc.cc.

1520 {
1521  outFile << "G4ComponentGGHadronNucleusXsc calculates total, inelastic and\n"
1522  << "elastic cross sections for hadron-nucleus cross sections using\n"
1523  << "the Glauber model with Gribov corrections. It is valid for all\n"
1524  << "targets except hydrogen, and for incident p, pbar, n, sigma-,\n"
1525  << "pi+, pi-, K+, K- and gammas with energies above 3 GeV. This is\n"
1526  << "a cross section component which is to be used to build a cross\n"
1527  << "data set.\n";
1528 }
static const char* G4ComponentGGHadronNucleusXsc::Default_Name ( )
inlinestatic

Definition at line 58 of file G4ComponentGGHadronNucleusXsc.hh.

58 { return "Glauber-Gribov"; }
G4double G4ComponentGGHadronNucleusXsc::GetAxsc2piR2 ( )
inline

Definition at line 143 of file G4ComponentGGHadronNucleusXsc.hh.

143 {return fAxsc2piR2;};
G4double G4ComponentGGHadronNucleusXsc::GetDiffractionGlauberGribovXsc ( )
inline

Definition at line 155 of file G4ComponentGGHadronNucleusXsc.hh.

155 { return fDiffractionXsc; };

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetElasticElementCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4double  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 184 of file G4ComponentGGHadronNucleusXsc.cc.

187 {
188  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
189  kinEnergy);
190  fTotalXsc = GetIsoCrossSection(aDP, Z, G4int(A));
191  delete aDP;
192 
193  return fElasticXsc;
194 }
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)
double A(double temperature)
G4ThreeVector G4ParticleMomentum

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetElasticGlauberGribov ( const G4DynamicParticle dp,
G4int  Z,
G4int  A 
)
inline

Definition at line 226 of file G4ComponentGGHadronNucleusXsc.hh.

228 {
229  GetIsoCrossSection(dp, Z, A);
230  return fElasticXsc;
231 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
double A(double temperature)

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetElasticGlauberGribovXsc ( )
inline

Definition at line 152 of file G4ComponentGGHadronNucleusXsc.hh.

152 { return fElasticXsc; };

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetElasticIsotopeCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4int  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 198 of file G4ComponentGGHadronNucleusXsc.cc.

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

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXsc ( const G4DynamicParticle aParticle,
const G4Element anElement 
)

Definition at line 517 of file G4ComponentGGHadronNucleusXsc.cc.

519 {
520  G4int At = G4lrint(anElement->GetN()); // number of nucleons
521  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
522 
523  return GetHadronNucleonXsc(aParticle, At, Zt);
524 }
G4double GetN() const
Definition: G4Element.hh:135
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

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXsc ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 534 of file G4ComponentGGHadronNucleusXsc.cc.

536 {
537  G4double xsection;
538 
539  //G4double targ_mass = G4NucleiProperties::GetNuclearMass(At, Zt);
540 
541  G4double targ_mass = 0.939*GeV; // ~mean neutron and proton ???
542 
543  G4double proj_mass = aParticle->GetMass();
544  G4double proj_momentum = aParticle->GetMomentum().mag();
545  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
546 
547  sMand /= GeV*GeV; // in GeV for parametrisation
548  proj_momentum /= GeV;
549 
550  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
551 
552  G4double aa = At;
553 
554  if(theParticle == theGamma)
555  {
556  xsection = aa*(0.0677*G4Pow::GetInstance()->powA(sMand,0.0808) + 0.129*G4Pow::GetInstance()->powA(sMand,-0.4525));
557  }
558  else if(theParticle == theNeutron) // as proton ???
559  {
560  xsection = aa*(21.70*G4Pow::GetInstance()->powA(sMand,0.0808) + 56.08*G4Pow::GetInstance()->powA(sMand,-0.4525));
561  }
562  else if(theParticle == theProton)
563  {
564  xsection = aa*(21.70*G4Pow::GetInstance()->powA(sMand,0.0808) + 56.08*G4Pow::GetInstance()->powA(sMand,-0.4525));
565  // xsection = At*( 49.51*G4Pow::GetInstance()->powA(sMand,-0.097) + 0.314*G4Log(sMand)*G4Log(sMand) );
566  // xsection = At*( 38.4 + 0.85*std::abs(G4Pow::GetInstance()->powA(log(sMand),1.47)) );
567  }
568  else if(theParticle == theAProton)
569  {
570  xsection = aa*( 21.70*G4Pow::GetInstance()->powA(sMand,0.0808) + 98.39*G4Pow::GetInstance()->powA(sMand,-0.4525));
571  }
572  else if(theParticle == thePiPlus)
573  {
574  xsection = aa*(13.63*G4Pow::GetInstance()->powA(sMand,0.0808) + 27.56*G4Pow::GetInstance()->powA(sMand,-0.4525));
575  }
576  else if(theParticle == thePiMinus)
577  {
578  // xsection = At*( 55.2*G4Pow::GetInstance()->powA(sMand,-0.255) + 0.346*G4Log(sMand)*G4Log(sMand) );
579  xsection = aa*(13.63*G4Pow::GetInstance()->powA(sMand,0.0808) + 36.02*G4Pow::GetInstance()->powA(sMand,-0.4525));
580  }
581  else if(theParticle == theKPlus)
582  {
583  xsection = aa*(11.82*G4Pow::GetInstance()->powA(sMand,0.0808) + 8.15*G4Pow::GetInstance()->powA(sMand,-0.4525));
584  }
585  else if(theParticle == theKMinus)
586  {
587  xsection = aa*(11.82*G4Pow::GetInstance()->powA(sMand,0.0808) + 26.36*G4Pow::GetInstance()->powA(sMand,-0.4525));
588  }
589  else // as proton ???
590  {
591  xsection = aa*(21.70*G4Pow::GetInstance()->powA(sMand,0.0808) + 56.08*G4Pow::GetInstance()->powA(sMand,-0.4525));
592  }
593  xsection *= millibarn;
594  return xsection;
595 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
G4double CalcMandelstamS(const G4double, const G4double, const G4double)
G4ParticleDefinition * GetDefinition() const
G4double GetMass() const
static constexpr double GeV
Definition: G4SIunits.hh:217
double G4double
Definition: G4Types.hh:76
double mag() const
static constexpr double millibarn
Definition: G4SIunits.hh:106
G4ThreeVector GetMomentum() const

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscNS ( const G4DynamicParticle aParticle,
const G4Element anElement 
)

Definition at line 738 of file G4ComponentGGHadronNucleusXsc.cc.

740 {
741  G4int At = G4lrint(anElement->GetN()); // number of nucleons
742  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
743 
744  return GetHadronNucleonXscNS(aParticle, At, Zt);
745 }
G4double GetN() const
Definition: G4Element.hh:135
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

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscNS ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 756 of file G4ComponentGGHadronNucleusXsc.cc.

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

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscPDG ( const G4DynamicParticle aParticle,
const G4Element anElement 
)

Definition at line 604 of file G4ComponentGGHadronNucleusXsc.cc.

606 {
607  G4int At = G4lrint(anElement->GetN()); // number of nucleons
608  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
609 
610  return GetHadronNucleonXscPDG(aParticle, At, Zt);
611 }
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4Element *)
G4double GetN() const
Definition: G4Element.hh:135
G4double GetZ() const
Definition: G4Element.hh:131
int G4int
Definition: G4Types.hh:78
int G4lrint(double ad)
Definition: templates.hh:163

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscPDG ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 623 of file G4ComponentGGHadronNucleusXsc.cc.

625 {
626  G4double xsection;
627 
628  G4int Nt = At-Zt; // number of neutrons
629  if (Nt < 0) Nt = 0;
630 
631  G4double zz = Zt;
632  G4double aa = At;
633  G4double nn = Nt;
634 
636  GetIonTable()->GetIonMass(Zt, At);
637 
638  targ_mass = 0.939*GeV; // ~mean neutron and proton ???
639 
640  G4double proj_mass = aParticle->GetMass();
641  G4double proj_momentum = aParticle->GetMomentum().mag();
642 
643  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
644 
645  sMand /= GeV*GeV; // in GeV for parametrisation
646 
647  // General PDG fit constants
648 
649  G4double s0 = 5.38*5.38; // in Gev^2
650  G4double eta1 = 0.458;
651  G4double eta2 = 0.458;
652  G4double B = 0.308;
653 
654 
655  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
656 
657 
658  if(theParticle == theNeutron) // proton-neutron fit
659  {
660  xsection = zz*( 35.80 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
661  + 40.15*G4Pow::GetInstance()->powA(sMand,-eta1) - 30.*G4Pow::GetInstance()->powA(sMand,-eta2));
662  xsection += nn*( 35.45 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
663  + 42.53*G4Pow::GetInstance()->powA(sMand,-eta1) - 33.34*G4Pow::GetInstance()->powA(sMand,-eta2)); // pp for nn
664  }
665  else if(theParticle == theProton)
666  {
667 
668  xsection = zz*( 35.45 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
669  + 42.53*G4Pow::GetInstance()->powA(sMand,-eta1) - 33.34*G4Pow::GetInstance()->powA(sMand,-eta2));
670 
671  xsection += nn*( 35.80 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
672  + 40.15*G4Pow::GetInstance()->powA(sMand,-eta1) - 30.*G4Pow::GetInstance()->powA(sMand,-eta2));
673  }
674  else if(theParticle == theAProton)
675  {
676  xsection = zz*( 35.45 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
677  + 42.53*G4Pow::GetInstance()->powA(sMand,-eta1) + 33.34*G4Pow::GetInstance()->powA(sMand,-eta2));
678 
679  xsection += nn*( 35.80 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
680  + 40.15*G4Pow::GetInstance()->powA(sMand,-eta1) + 30.*G4Pow::GetInstance()->powA(sMand,-eta2));
681  }
682  else if(theParticle == thePiPlus)
683  {
684  xsection = aa*( 20.86 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
685  + 19.24*G4Pow::GetInstance()->powA(sMand,-eta1) - 6.03*G4Pow::GetInstance()->powA(sMand,-eta2));
686  }
687  else if(theParticle == thePiMinus)
688  {
689  xsection = aa*( 20.86 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
690  + 19.24*G4Pow::GetInstance()->powA(sMand,-eta1) + 6.03*G4Pow::GetInstance()->powA(sMand,-eta2));
691  }
692  else if(theParticle == theKPlus || theParticle == theK0L )
693  {
694  xsection = zz*( 17.91 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
695  + 7.14*G4Pow::GetInstance()->powA(sMand,-eta1) - 13.45*G4Pow::GetInstance()->powA(sMand,-eta2));
696 
697  xsection += nn*( 17.87 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
698  + 5.17*G4Pow::GetInstance()->powA(sMand,-eta1) - 7.23*G4Pow::GetInstance()->powA(sMand,-eta2));
699  }
700  else if(theParticle == theKMinus || theParticle == theK0S )
701  {
702  xsection = zz*( 17.91 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
703  + 7.14*G4Pow::GetInstance()->powA(sMand,-eta1) + 13.45*G4Pow::GetInstance()->powA(sMand,-eta2));
704 
705  xsection += nn*( 17.87 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
706  + 5.17*G4Pow::GetInstance()->powA(sMand,-eta1) + 7.23*G4Pow::GetInstance()->powA(sMand,-eta2));
707  }
708  else if(theParticle == theSMinus)
709  {
710  xsection = aa*( 35.20 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
711  - 199.*G4Pow::GetInstance()->powA(sMand,-eta1) + 264.*G4Pow::GetInstance()->powA(sMand,-eta2));
712  }
713  else if(theParticle == theGamma) // modify later on
714  {
715  xsection = aa*( 0.0 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
716  + 0.032*G4Pow::GetInstance()->powA(sMand,-eta1) - 0.0*G4Pow::GetInstance()->powA(sMand,-eta2));
717 
718  }
719  else // as proton ???
720  {
721  xsection = zz*( 35.45 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
722  + 42.53*G4Pow::GetInstance()->powA(sMand,-eta1) - 33.34*G4Pow::GetInstance()->powA(sMand,-eta2));
723 
724  xsection += nn*( 35.80 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
725  + 40.15*G4Pow::GetInstance()->powA(sMand,-eta1) - 30.*G4Pow::GetInstance()->powA(sMand,-eta2));
726  }
727  xsection *= millibarn; // parametrised in mb
728  return xsection;
729 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
G4double CalcMandelstamS(const G4double, const G4double, const G4double)
G4ParticleDefinition * GetDefinition() const
double B(double temperature)
int G4int
Definition: G4Types.hh:78
G4double GetMass() const
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static G4ParticleTable * GetParticleTable()
static constexpr double GeV
Definition: G4SIunits.hh:217
double G4double
Definition: G4Types.hh:76
double mag() const
static constexpr double millibarn
Definition: G4SIunits.hh:106
G4ThreeVector GetMomentum() const

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHNinelasticXsc ( const G4DynamicParticle aParticle,
const G4Element anElement 
)

Definition at line 1135 of file G4ComponentGGHadronNucleusXsc.cc.

1137 {
1138  G4int At = G4lrint(anElement->GetN()); // number of nucleons
1139  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
1140 
1141  return GetHNinelasticXsc(aParticle, At, Zt);
1142 }
G4double GetN() const
Definition: G4Element.hh:135
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

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHNinelasticXsc ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 1149 of file G4ComponentGGHadronNucleusXsc.cc.

1151 {
1152  const G4ParticleDefinition* hadron = aParticle->GetDefinition();
1153  G4double sumInelastic;
1154  G4int Nt = At - Zt;
1155  if(Nt < 0) Nt = 0;
1156 
1157  if( hadron == theKPlus )
1158  {
1159  sumInelastic = GetHNinelasticXscVU(aParticle, At, Zt);
1160  }
1161  else
1162  {
1163  //sumInelastic = Zt*GetHadronNucleonXscMK(aParticle, theProton);
1164  // sumInelastic += Nt*GetHadronNucleonXscMK(aParticle, theNeutron);
1165  sumInelastic = G4double(Zt)*GetHadronNucleonXscNS(aParticle, 1, 1);
1166  sumInelastic += G4double(Nt)*GetHadronNucleonXscNS(aParticle, 1, 0);
1167  }
1168  return sumInelastic;
1169 }
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

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHNinelasticXscVU ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 1177 of file G4ComponentGGHadronNucleusXsc.cc.

1179 {
1180  G4int PDGcode = aParticle->GetDefinition()->GetPDGEncoding();
1181  G4int absPDGcode = std::abs(PDGcode);
1182 
1183  G4double Elab = aParticle->GetTotalEnergy();
1184  // (s - 2*0.88*GeV*GeV)/(2*0.939*GeV)/GeV;
1185  G4double Plab = aParticle->GetMomentum().mag();
1186  // std::sqrt(Elab * Elab - 0.88);
1187 
1188  Elab /= GeV;
1189  Plab /= GeV;
1190 
1191  G4double LogPlab = G4Log( Plab );
1192  G4double sqrLogPlab = LogPlab * LogPlab;
1193 
1194  //G4cout<<"Plab = "<<Plab<<G4endl;
1195 
1196  G4double NumberOfTargetProtons = G4double(Zt);
1197  G4double NumberOfTargetNucleons = G4double(At);
1198  G4double NumberOfTargetNeutrons = NumberOfTargetNucleons - NumberOfTargetProtons;
1199 
1200  if(NumberOfTargetNeutrons < 0.0) NumberOfTargetNeutrons = 0.0;
1201 
1202  G4double Xtotal, Xelastic, Xinelastic;
1203 
1204  if( absPDGcode > 1000 ) //------Projectile is baryon --------
1205  {
1206  G4double XtotPP = 48.0 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1207  0.522*sqrLogPlab - 4.51*LogPlab;
1208 
1209  G4double XtotPN = 47.3 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1210  0.513*sqrLogPlab - 4.27*LogPlab;
1211 
1212  G4double XelPP = 11.9 + 26.9*G4Pow::GetInstance()->powA(Plab,-1.21) +
1213  0.169*sqrLogPlab - 1.85*LogPlab;
1214 
1215  G4double XelPN = 11.9 + 26.9*G4Pow::GetInstance()->powA(Plab,-1.21) +
1216  0.169*sqrLogPlab - 1.85*LogPlab;
1217 
1218  Xtotal = (NumberOfTargetProtons * XtotPP +
1219  NumberOfTargetNeutrons * XtotPN);
1220 
1221  Xelastic = (NumberOfTargetProtons * XelPP +
1222  NumberOfTargetNeutrons * XelPN);
1223  }
1224  else if( PDGcode == 211 ) //------Projectile is PionPlus -------
1225  {
1226  G4double XtotPiP = 16.4 + 19.3 *G4Pow::GetInstance()->powA(Plab,-0.42) +
1227  0.19 *sqrLogPlab - 0.0 *LogPlab;
1228 
1229  G4double XtotPiN = 33.0 + 14.0 *G4Pow::GetInstance()->powA(Plab,-1.36) +
1230  0.456*sqrLogPlab - 4.03*LogPlab;
1231 
1232  G4double XelPiP = 0.0 + 11.4*G4Pow::GetInstance()->powA(Plab,-0.40) +
1233  0.079*sqrLogPlab - 0.0 *LogPlab;
1234 
1235  G4double XelPiN = 1.76 + 11.2*G4Pow::GetInstance()->powA(Plab,-0.64) +
1236  0.043*sqrLogPlab - 0.0 *LogPlab;
1237 
1238  Xtotal = ( NumberOfTargetProtons * XtotPiP +
1239  NumberOfTargetNeutrons * XtotPiN );
1240 
1241  Xelastic = ( NumberOfTargetProtons * XelPiP +
1242  NumberOfTargetNeutrons * XelPiN );
1243  }
1244  else if( PDGcode == -211 ) //------Projectile is PionMinus -------
1245  {
1246  G4double XtotPiP = 33.0 + 14.0 *G4Pow::GetInstance()->powA(Plab,-1.36) +
1247  0.456*sqrLogPlab - 4.03*LogPlab;
1248 
1249  G4double XtotPiN = 16.4 + 19.3 *G4Pow::GetInstance()->powA(Plab,-0.42) +
1250  0.19 *sqrLogPlab - 0.0 *LogPlab;
1251 
1252  G4double XelPiP = 1.76 + 11.2*G4Pow::GetInstance()->powA(Plab,-0.64) +
1253  0.043*sqrLogPlab - 0.0 *LogPlab;
1254 
1255  G4double XelPiN = 0.0 + 11.4*G4Pow::GetInstance()->powA(Plab,-0.40) +
1256  0.079*sqrLogPlab - 0.0 *LogPlab;
1257 
1258  Xtotal = ( NumberOfTargetProtons * XtotPiP +
1259  NumberOfTargetNeutrons * XtotPiN );
1260 
1261  Xelastic = ( NumberOfTargetProtons * XelPiP +
1262  NumberOfTargetNeutrons * XelPiN );
1263  }
1264  else if( PDGcode == 111 ) //------Projectile is PionZero -------
1265  {
1266  G4double XtotPiP =(16.4 + 19.3 *G4Pow::GetInstance()->powA(Plab,-0.42) +
1267  0.19 *sqrLogPlab - 0.0 *LogPlab + //Pi+
1268  33.0 + 14.0 *G4Pow::GetInstance()->powA(Plab,-1.36) +
1269  0.456*sqrLogPlab - 4.03*LogPlab)/2; //Pi-
1270 
1271  G4double XtotPiN =(33.0 + 14.0 *G4Pow::GetInstance()->powA(Plab,-1.36) +
1272  0.456*sqrLogPlab - 4.03*LogPlab + //Pi+
1273  16.4 + 19.3 *G4Pow::GetInstance()->powA(Plab,-0.42) +
1274  0.19 *sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1275 
1276  G4double XelPiP =( 0.0 + 11.4*G4Pow::GetInstance()->powA(Plab,-0.40) +
1277  0.079*sqrLogPlab - 0.0 *LogPlab + //Pi+
1278  1.76 + 11.2*G4Pow::GetInstance()->powA(Plab,-0.64) +
1279  0.043*sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1280 
1281  G4double XelPiN =( 1.76 + 11.2*G4Pow::GetInstance()->powA(Plab,-0.64) +
1282  0.043*sqrLogPlab - 0.0 *LogPlab + //Pi+
1283  0.0 + 11.4*G4Pow::GetInstance()->powA(Plab,-0.40) +
1284  0.079*sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1285 
1286  Xtotal = ( NumberOfTargetProtons * XtotPiP +
1287  NumberOfTargetNeutrons * XtotPiN );
1288 
1289  Xelastic = ( NumberOfTargetProtons * XelPiP +
1290  NumberOfTargetNeutrons * XelPiN );
1291  }
1292  else if( PDGcode == 321 ) //------Projectile is KaonPlus -------
1293  {
1294  G4double XtotKP = 18.1 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1295  0.26 *sqrLogPlab - 1.0 *LogPlab;
1296  G4double XtotKN = 18.7 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1297  0.21 *sqrLogPlab - 0.89*LogPlab;
1298 
1299  G4double XelKP = 5.0 + 8.1*G4Pow::GetInstance()->powA(Plab,-1.8 ) +
1300  0.16 *sqrLogPlab - 1.3 *LogPlab;
1301 
1302  G4double XelKN = 7.3 + 0. *G4Pow::GetInstance()->powA(Plab,-0. ) +
1303  0.29 *sqrLogPlab - 2.4 *LogPlab;
1304 
1305  Xtotal = ( NumberOfTargetProtons * XtotKP +
1306  NumberOfTargetNeutrons * XtotKN );
1307 
1308  Xelastic = ( NumberOfTargetProtons * XelKP +
1309  NumberOfTargetNeutrons * XelKN );
1310  }
1311  else if( PDGcode ==-321 ) //------Projectile is KaonMinus ------
1312  {
1313  G4double XtotKP = 32.1 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1314  0.66 *sqrLogPlab - 5.6 *LogPlab;
1315  G4double XtotKN = 25.2 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1316  0.38 *sqrLogPlab - 2.9 *LogPlab;
1317 
1318  G4double XelKP = 7.3 + 0. *G4Pow::GetInstance()->powA(Plab,-0. ) +
1319  0.29 *sqrLogPlab - 2.4 *LogPlab;
1320 
1321  G4double XelKN = 5.0 + 8.1*G4Pow::GetInstance()->powA(Plab,-1.8 ) +
1322  0.16 *sqrLogPlab - 1.3 *LogPlab;
1323 
1324  Xtotal = ( NumberOfTargetProtons * XtotKP +
1325  NumberOfTargetNeutrons * XtotKN );
1326 
1327  Xelastic = ( NumberOfTargetProtons * XelKP +
1328  NumberOfTargetNeutrons * XelKN );
1329  }
1330  else if( PDGcode == 311 ) //------Projectile is KaonZero ------
1331  {
1332  G4double XtotKP = ( 18.1 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1333  0.26 *sqrLogPlab - 1.0 *LogPlab + //K+
1334  32.1 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1335  0.66 *sqrLogPlab - 5.6 *LogPlab)/2; //K-
1336 
1337  G4double XtotKN = ( 18.7 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1338  0.21 *sqrLogPlab - 0.89*LogPlab + //K+
1339  25.2 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1340  0.38 *sqrLogPlab - 2.9 *LogPlab)/2; //K-
1341 
1342  G4double XelKP = ( 5.0 + 8.1*G4Pow::GetInstance()->powA(Plab,-1.8 )
1343  + 0.16 *sqrLogPlab - 1.3 *LogPlab + //K+
1344  7.3 + 0. *G4Pow::GetInstance()->powA(Plab,-0. ) +
1345  0.29 *sqrLogPlab - 2.4 *LogPlab)/2; //K-
1346 
1347  G4double XelKN = ( 7.3 + 0. *G4Pow::GetInstance()->powA(Plab,-0. ) +
1348  0.29 *sqrLogPlab - 2.4 *LogPlab + //K+
1349  5.0 + 8.1*G4Pow::GetInstance()->powA(Plab,-1.8 ) +
1350  0.16 *sqrLogPlab - 1.3 *LogPlab)/2; //K-
1351 
1352  Xtotal = ( NumberOfTargetProtons * XtotKP +
1353  NumberOfTargetNeutrons * XtotKN );
1354 
1355  Xelastic = ( NumberOfTargetProtons * XelKP +
1356  NumberOfTargetNeutrons * XelKN );
1357  }
1358  else //------Projectile is undefined, Nucleon assumed
1359  {
1360  G4double XtotPP = 48.0 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1361  0.522*sqrLogPlab - 4.51*LogPlab;
1362 
1363  G4double XtotPN = 47.3 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1364  0.513*sqrLogPlab - 4.27*LogPlab;
1365 
1366  G4double XelPP = 11.9 + 26.9*G4Pow::GetInstance()->powA(Plab,-1.21) +
1367  0.169*sqrLogPlab - 1.85*LogPlab;
1368  G4double XelPN = 11.9 + 26.9*G4Pow::GetInstance()->powA(Plab,-1.21) +
1369  0.169*sqrLogPlab - 1.85*LogPlab;
1370 
1371  Xtotal = ( NumberOfTargetProtons * XtotPP +
1372  NumberOfTargetNeutrons * XtotPN );
1373 
1374  Xelastic = ( NumberOfTargetProtons * XelPP +
1375  NumberOfTargetNeutrons * XelPN );
1376  }
1377  Xinelastic = Xtotal - Xelastic;
1378 
1379  if( Xinelastic < 0.) Xinelastic = 0.;
1380 
1381  return Xinelastic*= millibarn;
1382 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
G4double GetTotalEnergy() const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static constexpr double GeV
Definition: G4SIunits.hh:217
double G4double
Definition: G4Types.hh:76
double mag() const
static constexpr double millibarn
Definition: G4SIunits.hh:106
G4ThreeVector GetMomentum() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetInelasticElementCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4double  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 156 of file G4ComponentGGHadronNucleusXsc.cc.

159 {
160  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
161  kinEnergy);
162  fTotalXsc = GetIsoCrossSection(aDP, Z, G4int(A));
163  delete aDP;
164 
165  return fInelasticXsc;
166 }
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)
double A(double temperature)
G4ThreeVector G4ParticleMomentum

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetInelasticGlauberGribov ( const G4DynamicParticle dp,
G4int  Z,
G4int  A 
)
inline

Definition at line 237 of file G4ComponentGGHadronNucleusXsc.hh.

239 {
240  GetIsoCrossSection(dp, Z, A);
241  return fInelasticXsc;
242 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
double A(double temperature)

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetInelasticGlauberGribovXsc ( )
inline

Definition at line 153 of file G4ComponentGGHadronNucleusXsc.hh.

153 { return fInelasticXsc; };

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetInelasticIsotopeCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4int  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 128 of file G4ComponentGGHadronNucleusXsc.cc.

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

Here is the call graph for this function:

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 294 of file G4ComponentGGHadronNucleusXsc.cc.

299 {
300  G4double xsection, sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
301  G4double hpInXsc(0.), hnInXsc(0.);
303 
304  G4int N = A - Z; // number of neutrons
305  if (N < 0) N = 0;
306 
307  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
308 
309  if( theParticle == theProton ||
310  theParticle == theNeutron ||
311  theParticle == thePiPlus ||
312  theParticle == thePiMinus )
313  {
314  // sigma = GetHadronNucleonXscNS(aParticle, A, Z);
315 
316  sigma = Z*hnXsc->GetHadronNucleonXscNS(aParticle, theProton);
317 
318  hpInXsc = hnXsc->GetInelasticHadronNucleonXsc();
319 
320  sigma += N*hnXsc->GetHadronNucleonXscNS(aParticle, theNeutron);
321 
322  hnInXsc = hnXsc->GetInelasticHadronNucleonXsc();
323 
324  cofInelastic = 2.4;
325  cofTotal = 2.0;
326  }
327  else if( theParticle == theKPlus ||
328  theParticle == theKMinus ||
329  theParticle == theK0S ||
330  theParticle == theK0L )
331  {
332  // sigma = GetKaonNucleonXscVector(aParticle, A, Z);
333 
334  sigma = Z*hnXsc->GetKaonNucleonXscGG(aParticle, theProton);
335 
336  hpInXsc = hnXsc->GetInelasticHadronNucleonXsc();
337 
338  sigma += N*hnXsc->GetKaonNucleonXscGG(aParticle, theNeutron);
339 
340  hnInXsc = hnXsc->GetInelasticHadronNucleonXsc();
341 
342  cofInelastic = 2.2;
343  cofTotal = 2.0;
344  R = 1.3*fermi;
345  R *= G4Pow::GetInstance()->powA(G4double(A), 0.3333);
346  }
347  else
348  {
349  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
350  cofInelastic = 2.2;
351  cofTotal = 2.0;
352  }
353  // cofInelastic = 2.0;
354 
355  if( A > 1 )
356  {
357  nucleusSquare = cofTotal*pi*R*R; // basically 2piRR
358  ratio = sigma/nucleusSquare;
359 
360  xsection = nucleusSquare*G4Log( 1. + ratio );
361 
362  xsection *= GetParticleBarCorTot(theParticle, Z);
363 
364  fTotalXsc = xsection;
365 
366  // inelastic xsc
367 
368  fAxsc2piR2 = cofInelastic*ratio;
369 
370  fModelInLog = G4Log( 1. + fAxsc2piR2 );
371 
372  fInelasticXsc = nucleusSquare*fModelInLog/cofInelastic;
373 
374  fInelasticXsc *= GetParticleBarCorIn(theParticle, Z);
375 
376  fElasticXsc = fTotalXsc - fInelasticXsc;
377 
378  if( fElasticXsc < 0. ) fElasticXsc = 0.;
379 
380  G4double difratio = ratio/(1.+ratio);
381 
382  fDiffractionXsc = 0.5*nucleusSquare*( difratio - G4Log( 1. + difratio ) );
383 
384 
385  // sigma = GetHNinelasticXsc(aParticle, A, Z);
386 
387  sigma = Z*hpInXsc + N*hnInXsc;
388 
389  ratio = sigma/nucleusSquare;
390 
391  fProductionXsc = nucleusSquare*G4Log( 1. + cofInelastic*ratio )/cofInelastic;
392 
393  fProductionXsc *= GetParticleBarCorIn(theParticle, Z);
394 
395  if (fElasticXsc < 0.) fElasticXsc = 0.;
396  }
397  else // H
398  {
399  if( theParticle == theKPlus ||
400  theParticle == theKMinus ||
401  theParticle == theK0S ||
402  theParticle == theK0L )
403  {
404  fTotalXsc = hnXsc->GetHadronNucleonXscNS(aParticle, theProton);
405  xsection = fTotalXsc;
406  fInelasticXsc = hnXsc->GetInelasticHadronNucleonXsc();
407  fElasticXsc = hnXsc->GetElasticHadronNucleonXsc();
408  }
409  else
410  {
411  fTotalXsc = sigma;
412  xsection = sigma;
413 
414  fInelasticXsc = hpInXsc;
415  fElasticXsc = fTotalXsc - fInelasticXsc;
416 
417  if (fElasticXsc < 0.) fElasticXsc = 0.;
418  }
419  }
420  return xsection;
421 }
G4double GetElasticHadronNucleonXsc()
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double GetNucleusRadius(const G4DynamicParticle *, const G4Element *)
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4ParticleDefinition *)
double A(double temperature)
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
G4double GetKaonNucleonXscGG(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double GetParticleBarCorIn(const G4ParticleDefinition *theParticle, G4int Z)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
**D E S C R I P T I O N
Definition: HEPEvtcom.cc:77
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76
static constexpr double fermi
Definition: G4SIunits.hh:103
G4double GetParticleBarCorTot(const G4ParticleDefinition *theParticle, G4int Z)
G4double GetInelasticHadronNucleonXsc()

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetModelInLog ( )
inline

Definition at line 144 of file G4ComponentGGHadronNucleusXsc.hh.

144 {return fModelInLog;};
G4double G4ComponentGGHadronNucleusXsc::GetNucleusRadius ( const G4DynamicParticle ,
const G4Element anElement 
)

Definition at line 1389 of file G4ComponentGGHadronNucleusXsc.cc.

1391 {
1392  G4int At = G4lrint(anElement->GetN());
1393  G4double oneThird = 1.0/3.0;
1394  G4double cubicrAt = G4Pow::GetInstance()->powA(G4double(At), oneThird);
1395 
1396  G4double R; // = fRadiusConst*cubicrAt;
1397  /*
1398  G4double tmp = G4Pow::GetInstance()->powA( cubicrAt-1., 3.);
1399  tmp += At;
1400  tmp *= 0.5;
1401 
1402  if (At > 20.) // 20.
1403  {
1404  R = fRadiusConst*G4Pow::GetInstance()->powA (tmp, oneThird);
1405  }
1406  else
1407  {
1408  R = fRadiusConst*cubicrAt;
1409  }
1410  */
1411 
1412  R = fRadiusConst*cubicrAt;
1413 
1414  G4double meanA = 21.;
1415 
1416  G4double tauA1 = 40.;
1417  G4double tauA2 = 10.;
1418  G4double tauA3 = 5.;
1419 
1420  G4double a1 = 0.85;
1421  G4double b1 = 1. - a1;
1422 
1423  G4double b2 = 0.3;
1424  G4double b3 = 4.;
1425 
1426  if (At > 20) // 20.
1427  {
1428  R *= ( a1 + b1*G4Exp( -(At - meanA)/tauA1) );
1429  }
1430  else if (At > 3)
1431  {
1432  R *= ( 1.0 + b2*( 1. - G4Exp( (At - meanA)/tauA2) ) );
1433  }
1434  else
1435  {
1436  R *= ( 1.0 + b3*( 1. - G4Exp( (At - meanA)/tauA3) ) );
1437  }
1438  return R;
1439 
1440 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
G4double GetN() const
Definition: G4Element.hh:135
int G4int
Definition: G4Types.hh:78
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
int G4lrint(double ad)
Definition: templates.hh:163
double G4double
Definition: G4Types.hh:76
const G4double oneThird

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetNucleusRadius ( G4int  At)

Definition at line 1446 of file G4ComponentGGHadronNucleusXsc.cc.

1447 {
1448  G4double oneThird = 1.0/3.0;
1449  G4double cubicrAt = G4Pow::GetInstance()->powA(G4double(At), oneThird);
1450 
1451  G4double R; // = fRadiusConst*cubicrAt;
1452 
1453  /*
1454  G4double tmp = G4Pow::GetInstance()->powA( cubicrAt-1., 3.);
1455  tmp += At;
1456  tmp *= 0.5;
1457 
1458  if (At > 20.)
1459  {
1460  R = fRadiusConst*G4Pow::GetInstance()->powA (tmp, oneThird);
1461  }
1462  else
1463  {
1464  R = fRadiusConst*cubicrAt;
1465  }
1466  */
1467 
1468  R = fRadiusConst*cubicrAt;
1469 
1470  G4double meanA = 20.;
1471  G4double tauA = 20.;
1472 
1473  if (At > 20) // 20.
1474  {
1475  R *= ( 0.8 + 0.2*G4Exp( -(G4double(At) - meanA)/tauA) );
1476  }
1477  else
1478  {
1479  R *= ( 1.0 + 0.1*( 1. - G4Exp( (G4double(At) - meanA)/tauA) ) );
1480  }
1481 
1482  return R;
1483 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
double G4double
Definition: G4Types.hh:76
const G4double oneThird

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetParticleBarCorIn ( const G4ParticleDefinition theParticle,
G4int  Z 
)
inline

Definition at line 270 of file G4ComponentGGHadronNucleusXsc.hh.

272 {
273  if(Z >= 2 && Z <= 92)
274  {
275  if( theParticle == theProton ) return fProtonBarCorrectionIn[Z];
276  else if( theParticle == theNeutron) return fNeutronBarCorrectionIn[Z];
277  else if( theParticle == thePiPlus ) return fPionPlusBarCorrectionIn[Z];
278  else if( theParticle == thePiMinus) return fPionMinusBarCorrectionIn[Z];
279  else return 1.0;
280  }
281  else return 1.0;
282 }

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetParticleBarCorTot ( const G4ParticleDefinition theParticle,
G4int  Z 
)
inline

Definition at line 250 of file G4ComponentGGHadronNucleusXsc.hh.

252 {
253  if(Z >= 2 && Z <= 92)
254  {
255  if( theParticle == theProton ) return fProtonBarCorrectionTot[Z];
256  else if( theParticle == theNeutron) return fNeutronBarCorrectionTot[Z];
257  else if( theParticle == thePiPlus ) return fPionPlusBarCorrectionTot[Z];
258  else if( theParticle == thePiMinus) return fPionMinusBarCorrectionTot[Z];
259  else return 1.0;
260  }
261  else return 1.0;
262 }

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetProductionElementCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4double  A 
)
virtual

Definition at line 170 of file G4ComponentGGHadronNucleusXsc.cc.

173 {
174  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
175  kinEnergy);
176  fTotalXsc = GetIsoCrossSection(aDP, Z, G4int(A));
177  delete aDP;
178 
179  return fProductionXsc;
180 }
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)
double A(double temperature)
G4ThreeVector G4ParticleMomentum

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetProductionGlauberGribovXsc ( )
inline

Definition at line 154 of file G4ComponentGGHadronNucleusXsc.hh.

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

Definition at line 142 of file G4ComponentGGHadronNucleusXsc.cc.

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

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetRadiusConst ( )
inline

Definition at line 156 of file G4ComponentGGHadronNucleusXsc.hh.

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

Definition at line 470 of file G4ComponentGGHadronNucleusXsc.cc.

471 {
472  G4double sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
474 
475  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
476 
477  if( theParticle == theProton ||
478  theParticle == theNeutron ||
479  theParticle == thePiPlus ||
480  theParticle == thePiMinus )
481  {
482  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
483  cofInelastic = 2.4;
484  cofTotal = 2.0;
485  }
486  else
487  {
488  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
489  cofInelastic = 2.2;
490  cofTotal = 2.0;
491  }
492  nucleusSquare = cofTotal*pi*R*R; // basically 2piRR
493  ratio = sigma/nucleusSquare;
494 
495  fInelasticXsc = nucleusSquare*G4Log( 1. + cofInelastic*ratio )/cofInelastic;
496 
497  sigma = GetHNinelasticXsc(aParticle, A, Z);
498  ratio = sigma/nucleusSquare;
499 
500  fProductionXsc = nucleusSquare*G4Log( 1. + cofInelastic*ratio )/cofInelastic;
501 
502  if (fInelasticXsc > fProductionXsc) ratio = (fInelasticXsc-fProductionXsc)/fInelasticXsc;
503  else ratio = 0.;
504  if ( ratio < 0. ) ratio = 0.;
505 
506  return ratio;
507 }
G4double GetNucleusRadius(const G4DynamicParticle *, const G4Element *)
G4ParticleDefinition * GetDefinition() const
double A(double temperature)
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
G4double GetHNinelasticXsc(const G4DynamicParticle *, const G4Element *)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetRatioSD ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 428 of file G4ComponentGGHadronNucleusXsc.cc.

429 {
430  G4double sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
432 
433  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
434 
435  if( theParticle == theProton ||
436  theParticle == theNeutron ||
437  theParticle == thePiPlus ||
438  theParticle == thePiMinus )
439  {
440  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
441  cofInelastic = 2.4;
442  cofTotal = 2.0;
443  }
444  else
445  {
446  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
447  cofInelastic = 2.2;
448  cofTotal = 2.0;
449  }
450  nucleusSquare = cofTotal*pi*R*R; // basically 2piRR
451  ratio = sigma/nucleusSquare;
452 
453  fInelasticXsc = nucleusSquare*G4Log( 1. + cofInelastic*ratio )/cofInelastic;
454 
455  G4double difratio = ratio/(1.+ratio);
456 
457  fDiffractionXsc = 0.5*nucleusSquare*( difratio - G4Log( 1. + difratio ) );
458 
459  if (fInelasticXsc > 0.) ratio = fDiffractionXsc/fInelasticXsc;
460  else ratio = 0.;
461 
462  return ratio;
463 }
G4double GetNucleusRadius(const G4DynamicParticle *, const G4Element *)
G4ParticleDefinition * GetDefinition() const
double A(double temperature)
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetTotalElementCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4double  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 114 of file G4ComponentGGHadronNucleusXsc.cc.

117 {
118  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
119  kinEnergy);
120  fTotalXsc = GetIsoCrossSection(aDP, Z, G4int(A));
121  delete aDP;
122 
123  return fTotalXsc;
124 }
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)
double A(double temperature)
G4ThreeVector G4ParticleMomentum

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetTotalGlauberGribovXsc ( )
inline

Definition at line 151 of file G4ComponentGGHadronNucleusXsc.hh.

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

Implements G4VComponentCrossSection.

Definition at line 100 of file G4ComponentGGHadronNucleusXsc.cc.

103 {
104  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
105  kinEnergy);
106  fTotalXsc = GetIsoCrossSection(aDP, Z, A);
107  delete aDP;
108 
109  return fTotalXsc;
110 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
double A(double temperature)
G4ThreeVector G4ParticleMomentum

Here is the call graph for this function:

G4bool G4ComponentGGHadronNucleusXsc::IsIsoApplicable ( const G4DynamicParticle aDP,
G4int  Z,
G4int  A,
const G4Element elm = 0,
const G4Material mat = 0 
)

Definition at line 236 of file G4ComponentGGHadronNucleusXsc.cc.

240 {
241  G4bool applicable = false;
242  // G4int baryonNumber = aDP->GetDefinition()->GetBaryonNumber();
243  G4double kineticEnergy = aDP->GetKineticEnergy();
244 
245  const G4ParticleDefinition* theParticle = aDP->GetDefinition();
246 
247  if (
248  Z >= 1 // >= H for kaons
249  &&
250  (
251  kineticEnergy >= fLowerLimit
252  &&
253  // Z > 1 && // >= He
254  (
255  theParticle == theAProton ||
256  theParticle == theGamma ||
257  theParticle == theSMinus ||
258  theParticle == theProton ||
259  theParticle == theNeutron ||
260  theParticle == thePiPlus ||
261  theParticle == thePiMinus
262  )
263  )
264  )
265  applicable = true;
266 
267  if (
268  Z >= 1 // >= H for kaons
269  &&
270  (
271  kineticEnergy >= 0.01*fLowerLimit
272  &&
273  (
274  theParticle == theKPlus ||
275  theParticle == theKMinus ||
276  theParticle == theK0L ||
277  theParticle == theK0S
278  )
279  )
280  )
281  applicable = true;
282 
283  return applicable;
284 }
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
bool G4bool
Definition: G4Types.hh:79
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4ComponentGGHadronNucleusXsc::SetEnergyLowerLimit ( G4double  E)
inline

Definition at line 161 of file G4ComponentGGHadronNucleusXsc.hh.

161 {fLowerLimit=E;};

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