Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
G4AntiNuclElastic Class Reference

#include <G4AntiNuclElastic.hh>

Inheritance diagram for G4AntiNuclElastic:
G4HadronElastic G4HadronicInteraction

Public Member Functions

void ActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Material *aMaterial)
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) override
 
G4double BesselJone (G4double z)
 
G4double BesselJzero (G4double z)
 
G4double BesselOneByArg (G4double z)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4double CalculateAm (G4double momentum, G4double n, G4double Z)
 
G4double CalculateParticleBeta (const G4ParticleDefinition *particle, G4double momentum)
 
G4double CalculateZommerfeld (G4double beta, G4double Z1, G4double Z2)
 
G4double ComputeMomentumCMS (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
G4double DampFactor (G4double z)
 
void DeActivateFor (const G4Element *anElement)
 
void DeActivateFor (const G4Material *aMaterial)
 
 G4AntiNuclElastic ()
 
G4ComponentAntiNuclNuclearXSGetComponentCrossSection ()
 
G4double GetcosTeta1 (G4double plab, G4int A)
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
const G4StringGetModelName () const
 
G4double GetRecoilEnergyThreshold () const
 
G4double GetSlopeCof (const G4int pdg)
 
G4int GetVerboseLevel () const
 
virtual void InitialiseModel ()
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4bool IsBlocked (const G4Element *anElement) const
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4double LowestEnergyLimit () const
 
void ModelDescription (std::ostream &) const override
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A) override
 
