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

#include <G4UnknownDecay.hh>

+ Inheritance diagram for G4UnknownDecay:

Public Member Functions

 G4UnknownDecay (const G4String &processName="UnknownDecay")
 
virtual ~G4UnknownDecay ()
 
virtual G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) override
 
virtual G4bool IsApplicable (const G4ParticleDefinition &) override
 
virtual void ProcessDescription (std::ostream &outFile) const override
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (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 &)
 

Protected Member Functions

virtual G4VParticleChangeDecayIt (const G4Track &aTrack, const G4Step &aStep)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

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

Constructor & Destructor Documentation

◆ G4UnknownDecay()

G4UnknownDecay::G4UnknownDecay ( const G4String processName = "UnknownDecay")

Definition at line 46 of file G4UnknownDecay.cc.

47 :G4VDiscreteProcess(processName, fDecay),
48 verboseLevel(1),
49 HighestValue(20.0)
50{
51 // set Process Sub Type
52 SetProcessSubType(static_cast<int>(DECAY_Unknown));
53
54#ifdef G4VERBOSE
55 if (GetVerboseLevel()>1) {
56 G4cout << "G4UnknownDecay constructor " << " Name:" << processName << G4endl;
57 }
58#endif
59 pParticleChange = &fParticleChangeForDecay;
60}
@ DECAY_Unknown
@ fDecay
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4int GetVerboseLevel() const
Definition: G4VProcess.hh:418
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321

◆ ~G4UnknownDecay()

G4UnknownDecay::~G4UnknownDecay ( )
virtual

Definition at line 62 of file G4UnknownDecay.cc.

63{
64}

Member Function Documentation

◆ BuildPhysicsTable()

void G4UnknownDecay::BuildPhysicsTable ( const G4ParticleDefinition )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 77 of file G4UnknownDecay.cc.

78{
79 return;
80}

◆ DecayIt()

G4VParticleChange * G4UnknownDecay::DecayIt ( const G4Track aTrack,
const G4Step aStep 
)
protectedvirtual

Definition at line 82 of file G4UnknownDecay.cc.

