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

#include <G4HadronStoppingProcess.hh>

+ Inheritance diagram for G4HadronStoppingProcess:

Public Member Functions

 G4HadronStoppingProcess (const G4String &name="hadronCaptureAtRest")
 
virtual ~G4HadronStoppingProcess ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual void ProcessDescription (std::ostream &outFile) const
 
void SetElementSelector (G4ElementSelector *ptr)
 
void SetEmCascade (G4HadronicInteraction *ptr)
 
void SetBoundDecay (G4HadronicInteraction *ptr)
 
- Public Member Functions inherited from G4HadronicProcess
 G4HadronicProcess (const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)
 
 G4HadronicProcess (const G4String &processName, G4HadronicProcessType subType)
 
 ~G4HadronicProcess () override
 
void RegisterMe (G4HadronicInteraction *a)
 
G4double GetElementCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
 
G4double GetMicroscopicCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
 
void StartTracking (G4Track *track) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double, G4ForceCondition *) override
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
void PreparePhysicsTable (const G4ParticleDefinition &) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void DumpPhysicsTable (const G4ParticleDefinition &p)
 
void AddDataSet (G4VCrossSectionDataSet *aDataSet)
 
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList ()
 
G4HadronicInteractionGetHadronicModel (const G4String &)
 
G4HadronicInteractionGetHadronicInteraction () const
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *) override
 
const G4NucleusGetTargetNucleus () const
 
G4NucleusGetTargetNucleusPointer ()
 
const G4IsotopeGetTargetIsotope ()
 
G4double ComputeCrossSection (const G4ParticleDefinition *, const G4Material *, const G4double kinEnergy)
 
G4HadXSType CrossSectionType () const
 
void SetCrossSectionType (G4HadXSType val)
 
void ProcessDescription (std::ostream &outFile) const override
 
void BiasCrossSectionByFactor (G4double aScale)
 
void MultiplyCrossSectionBy (G4double factor)
 
G4double CrossSectionFactor () const
 
void SetIntegral (G4bool val)
 
void SetEpReportLevel (G4int level)
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
G4CrossSectionDataStoreGetCrossSectionDataStore ()
 
std::vector< G4TwoPeaksHadXS * > * TwoPeaksXS () const
 
std::vector< G4double > * EnergyOfCrossSectionMax () const
 
G4HadronicProcessoperator= (const G4HadronicProcess &right)=delete
 
 G4HadronicProcess (const G4HadronicProcess &)=delete
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

G4double GetMeanLifeTime (const G4Track &, G4ForceCondition *)
 
- Protected Member Functions inherited from G4HadronicProcess
G4HadronicInteractionChooseHadronicInteraction (const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, const G4Material *aMaterial, const G4Element *anElement)
 
G4double GetLastCrossSection ()
 
void FillResult (G4HadFinalState *aR, const G4Track &aT)
 
void DumpState (const G4Track &, const G4String &, G4ExceptionDescription &)
 
G4HadFinalStateCheckResult (const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result)
 
void CheckEnergyMomentumConservation (const G4Track &, const G4Nucleus &)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4HadronicProcess
G4HadProjectile thePro
 
G4ParticleChangetheTotalResult
 
G4CrossSectionDataStoretheCrossSectionDataStore
 
G4double fWeight = 1.0
 
G4double aScaleFactor = 1.0
 
G4double theLastCrossSection = 0.0
 
G4double mfpKinEnergy = DBL_MAX
 
G4long epReportLevel = 0
 
G4HadXSType fXSType = fHadNoIntegral
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 62 of file G4HadronStoppingProcess.hh.

Constructor & Destructor Documentation

◆ G4HadronStoppingProcess()

G4HadronStoppingProcess::G4HadronStoppingProcess ( const G4String name = "hadronCaptureAtRest")
explicit

Definition at line 63 of file G4HadronStoppingProcess.cc.

65 fElementSelector(new G4ElementSelector()),
66 fEmCascade(new G4EmCaptureCascade()), // Owned by InteractionRegistry
67 fBoundDecay(0),
68 emcID(-1),
69 ncID(-1),
70 dioID(-1)
71{
72 // Modify G4VProcess flags to emulate G4VRest instead of G4VDiscrete
73 enableAtRestDoIt = true;
74 enablePostStepDoIt = false;
75
77}
@ fHadronAtRest
void RegisterExtraProcess(G4VProcess *)
static G4HadronicProcessStore * Instance()
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:363
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:365

◆ ~G4HadronStoppingProcess()

G4HadronStoppingProcess::~G4HadronStoppingProcess ( )
virtual

Definition at line 81 of file G4HadronStoppingProcess.cc.

82{
83 //G4HadronicProcessStore::Instance()->DeRegisterExtraProcess(this);
84 delete fElementSelector;
85 // NOTE: fEmCascade and fEmBoundDecay owned by registry, not locally
86}

Member Function Documentation

◆ AtRestDoIt()

