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

#include <G4ECDecay.hh>

+ Inheritance diagram for G4ECDecay:

Public Member Functions

 G4ECDecay (const G4ParticleDefinition *theParentNucleus, const G4double &theBR, const G4double &Qvalue, const G4double &excitation, const G4Ions::G4FloatLevelBase &flb, const G4RadioactiveDecayMode &mode)
 
virtual ~G4ECDecay ()
 
virtual G4DecayProductsDecayIt (G4double)
 
virtual void DumpNuclearInfo ()
 
void SetARM (G4bool onoff)
 
- Public Member Functions inherited from G4NuclearDecay
 G4NuclearDecay (const G4String &channelName, const G4RadioactiveDecayMode &mode, const G4double &excitation, const G4Ions::G4FloatLevelBase &floatingLevel)
 
virtual ~G4NuclearDecay ()
 
G4RadioactiveDecayMode GetDecayMode ()
 
G4double GetDaughterExcitation ()
 
G4Ions::G4FloatLevelBase GetFloatingLevel ()
 
G4ParticleDefinitionGetDaughterNucleus ()
 
void SetHLThreshold (G4double HLT)
 
G4double GetHLThreshold ()
 
virtual void DumpNuclearInfo ()=0
 
- 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
 
virtual G4DecayProductsDecayIt (G4double parentMass=-1.0)=0
 
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)
 
virtual G4bool IsOKWithParentMass (G4double parentMass)
 
void SetPolarization (const G4ThreeVector &)
 
const G4ThreeVectorGetPolarization () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VDecayChannel
void ClearDaughtersName ()
 
void CheckAndFillDaughters ()
 
void CheckAndFillParent ()
 
G4double DynamicalMass (G4double massPDG, G4double width, G4double maxDev=1.0) const
 
 G4VDecayChannel ()
 
 G4VDecayChannel (const G4VDecayChannel &)
 
G4VDecayChanneloperator= (const G4VDecayChannel &)
 
- 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 42 of file G4ECDecay.hh.

Constructor & Destructor Documentation

◆ G4ECDecay()

G4ECDecay::G4ECDecay ( const G4ParticleDefinition theParentNucleus,
const G4double theBR,
const G4double Qvalue,
const G4double excitation,
const G4Ions::G4FloatLevelBase flb,
const G4RadioactiveDecayMode mode 
)

Definition at line 47 of file G4ECDecay.cc.

52 : G4NuclearDecay("electron capture", mode, excitationE, flb), transitionQ(Qvalue),
53 applyARM(true)
54{
55 SetParent(theParentNucleus); // Store name of parent nucleus, delete G4MT_parent
56 SetBR(branch);
57
59 G4IonTable* theIonTable =
61 G4int daughterZ = theParentNucleus->GetAtomicNumber() - 1;
62 G4int daughterA = theParentNucleus->GetAtomicMass();
63 SetDaughter(0, theIonTable->GetIon(daughterZ, daughterA, excitationE, flb) );
64 SetDaughter(1, "nu_e");
65 DefineSubshellProbabilities(daughterZ, daughterZ);
66}
int G4int
Definition: G4Types.hh:85
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
G4int GetAtomicNumber() const
G4int GetAtomicMass() const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
void SetBR(G4double value)
void SetNumberOfDaughters(G4int value)
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
void SetParent(const G4ParticleDefinition *particle_type)

◆ ~G4ECDecay()

G4ECDecay::~G4ECDecay ( )
virtual

Definition at line 69 of file G4ECDecay.cc.

70{}

Member Function Documentation

◆ DecayIt()

G4DecayProducts * G4ECDecay::DecayIt ( G4double  )
virtual

Implements G4VDecayChannel.

Definition at line 73 of file G4ECDecay.cc.