G4double SampleThetaCMS (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
G4double SampleThetaLab (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
void SetLowestEnergyLimit (G4double value)
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
void SetRecoilEnergyThreshold (G4double val)
 
void SetVerboseLevel (G4int value)
 
 ~G4AntiNuclElastic () override
 

Protected Member Functions

void Block ()
 
G4bool IsBlocked () const
 
void SetModelName (const G4String &nam)
 

Protected Attributes

G4bool isBlocked
 
G4double pLocalTmax
 
G4int secID
 
G4double theMaxEnergy
 
G4double theMinEnergy
 
G4HadFinalState theParticleChange
 
G4int verboseLevel
 

Private Member Functions

 G4AntiNuclElastic (const G4AntiNuclElastic &)
 
G4AntiNuclElasticoperator= (const G4AntiNuclElastic &right)
 

Private Attributes

G4ComponentAntiNuclNuclearXScs
 
std::pair< G4double, G4doubleepCheckLevels
 
G4double fAm
 
G4double fBeta
 
G4ThreeVector fbst
 
G4double fceff
 
const G4ParticleDefinitionfParticle
 
G4double fptot
 
G4double fRa
 
G4double fRef
 
G4double fTetaCMS
 
G4double fThetaLab
 
G4double fTmax
 
G4double fWaveVector
 
G4double fZommerfeld
 
G4double lowestEnergyLimit
 
G4int nwarn
 
G4double recoilEnergyThreshold
 
G4HadronicInteractionRegistryregistry
 
G4ParticleDefinitiontheAAlpha
 
G4ParticleDefinitiontheADeuteron
 
G4ParticleDefinitiontheAHe3
 
G4ParticleDefinitiontheAlpha
 
G4ParticleDefinitiontheANeutron
 
G4ParticleDefinitiontheAProton
 
G4ParticleDefinitiontheATriton
 
std::vector< const G4Material * > theBlockedList
 
std::vector< const G4Element * > theBlockedListElements
 
G4ParticleDefinitiontheDeuteron
 
std::vector< std::pair< G4double, const G4Material * > > theMaxEnergyList
 
std::vector< std::pair< G4double, const G4Element * > > theMaxEnergyListElements
 
std::vector< std::pair< G4double, const G4Material * > > theMinEnergyList
 
std::vector< std::pair< G4double, const G4Element * > > theMinEnergyListElements
 
G4String theModelName
 
G4ParticleDefinitiontheNeutron
 
G4ParticleDefinitiontheProton
 

Detailed Description

Definition at line 46 of file G4AntiNuclElastic.hh.

Constructor & Destructor Documentation

◆ G4AntiNuclElastic() [1/2]

G4AntiNuclElastic::G4AntiNuclElastic ( )
explicit

Definition at line 56 of file G4AntiNuclElastic.cc.

57 : G4HadronElastic("AntiAElastic")
58{
59 //V.Ivanchenko commented out
60 //SetMinEnergy( 0.1*GeV );
61 //SetMaxEnergy( 10.*TeV );
62
69
74
76 cs = static_cast<G4ComponentAntiNuclNuclearXS*>(reg->GetComponentCrossSection("AntiAGlauber"));
77 if(!cs) { cs = new G4ComponentAntiNuclNuclearXS(); }
78
79 fParticle = 0;
80 fWaveVector = 0.;
81 fBeta = 0.;
82 fZommerfeld = 0.;
83 fAm = 0.;
84 fTetaCMS = 0.;
85 fRa = 0.;
86 fRef = 0.;
87 fceff = 0.;
88 fptot = 0.;
89 fTmax = 0.;
90 fThetaLab = 0.;
91}
static const G4double reg
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
static G4AntiAlpha * AntiAlpha()
Definition: G4AntiAlpha.cc:88
static G4AntiDeuteron * AntiDeuteron()
static G4AntiHe3 * AntiHe3()
Definition: G4AntiHe3.cc:93
static G4AntiNeutron * AntiNeutron()
G4ParticleDefinition * theAProton
G4ParticleDefinition * theAAlpha
G4ParticleDefinition * theNeutron
G4ComponentAntiNuclNuclearXS * cs
G4ParticleDefinition * theProton
G4ParticleDefinition * theADeuteron
G4ParticleDefinition * theANeutron
G4ParticleDefinition * theATriton
G4ParticleDefinition * theAHe3
G4ParticleDefinition * theDeuteron
G4ParticleDefinition * theAlpha
const G4ParticleDefinition * fParticle
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:92
static G4AntiTriton * AntiTriton()
Definition: G4AntiTriton.cc:93
static G4CrossSectionDataSetRegistry * Instance()
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
G4HadronElastic(const G4String &name="hElasticLHEP")
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4Proton * Proton()
Definition: G4Proton.cc:92

References G4Alpha::Alpha(), G4AntiAlpha::AntiAlpha(), G4AntiDeuteron::AntiDeuteron(), G4AntiHe3::AntiHe3(), G4AntiNeutron::AntiNeutron(), G4AntiProton::AntiProton(), G4AntiTriton::AntiTriton(), cs, G4Deuteron::Deuteron(), fAm, fBeta, fceff, fParticle, fptot, fRa, fRef, fTetaCMS, fThetaLab, fTmax, fWaveVector, fZommerfeld, G4CrossSectionDataSetRegistry::Instance(), G4Neutron::Neutron(), G4Proton::Proton(), reg, theAAlpha, theADeuteron, theAHe3, theAlpha, theANeutron, theAProton, theATriton, theDeuteron, theNeutron, and theProton.

◆ ~G4AntiNuclElastic()

G4AntiNuclElastic::~G4AntiNuclElastic ( )
override

Definition at line 94 of file G4AntiNuclElastic.cc.

95{}

◆ G4AntiNuclElastic() [2/2]

G4AntiNuclElastic::G4AntiNuclElastic ( const G4AntiNuclElastic )
private

Member Function Documentation

◆ ActivateFor() [1/2]

void G4HadronicInteraction::ActivateFor ( const G4Element anElement)
inlineinherited

◆ ActivateFor() [2/2]

void G4HadronicInteraction::ActivateFor ( const G4Material aMaterial)
inlineinherited

◆ ApplyYourself()

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

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
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}
static const G4double e1[44]
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
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
G4ParticleDefinition * theAlpha
G4ParticleDefinition * theProton
G4double lowestEnergyLimit
G4ParticleDefinition * theDeuteron
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
static constexpr double twopi
Definition: SystemOfUnits.h:56
T max(const T t1, const T t2)
brief Return the largest of the two arguments

References A, G4HadFinalState::AddSecondary(), CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::boostVector(), G4HadFinalState::Clear(), CLHEP::HepLorentzVector::e(), e1, G4cout, G4endl, G4Exception(), G4UniformRand, G4Nucleus::GetA_asInt(), G4HadProjectile::GetDefinition(), G4IonTable::GetIon(), G4ParticleTable::GetIonTable(), G4HadProjectile::GetKineticEnergy(), G4HadronicInteraction::GetModelName(), G4NucleiProperties::GetNuclearMass(), G4ParticleDefinition::GetParticleName(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetPDGMass(), G4HadronicInteraction::GetRecoilEnergyThreshold(), G4Nucleus::GetZ_asInt(), GeV, G4He3::He3(), JustWarning, G4HadronElastic::lowestEnergyLimit, G4INCL::Math::max(), MeV, G4HadronElastic::nwarn, G4HadronElastic::pLocalTmax, G4HadronElastic::SampleInvariantT(), G4HadronElastic::secID, G4HadFinalState::SetEnergyChange(), G4HadFinalState::SetLocalEnergyDeposit(), G4HadFinalState::SetMomentumChange(), G4HadronElastic::theAlpha, G4HadronElastic::theDeuteron, G4HadronicInteraction::theParticleChange, G4HadronElastic::theProton, G4Triton::Triton(), CLHEP::twopi, CLHEP::Hep3Vector::unit(), CLHEP::HepLorentzVector::vect(), G4HadronicInteraction::verboseLevel, and Z.

◆ BesselJone()

G4double G4AntiNuclElastic::BesselJone ( G4double  z)

Definition at line 570 of file G4AntiNuclElastic.cc.

571{
572 G4double modvalue, value2, fact1, fact2, arg, shift, bessel;
573
574 modvalue = std::fabs(value);
575
576 if ( modvalue < 8.0 )
577 {
578 value2 = value*value;
579 fact1 = value*(72362614232.0 + value2*(-7895059235.0
580 + value2*( 242396853.1
581 + value2*(-2972611.439
582 + value2*( 15704.48260
583 + value2*(-30.16036606 ) ) ) ) ) );
584
585 fact2 = 144725228442.0 + value2*(2300535178.0
586 + value2*(18583304.74
587 + value2*(99447.43394
588 + value2*(376.9991397
589 + value2*1.0 ) ) ) );
590 bessel = fact1/fact2;
591 }
592 else
593 {
594 arg = 8.0/modvalue;
595 value2 = arg*arg;
596
597 shift = modvalue - 2.356194491;
598
599 fact1 = 1.0 + value2*( 0.183105e-2
600 + value2*(-0.3516396496e-4
601 + value2*(0.2457520174e-5
602 + value2*(-0.240337019e-6 ) ) ) );
603
604 fact2 = 0.04687499995 + value2*(-0.2002690873e-3
605 + value2*( 0.8449199096e-5
606 + value2*(-0.88228987e-6
607 + value2*0.105787412e-6 ) ) );
608
609 bessel = std::sqrt( 0.636619772/modvalue)*(std::cos(shift)*fact1 - arg*std::sin(shift)*fact2);
610 if (value < 0.0) bessel = -bessel;
611 }
612 return bessel;
613}

Referenced by BesselOneByArg().

◆ BesselJzero()

G4double G4AntiNuclElastic::BesselJzero ( G4double  z)

Definition at line 519 of file G4AntiNuclElastic.cc.

520{
521 G4double modvalue, value2, fact1, fact2, arg, shift, bessel;
522
523 modvalue = std::fabs(value);
524
525 if ( value < 8.0 && value > -8.0 )
526 {
527 value2 = value*value;
528
529 fact1 = 57568490574.0 + value2*(-13362590354.0
530 + value2*( 651619640.7
531 + value2*(-11214424.18
532 + value2*( 77392.33017
533 + value2*(-184.9052456 ) ) ) ) );
534
535 fact2 = 57568490411.0 + value2*( 1029532985.0
536 + value2*( 9494680.718
537 + value2*(59272.64853
538 + value2*(267.8532712
539 + value2*1.0 ) ) ) );
540
541 bessel = fact1/fact2;
542 }
543 else
544 {
545 arg = 8.0/modvalue;
546
547 value2 = arg*arg;
548
549 shift = modvalue-0.785398164;
550
551 fact1 = 1.0 + value2*(-0.1098628627e-2
552 + value2*(0.2734510407e-4
553 + value2*(-0.2073370639e-5
554 + value2*0.2093887211e-6 ) ) );
555 fact2 = -0.1562499995e-1 + value2*(0.1430488765e-3
556 + value2*(-0.6911147651e-5
557 + value2*(0.7621095161e-6
558 - value2*0.934945152e-7 ) ) );
559
560 bessel = std::sqrt(0.636619772/modvalue)*(std::cos(shift)*fact1 - arg*std::sin(shift)*fact2);
561 }
562 return bessel;
563}

Referenced by SampleInvariantT().

◆ BesselOneByArg()

G4double G4AntiNuclElastic::BesselOneByArg ( G4double  z)

Definition at line 617 of file G4AntiNuclElastic.cc.

618{
619 G4double x2, result;
620
621 if( std::fabs(x) < 0.01 )
622 {
623 x *= 0.5;
624 x2 = x*x;
625 result = (2.- x2 + x2*x2/6.)/4.;
626 }
627 else
628 {
629 result = BesselJone(x)/x;
630 }
631 return result;
632}
G4double BesselJone(G4double z)

References BesselJone().

Referenced by SampleInvariantT().

◆ Block()

void G4HadronicInteraction::Block ( )
inlineprotectedinherited

◆ BuildPhysicsTable()

void G4HadronicInteraction::BuildPhysicsTable ( const G4ParticleDefinition )
virtualinherited

◆ CalculateAm()

G4double G4AntiNuclElastic::CalculateAm ( G4double  momentum,
G4double  n,
G4double  Z 
)

Definition at line 503 of file G4AntiNuclElastic.cc.

504{
505 G4double k = momentum/hbarc;
506 G4double ch = 1.13 + 3.76*n*n;
508 G4double zn2 = zn*zn;
509 fAm = ch/zn2;
510
511 return fAm;
512}
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double A13(G4double A) const
Definition: G4Pow.cc:120
float Bohr_radius
Definition: hepunit.py:289
float hbarc
Definition: hepunit.py:264

References G4Pow::A13(), source.hepunit::Bohr_radius, fAm, G4Pow::GetInstance(), source.hepunit::hbarc, CLHEP::detail::n, and Z.

Referenced by SampleInvariantT().

◆ CalculateParticleBeta()

G4double G4AntiNuclElastic::CalculateParticleBeta ( const G4ParticleDefinition particle,
G4double  momentum 
)

Definition at line 480 of file G4AntiNuclElastic.cc.

482{
483 G4double mass = particle->GetPDGMass();
484 G4double a = momentum/mass;
485 fBeta = a/std::sqrt(1+a*a);
486
487 return fBeta;
488}

References fBeta, and G4ParticleDefinition::GetPDGMass().

Referenced by SampleInvariantT().

◆ CalculateZommerfeld()

G4double G4AntiNuclElastic::CalculateZommerfeld ( G4double  beta,
G4double  Z1,
G4double  Z2 
)

◆ ComputeMomentumCMS()

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

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}

References A, G4NucleiProperties::GetNuclearMass(), G4ParticleDefinition::GetPDGMass(), and Z.

◆ DampFactor()

G4double G4AntiNuclElastic::DampFactor ( G4double  z)

Definition at line 460 of file G4AntiNuclElastic.cc.

461{
462 G4double df;
463 G4double f3 = 6.; // first factorials
464
465 if( std::fabs(x) < 0.01 )
466 {
467 df=1./(1.+x*x/f3);
468 }
469 else
470 {
471 df = x/std::sinh(x);
472 }
473 return df;
474}

Referenced by SampleInvariantT().

◆ DeActivateFor() [1/2]

void G4HadronicInteraction::DeActivateFor ( const G4Element anElement)
inherited

Definition at line 186 of file G4HadronicInteraction.cc.

187{
188 Block();
189 theBlockedListElements.push_back(anElement);
190}
std::vector< const G4Element * > theBlockedListElements

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theBlockedListElements.

◆ DeActivateFor() [2/2]

void G4HadronicInteraction::DeActivateFor ( const G4Material aMaterial)
inherited

Definition at line 180 of file G4HadronicInteraction.cc.

181{
182 Block();
183 theBlockedList.push_back(aMaterial);
184}
std::vector< const G4Material * > theBlockedList

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theBlockedList.

Referenced by G4HadronHElasticPhysics::ConstructProcess().

◆ GetComponentCrossSection()

G4ComponentAntiNuclNuclearXS * G4AntiNuclElastic::GetComponentCrossSection ( )
inline

Definition at line 121 of file G4AntiNuclElastic.hh.

122{
123 return cs;
124}

References cs.

Referenced by LBE::ConstructHad().

◆ GetcosTeta1()

G4double G4AntiNuclElastic::GetcosTeta1 ( G4double  plab,
G4int  A 
)

Definition at line 636 of file G4AntiNuclElastic.cc.

637{
638
639// G4double p0 =G4LossTableManager::Instance()->FactorForAngleLimit()*CLHEP::hbarc/CLHEP::fermi;
640 G4double p0 = 1.*hbarc/fermi;
641//G4double cteta1 = 1.0 - p0*p0/2.0 * pow(A,2./3.)/(plab*plab);
642 G4double cteta1 = 1.0 - p0*p0/2.0 * G4Pow::GetInstance()->Z23(A)/(plab*plab);
644 if(cteta1 < -1.) cteta1 = -1.0;
645 return cteta1;
646}
static constexpr double fermi
Definition: G4SIunits.hh:83
G4double Z23(G4int Z) const
Definition: G4Pow.hh:125

References A, fermi, G4Pow::GetInstance(), source.hepunit::hbarc, and G4Pow::Z23().

Referenced by SampleInvariantT().

◆ GetEnergyMomentumCheckLevels()

std::pair< G4double, G4double > G4HadronicInteraction::GetEnergyMomentumCheckLevels ( ) const
virtualinherited

◆ GetFatalEnergyCheckLevels()

const std::pair< G4double, G4double > G4HadronicInteraction::GetFatalEnergyCheckLevels ( ) const
virtualinherited

Reimplemented in G4FissLib, G4LFission, G4LENDFission, G4ParticleHPCapture, G4ParticleHPElastic, G4ParticleHPFission, G4ParticleHPInelastic, and G4ParticleHPThermalScattering.

Definition at line 210 of file G4HadronicInteraction.cc.

211{
212 // default level of Check
213 return std::pair<G4double, G4double>(2.*perCent, 1. * GeV);
214}
static constexpr double perCent
Definition: G4SIunits.hh:325

References GeV, and perCent.

Referenced by G4HadronicProcess::CheckResult().

◆ GetMaxEnergy() [1/2]

G4double G4HadronicInteraction::GetMaxEnergy ( ) const
inlineinherited

◆ GetMaxEnergy() [2/2]

G4double G4HadronicInteraction::GetMaxEnergy ( const G4Material aMaterial,
const G4Element anElement 
) const
inherited

Definition at line 131 of file G4HadronicInteraction.cc.

133{
134 if(!IsBlocked()) { return theMaxEnergy; }
135 if( IsBlocked(aMaterial) || IsBlocked(anElement) ) { return 0.0; }
136 if(!theMaxEnergyListElements.empty()) {
137 for(auto const& elmlist : theMaxEnergyListElements) {
138 if( anElement == elmlist.second )
139 { return elmlist.first; }
140 }
141 }
142 if(!theMaxEnergyList.empty()) {
143 for(auto const& matlist : theMaxEnergyList) {
144 if( aMaterial == matlist.second )
145 { return matlist.first; }
146 }
147 }
148 return theMaxEnergy;
149}
std::vector< std::pair< G4double, const G4Material * > > theMaxEnergyList
std::vector< std::pair< G4double, const G4Element * > > theMaxEnergyListElements

References G4HadronicInteraction::IsBlocked(), G4HadronicInteraction::theMaxEnergy, G4HadronicInteraction::theMaxEnergyList, and G4HadronicInteraction::theMaxEnergyListElements.

◆ GetMinEnergy() [1/2]

G4double G4HadronicInteraction::GetMinEnergy ( ) const
inlineinherited

◆ GetMinEnergy() [2/2]

G4double G4HadronicInteraction::GetMinEnergy ( const G4Material aMaterial,
const G4Element anElement 
) const
inherited

Definition at line 81 of file G4HadronicInteraction.cc.

83{
84 if(!IsBlocked()) { return theMinEnergy; }
85 if( IsBlocked(aMaterial) || IsBlocked(anElement) ) { return DBL_MAX; }
86 if(!theMinEnergyListElements.empty()) {
87 for(auto const& elmlist : theMinEnergyListElements) {
88 if( anElement == elmlist.second )
89 { return elmlist.first; }
90 }
91 }
92 if(!theMinEnergyList.empty()) {
93 for(auto const & matlist : theMinEnergyList) {
94 if( aMaterial == matlist.second )
95 { return matlist.first; }
96 }
97 }
98 return theMinEnergy;
99}
std::vector< std::pair< G4double, const G4Element * > > theMinEnergyListElements
std::vector< std::pair< G4double, const G4Material * > > theMinEnergyList
#define DBL_MAX
Definition: templates.hh:62

References DBL_MAX, G4HadronicInteraction::IsBlocked(), G4HadronicInteraction::theMinEnergy, G4HadronicInteraction::theMinEnergyList, and G4HadronicInteraction::theMinEnergyListElements.

◆ GetModelName()

const G4String & G4HadronicInteraction::GetModelName ( ) const
inlineinherited

Definition at line 115 of file G4HadronicInteraction.hh.

116 { return theModelName; }

References G4HadronicInteraction::theModelName.

Referenced by G4MuMinusCapturePrecompound::ApplyYourself(), G4HadronElastic::ApplyYourself(), G4INCLXXInterface::ApplyYourself(), G4TheoFSGenerator::ApplyYourself(), G4HadronStoppingProcess::AtRestDoIt(), G4VHadronPhysics::BuildModel(), G4HadronicProcess::CheckEnergyMomentumConservation(), G4HadronicProcess::CheckResult(), G4ChargeExchangePhysics::ConstructProcess(), G4MuonicAtomDecay::DecayIt(), G4LENDModel::DumpLENDTargetInfo(), G4AblaInterface::G4AblaInterface(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4EMDissociation::G4EMDissociation(), G4ExcitedStringDecay::G4ExcitedStringDecay(), G4LEHadronProtonElastic::G4LEHadronProtonElastic(), G4LENDModel::G4LENDModel(), G4LENDorBERTModel::G4LENDorBERTModel(), G4LEnp::G4LEnp(), G4LEpp::G4LEpp(), G4LFission::G4LFission(), G4LowEGammaNuclearModel::G4LowEGammaNuclearModel(), G4LowEIonFragmentation::G4LowEIonFragmentation(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4NeutrinoElectronCcModel::G4NeutrinoElectronCcModel(), G4NeutrinoNucleusModel::G4NeutrinoNucleusModel(), G4WilsonAbrasionModel::G4WilsonAbrasionModel(), G4INCLXXInterface::GetDeExcitationModelName(), G4EnergyRangeManager::GetHadronicInteraction(), G4VHighEnergyGenerator::GetProjectileNucleus(), G4NeutronRadCapture::InitialiseModel(), G4BinaryCascade::ModelDescription(), G4LMsdGenerator::ModelDescription(), G4VPartonStringModel::ModelDescription(), G4TheoFSGenerator::ModelDescription(), G4VHadronPhysics::NewModel(), G4NeutrinoElectronProcess::PostStepDoIt(), G4HadronicProcess::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4HadronicProcessStore::PrintModelHtml(), G4BinaryCascade::PropagateModelDescription(), G4HadronicProcessStore::RegisterInteraction(), and G4LENDModel::returnUnchanged().

◆ GetRecoilEnergyThreshold()

G4double G4HadronicInteraction::GetRecoilEnergyThreshold ( ) const
inlineinherited

◆ GetSlopeCof()

G4double G4HadronElastic::GetSlopeCof ( const G4int  pdg)
inherited

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}

◆ GetVerboseLevel()

G4int G4HadronicInteraction::GetVerboseLevel ( ) const
inlineinherited

Definition at line 109 of file G4HadronicInteraction.hh.

110 { return verboseLevel; }

References G4HadronicInteraction::verboseLevel.

◆ InitialiseModel()

void G4HadronicInteraction::InitialiseModel ( )
virtualinherited

◆ IsApplicable()

G4bool G4HadronicInteraction::IsApplicable ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
)
virtualinherited

◆ IsBlocked() [1/3]

G4bool G4HadronicInteraction::IsBlocked ( ) const
inlineprotectedinherited

◆ IsBlocked() [2/3]

G4bool G4HadronicInteraction::IsBlocked ( const G4Element anElement) const
inherited

Definition at line 202 of file G4HadronicInteraction.cc.

203{
204 for (auto const& elm : theBlockedListElements) {
205 if (anElement == elm) return true;
206 }
207 return false;
208}

References G4HadronicInteraction::theBlockedListElements.

◆ IsBlocked() [3/3]

G4bool G4HadronicInteraction::IsBlocked ( const G4Material aMaterial) const
inherited

Definition at line 193 of file G4HadronicInteraction.cc.

194{
195 for (auto const& mat : theBlockedList) {
196 if (aMaterial == mat) return true;
197 }
198 return false;
199}

References G4HadronicInteraction::theBlockedList.

◆ LowestEnergyLimit()

G4double G4HadronElastic::LowestEnergyLimit ( ) const
inlineinherited

◆ ModelDescription()

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

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}

