Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4HadronElastic Class Reference

#include <G4HadronElastic.hh>

+ Inheritance diagram for G4HadronElastic:

Public Member Functions

 G4HadronElastic (const G4String &name="hElasticLHEP")
 
 ~G4HadronElastic () override
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) override
 
G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A) override
 
G4double GetSlopeCof (const G4int pdg)
 
void SetLowestEnergyLimit (G4double value)
 
G4double LowestEnergyLimit () const
 
G4double ComputeMomentumCMS (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
void ModelDescription (std::ostream &) const override
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void InitialiseModel ()
 
 G4HadronicInteraction (const G4HadronicInteraction &right)=delete
 
const G4HadronicInteractionoperator= (const G4HadronicInteraction &right)=delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 

Protected Attributes

G4double pLocalTmax
 
G4int secID
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 

Detailed Description

Definition at line 48 of file G4HadronElastic.hh.

Constructor & Destructor Documentation

◆ G4HadronElastic()

G4HadronElastic::G4HadronElastic ( const G4String name = "hElasticLHEP")
explicit

Definition at line 49 of file G4HadronElastic.cc.

50 : G4HadronicInteraction(name), secID(-1)
51{
52 SetMinEnergy( 0.0*GeV );
54 lowestEnergyLimit= 1.e-6*eV;
55 pLocalTmax = 0.0;
56 nwarn = 0;
57
58 theProton = G4Proton::Proton();
59 theNeutron = G4Neutron::Neutron();
60 theDeuteron = G4Deuteron::Deuteron();
61 theAlpha = G4Alpha::Alpha();
62
63 secID = G4PhysicsModelCatalog::GetModelID( "model_" + name );
64}
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
static G4HadronicParameters * Instance()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4int GetModelID(const G4int modelIndex)
static G4Proton * Proton()
Definition: G4Proton.cc:92

◆ ~G4HadronElastic()

G4HadronElastic::~G4HadronElastic ( )
override

Definition at line 66 of file G4HadronElastic.cc.

67{}

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4HadronElastic::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
)
overridevirtual

Reimplemented from G4HadronicInteraction.

Reimplemented in G4NeutrinoElectronNcModel, G4NeutronElectronElModel, G4LEHadronProtonElastic, G4LEnp, and G4LEpp.

Definition at line 81 of file G4HadronElastic.cc.

