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

#include <G4TheoFSGenerator.hh>

+ Inheritance diagram for G4TheoFSGenerator:

Public Member Functions

 G4TheoFSGenerator (const G4String &name="TheoFSGenerator")
 
 ~G4TheoFSGenerator () override
 
 G4TheoFSGenerator (const G4TheoFSGenerator &right)=delete
 
const G4TheoFSGeneratoroperator= (const G4TheoFSGenerator &right)=delete
 
G4bool operator== (const G4TheoFSGenerator &right) const =delete
 
G4bool operator!= (const G4TheoFSGenerator &right) const =delete
 
G4HadFinalStateApplyYourself (const G4HadProjectile &thePrimary, G4Nucleus &theNucleus) override
 
void SetTransport (G4VIntraNuclearTransportModel *const value)
 
void SetHighEnergyGenerator (G4VHighEnergyGenerator *const value)
 
void SetQuasiElasticChannel (G4QuasiElasticChannel *const value)
 
const G4VIntraNuclearTransportModelGetTransport () const
 
const G4VHighEnergyGeneratorGetHighEnergyGenerator () const
 
const G4QuasiElasticChannelGetQuasiElasticChannel () const
 
std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const override
 
void ModelDescription (std::ostream &outFile) const override
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void InitialiseModel ()
 
 G4HadronicInteraction (const G4HadronicInteraction &right)=delete
 
const G4HadronicInteractionoperator= (const G4HadronicInteraction &right)=delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 54 of file G4TheoFSGenerator.hh.

Constructor & Destructor Documentation

◆ G4TheoFSGenerator() [1/2]

G4TheoFSGenerator::G4TheoFSGenerator ( const G4String & name = "TheoFSGenerator")
explicit

Definition at line 42 of file G4TheoFSGenerator.cc.

44 , theTransport(nullptr), theHighEnergyGenerator(nullptr)
45 , theQuasielastic(nullptr)
46 , theCosmicCoalescence(nullptr)
47 , theStringModelID(-1)
48{
49 theParticleChange = new G4HadFinalState;
50 theStringModelID = G4PhysicsModelCatalog::GetModelID( "model_" + name );
51}
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static G4int GetModelID(const G4int modelIndex)

◆ ~G4TheoFSGenerator()

G4TheoFSGenerator::~G4TheoFSGenerator ( )
override

Definition at line 53 of file G4TheoFSGenerator.cc.

54{
55 delete theParticleChange;
56}

◆ G4TheoFSGenerator() [2/2]

G4TheoFSGenerator::G4TheoFSGenerator ( const G4TheoFSGenerator & right)
delete

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4TheoFSGenerator::ApplyYourself ( const G4HadProjectile & thePrimary,
G4Nucleus & theNucleus )
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 72 of file G4TheoFSGenerator.cc.

