Belle II Software development
G4LongLivedNeutralDecay Class Reference

This class is a decay process. More...

#include <G4LongLivedNeutralDecay.h>

Inheritance diagram for G4LongLivedNeutralDecay:

Public Member Functions

 G4LongLivedNeutralDecay (const G4String &processName="LongLivedNeutralDecay")
 Constructor.
 
virtual ~G4LongLivedNeutralDecay ()
 Destructor.
 
virtual G4VParticleChange * PostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 G4VProcess::PostStepDoIt() implemention.
 
virtual G4VParticleChange * AtRestDoIt (const G4Track &aTrack, const G4Step &aStep) override
 G4VProcess::AtRestDoIt() implementation for decay at rest.
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) override
 In G4Decay, thePhysicsTable stores values of beta * std::sqrt( 1 - beta*beta) as a function of normalized kinetic enregy (=Ekin/mass), becasuse this table is universal for all particle types.
 
virtual G4bool IsApplicable (const G4ParticleDefinition &) override
 returns "true" if the decay process can be applied to the particle type.
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition) override
 G4VProcess::AtRestGetPhysicalInteractionLength() implemention.
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 G4VProcess::PostStepGetPhysicalInteractionLength() implemention.
 
virtual void StartTracking (G4Track *) override
 inform Start of tracking for each track to the physics process
 
virtual void EndTracking () override
 inform End of tracking for each track to the physics process
 
G4double GetRemainderLifeTime () const
 Get Remainder of life time at rest decay.
 

Protected Member Functions

virtual G4VParticleChange * DecayIt (const G4Track &aTrack, const G4Step &aStep)
 G4Decay::DecayIt() implemention.
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
 GetMeanFreePath returns ctau*beta*gamma for decay in flight.
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition) override
 GetMeanLifeTime returns ctau for decay at rest.
 

Protected Attributes

G4int verboseLevel
 controle flag for output message 0: Silent 1: Warning message 2: More
 
const G4double HighestValue
 Remainder of life time at rest.
 
G4double fRemainderLifeTime
 ParticleChange for decay process.
 
G4ParticleChangeForDecay fParticleChangeForDecay
 ParticleChange for decay process.
 

Private Member Functions

 G4LongLivedNeutralDecay (const G4LongLivedNeutralDecay &right)
 Copy Constructor.
 
G4LongLivedNeutralDecayoperator= (const G4LongLivedNeutralDecay &right)
 Assignment Operator.
 

Detailed Description

This class is a decay process.

Definition at line 30 of file G4LongLivedNeutralDecay.h.

Constructor & Destructor Documentation

◆ G4LongLivedNeutralDecay() [1/2]

G4LongLivedNeutralDecay ( const G4String &  processName = "LongLivedNeutralDecay")

Constructor.

Parameters
processNamedecay process name

Definition at line 28 of file G4LongLivedNeutralDecay.cc.

29 : G4VRestDiscreteProcess(processName, fDecay),
30 verboseLevel(1),
31 HighestValue(20.0),
33
34{
35 // set Process Sub Type
36 SetProcessSubType(static_cast<int>(DECAY));
37
38#ifdef G4VERBOSE
39 if (GetVerboseLevel() > 1) {
40 G4cout << "G4LongLivedNeutralDecay constructor " << " Name:" << processName << G4endl;
41 }
42#endif
43
44 pParticleChange = &fParticleChangeForDecay;
45}
G4ParticleChangeForDecay fParticleChangeForDecay
ParticleChange for decay process.
G4double fRemainderLifeTime
ParticleChange for decay process.
const G4double HighestValue
Remainder of life time at rest.
G4int verboseLevel
controle flag for output message 0: Silent 1: Warning message 2: More

◆ ~G4LongLivedNeutralDecay()

Destructor.

Definition at line 47 of file G4LongLivedNeutralDecay.cc.

48{
49
50}

◆ G4LongLivedNeutralDecay() [2/2]

Copy Constructor.

Parameters
rightcopy reference

Member Function Documentation

◆ AtRestGetPhysicalInteractionLength()

G4double AtRestGetPhysicalInteractionLength ( const G4Track &  track,
G4ForceCondition *  condition 
)
overridevirtual

G4VProcess::AtRestGetPhysicalInteractionLength() implemention.

Parameters
track
condition

No operation in AtRestDoIt.

Definition at line 263 of file G4LongLivedNeutralDecay.cc.