◆ operator!=()

G4bool G4HadronicInteraction::operator!= ( const G4HadronicInteraction right) const
deleteinherited

◆ operator=()

G4AntiNuclElastic & G4AntiNuclElastic::operator= ( const G4AntiNuclElastic right)
private

◆ operator==()

G4bool G4HadronicInteraction::operator== ( const G4HadronicInteraction right) const
deleteinherited

◆ SampleInvariantT()

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

Reimplemented from G4HadronicInteraction.

Definition at line 99 of file G4AntiNuclElastic.cc.

101{
102 G4double T;
103 G4double Mproj = particle->GetPDGMass();
104 G4LorentzVector Pproj(0.,0.,Plab,std::sqrt(Plab*Plab+Mproj*Mproj));
105 G4double ctet1 = GetcosTeta1(Plab, A);
106
107 G4double energy=Pproj.e()-Mproj;
108
109 const G4ParticleDefinition* theParticle = particle;
110
111 G4ParticleDefinition * theTargetDef = 0;
112
113 if (Z == 1 && A == 1) theTargetDef = theProton;
114 else if (Z == 1 && A == 2) theTargetDef = theDeuteron;
115 else if (Z == 1 && A == 3) theTargetDef = G4Triton::Triton();
116 else if (Z == 2 && A == 3) theTargetDef = G4He3::He3();
117 else if (Z == 2 && A == 4) theTargetDef = theAlpha;
118
119
121
122 //transform to CMS
123
124 G4LorentzVector lv(0.0,0.0,0.0,TargMass);
125 lv += Pproj;
126 G4double S = lv.mag2()/(GeV*GeV);
127
128 G4ThreeVector bst = lv.boostVector();
129 Pproj.boost(-bst);
130
131 G4ThreeVector p1 = Pproj.vect();
132 G4double ptot = p1.mag();
133
134 fbst = bst;
135 fptot= ptot;
136 fTmax = 4.0*ptot*ptot; // In (MeV/c)^2
137
138 if(Plab < (std::abs(particle->GetBaryonNumber())*100)*MeV)
139 {return fTmax*G4UniformRand();}
140
141 // Calculation of NN collision properties
142 G4double PlabPerN = Plab/std::abs(theParticle->GetBaryonNumber());
143 G4double NucleonMass = 0.5*( theProton->GetPDGMass() + theNeutron->GetPDGMass() );
144 G4double PrNucleonMass(0.); // Projectile average nucleon mass
145 if( std::abs(theParticle->GetBaryonNumber()) == 1 ) { PrNucleonMass = theParticle->GetPDGMass(); }
146 else { PrNucleonMass = NucleonMass; }
147 G4double energyPerN = std::sqrt( sqr(PlabPerN) + sqr(PrNucleonMass));
148 energyPerN -= PrNucleonMass;
149 //---
150
151 G4double Z1 = particle->GetPDGCharge();
152 G4double Z2 = Z;
153
154 G4double beta = CalculateParticleBeta(particle, ptot);
156 G4double Am = CalculateAm( ptot, n, Z2 );
157 fWaveVector = ptot; // /hbarc;
158
159 G4LorentzVector Fproj(0.,0.,0.,0.);
160 const G4double mevToBarn = 0.38938e+6;
161 G4double XsCoulomb = mevToBarn*sqr(n/fWaveVector)*pi*(1+ctet1)/(1.+Am)/(1.+2.*Am-ctet1);
162
163 G4double XsElastHadronic =cs->GetElasticElementCrossSection(particle, energy, Z, (G4double)A);
164 G4double XsTotalHadronic =cs->GetTotalElementCrossSection(particle, energy, Z, (G4double)A);
165
166 XsElastHadronic/=millibarn; XsTotalHadronic/=millibarn;
167
168 G4double CoulombProb = XsCoulomb/(XsCoulomb+XsElastHadronic);
169
170 if(G4UniformRand() < CoulombProb)
171 { // Simulation of Coulomb scattering
172
173 G4double phi = twopi * G4UniformRand();
174 G4double Ksi = G4UniformRand();
175
176 G4double par1 = 2.*(1.+Am)/(1.+ctet1);
177
178 // ////sample ThetaCMS in Coulomb part
179
180 G4double cosThetaCMS = (par1*ctet1- Ksi*(1.+2.*Am))/(par1-Ksi);
181
182 G4double PtZ=ptot*cosThetaCMS;
183 Fproj.setPz(PtZ);
184 G4double PtProjCMS = ptot*std::sqrt(1.0 - cosThetaCMS*cosThetaCMS);
185 G4double PtX= PtProjCMS * std::cos(phi);
186 G4double PtY= PtProjCMS * std::sin(phi);
187 Fproj.setPx(PtX);
188 Fproj.setPy(PtY);
189 Fproj.setE(std::sqrt(PtX*PtX+PtY*PtY+PtZ*PtZ+Mproj*Mproj));
190 T = -(Pproj-Fproj).mag2();
191 }
192 else
193 {
194 // Simulation of strong interaction scattering
195
196 G4double Qmax = 2.*ptot/197.33; // in fm^-1
197
198 G4double Amag = 1.0; // A1 in Majorant funct:A1*exp(-q*A2)
199 G4double SlopeMag = 0.5; // A2 in Majorant funct:A1*exp(-q*A2)
200
201 G4double sig_pbarp = cs->GetAntiHadronNucleonTotCrSc(theAProton,energyPerN); //mb
202
203 fRa = 1.113*G4Pow::GetInstance()->Z13(A) -
204 0.227/G4Pow::GetInstance()->Z13(A);
205 if(A == 3) fRa=1.81;
206 if(A == 4) fRa=1.37;
207
208 if((A>=12.) && (A<27) ) fRa=fRa*0.85;
209 if((A>=27.) && (A<48) ) fRa=fRa*0.90;
210 if((A>=48.) && (A<65) ) fRa=fRa*0.95;
211
212 G4double Ref2 = XsTotalHadronic/10./2./pi; // in fm^2
213 G4double ceff2 =0;
214 G4double rho = 0;
215
216 if ((theParticle == theAProton) || (theParticle == theANeutron))
217 {
218 if(theTargetDef == theProton)
219 {
220 // Determination of the real part of Pbar+N amplitude
221 if(Plab < 610.)
222 { rho = 1.3347-10.342*Plab/1000.+22.277*Plab/1000.*Plab/1000.-
223 13.634*Plab/1000.*Plab/1000.*Plab/1000. ;}
224 if((Plab < 5500.)&&(Plab >= 610.) )
225 { rho = 0.22; }
226 if((Plab >= 5500.)&&(Plab < 12300.) )
227 { rho = -0.32; }
228 if( Plab >= 12300.)
229 { rho = 0.135-2.26/(std::sqrt(S)) ;}
230 Ref2 = 0.35 + 0.9/std::sqrt(std::sqrt(S-4.*0.88))+0.04*G4Log(S) ;
231 ceff2 = 0.375 - 2./S + 0.44/(sqr(S-4.)+1.5) ;
232 Ref2 =Ref2*Ref2;
233 ceff2 = ceff2*ceff2;
234 }
235
236 if( (Z==1)&&(A==2) )
237 {
238 Ref2 = fRa*fRa - 0.28 + 0.019 * sig_pbarp + 2.06e-6 * sig_pbarp*sig_pbarp;
239 ceff2 = 0.297 + 7.853e-04*sig_pbarp + 0.2899*G4Exp(-0.03*sig_pbarp);
240 }
241 if( (Z==1)&&(A==3) )
242 {
243 Ref2 = fRa*fRa - 1.36 + 0.025 * sig_pbarp - 3.69e-7 * sig_pbarp*sig_pbarp;
244 ceff2 = 0.149 + 7.091e-04*sig_pbarp + 0.3743*G4Exp(-0.03*sig_pbarp);
245 }
246 if( (Z==2)&&(A==3) )
247 {
248 Ref2 = fRa*fRa - 1.36 + 0.025 * sig_pbarp - 3.69e-7 * sig_pbarp*sig_pbarp;
249 ceff2 = 0.149 + 7.091e-04*sig_pbarp + 0.3743*G4Exp(-0.03*sig_pbarp);
250 }
251 if( (Z==2)&&(A==4) )
252 {
253 Ref2 = fRa*fRa -0.46 +0.03*sig_pbarp - 2.98e-6*sig_pbarp*sig_pbarp;
254 ceff2= 0.078 + 6.657e-4*sig_pbarp + 0.3359*G4Exp(-0.03*sig_pbarp);
255 }
256 if(Z>2)
257 {
258 Ref2 = fRa*fRa +2.48*0.01*sig_pbarp*fRa - 2.23e-6*sig_pbarp*sig_pbarp*fRa*fRa;
259 ceff2 = 0.16+3.3e-4*sig_pbarp+0.35*G4Exp(-0.03*sig_pbarp);
260 }
261 } // End of if ((theParticle == theAProton) || (theParticle == theANeutron))
262
263 if (theParticle == theADeuteron)
264 {
265 if(theTargetDef == theProton)
266 {
267 ceff2 = 0.297 + 7.853e-04*sig_pbarp + 0.2899*G4Exp(-0.03*sig_pbarp);
268 }
269 if(theTargetDef == theDeuteron)
270 {
271 ceff2 = 0.65 + 3.0e-4*sig_pbarp + 0.55 * G4Exp(-0.03*sig_pbarp);
272 }
273 if( (theTargetDef == G4Triton::Triton()) || (theTargetDef == G4He3::He3() ) )
274 {
275 ceff2 = 0.57 + 2.5e-4*sig_pbarp + 0.65 * G4Exp(-0.02*sig_pbarp);
276 }
277 if(theTargetDef == theAlpha)
278 {
279 ceff2 = 0.40 + 3.5e-4 *sig_pbarp + 0.45 * G4Exp(-0.02*sig_pbarp);
280 }
281 if(Z>2)
282 {
283 ceff2 = 0.38 + 2.0e-4 *sig_pbarp + 0.5 * G4Exp(-0.03*sig_pbarp);
284 }
285 }
286
287 if( (theParticle ==theAHe3) || (theParticle ==theATriton) )
288 {
289 if(theTargetDef == theProton)
290 {
291 ceff2 = 0.149 + 7.091e-04*sig_pbarp + 0.3743*G4Exp(-0.03*sig_pbarp);
292 }
293 if(theTargetDef == theDeuteron)
294 {
295 ceff2 = 0.57 + 2.5e-4*sig_pbarp + 0.65 * G4Exp(-0.02*sig_pbarp);
296 }
297 if( (theTargetDef == G4Triton::Triton()) || (theTargetDef == G4He3::He3() ) )
298 {
299 ceff2 = 0.39 + 2.7e-4*sig_pbarp + 0.7 * G4Exp(-0.02*sig_pbarp);
300 }
301 if(theTargetDef == theAlpha)
302 {
303 ceff2 = 0.24 + 3.5e-4*sig_pbarp + 0.75 * G4Exp(-0.03*sig_pbarp);
304 }
305 if(Z>2)
306 {
307 ceff2 = 0.26 + 2.2e-4*sig_pbarp + 0.33*G4Exp(-0.03*sig_pbarp);
308 }
309 }
310
311 if (theParticle == theAAlpha)
312 {
313 if(theTargetDef == theProton)
314 {
315 ceff2= 0.078 + 6.657e-4*sig_pbarp + 0.3359*G4Exp(-0.03*sig_pbarp);
316 }
317 if(theTargetDef == theDeuteron)
318 {
319 ceff2 = 0.40 + 3.5e-4 *sig_pbarp + 0.45 * G4Exp(-0.02*sig_pbarp);
320 }
321 if( (theTargetDef == G4Triton::Triton()) || (theTargetDef == G4He3::He3() ) )
322 {
323 ceff2 = 0.24 + 3.5e-4*sig_pbarp + 0.75 * G4Exp(-0.03*sig_pbarp);
324 }
325 if(theTargetDef == theAlpha)
326 {
327 ceff2 = 0.17 + 3.5e-4*sig_pbarp + 0.45 * G4Exp(-0.03*sig_pbarp);
328 }
329 if(Z>2)
330 {
331 ceff2 = 0.22 + 2.0e-4*sig_pbarp + 0.2 * G4Exp(-0.03*sig_pbarp);
332 }
333 }
334
335 fRef=std::sqrt(Ref2);
336 fceff = std::sqrt(ceff2);
337
338 G4double Q = 0.0 ;
339 G4double BracFunct;
340
341 const G4int maxNumberOfLoops = 10000;
342 G4int loopCounter = 0;
343 do
344 {
345 Q = -G4Log(1.-(1.- G4Exp(-SlopeMag * Qmax))* G4UniformRand() )/SlopeMag;
346 G4double x = fRef * Q;
347 BracFunct = ( ( sqr(BesselOneByArg(x))+sqr(rho/2. * BesselJzero(x)) )
348 * sqr(DampFactor(pi*fceff*Q))) /(Amag*G4Exp(-SlopeMag*Q));
349 BracFunct = BracFunct * Q;
350 }
351 while ( (G4UniformRand()>BracFunct) &&
352 ++loopCounter < maxNumberOfLoops ); /* Loop checking, 10.08.2015, A.Ribon */
353 if ( loopCounter >= maxNumberOfLoops ) {
354 fTetaCMS = 0.0;
355 return 0.0;
356 }
357
358 T= sqr(Q);
359 T*=3.893913e+4; // fm^(-2) -> MeV^2
360
361 } // End of simulation of strong interaction scattering
362
363 return T;
364}
G4double S(G4double temp)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double millibarn
Definition: G4SIunits.hh:86
static constexpr double pi
Definition: G4SIunits.hh:55
double mag() const
G4double BesselOneByArg(G4double z)
G4double CalculateParticleBeta(const G4ParticleDefinition *particle, G4double momentum)
G4double CalculateZommerfeld(G4double beta, G4double Z1, G4double Z2)
G4double DampFactor(G4double z)
G4double GetcosTeta1(G4double plab, G4int A)
G4double BesselJzero(G4double z)
G4double CalculateAm(G4double momentum, G4double n, G4double Z)
virtual G4double GetTotalElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
virtual G4double GetElasticElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
G4double GetAntiHadronNucleonTotCrSc(const G4ParticleDefinition *aParticle, G4double kinEnergy)
G4double Z13(G4int Z) const
Definition: G4Pow.hh:123
G4double energy(const ThreeVector &p, const G4double m)
static double Q[]
T sqr(const T &x)
Definition: templates.hh:128

