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

#include <G4DiffuseElasticV2.hh>

Inheritance diagram for G4DiffuseElasticV2:
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)
 
void BuildAngleTable ()
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4double CalculateAm (G4double momentum, G4double n, G4double Z)
 
G4double CalculateNuclearRad (G4double A)
 
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)
 
 G4DiffuseElasticV2 ()
 
G4double GetDiffElasticSumProbA (G4double alpha)
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
G4double GetIntegrandFunction (G4double theta)
 
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 GetNuclearRadius ()
 
G4double GetRecoilEnergyThreshold () const
 
G4double GetScatteringAngle (G4int iMomentum, unsigned long iAngle, G4double position)
 
G4double GetSlopeCof (const G4int pdg)
 
G4int GetVerboseLevel () const
 
void Initialise ()
 
virtual void InitialiseModel ()
 
void InitialiseOnFly (G4double Z, G4double A)
 
virtual G4bool IsApplicable (const G4HadProjectile &, G4Nucleus &)
 
G4bool IsBlocked (const G4Element *anElement) const
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4double LowestEnergyLimit () const
 
void ModelDescription (std::ostream &) const override
 
G4double NeutronTuniform (G4int Z)
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
G4double SampleTableT (const G4ParticleDefinition *aParticle, G4double p, G4double Z, G4double A)
 
G4double SampleTableThetaCMS (const G4ParticleDefinition *aParticle, G4double p, G4double Z, G4double A)
 
G4double SampleThetaCMS (const G4ParticleDefinition *aParticle, G4double p, G4double A)
 
G4double SampleThetaLab (const G4HadProjectile *aParticle, G4double tmass, G4double A)
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
void SetHEModelLowLimit (G4double value)
 
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 SetPlabLowLimit (G4double value)
 
void SetQModelLowLimit (G4double value)
 
void SetRecoilEnergyThreshold (G4double val)
 
void SetRecoilKinEnergyLimit (G4double value)
 
void SetVerboseLevel (G4int value)
 
G4double ThetaCMStoThetaLab (const G4DynamicParticle *aParticle, G4double tmass, G4double thetaCMS)
 
G4double ThetaLabToThetaCMS (const G4DynamicParticle *aParticle, G4double tmass, G4double thetaLab)
 
virtual ~G4DiffuseElasticV2 ()
 

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 Attributes

std::pair< G4double, G4doubleepCheckLevels
 
G4bool fAddCoulomb
 
G4double fAm
 
unsigned long fAngleBin
 
G4double fAtomicNumber
 
G4double fAtomicWeight
 
G4double fBeta
 
std::vector< G4StringfElementNameVector
 
std::vector< G4doublefElementNumberVector
 
std::vector< std::vector< double > * > * fEnergyAngleVector
 
std::vector< std::vector< std::vector< double > * > * > fEnergyAngleVectorBank
 
G4int fEnergyBin
 
std::vector< std::vector< double > * > * fEnergySumVector
 
std::vector< std::vector< std::vector< double > * > * > fEnergySumVectorBank
 
G4PhysicsLogVectorfEnergyVector
 
G4double fNuclearRadius
 
const G4ParticleDefinitionfParticle
 
G4double fWaveVector
 
G4double fZommerfeld
 
G4double lowEnergyLimitHE
 
G4double lowEnergyLimitQ
 
G4double lowEnergyRecoilLimit
 
G4double lowestEnergyLimit
 
G4int nwarn
 
G4double plabLowLimit
 
G4double recoilEnergyThreshold
 
G4HadronicInteractionRegistryregistry
 
G4ParticleDefinitiontheAlpha
 
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 61 of file G4DiffuseElasticV2.hh.

Constructor & Destructor Documentation

◆ G4DiffuseElasticV2()

G4DiffuseElasticV2::G4DiffuseElasticV2 ( )

Definition at line 76 of file G4DiffuseElasticV2.cc.

77 : G4HadronElastic("DiffuseElasticV2"), fParticle(0)
78{
79 SetMinEnergy( 0.01*MeV );
81
82 verboseLevel = 0;
84 lowEnergyLimitQ = 0.0*GeV;
85 lowEnergyLimitHE = 0.0*GeV;
87 plabLowLimit = 20.0*MeV;
88
91
92 fEnergyBin = 300; // Increased from the original 200 to have no wider log-energy-bins up to 10 PeV
93 fAngleBin = 200;
94
96
99
100 fParticle = 0;
101 fWaveVector = 0.;
102 fAtomicWeight = 0.;
103 fAtomicNumber = 0.;
104 fNuclearRadius = 0.;
105 fBeta = 0.;
106 fZommerfeld = 0.;
107 fAm = 0.;
108 fAddCoulomb = false;
109}
static constexpr double keV
Definition: G4SIunits.hh:202
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
G4ParticleDefinition * theProton
const G4ParticleDefinition * fParticle
G4ParticleDefinition * theNeutron
G4PhysicsLogVector * fEnergyVector
std::vector< std::vector< double > * > * fEnergyAngleVector
std::vector< std::vector< double > * > * fEnergySumVector
G4HadronElastic(const G4String &name="hElasticLHEP")
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
static G4HadronicParameters * Instance()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4Proton * Proton()
Definition: G4Proton.cc:92

References fAddCoulomb, fAm, fAngleBin, fAtomicNumber, fAtomicWeight, fBeta, fEnergyAngleVector, fEnergyBin, fEnergySumVector, fEnergyVector, fNuclearRadius, fParticle, fWaveVector, fZommerfeld, G4HadronicInteraction::GetMaxEnergy(), GeV, G4HadronicParameters::Instance(), keV, lowEnergyLimitHE, lowEnergyLimitQ, lowEnergyRecoilLimit, lowestEnergyLimit, MeV, G4Neutron::Neutron(), plabLowLimit, G4Proton::Proton(), G4HadronicInteraction::SetMaxEnergy(), G4HadronicInteraction::SetMinEnergy(), G4HadronicInteraction::theMaxEnergy, G4HadronicInteraction::theMinEnergy, theNeutron, theProton, and G4HadronicInteraction::verboseLevel.

