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

#include <G4AntiNeutronAnnihilationAtRest.hh>

+ Inheritance diagram for G4AntiNeutronAnnihilationAtRest:

Public Member Functions

 G4AntiNeutronAnnihilationAtRest (const G4String &processName="AntiNeutronAnnihilationAtRest", G4ProcessType aType=fHadronic)
 
 ~G4AntiNeutronAnnihilationAtRest ()
 
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 &aName, G4ProcessType aType=fNotDefined)
 
 G4VRestProcess (G4VRestProcess &)
 
virtual ~G4VRestProcess ()
 
G4VRestProcessoperator= (const G4VRestProcess &)=delete
 
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 ()
 
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 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 &)
 

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 prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- 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 46 of file G4AntiNeutronAnnihilationAtRest.hh.

Constructor & Destructor Documentation

◆ G4AntiNeutronAnnihilationAtRest()

G4AntiNeutronAnnihilationAtRest::G4AntiNeutronAnnihilationAtRest ( const G4String processName = "AntiNeutronAnnihilationAtRest",
G4ProcessType  aType = fHadronic 
)

Definition at line 45 of file G4AntiNeutronAnnihilationAtRest.cc.

46 :
47 G4VRestProcess (processName, aType), // initialization
48 massPionMinus(G4PionMinus::PionMinus()->GetPDGMass()/GeV),
49 massPionZero(G4PionZero::PionZero()->GetPDGMass()/GeV),
50 massPionPlus(G4PionPlus::PionPlus()->GetPDGMass()/GeV),
51 massGamma(G4Gamma::Gamma()->GetPDGMass()/GeV),
52 massAntiNeutron(G4AntiNeutron::AntiNeutron()->GetPDGMass()/GeV),
53 massNeutron(G4Neutron::Neutron()->GetPDGMass()/GeV),
54 pdefGamma(G4Gamma::Gamma()),
55 pdefPionPlus(G4PionPlus::PionPlus()),
56 pdefPionZero(G4PionZero::PionZero()),
57 pdefPionMinus(G4PionMinus::PionMinus()),
58 pdefProton(G4Proton::Proton()),
59 pdefNeutron(G4Neutron::Neutron()),
60 pdefAntiNeutron(G4AntiNeutron::AntiNeutron()),
61 pdefDeuteron(G4Deuteron::Deuteron()),
62 pdefTriton(G4Triton::Triton()),
63 pdefAlpha(G4Alpha::Alpha())
64{
65 G4HadronicDeprecate("G4AntiNeutronAnnihilationAtRest");
66 if (verboseLevel>0) {
67 G4cout << GetProcessName() << " is created "<< G4endl;
68 }
73
75 globalTime = targetAtomicMass = targetCharge = evapEnergy1
76 = evapEnergy3 = 0.0;
77 ngkine = ntot = 0;
78}
#define MAX_SECONDARIES
#define G4HadronicDeprecate(name)
@ fHadronAtRest
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
static G4AntiNeutron * AntiNeutron()
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void RegisterExtraProcess(G4VProcess *)
static G4HadronicProcessStore * Instance()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:97
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:97
static G4PionZero * PionZero()
Definition: G4PionZero.cc:107
static G4Proton * Proton()
Definition: G4Proton.cc:92
static G4Triton * Triton()
Definition: G4Triton.cc:94
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

◆ ~G4AntiNeutronAnnihilationAtRest()

G4AntiNeutronAnnihilationAtRest::~G4AntiNeutronAnnihilationAtRest ( )

Definition at line 82 of file G4AntiNeutronAnnihilationAtRest.cc.

83{
85 delete [] pv;
86 delete [] eve;
87 delete [] gkin;
88}
void DeRegisterExtraProcess(G4VProcess *)

Member Function Documentation

◆ AtRestDoIt()

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

Reimplemented from G4VRestProcess.

Definition at line 150 of file G4AntiNeutronAnnihilationAtRest.cc.