References A, BesselJzero(), BesselOneByArg(), anonymous_namespace{G4PionRadiativeDecayChannel.cc}::beta, CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::boostVector(), CalculateAm(), CalculateParticleBeta(), CalculateZommerfeld(), cs, DampFactor(), CLHEP::HepLorentzVector::e(), G4INCL::KinematicsUtils::energy(), fbst, fceff, fptot, fRa, fRef, fTetaCMS, fTmax, fWaveVector, G4Exp(), G4Log(), G4UniformRand, G4ComponentAntiNuclNuclearXS::GetAntiHadronNucleonTotCrSc(), G4ParticleDefinition::GetBaryonNumber(), GetcosTeta1(), G4ComponentAntiNuclNuclearXS::GetElasticElementCrossSection(), G4Pow::GetInstance(), G4NucleiProperties::GetNuclearMass(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), G4ComponentAntiNuclNuclearXS::GetTotalElementCrossSection(), GeV, G4He3::He3(), CLHEP::Hep3Vector::mag(), CLHEP::HepLorentzVector::mag2(), MeV, millibarn, CLHEP::detail::n, pi, Q, S(), CLHEP::HepLorentzVector::setE(), CLHEP::HepLorentzVector::setPx(), CLHEP::HepLorentzVector::setPy(), CLHEP::HepLorentzVector::setPz(), sqr(), theAAlpha, theADeuteron, theAHe3, theAlpha, theANeutron, theAProton, theATriton, theDeuteron, theNeutron, theProton, G4Triton::Triton(), twopi, CLHEP::HepLorentzVector::vect(), Z, anonymous_namespace{paraMaker.cc}::Z1, and G4Pow::Z13().