83{
85
86 const G4HadProjectile* aParticle = &aTrack;
87 G4double ekin = aParticle->GetKineticEnergy();
88
89 // no scattering below the limit
90 if(ekin <= lowestEnergyLimit) {
93 return &theParticleChange;
94 }
95
96 G4int A = targetNucleus.GetA_asInt();
97 G4int Z = targetNucleus.GetZ_asInt();
98
99 // Scattered particle referred to axis of incident particle
100 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
101 G4double m1 = theParticle->GetPDGMass();
102 G4double plab = std::sqrt(ekin*(ekin + 2.0*m1));
103
104 if (verboseLevel>1) {
105 G4cout << "G4HadronElastic: "
106 << aParticle->GetDefinition()->GetParticleName()
107 << " Plab(GeV/c)= " << plab/GeV
108 << " Ekin(MeV) = " << ekin/MeV
109 << " scattered off Z= " << Z
110 << " A= " << A
111 << G4endl;
112 }
113
115 G4double e1 = m1 + ekin;
116 G4LorentzVector lv(0.0,0.0,plab,e1+mass2);
117 G4ThreeVector bst = lv.boostVector();
118 G4double momentumCMS = plab*mass2/std::sqrt(m1*m1 + mass2*mass2 + 2.*mass2*e1);
119
120 pLocalTmax = 4.0*momentumCMS*momentumCMS;
121
122 // Sampling in CM system
123 G4double t = SampleInvariantT(theParticle, plab, Z, A);
124
125 if(t < 0.0 || t > pLocalTmax) {
126 // For the very rare cases where cos(theta) is greater than 1 or smaller than -1,
127 // print some debugging information via a "JustWarning" exception, and resample
128 // using the default algorithm
129#ifdef G4VERBOSE
130 if(nwarn < 2) {
132 ed << GetModelName() << " wrong sampling t= " << t << " tmax= " << pLocalTmax
133 << " for " << aParticle->GetDefinition()->GetParticleName()
134 << " ekin=" << ekin << " MeV"
135 << " off (Z,A)=(" << Z << "," << A << ") - will be resampled" << G4endl;
136 G4Exception( "G4HadronElastic::ApplyYourself", "hadEla001", JustWarning, ed);
137 ++nwarn;
138 }
139#endif
140 t = G4HadronElastic::SampleInvariantT(theParticle, plab, Z, A);
141 }
142
143 G4double phi = G4UniformRand()*CLHEP::twopi;
144 G4double cost = 1. - 2.0*t/pLocalTmax;
145
146 if (cost > 1.0) { cost = 1.0; }
147 else if(cost < -1.0) { cost = -1.0; }
148
149 G4double sint = std::sqrt((1.0-cost)*(1.0+cost));
150
151 if (verboseLevel>1) {
152 G4cout << " t= " << t << " tmax(GeV^2)= " << pLocalTmax/(GeV*GeV)
153 << " Pcms(GeV)= " << momentumCMS/GeV << " cos(t)=" << cost
154 << " sin(t)=" << sint << G4endl;
155 }
156 G4LorentzVector nlv1(momentumCMS*sint*std::cos(phi),
157 momentumCMS*sint*std::sin(phi),
158 momentumCMS*cost,
159 std::sqrt(momentumCMS*momentumCMS + m1*m1));
160
161 nlv1.boost(bst);
162
163 G4double eFinal = nlv1.e() - m1;
164 if (verboseLevel > 1) {
165 G4cout <<"G4HadronElastic: m= " << m1 << " Efin(MeV)= " << eFinal
166 << " 4-M Final: " << nlv1
167 << G4endl;
168 }
169
170 if(eFinal <= 0.0) {
173 } else {
174 theParticleChange.SetMomentumChange(nlv1.vect().unit());
176 }
177 lv -= nlv1;
178 G4double erec = std::max(lv.e() - mass2, 0.0);
179 if (verboseLevel > 1) {
180 G4cout << "Recoil: " <<" m= " << mass2 << " Erec(MeV)= " << erec
181 << " 4-mom: " << lv
182 << G4endl;
183 }
184
185 // the recoil is created if kinetic energy above the threshold
186 if(erec > GetRecoilEnergyThreshold()) {
187 G4ParticleDefinition * theDef = nullptr;
188 if(Z == 1 && A == 1) { theDef = theProton; }
189 else if (Z == 1 && A == 2) { theDef = theDeuteron; }
190 else if (Z == 1 && A == 3) { theDef = G4Triton::Triton(); }
191 else if (Z == 2 && A == 3) { theDef = G4He3::He3(); }
192 else if (Z == 2 && A == 4) { theDef = theAlpha; }
193 else {
194 theDef =
196 }
197 G4DynamicParticle * aSec = new G4DynamicParticle(theDef, lv.vect().unit(), erec);
199 } else {
201 }
202
203 return &theParticleChange;
204}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
void SetLocalEnergyDeposit(G4double aE)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double SampleInvariantT(const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A) override
G4double GetRecoilEnergyThreshold() const
const G4String & GetModelName() const
static G4He3 * He3()
Definition: G4He3.cc:93
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetA_asInt() const
Definition: G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
const G4String & GetParticleName() const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
static G4Triton * Triton()
Definition: G4Triton.cc:93

◆ ComputeMomentumCMS()

G4double G4HadronElastic::ComputeMomentumCMS ( const G4ParticleDefinition p,
G4double  plab,
G4int  Z,
G4int  A 
)
inline

Definition at line 102 of file G4HadronElastic.hh.

104{
105 G4double m1 = p->GetPDGMass();
106 G4double m12= m1*m1;
108 return plab*mass2/std::sqrt(m12 + mass2*mass2 + 2.*mass2*std::sqrt(m12 + plab*plab));
109}

◆ GetSlopeCof()

G4double G4HadronElastic::GetSlopeCof ( const G4int  pdg)

Definition at line 277 of file G4HadronElastic.cc.