267{
268 // condition is set to "Not Forced"
269 *condition = NotForced;
270 G4double pTime = track.GetDynamicParticle()->GetPrimaryParticle()->GetProperTime();
271 if (pTime >= 0.) {
272 // general case for LongLivedNeutral
273 fRemainderLifeTime = pTime - track.GetProperTime();
274 if (fRemainderLifeTime <= 0.0) fRemainderLifeTime = DBL_MIN;
275 } else {
276 // LongLivedNeutral case: set stable if there is no pre-assigned decay time
277 fRemainderLifeTime = 1e24 * s;
278 }
279 return fRemainderLifeTime;
280}

◆ BuildPhysicsTable()

void BuildPhysicsTable ( const G4ParticleDefinition &  )
overridevirtual

In G4Decay, thePhysicsTable stores values of beta * std::sqrt( 1 - beta*beta) as a function of normalized kinetic enregy (=Ekin/mass), becasuse this table is universal for all particle types.

Definition at line 83 of file G4LongLivedNeutralDecay.cc.

84{
85 return;
86}

◆ DecayIt()

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

G4Decay::DecayIt() implemention.

Parameters
aTrack
aStep

The DecayIt() method returns by pointer a particle-change object, which has information of daughter particles.

Definition at line 88 of file G4LongLivedNeutralDecay.cc.

89{
90 // The DecayIt() method returns by pointer a particle-change object.
91 // Units are expressed in GEANT4 internal units.
92 // Initialize ParticleChange
93 // all members of G4VParticleChange are set to equal to
94 // corresponding member in G4Track
95//
96
97 fParticleChangeForDecay.Initialize(aTrack);
98
99 // LongLivedNeutral: work only with DynamicParticle, not ParticleDefinition
100 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
101 const G4PrimaryParticle* aPrimaryParticle = aParticle->GetPrimaryParticle();
102
103 //check if thePreAssignedDecayProducts exists
104 const G4DecayProducts* o_products = (aParticle->GetPreAssignedDecayProducts());
105 G4bool isPreAssigned = (o_products != nullptr);
106 G4DecayProducts* products = nullptr;
107
108 if (isPreAssigned) {
109 // default: copy decay products
110 products = new G4DecayProducts(*o_products);
111 } else {
112 // Should not happen, since interaction and free path length are set to max
113 G4ExceptionDescription ed;
114 ed << "LongLivedNeutral particle:"
115 << " decay probability exist but no pre-assigned decay products defined "
116 << "- the particle will be killed;\n";
117 G4Exception("G4LongLivedNeutralDecay::DecayIt ",
118 "DECAY101", JustWarning, ed);
119
120 fParticleChangeForDecay.SetNumberOfSecondaries(0);
121 // Kill the parent particle
122 fParticleChangeForDecay.ProposeTrackStatus(fStopAndKill) ;
123 fParticleChangeForDecay.ProposeLocalEnergyDeposit(0.0);
124
125 ClearNumberOfInteractionLengthLeft();
127 }
128
129 // get parent particle information ...................................
130 G4double ParentEnergy = aPrimaryParticle->GetTotalEnergy();
131 G4double ParentMass = aPrimaryParticle->GetMass();
132 if (ParentEnergy < ParentMass) {
133 G4ExceptionDescription ed;
134 ed << "Total Energy is less than its mass - increased the energy"
135 << "\n Particle: " << "LongLivedNeutral"
136 << "\n Energy:" << ParentEnergy / MeV << "[MeV]"
137 << "\n Mass:" << ParentMass / MeV << "[MeV]";
138 G4Exception("G4LongLivedNeutralDecay::DecayIt ",
139 "DECAY102", JustWarning, ed);
140 ParentEnergy = ParentMass;
141 }
142
143 G4ThreeVector ParentDirection(aPrimaryParticle->GetMomentumDirection());
144
145 // assign a new parent to the products from info of the primary particle
146 const G4DynamicParticle dParticle(aPrimaryParticle->GetParticleDefinition(),
147 aPrimaryParticle->GetMomentumDirection(),
148 aPrimaryParticle->GetKineticEnergy(),
149 aPrimaryParticle->GetMass());
150
151 products->SetParentParticle(dParticle);
152 //boost all decay products to laboratory frame
153 G4double energyDeposit = 0.0;
154 G4double finalGlobalTime = aTrack.GetGlobalTime();
155 G4double finalLocalTime = aTrack.GetLocalTime();
156 if (aTrack.GetTrackStatus() == fStopButAlive) {
157 // AtRest case
158 finalGlobalTime += fRemainderLifeTime;
159 finalLocalTime += fRemainderLifeTime;
160 energyDeposit += aPrimaryParticle->GetKineticEnergy();
161 products->Boost(ParentEnergy, ParentDirection);
162 } else {
163 // default for LongLivedNeutral: PostStep case (decay in flight)
164 products->Boost(ParentEnergy, ParentDirection);
165 }
166
167 //add products in fParticleChangeForDecay
168 G4int numberOfSecondaries = products->entries();
169 fParticleChangeForDecay.SetNumberOfSecondaries(numberOfSecondaries);
170#ifdef G4VERBOSE
171 if (GetVerboseLevel() > 1) {
172 G4cout << "G4LongLivedNeutralDecay::DoIt : Decay vertex :";
173 G4cout << " Time: " << finalGlobalTime / ns << "[ns]";
174 G4cout << " proper time: " << finalLocalTime / ns << "[ns]";
175 G4cout << " X:" << (aTrack.GetPosition()).getX() / cm << "[cm]";
176 G4cout << " Y:" << (aTrack.GetPosition()).getY() / cm << "[cm]";
177 G4cout << " Z:" << (aTrack.GetPosition()).getZ() / cm << "[cm]";
178 G4cout << G4endl;
179 G4cout << "G4LongLivedNeutralDecay::DoIt : decay products in Lab. Frame" << G4endl;
180 products->DumpInfo();
181 }
182#endif
183
184 G4int index;
185 G4ThreeVector currentPosition;
186 const G4TouchableHandle thand = aTrack.GetTouchableHandle();
187 for (index = 0; index < numberOfSecondaries; index++) {
188 // get current position of the track
189 currentPosition = aTrack.GetPosition();
190 // create a new track object
191 G4Track* secondary = new G4Track(products->PopProducts(),
192 finalGlobalTime,
193 currentPosition);
194 // switch on good for tracking flag
195 secondary->SetGoodForTrackingFlag();
196 secondary->SetTouchableHandle(thand);
197 // add the secondary track in the List
198 fParticleChangeForDecay.AddSecondary(secondary);
199 }
200 delete products;
201
202 // Kill the parent particle
203 fParticleChangeForDecay.ProposeTrackStatus(fStopAndKill) ;
204 fParticleChangeForDecay.ProposeLocalEnergyDeposit(energyDeposit);
205 fParticleChangeForDecay.ProposeLocalTime(finalLocalTime);
206
207 // Clear NumberOfInteractionLengthLeft
208 ClearNumberOfInteractionLengthLeft();
209
211}

