Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ITDecay Class Reference

#include <G4ITDecay.hh>

Inheritance diagram for G4ITDecay:
Collaboration diagram for G4ITDecay:

Public Member Functions

 G4ITDecay (const G4ParticleDefinition *theParentNucleus, const G4double &theBR, const G4double &Qvalue, const G4double &excitation)
 
virtual ~G4ITDecay ()
 
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 ()
 
- 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="")
 
virtual ~G4VDecayChannel ()
 
G4int operator== (const G4VDecayChannel &right) const
 
G4int operator!= (const G4VDecayChannel &right) const
 
G4int 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)
 
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 ()
 
 G4VDecayChannel ()
 
 G4VDecayChannel (const G4VDecayChannel &)
 
G4VDecayChanneloperator= (const G4VDecayChannel &)
 
G4double DynamicalMass (G4double massPDG, G4double width, G4double maxDev=+1.) const
 
- Protected Attributes inherited from G4NuclearDecay
const G4RadioactiveDecayMode theMode
 
- Protected Attributes inherited from G4VDecayChannel
G4String kinematics_name
 
G4double rbranch
 
G4int numberOfDaughters
 
G4Stringparent_name
 
G4String ** daughters_name
 
G4double rangeMass
 
G4ThreeVector parent_polarization
 
G4ParticleTableparticletable
 
G4int verboseLevel
 
G4ParticleDefinitionG4MT_parent
 
G4ParticleDefinition ** G4MT_daughters
 
G4double G4MT_parent_mass
 
G4doubleG4MT_daughters_mass
 
G4doubleG4MT_daughters_width
 
G4Mutex daughtersMutex
 
G4Mutex parentMutex
 
- Static Protected Attributes inherited from G4VDecayChannel
static const G4String noName = " "
 

Detailed Description

Definition at line 46 of file G4ITDecay.hh.

Constructor & Destructor Documentation

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

Definition at line 51 of file G4ITDecay.cc.

54  : G4NuclearDecay("IT decay", IT, excitationE, noFloat), transitionQ(Qvalue),
55  applyARM(true), daughterNucleus(nullptr), photoEvap(0)
56 {
57  SetParent(theParentNucleus); // Store name of parent nucleus, delete G4MT_parent
58  SetBR(branch);
59 
60  parentZ = theParentNucleus->GetAtomicNumber();
61  parentA = theParentNucleus->GetAtomicMass();
62 
64  G4IonTable* theIonTable =
66  SetDaughter(0, theIonTable->GetIon(parentZ, parentA, excitationE, noFloat) );
67 
68  // Let G4PhotonEvaporation do the decay
69  photoEvap = new G4PhotonEvaporation;
70  photoEvap->RDMForced(true);
71  photoEvap->SetICM(true);
72 }
void SetBR(G4double value)
virtual void RDMForced(G4bool)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:503
G4int GetAtomicNumber() const
G4IonTable * GetIonTable() const
void SetNumberOfDaughters(G4int value)
G4int GetAtomicMass() const
G4NuclearDecay(const G4String &channelName, const G4RadioactiveDecayMode &mode, const G4double &excitation, const G4Ions::G4FloatLevelBase &floatingLevel)
static G4ParticleTable * GetParticleTable()
void SetParent(const G4ParticleDefinition *particle_type)
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
#define noFloat
Definition: G4Ions.hh:118

Here is the call graph for this function:

G4ITDecay::~G4ITDecay ( )
virtual

Definition at line 75 of file G4ITDecay.cc.

76 {
77  delete photoEvap;
78 }

Member Function Documentation

G4DecayProducts * G4ITDecay::DecayIt ( G4double  )
virtual

Implements G4VDecayChannel.

Definition at line 81 of file G4ITDecay.cc.

