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

#include <G4PiMinusAbsorptionAtRest.hh>

+ Inheritance diagram for G4PiMinusAbsorptionAtRest:

Public Member Functions

 G4PiMinusAbsorptionAtRest (const G4String &processName="PiMinusAbsorptionAtRest", G4ProcessType aType=fHadronic)
 
 ~G4PiMinusAbsorptionAtRest ()
 
G4bool IsApplicable (const G4ParticleDefinition &particle)
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4VParticleChangeAtRestDoIt (const G4Track &aTrack, const G4Step &aStep)
 
void SetDeexcitationAlgorithm (G4int index)
 
- Public Member Functions inherited from G4VRestProcess
 G4VRestProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestProcess (G4VRestProcess &)
 
virtual ~G4VRestProcess ()
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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
 
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 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
 

Protected Member Functions

G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *)
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)=0
 
- 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 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 54 of file G4PiMinusAbsorptionAtRest.hh.

Constructor & Destructor Documentation

◆ G4PiMinusAbsorptionAtRest()

G4PiMinusAbsorptionAtRest::G4PiMinusAbsorptionAtRest ( const G4String processName = "PiMinusAbsorptionAtRest",
G4ProcessType  aType = fHadronic 
)

Definition at line 65 of file G4PiMinusAbsorptionAtRest.cc.

66 :
67 G4VRestProcess (processName, aType)
68{
69 G4HadronicDeprecate("G4PiMinusAbsorptionAtRest");
70
72
73 _indexDeexcitation = 0;
74
75 if (verboseLevel>0)
76 { G4cout << GetProcessName() << " is created "<< G4endl; }
78}
#define G4HadronicDeprecate(name)
@ fHadronAtRest
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void RegisterExtraProcess(G4VProcess *)
static G4HadronicProcessStore * Instance()
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ ~G4PiMinusAbsorptionAtRest()

G4PiMinusAbsorptionAtRest::~G4PiMinusAbsorptionAtRest ( )

Definition at line 83 of file G4PiMinusAbsorptionAtRest.cc.

Member Function Documentation

◆ AtRestDoIt()

G4VParticleChange * G4PiMinusAbsorptionAtRest::AtRestDoIt ( const G4Track aTrack,
const G4Step aStep 
)
virtual

Reimplemented from G4VRestProcess.

Definition at line 98 of file G4PiMinusAbsorptionAtRest.cc.