◆ EndTracking()

void EndTracking ( )
overridevirtual

inform End of tracking for each track to the physics process

Definition at line 223 of file G4LongLivedNeutralDecay.cc.

224{
225 // Clear NumberOfInteractionLengthLeft
226 ClearNumberOfInteractionLengthLeft();
227
228 currentInteractionLength = -1.0;
229}

◆ GetMeanFreePath()

G4double GetMeanFreePath ( const G4Track &  aTrack,
G4double  previousStepSize,
G4ForceCondition *  condition 
)
overrideprotectedvirtual

GetMeanFreePath returns ctau*beta*gamma for decay in flight.

Definition at line 72 of file G4LongLivedNeutralDecay.cc.

73{
74
75 const G4DynamicParticle* dynamicParticle = aTrack.GetDynamicParticle();
76 const G4PrimaryParticle* primaryParticle = dynamicParticle->GetPrimaryParticle();
77 G4double pTime = primaryParticle->GetProperTime();
78 G4double aMass = primaryParticle->GetMass();
79 G4double aMomentum = primaryParticle->GetTotalMomentum();
80 return c_light * pTime * aMomentum / aMass;
81}

◆ GetMeanLifeTime()

G4double GetMeanLifeTime ( const G4Track &  aTrack,
G4ForceCondition *  condition 
)
overrideprotectedvirtual

GetMeanLifeTime returns ctau for decay at rest.

Definition at line 62 of file G4LongLivedNeutralDecay.cc.