Referenced by SampleThetaCMS(), and SampleThetaLab().

◆ SampleThetaCMS()

G4double G4AntiNuclElastic::SampleThetaCMS ( const G4ParticleDefinition p,
G4double  plab,
G4int  Z,
G4int  A 
)

Definition at line 368 of file G4AntiNuclElastic.cc.

370{
371 G4double T;
372 T = SampleInvariantT( p, plab, Z, A);
373
374 // NaN finder
375 if(!(T < 0.0 || T >= 0.0))
376 {
377 if (verboseLevel > 0)
378 {
379 G4cout << "G4DiffuseElastic:WARNING: A = " << A
380 << " mom(GeV)= " << plab/GeV
381 << " S-wave will be sampled"
382 << G4endl;
383 }
384 T = G4UniformRand()*fTmax;
385
386 }
387
388 if(fptot > 0.)
389 {
390 G4double cosTet=1.0-T/(2.*fptot*fptot);
391 if(cosTet > 1.0 ) cosTet= 1.;
392 if(cosTet < -1.0 ) cosTet=-1.;
393 fTetaCMS=std::acos(cosTet);
394 return fTetaCMS;
395 } else
396 {
397 return 2.*G4UniformRand()-1.;
398 }
399}
G4double SampleInvariantT(const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A) override

References A, fptot, fTetaCMS, fTmax, G4cout, G4endl, G4UniformRand, GeV, SampleInvariantT(), G4HadronicInteraction::verboseLevel, and Z.

◆ SampleThetaLab()

G4double G4AntiNuclElastic::SampleThetaLab ( const G4ParticleDefinition p,
G4double  plab,
G4int  Z,
G4int  A 
)

Definition at line 404 of file G4AntiNuclElastic.cc.

406{
407 G4double T;
408 T = SampleInvariantT( p, plab, Z, A);
409
410 // NaN finder
411 if(!(T < 0.0 || T >= 0.0))
412 {
413 if (verboseLevel > 0)
414 {
415 G4cout << "G4DiffuseElastic:WARNING: A = " << A
416 << " mom(GeV)= " << plab/GeV
417 << " S-wave will be sampled"
418 << G4endl;
419 }
420 T = G4UniformRand()*fTmax;
421 }
422
424
425 G4double cost(1.);
426 if(fTmax > 0.) {cost = 1. - 2.0*T/fTmax;}
427
428 G4double sint;
429 if( cost >= 1.0 )
430 {
431 cost = 1.0;
432 sint = 0.0;
433 }
434 else if( cost <= -1.0)
435 {
436 cost = -1.0;
437 sint = 0.0;
438 }
439 else
440 {
441 sint = std::sqrt((1.0-cost)*(1.0+cost));
442 }
443
444 G4double m1 = p->GetPDGMass();
445 G4ThreeVector v(sint*std::cos(phi),sint*std::sin(phi),cost);
446 v *= fptot;
447 G4LorentzVector nlv(v.x(),v.y(),v.z(),std::sqrt(fptot*fptot + m1*m1));
448
449 nlv.boost(fbst);
450
451 G4ThreeVector np = nlv.vect();
452 G4double theta = np.theta();
453 fThetaLab = theta;
454
455 return theta;
456}
double theta() const

References A, CLHEP::HepLorentzVector::boost(), fbst, fptot, fThetaLab, fTmax, G4cout, G4endl, G4UniformRand, G4ParticleDefinition::GetPDGMass(), GeV, SampleInvariantT(), CLHEP::Hep3Vector::theta(), twopi, CLHEP::HepLorentzVector::vect(), G4HadronicInteraction::verboseLevel, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), CLHEP::Hep3Vector::z(), and Z.

◆ SetEnergyMomentumCheckLevels()

void G4HadronicInteraction::SetEnergyMomentumCheckLevels ( G4double  relativeLevel,
G4double  absoluteLevel 
)
inlineinherited

Definition at line 149 of file G4HadronicInteraction.hh.

150 { epCheckLevels.first = relativeLevel;
151 epCheckLevels.second = absoluteLevel; }

References G4HadronicInteraction::epCheckLevels.

Referenced by G4BinaryCascade::G4BinaryCascade(), G4CascadeInterface::G4CascadeInterface(), and G4FTFModel::G4FTFModel().

◆ SetLowestEnergyLimit()

void G4HadronElastic::SetLowestEnergyLimit ( G4double  value)
inlineinherited

◆ SetMaxEnergy() [1/3]

void G4HadronicInteraction::SetMaxEnergy ( const G4double  anEnergy)
inlineinherited

Definition at line 102 of file G4HadronicInteraction.hh.

103 { theMaxEnergy = anEnergy; }

References G4HadronicInteraction::theMaxEnergy.