G4VParticleChange * G4HadronStoppingProcess::AtRestDoIt ( const G4Track track,
const G4Step  
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 133 of file G4HadronStoppingProcess.cc.

135{
136 // if primary is not Alive then do nothing
138
140 const G4Element* elm = fElementSelector->SelectZandA(track, nucleus);
141
142 G4HadFinalState* result = 0;
143 thePro.Initialise(track);
144
145 // save track time an dstart capture from zero time
147 G4double time0 = track.GetGlobalTime();
148
149 G4bool nuclearCapture = true;
150
151 // Do the electromagnetic cascade in the nuclear field.
152 // EM cascade should keep G4HadFinalState object,
153 // because it will not be deleted at the end of this method
154 //
155 result = fEmCascade->ApplyYourself(thePro, *nucleus);
156 G4double ebound = result->GetLocalEnergyDeposit();
157 G4double edep = 0.0;
158 G4int nSecondaries = (G4int)result->GetNumberOfSecondaries();
159 G4int nEmCascadeSec = nSecondaries;
160
161 // Try decay from bound level
162 // For mu- the time of projectile should be changed.
163 // Decay should keep G4HadFinalState object,
164 // because it will not be deleted at the end of this method.
165 //
166 thePro.SetBoundEnergy(ebound);
167 if(fBoundDecay) {
168 G4HadFinalState* resultDecay =
169 fBoundDecay->ApplyYourself(thePro, *nucleus);
170 G4int n = (G4int)resultDecay->GetNumberOfSecondaries();
171 if(0 < n) {
172 nSecondaries += n;
173 result->AddSecondaries(resultDecay);
174 }
175 if(resultDecay->GetStatusChange() == stopAndKill) {
176 nuclearCapture = false;
177 }
178 resultDecay->Clear();
179 }
180
181 if(nuclearCapture) {
182
183 // delay of capture
184 G4double capTime = thePro.GetGlobalTime();
186
187 // select model
188 G4HadronicInteraction* model = 0;
189 try {
190 model = ChooseHadronicInteraction(thePro, *nucleus,
191 track.GetMaterial(), elm);
192 }
193 catch(G4HadronicException & aE) {
195 ed << "Target element "<<elm->GetName()<<" Z= "
196 << nucleus->GetZ_asInt() << " A= "
197 << nucleus->GetA_asInt() << G4endl;
198 DumpState(track,"ChooseHadronicInteraction",ed);
199 ed << " No HadronicInteraction found out" << G4endl;
200 G4Exception("G4HadronStoppingProcess::AtRestDoIt", "had005",
201 FatalException, ed);
202 }
203
204 G4HadFinalState* resultNuc = 0;
205 G4int reentryCount = 0;
206 do {
207 // sample final state
208 // nuclear interaction should keep G4HadFinalState object
209 // model should define time of each secondary particle
210 try {
211 resultNuc = model->ApplyYourself(thePro, *nucleus);
212 ++reentryCount;
213 }
214 catch(G4HadronicException & aR) {
216 ed << "Call for " << model->GetModelName() << G4endl;
217 ed << "Target element "<<elm->GetName()<<" Z= "
218 << nucleus->GetZ_asInt()
219 << " A= " << nucleus->GetA_asInt() << G4endl;
220 DumpState(track,"ApplyYourself",ed);
221 ed << " ApplyYourself failed" << G4endl;
222 G4Exception("G4HadronStoppingProcess::AtRestDoIt", "had006",
223 FatalException, ed);
224 }
225
226 // Check the result for catastrophic energy non-conservation
227 resultNuc = CheckResult(thePro, *nucleus, resultNuc);
228
229 if(reentryCount>100) {
231 ed << "Call for " << model->GetModelName() << G4endl;
232 ed << "Target element "<<elm->GetName()<<" Z= "
233 << nucleus->GetZ_asInt()
234 << " A= " << nucleus->GetA_asInt() << G4endl;
235 DumpState(track,"ApplyYourself",ed);
236 ed << " ApplyYourself does not completed after 100 attempts" << G4endl;
237 G4Exception("G4HadronStoppingProcess::AtRestDoIt", "had006",
238 FatalException, ed);
239 }
240 // Loop checking, 06-Aug-2015, Vladimir Ivanchenko
241 } while(!resultNuc);
242
243 edep = resultNuc->GetLocalEnergyDeposit();
244 std::size_t nnuc = resultNuc->GetNumberOfSecondaries();
245
246 // add delay time of capture
247 for(std::size_t i=0; i<nnuc; ++i) {
248 G4HadSecondary* sec = resultNuc->GetSecondary(i);
249 sec->SetTime(capTime + sec->GetTime());
250 }
251
252 nSecondaries += nnuc;
253 result->AddSecondaries(resultNuc);
254 resultNuc->Clear();
255 }
256
257 // Fill results
258 //
262 G4double w = track.GetWeight();
264 for(G4int i=0; i<nSecondaries; ++i) {
265 G4HadSecondary* sec = result->GetSecondary(i);
266
267 // add track global time to the reaction time
268 G4double time = sec->GetTime();
269 if(time < 0.0) { time = 0.0; }
270 time += time0;
271
272 // create secondary track
273 G4Track* t = new G4Track(sec->GetParticle(),
274 time,
275 track.GetPosition());
276 t->SetWeight(w*sec->GetWeight());
277
278 // use SetCreatorModelID to "label" the track
279 if (i<nEmCascadeSec) {
280 t->SetCreatorModelID(emcID);
281 } else if (nuclearCapture) {
282 t->SetCreatorModelID(ncID);
283 } else {
284 t->SetCreatorModelID(dioID);
285 }
286
289 }
290 result->Clear();
291
292 if (epReportLevel != 0) {
293 CheckEnergyMomentumConservation(track, *nucleus);
294 }
295 return theTotalResult;
296}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
@ stopAndKill
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
virtual const G4Element * SelectZandA(const G4Track &track, G4Nucleus *)
const G4String & GetName() const
Definition: G4Element.hh:127
G4HadFinalStateStatus GetStatusChange() const
G4double GetLocalEnergyDeposit() const
void AddSecondaries(const std::vector< G4HadSecondary > &addSecs)
std::size_t GetNumberOfSecondaries() const
G4HadSecondary * GetSecondary(size_t i)
void Initialise(const G4Track &aT)
void SetGlobalTime(G4double t)
void SetBoundEnergy(G4double e)
G4double GetGlobalTime() const
G4DynamicParticle * GetParticle()
G4double GetWeight() const
void SetTime(G4double aT)
G4double GetTime() const
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
const G4String & GetModelName() const
G4HadProjectile thePro
G4Nucleus * GetTargetNucleusPointer()
G4HadFinalState * CheckResult(const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result)
G4ParticleChange * theTotalResult
G4HadronicInteraction * ChooseHadronicInteraction(const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, const G4Material *aMaterial, const G4Element *anElement)
void CheckEnergyMomentumConservation(const G4Track &, const G4Nucleus &)
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
G4int GetA_asInt() const
Definition: G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
G4double GetWeight() const
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4Material * GetMaterial() const
const G4TouchableHandle & GetTouchableHandle() const
void SetCreatorModelID(const G4int id)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeWeight(G4double finalWeight)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)

