Geant4-11
G4RadioactiveDecay.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//
27// //
28// File: G4RadioactiveDecay.cc //
29// Author: D.H. Wright (SLAC) //
30// Date: 9 August 2017 //
31// Description: version the G4RadioactiveDecay process by F. Lei and //
32// P.R. Truscott with biasing and activation calculations //
33// removed to a derived class. It performs alpha, beta, //
34// electron capture and isomeric transition decays of //
35// radioactive nuclei. //
36// //
38
39#include "G4RadioactiveDecay.hh"
41
42#include "G4SystemOfUnits.hh"
43#include "G4DynamicParticle.hh"
44#include "G4DecayProducts.hh"
45#include "G4DecayTable.hh"
47#include "G4ITDecay.hh"
48#include "G4BetaDecayType.hh"
49#include "G4BetaMinusDecay.hh"
50#include "G4BetaPlusDecay.hh"
51#include "G4ECDecay.hh"
52#include "G4AlphaDecay.hh"
53#include "G4TritonDecay.hh"
54#include "G4ProtonDecay.hh"
55#include "G4NeutronDecay.hh"
56#include "G4SFDecay.hh"
57#include "G4VDecayChannel.hh"
58#include "G4NuclearDecay.hh"
60#include "G4Fragment.hh"
61#include "G4Ions.hh"
62#include "G4IonTable.hh"
63#include "G4BetaDecayType.hh"
64#include "Randomize.hh"
66#include "G4NuclearLevelData.hh"
68#include "G4LevelManager.hh"
69#include "G4ThreeVector.hh"
70#include "G4Electron.hh"
71#include "G4Positron.hh"
72#include "G4Neutron.hh"
73#include "G4Gamma.hh"
74#include "G4Alpha.hh"
75#include "G4Triton.hh"
76#include "G4Proton.hh"
77
81#include "G4LossTableManager.hh"
86
87#include <vector>
88#include <sstream>
89#include <algorithm>
90#include <fstream>
91
93
94using namespace CLHEP;
95
98
99#ifdef G4MULTITHREADED
100#include "G4AutoLock.hh"
101G4Mutex G4RadioactiveDecay::radioactiveDecayMutex = G4MUTEX_INITIALIZER;
102DecayTableMap* G4RadioactiveDecay::master_dkmap = 0;
103
104G4int& G4RadioactiveDecay::NumberOfInstances()
105{
106 static G4int numberOfInstances = 0;
107 return numberOfInstances;
108}
109#endif
110
112 : G4VRestDiscreteProcess(processName, fDecay), isInitialised(false),
113 forceDecayDirection(0.,0.,0.), forceDecayHalfAngle(0.*deg), dirPath(""),
114 verboseLevel(1),
115 fThresholdForVeryLongDecayTime( 1.0e+27*CLHEP::nanosecond ) // Longer than twice Universe's age
116{
117#ifdef G4VERBOSE
118 if (GetVerboseLevel() > 1) {
119 G4cout << "G4RadioactiveDecay constructor: processName = " << processName
120 << G4endl;
121 }
122#endif
123
125
128
129 // Set up photon evaporation for use in G4ITDecay
133
134 // DHW G4DeexPrecoParameters* deex = G4NuclearLevelData::GetInstance()->GetParameters();
135 // DHW deex->SetCorrelatedGamma(true);
136
137 // Check data directory
138 char* path_var = std::getenv("G4RADIOACTIVEDATA");
139 if (!path_var) {
140 G4Exception("G4RadioactiveDecay()","HAD_RDM_200",FatalException,
141 "Environment variable G4RADIOACTIVEDATA is not set");
142 } else {
143 dirPath = path_var; // convert to string
144 std::ostringstream os;
145 os << dirPath << "/z1.a3"; // used as a dummy
146 std::ifstream testFile;
147 testFile.open(os.str() );
148 if (!testFile.is_open() )
149 G4Exception("G4RadioactiveDecay()","HAD_RDM_201",FatalException,
150 "Environment variable G4RADIOACTIVEDATA is set, but does not point to correct directory");
151 }
152
153 // Reset the list of user defined data files
155
156 // Instantiate the map of decay tables
157#ifdef G4MULTITHREADED
158 G4AutoLock lk(&G4RadioactiveDecay::radioactiveDecayMutex);
159 NumberOfInstances()++;
160 if(!master_dkmap) master_dkmap = new DecayTableMap;
161#endif
162 dkmap = new DecayTableMap;
163
164 // Apply default values
165 applyARM = true;
166 applyICM = true; // Always on; keep only for backward compatibility
167
168 // RDM applies to all logical volumes by default
169 isAllVolumesMode = true;
172}
173
174void G4RadioactiveDecay::ProcessDescription(std::ostream& outFile) const
175{
176 outFile << "The radioactive decay process (G4RadioactiveDecay) handles the\n"
177 << "alpha, beta+, beta-, electron capture and isomeric transition\n"
178 << "decays of nuclei (G4GenericIon) with masses A > 4.\n"
179 << "The required half-lives and decay schemes are retrieved from\n"
180 << "the RadioactiveDecay database which was derived from ENSDF.\n";
181}
182
183
185{
187 delete photonEvaporation;
188 for (DecayTableMap::iterator i = dkmap->begin(); i != dkmap->end(); i++) {
189 delete i->second;
190 }
191 dkmap->clear();
192 delete dkmap;
193#ifdef G4MULTITHREADED
194 G4AutoLock lk(&G4RadioactiveDecay::radioactiveDecayMutex);
195 --NumberOfInstances();
196 if(NumberOfInstances()==0)
197 {
198 for (DecayTableMap::iterator i = master_dkmap->begin(); i != master_dkmap->end(); i++) {
199 delete i->second;
200 }
201 master_dkmap->clear();
202 delete master_dkmap;
203 }
204#endif
205}
206
207
209{
210 // All particles other than G4Ions, are rejected by default
211 if (((const G4Ions*)(&aParticle))->GetExcitationEnergy() > 0.) {return true;}
212 if (aParticle.GetParticleName() == "GenericIon") {
213 return true;
214 } else if (!(aParticle.GetParticleType() == "nucleus")
215 || aParticle.GetPDGLifeTime() < 0. ) {
216 return false;
217 }
218
219 // Determine whether the nuclide falls into the correct A and Z range
220 G4int A = ((const G4Ions*) (&aParticle))->GetAtomicMass();
221 G4int Z = ((const G4Ions*) (&aParticle))->GetAtomicNumber();
222
224 {return false;}
226 {return false;}
227 return true;
228}
229
231{
232 G4String key = aNucleus->GetParticleName();
233 DecayTableMap::iterator table_ptr = dkmap->find(key);
234
235 G4DecayTable* theDecayTable = 0;
236 if (table_ptr == dkmap->end() ) { // If table not there,
237 theDecayTable = LoadDecayTable(*aNucleus); // load from file and
238 if(theDecayTable) (*dkmap)[key] = theDecayTable; // store in library
239 } else {
240 theDecayTable = table_ptr->second;
241 }
242 return theDecayTable;
243}
244
245
247{
249 G4LogicalVolume* volume = nullptr;
250 volume = theLogicalVolumes->GetVolume(aVolume);
251 if (volume != nullptr)
252 {
253 ValidVolumes.push_back(aVolume);
254 std::sort(ValidVolumes.begin(), ValidVolumes.end());
255 // sort need for performing binary_search
256
257 if (GetVerboseLevel() > 0)
258 G4cout << " Radioactive decay applied to " << aVolume << G4endl;
259 }
260 else
261 {
263 ed << aVolume << " is not a valid logical volume name."
264 << " Decay not activated for it."
265 << G4endl;
266 G4Exception("G4RadioactiveDecay::SelectAVolume()", "HAD_RDM_300",
267 JustWarning, ed);
268 }
269}
270
271
273{
275 G4LogicalVolume* volume = nullptr;
276 volume = theLogicalVolumes->GetVolume(aVolume);
277 if (volume != nullptr)
278 {
279 auto location= std::find(ValidVolumes.cbegin(),ValidVolumes.cend(),aVolume);
280 if (location != ValidVolumes.cend() )
281 {
282 ValidVolumes.erase(location);
283 std::sort(ValidVolumes.begin(), ValidVolumes.end());
284 isAllVolumesMode = false;
285 if (GetVerboseLevel() > 0)
286 G4cout << " G4RadioactiveDecay::DeselectAVolume: " << aVolume
287 << " is removed from list " << G4endl;
288 }
289 else
290 {
292 ed << aVolume << " is not in the list. No action taken." << G4endl;
293 G4Exception("G4RadioactiveDecay::DeselectAVolume()", "HAD_RDM_300",
294 JustWarning, ed);
295 }
296 }
297 else
298 {
300 ed << aVolume << " is not a valid logical volume name. No action taken."
301 << G4endl;
302 G4Exception("G4RadioactiveDecay::DeselectAVolume()", "HAD_RDM_300",
303 JustWarning, ed);
304 }
305}
306
307
309{
311 G4LogicalVolume* volume = nullptr;
312 ValidVolumes.clear();
313#ifdef G4VERBOSE
314 if (GetVerboseLevel()>1)
315 G4cout << " RDM Applies to all Volumes" << G4endl;
316#endif
317 for (std::size_t i = 0; i < theLogicalVolumes->size(); ++i){
318 volume = (*theLogicalVolumes)[i];
319 ValidVolumes.push_back(volume->GetName());
320#ifdef G4VERBOSE
321 if (GetVerboseLevel()>1)
322 G4cout << " RDM Applies to Volume " << volume->GetName() << G4endl;
323#endif
324 }
325 std::sort(ValidVolumes.begin(), ValidVolumes.end());
326 // sort needed in order to allow binary_search
327 isAllVolumesMode=true;
328}
329
330
332{
333 ValidVolumes.clear();
334 isAllVolumesMode=false;
335#ifdef G4VERBOSE
336 if (GetVerboseLevel() > 1) G4cout << "RDM removed from all volumes" << G4endl;
337#endif
338}
339
340
342// //
343// GetMeanLifeTime (required by the base class) //
344// //
346
349{
350 G4double meanlife = 0.;
351 const G4DynamicParticle* theParticle = theTrack.GetDynamicParticle();
352 const G4ParticleDefinition* theParticleDef = theParticle->GetDefinition();
353 G4double theLife = theParticleDef->GetPDGLifeTime();
354#ifdef G4VERBOSE
355 if (GetVerboseLevel() > 2) {
356 G4cout << "G4RadioactiveDecay::GetMeanLifeTime() " << G4endl;
357 G4cout << "KineticEnergy: " << theParticle->GetKineticEnergy()/GeV
358 << " GeV, Mass: " << theParticle->GetMass()/GeV
359 << " GeV, Life time: " << theLife/ns << " ns " << G4endl;
360 }
361#endif
362 if (theParticleDef->GetPDGStable()) {meanlife = DBL_MAX;}
363 else if (theLife < 0.0) {meanlife = DBL_MAX;}
364 else {meanlife = theLife;}
365 // Set meanlife to zero for excited istopes which are not in the
366 // RDM database
367 if (((const G4Ions*)(theParticleDef))->GetExcitationEnergy() > 0. &&
368 meanlife == DBL_MAX) {meanlife = 0.;}
369#ifdef G4VERBOSE
370 if (GetVerboseLevel() > 2)
371 G4cout << " mean life time: " << meanlife/s << " s " << G4endl;
372#endif
373
374 return meanlife;
375}
376
378// //
379// GetMeanFreePath for decay in flight //
380// //
382
385{
386 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
387 const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
388 G4double tau = aParticleDef->GetPDGLifeTime();
389 G4double aMass = aParticle->GetMass();
390
391#ifdef G4VERBOSE
392 if (GetVerboseLevel() > 2) {
393 G4cout << "G4RadioactiveDecay::GetMeanFreePath() " << G4endl;
394 G4cout << " KineticEnergy: " << aParticle->GetKineticEnergy()/GeV
395 << " GeV, Mass: " << aMass/GeV << " GeV, tau: " << tau << " ns "
396 << G4endl;
397 }
398#endif
399 G4double pathlength = DBL_MAX;
400 if (tau != -1) {
401 // Ion can decay
402
403 if (tau < -1000.0) {
404 pathlength = DBL_MIN; // nuclide had very short lifetime or wasn't in table
405
406 } else if (tau < 0.0) {
407 G4cout << aParticleDef->GetParticleName() << " has lifetime " << tau << G4endl;
409 ed << "Ion has negative lifetime " << tau
410 << " but is not stable. Setting mean free path to DBL_MAX" << G4endl;
411 G4Exception("G4RadioactiveDecay::GetMeanFreePath()", "HAD_RDM_011",
412 JustWarning, ed);
413 pathlength = DBL_MAX;
414
415 } else {
416 // Calculate mean free path
417 G4double betaGamma = aParticle->GetTotalMomentum()/aMass;
418 pathlength = c_light*tau*betaGamma;
419
420 if (pathlength < DBL_MIN) {
421 pathlength = DBL_MIN;
422#ifdef G4VERBOSE
423 if (GetVerboseLevel() > 2) {
424 G4cout << "G4Decay::GetMeanFreePath: "
425 << aParticleDef->GetParticleName()
426 << " stops, kinetic energy = "
427 << aParticle->GetKineticEnergy()/keV <<" keV " << G4endl;
428 }
429#endif
430 }
431 }
432 }
433
434#ifdef G4VERBOSE
435 if (GetVerboseLevel() > 2) {
436 G4cout << "mean free path: "<< pathlength/m << " m" << G4endl;
437 }
438#endif
439 return pathlength;
440}
441
443// //
444// BuildPhysicsTable - initialization of atomic de-excitation //
445// //
447
449{
450 if (!isInitialised) {
451 isInitialised = true;
452#ifdef G4VERBOSE
455#endif
456 }
459}
460
462// //
463// StreamInfo - stream out parameters //
464// //
466
467void
468G4RadioactiveDecay::StreamInfo(std::ostream& os, const G4String& endline)
469{
473
474 G4int prec = os.precision(5);
475 os << "======================================================================"
476 << endline;
477 os << "====== Radioactive Decay Physics Parameters ======="
478 << endline;
479 os << "======================================================================"
480 << endline;
481 os << "Max life time "
482 << deex->GetMaxLifeTime()/CLHEP::ps << " ps" << endline;
483 os << "Internal e- conversion flag "
484 << deex->GetInternalConversionFlag() << endline;
485 os << "Stored internal conversion coefficients "
486 << deex->StoreICLevelData() << endline;
487 os << "Enable correlated gamma emission "
488 << deex->CorrelatedGamma() << endline;
489 os << "Max 2J for sampling of angular correlations "
490 << deex->GetTwoJMAX() << endline;
491 os << "Atomic de-excitation enabled "
492 << emparam->Fluo() << endline;
493 os << "Auger electron emission enabled "
494 << emparam->Auger() << endline;
495 os << "Check EM cuts disabled for atomic de-excitation "
496 << emparam->DeexcitationIgnoreCut() << endline;
497 os << "Use Bearden atomic level energies "
498 << emparam->BeardenFluoDir() << endline;
499 os << "Use ANSTO fluorescence model "
500 << emparam->ANSTOFluoDir() << endline;
501 os << "Threshold for very long decay time at rest "
502 << fThresholdForVeryLongDecayTime/CLHEP::ns << " ns" << endline;
503 os << "======================================================================"
504 << G4endl;
505 os.precision(prec);
506}
507
508
510// //
511// LoadDecayTable loads the decay scheme from the RadioactiveDecay database //
512// for the parent nucleus. //
513// //
515
518{
519 // Generate input data file name using Z and A of the parent nucleus
520 // file containing radioactive decay data.
521 G4int A = ((const G4Ions*)(&theParentNucleus))->GetAtomicMass();
522 G4int Z = ((const G4Ions*)(&theParentNucleus))->GetAtomicNumber();
523
524 G4double levelEnergy = ((const G4Ions*)(&theParentNucleus))->GetExcitationEnergy();
525 G4Ions::G4FloatLevelBase floatingLevel =
526 ((const G4Ions*)(&theParentNucleus))->GetFloatLevelBase();
527
528#ifdef G4MULTITHREADED
529 G4AutoLock lk(&G4RadioactiveDecay::radioactiveDecayMutex);
530
531 G4String key = theParentNucleus.GetParticleName();
532 DecayTableMap::iterator master_table_ptr = master_dkmap->find(key);
533
534 if (master_table_ptr != master_dkmap->end() ) { // If table is there
535 return master_table_ptr->second;
536 }
537#endif
538
539 //Check if data have been provided by the user
541
542 if (file == "") {
543 std::ostringstream os;
544 os << dirPath << "/z" << Z << ".a" << A << '\0';
545 file = os.str();
546 }
547
548 G4DecayTable* theDecayTable = new G4DecayTable();
549 G4bool found(false); // True if energy level matches one in table
550
551 std::ifstream DecaySchemeFile;
552 DecaySchemeFile.open(file);
553
554 if (DecaySchemeFile.good()) {
555 // Initialize variables used for reading in radioactive decay data
556 G4bool floatMatch(false);
557 const G4int nMode = G4RadioactiveDecayModeSize;
558 G4double modeTotalBR[nMode] = {0.0};
559 G4double modeSumBR[nMode];
560 for (G4int i = 0; i < nMode; i++) {
561 modeSumBR[i] = 0.0;
562 }
563
564 char inputChars[120]={' '};
565 G4String inputLine;
566 G4String recordType("");
567 G4String floatingFlag("");
568 G4String daughterFloatFlag("");
569 G4Ions::G4FloatLevelBase daughterFloatLevel;
570 G4RadioactiveDecayMode theDecayMode;
571 G4double decayModeTotal(0.0);
572 G4double parentExcitation(0.0);
573 G4double a(0.0);
574 G4double b(0.0);
575 G4double c(0.0);
576 G4double dummy(0.0);
577 G4BetaDecayType betaType(allowed);
578
579 // Loop through each data file record until you identify the decay
580 // data relating to the nuclide of concern.
581
582 G4bool complete(false); // bool insures only one set of values read for any
583 // given parent energy level
584 G4int loop = 0;
585 while (!complete && !DecaySchemeFile.getline(inputChars, 120).eof()) { /* Loop checking, 01.09.2015, D.Wright */
586 loop++;
587 if (loop > 100000) {
588 G4Exception("G4RadioactiveDecay::LoadDecayTable()", "HAD_RDM_100",
589 JustWarning, "While loop count exceeded");
590 break;
591 }
592
593 inputLine = inputChars;
594 G4StrUtil::rstrip(inputLine);
595 if (inputChars[0] != '#' && inputLine.length() != 0) {
596 std::istringstream tmpStream(inputLine);
597
598 if (inputChars[0] == 'P') {
599 // Nucleus is a parent type. Check excitation level to see if it
600 // matches that of theParentNucleus
601 tmpStream >> recordType >> parentExcitation >> floatingFlag >> dummy;
602 // "dummy" takes the place of half-life
603 // Now read in from ENSDFSTATE in particle category
604
605 if (found) {
606 complete = true;
607 } else {
608 // Take first level which matches excitation energy regardless of floating level
609 found = (std::abs(parentExcitation*keV - levelEnergy) < levelTolerance);
610 if (floatingLevel != noFloat) {
611 // If floating level specificed, require match of both energy and floating level
612 floatMatch = (floatingLevel == G4Ions::FloatLevelBase(floatingFlag.back()) );
613 if (!floatMatch) found = false;
614 }
615 }
616
617 } else if (found) {
618 // The right part of the radioactive decay data file has been found. Search
619 // through it to determine the mode of decay of the subsequent records.
620
621 // Store for later the total decay probability for each decay mode
622 if (inputLine.length() < 72) {
623 tmpStream >> theDecayMode >> dummy >> decayModeTotal;
624 switch (theDecayMode) {
625 case IT:
626 {
627 G4ITDecay* anITChannel = new G4ITDecay(&theParentNucleus, decayModeTotal,
628 0.0, 0.0, photonEvaporation);
629// anITChannel->SetHLThreshold(halflifethreshold);
630 anITChannel->SetARM(applyARM);
631 theDecayTable->Insert(anITChannel);
632// anITChannel->DumpNuclearInfo();
633 }
634 break;
635 case BetaMinus:
636 modeTotalBR[BetaMinus] = decayModeTotal; break;
637 case BetaPlus:
638 modeTotalBR[BetaPlus] = decayModeTotal; break;
639 case KshellEC:
640 modeTotalBR[KshellEC] = decayModeTotal; break;
641 case LshellEC:
642 modeTotalBR[LshellEC] = decayModeTotal; break;
643 case MshellEC:
644 modeTotalBR[MshellEC] = decayModeTotal; break;
645 case NshellEC:
646 modeTotalBR[NshellEC] = decayModeTotal; break;
647 case Alpha:
648 modeTotalBR[Alpha] = decayModeTotal; break;
649 case Proton:
650 modeTotalBR[Proton] = decayModeTotal; break;
651 case Neutron:
652 modeTotalBR[Neutron] = decayModeTotal; break;
653 case SpFission:
654 modeTotalBR[SpFission] = decayModeTotal; break;
655 case BDProton:
656 /* Not yet implemented */ break;
657 case BDNeutron:
658 /* Not yet implemented */ break;
659 case Beta2Minus:
660 /* Not yet implemented */ break;
661 case Beta2Plus:
662 /* Not yet implemented */ break;
663 case Proton2:
664 /* Not yet implemented */ break;
665 case Neutron2:
666 /* Not yet implemented */ break;
667 case Triton:
668 modeTotalBR[Triton] = decayModeTotal; break;
669 case RDM_ERROR:
670
671 default:
672 G4Exception("G4RadioactiveDecay::LoadDecayTable()", "HAD_RDM_000",
673 FatalException, "Selected decay mode does not exist");
674 } // switch
675
676 } else {
677 if (inputLine.length() < 84) {
678 tmpStream >> theDecayMode >> a >> daughterFloatFlag >> b >> c;
679 betaType = allowed;
680 } else {
681 tmpStream >> theDecayMode >> a >> daughterFloatFlag >> b >> c >> betaType;
682 }
683
684 // Allowed transitions are the default. Forbidden transitions are
685 // indicated in the last column.
686 a /= 1000.;
687 c /= 1000.;
688 b /= 100.;
689 daughterFloatLevel = G4Ions::FloatLevelBase(daughterFloatFlag.back());
690
691 switch (theDecayMode) {
692 case BetaMinus:
693 {
694 G4BetaMinusDecay* aBetaMinusChannel =
695 new G4BetaMinusDecay(&theParentNucleus, b, c*MeV, a*MeV,
696 daughterFloatLevel, betaType);
697// aBetaMinusChannel->DumpNuclearInfo();
698// aBetaMinusChannel->SetHLThreshold(halflifethreshold);
699 theDecayTable->Insert(aBetaMinusChannel);
700 modeSumBR[BetaMinus] += b;
701 }
702 break;
703
704 case BetaPlus:
705 {
706 G4BetaPlusDecay* aBetaPlusChannel =
707 new G4BetaPlusDecay(&theParentNucleus, b, c*MeV, a*MeV,
708 daughterFloatLevel, betaType);
709// aBetaPlusChannel->DumpNuclearInfo();
710// aBetaPlusChannel->SetHLThreshold(halflifethreshold);
711 theDecayTable->Insert(aBetaPlusChannel);
712 modeSumBR[BetaPlus] += b;
713 }
714 break;
715
716 case KshellEC: // K-shell electron capture
717 {
718 G4ECDecay* aKECChannel =
719 new G4ECDecay(&theParentNucleus, b, c*MeV, a*MeV,
720 daughterFloatLevel, KshellEC);
721// aKECChannel->DumpNuclearInfo();
722// aKECChannel->SetHLThreshold(halflifethreshold);
723 aKECChannel->SetARM(applyARM);
724 theDecayTable->Insert(aKECChannel);
725 modeSumBR[KshellEC] += b;
726 }
727 break;
728
729 case LshellEC: // L-shell electron capture
730 {
731 G4ECDecay* aLECChannel =
732 new G4ECDecay(&theParentNucleus, b, c*MeV, a*MeV,
733 daughterFloatLevel, LshellEC);
734// aLECChannel->DumpNuclearInfo();
735// aLECChannel->SetHLThreshold(halflifethreshold);
736 aLECChannel->SetARM(applyARM);
737 theDecayTable->Insert(aLECChannel);
738 modeSumBR[LshellEC] += b;
739 }
740 break;
741
742 case MshellEC: // M-shell electron capture
743 {
744 G4ECDecay* aMECChannel =
745 new G4ECDecay(&theParentNucleus, b, c*MeV, a*MeV,
746 daughterFloatLevel, MshellEC);
747// aMECChannel->DumpNuclearInfo();
748// aMECChannel->SetHLThreshold(halflifethreshold);
749 aMECChannel->SetARM(applyARM);
750 theDecayTable->Insert(aMECChannel);
751 modeSumBR[MshellEC] += b;
752 }
753 break;
754
755 case NshellEC: // N-shell electron capture
756 {
757 G4ECDecay* aNECChannel =
758 new G4ECDecay(&theParentNucleus, b, c*MeV, a*MeV,
759 daughterFloatLevel, NshellEC);
760// aNECChannel->DumpNuclearInfo();
761// aNECChannel->SetHLThreshold(halflifethreshold);
762 aNECChannel->SetARM(applyARM);
763 theDecayTable->Insert(aNECChannel);
764 modeSumBR[NshellEC] += b;
765 }
766 break;
767
768 case Alpha:
769 {
770 G4AlphaDecay* anAlphaChannel =
771 new G4AlphaDecay(&theParentNucleus, b, c*MeV, a*MeV,
772 daughterFloatLevel);
773// anAlphaChannel->DumpNuclearInfo();
774// anAlphaChannel->SetHLThreshold(halflifethreshold);
775 theDecayTable->Insert(anAlphaChannel);
776 modeSumBR[Alpha] += b;
777 }
778 break;
779
780 case Proton:
781 {
782 G4ProtonDecay* aProtonChannel =
783 new G4ProtonDecay(&theParentNucleus, b, c*MeV, a*MeV,
784 daughterFloatLevel);
785// aProtonChannel->DumpNuclearInfo();
786// aProtonChannel->SetHLThreshold(halflifethreshold);
787 theDecayTable->Insert(aProtonChannel);
788 modeSumBR[Proton] += b;
789 }
790 break;
791
792 case Neutron:
793 {
794 G4NeutronDecay* aNeutronChannel =
795 new G4NeutronDecay(&theParentNucleus, b, c*MeV, a*MeV,
796 daughterFloatLevel);
797// aNeutronChannel->DumpNuclearInfo();
798// aNeutronChannel->SetHLThreshold(halflifethreshold);
799 theDecayTable->Insert(aNeutronChannel);
800 modeSumBR[Neutron] += b;
801 }
802 break;
803
804 case SpFission:
805 {
806 G4SFDecay* aSpontFissChannel =
807// new G4SFDecay(&theParentNucleus, decayModeTotal, 0.0, 0.0);
808 new G4SFDecay(&theParentNucleus, b, c*MeV, a*MeV,
809 daughterFloatLevel);
810 theDecayTable->Insert(aSpontFissChannel);
811 modeSumBR[SpFission] += b;
812 }
813 break;
814
815 case BDProton:
816 // Not yet implemented
817 // G4cout << " beta-delayed proton decay, a = " << a << ", b = " << b << ", c = " << c << G4endl;
818 break;
819
820 case BDNeutron:
821 // Not yet implemented
822 // G4cout << " beta-delayed neutron decay, a = " << a << ", b = " << b << ", c = " << c << G4endl;
823 break;
824
825 case Beta2Minus:
826 // Not yet implemented
827 // G4cout << " Double beta- decay, a = " << a << ", b = " << b << ", c = " << c << G4endl;
828 break;
829
830 case Beta2Plus:
831 // Not yet implemented
832 // G4cout << " Double beta+ decay, a = " << a << ", b = " << b << ", c = " << c << G4endl;
833 break;
834
835 case Proton2:
836 // Not yet implemented
837 // G4cout << " Double proton decay, a = " << a << ", b = " << b << ", c = " << c << G4endl;
838 break;
839
840 case Neutron2:
841 // Not yet implemented
842 // G4cout << " Double beta- decay, a = " << a << ", b = " << b << ", c = " << c << G4endl;
843 break;
844
845 case Triton:
846 {
847 G4TritonDecay* aTritonChannel =
848 new G4TritonDecay(&theParentNucleus, b, c*MeV, a*MeV,
849 daughterFloatLevel);
850 // anAlphaChannel->DumpNuclearInfo();
851 // anAlphaChannel->SetHLThreshold(halflifethreshold);
852 theDecayTable->Insert(aTritonChannel);
853 modeSumBR[Triton] += b;
854 }
855 break;
856
857 case RDM_ERROR:
858
859 default:
860 G4Exception("G4RadioactiveDecay::LoadDecayTable()", "HAD_RDM_000",
861 FatalException, "Selected decay mode does not exist");
862 } // switch
863 } // line < 72
864 } // if char == P
865 } // if char != #
866 } // While
867
868 // Go through the decay table and make sure that the branching ratios are
869 // correctly normalised.
870
871 G4VDecayChannel* theChannel = 0;
872 G4NuclearDecay* theNuclearDecayChannel = 0;
873 G4String mode = "";
874
875 G4double theBR = 0.0;
876 for (G4int i = 0; i < theDecayTable->entries(); i++) {
877 theChannel = theDecayTable->GetDecayChannel(i);
878 theNuclearDecayChannel = static_cast<G4NuclearDecay*>(theChannel);
879 theDecayMode = theNuclearDecayChannel->GetDecayMode();
880
881 if (theDecayMode != IT) {
882 theBR = theChannel->GetBR();
883 theChannel->SetBR(theBR*modeTotalBR[theDecayMode]/modeSumBR[theDecayMode]);
884 }
885 }
886 } // decay file exists
887
888 DecaySchemeFile.close();
889
890 if (!found && levelEnergy > 0) {
891 // Case where IT cascade for excited isotopes has no entries in RDM database
892 // Decay mode is isomeric transition.
893 G4ITDecay* anITChannel = new G4ITDecay(&theParentNucleus, 1.0, 0.0, 0.0,
895// anITChannel->SetHLThreshold(halflifethreshold);
896 anITChannel->SetARM(applyARM);
897 theDecayTable->Insert(anITChannel);
898 }
899
900 if (theDecayTable && GetVerboseLevel() > 1) {
901 theDecayTable->DumpInfo();
902 }
903
904#ifdef G4MULTITHREADED
905 //(*master_dkmap)[key] = theDecayTable; // store in master library
906#endif
907 return theDecayTable;
908}
909
910void
912{
913 if (Z < 1 || A < 2) G4cout << "Z and A not valid!" << G4endl;
914
915 std::ifstream DecaySchemeFile(filename);
916 if (DecaySchemeFile) {
917 G4int ID_ion = A*1000 + Z;
918 theUserRadioactiveDataFiles[ID_ion] = filename;
919 } else {
921 ed << filename << " does not exist! " << G4endl;
922 G4Exception("G4RadioactiveDecay::AddUserDecayDataFile()", "HAD_RDM_001",
923 FatalException, ed);
924 }
925}
926
928// //
929// DecayIt //
930// //
932
935{
936 // Initialize G4ParticleChange object, get particle details and decay table
939 const G4DynamicParticle* theParticle = theTrack.GetDynamicParticle();
940 const G4ParticleDefinition* theParticleDef = theParticle->GetDefinition();
941
942 // First check whether RDM applies to the current logical volume
943 if (!isAllVolumesMode) {
944 if (!std::binary_search(ValidVolumes.begin(), ValidVolumes.end(),
945 theTrack.GetVolume()->GetLogicalVolume()->GetName())) {
946#ifdef G4VERBOSE
947 if (GetVerboseLevel()>1) {
948 G4cout <<"G4RadioactiveDecay::DecayIt : "
949 << theTrack.GetVolume()->GetLogicalVolume()->GetName()
950 << " is not selected for the RDM"<< G4endl;
951 G4cout << " There are " << ValidVolumes.size() << " volumes" << G4endl;
952 G4cout << " The Valid volumes are " << G4endl;
953 for (size_t i = 0; i< ValidVolumes.size(); i++)
954 G4cout << ValidVolumes[i] << G4endl;
955 }
956#endif
958
959 // Kill the parent particle.
964 }
965 }
966
967 // Now check if particle is valid for RDM
968 if (!(IsApplicable(*theParticleDef) ) ) {
969 // Particle is not an ion or is outside the nucleuslimits for decay
970#ifdef G4VERBOSE
971 if (GetVerboseLevel() > 1) {
972 G4cout << "G4RadioactiveDecay::DecayIt : "
973 << theParticleDef->GetParticleName()
974 << " is not an ion or is outside (Z,A) limits set for the decay. "
975 << " Set particle change accordingly. "
976 << G4endl;
977 }
978#endif
980
981 // Kill the parent particle
986 }
987
988 G4DecayTable* theDecayTable = GetDecayTable(theParticleDef);
989
990 if (theDecayTable == 0 || theDecayTable->entries() == 0) {
991 // No data in the decay table. Set particle change parameters
992 // to indicate this.
993#ifdef G4VERBOSE
994 if (GetVerboseLevel() > 1) {
995 G4cout << "G4RadioactiveDecay::DecayIt : "
996 << "decay table not defined for "
997 << theParticleDef->GetParticleName()
998 << ". Set particle change accordingly. "
999 << G4endl;
1000 }
1001#endif
1003
1004 // Kill the parent particle.
1009
1010 } else {
1011 // Data found. Try to decay nucleus
1012
1013/*
1014 G4double energyDeposit = 0.0;
1015 G4double finalGlobalTime = theTrack.GetGlobalTime();
1016 G4double finalLocalTime = theTrack.GetLocalTime();
1017 G4int index;
1018 G4ThreeVector currentPosition;
1019 currentPosition = theTrack.GetPosition();
1020
1021 G4DecayProducts* products = DoDecay(*theParticleDef);
1022
1023 // If the product is the same as the input kill the track if
1024 // necessary to prevent infinite loop (11/05/10, F.Lei)
1025 if (products->entries() == 1) {
1026 fParticleChangeForRadDecay.SetNumberOfSecondaries(0);
1027 fParticleChangeForRadDecay.ProposeTrackStatus(fStopAndKill);
1028 fParticleChangeForRadDecay.ProposeLocalEnergyDeposit(0.0);
1029 ClearNumberOfInteractionLengthLeft();
1030 return &fParticleChangeForRadDecay;
1031 }
1032
1033 // Get parent particle information and boost the decay products to the
1034 // laboratory frame based on this information.
1035
1036 // The Parent Energy used for the boost should be the total energy of
1037 // the nucleus of the parent ion without the energy of the shell electrons
1038 // (correction for bug 1359 by L. Desorgher)
1039 G4double ParentEnergy = theParticle->GetKineticEnergy()
1040 + theParticle->GetParticleDefinition()->GetPDGMass();
1041 G4ThreeVector ParentDirection(theParticle->GetMomentumDirection());
1042
1043 if (theTrack.GetTrackStatus() == fStopButAlive) {
1044 // This condition seems to be always True, further investigation is needed
1045 // (L.Desorgher)
1046 // The particle is decayed at rest.
1047 // since the time is still for rest particle in G4 we need to add the
1048 // additional time lapsed between the particle come to rest and the
1049 // actual decay. This time is simply sampled with the mean-life of
1050 // the particle. But we need to protect the case PDGTime < 0.
1051 // (F.Lei 11/05/10)
1052 G4double temptime = -std::log( G4UniformRand())
1053 *theParticleDef->GetPDGLifeTime();
1054 if (temptime < 0.) temptime = 0.;
1055 finalGlobalTime += temptime;
1056 finalLocalTime += temptime;
1057 energyDeposit += theParticle->GetKineticEnergy();
1058 }
1059 products->Boost(ParentEnergy, ParentDirection);
1060
1061 // Add products in theParticleChangeForRadDecay.
1062 G4int numberOfSecondaries = products->entries();
1063 fParticleChangeForRadDecay.SetNumberOfSecondaries(numberOfSecondaries);
1064
1065#ifdef G4VERBOSE
1066 if (GetVerboseLevel()>1) {
1067 G4cout <<"G4RadioactiveDecay::DecayIt : Decay vertex :";
1068 G4cout <<" Time: " <<finalGlobalTime/ns <<"[ns]";
1069 G4cout <<" X:" <<(theTrack.GetPosition()).x() /cm <<"[cm]";
1070 G4cout <<" Y:" <<(theTrack.GetPosition()).y() /cm <<"[cm]";
1071 G4cout <<" Z:" <<(theTrack.GetPosition()).z() /cm <<"[cm]";
1072 G4cout << G4endl;
1073 G4cout <<"G4Decay::DecayIt : decay products in Lab. Frame" <<G4endl;
1074 products->DumpInfo();
1075 products->IsChecked();
1076 }
1077#endif
1078 for (index=0; index < numberOfSecondaries; index++) {
1079 G4Track* secondary = new G4Track(products->PopProducts(),
1080 finalGlobalTime, currentPosition);
1081 secondary->SetGoodForTrackingFlag();
1082 secondary->SetTouchableHandle(theTrack.GetTouchableHandle());
1083 fParticleChangeForRadDecay.AddSecondary(secondary);
1084 }
1085 delete products;
1086
1087 // Kill the parent particle
1088 fParticleChangeForRadDecay.ProposeTrackStatus(fStopAndKill) ;
1089 fParticleChangeForRadDecay.ProposeLocalEnergyDeposit(energyDeposit);
1090 fParticleChangeForRadDecay.ProposeLocalTime(finalLocalTime);
1091 // Reset NumberOfInteractionLengthLeft.
1092 ClearNumberOfInteractionLengthLeft();
1093*/
1094 // Decay without variance reduction
1095 DecayAnalog(theTrack);
1097 }
1098}
1099
1100
1102{
1103 const G4DynamicParticle* theParticle = theTrack.GetDynamicParticle();
1104 const G4ParticleDefinition* theParticleDef = theParticle->GetDefinition();
1105 G4DecayProducts* products = DoDecay(*theParticleDef);
1106
1107 // Check if the product is the same as input and kill the track if
1108 // necessary to prevent infinite loop (11/05/10, F.Lei)
1109 if (products->entries() == 1) {
1114 delete products;
1115 return;
1116 }
1117
1118 G4double energyDeposit = 0.0;
1119 G4double finalGlobalTime = theTrack.GetGlobalTime();
1120 G4double finalLocalTime = theTrack.GetLocalTime();
1121
1122 // Get parent particle information and boost the decay products to the
1123 // laboratory frame
1124
1125 // ParentEnergy used for the boost should be the total energy of the nucleus
1126 // of the parent ion without the energy of the shell electrons
1127 // (correction for bug 1359 by L. Desorgher)
1128 G4double ParentEnergy = theParticle->GetKineticEnergy()
1129 + theParticle->GetParticleDefinition()->GetPDGMass();
1130 G4ThreeVector ParentDirection(theParticle->GetMomentumDirection());
1131
1132 if (theTrack.GetTrackStatus() == fStopButAlive) {
1133 // this condition seems to be always True, further investigation is needed (L.Desorgher)
1134
1135 // The particle is decayed at rest
1136 // Since the time is for the particle at rest, need to add additional time
1137 // lapsed between particle coming to rest and the actual decay. This time
1138 // is sampled with the mean-life of the particle. Need to protect the case
1139 // PDGTime < 0. (F.Lei 11/05/10)
1140 G4double temptime = -std::log(G4UniformRand() ) *
1141 theParticleDef->GetPDGLifeTime();
1142 if (temptime < 0.) temptime = 0.;
1143 finalGlobalTime += temptime;
1144 finalLocalTime += temptime;
1145 energyDeposit += theParticle->GetKineticEnergy();
1146
1147 // Kill the parent particle, and ignore its decay, if it decays later than the
1148 // threshold fThresholdForVeryLongDecayTime (whose default value corresponds
1149 // to more than twice the age of the universe).
1150 // This kind of cut has been introduced (in April 2021) in order to avoid to
1151 // account energy depositions happening after many billions of years in
1152 // ordinary materials used in calorimetry, in particular Tungsten and Lead
1153 // (via their natural unstable, but very long lived, isotopes, such as
1154 // W183, W180 and Pb204).
1155 // Note that the cut is not on the average, mean lifetime, but on the actual
1156 // sampled global decay time.
1157 if ( finalGlobalTime > fThresholdForVeryLongDecayTime ) {
1162 delete products;
1163 return;
1164 }
1165 }
1166 products->Boost(ParentEnergy, ParentDirection);
1167
1168 // Add products in theParticleChangeForRadDecay.
1169 G4int numberOfSecondaries = products->entries();
1171
1172 if (GetVerboseLevel() > 1) {
1173 G4cout << "G4RadioactiveDecay::DecayAnalog: Decay vertex :";
1174 G4cout << " Time: " << finalGlobalTime/ns << "[ns]";
1175 G4cout << " X:" << (theTrack.GetPosition()).x() /cm << "[cm]";
1176 G4cout << " Y:" << (theTrack.GetPosition()).y() /cm << "[cm]";
1177 G4cout << " Z:" << (theTrack.GetPosition()).z() /cm << "[cm]";
1178 G4cout << G4endl;
1179 G4cout << "G4Decay::DecayIt : decay products in Lab. Frame" << G4endl;
1180 products->DumpInfo();
1181 products->IsChecked();
1182 }
1183
1184 const G4int modelID_forIT = G4PhysicsModelCatalog::GetModelID( "model_RDM_IT" );
1185 G4int modelID = modelID_forIT + 10*theRadDecayMode;
1186 const G4int modelID_forAtomicRelaxation =
1187 G4PhysicsModelCatalog::GetModelID( "model_RDM_AtomicRelaxation" );
1188 for ( G4int index = 0; index < numberOfSecondaries; ++index ) {
1189 G4Track* secondary = new G4Track( products->PopProducts(), finalGlobalTime,
1190 theTrack.GetPosition() );
1191 secondary->SetWeight( theTrack.GetWeight() );
1192 secondary->SetCreatorModelID( modelID );
1193 // Change for atomics relaxation
1194 if ( theRadDecayMode == IT && index > 0 ) {
1195 if ( index == numberOfSecondaries-1 ) {
1196 secondary->SetCreatorModelID( modelID_forIT );
1197 } else {
1198 secondary->SetCreatorModelID( modelID_forAtomicRelaxation) ;
1199 }
1200 } else if ( theRadDecayMode >= KshellEC && theRadDecayMode <= NshellEC &&
1201 index < numberOfSecondaries-1 ) {
1202 secondary->SetCreatorModelID( modelID_forAtomicRelaxation );
1203 }
1204 secondary->SetGoodForTrackingFlag();
1205 secondary->SetTouchableHandle( theTrack.GetTouchableHandle() );
1207 }
1208
1209 delete products;
1210
1211 // Kill the parent particle
1215
1216 // Reset NumberOfInteractionLengthLeft.
1218}
1219
1220
1223{
1224 G4DecayProducts* products = 0;
1225 G4DecayTable* theDecayTable = GetDecayTable(&theParticleDef);
1226
1227 // Choose a decay channel.
1228 // G4DecayTable::SelectADecayChannel checks to see if sum of daughter masses
1229 // exceeds parent mass. Pass it the parent mass + maximum Q value to account
1230 // for difference in mass defect.
1231 G4double parentPlusQ = theParticleDef.GetPDGMass() + 30.*MeV;
1232 G4VDecayChannel* theDecayChannel = theDecayTable->SelectADecayChannel(parentPlusQ);
1233
1234 if (theDecayChannel == 0) {
1235 // Decay channel not found.
1237 ed << " Cannot determine decay channel for " << theParticleDef.GetParticleName() << G4endl;
1238 G4Exception("G4RadioactiveDecay::DoDecay", "HAD_RDM_013",
1239 FatalException, ed);
1240 } else {
1241 // A decay channel has been identified, so execute the DecayIt.
1242#ifdef G4VERBOSE
1243 if (GetVerboseLevel() > 1) {
1244 G4cout << "G4RadioactiveDecay::DoIt : selected decay channel addr: "
1245 << theDecayChannel << G4endl;
1246 }
1247#endif
1248 theRadDecayMode = (static_cast<G4NuclearDecay*>(theDecayChannel))->GetDecayMode();
1249 products = theDecayChannel->DecayIt(theParticleDef.GetPDGMass() );
1250
1251 // Apply directional bias if requested by user
1252 CollimateDecay(products);
1253 }
1254
1255 return products;
1256}
1257
1258
1259// Apply directional bias for "visible" daughters (e+-, gamma, n, p, alpha)
1260
1262
1263 if (origin == forceDecayDirection) return; // No collimation requested
1264 if (180.*deg == forceDecayHalfAngle) return;
1265 if (0 == products || 0 == products->entries()) return;
1266
1267#ifdef G4VERBOSE
1268 if (GetVerboseLevel() > 1) G4cout << "Begin of CollimateDecay..." << G4endl;
1269#endif
1270
1271 // Particles suitable for directional biasing (for if-blocks below)
1275 static const G4ParticleDefinition* gamma = G4Gamma::Definition();
1279
1280 G4ThreeVector newDirection; // Re-use to avoid memory churn
1281 for (G4int i=0; i<products->entries(); i++) {
1282 G4DynamicParticle* daughter = (*products)[i];
1283 const G4ParticleDefinition* daughterType =
1284 daughter->GetParticleDefinition();
1285 if (daughterType == electron || daughterType == positron ||
1286 daughterType == neutron || daughterType == gamma ||
1287 daughterType == alpha || daughterType == triton || daughterType == proton) CollimateDecayProduct(daughter);
1288 }
1289}
1290
1292#ifdef G4VERBOSE
1293 if (GetVerboseLevel() > 1) {
1294 G4cout << "CollimateDecayProduct for daughter "
1295 << daughter->GetParticleDefinition()->GetParticleName() << G4endl;
1296 }
1297#endif
1298
1300 if (origin != collimate) daughter->SetMomentumDirection(collimate);
1301}
1302
1303
1304// Choose random direction within collimation cone
1305
1307 if (origin == forceDecayDirection) return origin; // Don't do collimation
1308 if (forceDecayHalfAngle == 180.*deg) return origin;
1309
1311
1312 // Return direction offset by random throw
1313 if (forceDecayHalfAngle > 0.) {
1314 // Generate uniform direction around central axis
1315 G4double phi = 2.*pi*G4UniformRand();
1316 G4double cosMin = std::cos(forceDecayHalfAngle);
1317 G4double cosTheta = (1.-cosMin)*G4UniformRand() + cosMin; // [cosMin,1.)
1318
1319 dir.setPhi(dir.phi()+phi);
1320 dir.setTheta(dir.theta()+std::acos(cosTheta));
1321 }
1322
1323#ifdef G4VERBOSE
1324 if (GetVerboseLevel()>1)
1325 G4cout << " ChooseCollimationDirection returns " << dir << G4endl;
1326#endif
1327
1328 return dir;
1329}
1330
G4BetaDecayType
@ allowed
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
@ fRadioactiveDecay
#define noFloat
Definition: G4Ions.hh:112
static const G4double alpha
@ fDecay
std::map< G4String, G4DecayTable * > DecayTableMap
G4RadioactiveDecayMode
@ G4RadioactiveDecayModeSize
static constexpr double m
Definition: G4SIunits.hh:109
static constexpr double s
Definition: G4SIunits.hh:154
static constexpr double keV
Definition: G4SIunits.hh:202
static constexpr double eV
Definition: G4SIunits.hh:201
static constexpr double nanosecond
Definition: G4SIunits.hh:138
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
static constexpr double pi
Definition: G4SIunits.hh:55
static constexpr double cm
Definition: G4SIunits.hh:99
static constexpr double deg
Definition: G4SIunits.hh:132
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
std::mutex G4Mutex
Definition: G4Threading.hh:81
@ fStopAndKill
@ fStopButAlive
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
double phi() const
double theta() const
void setTheta(double)
void setPhi(double)
static G4Alpha * Definition()
Definition: G4Alpha.cc:48
void DumpInfo() const
G4int entries() const
G4DynamicParticle * PopProducts()
G4bool IsChecked() const
void Boost(G4double totalEnergy, const G4ThreeVector &momentumDirection)
G4VDecayChannel * GetDecayChannel(G4int index) const
void Insert(G4VDecayChannel *aChannel)
Definition: G4DecayTable.cc:53
G4VDecayChannel * SelectADecayChannel(G4double parentMass=-1.)
Definition: G4DecayTable.cc:82
void DumpInfo() const
G4int entries() const
G4bool GetInternalConversionFlag() const
G4double GetMass() const
void SetMomentumDirection(const G4ThreeVector &aDirection)
const G4ThreeVector & GetMomentumDirection() const
const G4ParticleDefinition * GetParticleDefinition() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetTotalMomentum() const
void SetARM(G4bool onoff)
Definition: G4ECDecay.hh:56
static G4Electron * Definition()
Definition: G4Electron.cc:48
static G4EmParameters * Instance()
G4bool BeardenFluoDir() const
G4bool ANSTOFluoDir() const
G4bool Fluo() const
G4bool DeexcitationIgnoreCut() const
G4bool Auger() const
static G4Gamma * Definition()
Definition: G4Gamma.cc:48
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
static G4HadronicParameters * Instance()
void RegisterExtraProcess(G4VProcess *)
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)
static G4HadronicProcessStore * Instance()
void SetARM(G4bool onoff)
Definition: G4ITDecay.hh:59
Definition: G4Ions.hh:52
static G4Ions::G4FloatLevelBase FloatLevelBase(char flbChar)
Definition: G4Ions.cc:103
G4FloatLevelBase
Definition: G4Ions.hh:83
G4LogicalVolume * GetVolume(const G4String &name, G4bool verbose=true, G4bool reverseSearch=false) const
static G4LogicalVolumeStore * GetInstance()
const G4String & GetName() const
static G4Neutron * Definition()
Definition: G4Neutron.cc:53
G4RadioactiveDecayMode GetDecayMode()
G4DeexPrecoParameters * GetParameters()
static G4NuclearLevelData * GetInstance()
G4int GetZMax() const
G4int GetZMin() const
G4int GetAMin() const
G4int GetAMax() const
virtual void Initialize(const G4Track &)
void AddSecondary(G4Track *aSecondary)
G4bool GetPDGStable() const
const G4String & GetParticleType() const
G4double GetPDGLifeTime() const
const G4String & GetParticleName() const
virtual void SetICM(G4bool)
virtual void RDMForced(G4bool)
static G4int GetModelID(const G4int modelIndex)
static G4Positron * Definition()
Definition: G4Positron.cc:48
static G4Proton * Definition()
Definition: G4Proton.cc:48
static const G4ThreeVector origin
std::map< G4int, G4String > theUserRadioactiveDataFiles
void SelectAVolume(const G4String &aVolume)
void DecayAnalog(const G4Track &theTrack)
G4RadioactiveDecayMessenger * theRadioactiveDecayMessenger
void StreamInfo(std::ostream &os, const G4String &endline)
void BuildPhysicsTable(const G4ParticleDefinition &)
std::vector< G4String > ValidVolumes
G4DecayTable * LoadDecayTable(const G4ParticleDefinition &theParentNucleus)
G4RadioactiveDecay(const G4String &processName="RadioactiveDecay")
G4double GetMeanLifeTime(const G4Track &theTrack, G4ForceCondition *condition)
G4DecayProducts * DoDecay(const G4ParticleDefinition &theParticleDef)
void CollimateDecayProduct(G4DynamicParticle *product)
G4NucleusLimits theNucleusLimits
G4DecayTable * GetDecayTable(const G4ParticleDefinition *)
static const G4double levelTolerance
G4ParticleChangeForRadDecay fParticleChangeForRadDecay
G4ThreeVector forceDecayDirection
G4int GetVerboseLevel() const
G4double GetMeanFreePath(const G4Track &theTrack, G4double previousStepSize, G4ForceCondition *condition)
G4ThreeVector ChooseCollimationDirection() const
G4VParticleChange * DecayIt(const G4Track &theTrack, const G4Step &theStep)
G4bool IsApplicable(const G4ParticleDefinition &)
G4RadioactiveDecayMode theRadDecayMode
virtual void ProcessDescription(std::ostream &outFile) const
void AddUserDecayDataFile(G4int Z, G4int A, G4String filename)
G4double fThresholdForVeryLongDecayTime
void DeselectAVolume(const G4String &aVolume)
void CollimateDecay(G4DecayProducts *products)
G4PhotonEvaporation * photonEvaporation
Definition: G4Step.hh:62
G4TrackStatus GetTrackStatus() const
G4VPhysicalVolume * GetVolume() const
G4double GetWeight() const
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4double GetLocalTime() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
void SetCreatorModelID(const G4int id)
void SetGoodForTrackingFlag(G4bool value=true)
static G4Triton * Definition()
Definition: G4Triton.cc:49
G4double GetBR() const
void SetBR(G4double value)
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0
void ProposeTrackStatus(G4TrackStatus status)
void ProposeWeight(G4double finalWeight)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4LogicalVolume * GetLogicalVolume() const
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:424
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
Definition: DoubConv.h:17
static const double prec
Definition: RanecuEngine.cc:61
static constexpr double ns
static constexpr double ps
void rstrip(G4String &str, char c=' ')
Remove trailing characters from string.
G4bool IsMasterThread()
Definition: G4Threading.cc:124
float c_light
Definition: hepunit.py:256
#define DBL_MIN
Definition: templates.hh:54
#define DBL_MAX
Definition: templates.hh:62
#define ns
Definition: xmlparse.cc:614