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

#include <G4NeutronCaptureAtRest.hh>

+ Inheritance diagram for G4NeutronCaptureAtRest:

Public Member Functions

 G4NeutronCaptureAtRest (const G4String &processName="NeutronCaptureAtRest", G4ProcessType aType=fHadronic)
 
 ~G4NeutronCaptureAtRest ()
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
G4double GetMeanLifeTime (const G4Track &, G4ForceCondition *)
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
G4int GetNumberOfSecondaries ()
 
G4GHEKinematicsVectorGetSecondaryKinematics ()
 
- 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
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- 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 47 of file G4NeutronCaptureAtRest.hh.

Constructor & Destructor Documentation

◆ G4NeutronCaptureAtRest()

G4NeutronCaptureAtRest::G4NeutronCaptureAtRest ( const G4String processName = "NeutronCaptureAtRest",
G4ProcessType  aType = fHadronic 
)

Definition at line 46 of file G4NeutronCaptureAtRest.cc.

47 :
48 G4VRestProcess (processName, aType), // initialization
49 massProton(G4Proton::Proton()->GetPDGMass()/GeV),
50 massNeutron(G4Neutron::Neutron()->GetPDGMass()/GeV),
51 massElectron(G4Electron::Electron()->GetPDGMass()/GeV),
52 massDeuteron(G4Deuteron::Deuteron()->GetPDGMass()/GeV),
53 massAlpha(G4Alpha::Alpha()->GetPDGMass()/GeV),
54 pdefGamma(G4Gamma::Gamma()),
55 pdefNeutron(G4Neutron::Neutron())
56{
57 G4HadronicDeprecate("G4NeutronCaptureAtRest");
58 if (verboseLevel>0) {
59 G4cout << GetProcessName() << " is created "<< G4endl;
60 }
65
67}
#define MAX_SECONDARIES
#define G4HadronicDeprecate(name)
@ fHadronAtRest
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void RegisterExtraProcess(G4VProcess *)
static G4HadronicProcessStore * Instance()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ ~G4NeutronCaptureAtRest()

G4NeutronCaptureAtRest::~G4NeutronCaptureAtRest ( )

Definition at line 71 of file G4NeutronCaptureAtRest.cc.

72{
74 delete [] pv;
75 delete [] eve;
76 delete [] gkin;
77}
void DeRegisterExtraProcess(G4VProcess *)

Member Function Documentation

◆ AtRestDoIt()

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

Reimplemented from G4VRestProcess.

Definition at line 139 of file G4NeutronCaptureAtRest.cc.