Referenced by BuildAngleTable().

◆ ~G4DiffuseElasticV2()

G4DiffuseElasticV2::~G4DiffuseElasticV2 ( )
virtual

Definition at line 115 of file G4DiffuseElasticV2.cc.

116{
117 if ( fEnergyVector )
118 {
119 delete fEnergyVector;
120 fEnergyVector = 0;
121 }
122}

References fEnergyVector.

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
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 G4DiffuseElasticV2::BesselJone ( G4double  z)
inline

Definition at line 267 of file G4DiffuseElasticV2.hh.

268{
269 G4double modvalue, value2, fact1, fact2, arg, shift, bessel;
270
271 modvalue = std::fabs(value);
272
273 if ( modvalue < 8.0 )
274 {
275 value2 = value*value;
276
277 fact1 = value*(72362614232.0 + value2*(-7895059235.0
278 + value2*( 242396853.1
279 + value2*(-2972611.439
280 + value2*( 15704.48260
281 + value2*(-30.16036606 ) ) ) ) ) );
282
283 fact2 = 144725228442.0 + value2*(2300535178.0
284 + value2*(18583304.74
285 + value2*(99447.43394
286 + value2*(376.9991397
287 + value2*1.0 ) ) ) );
288 bessel = fact1/fact2;
289 }
290 else
291 {
292 arg = 8.0/modvalue;
293
294 value2 = arg*arg;
295
296 shift = modvalue - 2.356194491;
297
298 fact1 = 1.0 + value2*( 0.183105e-2
299 + value2*(-0.3516396496e-4
300 + value2*(0.2457520174e-5
301 + value2*(-0.240337019e-6 ) ) ) );
302
303 fact2 = 0.04687499995 + value2*(-0.2002690873e-3
304 + value2*( 0.8449199096e-5
305 + value2*(-0.88228987e-6
306 + value2*0.105787412e-6 ) ) );
307
308 bessel = std::sqrt( 0.636619772/modvalue)*(std::cos(shift)*fact1 - arg*std::sin(shift)*fact2);
309
310 if (value < 0.0) bessel = -bessel;
311 }
312 return bessel;
313}

Referenced by BesselOneByArg(), and GetDiffElasticSumProbA().

◆ BesselJzero()

G4double G4DiffuseElasticV2::BesselJzero ( G4double  z)
inline

Definition at line 215 of file G4DiffuseElasticV2.hh.

216{
217 G4double modvalue, value2, fact1, fact2, arg, shift, bessel;
218
219 modvalue = std::fabs(value);
220
221 if ( value < 8.0 && value > -8.0 )
222 {
223 value2 = value*value;
224
225 fact1 = 57568490574.0 + value2*(-13362590354.0
226 + value2*( 651619640.7
227 + value2*(-11214424.18
228 + value2*( 77392.33017
229 + value2*(-184.9052456 ) ) ) ) );
230
231 fact2 = 57568490411.0 + value2*( 1029532985.0
232 + value2*( 9494680.718
233 + value2*(59272.64853
234 + value2*(267.8532712
235 + value2*1.0 ) ) ) );
236
237 bessel = fact1/fact2;
238 }
239 else
240 {
241 arg = 8.0/modvalue;
242
243 value2 = arg*arg;
244
245 shift = modvalue-0.785398164;
246
247 fact1 = 1.0 + value2*(-0.1098628627e-2
248 + value2*(0.2734510407e-4
249 + value2*(-0.2073370639e-5
250 + value2*0.2093887211e-6 ) ) );
251
252 fact2 = -0.1562499995e-1 + value2*(0.1430488765e-3
253 + value2*(-0.6911147651e-5
254 + value2*(0.7621095161e-6
255 - value2*0.934945152e-7 ) ) );
256
257 bessel = std::sqrt(0.636619772/modvalue)*(std::cos(shift)*fact1 - arg*std::sin(shift)*fact2 );
258 }
259 return bessel;
260}

Referenced by GetDiffElasticSumProbA().

◆ BesselOneByArg()

G4double G4DiffuseElasticV2::BesselOneByArg ( G4double  z)
inline

Definition at line 342 of file G4DiffuseElasticV2.hh.

343{
344 G4double x2, result;
345
346 if( std::fabs(x) < 0.01 )
347 {
348 x *= 0.5;
349 x2 = x*x;
350 result = 2. - x2 + x2*x2/6.;
351 }
352 else
353 {
354 result = BesselJone(x)/x;
355 }
356 return result;
357}
G4double BesselJone(G4double z)

References BesselJone().

Referenced by GetDiffElasticSumProbA().

◆ Block()

void G4HadronicInteraction::Block ( )
inlineprotectedinherited

◆ BuildAngleTable()

void G4DiffuseElasticV2::BuildAngleTable ( )

Definition at line 435 of file G4DiffuseElasticV2.cc.

