Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ITDecay Class Reference

#include <G4ITDecay.hh>

+ Inheritance diagram for G4ITDecay:

Public Member Functions

 G4ITDecay (G4PhotonEvaporation *aPhotonEvap)
 
 G4ITDecay (const G4ParticleDefinition *theParentNucleus, const G4double &theBR, const G4double &Qvalue, const G4double &excitation)
 
 ~G4ITDecay () override=default
 
void SetupDecay (const G4ParticleDefinition *)
 
G4DecayProductsDecayIt (G4double) override
 
void DumpNuclearInfo () override
 
void SetARM (G4bool onoff)
 
- Public Member Functions inherited from G4NuclearDecay
 G4NuclearDecay (const G4String &channelName, const G4RadioactiveDecayMode &mode, const G4double &excitation, const G4Ions::G4FloatLevelBase &floatingLevel)
 
 ~G4NuclearDecay () override=default
 
G4bool IsOKWithParentMass (G4double parentMass) override
 
G4RadioactiveDecayMode GetDecayMode () const
 
G4double GetDaughterExcitation () const
 
G4Ions::G4FloatLevelBase GetFloatingLevel () const
 
G4ParticleDefinitionGetDaughterNucleus ()
 
- Public Member Functions inherited from G4VDecayChannel
 G4VDecayChannel (const G4String &aName, G4int Verbose=1)
 
 G4VDecayChannel (const G4String &aName, const G4String &theParentName, G4double theBR, G4int theNumberOfDaughters, const G4String &theDaughterName1, const G4String &theDaughterName2="", const G4String &theDaughterName3="", const G4String &theDaughterName4="", const G4String &theDaughterName5="")
 
virtual ~G4VDecayChannel ()
 
G4bool operator== (const G4VDecayChannel &r) const
 
G4bool operator!= (const G4VDecayChannel &r) const
 
G4bool operator< (const G4VDecayChannel &right) const
 
const G4StringGetKinematicsName () const
 
G4double GetBR () const
 
G4int GetNumberOfDaughters () const
 
G4ParticleDefinitionGetParent ()
 
G4ParticleDefinitionGetDaughter (G4int anIndex)
 
G4int GetAngularMomentum ()
 
const G4StringGetParentName () const
 
const G4StringGetDaughterName (G4int anIndex) const
 
G4double GetParentMass () const
 
G4double GetDaughterMass (G4int anIndex) const
 
void SetParent (const G4ParticleDefinition *particle_type)
 
void SetParent (const G4String &particle_name)
 
void SetBR (G4double value)
 
void SetNumberOfDaughters (G4int value)
 
void SetDaughter (G4int anIndex, const G4ParticleDefinition *particle_type)
 
void SetDaughter (G4int anIndex, const G4String &particle_name)
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
void DumpInfo ()
 
G4double GetRangeMass () const
 
void SetRangeMass (G4double val)
 
void SetPolarization (const G4ThreeVector &)
 
const G4ThreeVectorGetPolarization () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VDecayChannel
 G4VDecayChannel ()
 
 G4VDecayChannel (const G4VDecayChannel &)
 
G4VDecayChanneloperator= (const G4VDecayChannel &)
 
void ClearDaughtersName ()
 
void CheckAndFillDaughters ()
 
void CheckAndFillParent ()
 
G4double DynamicalMass (G4double massPDG, G4double width, G4double maxDev=1.0) const
 
- Protected Attributes inherited from G4NuclearDecay
const G4RadioactiveDecayMode theMode
 
- Protected Attributes inherited from G4VDecayChannel
G4String kinematics_name = ""
 
G4double rbranch = 0.0
 
G4Stringparent_name = nullptr
 
G4String ** daughters_name = nullptr
 
G4double rangeMass = 2.5
 
G4ThreeVector parent_polarization
 
G4ParticleTableparticletable = nullptr
 
G4ParticleDefinitionG4MT_parent = nullptr
 
G4ParticleDefinition ** G4MT_daughters = nullptr
 
G4double G4MT_parent_mass = 0.0
 
G4doubleG4MT_daughters_mass = nullptr
 
G4doubleG4MT_daughters_width = nullptr
 
G4Mutex daughtersMutex
 
G4Mutex parentMutex
 
G4int numberOfDaughters = 0
 