83{
84 // The DecayIt() method returns by pointer a particle-change object.
85 // Units are expressed in GEANT4 internal units.
86
87 // Initialize ParticleChange
88 // all members of G4VParticleChange are set to equal to
89 // corresponding member in G4Track
90 fParticleChangeForDecay.Initialize(aTrack);
91
92 // get particle
93 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
94
95 //check if thePreAssignedDecayProducts exists
96 const G4DecayProducts* o_products = (aParticle->GetPreAssignedDecayProducts());
97 G4bool isPreAssigned = (o_products != nullptr);
98 G4DecayProducts* products = nullptr;
99
100 if (!isPreAssigned ){
101 fParticleChangeForDecay.SetNumberOfSecondaries(0);
102 // Kill the parent particle
103 fParticleChangeForDecay.ProposeTrackStatus( fStopAndKill ) ;
104 fParticleChangeForDecay.ProposeLocalEnergyDeposit(0.0);
105
107 return &fParticleChangeForDecay ;
108 }
109
110 // copy decay products
111 products = new G4DecayProducts(*o_products);
112
113 // get parent particle information ...................................
114 G4double ParentEnergy = aParticle->GetTotalEnergy();
115 G4double ParentMass = aParticle->GetMass();
116 if (ParentEnergy < ParentMass) {
117 ParentEnergy = ParentMass;
118#ifdef G4VERBOSE
119 if (GetVerboseLevel()>1) {
120 G4cout << "G4UnknownDecay::DoIt : Total Energy is less than its mass" << G4endl;
121 G4cout << " Particle: " << aParticle->GetDefinition()->GetParticleName();
122 G4cout << " Energy:" << ParentEnergy/MeV << "[MeV]";
123 G4cout << " Mass:" << ParentMass/MeV << "[MeV]";
124 G4cout << G4endl;
125 }
126#endif
127 }
128 G4ThreeVector ParentDirection(aParticle->GetMomentumDirection());
129
130 G4double energyDeposit = 0.0;
131 G4double finalGlobalTime = aTrack.GetGlobalTime();
132 //boost all decay products to laboratory frame
133 //if the particle has traveled
134 if(aParticle->GetPreAssignedDecayProperTime()>=0.) {
135 products->Boost( ParentEnergy, ParentDirection);
136 }
137
138 //add products in fParticleChangeForDecay
139 G4int numberOfSecondaries = products->entries();
140 fParticleChangeForDecay.SetNumberOfSecondaries(numberOfSecondaries);
141#ifdef G4VERBOSE
142 if (GetVerboseLevel()>1) {
143 G4cout << "G4UnknownDecay::DoIt : Decay vertex :";
144 G4cout << " Time: " << finalGlobalTime/ns << "[ns]";
145 G4cout << " X:" << (aTrack.GetPosition()).x() /cm << "[cm]";
146 G4cout << " Y:" << (aTrack.GetPosition()).y() /cm << "[cm]";
147 G4cout << " Z:" << (aTrack.GetPosition()).z() /cm << "[cm]";
148 G4cout << G4endl;
149 G4cout << "G4UnknownDecay::DoIt : decay products in Lab. Frame" << G4endl;
150 products->DumpInfo();
151 }
152#endif
153 G4int index;
154 G4ThreeVector currentPosition;
155 const G4TouchableHandle thand = aTrack.GetTouchableHandle();
156 for (index=0; index < numberOfSecondaries; index++){
157 // get current position of the track
158 currentPosition = aTrack.GetPosition();
159 // create a new track object
160 G4Track* secondary = new G4Track( products->PopProducts(),
161 finalGlobalTime ,
162 currentPosition );
163 // switch on good for tracking flag
164 secondary->SetGoodForTrackingFlag();
165 secondary->SetTouchableHandle(thand);
166 // add the secondary track in the List
167 fParticleChangeForDecay.AddSecondary(secondary);
168 }
169 delete products;
170
171 // Kill the parent particle
172 fParticleChangeForDecay.ProposeTrackStatus( fStopAndKill ) ;
173 fParticleChangeForDecay.ProposeLocalEnergyDeposit(energyDeposit);
174 fParticleChangeForDecay.ProposeGlobalTime( finalGlobalTime );
175 // reset NumberOfInteractionLengthLeft
177
178 return &fParticleChangeForDecay ;
179}
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
void DumpInfo() const
G4int entries() const
G4DynamicParticle * PopProducts()
void Boost(G4double totalEnergy, const G4ThreeVector &momentumDirection)
G4double GetMass() const
const G4ThreeVector & GetMomentumDirection() const
const G4DecayProducts * GetPreAssignedDecayProducts() const
G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
G4double GetPreAssignedDecayProperTime() const
virtual void Initialize(const G4Track &)
const G4String & GetParticleName() const
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
void SetGoodForTrackingFlag(G4bool value=true)
void ProposeTrackStatus(G4TrackStatus status)
void AddSecondary(G4Track *aSecondary)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:424
#define ns
Definition: xmlparse.cc:614

Referenced by PostStepDoIt().

◆ GetMeanFreePath()

G4double G4UnknownDecay::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
overrideprotectedvirtual

Implements G4VDiscreteProcess.

Definition at line 72 of file G4UnknownDecay.cc.

73{
74 return 0.0;
75}

◆ IsApplicable()

G4bool G4UnknownDecay::IsApplicable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 66 of file G4UnknownDecay.cc.

67{
68 if(aParticleType.GetParticleName()=="unknown") return true;
69 return false;
70}

Referenced by G4UnknownDecayPhysics::ConstructProcess().

◆ PostStepDoIt()

G4VParticleChange * G4UnknownDecay::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
inlineoverridevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 149 of file G4UnknownDecay.hh.

153{
154 return DecayIt(aTrack, aStep);
155}
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)

◆ PostStepGetPhysicalInteractionLength()

G4double G4UnknownDecay::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
inlineoverridevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 124 of file G4UnknownDecay.hh.

129{
130 // pre-assigned UnknownDecay time
132
133 if (pTime < 0.) pTime = DBL_MIN;
134
135 // condition is set to "Not Forced"
137
138 // reminder proper time
139 G4double remainder = pTime - track.GetProperTime();
140 if (remainder <= 0.0) remainder = DBL_MIN;
141
142 // use pre-assigned Decay time to determine PIL
143 //return GetMeanFreePath(track, previousStepSize, condition);
144 return remainder*CLHEP::c_light;
145
146}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
G4double GetProperTime() const
#define DBL_MIN
Definition: templates.hh:54

◆ ProcessDescription()

void G4UnknownDecay::ProcessDescription ( std::ostream &  outFile) const
overridevirtual

Reimplemented from G4VProcess.

Definition at line 181 of file G4UnknownDecay.cc.

182{
183 outFile << GetProcessName()
184 << ": Decay of 'unknown' particles. \n"
185 << "kinematics of daughters are dertermined "
186 << "by PreAssignedDecayProducts. \n";
187}
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

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