436{
437 G4int i, j;
438 G4double partMom, kinE, a = 0., z = fParticle->GetPDGCharge(), m1 = fParticle->GetPDGMass();
439 G4double alpha1, alpha2, alphaMax, alphaCoulomb, delta = 0., sum = 0.;
440
442
443 fEnergyAngleVector = new std::vector<std::vector<double>*>;
444 fEnergySumVector = new std::vector<std::vector<double>*>;
445
446 for( i = 0; i < fEnergyBin; i++)
447 {
448 kinE = fEnergyVector->Energy(i);
449 partMom = std::sqrt( kinE*(kinE + 2*m1) );
450
451 fWaveVector = partMom/hbarc;
452
454 G4double kRmax = 18.6; // 10.6; 10.6, 18, 10.174; ~ 3 maxima of J1 or 15., 25.
455 G4double kRcoul = 1.9; // 1.2; 1.4, 2.5; // on the first slope of J1
456
457 alphaMax = kRmax/kR;
458
459 if ( alphaMax >= CLHEP::pi ) alphaMax = CLHEP::pi; // vmg21.10.15
460
461 alphaCoulomb = kRcoul/kR;
462
463 if( z )
464 {
465 a = partMom/m1; // beta*gamma for m1
466 fBeta = a/std::sqrt(1+a*a);
469 fAddCoulomb = true;
470 }
471
472 std::vector<double>* angleVector = new std::vector<double>(fAngleBin);
473 std::vector<double>* sumVector = new std::vector<double>(fAngleBin);
474
475
476 G4double delth = alphaMax/fAngleBin;
477
478 sum = 0.;
479
480 for(j = fAngleBin-1; j >= 0; j--)
481 {
482 alpha1 = delth*j;
483 alpha2 = alpha1 + delth;
484
485 if( fAddCoulomb && ( alpha2 < alphaCoulomb)) fAddCoulomb = false;
486
487 delta = integral.Legendre10(this, &G4DiffuseElasticV2::GetIntegrandFunction, alpha1, alpha2);
488
489 sum += delta;
490
491 (*angleVector)[j] = alpha1;
492 (*sumVector)[j] = sum;
493
494 }
495 fEnergyAngleVector->push_back(angleVector);
496 fEnergySumVector->push_back(sumVector);
497
498 }
499 return;
500}
const G4double alpha2
G4double CalculateZommerfeld(G4double beta, G4double Z1, G4double Z2)
G4double CalculateAm(G4double momentum, G4double n, G4double Z)
G4double GetIntegrandFunction(G4double theta)
G4double GetPDGCharge() const
G4double Energy(const std::size_t index) const
static constexpr double pi
Definition: SystemOfUnits.h:55
float hbarc
Definition: hepunit.py:264

References alpha2, CalculateAm(), CalculateZommerfeld(), G4PhysicsVector::Energy(), fAddCoulomb, fAm, fAngleBin, fAtomicNumber, fBeta, fEnergyAngleVector, fEnergyBin, fEnergySumVector, fEnergyVector, fNuclearRadius, fParticle, fWaveVector, fZommerfeld, G4DiffuseElasticV2(), GetIntegrandFunction(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), source.hepunit::hbarc, and CLHEP::pi.

Referenced by Initialise(), and InitialiseOnFly().

◆ BuildPhysicsTable()

void G4HadronicInteraction::BuildPhysicsTable ( const G4ParticleDefinition )
virtualinherited

◆ CalculateAm()

G4double G4DiffuseElasticV2::CalculateAm ( G4double  momentum,
G4double  n,
G4double  Z 
)
inline

Definition at line 375 of file G4DiffuseElasticV2.hh.

376{
377 G4double k = momentum/CLHEP::hbarc;
378 G4double ch = 1.13 + 3.76*n*n;
379 G4double zn = 1.77*k*(1.0/G4Pow::GetInstance()->A13(Z))*CLHEP::Bohr_radius;
380 G4double zn2 = zn*zn;
381 fAm = ch/zn2;
382
383 return fAm;
384}
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double A13(G4double A) const
Definition: G4Pow.cc:120
static constexpr double Bohr_radius
static constexpr double hbarc

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

Referenced by BuildAngleTable().

◆ CalculateNuclearRad()

G4double G4DiffuseElasticV2::CalculateNuclearRad ( G4double  A)
inline

Definition at line 390 of file G4DiffuseElasticV2.hh.

391{
392 G4double R, r0, a11, a12, a13, a2, a3;
393
394 a11 = 1.26; // 1.08, 1.16
395 a12 = 1.; // 1.08, 1.16
396 a13 = 1.12; // 1.08, 1.16
397 a2 = 1.1;
398 a3 = 1.;
399
400 // Special rms radii for light nucleii
401
402 if (A < 50.)
403 {
404 if (std::abs(A-1.) < 0.5) return 0.89*CLHEP::fermi; // p
405 else if(std::abs(A-2.) < 0.5) return 2.13*CLHEP::fermi; // d
406 else if( // std::abs(Z-1.) < 0.5 &&
407 std::abs(A-3.) < 0.5) return 1.80*CLHEP::fermi; // t
408
409 // else if(std::abs(Z-2.) < 0.5 && std::abs(A-3.) < 0.5) return 1.96CLHEP::fermi; // He3
410 else if( // std::abs(Z-2.) < 0.5 &&
411 std::abs(A-4.) < 0.5) return 1.68*CLHEP::fermi; // He4
412
413 else if( // std::abs(Z-3.) < 0.5
414 std::abs(A-7.) < 0.5 ) return 2.40*CLHEP::fermi; // Li7
415 else if( // std::abs(Z-4.) < 0.5
416 std::abs(A-9.) < 0.5) return 2.51*CLHEP::fermi; // Be9
417
418 else if( 10. < A && A <= 16. ) r0 = a11*( 1 - (1.0/G4Pow::GetInstance()->A23(A)) )*CLHEP::fermi; // 1.08CLHEP::fermi;
419 else if( 15. < A && A <= 20. ) r0 = a12*( 1 - (1.0/G4Pow::GetInstance()->A23(A)) )*CLHEP::fermi;
420 else if( 20. < A && A <= 30. ) r0 = a13*( 1 - (1.0/G4Pow::GetInstance()->A23(A)) )*CLHEP::fermi;
421 else r0 = a2*CLHEP::fermi;
422
423 R = r0*G4Pow::GetInstance()->A13(A);
424 }
425 else
426 {
427 r0 = a3*CLHEP::fermi;
428
429 R = r0*G4Pow::GetInstance()->powA(A, 0.27);
430 }
431 fNuclearRadius = R;
432
433 return R;
434}
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:230
G4double A23(G4double A) const
Definition: G4Pow.hh:131
static constexpr double fermi
Definition: SystemOfUnits.h:84

References A, G4Pow::A13(), G4Pow::A23(), CLHEP::fermi, fNuclearRadius, G4Pow::GetInstance(), and G4Pow::powA().

Referenced by Initialise(), and InitialiseOnFly().

◆ CalculateZommerfeld()

G4double G4DiffuseElasticV2::CalculateZommerfeld ( G4double  beta,
G4double  Z1,
G4double  Z2 
)
inline

◆ 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 G4DiffuseElasticV2::DampFactor ( G4double  z)
inline