82 {
83  // Fill G4MT_parent with theParentNucleus (stored by SetParent in ctor)
85 
86  // Set up final state
87  // parentParticle is set at rest here because boost with correct momentum
88  // is done later
90  G4ThreeVector(0.,0.,0.) );
91  G4DynamicParticle parentParticle(G4MT_parent, atRest);
92  G4DecayProducts* products = new G4DecayProducts(parentParticle);
93 
94  // Let G4PhotonEvaporation do the decay
95  G4Fragment parentNucleus(parentA, parentZ, atRest);
96  // G4cout << " START of G4ITDecay::DecayIt: " << G4endl;
97 
98  // Check if the old nuclear polarization can be used
99  G4NuclearPolarization* nucPol(nullptr);
100  G4int oldZ(0), oldA(0);
101  G4double oldMass(0.0);
102  if(daughterNucleus) {
103  nucPol = daughterNucleus->GetNuclearPolarization();
104  if(nucPol) {
105  oldZ = daughterNucleus->GetZ_asInt();
106  oldA = daughterNucleus->GetA_asInt();
107  oldMass = daughterNucleus->GetGroundStateMass() +
108  daughterNucleus->GetExcitationEnergy();
109  }
110  }
111  static const G4double mlimit = 10*CLHEP::eV;
112  if (nucPol && oldZ == parentZ && oldA == parentA
113  && std::abs(atRest.e() - oldMass) < mlimit) {
114  // Continue with existing chain
115  parentNucleus = *daughterNucleus;
116  parentNucleus.SetMomentum(atRest);
117  }
118 
119  // G4cout << " BEFORE TRANSITION fragment = " << G4endl;
120  // G4cout << parentNucleus << G4endl;
121 
122  G4Fragment* eOrGamma = photoEvap->EmittedFragment(&parentNucleus);
123 
124  // G4cout << " AFTER TRANSITION " << G4endl;
125  // G4cout << parentNucleus << G4endl;
126 
127  // Check if IT chain has ended with excited isomere
128  // Take care for deletion of cached G4Fragment
129  if (parentNucleus.GetExcitationEnergy() > mlimit) {
130  // G4cout << " IT chain is continue " << G4endl;
131  delete daughterNucleus;
132  daughterNucleus = new G4Fragment(parentNucleus);
133  } else if(daughterNucleus) {
134  // G4cout << " End of IT chain " << G4endl;
135  delete daughterNucleus;
136  daughterNucleus = nullptr;
137  }
138 
139  // Modified nuclide is returned as dynDaughter
140  G4IonTable* theIonTable =
142  G4ParticleDefinition* daughterIon =
143  theIonTable->GetIon(parentZ, parentA, parentNucleus.GetExcitationEnergy(),
144  G4Ions::FloatLevelBase(parentNucleus.GetFloatingLevelNumber()));
145  G4DynamicParticle* dynDaughter = new G4DynamicParticle(daughterIon,
146  parentNucleus.GetMomentum());
147 
148  if (eOrGamma) {
149  G4DynamicParticle* eOrGammaDyn =
151  eOrGamma->GetMomentum() );
152  eOrGammaDyn->SetProperTime(eOrGamma->GetCreationTime() );
153  products->PushProducts(eOrGammaDyn);
154  delete eOrGamma;
155 
156  // Now do atomic relaxation if e- is emitted
157  if (applyARM) {
158  G4int shellIndex = photoEvap->GetVacantShellNumber();
159  if (shellIndex > -1) {
160  G4VAtomDeexcitation* atomDeex =
162  if (atomDeex->IsFluoActive() && parentZ > 5 && parentZ < 100) {
163  G4int nShells = G4AtomicShells::GetNumberOfShells(parentZ);
164  if (shellIndex >= nShells) shellIndex = nShells;
166  const G4AtomicShell* shell = atomDeex->GetAtomicShell(parentZ, as);
167  std::vector<G4DynamicParticle*> armProducts;
168 
169  // VI, SI
170  // Allows fixing of Bugzilla 1727
171  G4double deexLimit = 0.1*keV;
172  if (G4EmParameters::Instance()->DeexcitationIgnoreCut()) deexLimit =0.;
173  //
174 
175  atomDeex->GenerateParticles(&armProducts, shell, parentZ, deexLimit,
176  deexLimit);
177  G4double productEnergy = 0.;
178  for (G4int i = 0; i < G4int(armProducts.size()); i++)
179  productEnergy += armProducts[i]->GetKineticEnergy();
180 
181  G4double deficit = shell->BindingEnergy() - productEnergy;
182  if (deficit > 0.0) {
183  // Add a dummy electron to make up extra energy
184  G4double cosTh = 1.-2.*G4UniformRand();
185  G4double sinTh = std::sqrt(1.- cosTh*cosTh);
186  G4double phi = twopi*G4UniformRand();
187 
188  G4ThreeVector electronDirection(sinTh*std::sin(phi),
189  sinTh*std::cos(phi), cosTh);
190  G4DynamicParticle* extra =
191  new G4DynamicParticle(G4Electron::Electron(), electronDirection,
192  deficit);
193  armProducts.push_back(extra);
194  }
195 
196  G4int nArm = armProducts.size();
197  if (nArm > 0) {
198  G4ThreeVector bst = dynDaughter->Get4Momentum().boostVector();
199  for (G4int i = 0; i < nArm; ++i) {
200  G4DynamicParticle* dp = armProducts[i];
201  G4LorentzVector lv = dp->Get4Momentum().boost(bst);
202  dp->Set4Momentum(lv);
203  products->PushProducts(dp);
204  }
205  }
206  }
207  }
208  } // if ARM on
209  } // eOrGamma
210 
211  products->PushProducts(dynDaughter);
212 
213  // Energy conservation check
214  /*
215  G4int newSize = products->entries();
216  G4DynamicParticle* temp = 0;
217  G4double KEsum = 0.0;
218  for (G4int i = 0; i < newSize; i++) {
219  temp = products->operator[](i);
220  KEsum += temp->GetKineticEnergy();
221  }
222  G4double eCons = G4MT_parent->GetPDGMass() - dynDaughter->GetMass() - KEsum;
223  G4cout << " IT check: Ediff (keV) = " << eCons/keV << G4endl;
224  */
225 // delete photoEvap;
226  // G4cout << " END G4ITDecay::DecayIt " << G4endl;
227  return products;
228 }
Hep3Vector boostVector() const
static G4LossTableManager * Instance()
G4bool IsFluoActive() const
CLHEP::Hep3Vector G4ThreeVector
G4int PushProducts(G4DynamicParticle *aParticle)
G4ParticleDefinition * G4MT_parent
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:503
virtual G4Fragment * EmittedFragment(G4Fragment *theNucleus) final
const G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Fragment.hh:438
G4double BindingEnergy() const
int G4int
Definition: G4Types.hh:78
static constexpr double twopi
Definition: G4SIunits.hh:76
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4IonTable * GetIonTable() const
static G4Ions::G4FloatLevelBase FloatLevelBase(char flbChar)
Definition: G4Ions.hh:189
#define G4UniformRand()
Definition: Randomize.hh:97
G4int GetA_asInt() const
Definition: G4Fragment.hh:266
G4double GetCreationTime() const
Definition: G4Fragment.hh:448
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:307
HepLorentzVector & boost(double, double, double)
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:312
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:288
G4LorentzVector Get4Momentum() const
G4int GetVacantShellNumber() const
static constexpr double eV
void Set4Momentum(const G4LorentzVector &momentum)
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
static G4EmParameters * Instance()
G4int GetZ_asInt() const
Definition: G4Fragment.hh:271
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4VAtomDeexcitation * AtomDeexcitation()
double G4double
Definition: G4Types.hh:76
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetProperTime(G4double)
static constexpr double keV
Definition: G4SIunits.hh:216
G4AtomicShellEnumerator
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:283
G4NuclearPolarization * GetNuclearPolarization() const
Definition: G4Fragment.hh:458
static G4int GetNumberOfShells(G4int Z)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ITDecay::DumpNuclearInfo ( )
virtual

Implements G4NuclearDecay.

Definition at line 231 of file G4ITDecay.cc.

232 {
233  G4cout << " G4ITDecay for parent nucleus " << GetParentName() << G4endl;
234  G4cout << " decays to " << GetDaughterName(0)
235  << " + gammas (or electrons), with branching ratio " << GetBR()
236  << "% and Q value " << transitionQ << G4endl;
237 }
G4double GetBR() const
const G4String & GetParentName() const
G4GLOB_DLL std::ostream G4cout
const G4String & GetDaughterName(G4int anIndex) const
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void G4ITDecay::SetARM ( G4bool  onoff)
inline

Definition at line 59 of file G4ITDecay.hh.

59 {applyARM = onoff;}

Here is the caller graph for this function:


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