74{
75 // Fill G4MT_parent with theParentNucleus (stored by SetParent in ctor)
77
78 // Fill G4MT_daughters with alpha and residual nucleus (stored by SetDaughter)
80
81 // Get shell number of captured electron
82 G4int shellIndex = -1;
83 G4double ran;
84 switch (theMode)
85 {
86 case KshellEC:
87 shellIndex = 0;
88 break;
89 case LshellEC: // PL1+PL2+PL3=1
90 ran=G4UniformRand();
91 if (ran <= PL1) shellIndex =1;
92 else if (ran<= (PL1+PL2)) shellIndex =2;
93 else shellIndex =3;
94 break;
95 case MshellEC: // PM1+PM2+PM3=1
96 ran=G4UniformRand();
97 if (ran < PM1) shellIndex =4;
98 else if (ran< (PM1+PM2)) shellIndex =5;
99 else shellIndex = 6;
100 break;
101 case NshellEC: // PN1+PN2+PN3=1
102 ran=G4UniformRand();
103 if (ran < PN1) shellIndex = 9;
104 else if (ran<= (PN1+PN2)) shellIndex =2;
105 else shellIndex =10;
106 break;
107 default:
108 G4Exception("G4ECDecay::DecayIt()", "HAD_RDM_009",
109 FatalException, "Invalid electron shell selected");
110 }
111
112 // Initialize decay products with parent nucleus at rest
113 G4DynamicParticle parentParticle(G4MT_parent, G4ThreeVector(0,0,0), 0.0);
114 G4DecayProducts* products = new G4DecayProducts(parentParticle);
115 G4double eBind = 0.0;
116
117 // G4LossTableManager must already be initialized with G4UAtomicDeexcitation
118 // This is currently done in G4RadioactiveDecay::BuildPhysicsTable
119 G4VAtomDeexcitation* atomDeex =
121 std::vector<G4DynamicParticle*> armProducts;
122
123 if (applyARM) {
124 if (nullptr != atomDeex) {
127 if (shellIndex >= nShells) shellIndex = nShells;
129 const G4AtomicShell* shell = atomDeex->GetAtomicShell(aZ, as);
130 eBind = shell->BindingEnergy();
131 if (atomDeex->IsFluoActive() && aZ > 5 && aZ < 105) {
132 // Do atomic relaxation
133 // VI, SI
134 // Allows fixing of Bugzilla 1727
135 //const G4double deexLimit = 0.1*keV;
136 G4double deexLimit = 0.1*keV;
137 if (G4EmParameters::Instance()->DeexcitationIgnoreCut()) deexLimit = 0.;
138
139 atomDeex->GenerateParticles(&armProducts, shell, aZ, deexLimit, deexLimit);
140 }
141
142 G4double productEnergy = 0.;
143 for (std::size_t i = 0; i < armProducts.size(); ++i) {
144 productEnergy += armProducts[i]->GetKineticEnergy();
145 }
146 G4double deficit = shell->BindingEnergy() - productEnergy;
147 if (deficit > 0.0) {
148 // Add a dummy electron to make up extra energy
149 G4double cosTh = 1.-2.*G4UniformRand();
150 G4double sinTh = std::sqrt(1.- cosTh*cosTh);
151 G4double phi = twopi*G4UniformRand();
152
153 G4ThreeVector electronDirection(sinTh*std::sin(phi),
154 sinTh*std::cos(phi), cosTh);
155 G4DynamicParticle* extra =
156 new G4DynamicParticle(G4Electron::Electron(), electronDirection,
157 deficit);
158 armProducts.push_back(extra);
159 }
160 } // atomDeex
161 } // applyARM
162
163 G4double daughterMass = G4MT_daughters[0]->GetPDGMass();
164
165 // CM momentum using Q value corrected for binding energy of captured electron
166 G4double Q = transitionQ - eBind;
167
168 // Negative transitionQ values for some rare nuclides require this
169 // Absolute values in these cases are small
170 if (Q < 0.0) Q = 0.0;
171
172 G4double cmMomentum = Q*(Q + 2.*daughterMass)/(Q + daughterMass)/2.;
173
174 G4double costheta = 2.*G4UniformRand() - 1.0;
175 G4double sintheta = std::sqrt(1.0 - costheta*costheta);
176 G4double phi = twopi*G4UniformRand()*rad;
177 G4ThreeVector direction(sintheta*std::cos(phi),sintheta*std::sin(phi),
178 costheta);
179 G4double KE = cmMomentum;
180 G4DynamicParticle* daughterParticle =
181 new G4DynamicParticle(G4MT_daughters[1], direction, KE, 0.0);
182 products->PushProducts(daughterParticle);
183
184 KE = std::sqrt(cmMomentum*cmMomentum + daughterMass*daughterMass) - daughterMass;
185 daughterParticle =
186 new G4DynamicParticle(G4MT_daughters[0], -1.0*direction, KE, daughterMass);
187 products->PushProducts(daughterParticle);
188
189 std::size_t nArm = armProducts.size();
190 if (nArm > 0) {
191 G4ThreeVector bst = daughterParticle->Get4Momentum().boostVector();
192 for (std::size_t i = 0; i < nArm; ++i) {
193 G4DynamicParticle* dp = armProducts[i];
194 G4LorentzVector lv = dp->Get4Momentum().boost(bst);
195 dp->Set4Momentum(lv);
196 products->PushProducts(dp);
197 }
198 }
199
200 // Energy conservation check
201 /*
202 G4int newSize = products->entries();
203 G4DynamicParticle* temp = 0;
204 G4double KEsum = 0.0;
205 for (G4int i = 0; i < newSize; i++) {
206 temp = products->operator[](i);
207 KEsum += temp->GetKineticEnergy();
208 }
209
210 G4double eCons = (transitionQ - KEsum)/keV;
211 G4cout << " EC check: Ediff (keV) = " << eCons << G4endl;
212 */
213 return products;
214}
G4AtomicShellEnumerator
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
#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)
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4EmParameters * Instance()
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4RadioactiveDecayMode theMode
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
G4ParticleDefinition ** G4MT_daughters
G4ParticleDefinition * G4MT_parent
void CheckAndFillDaughters()

◆ DumpNuclearInfo()

void G4ECDecay::DumpNuclearInfo ( )
virtual

Implements G4NuclearDecay.

Definition at line 217 of file G4ECDecay.cc.

218{
219 G4cout << " G4ECDecay of parent nucleus " << GetParentName() << " from ";
220 if (theMode == 3) {
221 G4cout << "K shell";
222 } else if (theMode == 4) {
223 G4cout << "L shell";
224 } else if (theMode == 5) {
225 G4cout << "M shell";
226 }
227 else if (theMode == 6) {
228 G4cout << "N shell";
229 }
230 G4cout << G4endl;
231 G4cout << " to " << GetDaughterName(0) << " + " << GetDaughterName(1)
232 << " with branching ratio " << GetBR() << "% and Q value "
233 << transitionQ << G4endl;
234}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetBR() const
const G4String & GetParentName() const
const G4String & GetDaughterName(G4int anIndex) const

◆ SetARM()

void G4ECDecay::SetARM ( G4bool  onoff)
inline

Definition at line 56 of file G4ECDecay.hh.

56{applyARM = onoff;}

Referenced by G4RadioactiveDecay::LoadDecayTable().


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