Definition at line 319 of file G4DiffuseElasticV2.hh.

320{
321 G4double df;
322 G4double f2 = 2., f3 = 6., f4 = 24.; // first factorials
323
324 // x *= pi;
325
326 if( std::fabs(x) < 0.01 )
327 {
328 df = 1./(1. + x/f2 + x*x/f3 + x*x*x/f4);
329 }
330 else
331 {
332 df = x/std::sinh(x);
333 }
334 return df;
335}

Referenced by GetDiffElasticSumProbA().

◆ 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().

◆ GetDiffElasticSumProbA()

G4double G4DiffuseElasticV2::GetDiffElasticSumProbA ( G4double  alpha)

Definition at line 164 of file G4DiffuseElasticV2.cc.

165{
166
167 G4double sigma, bzero, bzero2, bonebyarg, bonebyarg2, damp, damp2;
168 G4double delta, diffuse, gamma;
169 G4double e1, e2, bone, bone2;
170
171 // G4double wavek = momentum/hbarc; // wave vector
172 // G4double r0 = 1.08*fermi;
173 // G4double rad = r0*G4Pow::GetInstance()->A13(A);
174
175 G4double kr = fWaveVector*fNuclearRadius; // wavek*rad;
176 G4double kr2 = kr*kr;
177 G4double krt = kr*theta;
178
179 bzero = BesselJzero(krt);
180 bzero2 = bzero*bzero;
181 bone = BesselJone(krt);
182 bone2 = bone*bone;
183 bonebyarg = BesselOneByArg(krt);
184 bonebyarg2 = bonebyarg*bonebyarg;
185
186 if ( fParticle == theProton )
187 {
188 diffuse = 0.63*fermi;
189 gamma = 0.3*fermi;
190 delta = 0.1*fermi*fermi;
191 e1 = 0.3*fermi;
192 e2 = 0.35*fermi;
193 }
194 else if ( fParticle == theNeutron )
195 {
196 diffuse = 0.63*fermi;
197 gamma = 0.3*fermi;
198 delta = 0.1*fermi*fermi;
199 e1 = 0.3*fermi;
200 e2 = 0.35*fermi;
201 }
202 else // as proton, if were not defined
203 {
204 diffuse = 0.63*fermi;
205 gamma = 0.3*fermi;
206 delta = 0.1*fermi*fermi;
207 e1 = 0.3*fermi;
208 e2 = 0.35*fermi;
209 }
210
211 G4double lambda = 15; // 15 ok
212 // G4double kgamma = fWaveVector*gamma; // wavek*delta;
213 G4double kgamma = lambda*(1.-G4Exp(-fWaveVector*gamma/lambda)); // wavek*delta;
214
215 if( fAddCoulomb ) // add Coulomb correction
216 {
217 G4double sinHalfTheta = std::sin(0.5*theta);
218 G4double sinHalfTheta2 = sinHalfTheta*sinHalfTheta;
219
220 kgamma += 0.5*fZommerfeld/kr/(sinHalfTheta2+fAm); // correction at J0()
221 }
222
223 G4double kgamma2 = kgamma*kgamma;
224
225 // G4double dk2t = delta*fWaveVector*fWaveVector*theta; // delta*wavek*wavek*theta;
226 // G4double dk2t2 = dk2t*dk2t;
227
228 // G4double pikdt = pi*fWaveVector*diffuse*theta;// pi*wavek*diffuse*theta;
229 G4double pikdt = lambda*(1. - G4Exp( -pi*fWaveVector*diffuse*theta/lambda ) ); // wavek*delta;
230
231 damp = DampFactor( pikdt );
232 damp2 = damp*damp;
233
234 G4double mode2k2 = ( e1*e1 + e2*e2 )*fWaveVector*fWaveVector;
235 G4double e2dk3t = -2.*e2*delta*fWaveVector*fWaveVector*fWaveVector*theta;
236
237 sigma = kgamma2;
238 // sigma += dk2t2;
239 sigma *= bzero2;
240 sigma += mode2k2*bone2;
241 sigma += e2dk3t*bzero*bone;
242
243 // sigma += kr2*(1 + 8.*fZommerfeld*fZommerfeld/kr2)*bonebyarg2; // correction at J1()/()
244 sigma += kr2*bonebyarg2; // correction at J1()/()
245
246 sigma *= damp2; // *rad*rad;
247
248 return sigma;
249}
static const G4double e2[44]
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
static constexpr double fermi
Definition: G4SIunits.hh:83
static constexpr double pi
Definition: G4SIunits.hh:55
G4double BesselJzero(G4double z)
G4double DampFactor(G4double z)
G4double BesselOneByArg(G4double z)

References BesselJone(), BesselJzero(), BesselOneByArg(), DampFactor(), e1, e2, fAddCoulomb, fAm, fermi, fNuclearRadius, fParticle, fWaveVector, fZommerfeld, G4Exp(), G4InuclParticleNames::lambda, pi, theNeutron, and theProton.

Referenced by GetIntegrandFunction().

◆ 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().

◆ GetIntegrandFunction()

G4double G4DiffuseElasticV2::GetIntegrandFunction ( G4double  theta)

Definition at line 257 of file G4DiffuseElasticV2.cc.

258{
259 G4double result;
260
261 result = GetDiffElasticSumProbA(alpha) * 2 * CLHEP::pi * std::sin(alpha);
262
263 return result;
264}
static const G4double alpha
G4double GetDiffElasticSumProbA(G4double alpha)

References alpha, GetDiffElasticSumProbA(), and CLHEP::pi.

Referenced by BuildAngleTable().

◆ 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().

◆ GetNuclearRadius()

G4double G4DiffuseElasticV2::GetNuclearRadius ( )
inline

Definition at line 129 of file G4DiffuseElasticV2.hh.

129{return fNuclearRadius;};

References fNuclearRadius.

◆ GetRecoilEnergyThreshold()

G4double G4HadronicInteraction::GetRecoilEnergyThreshold ( ) const
inlineinherited

◆ GetScatteringAngle()

