Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ChargeExchangeXS Class Referencefinal

#include <G4ChargeExchangeXS.hh>

+ Inheritance diagram for G4ChargeExchangeXS:

Public Member Functions

 G4ChargeExchangeXS ()
 
 ~G4ChargeExchangeXS () override=default
 
G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *) final
 
G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *) final
 
void CrossSectionDescription (std::ostream &) const final
 
const G4ParticleDefinitionSampleSecondaryType (const G4ParticleDefinition *, const G4int Z, const G4int A)
 
void SetEnergyLimit (G4double val)
 
void SetCrossSectionFactor (G4double val)
 
G4double GetCrossSectionFactor () const
 
G4ChargeExchangeXSoperator= (const G4ChargeExchangeXS &right)=delete
 
 G4ChargeExchangeXS (const G4ChargeExchangeXS &)=delete
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double ComputeCrossSectionPerElement (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
virtual G4double ComputeIsoCrossSection (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
virtual const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
bool ForAllAtomsAndEnergies () const
 
void SetForAllAtomsAndEnergies (G4bool val)
 
const G4StringGetName () const
 
void SetName (const G4String &nam)
 
G4VCrossSectionDataSetoperator= (const G4VCrossSectionDataSet &right)=delete
 
 G4VCrossSectionDataSet (const G4VCrossSectionDataSet &)=delete
 

Additional Inherited Members

- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel {0}
 
G4String name
 

Detailed Description

Definition at line 65 of file G4ChargeExchangeXS.hh.

Constructor & Destructor Documentation

◆ G4ChargeExchangeXS() [1/2]

G4ChargeExchangeXS::G4ChargeExchangeXS ( )

Definition at line 66 of file G4ChargeExchangeXS.cc.

67{
68 if (verboseLevel > 1) {
69 G4cout << "G4ChargeExchangeXS::G4ChargeExchangeXS" << G4endl;
70 }
71 g4calc = G4Pow::GetInstance();
73 const G4String nam[5] = {"pi0", "eta", "eta_prime", "omega", "f2(1270)"};
74 for (G4int i=0; i<5; ++i) {
75 fPionSecPD[i] = table->FindParticle(nam[i]);
76 if (nullptr == fPionSecPD[i]) {
78 ed << "### meson " << nam[i] << " is not found out in the particle table";
79 G4Exception("G4ChargeExchangeXS::G4ChargeExchangeXS()","had044",
80 FatalException, ed,"");
81 }
82 }
83}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4ParticleTable * GetParticleTable()
static G4Pow * GetInstance()
Definition G4Pow.cc:41

Referenced by G4ChargeExchangeXS(), and operator=().

◆ ~G4ChargeExchangeXS()

G4ChargeExchangeXS::~G4ChargeExchangeXS ( )
overridedefault

◆ G4ChargeExchangeXS() [2/2]

G4ChargeExchangeXS::G4ChargeExchangeXS ( const G4ChargeExchangeXS & )
delete

Member Function Documentation

◆ CrossSectionDescription()

void G4ChargeExchangeXS::CrossSectionDescription ( std::ostream & outFile) const
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 86 of file G4ChargeExchangeXS.cc.

87{
88 outFile << "G4ChargeExchangeXS calculates charge exchange cross section for "
89 << "pi+, pi-, K+, K-, KL\n";
90}

◆ GetCrossSectionFactor()

G4double G4ChargeExchangeXS::GetCrossSectionFactor ( ) const
inline

Definition at line 89 of file G4ChargeExchangeXS.hh.

89{ return fFactor; };

◆ GetElementCrossSection()

G4double G4ChargeExchangeXS::GetElementCrossSection ( const G4DynamicParticle * aParticle,
G4int Z,
const G4Material * mat )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 99 of file G4ChargeExchangeXS.cc.

101{
102 G4double result = 0.0;
103 const G4double pE = aParticle->GetTotalEnergy();
104 if (pE <= fEnergyLimit) { return result; }
105
106 auto part = aParticle->GetDefinition();
107 G4int pdg = part->GetPDGEncoding();
108
109 // Get or calculate the proton mass, particle mass, and s(Lorentz invariant)
110 G4double tM = CLHEP::proton_mass_c2;
111 G4double pM = part->GetPDGMass();
112 G4double lorentz_s = tM*tM + 2*tM*pE + pM*pM;
113 if (lorentz_s <= (tM + pM)*(tM + pM)) { return result; }
114
115 const G4int Z = std::min(ZZ, ZMAXNUCLEARDATA);
116 const G4int A = G4lrint(aeff[Z]);
117
118 if (verboseLevel > 1) {
119 G4cout << "### G4ChargeExchangeXS: " << part->GetParticleName()
120 << " Z=" << Z << " A=" << A << " Etot(GeV)=" << pE/CLHEP::GeV
121 << " s(GeV^2)=" << lorentz_s/(CLHEP::GeV*CLHEP::GeV) << G4endl;
122 }
123
124 // For unit conversion
125 const G4double inv1e7 = 0.1/(CLHEP::GeV*CLHEP::GeV);
126 const G4double fact = 1e-30*CLHEP::cm2;
127 const G4double pfact = 0.1/CLHEP::GeV;
128 const G4double kfact = 56.3*fact;
129 const G4double csmax = 1e-16;
130
131 // The approximation of Glauber-Gribov formula -> extend it from interaction with
132 // proton to nuclei Z^(2/3). The factor g4calc->powA(A,-beta_prime_pi*G4Log(A))
133 // takes into account absorption of pi0 and eta
134
135 // pi- + p -> n + meson (0- pi0, 1- eta, 2- eta', 3- omega, 4- f2(1270))
136 if (pdg == -211) {
137 G4double z23 = g4calc->Z23(Z);
138 G4double x = lorentz_s*inv1e7;
139 G4double logX = G4Log(x);
140 G4double logA = g4calc->logZ(A);
141 G4double xf = g4calc->powZ(A, -beta_prime_pi*logA);
142 G4double sum = 0.0;
143 for (G4int i=0; i<5; ++i) {
144 G4double xg = std::max(1.0 + pG0[i] + pG1[i]*logX, 0.0);
145 G4double xc = std::max(pC0[i] + pC1[i]*logX, csmax);
146 G4double xs = z23*piA[i]*g4calc->powA(x, -pAP[i])*xf*xg/xc;
147 sum += xs;
148 fXSecPion[i] = sum;
149 }
150 result = sum*fact;
151 }
152
153 // pi+ + n -> p + meson (0- pi0, 1- eta, 2- eta', 3- omega, 4- f2(1270))
154 else if (pdg == 211) {
155 G4double n23 = g4calc->Z23(A - Z);
156 G4double x = lorentz_s*inv1e7;
157 G4double logX = G4Log(x);
158 G4double logA = g4calc->logZ(A);
159 G4double xf = g4calc->powZ(A, -beta_prime_pi*logA);
160
161 // hydrogen target case Z = A = 1
162 // the cross section is defined by fraction of deuteron and tritium
163 if (1 == Z) { n23 = ComputeDeuteronFraction(mat); }
164 G4double sum = 0.0;
165 for (G4int i=0; i<5; ++i) {
166 G4double xg = std::max(1.0 + pG0[i] + pG1[i]*logX, 0.0);
167 G4double xc = std::max(pC0[i] + pC1[i]*logX, csmax);
168 G4double xs = n23*piA[i]*g4calc->powA(x, -pAP[i])*xf*xg/xc;
169 sum += xs;
170 fXSecPion[i] = sum;
171 }
172 result = sum*fact;
173 }
174
175 // Kaon x-sections depend on the primary particles momentum
176 // K- + p -> Kbar + n
177 else if (pdg == -321) {
178 G4double p_momentum = std::sqrt(pE*pE - pM*pM)*pfact;
179 result = g4calc->Z23(Z)*g4calc->powA(p_momentum, -1.60)*kfact;
180 }
181
182 // K+ + n -> Kbar + p
183 else if (pdg == 321) {
184 G4double p_momentum = std::sqrt(pE*pE - pM*pM)*pfact;
185 G4double n23 = g4calc->Z23(A-Z);
186 // hydrogen target case Z = A = 1
187 // the cross section is defined by fraction of deuteron and tritium
188 if (1 == Z) { n23 = ComputeDeuteronFraction(mat); }
189 result = n23*g4calc->powA(p_momentum, -1.60)*kfact;
190 }
191
192 // KL
193 else if (pdg == 130) {
194 // Cross section of KL = 0.5*(Cross section of K+ + Cross section of K-)
195 const G4double p_momentum = std::sqrt(pE*pE - pM*pM)*pfact;
196 result = 0.5*(g4calc->Z23(Z) + g4calc->Z23(A-Z))*
197 g4calc->powA(p_momentum, -1.60)*kfact;
198 }
199 result *= fFactor;
200 if (verboseLevel > 1) {
201 G4cout << " Done for " << part->GetParticleName() << " Etot(GeV)=" << pE/CLHEP::GeV
202 << " res(mb)=" << result/CLHEP::millibarn << G4endl;
203 }
204 return result;
205}
G4double G4Log(G4double x)
Definition G4Log.hh:227
double G4double
Definition G4Types.hh:83
const G4double A[17]
G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
int G4lrint(double ad)
Definition templates.hh:134

◆ IsElementApplicable()

G4bool G4ChargeExchangeXS::IsElementApplicable ( const G4DynamicParticle * ,
G4int Z,
const G4Material *  )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 92 of file G4ChargeExchangeXS.cc.

94{
95 return true;
96}

◆ operator=()

G4ChargeExchangeXS & G4ChargeExchangeXS::operator= ( const G4ChargeExchangeXS & right)
delete

◆ SampleSecondaryType()

const G4ParticleDefinition * G4ChargeExchangeXS::SampleSecondaryType ( const G4ParticleDefinition * part,
const G4int Z,
const G4int A )

Definition at line 208 of file G4ChargeExchangeXS.cc.

210{
211 const G4ParticleDefinition* pd = nullptr;
212 G4int pdg = std::abs(part->GetPDGEncoding());
213
214 // pi- + p / pi+ + n
215 if (pdg == 211) {
216 pd = fPionSecPD[0];
217 G4double x = fXSecPion[4]*G4UniformRand();
218 for (G4int i=0; i<5; ++i) {
219 if (x <= fXSecPion[i]) {
220 pd = fPionSecPD[i];
221 break;
222 }
223 }
224 }
225
226 // K- + p / K+ + n
227 // Equal opportunity of producing k-short and k-long
228 else if (pdg == 321) {
229 if (G4UniformRand() >= 0.5) {
231 }
232 else {
234 }
235 }
236
237 // KL + nucleus
238 else if (pdg == 130) {
239 G4double prob = (G4double)Z/(G4double)A;
240 if (G4UniformRand() >= prob) {
242 }
243 else {
245 }
246 }
247
248 return pd;
249}
#define G4UniformRand()
Definition Randomize.hh:52
static G4KaonMinus * KaonMinus()
static G4KaonPlus * KaonPlus()
static G4KaonZeroLong * KaonZeroLong()
static G4KaonZeroShort * KaonZeroShort()

◆ SetCrossSectionFactor()

void G4ChargeExchangeXS::SetCrossSectionFactor ( G4double val)
inline

Definition at line 87 of file G4ChargeExchangeXS.hh.

87{ fFactor = val; };

◆ SetEnergyLimit()

void G4ChargeExchangeXS::SetEnergyLimit ( G4double val)
inline

Definition at line 85 of file G4ChargeExchangeXS.hh.

85{ fEnergyLimit = val; };

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