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