G4double G4DiffuseElasticV2::GetScatteringAngle ( G4int  iMomentum,
unsigned long  iAngle,
G4double  position 
)

Definition at line 507 of file G4DiffuseElasticV2.cc.

508{
509 G4double x1, x2, y1, y2, randAngle = 0;
510
511 if( iAngle == 0 )
512 {
513 randAngle = (*(*fEnergyAngleVector)[iMomentum])[iAngle];
514 }
515 else
516 {
517 if ( iAngle >= (*fEnergyAngleVector)[iMomentum]->size() )
518 {
519 iAngle = (*fEnergyAngleVector)[iMomentum]->size() - 1;
520 }
521
522 y1 = (*(*fEnergySumVector)[iMomentum])[iAngle-1];
523 y2 = (*(*fEnergySumVector)[iMomentum])[iAngle];
524
525 x1 = (*(*fEnergyAngleVector)[iMomentum])[iAngle-1];
526 x2 = (*(*fEnergyAngleVector)[iMomentum])[iAngle];
527
528 if ( x1 == x2 ) randAngle = x2;
529 else
530 {
531 if ( y1 == y2 ) randAngle = x1 + ( x2 - x1 )*G4UniformRand();
532 else
533 {
534 randAngle = x1 + ( position - y1 )*( x2 - x1 )/( y2 - y1 );
535 }
536 }
537 }
538
539 return randAngle;
540}

References fEnergyAngleVector, and G4UniformRand.

Referenced by SampleTableThetaCMS().

◆ 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.

◆ Initialise()

void G4DiffuseElasticV2::Initialise ( )

Definition at line 128 of file G4DiffuseElasticV2.cc.

129{
130
131 const G4ElementTable* theElementTable = G4Element::GetElementTable();
132
133 size_t jEl, numOfEl = G4Element::GetNumberOfElements();
134
135 for( jEl = 0; jEl < numOfEl; ++jEl) // application element loop
136 {
137 fAtomicNumber = (*theElementTable)[jEl]->GetZ(); // atomic number
140
141 if( verboseLevel > 0 )
142 {
143 G4cout<<"G4DiffuseElasticV2::Initialise() the element: "
144 <<(*theElementTable)[jEl]->GetName()<<G4endl;
145 }
147 fElementNameVector.push_back((*theElementTable)[jEl]->GetName());
148
150
153
154 }
155 return;
156}
std::vector< G4Element * > G4ElementTable
std::vector< G4String > fElementNameVector
std::vector< std::vector< std::vector< double > * > * > fEnergyAngleVectorBank
G4double CalculateNuclearRad(G4double A)
std::vector< std::vector< std::vector< double > * > * > fEnergySumVectorBank
std::vector< G4double > fElementNumberVector
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:397
static size_t GetNumberOfElements()
Definition: G4Element.cc:404
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const

References BuildAngleTable(), CalculateNuclearRad(), fAtomicNumber, fAtomicWeight, fElementNameVector, fElementNumberVector, fEnergyAngleVector, fEnergyAngleVectorBank, fEnergySumVector, fEnergySumVectorBank, fNuclearRadius, G4cout, G4endl, G4NistManager::GetAtomicMassAmu(), G4Element::GetElementTable(), G4Element::GetNumberOfElements(), G4NistManager::Instance(), and G4HadronicInteraction::verboseLevel.

◆ InitialiseModel()

void G4HadronicInteraction::InitialiseModel ( )
virtualinherited

◆ InitialiseOnFly()

void G4DiffuseElasticV2::InitialiseOnFly ( G4double  Z,
G4double  A 
)

◆ IsApplicable()

G4bool G4DiffuseElasticV2::IsApplicable ( const G4HadProjectile projectile,
G4Nucleus nucleus 
)
inlinevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 170 of file G4DiffuseElasticV2.hh.

172{
173 if( ( projectile.GetDefinition() == G4Proton::Proton() ||
174 projectile.GetDefinition() == G4Neutron::Neutron() ||
175 projectile.GetDefinition() == G4PionPlus::PionPlus() ||
176 projectile.GetDefinition() == G4PionMinus::PionMinus() ||
177 projectile.GetDefinition() == G4KaonPlus::KaonPlus() ||
178 projectile.GetDefinition() == G4KaonMinus::KaonMinus() ) &&
179
180 nucleus.GetZ_asInt() >= 2 ) return true;
181 else return false;
182}
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:112
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:112
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:97
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:97

References G4HadProjectile::GetDefinition(), G4Nucleus::GetZ_asInt(), G4KaonMinus::KaonMinus(), G4KaonPlus::KaonPlus(), G4Neutron::Neutron(), G4PionMinus::PionMinus(), G4PionPlus::PionPlus(), and G4Proton::Proton().

◆ 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}

◆ NeutronTuniform()

G4double G4DiffuseElasticV2::NeutronTuniform ( G4int  Z)

Definition at line 310 of file G4DiffuseElasticV2.cc.

311{
312 G4double elZ = G4double(Z);
313 elZ -= 1.;
314 G4double Tkin = 12.*G4Exp(-elZ/10.) + 1.;
315
316 return Tkin;
317}

References G4Exp(), and Z.

Referenced by SampleInvariantT().

◆ operator!=()

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

◆ operator==()

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

◆ SampleInvariantT()

G4double G4DiffuseElasticV2::SampleInvariantT ( const G4ParticleDefinition p,
G4double  plab,
G4int  Z,
G4int  A 
)
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 273 of file G4DiffuseElasticV2.cc.

