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