159{
160
161// Initialize ParticleChange
162// all members of G4VParticleChange are set to equal to
163// corresponding member in G4Track
164
166
167// Store some global quantities that depend on current material and particle
168
169 globalTime = track.GetGlobalTime()/s;
170 G4Material * aMaterial = track.GetMaterial();
171 const G4int numberOfElements = aMaterial->GetNumberOfElements();
172 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
173
174 const G4double* theAtomicNumberDensity = aMaterial->GetAtomicNumDensityVector();
175 G4double normalization = 0;
176 for ( G4int i1=0; i1 < numberOfElements; i1++ )
177 {
178 normalization += theAtomicNumberDensity[i1] ; // change when nucleon specific
179 // probabilities are included.
180 }
181 G4double runningSum= 0.;
182 G4double random = G4UniformRand()*normalization;
183 for ( G4int i2=0; i2 < numberOfElements; i2++ )
184 {
185 runningSum += theAtomicNumberDensity[i2]; // change when nucleon specific
186 // probabilities are included.
187 if (random<=runningSum)
188 {
189 targetCharge = G4double( ((*theElementVector)[i2])->GetZ());
190 targetAtomicMass = (*theElementVector)[i2]->GetN();
191 }
192 }
193 if (random>runningSum)
194 {
195 targetCharge = G4double( ((*theElementVector)[numberOfElements-1])->GetZ());
196 targetAtomicMass = (*theElementVector)[numberOfElements-1]->GetN();
197 }
198
199 if (verboseLevel>1) {
200 G4cout << "G4AntiNeutronAnnihilationAtRest::AtRestDoIt is invoked " <<G4endl;
201 }
202
203 G4ParticleMomentum momentum;
204 G4float localtime;
205
207
208 GenerateSecondaries(); // Generate secondaries
209
211
212 for ( G4int isec = 0; isec < ngkine; isec++ ) {
213 G4DynamicParticle* aNewParticle = new G4DynamicParticle;
214 aNewParticle->SetDefinition( gkin[isec].GetParticleDef() );
215 aNewParticle->SetMomentum( gkin[isec].GetMomentum() * GeV );
216
217 localtime = globalTime + gkin[isec].GetTOF();
218
219 G4Track* aNewTrack = new G4Track( aNewParticle, localtime*s, position );
220 aNewTrack->SetTouchableHandle(track.GetTouchableHandle());
221 aParticleChange.AddSecondary( aNewTrack );
222
223 }
224
226
227 aParticleChange.ProposeTrackStatus(fStopAndKill); // Kill the incident AntiNeutron
228
229// clear InteractionLengthLeft
230
232
233 return &aParticleChange;
234
235}
std::vector< G4Element * > G4ElementVector
@ fStopAndKill
float G4float
Definition: G4Types.hh:84
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetMomentum(const G4ThreeVector &momentum)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:214
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:80
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:327

◆ AtRestGetPhysicalInteractionLength()

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

Reimplemented from G4VRestProcess.

Definition at line 124 of file G4AntiNeutronAnnihilationAtRest.cc.

128{
129 // beggining of tracking
131
132 // condition is set to "Not Forced"
134
135 // get mean life time
137
138 if ((currentInteractionLength <0.0) || (verboseLevel>2)){
139 G4cout << "G4AntiNeutronAnnihilationAtRestProcess::AtRestGetPhysicalInteractionLength ";
140 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
141 track.GetDynamicParticle()->DumpInfo();
142 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
143 G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl;
144 }
145
147
148}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
G4double GetMeanLifeTime(const G4Track &, G4ForceCondition *)
void DumpInfo(G4int mode=0) const
const G4String & GetName() const
Definition: G4Material.hh:175
const G4DynamicParticle * GetDynamicParticle() const
G4double currentInteractionLength
Definition: G4VProcess.hh:335
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:331
#define ns
Definition: xmlparse.cc:614

◆ BuildPhysicsTable()

void G4AntiNeutronAnnihilationAtRest::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 95 of file G4AntiNeutronAnnihilationAtRest.cc.

96{
98}
void PrintInfo(const G4ParticleDefinition *)

◆ GetMeanLifeTime()

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

Implements G4VRestProcess.

Definition at line 71 of file G4AntiNeutronAnnihilationAtRest.hh.

72 {return 0.0;}

Referenced by AtRestGetPhysicalInteractionLength().

◆ GetNumberOfSecondaries()

G4int G4AntiNeutronAnnihilationAtRest::GetNumberOfSecondaries ( )

Definition at line 111 of file G4AntiNeutronAnnihilationAtRest.cc.

112{
113 return ( ngkine );
114
115}

◆ GetSecondaryKinematics()

G4GHEKinematicsVector * G4AntiNeutronAnnihilationAtRest::GetSecondaryKinematics ( )

Definition at line 118 of file G4AntiNeutronAnnihilationAtRest.cc.

119{
120 return ( &gkin[0] );
121
122}

◆ IsApplicable()

G4bool G4AntiNeutronAnnihilationAtRest::IsApplicable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4VProcess.

Definition at line 102 of file G4AntiNeutronAnnihilationAtRest.cc.

105{
106 return ( &particle == pdefAntiNeutron );
107
108}

◆ PreparePhysicsTable()

void G4AntiNeutronAnnihilationAtRest::PreparePhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 90 of file G4AntiNeutronAnnihilationAtRest.cc.

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

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