73{
74 // init particle change
75 theParticleChange->Clear();
76 theParticleChange->SetStatusChange(stopAndKill);
77 G4double timePrimary=thePrimary.GetGlobalTime();
78
79 // Temporarily dummy treatment of heavy (charm and bottom) hadron projectiles at low energies.
80 // Cascade models are currently not applicable for heavy hadrons and string models cannot
81 // handle them properly at low energies - let's say safely below ~100 MeV.
82 // In these cases, we return as final state the initial state unchanged.
83 // For most applications, this is a safe simplification, giving that the nearly all
84 // slowly moving charm and bottom hadrons decay before any hadronic interaction can occur.
85 // Note that we prefer not to use G4HadronicParameters::GetMinEnergyTransitionFTF_Cascade()
86 // (typically ~3 GeV) because FTFP works reasonably well below such a value.
87 const G4double energyThresholdForCharmAndBottomHadrons = 100.0*CLHEP::MeV;
88 if ( thePrimary.GetKineticEnergy() < energyThresholdForCharmAndBottomHadrons &&
89 ( thePrimary.GetDefinition()->GetQuarkContent( 4 ) != 0 || // Has charm constituent quark
90 thePrimary.GetDefinition()->GetAntiQuarkContent( 4 ) != 0 || // Has anti-charm constituent anti-quark
91 thePrimary.GetDefinition()->GetQuarkContent( 5 ) != 0 || // Has bottom constituent quark
92 thePrimary.GetDefinition()->GetAntiQuarkContent( 5 ) != 0 ) ) { // Has anti-bottom constituent anti-quark
93 theParticleChange->SetStatusChange( isAlive );
94 theParticleChange->SetEnergyChange( thePrimary.GetKineticEnergy() );
95 theParticleChange->SetMomentumChange( thePrimary.Get4Momentum().vect().unit() );
96 return theParticleChange;
97 }
98
99 // Temporarily dummy treatment of light hypernuclei projectiles at low energies.
100 // Bertini and Binary intra-nuclear cascade models are currently not applicable
101 // for hypernuclei and string models cannot handle them properly at low energies -
102 // let's say safely below ~100 MeV.
103 // In these cases, we return as final state the initial state unchanged.
104 // For most applications, this is a safe simplification, giving that the nearly all
105 // slowly moving hypernuclei decay before any hadronic interaction can occur.
106 // Note that we prefer not to use G4HadronicParameters::GetMinEnergyTransitionFTF_Cascade()
107 // (typically ~3 GeV) because FTFP works reasonably well below such a value.
108 // Note that for anti-hypernuclei, FTF model works fine down to zero kinetic energy,
109 // so there is no need of any extra, dummy treatment.
110 const G4double energyThresholdForHyperNuclei = 100.0*CLHEP::MeV;
111 if ( thePrimary.GetKineticEnergy() < energyThresholdForHyperNuclei &&
112 thePrimary.GetDefinition()->IsHypernucleus() ) {
113 theParticleChange->SetStatusChange( isAlive );
114 theParticleChange->SetEnergyChange( thePrimary.GetKineticEnergy() );
115 theParticleChange->SetMomentumChange( thePrimary.Get4Momentum().vect().unit() );
116 return theParticleChange;
117 }
118
119 // check if models have been registered, and use default, in case this is not true @@
120
121 const G4DynamicParticle aPart(thePrimary.GetDefinition(),thePrimary.Get4Momentum().vect());
122
123 if ( theQuasielastic )
124 {
125 if ( theQuasielastic->GetFraction(theNucleus, aPart) > G4UniformRand() )
126 {
127 //G4cout<<"___G4TheoFSGenerator: before Scatter (1) QE=" << theQuasielastic<<G4endl;
128 G4KineticTrackVector *result= theQuasielastic->Scatter(theNucleus, aPart);
129 //G4cout << "^^G4TheoFSGenerator: after Scatter (1) " << G4endl;
130 if (result)
131 {
132 for(auto & ptr : *result)
133 {
134 G4DynamicParticle * aNew =
135 new G4DynamicParticle(ptr->GetDefinition(),
136 ptr->Get4Momentum().e(),
137 ptr->Get4Momentum().vect());
138 theParticleChange->AddSecondary(aNew, ptr->GetCreatorModelID());
139 delete ptr;
140 }
141 delete result;
142 }
143 else
144 {
145 theParticleChange->SetStatusChange(isAlive);
146 theParticleChange->SetEnergyChange(thePrimary.GetKineticEnergy());
147 theParticleChange->SetMomentumChange(thePrimary.Get4Momentum().vect().unit());
148 }
149 return theParticleChange;
150 }
151 }
152
153 // get result from high energy model
154 G4KineticTrackVector * theInitialResult =
155 theHighEnergyGenerator->Scatter(theNucleus, aPart);
156
157 // Assign the creator model ID
158 for ( auto & ptr : *theInitialResult ) {
159 ptr->SetCreatorModelID( theStringModelID );
160 }
161
162 //#define DEBUG_initial_result
163 #ifdef DEBUG_initial_result
164 G4double E_out(0);
166 for(auto & ptr : *theInitialResult)
167 {
168 //G4cout << "TheoFS secondary, mom " << ptr->GetDefinition()->GetParticleName()
169 // << " " << ptr->Get4Momentum() << G4endl;
170 E_out += ptr->Get4Momentum().e();
171 }
172 G4double init_mass= ionTable->GetIonMass(theNucleus.GetZ_asInt(),theNucleus.GetA_asInt());
173 G4double init_E=aPart.Get4Momentum().e();
174 // residual nucleus
175
176 const std::vector<G4Nucleon> & thy = theHighEnergyGenerator->GetWoundedNucleus()->GetNucleons();
177
178 G4int resZ(0),resA(0);
179 G4double delta_m(0);
180 for(auto & nuc : thy)
181 {
182 if(nuc.AreYouHit()) {
183 ++resA;
184 if ( nuc.GetDefinition() == G4Proton::Proton() ) {
185 ++resZ;
186 delta_m += CLHEP::proton_mass_c2;
187 } else {
188 delta_m += CLHEP::neutron_mass_c2;
189 }
190 }
191 }
192
193 G4double final_mass(0);
194 if ( theNucleus.GetA_asInt() ) {
195 final_mass=ionTable->GetIonMass(theNucleus.GetZ_asInt()-resZ,theNucleus.GetA_asInt()- resA);
196 }
197 G4double E_excit=init_mass + init_E - final_mass - E_out;
198 G4cout << " Corrected delta mass " << init_mass - final_mass - delta_m << G4endl;
199 G4cout << "initial E, mass = " << init_E << ", " << init_mass << G4endl;
200 G4cout << " final E, mass = " << E_out <<", " << final_mass << " excitation_E " << E_excit << G4endl;
201 #endif
202
203 G4ReactionProductVector * theTransportResult = nullptr;
204
205 G4V3DNucleus* theProjectileNucleus = theHighEnergyGenerator->GetProjectileNucleus();
206 if(theProjectileNucleus == nullptr)
207 {
208 G4int hitCount = 0;
209 const std::vector<G4Nucleon>& they = theHighEnergyGenerator->GetWoundedNucleus()->GetNucleons();
210 for(auto & nuc : they)
211 {
212 if(nuc.AreYouHit()) ++hitCount;
213 }
214 if(hitCount != theHighEnergyGenerator->GetWoundedNucleus()->GetMassNumber() )
215 {
216 theTransport->SetPrimaryProjectile(thePrimary);
217 theTransportResult =
218 theTransport->Propagate(theInitialResult,
219 theHighEnergyGenerator->GetWoundedNucleus());
220 if ( !theTransportResult ) {
221 G4cout << "G4TheoFSGenerator: null ptr from transport propagate " << G4endl;
222 throw G4HadronicException(__FILE__, __LINE__, "Null ptr from transport propagate");
223 }
224 }
225 else
226 {
227 theTransportResult = theDecay.Propagate(theInitialResult,
228 theHighEnergyGenerator->GetWoundedNucleus());
229 if ( theTransportResult == nullptr ) {
230 G4cout << "G4TheoFSGenerator: null ptr from decay propagate " << G4endl;
231 throw G4HadronicException(__FILE__, __LINE__, "Null ptr from decay propagate");
232 }
233 }
234 }
235 else
236 {
237 theTransport->SetPrimaryProjectile(thePrimary);
238 theTransportResult = theTransport->PropagateNuclNucl(theInitialResult,
239 theHighEnergyGenerator->GetWoundedNucleus(),
240 theProjectileNucleus);
241 if ( !theTransportResult ) {
242 G4cout << "G4TheoFSGenerator: null ptr from transport propagate " << G4endl;
243 throw G4HadronicException(__FILE__, __LINE__, "Null ptr from transport propagate");
244 }
245 }
246
247 // If enabled, apply the Cosmic Rays (CR) coalescence to the list of secondaries produced so far.
248 // This algorithm can form deuterons and antideuterons by coalescence of, respectively,
249 // proton-neutron and antiproton-antineutron pairs close in momentum space.
250 // This can be useful in particular for Cosmic Ray applications.
251 if ( G4HadronicParameters::Instance()->EnableCRCoalescence() ) {
252 if(nullptr == theCosmicCoalescence) {
253 theCosmicCoalescence = (G4CRCoalescence*)
255 if(nullptr == theCosmicCoalescence) {
256 theCosmicCoalescence = new G4CRCoalescence();
257 }
258 }
259 theCosmicCoalescence->SetP0Coalescence( thePrimary, theHighEnergyGenerator->GetModelName() );
260 theCosmicCoalescence->GenerateDeuterons( theTransportResult );
261 }
262
263 // Fill particle change
264 for(auto & ptr : *theTransportResult)
265 {
266 G4DynamicParticle * aNewDP =
267 new G4DynamicParticle(ptr->GetDefinition(),
268 ptr->GetTotalEnergy(),
269 ptr->GetMomentum());
270 G4HadSecondary aNew = G4HadSecondary(aNewDP);
271 G4double time = std::max(ptr->GetFormationTime(), 0.0);
272 aNew.SetTime(timePrimary + time);
273 aNew.SetCreatorModelID(ptr->GetCreatorModelID());
274 aNew.SetParentResonanceDef(ptr->GetParentResonanceDef());
275 aNew.SetParentResonanceID(ptr->GetParentResonanceID());
276 theParticleChange->AddSecondary(aNew);
277 delete ptr;
278 }
279
280 // some garbage collection
281 delete theTransportResult;
282
283 // Done
284 return theParticleChange;
285}
@ isAlive
@ stopAndKill
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
Hep3Vector vect() const
void GenerateDeuterons(G4ReactionProductVector *result)
void SetP0Coalescence(const G4HadProjectile &thePrimary, G4String)
G4ReactionProductVector * Propagate(G4KineticTrackVector *theSecondaries, G4V3DNucleus *)
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetGlobalTime() const
void SetParentResonanceDef(const G4ParticleDefinition *parentDef)
void SetTime(G4double aT)
void SetCreatorModelID(G4int id)
void SetParentResonanceID(const G4int parentID)
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
const G4String & GetModelName() const
static G4HadronicParameters * Instance()
G4double GetIonMass(G4int Z, G4int A, G4int nL=0, G4int lvl=0) const
G4int GetA_asInt() const
Definition G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition G4Nucleus.hh:105
G4int GetQuarkContent(G4int flavor) const
G4bool IsHypernucleus() const
G4int GetAntiQuarkContent(G4int flavor) const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
static G4Proton * Proton()
Definition G4Proton.cc:90
G4KineticTrackVector * Scatter(G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
G4double GetFraction(G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
virtual const std::vector< G4Nucleon > & GetNucleons()=0
virtual G4int GetMassNumber()=0
virtual G4V3DNucleus * GetProjectileNucleus() const
virtual G4V3DNucleus * GetWoundedNucleus() const =0
virtual G4KineticTrackVector * Scatter(const G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)=0
virtual G4ReactionProductVector * PropagateNuclNucl(G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus, G4V3DNucleus *theProjectileNucleus)
void SetPrimaryProjectile(const G4HadProjectile &aPrimary)
virtual G4ReactionProductVector * Propagate(G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus)=0

◆ GetEnergyMomentumCheckLevels()

std::pair< G4double, G4double > G4TheoFSGenerator::GetEnergyMomentumCheckLevels ( ) const
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 287 of file G4TheoFSGenerator.cc.

288{
289 if ( theHighEnergyGenerator ) {
290 return theHighEnergyGenerator->GetEnergyMomentumCheckLevels();
291 } else {
292 return std::pair<G4double, G4double>(DBL_MAX, DBL_MAX);
293 }
294}
virtual std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
#define DBL_MAX
Definition templates.hh:62

◆ GetHighEnergyGenerator()

const G4VHighEnergyGenerator * G4TheoFSGenerator::GetHighEnergyGenerator ( ) const
inline

Definition at line 111 of file G4TheoFSGenerator.hh.

112{
113 return theHighEnergyGenerator;
114}

◆ GetQuasiElasticChannel()

const G4QuasiElasticChannel * G4TheoFSGenerator::GetQuasiElasticChannel ( ) const
inline

Definition at line 116 of file G4TheoFSGenerator.hh.

117{
118 return theQuasielastic;
119}

◆ GetTransport()

const G4VIntraNuclearTransportModel * G4TheoFSGenerator::GetTransport ( ) const
inline

Definition at line 106 of file G4TheoFSGenerator.hh.

107{
108 return theTransport;
109}

◆ ModelDescription()

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

Reimplemented from G4HadronicInteraction.

Definition at line 58 of file G4TheoFSGenerator.cc.

59{
60 outFile << GetModelName() <<" consists of a " << theHighEnergyGenerator->GetModelName()
61 << " string model and a stage to de-excite the excited nuclear fragment.\n<p>"
62 << "The string model simulates the interaction of\n"
63 << "an incident hadron with a nucleus, forming \n"
64 << "excited strings, decays these strings into hadrons,\n"
65 << "and leaves an excited nucleus. \n"
66 << "<p>The string model:\n";
67 theHighEnergyGenerator->ModelDescription(outFile);
68 outFile <<"\n<p>";
69 theTransport->PropagateModelDescription(outFile);
70}
void ModelDescription(std::ostream &) const override
virtual void PropagateModelDescription(std::ostream &outFile) const

◆ operator!=()

G4bool G4TheoFSGenerator::operator!= ( const G4TheoFSGenerator & right) const
delete

◆ operator=()

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

◆ operator==()

G4bool G4TheoFSGenerator::operator== ( const G4TheoFSGenerator & right) const
delete

◆ SetHighEnergyGenerator()

void G4TheoFSGenerator::SetHighEnergyGenerator ( G4VHighEnergyGenerator *const value)
inline

Definition at line 96 of file G4TheoFSGenerator.hh.

97{
98 theHighEnergyGenerator= value;
99}

Referenced by G4FTFBuilder::BuildModel(), G4QGSBuilder::BuildModel(), LBE::ConstructHad(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4FTFBinaryKaonBuilder::G4FTFBinaryKaonBuilder(), G4FTFBinaryNeutronBuilder::G4FTFBinaryNeutronBuilder(), G4FTFBinaryPiKBuilder::G4FTFBinaryPiKBuilder(), G4FTFBinaryPionBuilder::G4FTFBinaryPionBuilder(), G4FTFBinaryProtonBuilder::G4FTFBinaryProtonBuilder(), G4FTFPAntiBarionBuilder::G4FTFPAntiBarionBuilder(), G4FTFPKaonBuilder::G4FTFPKaonBuilder(), G4FTFPNeutronBuilder::G4FTFPNeutronBuilder(), G4FTFPPiKBuilder::G4FTFPPiKBuilder(), G4FTFPPionBuilder::G4FTFPPionBuilder(), G4FTFPProtonBuilder::G4FTFPProtonBuilder(), G4HadronicAbsorptionFritiof::G4HadronicAbsorptionFritiof(), G4HadronicAbsorptionFritiofWithBinaryCascade::G4HadronicAbsorptionFritiofWithBinaryCascade(), G4HyperonFTFPBuilder::G4HyperonFTFPBuilder(), G4HyperonQGSPBuilder::G4HyperonQGSPBuilder(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4QGSBinaryKaonBuilder::G4QGSBinaryKaonBuilder(), G4QGSBinaryNeutronBuilder::G4QGSBinaryNeutronBuilder(), G4QGSBinaryPiKBuilder::G4QGSBinaryPiKBuilder(), G4QGSBinaryPionBuilder::G4QGSBinaryPionBuilder(), G4QGSBinaryProtonBuilder::G4QGSBinaryProtonBuilder(), G4QGSPAntiBarionBuilder::G4QGSPAntiBarionBuilder(), G4QGSPKaonBuilder::G4QGSPKaonBuilder(), G4QGSPLundStrFragmProtonBuilder::G4QGSPLundStrFragmProtonBuilder(), G4QGSPNeutronBuilder::G4QGSPNeutronBuilder(), G4QGSPPiKBuilder::G4QGSPPiKBuilder(), G4QGSPPionBuilder::G4QGSPPionBuilder(), and G4QGSPProtonBuilder::G4QGSPProtonBuilder().

◆ SetQuasiElasticChannel()

◆ SetTransport()

void G4TheoFSGenerator::SetTransport ( G4VIntraNuclearTransportModel *const value)
inline

Definition at line 91 of file G4TheoFSGenerator.hh.

92{
93 theTransport = value;
94}

Referenced by G4FTFBuilder::BuildModel(), G4QGSBuilder::BuildModel(), LBE::ConstructHad(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4FTFBinaryKaonBuilder::G4FTFBinaryKaonBuilder(), G4FTFBinaryNeutronBuilder::G4FTFBinaryNeutronBuilder(), G4FTFBinaryPiKBuilder::G4FTFBinaryPiKBuilder(), G4FTFBinaryPionBuilder::G4FTFBinaryPionBuilder(), G4FTFBinaryProtonBuilder::G4FTFBinaryProtonBuilder(), G4FTFPAntiBarionBuilder::G4FTFPAntiBarionBuilder(), G4FTFPKaonBuilder::G4FTFPKaonBuilder(), G4FTFPNeutronBuilder::G4FTFPNeutronBuilder(), G4FTFPPiKBuilder::G4FTFPPiKBuilder(), G4FTFPPionBuilder::G4FTFPPionBuilder(), G4FTFPProtonBuilder::G4FTFPProtonBuilder(), G4HadronicAbsorptionFritiof::G4HadronicAbsorptionFritiof(), G4HadronicAbsorptionFritiofWithBinaryCascade::G4HadronicAbsorptionFritiofWithBinaryCascade(), G4HyperonFTFPBuilder::G4HyperonFTFPBuilder(), G4HyperonQGSPBuilder::G4HyperonQGSPBuilder(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4QGSBinaryKaonBuilder::G4QGSBinaryKaonBuilder(), G4QGSBinaryNeutronBuilder::G4QGSBinaryNeutronBuilder(), G4QGSBinaryPiKBuilder::G4QGSBinaryPiKBuilder(), G4QGSBinaryPionBuilder::G4QGSBinaryPionBuilder(), G4QGSBinaryProtonBuilder::G4QGSBinaryProtonBuilder(), G4QGSPAntiBarionBuilder::G4QGSPAntiBarionBuilder(), G4QGSPKaonBuilder::G4QGSPKaonBuilder(), G4QGSPLundStrFragmProtonBuilder::G4QGSPLundStrFragmProtonBuilder(), G4QGSPNeutronBuilder::G4QGSPNeutronBuilder(), G4QGSPPiKBuilder::G4QGSPPiKBuilder(), G4QGSPPionBuilder::G4QGSPPionBuilder(), and G4QGSPProtonBuilder::G4QGSPProtonBuilder().


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