275{
276 fParticle = aParticle;
277 G4double m1 = fParticle->GetPDGMass(), t;
278 G4double totElab = std::sqrt(m1*m1+p*p);
280 G4LorentzVector lv1(p,0.0,0.0,totElab);
281 G4LorentzVector lv(0.0,0.0,0.0,mass2);
282 lv += lv1;
283
284 G4ThreeVector bst = lv.boostVector();
285 lv1.boost(-bst);
286
287 G4ThreeVector p1 = lv1.vect();
288 G4double momentumCMS = p1.mag();
289
290 if( aParticle == theNeutron)
291 {
292 G4double Tmax = NeutronTuniform( Z );
293 G4double pCMS2 = momentumCMS*momentumCMS;
294 G4double Tkin = std::sqrt(pCMS2+m1*m1)-m1;
295
296 if( Tkin <= Tmax )
297 {
298 t = 4.*pCMS2*G4UniformRand();
299 return t;
300 }
301 }
302
303 t = SampleTableT( aParticle, momentumCMS, G4double(Z), G4double(A) ); // sample theta in cms
304
305 return t;
306}
double mag() const
G4double SampleTableT(const G4ParticleDefinition *aParticle, G4double p, G4double Z, G4double A)
G4double NeutronTuniform(G4int Z)

References A, CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::boostVector(), fParticle, G4UniformRand, G4NucleiProperties::GetNuclearMass(), G4ParticleDefinition::GetPDGMass(), CLHEP::Hep3Vector::mag(), NeutronTuniform(), SampleTableT(), theNeutron, CLHEP::HepLorentzVector::vect(), and Z.

◆ SampleTableT()

G4double G4DiffuseElasticV2::SampleTableT ( const G4ParticleDefinition aParticle,
G4double  p,
G4double  Z,
G4double  A 
)

Definition at line 324 of file G4DiffuseElasticV2.cc.

326{
327 G4double alpha = SampleTableThetaCMS( aParticle, p, Z, A); // sample theta in cms
328 G4double t = 2*p*p*( 1 - std::cos(alpha) ); // -t !!!
329
330 return t;
331}
G4double SampleTableThetaCMS(const G4ParticleDefinition *aParticle, G4double p, G4double Z, G4double A)

References A, alpha, SampleTableThetaCMS(), and Z.

Referenced by SampleInvariantT().

◆ SampleTableThetaCMS()

G4double G4DiffuseElasticV2::SampleTableThetaCMS ( const G4ParticleDefinition aParticle,
G4double  p,
G4double  Z,
G4double  A 
)

Definition at line 339 of file G4DiffuseElasticV2.cc.

341{
342 size_t iElement;
343 G4int iMomentum;
344 unsigned long iAngle = 0;
345 G4double randAngle, position, theta1, theta2, E1, E2, W1, W2, W;
346 G4double m1 = particle->GetPDGMass();
347
348 for(iElement = 0; iElement < fElementNumberVector.size(); iElement++)
349 {
350 if( std::fabs(Z - fElementNumberVector[iElement]) < 0.5) break;
351 }
352
353 if ( iElement == fElementNumberVector.size() )
354 {
355 InitialiseOnFly(Z,A); // table preparation, if needed
356 }
357
360
361
362 G4double kinE = std::sqrt(momentum*momentum + m1*m1) - m1;
363
364 iMomentum = fEnergyVector->FindBin(kinE,1000) + 1;
365
366 position = (*(*fEnergySumVector)[iMomentum])[0]*G4UniformRand();
367
368 for(iAngle = 0; iAngle < fAngleBin; iAngle++)
369 {
370 if (position > (*(*fEnergySumVector)[iMomentum])[iAngle]) break;
371 }
372
373
374 if (iMomentum == fEnergyBin -1 || iMomentum == 0 ) // the table edges
375 {
376 randAngle = GetScatteringAngle(iMomentum, iAngle, position);
377 }
378 else // kinE inside between energy table edges
379 {
380 theta2 = GetScatteringAngle(iMomentum, iAngle, position);
381
382 E2 = fEnergyVector->Energy(iMomentum);
383
384 iMomentum--;
385
386 theta1 = GetScatteringAngle(iMomentum, iAngle, position);
387
388 E1 = fEnergyVector->Energy(iMomentum);
389
390 W = 1.0/(E2 - E1);
391 W1 = (E2 - kinE)*W;
392 W2 = (kinE - E1)*W;
393
394 randAngle = W1*theta1 + W2*theta2;
395 }
396
397
398
399 if(randAngle < 0.) randAngle = 0.;
400
401 return randAngle;
402}
G4double GetScatteringAngle(G4int iMomentum, unsigned long iAngle, G4double position)
void InitialiseOnFly(G4double Z, G4double A)
std::size_t FindBin(const G4double energy, std::size_t idx) const
#define position
Definition: xmlparse.cc:622

References A, G4PhysicsVector::Energy(), fAngleBin, fElementNumberVector, fEnergyAngleVector, fEnergyAngleVectorBank, fEnergyBin, fEnergySumVector, fEnergySumVectorBank, fEnergyVector, G4PhysicsVector::FindBin(), G4UniformRand, G4ParticleDefinition::GetPDGMass(), GetScatteringAngle(), InitialiseOnFly(), position, and Z.

Referenced by SampleTableT().

◆ SampleThetaCMS()

G4double G4DiffuseElasticV2::SampleThetaCMS ( const G4ParticleDefinition aParticle,
G4double  p,
G4double  A 
)

◆ SampleThetaLab()

G4double G4DiffuseElasticV2::SampleThetaLab ( const G4HadProjectile aParticle,
G4double  tmass,
G4double  A 
)

◆ 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().

◆ SetHEModelLowLimit()

void G4DiffuseElasticV2::SetHEModelLowLimit ( G4double  value)
inline

Definition at line 194 of file G4DiffuseElasticV2.hh.

195{
196 lowEnergyLimitHE = value;
197}

References lowEnergyLimitHE.

◆ SetLowestEnergyLimit()

void G4DiffuseElasticV2::SetLowestEnergyLimit ( G4double  value)
inline

Definition at line 204 of file G4DiffuseElasticV2.hh.

205{
206 lowestEnergyLimit = value;
207}

References lowestEnergyLimit.

◆ 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(), 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(), 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

◆ SetPlabLowLimit()

void G4DiffuseElasticV2::SetPlabLowLimit ( G4double  value)
inline

Definition at line 189 of file G4DiffuseElasticV2.hh.

190{
191 plabLowLimit = value;
192}

References plabLowLimit.

◆ SetQModelLowLimit()

void G4DiffuseElasticV2::SetQModelLowLimit ( G4double  value)
inline

