Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4Decay.cc
이 파일의 문서화 페이지로 가기
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 // $Id: G4Decay.cc 105727 2017-08-16 12:47:05Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------
31 // GEANT 4 class implementation file
32 //
33 // History: first implementation, based on object model of
34 // 2nd December 1995, G.Cosmo
35 // 7 July 1996 H.Kurashige
36 // ------------------------------------------------------------
37 // remove BuildPhysicsTable() 28 Nov. 1997 H.Kurashige
38 // change DBL_EPSIRON to DBL_MIN 14 Dec. 1997 H.Kurashige
39 // modified for new ParticleChange 12 Mar. 1998 H.Kurashige
40 // modified for "GoodForTrackingFlag" 19 June 1998 H.Kurashige
41 // rename thePhysicsTable to aPhyscisTable 2 Aug. 1998 H.Kurashige
42 // modified IsApplicable in order to protect the decay from registered
43 // to resonances 12 Dec. 1998 H.Kurashige
44 // remove G4ParticleMomentum 6 Feb. 99 H.Kurashige
45 // modified IsApplicable to activate G4Decay for resonances 1 Mar. 00 H.Kurashige
46 // Add External Decayer 23 Feb. 2001 H.Kurashige
47 // change LowestBinValue,HighestBinValue and TotBin(200) 9 Feb. 2002
48 //
49 
50 #include "G4Decay.hh"
51 
52 #include "G4PhysicalConstants.hh"
53 #include "G4SystemOfUnits.hh"
54 #include "G4DynamicParticle.hh"
55 #include "G4DecayProducts.hh"
56 #include "G4DecayTable.hh"
57 #include "G4VDecayChannel.hh"
58 #include "G4PhysicsLogVector.hh"
60 #include "G4VExtDecayer.hh"
61 
62 // constructor
63 G4Decay::G4Decay(const G4String& processName)
64  :G4VRestDiscreteProcess(processName, fDecay),
65  verboseLevel(1),
66  HighestValue(20.0),
67  fRemainderLifeTime(-1.0),
68  pExtDecayer(nullptr)
69 {
70  // set Process Sub Type
71  SetProcessSubType(static_cast<int>(DECAY));
72 
73 #ifdef G4VERBOSE
74  if (GetVerboseLevel()>1) {
75  G4cout << "G4Decay constructor " << " Name:" << processName << G4endl;
76  }
77 #endif
78 
80 }
81 
83 {
84  if (pExtDecayer != nullptr) {
85  delete pExtDecayer;
86  }
87 }
88 
90 {
91  // check if the particle is stable?
92  if (aParticleType.GetPDGLifeTime() <0.0) {
93  return false;
94  } else if (aParticleType.GetPDGMass() <= 0.0*MeV) {
95  return false;
96  } else {
97  return true;
98  }
99 }
100 
103 {
104  // returns the mean free path in GEANT4 internal units
105  G4double meanlife;
106 
107  // get particle
108  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
109  const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
110  G4double aLife = aParticleDef->GetPDGLifeTime();
111 
112  // check if the particle is stable?
113  if (aParticleDef->GetPDGStable()) {
114  //1000000 times the life time of the universe
115  meanlife = 1e24 * s;
116 
117  } else {
118  meanlife = aLife;
119  }
120 
121 #ifdef G4VERBOSE
122  if (GetVerboseLevel()>1) {
123  G4cout << "mean life time: "<< meanlife/ns << "[ns]" << G4endl;
124  }
125 #endif
126 
127  return meanlife;
128 }
129 
131 {
132  // get particle
133  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
134  const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
135  G4double aMass = aParticle->GetMass();
136  G4double aLife = aParticleDef->GetPDGLifeTime();
137 
138 
139  // returns the mean free path in GEANT4 internal units
140  G4double pathlength;
141  G4double aCtau = c_light * aLife;
142 
143  // check if the particle is stable?
144  if (aParticleDef->GetPDGStable()) {
145  pathlength = DBL_MAX;
146 
147  //check if the particle has very short life time ?
148  } else if (aCtau < DBL_MIN) {
149  pathlength = DBL_MIN;
150 
151  } else {
152  //calculate the mean free path
153  // by using normalized kinetic energy (= Ekin/mass)
154  G4double rKineticEnergy = aParticle->GetKineticEnergy()/aMass;
155  if ( rKineticEnergy > HighestValue) {
156  // gamma >> 1
157  pathlength = ( rKineticEnergy + 1.0)* aCtau;
158  } else if ( rKineticEnergy < DBL_MIN ) {
159  // too slow particle
160 #ifdef G4VERBOSE
161  if (GetVerboseLevel()>1) {
162  G4cout << "G4Decay::GetMeanFreePath() !!particle stops!!";
163  G4cout << aParticleDef->GetParticleName() << G4endl;
164  G4cout << "KineticEnergy:" << aParticle->GetKineticEnergy()/GeV <<"[GeV]";
165  }
166 #endif
167  pathlength = DBL_MIN;
168  } else {
169  // beta <1
170  pathlength = (aParticle->GetTotalMomentum())/aMass*aCtau ;
171  }
172  }
173  return pathlength;
174 }
175 
177 {
178  return;
179 }
180 
182 {
183  // The DecayIt() method returns by pointer a particle-change object.
184  // Units are expressed in GEANT4 internal units.
185 
186  // Initialize ParticleChange
187  // all members of G4VParticleChange are set to equal to
188  // corresponding member in G4Track
190 
191  // get particle
192  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
193  const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
194 
195  // check if the particle is stable
196  if (aParticleDef->GetPDGStable()) return &fParticleChangeForDecay ;
197 
198 
199  //check if thePreAssignedDecayProducts exists
200  const G4DecayProducts* o_products = (aParticle->GetPreAssignedDecayProducts());
201  G4bool isPreAssigned = (o_products != nullptr);
202  G4DecayProducts* products = nullptr;
203 
204  // decay table
205  G4DecayTable *decaytable = aParticleDef->GetDecayTable();
206 
207  // check if external decayer exists
208  G4bool isExtDecayer = (decaytable == nullptr) && (pExtDecayer != nullptr);
209 
210  // Error due to NO Decay Table
211  if ( (decaytable == nullptr) && !isExtDecayer && !isPreAssigned ){
212  if (GetVerboseLevel()>0) {
213  G4cout << "G4Decay::DoIt : decay table not defined for ";
214  G4cout << aParticle->GetDefinition()->GetParticleName()<< G4endl;
215  }
216  G4Exception( "G4Decay::DecayIt ",
217  "DECAY101",JustWarning,
218  "Decay table is not defined");
219 
221  // Kill the parent particle
224 
226  return &fParticleChangeForDecay ;
227  }
228 
229  if (isPreAssigned) {
230  // copy decay products
231  products = new G4DecayProducts(*o_products);
232  } else if ( isExtDecayer ) {
233  // decay according to external decayer
234  products = pExtDecayer->ImportDecayProducts(aTrack);
235  } else {
236  // Decay according to decay table.
237  // Keep trying to choose a candidate decay channel if the dynamic mass
238  // of the decaying particle is below the sum of the PDG masses of the
239  // candidate daughter particles.
240  // This is needed because the decay table used in Geant4 is based on
241  // the assumption of nominal PDG masses, but a wide resonance can have
242  // a dynamic masses well below its nominal PDG masses, and therefore
243  // some of its decay channels can be below the kinematical threshold.
244  // Note that, for simplicity, we ignore here the possibility that
245  // one or more of the candidate daughter particles can be, in turn,
246  // wide resonance. However, if this is the case, and the channel is
247  // accepted, then the masses of the resonance daughter particles will
248  // be sampled by taking into account their widths.
249  G4VDecayChannel* decaychannel = nullptr;
250  G4double massParent = aParticle->GetMass();
251  decaychannel = decaytable->SelectADecayChannel(massParent);
252  if ( decaychannel == nullptr) {
253  // decay channel not found
255  ed << "Can not determine decay channel for "
256  << aParticleDef->GetParticleName() << G4endl
257  << " mass of dynamic particle: "
258  << massParent/GeV << " (GEV)" << G4endl
259  << " dacay table has " << decaytable->entries()
260  << " entries" << G4endl;
261  G4double checkedmass=massParent;
262  if (massParent < 0.) {
263  checkedmass=aParticleDef->GetPDGMass();
264  ed << "Using PDG mass ("<<checkedmass/GeV
265  << "(GeV)) in IsOKWithParentMass" << G4endl;
266  }
267  for (G4int ic =0;ic <decaytable->entries();++ic) {
268  G4VDecayChannel * dc= decaytable->GetDecayChannel(ic);
269  ed << ic << ": BR " << dc->GetBR() << ", IsOK? "
270  << dc->IsOKWithParentMass(checkedmass)
271  << ", --> ";
272  G4int ndaughters=dc->GetNumberOfDaughters();
273  for (G4int id=0;id<ndaughters;++id) {
274  if (id>0) ed << " + "; // seperator, except for first
275  ed << dc->GetDaughterName(id);
276  }
277  ed << G4endl;
278  }
279  G4Exception("G4Decay::DoIt", "DECAY003", FatalException,ed);
280  } else {
281  // execute DecayIt()
282 #ifdef G4VERBOSE
283  G4int temp = decaychannel->GetVerboseLevel();
284  if (GetVerboseLevel()>1) {
285  G4cout << "G4Decay::DoIt : selected decay channel addr:"
286  << decaychannel <<G4endl;
287  decaychannel->SetVerboseLevel(GetVerboseLevel());
288  }
289 #endif
290  products = decaychannel->DecayIt(aParticle->GetMass());
291 #ifdef G4VERBOSE
292  if (GetVerboseLevel()>1) {
293  decaychannel->SetVerboseLevel(temp);
294  }
295 #endif
296 #ifdef G4VERBOSE
297  if (GetVerboseLevel()>2) {
298  if (! products->IsChecked() ) products->DumpInfo();
299  }
300 #endif
301  }
302  }
303 
304  // get parent particle information ...................................
305  G4double ParentEnergy = aParticle->GetTotalEnergy();
306  G4double ParentMass = aParticle->GetMass();
307  if (ParentEnergy < ParentMass) {
308  if (GetVerboseLevel()>0) {
309  G4cout << "G4Decay::DoIt : Total Energy is less than its mass" << G4endl;
310  G4cout << " Particle: " << aParticle->GetDefinition()->GetParticleName();
311  G4cout << " Energy:" << ParentEnergy/MeV << "[MeV]";
312  G4cout << " Mass:" << ParentMass/MeV << "[MeV]";
313  G4cout << G4endl;
314  }
315  G4Exception( "G4Decay::DecayIt ",
316  "DECAY102",JustWarning,
317  "Total Energy is less than its mass");
318  ParentEnergy = ParentMass;
319  }
320 
321  G4ThreeVector ParentDirection(aParticle->GetMomentumDirection());
322 
323  //boost all decay products to laboratory frame
324  G4double energyDeposit = 0.0;
325  G4double finalGlobalTime = aTrack.GetGlobalTime();
326  G4double finalLocalTime = aTrack.GetLocalTime();
327  if (aTrack.GetTrackStatus() == fStopButAlive ){
328  // AtRest case
329  finalGlobalTime += fRemainderLifeTime;
330  finalLocalTime += fRemainderLifeTime;
331  energyDeposit += aParticle->GetKineticEnergy();
332  if (isPreAssigned) products->Boost( ParentEnergy, ParentDirection);
333  } else {
334  // PostStep case
335  if (!isExtDecayer) products->Boost( ParentEnergy, ParentDirection);
336  }
337  // set polarization for daughter particles
338  DaughterPolarization(aTrack, products);
339 
340 
341  //add products in fParticleChangeForDecay
342  G4int numberOfSecondaries = products->entries();
343  fParticleChangeForDecay.SetNumberOfSecondaries(numberOfSecondaries);
344 #ifdef G4VERBOSE
345  if (GetVerboseLevel()>1) {
346  G4cout << "G4Decay::DoIt : Decay vertex :";
347  G4cout << " Time: " << finalGlobalTime/ns << "[ns]";
348  G4cout << " X:" << (aTrack.GetPosition()).x() /cm << "[cm]";
349  G4cout << " Y:" << (aTrack.GetPosition()).y() /cm << "[cm]";
350  G4cout << " Z:" << (aTrack.GetPosition()).z() /cm << "[cm]";
351  G4cout << G4endl;
352  G4cout << "G4Decay::DoIt : decay products in Lab. Frame" << G4endl;
353  products->DumpInfo();
354  }
355 #endif
356  G4int index;
357  G4ThreeVector currentPosition;
358  const G4TouchableHandle thand = aTrack.GetTouchableHandle();
359  for (index=0; index < numberOfSecondaries; index++){
360  // get current position of the track
361  currentPosition = aTrack.GetPosition();
362  // create a new track object
363  G4Track* secondary = new G4Track( products->PopProducts(),
364  finalGlobalTime ,
365  currentPosition );
366  // switch on good for tracking flag
367  secondary->SetGoodForTrackingFlag();
368  secondary->SetTouchableHandle(thand);
369  // add the secondary track in the List
371  }
372  delete products;
373 
374  // Kill the parent particle
377  fParticleChangeForDecay.ProposeLocalTime( finalLocalTime );
378 
379  // Clear NumberOfInteractionLengthLeft
381 
382  return &fParticleChangeForDecay ;
383 }
384 
386 {
387  // empty implementation
388 }
389 
390 
391 
393 {
396 
397  fRemainderLifeTime = -1.0;
398 }
399 
401 {
402  // Clear NumberOfInteractionLengthLeft
404 
406 }
407 
408 
410  const G4Track& track,
411  G4double previousStepSize,
413  )
414 {
415  // condition is set to "Not Forced"
416  *condition = NotForced;
417 
418  // pre-assigned Decay time
421 
422  if (pTime < 0.) {
423  // normal case
424  if ( previousStepSize > 0.0){
425  // subtract NumberOfInteractionLengthLeft
426  SubtractNumberOfInteractionLengthLeft(previousStepSize);
429  }
431  }
432  // get mean free path
433  currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
434 
435 #ifdef G4VERBOSE
436  if ((currentInteractionLength <=0.0) || (verboseLevel>2)){
437  G4cout << "G4Decay::PostStepGetPhysicalInteractionLength " << G4endl;
438  track.GetDynamicParticle()->DumpInfo();
439  G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
440  G4cout << "MeanFreePath = " << currentInteractionLength/cm << "[cm]" <<G4endl;
441  }
442 #endif
443 
444  G4double value;
447  //fRemainderLifeTime = theNumberOfInteractionLengthLeft*aLife;
448  } else {
449  value = DBL_MAX;
450  }
451 
452  return value;
453 
454  } else {
455  //pre-assigned Decay time case
456  // reminder proper time
457  fRemainderLifeTime = pTime - track.GetProperTime();
459 
460  G4double rvalue=0.0;
461  // use pre-assigned Decay time to determine PIL
462  if (aLife>0.0) {
463  // ordinary particle
464  rvalue = (fRemainderLifeTime/aLife)*GetMeanFreePath(track, previousStepSize, condition);
465  } else {
466  // shortlived particle
467  rvalue = c_light * fRemainderLifeTime;
468  // by using normalized kinetic energy (= Ekin/mass)
469  G4double aMass = track.GetDynamicParticle()->GetMass();
470  rvalue *= track.GetDynamicParticle()->GetTotalMomentum()/aMass;
471  }
472  return rvalue;
473  }
474 }
475 
477  const G4Track& track,
479  )
480 {
481  // condition is set to "Not Forced"
482  *condition = NotForced;
483 
485  if (pTime >= 0.) {
486  fRemainderLifeTime = pTime - track.GetProperTime();
488  } else {
491  }
492  return fRemainderLifeTime;
493 }
494 
495 
497 {
498  pExtDecayer = val;
499 
500  // set Process Sub Type
501  if ( pExtDecayer !=0 ) {
502  SetProcessSubType(static_cast<int>(DECAY_External));
503  }
504 }
505 
507  const G4Track& aTrack,
508  const G4Step& aStep
509  )
510 {
511  if ( (aTrack.GetTrackStatus() == fStopButAlive ) ||
512  (aTrack.GetTrackStatus() == fStopAndKill ) ){
514  return &fParticleChangeForDecay;
515  } else {
516  return DecayIt(aTrack, aStep);
517  }
518 }
519 
520 void G4Decay::ProcessDescription(std::ostream& outFile) const
521 {
522  outFile << GetProcessName() << ": Decay of particles. \n"
523  << "kinematics of daughters are dertermined by DecayChannels "
524  << " or by PreAssignedDecayProducts\n";
525 }
Float_t x
Definition: compare.C:6
G4double GetLocalTime() const
void AddSecondary(G4Track *aSecondary)
G4double GetBR() const
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:45
G4double fRemainderLifeTime
Definition: G4Decay.hh:174
G4DecayTable * GetDecayTable() const
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:96
G4int verboseLevel
Definition: G4Decay.hh:163
static constexpr double MeV
Definition: G4SIunits.hh:214
G4double GetPDGLifeTime() const
void SetExtDecayer(G4VExtDecayer *)
Definition: G4Decay.cc:496
void SetVerboseLevel(G4int value)
G4int GetNumberOfDaughters() const
#define G4endl
Definition: G4ios.hh:61
Float_t y
Definition: compare.C:6
const G4ThreeVector & GetMomentumDirection() const
Double_t z
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
Definition: G4Decay.cc:506
const G4String & GetParticleName() const
void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
Definition: G4VProcess.hh:547
const G4TouchableHandle & GetTouchableHandle() const
virtual ~G4Decay()
Definition: G4Decay.cc:82
static constexpr double perMillion
Definition: G4SIunits.hh:334
virtual G4bool IsApplicable(const G4ParticleDefinition &) override
Definition: G4Decay.cc:89
G4double condition(const G4ErrorSymMatrix &m)
G4bool IsChecked() const
G4DynamicParticle * PopProducts()
G4double GetPDGMass() const
void Boost(G4double totalEnergy, const G4ThreeVector &momentumDirection)
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition) override
Definition: G4Decay.cc:101
G4double GetProperTime() const
G4double currentInteractionLength
Definition: G4VProcess.hh:297
G4int entries() const
const G4String & GetName() const
Definition: G4Material.hh:179
const G4double HighestValue
Definition: G4Decay.hh:171
const XML_Char * s
Definition: expat.h:262
virtual void ProcessDescription(std::ostream &outFile) const override
Definition: G4Decay.cc:520
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4ParticleDefinition * GetDefinition() const
G4ParticleChangeForDecay fParticleChangeForDecay
Definition: G4Decay.hh:177
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
const XML_Char int const XML_Char * value
Definition: expat.h:331
G4VDecayChannel * GetDecayChannel(G4int index) const
void DumpInfo() const
G4double GetGlobalTime() const
G4double GetPreAssignedDecayProperTime() const
G4int GetVerboseLevel() const
Definition: G4VProcess.hh:445
void SetTouchableHandle(const G4TouchableHandle &apValue)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:411
const G4ThreeVector & GetPosition() const
Definition: G4Step.hh:76
G4Material * GetMaterial() const
G4int GetVerboseLevel() const
const G4String & GetDaughterName(G4int anIndex) const
G4bool GetPDGStable() const
virtual void EndTracking() override
Definition: G4Decay.cc:400
void SetGoodForTrackingFlag(G4bool value=true)
const G4DecayProducts * GetPreAssignedDecayProducts() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition) override
Definition: G4Decay.cc:476
int G4int
Definition: G4Types.hh:78
static constexpr double c_light
#define DBL_MIN
Definition: templates.hh:75
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:435
G4ForceCondition
G4double GetKineticEnergy() const
static constexpr double cm
Definition: G4SIunits.hh:119
G4GLOB_DLL std::ostream G4cout
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:450
virtual void StartTracking(G4Track *) override
Definition: G4Decay.cc:392
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0
virtual G4DecayProducts * ImportDecayProducts(const G4Track &aTrack)=0
G4VExtDecayer * pExtDecayer
Definition: G4Decay.hh:180
virtual G4bool IsOKWithParentMass(G4double parentMass)
G4VDecayChannel * SelectADecayChannel(G4double parentMass=-1.)
Definition: G4DecayTable.cc:85
G4double GetMass() const
G4double GetTotalMomentum() const
G4double GetTotalEnergy() const
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
Definition: G4Decay.cc:176
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4TrackStatus GetTrackStatus() const
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
Definition: G4Decay.cc:130
G4Decay(const G4String &processName="Decay")
Definition: G4Decay.cc:63
virtual void Initialize(const G4Track &)
#define DBL_MAX
Definition: templates.hh:83
void DumpInfo(G4int mode=0) const
void ProposeTrackStatus(G4TrackStatus status)
static constexpr double GeV
Definition: G4SIunits.hh:217
#define ns
Definition: xmlparse.cc:614
const G4DynamicParticle * GetDynamicParticle() const
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
Definition: G4Decay.cc:409
virtual void DaughterPolarization(const G4Track &aTrack, G4DecayProducts *products)
Definition: G4Decay.cc:385
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4Decay.cc:181
void SetNumberOfSecondaries(G4int totSecondaries)
G4int entries() const