Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions
G4HadronStoppingProcess Class Reference

#include <G4HadronStoppingProcess.hh>

Inheritance diagram for G4HadronStoppingProcess:
G4HadronicProcess G4VDiscreteProcess G4VProcess G4HadronicAbsorptionBertini G4HadronicAbsorptionFritiof G4MuonMinusCapture G4KaonMinusAbsorptionBertini G4PiMinusAbsorptionBertini G4SigmaMinusAbsorptionBertini G4AntiProtonAbsorptionFritiof G4AntiSigmaPlusAbsorptionFritiof

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)
 
virtual ~G4HadronicProcess ()
 
void RegisterMe (G4HadronicInteraction *a)
 
G4double GetElementCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=0)
 
G4double GetMicroscopicCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=0)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
void DumpPhysicsTable (const G4ParticleDefinition &p)
 
void AddDataSet (G4VCrossSectionDataSet *aDataSet)
 
G4EnergyRangeManagerGetManagerPointer ()
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *)
 
const G4NucleusGetTargetNucleus () const
 
const G4IsotopeGetTargetIsotope ()
 
void BiasCrossSectionByFactor (G4double aScale)
 
void SetEpReportLevel (G4int level)
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
G4CrossSectionDataStoreGetCrossSectionDataStore ()
 
void MultiplyCrossSectionBy (G4double factor)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
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 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 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
 
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 (G4double kineticEnergy, G4Material *aMaterial, G4Element *anElement)
 
G4NucleusGetTargetNucleusPointer ()
 
void DumpState (const G4Track &, const G4String &, G4ExceptionDescription &)
 
const G4EnergyRangeManagerGetEnergyRangeManager () const
 
void SetEnergyRangeManager (const G4EnergyRangeManager &value)
 
G4HadronicInteractionGetHadronicInteraction () const
 
G4double GetLastCrossSection ()
 
void FillResult (G4HadFinalState *aR, const G4Track &aT)
 
G4HadFinalStateCheckResult (const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result) const
 
void CheckEnergyMomentumConservation (const G4Track &, const G4Nucleus &)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4HadronicProcess
G4HadProjectile thePro
 
G4ParticleChangetheTotalResult
 
G4int epReportLevel
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 62 of file G4HadronStoppingProcess.hh.

Constructor & Destructor Documentation

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

Definition at line 63 of file G4HadronStoppingProcess.cc.

References G4VProcess::enableAtRestDoIt, G4VProcess::enablePostStepDoIt, G4HadronicProcessStore::Instance(), and G4HadronicProcessStore::RegisterExtraProcess().

65 {
66  // Modify G4VProcess flags to emulate G4VRest instead of G4VDiscrete
67  enableAtRestDoIt = true;
68  enablePostStepDoIt = false;
69 
70  fElementSelector = new G4ElementSelector();
71  fEmCascade = new G4EmCaptureCascade(); // Owned by InteractionRegistry
72  fBoundDecay = 0;
74 }
static G4HadronicProcessStore * Instance()
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:352
void RegisterExtraProcess(G4VProcess *)
G4HadronicProcess(const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)
G4HadronStoppingProcess::~G4HadronStoppingProcess ( )
virtual

Definition at line 78 of file G4HadronStoppingProcess.cc.

References G4HadronicProcessStore::DeRegisterExtraProcess(), and G4HadronicProcessStore::Instance().

79 {
81  delete fElementSelector;
82  // NOTE: fEmCascade and fEmBoundDecay owned by registry, not locally
83 }
void DeRegisterExtraProcess(G4VProcess *)
static G4HadronicProcessStore * Instance()

Member Function Documentation

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

Reimplemented from G4VDiscreteProcess.

Definition at line 126 of file G4HadronStoppingProcess.cc.

References G4HadFinalState::AddSecondaries(), G4ParticleChange::AddSecondary(), G4HadronicInteraction::ApplyYourself(), G4HadronicProcess::CheckEnergyMomentumConservation(), G4HadronicProcess::CheckResult(), G4HadronicProcess::ChooseHadronicInteraction(), G4HadFinalState::Clear(), G4HadronicProcess::DumpState(), G4HadronicProcess::epReportLevel, FatalException, fStopAndKill, G4endl, G4Exception(), G4Nucleus::GetA_asInt(), G4Track::GetGlobalTime(), G4HadFinalState::GetLocalEnergyDeposit(), G4Track::GetMaterial(), G4HadronicInteraction::GetModelName(), G4Element::GetName(), G4HadFinalState::GetNumberOfSecondaries(), G4HadSecondary::GetParticle(), G4Track::GetPosition(), G4HadFinalState::GetSecondary(), G4HadFinalState::GetStatusChange(), G4HadronicProcess::GetTargetNucleusPointer(), G4HadSecondary::GetTime(), G4Track::GetTouchableHandle(), G4HadSecondary::GetWeight(), G4Track::GetWeight(), G4Nucleus::GetZ_asInt(), G4HadProjectile::Initialise(), G4ParticleChange::Initialize(), n, G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::ProposeTrackStatus(), G4VParticleChange::ProposeWeight(), G4ElementSelector::SelectZandA(), G4HadProjectile::SetBoundEnergy(), G4VParticleChange::SetNumberOfSecondaries(), G4Track::SetTouchableHandle(), G4Track::SetWeight(), stopAndKill, G4HadronicProcess::thePro, and G4HadronicProcess::theTotalResult.

