Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4AntiNeutronAnnihilationAtRest Class Reference

#include <G4AntiNeutronAnnihilationAtRest.hh>

Inheritance diagram for G4AntiNeutronAnnihilationAtRest:
Collaboration 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 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
 
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 &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- 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 G4AntiNeutronAnnihilationAtRest.hh.

Constructor & Destructor Documentation

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  globalTime = targetAtomicMass = targetCharge = evapEnergy1
77  = evapEnergy3 = 0.0;
78  ngkine = ntot = 0;
79 }
G4int verboseLevel
Definition: G4VProcess.hh:368
static G4HadronicProcessStore * Instance()
#define G4HadronicDeprecate(name)
G4GLOB_DLL std::ostream G4cout
#define MAX_SECONDARIES
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
static G4PionZero * PionZero()
Definition: G4PionZero.cc:108
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
void RegisterExtraProcess(G4VProcess *)
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static constexpr double GeV
Definition: G4SIunits.hh:217
#define G4endl
Definition: G4ios.hh:61
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4AntiNeutron * AntiNeutron()

Here is the call graph for this function:

G4AntiNeutronAnnihilationAtRest::~G4AntiNeutronAnnihilationAtRest ( )

Definition at line 83 of file G4AntiNeutronAnnihilationAtRest.cc.

84 {
86  delete [] pv;
87  delete [] eve;
88  delete [] gkin;
89 }
void DeRegisterExtraProcess(G4VProcess *)
static G4HadronicProcessStore * Instance()

Here is the call graph for this function:

Member Function Documentation

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

Reimplemented from G4VRestProcess.

Definition at line 151 of file G4AntiNeutronAnnihilationAtRest.cc.

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

Here is the call graph for this function:

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

Reimplemented from G4VRestProcess.

Definition at line 125 of file G4AntiNeutronAnnihilationAtRest.cc.

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

Here is the call graph for this function:

void G4AntiNeutronAnnihilationAtRest::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 96 of file G4AntiNeutronAnnihilationAtRest.cc.

97 {
99 }
static G4HadronicProcessStore * Instance()
void PrintInfo(const G4ParticleDefinition *)

Here is the call graph for this function:

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

Implements G4VRestProcess.

Definition at line 72 of file G4AntiNeutronAnnihilationAtRest.hh.

73  {return 0.0;}

Here is the caller graph for this function:

G4int G4AntiNeutronAnnihilationAtRest::GetNumberOfSecondaries ( )

Definition at line 112 of file G4AntiNeutronAnnihilationAtRest.cc.

113 {
114  return ( ngkine );
115 
116 }
G4GHEKinematicsVector * G4AntiNeutronAnnihilationAtRest::GetSecondaryKinematics ( )

Definition at line 119 of file G4AntiNeutronAnnihilationAtRest.cc.

120 {
121  return ( &gkin[0] );
122 
123 }
G4bool G4AntiNeutronAnnihilationAtRest::IsApplicable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4VProcess.

Definition at line 103 of file G4AntiNeutronAnnihilationAtRest.cc.

106 {
107  return ( &particle == pdefAntiNeutron );
108 
109 }
void G4AntiNeutronAnnihilationAtRest::PreparePhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 91 of file G4AntiNeutronAnnihilationAtRest.cc.

92 {
94 }
static G4HadronicProcessStore * Instance()
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)

Here is the call graph for this function:


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