Definition at line 199 of file G4DiffuseElasticV2.hh.

200{
201 lowEnergyLimitQ = value;
202}

References lowEnergyLimitQ.

◆ SetRecoilEnergyThreshold()

void G4HadronicInteraction::SetRecoilEnergyThreshold ( G4double  val)
inlineinherited

◆ SetRecoilKinEnergyLimit()

void G4DiffuseElasticV2::SetRecoilKinEnergyLimit ( G4double  value)
inline

Definition at line 184 of file G4DiffuseElasticV2.hh.

185{
186 lowEnergyRecoilLimit = value;
187}

References lowEnergyRecoilLimit.

◆ SetVerboseLevel()

void G4HadronicInteraction::SetVerboseLevel ( G4int  value)
inlineinherited

◆ ThetaCMStoThetaLab()

G4double G4DiffuseElasticV2::ThetaCMStoThetaLab ( const G4DynamicParticle aParticle,
G4double  tmass,
G4double  thetaCMS 
)

Definition at line 552 of file G4DiffuseElasticV2.cc.

554{
555 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
556 G4double m1 = theParticle->GetPDGMass();
557 G4LorentzVector lv1 = aParticle->Get4Momentum();
558 G4LorentzVector lv(0.0,0.0,0.0,tmass);
559
560 lv += lv1;
561
562 G4ThreeVector bst = lv.boostVector();
563
564 lv1.boost(-bst);
565
566 G4ThreeVector p1 = lv1.vect();
567 G4double ptot = p1.mag();
568
570 G4double cost = std::cos(thetaCMS);
571 G4double sint;
572
573 if( cost >= 1.0 )
574 {
575 cost = 1.0;
576 sint = 0.0;
577 }
578 else if( cost <= -1.0)
579 {
580 cost = -1.0;
581 sint = 0.0;
582 }
583 else
584 {
585 sint = std::sqrt((1.0-cost)*(1.0+cost));
586 }
587 if (verboseLevel>1)
588 {
589 G4cout << "cos(tcms)=" << cost << " std::sin(tcms)=" << sint << G4endl;
590 }
591 G4ThreeVector v1(sint*std::cos(phi),sint*std::sin(phi),cost);
592 v1 *= ptot;
593 G4LorentzVector nlv1(v1.x(),v1.y(),v1.z(),std::sqrt(ptot*ptot + m1*m1));
594
595 nlv1.boost(bst);
596
597 G4ThreeVector np1 = nlv1.vect();
598
599 G4double thetaLab = np1.theta();
600
601 return thetaLab;
602}
static constexpr double twopi
Definition: G4SIunits.hh:56
double theta() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
G4ParticleDefinition * GetDefinition() const
G4LorentzVector Get4Momentum() const

References CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::boostVector(), G4cout, G4endl, G4UniformRand, G4DynamicParticle::Get4Momentum(), G4DynamicParticle::GetDefinition(), G4ParticleDefinition::GetPDGMass(), CLHEP::Hep3Vector::mag(), CLHEP::Hep3Vector::theta(), twopi, CLHEP::HepLorentzVector::vect(), G4HadronicInteraction::verboseLevel, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ ThetaLabToThetaCMS()

G4double G4DiffuseElasticV2::ThetaLabToThetaCMS ( const G4DynamicParticle aParticle,
G4double  tmass,
G4double  thetaLab 
)

Definition at line 610 of file G4DiffuseElasticV2.cc.

612{
613 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
614 G4double m1 = theParticle->GetPDGMass();
615 G4double plab = aParticle->GetTotalMomentum();
616 G4LorentzVector lv1 = aParticle->Get4Momentum();
617 G4LorentzVector lv(0.0,0.0,0.0,tmass);
618
619 lv += lv1;
620
621 G4ThreeVector bst = lv.boostVector();
622
624 G4double cost = std::cos(thetaLab);
625 G4double sint;
626
627 if( cost >= 1.0 )
628 {
629 cost = 1.0;
630 sint = 0.0;
631 }
632 else if( cost <= -1.0)
633 {
634 cost = -1.0;
635 sint = 0.0;
636 }
637 else
638 {
639 sint = std::sqrt((1.0-cost)*(1.0+cost));
640 }
641 if (verboseLevel>1)
642 {
643 G4cout << "cos(tlab)=" << cost << " std::sin(tlab)=" << sint << G4endl;
644 }
645 G4ThreeVector v1(sint*std::cos(phi),sint*std::sin(phi),cost);
646 v1 *= plab;
647 G4LorentzVector nlv1(v1.x(),v1.y(),v1.z(),std::sqrt(plab*plab + m1*m1));
648
649 nlv1.boost(-bst);
650
651 G4ThreeVector np1 = nlv1.vect();
652 G4double thetaCMS = np1.theta();
653
654 return thetaCMS;
655}
G4double GetTotalMomentum() const

References CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::boostVector(), G4cout, G4endl, G4UniformRand, G4DynamicParticle::Get4Momentum(), G4DynamicParticle::GetDefinition(), G4ParticleDefinition::GetPDGMass(), G4DynamicParticle::GetTotalMomentum(), CLHEP::Hep3Vector::theta(), twopi, CLHEP::HepLorentzVector::vect(), G4HadronicInteraction::verboseLevel, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Field Documentation

◆ epCheckLevels

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

◆ fAddCoulomb

G4bool G4DiffuseElasticV2::fAddCoulomb
private

◆ fAm

G4double G4DiffuseElasticV2::fAm
private

◆ fAngleBin

unsigned long G4DiffuseElasticV2::fAngleBin
private

Definition at line 144 of file G4DiffuseElasticV2.hh.

Referenced by BuildAngleTable(), G4DiffuseElasticV2(), and SampleTableThetaCMS().

◆ fAtomicNumber

G4double G4DiffuseElasticV2::fAtomicNumber
private

◆ fAtomicWeight

G4double G4DiffuseElasticV2::fAtomicWeight
private

Definition at line 160 of file G4DiffuseElasticV2.hh.