99{
100 const G4DynamicParticle* stoppedHadron = track.GetDynamicParticle();
101
102 // Check applicability
103 if (! IsApplicable(*(stoppedHadron->GetDefinition())))
104 {
105 G4cerr
106 << "G4PiMinusAbsorptionAtRest: ERROR, particle must be a pion minus!"
107 << G4endl;
108 return NULL;
109 }
110
111 // Get the current material
112 const G4Material* material = track.GetMaterial();
113
114 G4double A=-1;
115 G4double Z=-1;
116 G4double random = G4UniformRand();
117 const G4ElementVector* theElementVector = material->GetElementVector();
118 unsigned int i;
119 G4double sum = 0;
120 G4double totalsum=0;
121 for(i=0; i<material->GetNumberOfElements(); ++i)
122 {
123 if((*theElementVector)[i]->GetZ()!=1) totalsum+=material->GetFractionVector()[i];
124 }
125 for (i = 0; i<material->GetNumberOfElements(); ++i)
126 {
127 if((*theElementVector)[i]->GetZ()!=1) sum += material->GetFractionVector()[i];
128 if ( sum/totalsum > random )
129 {
130 A = (*theElementVector)[i]->GetA()*mole/g;
131 Z = (*theElementVector)[i]->GetZ();
132 break;
133 }
134 }
135
136 // Do the interaction with the nucleon cluster
137
138 G4PiMinusStopMaterial* algorithm = LoadAlgorithm(static_cast<G4int>(Z));
139 G4PiMinusStopAbsorption stopAbsorption(algorithm,Z,A);
140 stopAbsorption.SetVerboseLevel(verboseLevel);
141
142 G4DynamicParticleVector* absorptionProducts = stopAbsorption.DoAbsorption();
143
144 // Deal with the leftover nucleus
145
146 G4double pionEnergy = stoppedHadron->GetTotalEnergy();
147 G4double excitation = pionEnergy - stopAbsorption.Energy();
148 if (excitation < 0.)
149 {
150 G4Exception("G4PiMinusAbsorptionAtRest::AtRestDoIt()", "HAD_STOP_0000",
151 FatalException, "Excitation energy < 0");
152 }
153 if (verboseLevel>0) { G4cout << " excitation " << excitation << G4endl; }
154
155 G4StopDeexcitationAlgorithm* nucleusAlgorithm = LoadNucleusAlgorithm();
156 G4StopDeexcitation stopDeexcitation(nucleusAlgorithm);
157
158 G4double newZ = Z - stopAbsorption.NProtons();
159 G4double newN = A - Z - stopAbsorption.NNeutrons();
160 G4double newA = newZ + newN;
161 G4ReactionProductVector* fragmentationProducts = stopDeexcitation.DoBreakUp(newA,newZ,excitation,stopAbsorption.RecoilMomentum());
162
163 unsigned int nAbsorptionProducts = 0;
164 if (absorptionProducts != 0)
165 { nAbsorptionProducts = absorptionProducts->size(); }
166
167 unsigned int nFragmentationProducts = 0;
168 if (fragmentationProducts != 0)
169 { nFragmentationProducts = fragmentationProducts->size(); }
170
171 if (verboseLevel>0)
172 {
173 G4cout << "nAbsorptionProducts = " << nAbsorptionProducts
174 << " nFragmentationProducts = " << nFragmentationProducts
175 << G4endl;
176 }
177
178 // Deal with ParticleChange final state
179
181 aParticleChange.SetNumberOfSecondaries(G4int(nAbsorptionProducts + nFragmentationProducts));
182
183 for (i = 0; i<nAbsorptionProducts; i++)
184 { aParticleChange.AddSecondary((*absorptionProducts)[i]); }
185
186// for (i = 0; i<nFragmentationProducts; i++)
187// { aParticleChange.AddSecondary(fragmentationProducts->at(i)); }
188 for(i=0; i<nFragmentationProducts; i++)
189 {
190 G4DynamicParticle * aNew =
191 new G4DynamicParticle((*fragmentationProducts)[i]->GetDefinition(),
192 (*fragmentationProducts)[i]->GetMomentum());
193 G4double newTime = aParticleChange.GetGlobalTime((*fragmentationProducts)[i]->GetFormationTime());
194 aParticleChange.AddSecondary(aNew, newTime);
195 delete (*fragmentationProducts)[i];
196 }
197
198 if (fragmentationProducts != 0) delete fragmentationProducts;
199
200 if (_indexDeexcitation == 1) aParticleChange.ProposeLocalEnergyDeposit(excitation);
201
202 // Kill the absorbed pion
204
205 return &aParticleChange;
206
207}
std::vector< G4DynamicParticle * > G4DynamicParticleVector
std::vector< G4Element * > G4ElementVector
@ FatalException
std::vector< G4ReactionProduct * > G4ReactionProductVector
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4DLLIMPORT std::ostream G4cerr
#define G4UniformRand()
Definition: Randomize.hh:53
G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
const G4double * GetFractionVector() const
Definition: G4Material.hh:193
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:576
void AddSecondary(G4Track *aSecondary)
virtual void Initialize(const G4Track &)
G4double GetGlobalTime(G4double timeDelay=0.0) const
G4bool IsApplicable(const G4ParticleDefinition &particle)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ BuildPhysicsTable()

void G4PiMinusAbsorptionAtRest::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 93 of file G4PiMinusAbsorptionAtRest.cc.

94{
96}
void PrintInfo(const G4ParticleDefinition *)

◆ GetMeanLifeTime()

G4double G4PiMinusAbsorptionAtRest::GetMeanLifeTime ( const G4Track aTrack,
G4ForceCondition  
)
inlineprotectedvirtual

Implements G4VRestProcess.

Definition at line 88 of file G4PiMinusAbsorptionAtRest.hh.

90 {
91 G4double result = 0;
92 if(aTrack.GetMaterial()->GetNumberOfElements() == 1)
93 if(aTrack.GetMaterial()->GetZ()<1.5) result = DBL_MAX;
94 return result;
95 }
G4double GetZ() const
Definition: G4Material.cc:604
G4Material * GetMaterial() const
#define DBL_MAX
Definition: templates.hh:83

◆ IsApplicable()

G4bool G4PiMinusAbsorptionAtRest::IsApplicable ( const G4ParticleDefinition particle)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 74 of file G4PiMinusAbsorptionAtRest.hh.

75 { return ( particle == *(G4PionMinus::PionMinus()) ); }
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98

Referenced by AtRestDoIt().

◆ PreparePhysicsTable()

void G4PiMinusAbsorptionAtRest::PreparePhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 88 of file G4PiMinusAbsorptionAtRest.cc.

89{
91}
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)

◆ SetDeexcitationAlgorithm()

void G4PiMinusAbsorptionAtRest::SetDeexcitationAlgorithm ( G4int  index)

Definition at line 287 of file G4PiMinusAbsorptionAtRest.cc.

288{
289 _indexDeexcitation = index;
290}

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