G4int verboseLevel = 1
 
- Static Protected Attributes inherited from G4VDecayChannel
static const G4String noName = " "
 

Detailed Description

Definition at line 45 of file G4ITDecay.hh.

Constructor & Destructor Documentation

◆ G4ITDecay() [1/2]

G4ITDecay::G4ITDecay ( G4PhotonEvaporation * aPhotonEvap)

Definition at line 50 of file G4ITDecay.cc.

51 : G4NuclearDecay("IT decay", IT, 0.0, noFloat), photonEvaporation(ptr)
52{}
#define noFloat
Definition G4Ions.hh:119
G4NuclearDecay(const G4String &channelName, const G4RadioactiveDecayMode &mode, const G4double &excitation, const G4Ions::G4FloatLevelBase &floatingLevel)

◆ G4ITDecay() [2/2]

G4ITDecay::G4ITDecay ( const G4ParticleDefinition * theParentNucleus,
const G4double & theBR,
const G4double & Qvalue,
const G4double & excitation )

Definition at line 54 of file G4ITDecay.cc.

57 : G4NuclearDecay("IT decay", IT, excitationE, noFloat)
58{
59 SetParent(theParentNucleus); // Store name of parent nucleus, delete G4MT_parent
60 SetBR(branch);
62 SetDaughter(0, theParentNucleus);
63
64 SetupDecay(theParentNucleus);
65}
void SetupDecay(const G4ParticleDefinition *)
Definition G4ITDecay.cc:67
void SetBR(G4double value)
void SetNumberOfDaughters(G4int value)
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
void SetParent(const G4ParticleDefinition *particle_type)

◆ ~G4ITDecay()

G4ITDecay::~G4ITDecay ( )
overridedefault

Member Function Documentation

◆ DecayIt()

G4DecayProducts * G4ITDecay::DecayIt ( G4double )
overridevirtual

Implements G4VDecayChannel.

Definition at line 74 of file G4ITDecay.cc.