Referenced by G4DiffuseElasticV2(), Initialise(), and InitialiseOnFly().

◆ fBeta

G4double G4DiffuseElasticV2::fBeta
private

Definition at line 163 of file G4DiffuseElasticV2.hh.

Referenced by BuildAngleTable(), and G4DiffuseElasticV2().

◆ fElementNameVector

std::vector<G4String> G4DiffuseElasticV2::fElementNameVector
private

Definition at line 156 of file G4DiffuseElasticV2.hh.

Referenced by Initialise().

◆ fElementNumberVector

std::vector<G4double> G4DiffuseElasticV2::fElementNumberVector
private

Definition at line 155 of file G4DiffuseElasticV2.hh.

Referenced by Initialise(), InitialiseOnFly(), and SampleTableThetaCMS().

◆ fEnergyAngleVector

std::vector<std::vector<double>*>* G4DiffuseElasticV2::fEnergyAngleVector
private

◆ fEnergyAngleVectorBank

std::vector<std::vector<std::vector<double>*>*> G4DiffuseElasticV2::fEnergyAngleVectorBank
private

Definition at line 148 of file G4DiffuseElasticV2.hh.

Referenced by Initialise(), InitialiseOnFly(), and SampleTableThetaCMS().

◆ fEnergyBin

G4int G4DiffuseElasticV2::fEnergyBin
private

Definition at line 143 of file G4DiffuseElasticV2.hh.

Referenced by BuildAngleTable(), G4DiffuseElasticV2(), and SampleTableThetaCMS().

◆ fEnergySumVector

std::vector<std::vector<double>*>* G4DiffuseElasticV2::fEnergySumVector
private

◆ fEnergySumVectorBank

std::vector<std::vector<std::vector<double>*>*> G4DiffuseElasticV2::fEnergySumVectorBank
private

Definition at line 149 of file G4DiffuseElasticV2.hh.

Referenced by Initialise(), InitialiseOnFly(), and SampleTableThetaCMS().

◆ fEnergyVector

G4PhysicsLogVector* G4DiffuseElasticV2::fEnergyVector
private

◆ fNuclearRadius

G4double G4DiffuseElasticV2::fNuclearRadius
private

◆ fParticle

const G4ParticleDefinition* G4DiffuseElasticV2::fParticle
private

◆ fWaveVector

G4double G4DiffuseElasticV2::fWaveVector
private

◆ fZommerfeld

G4double G4DiffuseElasticV2::fZommerfeld
private

◆ isBlocked

G4bool G4HadronicInteraction::isBlocked
protectedinherited

◆ lowEnergyLimitHE

G4double G4DiffuseElasticV2::lowEnergyLimitHE
private

Definition at line 138 of file G4DiffuseElasticV2.hh.

Referenced by G4DiffuseElasticV2(), and SetHEModelLowLimit().

◆ lowEnergyLimitQ

G4double G4DiffuseElasticV2::lowEnergyLimitQ
private

Definition at line 139 of file G4DiffuseElasticV2.hh.

Referenced by G4DiffuseElasticV2(), and SetQModelLowLimit().

◆ lowEnergyRecoilLimit

G4double G4DiffuseElasticV2::lowEnergyRecoilLimit
private

Definition at line 137 of file G4DiffuseElasticV2.hh.

Referenced by G4DiffuseElasticV2(), and SetRecoilKinEnergyLimit().

◆ lowestEnergyLimit

G4double G4DiffuseElasticV2::lowestEnergyLimit
private

Definition at line 140 of file G4DiffuseElasticV2.hh.

Referenced by G4DiffuseElasticV2(), and SetLowestEnergyLimit().

◆ nwarn

G4int G4HadronElastic::nwarn
privateinherited

◆ plabLowLimit

G4double G4DiffuseElasticV2::plabLowLimit
private

Definition at line 141 of file G4DiffuseElasticV2.hh.

Referenced by G4DiffuseElasticV2(), and SetPlabLowLimit().

◆ pLocalTmax

G4double G4HadronElastic::pLocalTmax
protectedinherited

◆ recoilEnergyThreshold

G4double G4HadronicInteraction::recoilEnergyThreshold
privateinherited

◆ registry

G4HadronicInteractionRegistry* G4HadronicInteraction::registry
privateinherited

◆ secID

G4int G4HadronElastic::secID
protectedinherited

◆ theAlpha

G4ParticleDefinition* G4HadronElastic::theAlpha
privateinherited

◆ theBlockedList

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

◆ theBlockedListElements

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

◆ theDeuteron

G4ParticleDefinition* G4HadronElastic::theDeuteron
privateinherited

◆ 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* G4DiffuseElasticV2::theNeutron
private

◆ 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* G4DiffuseElasticV2::theProton
private

Definition at line 134 of file G4DiffuseElasticV2.hh.

Referenced by G4DiffuseElasticV2(), and GetDiffElasticSumProbA().

◆ 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(), 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(), Initialise(), G4NuclNuclDiffuseElastic::Initialise(), G4DiffuseElastic::InitialiseOnFly(), InitialiseOnFly(), G4NuclNuclDiffuseElastic::InitialiseOnFly(), G4CascadeInterface::makeDynamicParticle(), G4CascadeInterface::NoInteraction(), G4CascadeInterface::Propagate(), G4ElasticHadrNucleusHE::SampleInvariantT(), G4AntiNuclElastic::SampleThetaCMS(), G4DiffuseElastic::SampleThetaLab(), G4NuclNuclDiffuseElastic::SampleThetaLab(), G4AntiNuclElastic::SampleThetaLab(), G4WilsonAbrasionModel::SetUseAblation(), G4HadronicInteraction::SetVerboseLevel(), G4WilsonAbrasionModel::SetVerboseLevel(), G4DiffuseElastic::ThetaCMStoThetaLab(), ThetaCMStoThetaLab(), G4NuclNuclDiffuseElastic::ThetaCMStoThetaLab(), G4DiffuseElastic::ThetaLabToThetaCMS(), ThetaLabToThetaCMS(), and G4NuclNuclDiffuseElastic::ThetaLabToThetaCMS().


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