Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4Decay.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//
28//
29// --------------------------------------------------------------
30// GEANT 4 class implementation file
31//
32// History: first implementation, based on object model of
33// 2nd December 1995, G.Cosmo
34// 7 July 1996 H.Kurashige
35// ------------------------------------------------------------
36// remove BuildPhysicsTable() 28 Nov. 1997 H.Kurashige
37// change DBL_EPSIRON to DBL_MIN 14 Dec. 1997 H.Kurashige
38// modified for new ParticleChange 12 Mar. 1998 H.Kurashige
39// modified for "GoodForTrackingFlag" 19 June 1998 H.Kurashige
40// rename thePhysicsTable to aPhyscisTable 2 Aug. 1998 H.Kurashige
41// modified IsApplicable in order to protect the decay from registered
42// to resonances 12 Dec. 1998 H.Kurashige
43// remove G4ParticleMomentum 6 Feb. 99 H.Kurashige
44// modified IsApplicable to activate G4Decay for resonances 1 Mar. 00 H.Kurashige
45// Add External Decayer 23 Feb. 2001 H.Kurashige
46// change LowestBinValue,HighestBinValue and TotBin(200) 9 Feb. 2002
47//
48
49#include "G4Decay.hh"
50
52#include "G4SystemOfUnits.hh"
53#include "G4DynamicParticle.hh"
54#include "G4DecayProducts.hh"
55#include "G4DecayTable.hh"
56#include "G4VDecayChannel.hh"
57#include "G4PhysicsLogVector.hh"
59#include "G4VExtDecayer.hh"
60
61// constructor
62G4Decay::G4Decay(const G4String& processName)
63 :G4VRestDiscreteProcess(processName, fDecay),
64 verboseLevel(1),
65 HighestValue(20.0),
66 fRemainderLifeTime(-1.0),
67 pExtDecayer(nullptr)
68{
69 // set Process Sub Type
70 SetProcessSubType(static_cast<int>(DECAY));
71
72#ifdef G4VERBOSE
73 if (GetVerboseLevel()>1) {
74 G4cout << "G4Decay constructor " << " Name:" << processName << G4endl;
75 }
76#endif
77
79}
80
82{
83 if (pExtDecayer != nullptr) {
84 delete pExtDecayer;
85 }
86}
87
89{
90 // check if the particle is stable?
91 if (aParticleType.GetPDGLifeTime() <0.0) {
92 return false;
93 } else if (aParticleType.GetPDGMass() <= 0.0*MeV) {
94 return false;
95 } else {
96 return true;
97 }
98}
99
102{
103 // returns the mean free path in GEANT4 internal units
104 G4double meanlife;
105
106 // get particle
107 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
108 const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
109 G4double aLife = aParticleDef->GetPDGLifeTime();
110
111 // check if the particle is stable?
112 if (aParticleDef->GetPDGStable()) {
113 //1000000 times the life time of the universe
114 meanlife = 1e24 * s;
115
116 } else {
117 meanlife = aLife;
118 }
119
120#ifdef G4VERBOSE
121 if (GetVerboseLevel()>1) {
122 G4cout << "mean life time: "<< meanlife/ns << "[ns]" << G4endl;
123 }
124#endif
125
126 return meanlife;
127}
128
130{
131 // get particle
132 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
133 const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
134 G4double aMass = aParticle->GetMass();
135 G4double aLife = aParticleDef->GetPDGLifeTime();
136
137
138 // returns the mean free path in GEANT4 internal units
139 G4double pathlength;
140 G4double aCtau = c_light * aLife;
141
142 // check if the particle is stable?
143 if (aParticleDef->GetPDGStable()) {
144 pathlength = DBL_MAX;
145
146 //check if the particle has very short life time ?
147 } else if (aCtau < DBL_MIN) {
148 pathlength = DBL_MIN;
149
150 } else {
151 //calculate the mean free path
152 // by using normalized kinetic energy (= Ekin/mass)
153 G4double rKineticEnergy = aParticle->GetKineticEnergy()/aMass;
154 if ( rKineticEnergy > HighestValue) {
155 // gamma >> 1
156 pathlength = ( rKineticEnergy + 1.0)* aCtau;
157 } else if ( rKineticEnergy < DBL_MIN ) {
158 // too slow particle
159#ifdef G4VERBOSE
160 if (GetVerboseLevel()>1) {
161 G4cout << "G4Decay::GetMeanFreePath() !!particle stops!!";
162 G4cout << aParticleDef->GetParticleName() << G4endl;
163 G4cout << "KineticEnergy:" << aParticle->GetKineticEnergy()/GeV <<"[GeV]";
164 }
165#endif
166 pathlength = DBL_MIN;
167 } else {
168 // beta <1
169 pathlength = (aParticle->GetTotalMomentum())/aMass*aCtau ;
170 }
171 }
172 return pathlength;
173}
174
176{
177 return;
178}
179
181{
182 // The DecayIt() method returns by pointer a particle-change object.
183 // Units are expressed in GEANT4 internal units.
184
185 // Initialize ParticleChange
186 // all members of G4VParticleChange are set to equal to
187 // corresponding member in G4Track
189
190 // get particle
191 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
192 const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
193
194 // check if the particle is stable
195 if (aParticleDef->GetPDGStable()) return &fParticleChangeForDecay ;
196
197
198 //check if thePreAssignedDecayProducts exists
199 const G4DecayProducts* o_products = (aParticle->GetPreAssignedDecayProducts());
200 G4bool isPreAssigned = (o_products != nullptr);
201 G4DecayProducts* products = nullptr;
202
203 // decay table
204 G4DecayTable *decaytable = aParticleDef->GetDecayTable();
205
206 // check if external decayer exists
207 G4bool isExtDecayer = (decaytable == nullptr) && (pExtDecayer != nullptr);
208
209 // Error due to NO Decay Table
210 if ( (decaytable == nullptr) && !isExtDecayer && !isPreAssigned ){
211 if (GetVerboseLevel()>0) {
212 G4cout << "G4Decay::DoIt : decay table not defined for ";
213 G4cout << aParticle->GetDefinition()->GetParticleName()<< G4endl;
214 }
216 ed << "For " << aParticle->GetDefinition()->GetParticleName()
217 << " decay probability exist but decay table is not defined "
218 << "- the particle will be killed;\n"
219 << " isExtDecayer: " << isExtDecayer
220 << "; isPreAssigned: " << isPreAssigned;
221 G4Exception( "G4Decay::DecayIt ",
222 "DECAY101",JustWarning, ed);
223
225 // Kill the parent particle
228
231 }
232
233 if (isPreAssigned) {
234 // copy decay products
235 products = new G4DecayProducts(*o_products);
236 } else if ( isExtDecayer ) {
237 // decay according to external decayer
238 products = pExtDecayer->ImportDecayProducts(aTrack);
239 } else {
240 // Decay according to decay table.
241 // Keep trying to choose a candidate decay channel if the dynamic mass
242 // of the decaying particle is below the sum of the PDG masses of the
243 // candidate daughter particles.
244 // This is needed because the decay table used in Geant4 is based on
245 // the assumption of nominal PDG masses, but a wide resonance can have
246 // a dynamic masses well below its nominal PDG masses, and therefore
247 // some of its decay channels can be below the kinematical threshold.
248 // Note that, for simplicity, we ignore here the possibility that
249 // one or more of the candidate daughter particles can be, in turn,
250 // wide resonance. However, if this is the case, and the channel is
251 // accepted, then the masses of the resonance daughter particles will
252 // be sampled by taking into account their widths.
253 G4VDecayChannel* decaychannel = nullptr;
254 G4double massParent = aParticle->GetMass();
255 decaychannel = decaytable->SelectADecayChannel(massParent);
256 if ( decaychannel == nullptr) {
257 // decay channel not found
259 ed << "Can not determine decay channel for "
260 << aParticleDef->GetParticleName() << G4endl
261 << " mass of dynamic particle: "
262 << massParent/GeV << " (GEV)" << G4endl
263 << " dacay table has " << decaytable->entries()
264 << " entries" << G4endl;
265 G4double checkedmass=massParent;
266 if (massParent < 0.) {
267 checkedmass=aParticleDef->GetPDGMass();
268 ed << "Using PDG mass ("<<checkedmass/GeV
269 << "(GeV)) in IsOKWithParentMass" << G4endl;
270 }
271 for (G4int ic =0;ic <decaytable->entries();++ic) {
272 G4VDecayChannel * dc= decaytable->GetDecayChannel(ic);
273 ed << ic << ": BR " << dc->GetBR() << ", IsOK? "
274 << dc->IsOKWithParentMass(checkedmass)
275 << ", --> ";
276 G4int ndaughters=dc->GetNumberOfDaughters();
277 for (G4int id=0;id<ndaughters;++id) {
278 if (id>0) ed << " + "; // seperator, except for first
279 ed << dc->GetDaughterName(id);
280 }
281 ed << G4endl;
282 }
283 G4Exception("G4Decay::DoIt", "DECAY003", FatalException,ed);
284 } else {
285 // execute DecayIt()
286#ifdef G4VERBOSE
287 G4int temp = decaychannel->GetVerboseLevel();
288 if (GetVerboseLevel()>1) {
289 G4cout << "G4Decay::DoIt : selected decay channel addr:"
290 << decaychannel <<G4endl;
291 decaychannel->SetVerboseLevel(GetVerboseLevel());
292 }
293#endif
294 products = decaychannel->DecayIt(aParticle->GetMass());
295#ifdef G4VERBOSE
296 if (GetVerboseLevel()>1) {
297 decaychannel->SetVerboseLevel(temp);
298 }
299#endif
300#ifdef G4VERBOSE
301 if (GetVerboseLevel()>2) {
302 if (! products->IsChecked() ) products->DumpInfo();
303 }
304#endif
305 }
306 }
307
308 // get parent particle information ...................................
309 G4double ParentEnergy = aParticle->GetTotalEnergy();
310 G4double ParentMass = aParticle->GetMass();
311 if (ParentEnergy < ParentMass) {
313 ed << "Total Energy is less than its mass - increased the energy"
314 << "\n Particle: " << aParticle->GetDefinition()->GetParticleName()
315 << "\n Energy:" << ParentEnergy/MeV << "[MeV]"
316 << "\n Mass:" << ParentMass/MeV << "[MeV]";
317 G4Exception( "G4Decay::DecayIt ",
318 "DECAY102",JustWarning, ed);
319 ParentEnergy = ParentMass;
320 }
321
322 G4ThreeVector ParentDirection(aParticle->GetMomentumDirection());
323
324 //boost all decay products to laboratory frame
325 G4double energyDeposit = 0.0;
326 G4double finalGlobalTime = aTrack.GetGlobalTime();
327 G4double finalLocalTime = aTrack.GetLocalTime();
328 if (aTrack.GetTrackStatus() == fStopButAlive ){
329 // AtRest case
330 finalGlobalTime += fRemainderLifeTime;
331 finalLocalTime += fRemainderLifeTime;
332 energyDeposit += aParticle->GetKineticEnergy();
333 if (isPreAssigned) products->Boost( ParentEnergy, ParentDirection);
334 } else {
335 // PostStep case
336 if (!isExtDecayer) products->Boost( ParentEnergy, ParentDirection);
337 }
338 // set polarization for daughter particles
339 DaughterPolarization(aTrack, products);
340
341
342 //add products in fParticleChangeForDecay
343 G4int numberOfSecondaries = products->entries();
345#ifdef G4VERBOSE
346 if (GetVerboseLevel()>1) {
347 G4cout << "G4Decay::DoIt : Decay vertex :";
348 G4cout << " Time: " << finalGlobalTime/ns << "[ns]";
349 G4cout << " X:" << (aTrack.GetPosition()).x() /cm << "[cm]";
350 G4cout << " Y:" << (aTrack.GetPosition()).y() /cm << "[cm]";
351 G4cout << " Z:" << (aTrack.GetPosition()).z() /cm << "[cm]";
352 G4cout << G4endl;
353 G4cout << "G4Decay::DoIt : decay products in Lab. Frame" << G4endl;
354 products->DumpInfo();
355 }
356#endif
357 G4int index;
358 G4ThreeVector currentPosition;
359 const G4TouchableHandle thand = aTrack.GetTouchableHandle();
360 for (index=0; index < numberOfSecondaries; index++){
361 // get current position of the track
362 currentPosition = aTrack.GetPosition();
363 // create a new track object
364 G4Track* secondary = new G4Track( products->PopProducts(),
365 finalGlobalTime ,
366 currentPosition );
367 // switch on good for tracking flag
368 secondary->SetGoodForTrackingFlag();
369 secondary->SetTouchableHandle(thand);
370 // add the secondary track in the List
372 }
373 delete products;
374
375 // Kill the parent particle
379
380 // Clear NumberOfInteractionLengthLeft
382
384}
385
387{
388 // empty implementation
389}
390
391
392
394{
397
398 fRemainderLifeTime = -1.0;
399}
400
402{
403 // Clear NumberOfInteractionLengthLeft
405
407}
408
409
411 const G4Track& track,
412 G4double previousStepSize,
414 )
415{
416 // condition is set to "Not Forced"
418
419 // pre-assigned Decay time
422
423 if (pTime < 0.) {
424 // normal case
425 if ( previousStepSize > 0.0){
426 // subtract NumberOfInteractionLengthLeft
430 }
432 }
433 // get mean free path
434 currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
435
436#ifdef G4VERBOSE
437 if ((currentInteractionLength <=0.0) || (verboseLevel>2)){
438 G4cout << "G4Decay::PostStepGetPhysicalInteractionLength " << G4endl;
439 track.GetDynamicParticle()->DumpInfo();
440 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
441 G4cout << "MeanFreePath = " << currentInteractionLength/cm << "[cm]" <<G4endl;
442 }
443#endif
444
445 G4double value;
448 //fRemainderLifeTime = theNumberOfInteractionLengthLeft*aLife;
449 } else {
450 value = DBL_MAX;
451 }
452
453 return value;
454
455 } else {
456 //pre-assigned Decay time case
457 // reminder proper time
458 fRemainderLifeTime = pTime - track.GetProperTime();
459 if (fRemainderLifeTime <= 0.0) fRemainderLifeTime = 0.0;
460
461 G4double rvalue=0.0;
462 // use pre-assigned Decay time to determine PIL
463 if (aLife>0.0) {
464 // ordinary particle
465 rvalue = (fRemainderLifeTime/aLife)*GetMeanFreePath(track, previousStepSize, condition);
466 } else {
467 // shortlived particle
468 rvalue = c_light * fRemainderLifeTime;
469 // by using normalized kinetic energy (= Ekin/mass)
470 G4double aMass = track.GetDynamicParticle()->GetMass();
471 rvalue *= track.GetDynamicParticle()->GetTotalMomentum()/aMass;
472 }
473 return rvalue;
474 }
475}
476
478 const G4Track& track,
480 )
481{
482 // condition is set to "Not Forced"
484
486 if (pTime >= 0.) {
487 fRemainderLifeTime = pTime - track.GetProperTime();
489 } else {
492 }
493 return fRemainderLifeTime;
494}
495
496
498{
499 pExtDecayer = val;
500
501 // set Process Sub Type
502 if ( pExtDecayer !=0 ) {
503 SetProcessSubType(static_cast<int>(DECAY_External));
504 }
505}
506
508 const G4Track& aTrack,
509 const G4Step& aStep
510 )
511{
512 if ( (aTrack.GetTrackStatus() == fStopButAlive ) ||
513 (aTrack.GetTrackStatus() == fStopAndKill ) ){
516 } else {
517 return DecayIt(aTrack, aStep);
518 }
519}
520
521void G4Decay::ProcessDescription(std::ostream& outFile) const
522{
523 outFile << GetProcessName() << ": Decay of particles. \n"
524 << "kinematics of daughters are dertermined by DecayChannels "
525 << " or by PreAssignedDecayProducts\n";
526}
@ DECAY_External
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
@ NotForced
@ fDecay
@ fStopAndKill
@ fStopButAlive
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void DumpInfo() const
G4int entries() const
G4DynamicParticle * PopProducts()
G4bool IsChecked() const
void Boost(G4double totalEnergy, const G4ThreeVector &momentumDirection)
G4VDecayChannel * GetDecayChannel(G4int index) const
G4VDecayChannel * SelectADecayChannel(G4double parentMass=-1.)
Definition: G4DecayTable.cc:82
G4int entries() const
virtual void ProcessDescription(std::ostream &outFile) const override
Definition: G4Decay.cc:521
virtual ~G4Decay()
Definition: G4Decay.cc:81
G4VExtDecayer * pExtDecayer
Definition: G4Decay.hh:179
virtual void DaughterPolarization(const G4Track &aTrack, G4DecayProducts *products)
Definition: G4Decay.cc:386
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4Decay.cc:180
G4double fRemainderLifeTime
Definition: G4Decay.hh:173
G4ParticleChangeForDecay fParticleChangeForDecay
Definition: G4Decay.hh:176
virtual G4bool IsApplicable(const G4ParticleDefinition &) override
Definition: G4Decay.cc:88
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
Definition: G4Decay.cc:129
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
Definition: G4Decay.cc:410
virtual void EndTracking() override
Definition: G4Decay.cc:401
virtual void StartTracking(G4Track *) override
Definition: G4Decay.cc:393
void SetExtDecayer(G4VExtDecayer *)
Definition: G4Decay.cc:497
const G4double HighestValue
Definition: G4Decay.hh:170
G4Decay(const G4String &processName="Decay")
Definition: G4Decay.cc:62
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition) override
Definition: G4Decay.cc:477
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition) override
Definition: G4Decay.cc:100
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
Definition: G4Decay.cc:175
G4int verboseLevel
Definition: G4Decay.hh:162
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
Definition: G4Decay.cc:507
G4double GetMass() const
void DumpInfo(G4int mode=0) const
const G4ThreeVector & GetMomentumDirection() const
const G4DecayProducts * GetPreAssignedDecayProducts() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4double GetPreAssignedDecayProperTime() const
G4double GetTotalMomentum() const
const G4String & GetName() const
Definition: G4Material.hh:172
void Initialize(const G4Track &) final
G4bool GetPDGStable() const
G4DecayTable * GetDecayTable() const
G4double GetPDGLifeTime() const
const G4String & GetParticleName() const
Definition: G4Step.hh:62
G4TrackStatus GetTrackStatus() const
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4double GetProperTime() const
G4double GetLocalTime() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
void SetGoodForTrackingFlag(G4bool value=true)
G4double GetBR() const
void SetVerboseLevel(G4int value)
G4int GetVerboseLevel() const
G4int GetNumberOfDaughters() const
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0
virtual G4bool IsOKWithParentMass(G4double parentMass)
const G4String & GetDaughterName(G4int anIndex) const
virtual G4DecayProducts * ImportDecayProducts(const G4Track &aTrack)=0
void ProposeTrackStatus(G4TrackStatus status)
void AddSecondary(G4Track *aSecondary)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4double currentInteractionLength
Definition: G4VProcess.hh:339
void SubtractNumberOfInteractionLengthLeft(G4double prevStepSize)
Definition: G4VProcess.hh:528
G4int GetVerboseLevel() const
Definition: G4VProcess.hh:422
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:80
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:428
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:335
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
#define DBL_MIN
Definition: templates.hh:54
#define DBL_MAX
Definition: templates.hh:62
#define ns(x)
Definition: xmltok.c:1649