128 {
129  // if primary is not Alive then do nothing
130  theTotalResult->Initialize(track);
131 
132  G4Nucleus* nucleus = GetTargetNucleusPointer();
133  G4Element* elm = fElementSelector->SelectZandA(track, nucleus);
134 
135  G4HadFinalState* result = 0;
136  thePro.Initialise(track);
137  G4double time0 = track.GetGlobalTime();
138  G4bool nuclearCapture = true;
139 
140  // Do the electromagnetic cascade in the nuclear field.
141  // EM cascade should keep G4HadFinalState object,
142  // because it will not be deleted at the end of this method
143  //
144  result = fEmCascade->ApplyYourself(thePro, *nucleus);
145  G4double ebound = result->GetLocalEnergyDeposit();
146  G4double edep = 0.0;
147  G4int nSecondaries = result->GetNumberOfSecondaries();
148 
149  // Try decay from bound level
150  // For mu- the time of projectile should be changed.
151  // Decay should keep G4HadFinalState object,
152  // because it will not be deleted at the end of this method
153  //
154  thePro.SetBoundEnergy(ebound);
155  if(fBoundDecay) {
156  G4HadFinalState* resultDecay =
157  fBoundDecay->ApplyYourself(thePro, *nucleus);
158  G4int n = resultDecay->GetNumberOfSecondaries();
159  if(0 < n) {
160  nSecondaries += n;
161  result->AddSecondaries(resultDecay);
162  }
163  if(resultDecay->GetStatusChange() == stopAndKill) {
164  nuclearCapture = false;
165  }
166  resultDecay->Clear();
167  }
168 
169  if(nuclearCapture) {
170 
171  // select model
173  try {
174  model = ChooseHadronicInteraction(0.0, track.GetMaterial(), elm);
175  }
176  catch(G4HadronicException & aE) {
178  ed << "Target element "<<elm->GetName()<<" Z= "
179  << nucleus->GetZ_asInt() << " A= "
180  << nucleus->GetA_asInt() << G4endl;
181  DumpState(track,"ChooseHadronicInteraction",ed);
182  ed << " No HadronicInteraction found out" << G4endl;
183  G4Exception("G4HadronStoppingProcess::AtRestDoIt", "had005",
184  FatalException, ed);
185  }
186 
187  G4HadFinalState* resultNuc = 0;
188  G4int reentryCount = 0;
189  do {
190  // sample final state
191  // nuclear interaction should keep G4HadFinalState object
192  // model should define time of each secondary particle
193  try {
194  resultNuc = model->ApplyYourself(thePro, *nucleus);
195  ++reentryCount;
196  }
197  catch(G4HadronicException aR) {
199  ed << "Call for " << model->GetModelName() << G4endl;
200  ed << "Target element "<<elm->GetName()<<" Z= "
201  << nucleus->GetZ_asInt()
202  << " A= " << nucleus->GetA_asInt() << G4endl;
203  DumpState(track,"ApplyYourself",ed);
204  ed << " ApplyYourself failed" << G4endl;
205  G4Exception("G4HadronStoppingProcess::AtRestDoIt", "had006",
206  FatalException, ed);
207  }
208 
209  // Check the result for catastrophic energy non-conservation
210  resultNuc = CheckResult(thePro, *nucleus, resultNuc);
211 
212  if(reentryCount>100) {
214  ed << "Call for " << model->GetModelName() << G4endl;
215  ed << "Target element "<<elm->GetName()<<" Z= "
216  << nucleus->GetZ_asInt()
217  << " A= " << nucleus->GetA_asInt() << G4endl;
218  DumpState(track,"ApplyYourself",ed);
219  ed << " ApplyYourself does not completed after 100 attempts" << G4endl;
220  G4Exception("G4HadronStoppingProcess::AtRestDoIt", "had006",
221  FatalException, ed);
222  }
223  }
224  while(!resultNuc);
225 
226  edep = resultNuc->GetLocalEnergyDeposit();
227  nSecondaries += resultNuc->GetNumberOfSecondaries();
228  result->AddSecondaries(resultNuc);
229  resultNuc->Clear();
230  }
231 
232  // Fill results
233  //
236  theTotalResult->SetNumberOfSecondaries(nSecondaries);
237  G4double w = track.GetWeight();
239  for(G4int i=0; i<nSecondaries; ++i) {
240  G4HadSecondary* sec = result->GetSecondary(i);
241 
242  // add track global time to the reaction time
243  G4double time = sec->GetTime();
244  if(time < 0.0) { time = 0.0; }
245  time += time0;
246 
247  // create secondary track
248  G4Track* t = new G4Track(sec->GetParticle(),
249  time,
250  track.GetPosition());
251  t->SetWeight(w*sec->GetWeight());
254  }
255  result->Clear();
256 
257  if (epReportLevel != 0) {
258  CheckEnergyMomentumConservation(track, *nucleus);
259  }
260  return theTotalResult;
261 }
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
virtual G4Element * SelectZandA(const G4Track &track, G4Nucleus *)
void AddSecondaries(const std::vector< G4HadSecondary > &addSecs)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4HadSecondary * GetSecondary(size_t i)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4HadProjectile thePro
void CheckEnergyMomentumConservation(const G4Track &, const G4Nucleus &)
const G4String & GetModelName() const
G4double GetTime() const
int G4int
Definition: G4Types.hh:78
void SetWeight(G4double aValue)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void ProposeWeight(G4double finalWeight)
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
G4ParticleChange * theTotalResult
bool G4bool
Definition: G4Types.hh:79
G4HadronicInteraction * ChooseHadronicInteraction(G4double kineticEnergy, G4Material *aMaterial, G4Element *anElement)
G4HadFinalState * CheckResult(const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result) const
const XML_Char XML_Content * model
const G4int n
G4double GetGlobalTime() const
const G4TouchableHandle & GetTouchableHandle() const
G4Material * GetMaterial() const
void Initialise(const G4Track &aT)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual void Initialize(const G4Track &)
G4Nucleus * GetTargetNucleusPointer()
void SetNumberOfSecondaries(G4int totSecondaries)
G4DynamicParticle * GetParticle()
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
void SetBoundEnergy(G4double e)
void AddSecondary(G4Track *aSecondary)
G4double GetWeight() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
const G4String & GetName() const
Definition: G4Element.hh:127
G4int GetNumberOfSecondaries() const
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
G4double GetLocalEnergyDeposit() const
G4HadFinalStateStatus GetStatusChange() const
G4double GetWeight() const
G4double G4HadronStoppingProcess::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 108 of file G4HadronStoppingProcess.cc.