278{
279 // The input parameter "pdg" should be the absolute value of the PDG code
280 // (i.e. the same value for a particle and its antiparticle).
281
282 G4double coeff = 1.0;
283
284 // heavy barions
285
286 static const G4double lBarCof1S = 0.88;
287 static const G4double lBarCof2S = 0.76;
288 static const G4double lBarCof3S = 0.64;
289 static const G4double lBarCof1C = 0.784378;
290 static const G4double lBarCofSC = 0.664378;
291 static const G4double lBarCof2SC = 0.544378;
292 static const G4double lBarCof1B = 0.740659;
293 static const G4double lBarCofSB = 0.620659;
294 static const G4double lBarCof2SB = 0.500659;
295
296 if( pdg == 3122 || pdg == 3222 || pdg == 3112 || pdg == 3212 )
297 {
298 coeff = lBarCof1S; // Lambda, Sigma+, Sigma-, Sigma0
299
300 } else if( pdg == 3322 || pdg == 3312 )
301 {
302 coeff = lBarCof2S; // Xi-, Xi0
303 }
304 else if( pdg == 3324)
305 {
306 coeff = lBarCof3S; // Omega
307 }
308 else if( pdg == 4122 || pdg == 4212 || pdg == 4222 || pdg == 4112 )
309 {
310 coeff = lBarCof1C; // LambdaC+, SigmaC+, SigmaC++, SigmaC0
311 }
312 else if( pdg == 4332 )
313 {
314 coeff = lBarCof2SC; // OmegaC
315 }
316 else if( pdg == 4232 || pdg == 4132 )
317 {
318 coeff = lBarCofSC; // XiC+, XiC0
319 }
320 else if( pdg == 5122 || pdg == 5222 || pdg == 5112 || pdg == 5212 )
321 {
322 coeff = lBarCof1B; // LambdaB, SigmaB+, SigmaB-, SigmaB0
323 }
324 else if( pdg == 5332 )
325 {
326 coeff = lBarCof2SB; // OmegaB-
327 }
328 else if( pdg == 5132 || pdg == 5232 ) // XiB-, XiB0
329 {
330 coeff = lBarCofSB;
331 }
332 // heavy mesons Kaons?
333 static const G4double lMesCof1S = 0.82; // Kp/piP kaons?
334 static const G4double llMesCof1C = 0.676568;
335 static const G4double llMesCof1B = 0.610989;
336 static const G4double llMesCof2C = 0.353135;
337 static const G4double llMesCof2B = 0.221978;
338 static const G4double llMesCofSC = 0.496568;
339 static const G4double llMesCofSB = 0.430989;
340 static const G4double llMesCofCB = 0.287557;
341 static const G4double llMesCofEtaP = 0.88;
342 static const G4double llMesCofEta = 0.76;
343
344 if( pdg == 321 || pdg == 311 || pdg == 310 )
345 {
346 coeff = lMesCof1S; //K+-0
347 }
348 else if( pdg == 511 || pdg == 521 )
349 {
350 coeff = llMesCof1B; // BMeson0, BMeson+
351 }
352 else if(pdg == 421 || pdg == 411 )
353 {
354 coeff = llMesCof1C; // DMeson+, DMeson0
355 }
356 else if( pdg == 531 )
357 {
358 coeff = llMesCofSB; // BSMeson0
359 }
360 else if( pdg == 541 )
361 {
362 coeff = llMesCofCB; // BCMeson+-
363 }
364 else if(pdg == 431 )
365 {
366 coeff = llMesCofSC; // DSMeson+-
367 }
368 else if(pdg == 441 || pdg == 443 )
369 {
370 coeff = llMesCof2C; // Etac, JPsi
371 }
372 else if(pdg == 553 )
373 {
374 coeff = llMesCof2B; // Upsilon
375 }
376 else if(pdg == 221 )
377 {
378 coeff = llMesCofEta; // Eta
379 }
380 else if(pdg == 331 )
381 {
382 coeff = llMesCofEtaP; // Eta'
383 }
384 return coeff;
385}

◆ LowestEnergyLimit()

G4double G4HadronElastic::LowestEnergyLimit ( ) const
inline

Definition at line 96 of file G4HadronElastic.hh.

97{
98 return lowestEnergyLimit;
99}

Referenced by G4NeutrinoElectronNcModel::ApplyYourself(), and G4NeutronElectronElModel::ApplyYourself().

◆ ModelDescription()

void G4HadronElastic::ModelDescription ( std::ostream &  outFile) const
overridevirtual

Reimplemented from G4HadronicInteraction.

Reimplemented in G4NeutrinoElectronNcModel, and G4NeutronElectronElModel.

Definition at line 70 of file G4HadronElastic.cc.

71{
72 outFile << "G4HadronElastic is the base class for all hadron-nucleus\n"
73 << "elastic scattering models except HP.\n"
74 << "By default it uses the Gheisha two-exponential momentum\n"
75 << "transfer parameterization. The model is fully relativistic\n"
76 << "as opposed to the original Gheisha model which was not.\n"
77 << "This model may be used for all long-lived hadrons at all\n"
78 << "incident energies but fit the data only for relativistic scattering.\n";
79}

◆ SampleInvariantT()

