Geant4 9.6.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 &, 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 46 of file G4AntiNeutronAnnihilationAtRest.hh.

Constructor & Destructor Documentation

◆ G4AntiNeutronAnnihilationAtRest()

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

Definition at line 46 of file G4AntiNeutronAnnihilationAtRest.cc.

47 :
48 G4VRestProcess (processName, aType), // initialization
49 massPionMinus(G4PionMinus::PionMinus()->GetPDGMass()/GeV),
50 massPionZero(G4PionZero::PionZero()->GetPDGMass()/GeV),
51 massPionPlus(G4PionPlus::PionPlus()->GetPDGMass()/GeV),
52 massGamma(G4Gamma::Gamma()->GetPDGMass()/GeV),
53 massAntiNeutron(G4AntiNeutron::AntiNeutron()->GetPDGMass()/GeV),
54 massNeutron(G4Neutron::Neutron()->GetPDGMass()/GeV),
55 pdefGamma(G4Gamma::Gamma()),
56 pdefPionPlus(G4PionPlus::PionPlus()),
57 pdefPionZero(G4PionZero::PionZero()),
58 pdefPionMinus(G4PionMinus::PionMinus()),
59 pdefProton(G4Proton::Proton()),
60 pdefNeutron(G4Neutron::Neutron()),
61 pdefAntiNeutron(G4AntiNeutron::AntiNeutron()),
62 pdefDeuteron(G4Deuteron::Deuteron()),
63 pdefTriton(G4Triton::Triton()),
64 pdefAlpha(G4Alpha::Alpha())
65{
66 G4HadronicDeprecate("G4AntiNeutronAnnihilationAtRest");
67 if (verboseLevel>0) {
68 G4cout << GetProcessName() << " is created "<< G4endl;
69 }
74
76}
#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 G4AntiNeutron * AntiNeutron()
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void RegisterExtraProcess(G4VProcess *)
static G4HadronicProcessStore * Instance()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4PionZero * PionZero()
Definition: G4PionZero.cc:104
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Triton * Triton()
Definition: G4Triton.cc:95
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ ~G4AntiNeutronAnnihilationAtRest()

G4AntiNeutronAnnihilationAtRest::~G4AntiNeutronAnnihilationAtRest ( )

Definition at line 80 of file G4AntiNeutronAnnihilationAtRest.cc.

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

Member Function Documentation

◆ AtRestDoIt()

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

Reimplemented from G4VRestProcess.

Definition at line 148 of file G4AntiNeutronAnnihilationAtRest.cc.

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

Reimplemented from G4VRestProcess.

Definition at line 122 of file G4AntiNeutronAnnihilationAtRest.cc.

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

◆ BuildPhysicsTable()

void G4AntiNeutronAnnihilationAtRest::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 93 of file G4AntiNeutronAnnihilationAtRest.cc.

94{
96}
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 109 of file G4AntiNeutronAnnihilationAtRest.cc.

110{
111 return ( ngkine );
112
113}

◆ GetSecondaryKinematics()

G4GHEKinematicsVector * G4AntiNeutronAnnihilationAtRest::GetSecondaryKinematics ( )

Definition at line 116 of file G4AntiNeutronAnnihilationAtRest.cc.

117{
118 return ( &gkin[0] );
119
120}

◆ IsApplicable()

G4bool G4AntiNeutronAnnihilationAtRest::IsApplicable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4VProcess.

Definition at line 100 of file G4AntiNeutronAnnihilationAtRest.cc.

103{
104 return ( &particle == pdefAntiNeutron );
105
106}

◆ PreparePhysicsTable()

void G4AntiNeutronAnnihilationAtRest::PreparePhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 88 of file G4AntiNeutronAnnihilationAtRest.cc.

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

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