References NotForced.

110 {
111  *condition = NotForced;
112  return 0.0;
113 }
G4double condition(const G4ErrorSymMatrix &m)
void G4HadronStoppingProcess::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4HadronicProcess.

Definition at line 101 of file G4HadronStoppingProcess.cc.

References G4HadronicProcessStore::Instance(), and G4HadronicProcessStore::PrintInfo().

102 {
104 }
static G4HadronicProcessStore * Instance()
void PrintInfo(const G4ParticleDefinition *)
G4double G4HadronStoppingProcess::GetMeanLifeTime ( const G4Track ,
G4ForceCondition  
)
inlineprotected

Definition at line 98 of file G4HadronStoppingProcess.hh.

99  { return -1.0; }
G4bool G4HadronStoppingProcess::IsApplicable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Reimplemented in G4MuonMinusCapture, G4HadronicAbsorptionFritiof, and G4HadronicAbsorptionBertini.

Definition at line 87 of file G4HadronStoppingProcess.cc.

References G4ParticleDefinition::GetPDGCharge().

Referenced by G4HadronicAbsorptionBertini::IsApplicable().

88 {
89  return (p.GetPDGCharge() < 0.0);
90 }
G4double GetPDGCharge() const
G4double G4HadronStoppingProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 117 of file G4HadronStoppingProcess.cc.

References DBL_MAX, and NotForced.

119 {
120  *condition = NotForced;
121  return DBL_MAX;
122 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
void G4HadronStoppingProcess::PreparePhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4HadronicProcess.

Definition at line 94 of file G4HadronStoppingProcess.cc.

References G4HadronicProcessStore::Instance(), and G4HadronicProcessStore::RegisterParticleForExtraProcess().

95 {
97 }
static G4HadronicProcessStore * Instance()
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)
void G4HadronStoppingProcess::ProcessDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4HadronicProcess.

Reimplemented in G4MuonMinusCapture, G4HadronicAbsorptionFritiof, and G4HadronicAbsorptionBertini.

Definition at line 265 of file G4HadronStoppingProcess.cc.

266 {
267  outFile << "Base process for negatively charged particle capture at rest.\n";
268 }
std::ofstream outFile
Definition: GammaRayTel.cc:68
void G4HadronStoppingProcess::SetBoundDecay ( G4HadronicInteraction ptr)
inline

Definition at line 133 of file G4HadronStoppingProcess.hh.

Referenced by G4MuonMinusCapture::G4MuonMinusCapture().

134 {
135  fBoundDecay = ptr;
136 }
void G4HadronStoppingProcess::SetElementSelector ( G4ElementSelector ptr)
inline

Definition at line 118 of file G4HadronStoppingProcess.hh.

119 {
120  if(fElementSelector != ptr) {
121  delete fElementSelector;
122  fElementSelector = ptr;
123  }
124 }
void G4HadronStoppingProcess::SetEmCascade ( G4HadronicInteraction ptr)
inline

Definition at line 127 of file G4HadronStoppingProcess.hh.

128 {
129  fEmCascade = ptr;
130 }

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