148{
149
150// Initialize ParticleChange
151// all members of G4VParticleChange are set to equal to
152// corresponding member in G4Track
153
155
156// Store some global quantities that depend on current material and particle
157
158 globalTime = track.GetGlobalTime()/s;
159 G4Material * aMaterial = track.GetMaterial();
160 const G4int numberOfElements = aMaterial->GetNumberOfElements();
161 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
162
163 const G4double* theAtomicNumberDensity = aMaterial->GetAtomicNumDensityVector();
164 G4double normalization = 0;
165 for ( G4int i1=0; i1 < numberOfElements; i1++ )
166 {
167 normalization += theAtomicNumberDensity[i1] ; // change when nucleon specific
168 // probabilities are included.
169 }
170 G4double runningSum= 0.;
171 G4double random = G4UniformRand()*normalization;
172 for ( G4int i2=0; i2 < numberOfElements; i2++ )
173 {
174 runningSum += theAtomicNumberDensity[i2]; // change when nucleon specific
175 // probabilities are included.
176 if (random<=runningSum)
177 {
178 targetCharge = G4double((*theElementVector)[i2]->GetZ());
179 targetAtomicMass = (*theElementVector)[i2]->GetN();
180 }
181 }
182 if (random>runningSum)
183 {
184 targetCharge = G4double((*theElementVector)[numberOfElements-1]->GetZ());
185 targetAtomicMass = (*theElementVector)[numberOfElements-1]->GetN();
186
187 }
188
189 if (verboseLevel>1) {
190 G4cout << "G4NeutronCaptureAtRest::AtRestDoIt is invoked " <<G4endl;
191 }
192
193 G4ParticleMomentum momentum;
194 G4float localtime;
195
197
198 GenerateSecondaries(); // Generate secondaries
199
201
202 for ( G4int isec = 0; isec < ngkine; isec++ ) {
203 G4DynamicParticle* aNewParticle = new G4DynamicParticle;
204 aNewParticle->SetDefinition( gkin[isec].GetParticleDef() );
205 aNewParticle->SetMomentum( gkin[isec].GetMomentum() * GeV );
206
207 localtime = globalTime + gkin[isec].GetTOF();
208
209 G4Track* aNewTrack = new G4Track( aNewParticle, localtime*s, position );
210 aNewTrack->SetTouchableHandle(track.GetTouchableHandle());
211 aParticleChange.AddSecondary( aNewTrack );
212
213 }
214
216
217 aParticleChange.ProposeTrackStatus(fStopAndKill); // Kill the incident Neutron
218
219// clear InteractionLengthLeft
220
222
223 return &aParticleChange;
224
225}
std::vector< G4Element * > G4ElementVector
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
float G4float
Definition: G4Types.hh:65
int G4int
Definition: G4Types.hh:66
#define G4UniformRand()
Definition: Randomize.hh:53
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetMomentum(const G4ThreeVector &momentum)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:215
void AddSecondary(G4Track *aSecondary)
virtual void Initialize(const G4Track &)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4Material * GetMaterial() const
const G4TouchableHandle & GetTouchableHandle() const
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:92
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289

◆ AtRestGetPhysicalInteractionLength()

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

Reimplemented from G4VRestProcess.

Definition at line 113 of file G4NeutronCaptureAtRest.cc.

117{
118 // beggining of tracking
120
121 // condition is set to "Not Forced"
123
124 // get mean life time
126
127 if ((currentInteractionLength <0.0) || (verboseLevel>2)){
128 G4cout << "G4NeutronCaptureAtRestProcess::AtRestGetPhysicalInteractionLength ";
129 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
130 track.GetDynamicParticle()->DumpInfo();
131 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
132 G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl;
133 }
134
136
137}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
void DumpInfo(G4int mode=0) const
const G4String & GetName() const
Definition: G4Material.hh:177
G4double GetMeanLifeTime(const G4Track &, G4ForceCondition *)
const G4DynamicParticle * GetDynamicParticle() const
G4double currentInteractionLength
Definition: G4VProcess.hh:297
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
#define ns
Definition: xmlparse.cc:597

◆ BuildPhysicsTable()

void G4NeutronCaptureAtRest::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 84 of file G4NeutronCaptureAtRest.cc.

85{
87}
void PrintInfo(const G4ParticleDefinition *)

◆ GetMeanLifeTime()

G4double G4NeutronCaptureAtRest::GetMeanLifeTime ( const G4Track ,
G4ForceCondition  
)
inlinevirtual

Implements G4VRestProcess.

Definition at line 73 of file G4NeutronCaptureAtRest.hh.

74 {return 0.0;}

Referenced by AtRestGetPhysicalInteractionLength().

◆ GetNumberOfSecondaries()

G4int G4NeutronCaptureAtRest::GetNumberOfSecondaries ( )

Definition at line 100 of file G4NeutronCaptureAtRest.cc.

101{
102 return ( ngkine );
103
104}

◆ GetSecondaryKinematics()

G4GHEKinematicsVector * G4NeutronCaptureAtRest::GetSecondaryKinematics ( )

Definition at line 107 of file G4NeutronCaptureAtRest.cc.

108{
109 return ( &gkin[0] );
110
111}

◆ IsApplicable()

G4bool G4NeutronCaptureAtRest::IsApplicable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4VProcess.

Definition at line 91 of file G4NeutronCaptureAtRest.cc.

94{
95 return ( &particle == pdefNeutron );
96
97}

◆ PreparePhysicsTable()

void G4NeutronCaptureAtRest::PreparePhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 79 of file G4NeutronCaptureAtRest.cc.

80{
82}
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)

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