64{
65 const G4DynamicParticle* dynamicParticle = aTrack.GetDynamicParticle();
66 const G4PrimaryParticle* primaryParticle = dynamicParticle->GetPrimaryParticle();
67 G4double pTime = primaryParticle->GetProperTime();
68
69 return pTime;
70}

◆ IsApplicable()

G4bool IsApplicable ( const G4ParticleDefinition &  aParticleType)
overridevirtual

returns "true" if the decay process can be applied to the particle type.

Definition at line 52 of file G4LongLivedNeutralDecay.cc.

53{
54 // original G4Decay checks for particle definition mass and lifetime, which are zero for generality
55 // the correct mass and lifetime are set for the G4PrimaryParticle
56 if (!aParticleType.GetPDGStable()) {
57 return true;
58 } else return false;
59}

◆ operator=()

G4LongLivedNeutralDecay & operator= ( const G4LongLivedNeutralDecay right)
private

Assignment Operator.

Parameters
rightassign reference

◆ PostStepDoIt()

G4VParticleChange * PostStepDoIt ( const G4Track &  aTrack,
const G4Step &  aStep 
)
overridevirtual

G4VProcess::PostStepDoIt() implemention.

Parameters
aTrack
aStep

for decay in flight

Definition at line 283 of file G4LongLivedNeutralDecay.cc.

287{
288 if ((aTrack.GetTrackStatus() == fStopButAlive) ||
289 (aTrack.GetTrackStatus() == fStopAndKill)) {
290 fParticleChangeForDecay.Initialize(aTrack);
292 } else {
293 return DecayIt(aTrack, aStep);
294 }
295}
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
G4Decay::DecayIt() implemention.

◆ PostStepGetPhysicalInteractionLength()

G4double PostStepGetPhysicalInteractionLength ( const G4Track &  track,
G4double  previousStepSize,
G4ForceCondition *  condition 
)
overridevirtual

G4VProcess::PostStepGetPhysicalInteractionLength() implemention.

Parameters
trackThis argument of base function is ignored
previousStepSizeThis argument of base function is ignored
condition

Forces the PostStepDoIt action to be called, so that it can do the relocation if it is needed, but does not limit the step.

Definition at line 232 of file G4LongLivedNeutralDecay.cc.

237{
238 // condition is set to "Not Forced"
239 *condition = NotForced;
240
241 // pre-assigned Decay time
242 G4double pTime = track.GetDynamicParticle()->GetPrimaryParticle()->GetProperTime();
243 if (pTime >= 0.) {
244 // default for LongLived Particle: pre-assigned decay time case
245 fRemainderLifeTime = pTime - track.GetProperTime();
246 if (fRemainderLifeTime <= 0.0) fRemainderLifeTime = 0.0;
247
248 // use normalized kinetic energy (= Ekin/mass)
249 G4double rvalue = c_light * fRemainderLifeTime;
250 const G4DynamicParticle* dynamicParticle = track.GetDynamicParticle();
251 const G4PrimaryParticle* primaryParticle = dynamicParticle->GetPrimaryParticle();
252 G4double aMass = primaryParticle->GetMass();
253 G4double aMomentum = primaryParticle->GetTotalMomentum();
254 rvalue *= aMomentum / aMass;
255 return rvalue;
256
257 } else {
258 // no pre-assigned proper decay time --> stable
259 return previousStepSize * DBL_MAX;
260 }
261}

◆ StartTracking()

void StartTracking ( G4Track *  )
overridevirtual

inform Start of tracking for each track to the physics process

Definition at line 214 of file G4LongLivedNeutralDecay.cc.

215{
216 currentInteractionLength = -1.0;
217 ResetNumberOfInteractionLengthLeft();
218
219
220 fRemainderLifeTime = -1.0;
221}

Member Data Documentation

◆ fParticleChangeForDecay

G4ParticleChangeForDecay fParticleChangeForDecay
protected

ParticleChange for decay process.

Definition at line 180 of file G4LongLivedNeutralDecay.h.

◆ fRemainderLifeTime

G4double fRemainderLifeTime
protected

ParticleChange for decay process.

Definition at line 178 of file G4LongLivedNeutralDecay.h.

◆ HighestValue

const G4double HighestValue
protected

Remainder of life time at rest.

Definition at line 176 of file G4LongLivedNeutralDecay.h.

◆ verboseLevel

G4int verboseLevel
protected

controle flag for output message 0: Silent 1: Warning message 2: More

Definition at line 168 of file G4LongLivedNeutralDecay.h.


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