75{
76 // Set up final state
77 // parentParticle is set at rest here because boost with correct momentum
78 // is done later
79 G4LorentzVector atRest(theParent->GetPDGMass(), G4ThreeVector(0.,0.,0.) );
80 G4DynamicParticle parentParticle(theParent, atRest);
81 G4DecayProducts* products = new G4DecayProducts(parentParticle);
82
83 // Let G4PhotonEvaporation do the decay
84 G4Fragment parentNucleus(parentA, parentZ, atRest);
85
86 // one emission, parent nucleaus become less excited
87 G4Fragment* eOrGamma = photonEvaporation->EmittedFragment(&parentNucleus);
88
89 // Modified nuclide is returned as dynDaughter
90 auto theIonTable = G4ParticleTable::GetParticleTable()->GetIonTable();
91 G4ParticleDefinition* daughterIon =
92 theIonTable->GetIon(parentZ, parentA, parentNucleus.GetExcitationEnergy(),
93 G4Ions::FloatLevelBase(parentNucleus.GetFloatingLevelNumber()));
94 G4DynamicParticle* dynDaughter = new G4DynamicParticle(daughterIon,
95 parentNucleus.GetMomentum());
96
97 if (nullptr != eOrGamma) {
98 G4DynamicParticle* eOrGammaDyn =
100 eOrGamma->GetMomentum() );
101 eOrGammaDyn->SetProperTime(eOrGamma->GetCreationTime() );
102 products->PushProducts(eOrGammaDyn);
103 delete eOrGamma;
104
105 // Now do atomic relaxation if e- is emitted
106 if (applyARM) {
107 G4int shellIndex = photonEvaporation->GetVacantShellNumber();
108 if (shellIndex > -1) {
109 G4VAtomDeexcitation* atomDeex =
111 if (atomDeex->IsFluoActive() && parentZ > 5 && parentZ < 105) {
112 G4int nShells = G4AtomicShells::GetNumberOfShells(parentZ);
113 if (shellIndex >= nShells) shellIndex = nShells;
115 const G4AtomicShell* shell = atomDeex->GetAtomicShell(parentZ, as);
116 std::vector<G4DynamicParticle*> armProducts;
117
118 // VI, SI
119 // Allows fixing of Bugzilla 1727
120 G4double deexLimit = 0.1*keV;
121 if (G4EmParameters::Instance()->DeexcitationIgnoreCut()) deexLimit =0.;
122 //
123
124 atomDeex->GenerateParticles(&armProducts, shell, parentZ, deexLimit,
125 deexLimit);
126 G4double productEnergy = 0.;
127 for (G4int i = 0; i < G4int(armProducts.size()); i++)
128 productEnergy += armProducts[i]->GetKineticEnergy();
129
130 G4double deficit = shell->BindingEnergy() - productEnergy;
131 if (deficit > 0.0) {
132 // Add a dummy electron to make up extra energy
133 G4double cosTh = 1.-2.*G4UniformRand();
134 G4double sinTh = std::sqrt(1.- cosTh*cosTh);
135 G4double phi = twopi*G4UniformRand();
136
137 G4ThreeVector electronDirection(sinTh*std::sin(phi),
138 sinTh*std::cos(phi), cosTh);
139 G4DynamicParticle* extra =
140 new G4DynamicParticle(G4Electron::Electron(), electronDirection,
141 deficit);
142 armProducts.push_back(extra);
143 }
144
145 std::size_t nArm = armProducts.size();
146 if (nArm > 0) {
147 G4ThreeVector bst = dynDaughter->Get4Momentum().boostVector();
148 for (std::size_t i = 0; i < nArm; ++i) {
149 G4DynamicParticle* dp = armProducts[i];
150 G4LorentzVector lv = dp->Get4Momentum().boost(bst);
151 dp->Set4Momentum(lv);
152 products->PushProducts(dp);
153 }
154 }
155 }
156 }
157 } // if ARM on
158 } // eOrGamma
159
160 products->PushProducts(dynDaughter);
161
162 // Energy conservation check
163 /*
164 G4int newSize = products->entries();
165 G4DynamicParticle* temp = 0;
166 G4double KEsum = 0.0;
167 for (G4int i = 0; i < newSize; i++) {
168 temp = products->operator[](i);
169 KEsum += temp->GetKineticEnergy();
170 }
171 G4double eCons = G4MT_parent->GetPDGMass() - dynDaughter->GetMass() - KEsum;
172 G4cout << " IT check: Ediff (keV) = " << eCons/keV << G4endl;
173 */
174 return products;
175}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
G4double BindingEnergy() const
static G4int GetNumberOfShells(G4int Z)
G4int PushProducts(G4DynamicParticle *aParticle)
void SetProperTime(G4double)
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4EmParameters * Instance()
const G4LorentzVector & GetMomentum() const
G4double GetCreationTime() const
const G4ParticleDefinition * GetParticleDefinition() const
static G4Ions::G4FloatLevelBase FloatLevelBase(char flbChar)
Definition G4Ions.cc:107
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
G4int GetVacantShellNumber() const
G4Fragment * EmittedFragment(G4Fragment *theNucleus) override
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)

Referenced by G4Radioactivation::AddDeexcitationSpectrumForBiasMode(), and G4RadioactiveDecay::DoDecay().

◆ DumpNuclearInfo()

void G4ITDecay::DumpNuclearInfo ( )
overridevirtual

Implements G4NuclearDecay.

Definition at line 178 of file G4ITDecay.cc.

179{
180 if (theParent != nullptr) {
181 G4cout << " G4ITDecay for parent nucleus " << theParent->GetParticleName() << G4endl;
182 }
183}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
const G4String & GetParticleName() const

◆ SetARM()

void G4ITDecay::SetARM ( G4bool onoff)
inline

Definition at line 65 of file G4ITDecay.hh.

65{applyARM = onoff;}

Referenced by G4RadioactiveDecay::BuildPhysicsTable().

◆ SetupDecay()

void G4ITDecay::SetupDecay ( const G4ParticleDefinition * theParentNucleus)

Definition at line 67 of file G4ITDecay.cc.

68{
69 theParent = theParentNucleus;
70 parentZ = theParentNucleus->GetAtomicNumber();
71 parentA = theParentNucleus->GetAtomicMass();
72}
G4int GetAtomicNumber() const
G4int GetAtomicMass() const

Referenced by G4Radioactivation::AddDeexcitationSpectrumForBiasMode(), G4RadioactiveDecay::DoDecay(), and G4ITDecay().


The documentation for this class was generated from the following files: