Geant4 9.6.0
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 ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &thePrimary, G4Nucleus &theNucleus)
 
void SetTransport (G4VIntraNuclearTransportModel *const value)
 
void SetHighEnergyGenerator (G4VHighEnergyGenerator *const value)
 
void SetQuasiElasticChannel (G4QuasiElasticChannel *const value)
 
void SetProjectileDiffraction (G4ProjectileDiffractiveChannel *const value)
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void ModelDescription (std::ostream &outFile) const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &, G4Nucleus &)
 
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)
 
const G4HadronicInteractionGetMyPointer () const
 
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
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 

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

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

Definition at line 40 of file G4TheoFSGenerator.cc.

42 , theTransport(0), theHighEnergyGenerator(0)
43 , theQuasielastic(0), theProjectileDiffraction(0)
44 {
45 theParticleChange = new G4HadFinalState;
46}

◆ ~G4TheoFSGenerator()

G4TheoFSGenerator::~G4TheoFSGenerator ( )

Definition at line 48 of file G4TheoFSGenerator.cc.

49{
50 delete theParticleChange;
51}

Member Function Documentation

◆ ApplyYourself()

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

Implements G4HadronicInteraction.

Definition at line 68 of file G4TheoFSGenerator.cc.

69{
70 // init particle change
71 theParticleChange->Clear();
72 theParticleChange->SetStatusChange(stopAndKill);
73
74 // check if models have been registered, and use default, in case this is not true @@
75
76 // get result from high energy model
77 G4DynamicParticle aTemp(const_cast<G4ParticleDefinition *>(thePrimary.GetDefinition()),
78 thePrimary.Get4Momentum().vect());
79 const G4DynamicParticle * aPart = &aTemp;
80
81 if ( theQuasielastic ) {
82
83 if ( theQuasielastic->GetFraction(theNucleus, *aPart) > G4UniformRand() )
84 {
85 //G4cout<<"___G4TheoFSGenerator: before Scatter (1) QE=" << theQuasielastic<<G4endl;
86 G4KineticTrackVector *result= theQuasielastic->Scatter(theNucleus, *aPart);
87 //G4cout << "^^G4TheoFSGenerator: after Scatter (1) " << G4endl;
88 if (result)
89 {
90 for(unsigned int i=0; i<result->size(); i++)
91 {
92 G4DynamicParticle * aNew =
93 new G4DynamicParticle(result->operator[](i)->GetDefinition(),
94 result->operator[](i)->Get4Momentum().e(),
95 result->operator[](i)->Get4Momentum().vect());
96 theParticleChange->AddSecondary(aNew);
97 delete result->operator[](i);
98 }
99 delete result;
100
101 } else
102 {
103 theParticleChange->SetStatusChange(isAlive);
104 theParticleChange->SetEnergyChange(thePrimary.GetKineticEnergy());
105 theParticleChange->SetMomentumChange(thePrimary.Get4Momentum().vect().unit());
106
107 }
108 return theParticleChange;
109 }
110 }
111 if ( theProjectileDiffraction) {
112
113 if ( theProjectileDiffraction->GetFraction(theNucleus, *aPart) > G4UniformRand() )
114 // strictly, this returns fraction on inelastic, so quasielastic should
115 // already be substracted, ie. quasielastic should always be used
116 // before diffractive
117 {
118 //G4cout << "____G4TheoFSGenerator: before Scatter (2) " << G4endl;
119 G4KineticTrackVector *result= theProjectileDiffraction->Scatter(theNucleus, *aPart);
120 //G4cout << "^^^^G4TheoFSGenerator: after Scatter (2) " << G4endl;
121 if (result)
122 {
123 for(unsigned int i=0; i<result->size(); i++)
124 {
125 G4DynamicParticle * aNew =
126 new G4DynamicParticle(result->operator[](i)->GetDefinition(),
127 result->operator[](i)->Get4Momentum().e(),
128 result->operator[](i)->Get4Momentum().vect());
129 theParticleChange->AddSecondary(aNew);
130 delete result->operator[](i);
131 }
132 delete result;
133
134 } else
135 {
136 theParticleChange->SetStatusChange(isAlive);
137 theParticleChange->SetEnergyChange(thePrimary.GetKineticEnergy());
138 theParticleChange->SetMomentumChange(thePrimary.Get4Momentum().vect().unit());
139
140 }
141 return theParticleChange;
142 }
143 }
144 G4KineticTrackVector * theInitialResult =
145 theHighEnergyGenerator->Scatter(theNucleus, *aPart);
146
147//#define DEBUG_initial_result
148 #ifdef DEBUG_initial_result
149 G4double E_out(0);
151 std::vector<G4KineticTrack *>::iterator ir_iter;
152 for(ir_iter=theInitialResult->begin(); ir_iter!=theInitialResult->end(); ir_iter++)
153 {
154 //G4cout << "TheoFS secondary, mom " << (*ir_iter)->GetDefinition()->GetParticleName() << " " << (*ir_iter)->Get4Momentum() << G4endl;
155 E_out += (*ir_iter)->Get4Momentum().e();
156 }
157 G4double init_mass= ionTable->GetIonMass(theNucleus.GetZ_asInt(),theNucleus.GetA_asInt());
158 G4double init_E=aPart->Get4Momentum().e();
159 // residual nucleus
160 const std::vector<G4Nucleon> & thy = theHighEnergyGenerator->GetWoundedNucleus()->GetNucleons();
161 G4int resZ(0),resA(0);
162 G4double delta_m(0);
163 for(size_t them=0; them<thy.size(); them++)
164 {
165 if(thy[them].AreYouHit()) {
166 ++resA;
167 if ( thy[them].GetDefinition() == G4Proton::Proton() ) {
168 ++resZ;
169 delta_m +=G4Proton::Proton()->GetPDGMass();
170 } else {
171 delta_m +=G4Neutron::Neutron()->GetPDGMass();
172 }
173 }
174 }
175 G4double final_mass(0);
176 if ( theNucleus.GetA_asInt() ) {
177 final_mass=ionTable->GetIonMass(theNucleus.GetZ_asInt()-resZ,theNucleus.GetA_asInt()- resA);
178 }
179 G4double E_excit=init_mass + init_E - final_mass - E_out;
180 G4cout << " Corrected delta mass " << init_mass - final_mass - delta_m << G4endl;
181 G4cout << "initial E, mass = " << init_E << ", " << init_mass << G4endl;
182 G4cout << " final E, mass = " << E_out <<", " << final_mass << " excitation_E " << E_excit << G4endl;
183 #endif
184
185 G4ReactionProductVector * theTransportResult = NULL;
186 G4int hitCount = 0;
187 const std::vector<G4Nucleon>& they = theHighEnergyGenerator->GetWoundedNucleus()->GetNucleons();
188 for(size_t them=0; them<they.size(); them++)
189 {
190 if(they[them].AreYouHit()) hitCount ++;
191 }
192 if(hitCount != theHighEnergyGenerator->GetWoundedNucleus()->GetMassNumber() )
193 {
194 theTransport->SetPrimaryProjectile(thePrimary); // For Bertini Cascade
195 theTransportResult =
196 theTransport->Propagate(theInitialResult, theHighEnergyGenerator->GetWoundedNucleus());
197 if ( !theTransportResult ) {
198 G4cout << "G4TheoFSGenerator: null ptr from transport propagate " << G4endl;
199 throw G4HadronicException(__FILE__, __LINE__, "Null ptr from transport propagate");
200 }
201 }
202 else
203 {
204 theTransportResult = theDecay.Propagate(theInitialResult, theHighEnergyGenerator->GetWoundedNucleus());
205 if ( !theTransportResult ) {
206 G4cout << "G4TheoFSGenerator: null ptr from decay propagate " << G4endl;
207 throw G4HadronicException(__FILE__, __LINE__, "Null ptr from decay propagate");
208 }
209 }
210
211 // Fill particle change
212 unsigned int i;
213 for(i=0; i<theTransportResult->size(); i++)
214 {
215 G4DynamicParticle * aNew =
216 new G4DynamicParticle(theTransportResult->operator[](i)->GetDefinition(),
217 theTransportResult->operator[](i)->GetTotalEnergy(),
218 theTransportResult->operator[](i)->GetMomentum());
219 // @@@ - overkill? G4double newTime = theParticleChange->GetGlobalTime(theTransportResult->operator[](i)->GetFormationTime());
220 theParticleChange->AddSecondary(aNew);
221 delete theTransportResult->operator[](i);
222 }
223
224 // some garbage collection
225 delete theTransportResult;
226
227 // Done
228 return theParticleChange;
229}
@ isAlive
@ stopAndKill
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
Hep3Vector vect() const
G4ReactionProductVector * Propagate(G4KineticTrackVector *theSecondaries, G4V3DNucleus *)
G4LorentzVector Get4Momentum() const
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetIonMass(G4int Z, G4int A, G4int L=0) const
!! Only ground states are supported now
Definition: G4IonTable.cc:774
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
static G4ParticleTable * GetParticleTable()
G4IonTable * GetIonTable()
G4double GetFraction(G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
G4KineticTrackVector * Scatter(G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
static G4Proton * Proton()
Definition: G4Proton.cc:93
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 * GetWoundedNucleus() const =0
virtual G4KineticTrackVector * Scatter(const G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)=0
void SetPrimaryProjectile(const G4HadProjectile &aPrimary)
virtual G4ReactionProductVector * Propagate(G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus)=0

Referenced by G4ElectroNuclearReaction::ApplyYourself().

◆ GetEnergyMomentumCheckLevels()

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

Reimplemented from G4HadronicInteraction.

Definition at line 231 of file G4TheoFSGenerator.cc.

232{
233 if ( theHighEnergyGenerator ) {
234 return theHighEnergyGenerator->GetEnergyMomentumCheckLevels();
235 } else {
236 return std::pair<G4double, G4double>(DBL_MAX, DBL_MAX);
237 }
238}
std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
#define DBL_MAX
Definition: templates.hh:83

◆ ModelDescription()

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

Reimplemented from G4HadronicInteraction.

Definition at line 53 of file G4TheoFSGenerator.cc.

54{
55 outFile << GetModelName() <<" consists of a " << theHighEnergyGenerator->GetModelName()
56 << " string model and of "
57 << ".\n"
58 << "The string model simulates the interaction of\n"
59 << "an incident hadron with a nucleus, forming \n"
60 << "excited strings, decays these strings into hadrons,\n"
61 << "and leaves an excited nucleus.\n"
62 << "The string model:\n";
63 theHighEnergyGenerator->ModelDescription(outFile);
64//theTransport->IntraNuclearTransportDescription(outFile)
65}
const G4String & GetModelName() const
virtual G4String GetModelName() const
virtual void ModelDescription(std::ostream &) const

◆ SetHighEnergyGenerator()

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

Definition at line 106 of file G4TheoFSGenerator.hh.

107{
108 theHighEnergyGenerator= value;
109}

Referenced by G4BertiniElectroNuclearBuilder::Build(), G4ElectroNuclearBuilder::Build(), G4MiscQGSCBuilder::Build(), G4FTFBuilder::BuildModel(), G4QGSBuilder::BuildModel(), G4ElectroNuclearReaction::G4ElectroNuclearReaction(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4FTFBinaryKaonBuilder::G4FTFBinaryKaonBuilder(), G4FTFBinaryNeutronBuilder::G4FTFBinaryNeutronBuilder(), G4FTFBinaryPiKBuilder::G4FTFBinaryPiKBuilder(), G4FTFBinaryPionBuilder::G4FTFBinaryPionBuilder(), G4FTFBinaryProtonBuilder::G4FTFBinaryProtonBuilder(), G4FTFCNeutronBuilder::G4FTFCNeutronBuilder(), G4FTFCPiKBuilder::G4FTFCPiKBuilder(), G4FTFCProtonBuilder::G4FTFCProtonBuilder(), G4FTFPAntiBarionBuilder::G4FTFPAntiBarionBuilder(), G4FTFPNeutronBuilder::G4FTFPNeutronBuilder(), G4FTFPPiKBuilder::G4FTFPPiKBuilder(), G4FTFPProtonBuilder::G4FTFPProtonBuilder(), G4HadronicAbsorptionFritiof::G4HadronicAbsorptionFritiof(), G4HyperonFTFPBuilder::G4HyperonFTFPBuilder(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4QGSBinaryNeutronBuilder::G4QGSBinaryNeutronBuilder(), G4QGSBinaryPiKBuilder::G4QGSBinaryPiKBuilder(), G4QGSBinaryProtonBuilder::G4QGSBinaryProtonBuilder(), G4QGSC_CHIPSNeutronBuilder::G4QGSC_CHIPSNeutronBuilder(), G4QGSC_CHIPSPiKBuilder::G4QGSC_CHIPSPiKBuilder(), G4QGSC_CHIPSProtonBuilder::G4QGSC_CHIPSProtonBuilder(), G4QGSC_QGSCNeutronBuilder::G4QGSC_QGSCNeutronBuilder(), G4QGSC_QGSCPiKBuilder::G4QGSC_QGSCPiKBuilder(), G4QGSC_QGSCProtonBuilder::G4QGSC_QGSCProtonBuilder(), G4QGSCEflowNeutronBuilder::G4QGSCEflowNeutronBuilder(), G4QGSCEflowPiKBuilder::G4QGSCEflowPiKBuilder(), G4QGSCEflowProtonBuilder::G4QGSCEflowProtonBuilder(), G4QGSCNeutronBuilder::G4QGSCNeutronBuilder(), G4QGSCPiKBuilder::G4QGSCPiKBuilder(), G4QGSCProtonBuilder::G4QGSCProtonBuilder(), G4QGSPNeutronBuilder::G4QGSPNeutronBuilder(), G4QGSPPiKBuilder::G4QGSPPiKBuilder(), G4QGSPPionBuilder::G4QGSPPionBuilder(), G4QGSPProtonBuilder::G4QGSPProtonBuilder(), and TheoModelFactory< C, S, F >::New().

◆ SetProjectileDiffraction()

void G4TheoFSGenerator::SetProjectileDiffraction ( G4ProjectileDiffractiveChannel *const  value)
inline

◆ SetQuasiElasticChannel()

void G4TheoFSGenerator::SetQuasiElasticChannel ( G4QuasiElasticChannel *const  value)
inline

Definition at line 111 of file G4TheoFSGenerator.hh.

112{
113 theQuasielastic = value;
114}

Referenced by G4MiscQGSCBuilder::Build(), G4QGSBuilder::BuildModel(), G4FTFBinaryKaonBuilder::G4FTFBinaryKaonBuilder(), G4FTFBinaryNeutronBuilder::G4FTFBinaryNeutronBuilder(), G4FTFBinaryPiKBuilder::G4FTFBinaryPiKBuilder(), G4FTFBinaryPionBuilder::G4FTFBinaryPionBuilder(), G4FTFBinaryProtonBuilder::G4FTFBinaryProtonBuilder(), G4FTFCNeutronBuilder::G4FTFCNeutronBuilder(), G4FTFCPiKBuilder::G4FTFCPiKBuilder(), G4FTFCProtonBuilder::G4FTFCProtonBuilder(), G4FTFPAntiBarionBuilder::G4FTFPAntiBarionBuilder(), G4FTFPNeutronBuilder::G4FTFPNeutronBuilder(), G4FTFPPiKBuilder::G4FTFPPiKBuilder(), G4FTFPProtonBuilder::G4FTFPProtonBuilder(), G4QGSBinaryNeutronBuilder::G4QGSBinaryNeutronBuilder(), G4QGSBinaryPiKBuilder::G4QGSBinaryPiKBuilder(), G4QGSBinaryProtonBuilder::G4QGSBinaryProtonBuilder(), G4QGSC_CHIPSNeutronBuilder::G4QGSC_CHIPSNeutronBuilder(), G4QGSC_CHIPSPiKBuilder::G4QGSC_CHIPSPiKBuilder(), G4QGSC_CHIPSProtonBuilder::G4QGSC_CHIPSProtonBuilder(), G4QGSC_QGSCNeutronBuilder::G4QGSC_QGSCNeutronBuilder(), G4QGSC_QGSCPiKBuilder::G4QGSC_QGSCPiKBuilder(), G4QGSC_QGSCProtonBuilder::G4QGSC_QGSCProtonBuilder(), G4QGSCEflowNeutronBuilder::G4QGSCEflowNeutronBuilder(), G4QGSCEflowPiKBuilder::G4QGSCEflowPiKBuilder(), G4QGSCEflowProtonBuilder::G4QGSCEflowProtonBuilder(), G4QGSCNeutronBuilder::G4QGSCNeutronBuilder(), G4QGSCPiKBuilder::G4QGSCPiKBuilder(), G4QGSCProtonBuilder::G4QGSCProtonBuilder(), G4QGSPNeutronBuilder::G4QGSPNeutronBuilder(), G4QGSPPiKBuilder::G4QGSPPiKBuilder(), G4QGSPPionBuilder::G4QGSPPionBuilder(), and G4QGSPProtonBuilder::G4QGSPProtonBuilder().

◆ SetTransport()

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

Definition at line 96 of file G4TheoFSGenerator.hh.

97{
98 theTransport = value;
99}

Referenced by G4BertiniElectroNuclearBuilder::Build(), G4ElectroNuclearBuilder::Build(), G4MiscQGSCBuilder::Build(), G4FTFBuilder::BuildModel(), G4QGSBuilder::BuildModel(), G4ElectroNuclearReaction::G4ElectroNuclearReaction(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4FTFBinaryKaonBuilder::G4FTFBinaryKaonBuilder(), G4FTFBinaryNeutronBuilder::G4FTFBinaryNeutronBuilder(), G4FTFBinaryPiKBuilder::G4FTFBinaryPiKBuilder(), G4FTFBinaryPionBuilder::G4FTFBinaryPionBuilder(), G4FTFBinaryProtonBuilder::G4FTFBinaryProtonBuilder(), G4FTFCNeutronBuilder::G4FTFCNeutronBuilder(), G4FTFCPiKBuilder::G4FTFCPiKBuilder(), G4FTFCProtonBuilder::G4FTFCProtonBuilder(), G4FTFPAntiBarionBuilder::G4FTFPAntiBarionBuilder(), G4FTFPNeutronBuilder::G4FTFPNeutronBuilder(), G4FTFPPiKBuilder::G4FTFPPiKBuilder(), G4FTFPProtonBuilder::G4FTFPProtonBuilder(), G4HadronicAbsorptionFritiof::G4HadronicAbsorptionFritiof(), G4HyperonFTFPBuilder::G4HyperonFTFPBuilder(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4QGSBinaryNeutronBuilder::G4QGSBinaryNeutronBuilder(), G4QGSBinaryPiKBuilder::G4QGSBinaryPiKBuilder(), G4QGSBinaryProtonBuilder::G4QGSBinaryProtonBuilder(), G4QGSC_CHIPSNeutronBuilder::G4QGSC_CHIPSNeutronBuilder(), G4QGSC_CHIPSPiKBuilder::G4QGSC_CHIPSPiKBuilder(), G4QGSC_CHIPSProtonBuilder::G4QGSC_CHIPSProtonBuilder(), G4QGSC_QGSCNeutronBuilder::G4QGSC_QGSCNeutronBuilder(), G4QGSC_QGSCPiKBuilder::G4QGSC_QGSCPiKBuilder(), G4QGSC_QGSCProtonBuilder::G4QGSC_QGSCProtonBuilder(), G4QGSCEflowNeutronBuilder::G4QGSCEflowNeutronBuilder(), G4QGSCEflowPiKBuilder::G4QGSCEflowPiKBuilder(), G4QGSCEflowProtonBuilder::G4QGSCEflowProtonBuilder(), G4QGSCNeutronBuilder::G4QGSCNeutronBuilder(), G4QGSCPiKBuilder::G4QGSCPiKBuilder(), G4QGSCProtonBuilder::G4QGSCProtonBuilder(), G4QGSPNeutronBuilder::G4QGSPNeutronBuilder(), G4QGSPPiKBuilder::G4QGSPPiKBuilder(), G4QGSPPionBuilder::G4QGSPPionBuilder(), G4QGSPProtonBuilder::G4QGSPProtonBuilder(), and TheoModelFactory< C, S, F >::New().


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