Referenced by G4HadronicInteraction::ActivateFor(), G4IonINCLXXPhysics::AddProcess(), G4BertiniElectroNuclearBuilder::Build(), G4LENDBertiniGammaElectroNuclearBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronPHPBuilder::Build(), G4AlphaPHPBuilder::Build(), G4BertiniKaonBuilder::Build(), G4BertiniNeutronBuilder::Build(), G4BertiniPiKBuilder::Build(), G4BertiniPionBuilder::Build(), G4BertiniProtonBuilder::Build(), G4BinaryAlphaBuilder::Build(), G4BinaryDeuteronBuilder::Build(), G4BinaryHe3Builder::Build(), G4BinaryNeutronBuilder::Build(), G4BinaryPiKBuilder::Build(), G4BinaryPionBuilder::Build(), G4BinaryProtonBuilder::Build(), G4BinaryTritonBuilder::Build(), G4DeuteronPHPBuilder::Build(), G4FTFBinaryKaonBuilder::Build(), G4FTFBinaryNeutronBuilder::Build(), G4FTFBinaryPionBuilder::Build(), G4FTFBinaryProtonBuilder::Build(), G4FTFPAntiBarionBuilder::Build(), G4FTFPKaonBuilder::Build(), G4FTFPNeutronBuilder::Build(), G4FTFPPiKBuilder::Build(), G4FTFPPionBuilder::Build(), G4FTFPProtonBuilder::Build(), G4He3PHPBuilder::Build(), G4HyperonFTFPBuilder::Build(), G4HyperonQGSPBuilder::Build(), G4INCLXXNeutronBuilder::Build(), G4INCLXXPionBuilder::Build(), G4INCLXXProtonBuilder::Build(), G4PrecoNeutronBuilder::Build(), G4PrecoProtonBuilder::Build(), G4ProtonPHPBuilder::Build(), G4QGSBinaryKaonBuilder::Build(), G4QGSBinaryNeutronBuilder::Build(), G4QGSBinaryPiKBuilder::Build(), G4QGSBinaryPionBuilder::Build(), G4QGSBinaryProtonBuilder::Build(), G4QGSPAntiBarionBuilder::Build(), G4QGSPKaonBuilder::Build(), G4QGSPLundStrFragmProtonBuilder::Build(), G4QGSPNeutronBuilder::Build(), G4QGSPPiKBuilder::Build(), G4QGSPPionBuilder::Build(), G4TritonPHPBuilder::Build(), G4QGSPProtonBuilder::Build(), G4HadronicBuilder::BuildFTFP_BERT(), G4HadronicBuilder::BuildFTFQGSP_BERT(), G4QGSBuilder::BuildModel(), G4VHadronPhysics::BuildModel(), G4HadronicBuilder::BuildQGSP_FTFP_BERT(), G4EmExtraPhysics::ConstructGammaElectroNuclear(), LBE::ConstructHad(), G4EmExtraPhysics::ConstructLENDGammaNuclear(), G4HadronDElasticPhysics::ConstructProcess(), G4HadronElasticPhysics::ConstructProcess(), G4HadronHElasticPhysics::ConstructProcess(), G4IonINCLXXPhysics::ConstructProcess(), G4IonPhysics::ConstructProcess(), G4IonPhysicsPHP::ConstructProcess(), G4IonQMDPhysics::ConstructProcess(), G4ANuElNucleusNcModel::G4ANuElNucleusNcModel(), G4ANuMuNucleusNcModel::G4ANuMuNucleusNcModel(), G4BertiniKaonBuilder::G4BertiniKaonBuilder(), G4BertiniPiKBuilder::G4BertiniPiKBuilder(), G4BertiniPionBuilder::G4BertiniPionBuilder(), G4BinaryCascade::G4BinaryCascade(), G4BinaryPiKBuilder::G4BinaryPiKBuilder(), G4BinaryPionBuilder::G4BinaryPionBuilder(), G4ChargeExchange::G4ChargeExchange(), G4DiffuseElastic::G4DiffuseElastic(), G4DiffuseElasticV2::G4DiffuseElasticV2(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4EMDissociation::G4EMDissociation(), G4FissLib::G4FissLib(), G4FTFBinaryKaonBuilder::G4FTFBinaryKaonBuilder(), G4FTFBinaryNeutronBuilder::G4FTFBinaryNeutronBuilder(), G4FTFBinaryPiKBuilder::G4FTFBinaryPiKBuilder(), G4FTFBinaryPionBuilder::G4FTFBinaryPionBuilder(), G4FTFBinaryProtonBuilder::G4FTFBinaryProtonBuilder(), G4FTFPAntiBarionBuilder::G4FTFPAntiBarionBuilder(), G4FTFPKaonBuilder::G4FTFPKaonBuilder(), G4FTFPNeutronBuilder::G4FTFPNeutronBuilder(), G4FTFPPiKBuilder::G4FTFPPiKBuilder(), G4FTFPPionBuilder::G4FTFPPionBuilder(), G4FTFPProtonBuilder::G4FTFPProtonBuilder(), G4HadronElastic::G4HadronElastic(), G4HadronicAbsorptionFritiof::G4HadronicAbsorptionFritiof(), G4HadronicAbsorptionFritiofWithBinaryCascade::G4HadronicAbsorptionFritiofWithBinaryCascade(), G4hhElastic::G4hhElastic(), G4HyperonFTFPBuilder::G4HyperonFTFPBuilder(), G4HyperonQGSPBuilder::G4HyperonQGSPBuilder(), G4INCLXXPionBuilder::G4INCLXXPionBuilder(), G4LEHadronProtonElastic::G4LEHadronProtonElastic(), G4LENDModel::G4LENDModel(), G4LEnp::G4LEnp(), G4LEpp::G4LEpp(), G4LFission::G4LFission(), G4LowEGammaNuclearModel::G4LowEGammaNuclearModel(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4NeutrinoElectronCcModel::G4NeutrinoElectronCcModel(), G4NeutrinoElectronNcModel::G4NeutrinoElectronNcModel(), G4NeutrinoNucleusModel::G4NeutrinoNucleusModel(), G4NeutronElectronElModel::G4NeutronElectronElModel(), G4NeutronRadCapture::G4NeutronRadCapture(), G4NuclNuclDiffuseElastic::G4NuclNuclDiffuseElastic(), G4NuElNucleusNcModel::G4NuElNucleusNcModel(), G4NuMuNucleusNcModel::G4NuMuNucleusNcModel(), G4ParticleHPCapture::G4ParticleHPCapture(), G4ParticleHPElastic::G4ParticleHPElastic(), G4ParticleHPFission::G4ParticleHPFission(), G4ParticleHPInelastic::G4ParticleHPInelastic(), G4ParticleHPThermalScattering::G4ParticleHPThermalScattering(), G4QGSPAntiBarionBuilder::G4QGSPAntiBarionBuilder(), G4WilsonAbrasionModel::G4WilsonAbrasionModel(), G4HadronPhysicsFTFP_BERT_HP::Neutron(), G4HadronPhysicsINCLXX::Neutron(), G4HadronPhysicsQGSP_BERT_HP::Neutron(), G4HadronPhysicsQGSP_BIC_HP::Neutron(), G4HadronPhysicsShielding::Neutron(), and G4VHadronPhysics::NewModel().

◆ SetMaxEnergy() [2/3]

void G4HadronicInteraction::SetMaxEnergy ( G4double  anEnergy,
const G4Element anElement 
)
inherited

Definition at line 151 of file G4HadronicInteraction.cc.

153{
154 Block();
155 if(!theMaxEnergyListElements.empty()) {
156 for(auto & elmlist : theMaxEnergyListElements) {
157 if( anElement == elmlist.second ) {
158 elmlist.first = anEnergy;
159 return;
160 }
161 }
162 }
163 theMaxEnergyListElements.push_back(std::pair<G4double, const G4Element *>(anEnergy, anElement));
164}

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theMaxEnergyListElements.

◆ SetMaxEnergy() [3/3]

void G4HadronicInteraction::SetMaxEnergy ( G4double  anEnergy,
const G4Material aMaterial 
)
inherited

Definition at line 166 of file G4HadronicInteraction.cc.

167{
168 Block();
169 if(!theMaxEnergyList.empty()) {
170 for(auto & matlist: theMaxEnergyList) {
171 if( aMaterial == matlist.second ) {
172 matlist.first = anEnergy;
173 return;
174 }
175 }
176 }
177 theMaxEnergyList.push_back(std::pair<G4double, const G4Material *>(anEnergy, aMaterial));
178}

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theMaxEnergyList.

◆ SetMinEnergy() [1/3]

void G4HadronicInteraction::SetMinEnergy ( G4double  anEnergy)
inlineinherited

Definition at line 89 of file G4HadronicInteraction.hh.

90 { theMinEnergy = anEnergy; }

References G4HadronicInteraction::theMinEnergy.

Referenced by G4HadronicInteraction::ActivateFor(), G4BertiniElectroNuclearBuilder::Build(), G4LENDBertiniGammaElectroNuclearBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronPHPBuilder::Build(), G4AlphaPHPBuilder::Build(), G4BertiniKaonBuilder::Build(), G4BertiniNeutronBuilder::Build(), G4BertiniPiKBuilder::Build(), G4BertiniPionBuilder::Build(), G4BertiniProtonBuilder::Build(), G4BinaryAlphaBuilder::Build(), G4BinaryDeuteronBuilder::Build(), G4BinaryHe3Builder::Build(), G4BinaryNeutronBuilder::Build(), G4BinaryPiKBuilder::Build(), G4BinaryPionBuilder::Build(), G4BinaryProtonBuilder::Build(), G4BinaryTritonBuilder::Build(), G4DeuteronPHPBuilder::Build(), G4FTFBinaryKaonBuilder::Build(), G4FTFBinaryNeutronBuilder::Build(), G4FTFBinaryPiKBuilder::Build(), G4FTFBinaryPionBuilder::Build(), G4FTFBinaryProtonBuilder::Build(), G4FTFPAntiBarionBuilder::Build(), G4FTFPKaonBuilder::Build(), G4FTFPNeutronBuilder::Build(), G4FTFPPiKBuilder::Build(), G4FTFPPionBuilder::Build(), G4FTFPProtonBuilder::Build(), G4He3PHPBuilder::Build(), G4HyperonFTFPBuilder::Build(), G4HyperonQGSPBuilder::Build(), G4INCLXXNeutronBuilder::Build(), G4INCLXXPionBuilder::Build(), G4INCLXXProtonBuilder::Build(), G4PrecoNeutronBuilder::Build(), G4PrecoProtonBuilder::Build(), G4ProtonPHPBuilder::Build(), G4QGSBinaryKaonBuilder::Build(), G4QGSBinaryNeutronBuilder::Build(), G4QGSBinaryPiKBuilder::Build(), G4QGSBinaryPionBuilder::Build(), G4QGSBinaryProtonBuilder::Build(), G4QGSPAntiBarionBuilder::Build(), G4QGSPKaonBuilder::Build(), G4QGSPLundStrFragmProtonBuilder::Build(), G4QGSPNeutronBuilder::Build(), G4QGSPPiKBuilder::Build(), G4QGSPPionBuilder::Build(), G4TritonPHPBuilder::Build(), G4QGSPProtonBuilder::Build(), G4QGSBuilder::BuildModel(), G4VHadronPhysics::BuildModel(), G4EmExtraPhysics::ConstructGammaElectroNuclear(), LBE::ConstructHad(), G4EmExtraPhysics::ConstructLENDGammaNuclear(), G4HadronElasticPhysicsHP::ConstructProcess(), G4HadronElasticPhysicsLEND::ConstructProcess(), G4HadronElasticPhysicsPHP::ConstructProcess(), G4HadronDElasticPhysics::ConstructProcess(), G4HadronElasticPhysics::ConstructProcess(), G4HadronHElasticPhysics::ConstructProcess(), G4IonElasticPhysics::ConstructProcess(), G4IonINCLXXPhysics::ConstructProcess(), G4IonPhysics::ConstructProcess(), G4IonPhysicsPHP::ConstructProcess(), G4IonQMDPhysics::ConstructProcess(), G4ANuElNucleusNcModel::G4ANuElNucleusNcModel(), G4ANuMuNucleusNcModel::G4ANuMuNucleusNcModel(), G4BertiniKaonBuilder::G4BertiniKaonBuilder(), G4BertiniPiKBuilder::G4BertiniPiKBuilder(), G4BertiniPionBuilder::G4BertiniPionBuilder(), G4BinaryCascade::G4BinaryCascade(), G4BinaryPiKBuilder::G4BinaryPiKBuilder(), G4BinaryPionBuilder::G4BinaryPionBuilder(), G4ChargeExchange::G4ChargeExchange(), G4DiffuseElastic::G4DiffuseElastic(), G4DiffuseElasticV2::G4DiffuseElasticV2(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4EMDissociation::G4EMDissociation(), G4FissLib::G4FissLib(), G4FTFBinaryKaonBuilder::G4FTFBinaryKaonBuilder(), G4FTFBinaryNeutronBuilder::G4FTFBinaryNeutronBuilder(), G4FTFBinaryPiKBuilder::G4FTFBinaryPiKBuilder(), G4FTFBinaryPionBuilder::G4FTFBinaryPionBuilder(), G4FTFBinaryProtonBuilder::G4FTFBinaryProtonBuilder(), G4FTFPAntiBarionBuilder::G4FTFPAntiBarionBuilder(), G4FTFPKaonBuilder::G4FTFPKaonBuilder(), G4FTFPNeutronBuilder::G4FTFPNeutronBuilder(), G4FTFPPiKBuilder::G4FTFPPiKBuilder(), G4FTFPPionBuilder::G4FTFPPionBuilder(), G4FTFPProtonBuilder::G4FTFPProtonBuilder(), G4HadronElastic::G4HadronElastic(), G4HadronicAbsorptionBertini::G4HadronicAbsorptionBertini(), G4HadronicAbsorptionFritiof::G4HadronicAbsorptionFritiof(), G4HadronicAbsorptionFritiofWithBinaryCascade::G4HadronicAbsorptionFritiofWithBinaryCascade(), G4hhElastic::G4hhElastic(), G4HyperonFTFPBuilder::G4HyperonFTFPBuilder(), G4HyperonQGSPBuilder::G4HyperonQGSPBuilder(), G4INCLXXPionBuilder::G4INCLXXPionBuilder(), G4LEHadronProtonElastic::G4LEHadronProtonElastic(), G4LENDModel::G4LENDModel(), G4LEnp::G4LEnp(), G4LEpp::G4LEpp(), G4LFission::G4LFission(), G4LowEGammaNuclearModel::G4LowEGammaNuclearModel(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4NeutrinoElectronCcModel::G4NeutrinoElectronCcModel(), G4NeutrinoElectronNcModel::G4NeutrinoElectronNcModel(), G4NeutrinoNucleusModel::G4NeutrinoNucleusModel(), G4NeutronElectronElModel::G4NeutronElectronElModel(), G4NeutronRadCapture::G4NeutronRadCapture(), G4NuclNuclDiffuseElastic::G4NuclNuclDiffuseElastic(), G4NuElNucleusNcModel::G4NuElNucleusNcModel(), G4NuMuNucleusNcModel::G4NuMuNucleusNcModel(), G4ParticleHPCapture::G4ParticleHPCapture(), G4ParticleHPElastic::G4ParticleHPElastic(), G4ParticleHPFission::G4ParticleHPFission(), G4ParticleHPInelastic::G4ParticleHPInelastic(), G4ParticleHPThermalScattering::G4ParticleHPThermalScattering(), G4QGSPAntiBarionBuilder::G4QGSPAntiBarionBuilder(), G4WilsonAbrasionModel::G4WilsonAbrasionModel(), G4NeutrinoElectronCcModel::IsApplicable(), G4HadronPhysicsFTFP_BERT_HP::Neutron(), G4HadronPhysicsINCLXX::Neutron(), G4HadronPhysicsQGSP_BERT_HP::Neutron(), G4HadronPhysicsQGSP_BIC_HP::Neutron(), G4HadronPhysicsShielding::Neutron(), and G4VHadronPhysics::NewModel().

◆ SetMinEnergy() [2/3]

void G4HadronicInteraction::SetMinEnergy ( G4double  anEnergy,
const G4Element anElement 
)
inherited

Definition at line 101 of file G4HadronicInteraction.cc.

103{
104 Block();
105 if(!theMinEnergyListElements.empty()) {
106 for(auto & elmlist : theMinEnergyListElements) {
107 if( anElement == elmlist.second ) {
108 elmlist.first = anEnergy;
109 return;
110 }
111 }
112 }
113 theMinEnergyListElements.push_back(std::pair<G4double, const G4Element *>(anEnergy, anElement));
114}

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theMinEnergyListElements.

◆ SetMinEnergy() [3/3]

void G4HadronicInteraction::SetMinEnergy ( G4double  anEnergy,
const G4Material aMaterial 
)
inherited

Definition at line 116 of file G4HadronicInteraction.cc.

118{
119 Block();
120 if(!theMinEnergyList.empty()) {
121 for(auto & matlist : theMinEnergyList) {
122 if( aMaterial == matlist.second ) {
123 matlist.first = anEnergy;
124 return;
125 }
126 }
127 }
128 theMinEnergyList.push_back(std::pair<G4double, const G4Material *>(anEnergy, aMaterial));
129}

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theMinEnergyList.

◆ SetModelName()

void G4HadronicInteraction::SetModelName ( const G4String nam)
inlineprotectedinherited

◆ SetRecoilEnergyThreshold()

void G4HadronicInteraction::SetRecoilEnergyThreshold ( G4double  val)
inlineinherited

◆ SetVerboseLevel()

void G4HadronicInteraction::SetVerboseLevel ( G4int  value)
inlineinherited

Field Documentation

◆ cs

G4ComponentAntiNuclNuclearXS* G4AntiNuclElastic::cs
private

◆ epCheckLevels

std::pair<G4double, G4double> G4HadronicInteraction::epCheckLevels
privateinherited

◆ fAm

G4double G4AntiNuclElastic::fAm
private

Definition at line 97 of file G4AntiNuclElastic.hh.

Referenced by CalculateAm(), and G4AntiNuclElastic().

◆ fBeta

G4double G4AntiNuclElastic::fBeta
private

Definition at line 95 of file G4AntiNuclElastic.hh.

Referenced by CalculateParticleBeta(), and G4AntiNuclElastic().

◆ fbst

G4ThreeVector G4AntiNuclElastic::fbst
private

Definition at line 102 of file G4AntiNuclElastic.hh.

Referenced by SampleInvariantT(), and SampleThetaLab().

◆ fceff

G4double G4AntiNuclElastic::fceff
private

Definition at line 100 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic(), and SampleInvariantT().

◆ fParticle

const G4ParticleDefinition* G4AntiNuclElastic::fParticle
private

Definition at line 90 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic().

◆ fptot

G4double G4AntiNuclElastic::fptot
private

◆ fRa

G4double G4AntiNuclElastic::fRa
private

Definition at line 98 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic(), and SampleInvariantT().

◆ fRef

G4double G4AntiNuclElastic::fRef
private

Definition at line 99 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic(), and SampleInvariantT().

◆ fTetaCMS

G4double G4AntiNuclElastic::fTetaCMS
private

Definition at line 92 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic(), SampleInvariantT(), and SampleThetaCMS().

◆ fThetaLab

G4double G4AntiNuclElastic::fThetaLab
private

Definition at line 93 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic(), and SampleThetaLab().

◆ fTmax

G4double G4AntiNuclElastic::fTmax
private

◆ fWaveVector

G4double G4AntiNuclElastic::fWaveVector
private

Definition at line 94 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic(), and SampleInvariantT().

◆ fZommerfeld

G4double G4AntiNuclElastic::fZommerfeld
private

Definition at line 96 of file G4AntiNuclElastic.hh.

Referenced by CalculateZommerfeld(), and G4AntiNuclElastic().

◆ isBlocked

G4bool G4HadronicInteraction::isBlocked
protectedinherited

◆ lowestEnergyLimit

G4double G4HadronElastic::lowestEnergyLimit
privateinherited

◆ nwarn

G4int G4HadronElastic::nwarn
privateinherited

◆ pLocalTmax

G4double G4HadronElastic::pLocalTmax
protectedinherited

◆ recoilEnergyThreshold

G4double G4HadronicInteraction::recoilEnergyThreshold
privateinherited

◆ registry

G4HadronicInteractionRegistry* G4HadronicInteraction::registry
privateinherited

◆ secID

G4int G4HadronElastic::secID
protectedinherited

◆ theAAlpha

G4ParticleDefinition* G4AntiNuclElastic::theAAlpha
private

Definition at line 110 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic(), and SampleInvariantT().

◆ theADeuteron

G4ParticleDefinition* G4AntiNuclElastic::theADeuteron
private

Definition at line 108 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic(), and SampleInvariantT().

◆ theAHe3

G4ParticleDefinition* G4AntiNuclElastic::theAHe3
private

Definition at line 111 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic(), and SampleInvariantT().

◆ theAlpha

G4ParticleDefinition* G4AntiNuclElastic::theAlpha
private

Definition at line 116 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic(), and SampleInvariantT().

◆ theANeutron

G4ParticleDefinition* G4AntiNuclElastic::theANeutron
private

Definition at line 107 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic(), and SampleInvariantT().

◆ theAProton

G4ParticleDefinition* G4AntiNuclElastic::theAProton
private

Definition at line 106 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic(), and SampleInvariantT().

◆ theATriton

G4ParticleDefinition* G4AntiNuclElastic::theATriton
private

Definition at line 109 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic(), and SampleInvariantT().

◆ theBlockedList

std::vector<const G4Material *> G4HadronicInteraction::theBlockedList
privateinherited

◆ theBlockedListElements

std::vector<const G4Element *> G4HadronicInteraction::theBlockedListElements
privateinherited

◆ theDeuteron

G4ParticleDefinition* G4AntiNuclElastic::theDeuteron
private

Definition at line 115 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic(), and SampleInvariantT().

◆ theMaxEnergy

G4double G4HadronicInteraction::theMaxEnergy
protectedinherited

◆ theMaxEnergyList

std::vector<std::pair<G4double, const G4Material *> > G4HadronicInteraction::theMaxEnergyList
privateinherited

◆ theMaxEnergyListElements

std::vector<std::pair<G4double, const G4Element *> > G4HadronicInteraction::theMaxEnergyListElements
privateinherited

◆ theMinEnergy

G4double G4HadronicInteraction::theMinEnergy
protectedinherited

◆ theMinEnergyList

std::vector<std::pair<G4double, const G4Material *> > G4HadronicInteraction::theMinEnergyList
privateinherited

◆ theMinEnergyListElements

std::vector<std::pair<G4double, const G4Element *> > G4HadronicInteraction::theMinEnergyListElements
privateinherited

◆ theModelName

G4String G4HadronicInteraction::theModelName
privateinherited

◆ theNeutron

G4ParticleDefinition* G4AntiNuclElastic::theNeutron
private

Definition at line 114 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic(), and SampleInvariantT().

◆ theParticleChange

G4HadFinalState G4HadronicInteraction::theParticleChange
protectedinherited

Definition at line 172 of file G4HadronicInteraction.hh.

Referenced by G4WilsonAbrasionModel::ApplyYourself(), G4EMDissociation::ApplyYourself(), G4LENDCapture::ApplyYourself(), G4LENDElastic::ApplyYourself(), G4LENDFission::ApplyYourself(), G4LENDInelastic::ApplyYourself(), G4ElectroVDNuclearModel::ApplyYourself(), G4ParticleHPThermalScattering::ApplyYourself(), G4NeutrinoElectronNcModel::ApplyYourself(), G4NeutronElectronElModel::ApplyYourself(), G4LFission::ApplyYourself(), G4ANuElNucleusCcModel::ApplyYourself(), G4ANuElNucleusNcModel::ApplyYourself(), G4ANuMuNucleusCcModel::ApplyYourself(), G4ANuMuNucleusNcModel::ApplyYourself(), G4MuonVDNuclearModel::ApplyYourself(), G4NeutrinoElectronCcModel::ApplyYourself(), G4NuElNucleusCcModel::ApplyYourself(), G4NuElNucleusNcModel::ApplyYourself(), G4NuMuNucleusCcModel::ApplyYourself(), G4NuMuNucleusNcModel::ApplyYourself(), G4QMDReaction::ApplyYourself(), G4NeutronRadCapture::ApplyYourself(), G4LowEGammaNuclearModel::ApplyYourself(), G4ChargeExchange::ApplyYourself(), G4HadronElastic::ApplyYourself(), G4LEHadronProtonElastic::ApplyYourself(), G4LEnp::ApplyYourself(), G4LEpp::ApplyYourself(), G4BinaryCascade::ApplyYourself(), G4CascadeInterface::ApplyYourself(), G4LMsdGenerator::ApplyYourself(), G4ElectroVDNuclearModel::CalculateEMVertex(), G4MuonVDNuclearModel::CalculateEMVertex(), G4ElectroVDNuclearModel::CalculateHadronicVertex(), G4MuonVDNuclearModel::CalculateHadronicVertex(), G4NeutrinoNucleusModel::CoherentPion(), G4CascadeInterface::copyOutputToHadronicResult(), G4BinaryCascade::DebugEpConservation(), G4BinaryCascade::DebugFinalEpConservation(), G4NeutrinoNucleusModel::FinalBarion(), G4NeutrinoNucleusModel::FinalMeson(), G4WilsonAbrasionModel::GetAbradedNucleons(), G4CascadeInterface::NoInteraction(), G4CascadeInterface::Propagate(), G4NeutrinoNucleusModel::RecoilDeexcitation(), G4LEHadronProtonElastic::~G4LEHadronProtonElastic(), G4LEnp::~G4LEnp(), and G4LFission::~G4LFission().

◆ theProton

G4ParticleDefinition* G4AntiNuclElastic::theProton
private

Definition at line 113 of file G4AntiNuclElastic.hh.

Referenced by G4AntiNuclElastic(), and SampleInvariantT().

◆ verboseLevel

G4int G4HadronicInteraction::verboseLevel
protectedinherited

Definition at line 177 of file G4HadronicInteraction.hh.

Referenced by G4WilsonAbrasionModel::ApplyYourself(), G4EMDissociation::ApplyYourself(), G4LFission::ApplyYourself(), G4MuMinusCapturePrecompound::ApplyYourself(), G4NeutronRadCapture::ApplyYourself(), G4LowEGammaNuclearModel::ApplyYourself(), G4ChargeExchange::ApplyYourself(), G4HadronElastic::ApplyYourself(), G4LEHadronProtonElastic::ApplyYourself(), G4LEnp::ApplyYourself(), G4LEpp::ApplyYourself(), G4CascadeInterface::ApplyYourself(), G4CascadeInterface::checkFinalResult(), G4CascadeInterface::copyOutputToHadronicResult(), G4CascadeInterface::copyOutputToReactionProducts(), G4LENDModel::create_used_target_map(), G4CascadeInterface::createBullet(), G4CascadeInterface::createTarget(), G4ElasticHadrNucleusHE::DefineHadronValues(), G4ElasticHadrNucleusHE::FillData(), G4ElasticHadrNucleusHE::FillFq2(), G4DiffuseElastic::G4DiffuseElastic(), G4DiffuseElasticV2::G4DiffuseElasticV2(), G4ElasticHadrNucleusHE::G4ElasticHadrNucleusHE(), G4EMDissociation::G4EMDissociation(), G4hhElastic::G4hhElastic(), G4NuclNuclDiffuseElastic::G4NuclNuclDiffuseElastic(), G4WilsonAbrasionModel::G4WilsonAbrasionModel(), G4ElasticHadrNucleusHE::GetFt(), G4ElasticHadrNucleusHE::GetLightFq2(), G4ElasticHadrNucleusHE::GetQ2_2(), G4HadronicInteraction::GetVerboseLevel(), G4ElasticHadrNucleusHE::HadronNucleusQ2_2(), G4ElasticHadrNucleusHE::HadronProtonQ2(), G4LFission::init(), G4DiffuseElastic::Initialise(), G4DiffuseElasticV2::Initialise(), G4NuclNuclDiffuseElastic::Initialise(), G4DiffuseElastic::InitialiseOnFly(), G4DiffuseElasticV2::InitialiseOnFly(), G4NuclNuclDiffuseElastic::InitialiseOnFly(), G4CascadeInterface::makeDynamicParticle(), G4CascadeInterface::NoInteraction(), G4CascadeInterface::Propagate(), G4ElasticHadrNucleusHE::SampleInvariantT(), SampleThetaCMS(), G4DiffuseElastic::SampleThetaLab(), G4NuclNuclDiffuseElastic::SampleThetaLab(), SampleThetaLab(), G4WilsonAbrasionModel::SetUseAblation(), G4HadronicInteraction::SetVerboseLevel(), G4WilsonAbrasionModel::SetVerboseLevel(), G4DiffuseElastic::ThetaCMStoThetaLab(), G4DiffuseElasticV2::ThetaCMStoThetaLab(), G4NuclNuclDiffuseElastic::ThetaCMStoThetaLab(), G4DiffuseElastic::ThetaLabToThetaCMS(), G4DiffuseElasticV2::ThetaLabToThetaCMS(), and G4NuclNuclDiffuseElastic::ThetaLabToThetaCMS().


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