◆ AtRestGetPhysicalInteractionLength()

G4double G4HadronStoppingProcess::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 115 of file G4HadronStoppingProcess.cc.

117{
119 return 0.0;
120}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced

◆ BuildPhysicsTable()

void G4HadronStoppingProcess::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4HadronicProcess.

Definition at line 108 of file G4HadronStoppingProcess.cc.

109{
111}
void PrintInfo(const G4ParticleDefinition *)

◆ GetMeanLifeTime()

G4double G4HadronStoppingProcess::GetMeanLifeTime ( const G4Track ,
G4ForceCondition  
)
inlineprotected

Definition at line 98 of file G4HadronStoppingProcess.hh.

99 { return -1.0; }

◆ IsApplicable()

G4bool G4HadronStoppingProcess::IsApplicable ( const G4ParticleDefinition p)
virtual

◆ PostStepGetPhysicalInteractionLength()

G4double G4HadronStoppingProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Reimplemented from G4HadronicProcess.

Definition at line 124 of file G4HadronStoppingProcess.cc.

126{
128 return DBL_MAX;
129}
#define DBL_MAX
Definition: templates.hh:62

◆ PreparePhysicsTable()

void G4HadronStoppingProcess::PreparePhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4HadronicProcess.

Definition at line 98 of file G4HadronStoppingProcess.cc.

99{
101 emcID = G4PhysicsModelCatalog::GetModelID(G4String("model_" + (GetProcessName() + "_EMCascade")));
102 ncID = G4PhysicsModelCatalog::GetModelID(G4String("model_" + (GetProcessName() + "_NuclearCapture")));
103 dioID = G4PhysicsModelCatalog::GetModelID(G4String("model_" + (GetProcessName() + "_DIO")));
104}
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)
static G4int GetModelID(const G4int modelIndex)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386

◆ ProcessDescription()

void G4HadronStoppingProcess::ProcessDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4HadronicProcess.

Reimplemented in G4HadronicAbsorptionBertini, G4HadronicAbsorptionFritiof, G4HadronicAbsorptionFritiofWithBinaryCascade, and G4MuonMinusCapture.

Definition at line 300 of file G4HadronStoppingProcess.cc.

301{
302 outFile << "Base process for negatively charged particle capture at rest.\n";
303}

◆ SetBoundDecay()

void G4HadronStoppingProcess::SetBoundDecay ( G4HadronicInteraction ptr)
inline

Definition at line 137 of file G4HadronStoppingProcess.hh.

138{
139 fBoundDecay = ptr;
140}

Referenced by G4MuonMinusCapture::G4MuonMinusCapture().

◆ SetElementSelector()

void G4HadronStoppingProcess::SetElementSelector ( G4ElementSelector ptr)
inline

Definition at line 122 of file G4HadronStoppingProcess.hh.

123{
124 if(fElementSelector != ptr) {
125 delete fElementSelector;
126 fElementSelector = ptr;
127 }
128}

◆ SetEmCascade()

void G4HadronStoppingProcess::SetEmCascade ( G4HadronicInteraction ptr)
inline

Definition at line 131 of file G4HadronStoppingProcess.hh.

132{
133 fEmCascade = ptr;
134}

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