G4double G4HadronElastic::SampleInvariantT ( const G4ParticleDefinition p,
G4double  plab,
G4int  Z,
G4int  A 
)
overridevirtual

Reimplemented from G4HadronicInteraction.

Reimplemented in G4NuclNuclDiffuseElastic, G4LEHadronProtonElastic, G4LEnp, G4LEpp, G4LowEHadronElastic, and G4hhElastic.

Definition at line 208 of file G4HadronElastic.cc.

210{
211 const G4double plabLowLimit = 400.0*CLHEP::MeV;
212 const G4double GeV2 = GeV*GeV;
213 const G4double z07in13 = std::pow(0.7, 0.3333333333);
214 const G4double numLimit = 18.;
215
216 G4int pdg = std::abs(part->GetPDGEncoding());
217 G4double tmax = pLocalTmax/GeV2;
218
219 G4double aa, bb, cc, dd;
220 G4Pow* g4pow = G4Pow::GetInstance();
221 if (A <= 62) {
222 if (pdg == 211){ //Pions
223 if(mom >= plabLowLimit){ //High energy
224 bb = 14.5*g4pow->Z23(A);/*14.5*/
225 dd = 10.;
226 cc = 0.075*g4pow->Z13(A)/dd;//1.4
227 //aa = g4pow->powZ(A, 1.93)/bb;//1.63
228 aa = (A*A)/bb;//1.63
229 } else { //Low energy
230 bb = 29.*z07in13*z07in13*g4pow->Z23(A);
231 dd = 15.;
232 cc = 0.04*g4pow->Z13(A)/dd;//1.4
233 aa = g4pow->powZ(A, 1.63)/bb;//1.63
234 }
235 } else { //Other particles
236 bb = 14.5*g4pow->Z23(A);
237 dd = 20.;
238 aa = (A*A)/bb;//1.63
239 cc = 1.4*g4pow->Z13(A)/dd;
240 }
241 //===========================
242 } else { //(A>62)
243 if (pdg == 211) {
244 if(mom >= plabLowLimit){ //high
245 bb = 60.*z07in13*g4pow->Z13(A);//60
246 dd = 30.;
247 aa = 0.5*(A*A)/bb;//1.33
248 cc = 4.*g4pow->powZ(A,0.4)/dd;//1:0.4 --- 2: 0.4
249 } else { //low
250 bb = 120.*z07in13*g4pow->Z13(A);//60
251 dd = 30.;
252 aa = 2.*g4pow->powZ(A,1.33)/bb;
253 cc = 4.*g4pow->powZ(A,0.4)/dd;//1:0.4 --- 2: 0.4
254 }
255 } else {
256 bb = 60.*g4pow->Z13(A);
257 dd = 25.;
258 aa = g4pow->powZ(A,1.33)/bb;//1.33
259 cc = 0.2*g4pow->powZ(A,0.4)/dd;//1:0.4 --- 2: 0.4
260 }
261 }
262 G4double q1 = 1.0 - G4Exp(-std::min(bb*tmax, numLimit));
263 G4double q2 = 1.0 - G4Exp(-std::min(dd*tmax, numLimit));
264 G4double s1 = q1*aa;
265 G4double s2 = q2*cc;
266 if((s1 + s2)*G4UniformRand() < s2) {
267 q1 = q2;
268 bb = dd;
269 }
270 return -GeV2*G4Log(1.0 - G4UniformRand()*q1)/bb;
271}
const G4double GeV2
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4double G4Log(G4double x)
Definition: G4Log.hh:227
Definition: G4Pow.hh:49
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powZ(G4int Z, G4double y) const
Definition: G4Pow.hh:225
G4double Z13(G4int Z) const
Definition: G4Pow.hh:123
G4double Z23(G4int Z) const
Definition: G4Pow.hh:125

Referenced by ApplyYourself(), G4ChipsElasticModel::SampleInvariantT(), G4ElasticHadrNucleusHE::SampleInvariantT(), and G4LowEHadronElastic::SampleInvariantT().

◆ SetLowestEnergyLimit()

void G4HadronElastic::SetLowestEnergyLimit ( G4double  value)
inline

Definition at line 91 of file G4HadronElastic.hh.

92{
93 lowestEnergyLimit = value;
94}

Referenced by G4NeutrinoElectronNcModel::G4NeutrinoElectronNcModel(), and G4NeutronElectronElModel::G4NeutronElectronElModel().

Member Data Documentation

◆ pLocalTmax

◆ secID


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