Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4NuclearDecayChannel.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 // MODULES: G4NuclearDecayChannel.cc
29 //
30 // Version: 0.b.4
31 // Date: 14/04/00
32 // Author: F Lei & P R Truscott
33 // Organisation: DERA UK
34 // Customer: ESA/ESTEC, NOORDWIJK
35 // Contract: 12115/96/JG/NL Work Order No. 3
36 //
37 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
38 //
39 // CHANGE HISTORY
40 // --------------
41 //
42 // 29 February 2000, P R Truscott, DERA UK
43 // 0.b.3 release.
44 //
45 // 18 October 2002, F Lei
46 // modified link metheds in DecayIt() to G4PhotoEvaporation() in order to
47 // use the new Internal Coversion feature.
48 // 13 April 2000, F Lei, DERA UK
49 // Changes made are:
50 // 1) Use PhotonEvaporation instead of DiscreteGammaDeexcitation
51 // 2) verbose control
52 //
53 // 17 October 2011, L. Desorgher
54 // -Allow the atomic relaxation after de-excitation of exited
55 // nuclei even for beta and alpha
56 // decay. Bug found and solution proposed by Ko Abe.
57 // -Set halflifethreshold by default to a negative value
58 //
59 // 20 November 2011, V.Ivanchenko
60 // - Migration to new design of atomic deexcitation
61 //
62 ///////////////////////////////////////////////////////////////////////////////
63 
64 #include "G4PhysicalConstants.hh"
65 #include "G4SystemOfUnits.hh"
66 #include "G4NuclearLevelManager.hh"
67 #include "G4NuclearLevelStore.hh"
68 #include "G4NuclearDecayChannel.hh"
69 #include "G4DynamicParticle.hh"
70 #include "G4DecayProducts.hh"
71 #include "G4DecayTable.hh"
72 #include "G4PhysicsLogVector.hh"
74 #include "G4IonTable.hh"
75 
76 #include "G4BetaFermiFunction.hh"
77 #include "G4PhotonEvaporation.hh"
78 
79 #include "G4VAtomDeexcitation.hh"
80 #include "G4AtomicShells.hh"
81 #include "G4LossTableManager.hh"
82 
83 //#include "G4AtomicTransitionManager.hh"
84 //#include "G4AtomicShell.hh"
85 //#include "G4AtomicDeexcitation.hh"
86 
87 //Model const parameters
90 //const G4bool G4NuclearDecayChannel:: FermiOn = true;
91 
92 //These are a kind of "cache"
94 
95 // Constructor for one decay product (the nucleus).
96 //
99  G4int Verbose,
100  const G4ParticleDefinition* theParentNucleus,
101  G4double theBR,
102  G4double theQtransition,
103  G4int A,
104  G4int Z,
105  G4double theDaughterExcitation)
106  :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode),
107  Qtransition(theQtransition), RandomEnergy(0)
108 {
109 #ifdef G4VERBOSE
110  if (GetVerboseLevel() > 1) {
111  G4cout << "G4NuclearDecayChannel constructor for " << G4int(theMode)
112  << G4endl;
113  }
114 #endif
115  SetParent(theParentNucleus);
116  FillParent();
117  G4MT_parent_mass = theParentNucleus->GetPDGMass();
118  SetBR(theBR);
120  FillDaughterNucleus(0, A, Z, theDaughterExcitation);
122  applyICM = true;
123  applyARM = true;
124  FillDaughters();
125 }
126 
127 // Constructor for a daughter nucleus and one other particle.
128 //
131  G4int Verbose,
132  const G4ParticleDefinition *theParentNucleus,
133  G4double theBR,
134  G4double theQtransition,
135  G4int A,
136  G4int Z,
137  G4double theDaughterExcitation,
138  const G4String theDaughterName1)
139  :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode),
140  Qtransition(theQtransition), RandomEnergy(0)
141 {
142 #ifdef G4VERBOSE
143  if (GetVerboseLevel() > 1) {
144  G4cout << "G4NuclearDecayChannel constructor for " << G4int(theMode)
145  << G4endl;
146  }
147 #endif
148  SetParent(theParentNucleus);
149  FillParent();
150  G4MT_parent_mass = theParentNucleus->GetPDGMass();
151  SetBR(theBR);
153  SetDaughter(0, theDaughterName1);
154  FillDaughterNucleus(1, A, Z, theDaughterExcitation);
156  applyICM = true;
157  applyARM = true;
158  FillDaughters();
159 }
160 
161 // Constructor for a daughter nucleus and two other particles
162 //
165  G4int Verbose,
166  const G4ParticleDefinition *theParentNucleus,
167  G4double theBR,
168  G4double /* theFFN */,
169  G4bool /* betaS */,
170  G4RandGeneral* randBeta,
171  G4double theQtransition,
172  G4int A,
173  G4int Z,
174  G4double theDaughterExcitation,
175  const G4String theDaughterName1,
176  const G4String theDaughterName2)
177  :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode),
178  Qtransition(theQtransition)
179 {
180 #ifdef G4VERBOSE
181  if (GetVerboseLevel() > 1) {
182  G4cout << "G4NuclearDecayChannel constructor for " << G4int(theMode)
183  << G4endl;
184  }
185 #endif
186  SetParent(theParentNucleus);
187  FillParent();
188  G4MT_parent_mass = theParentNucleus->GetPDGMass();
189  SetBR (theBR);
191  SetDaughter(0, theDaughterName1);
192  SetDaughter(2, theDaughterName2);
193  FillDaughterNucleus(1, A, Z, theDaughterExcitation);
194  RandomEnergy = randBeta;
195 
197  applyICM = true;
198  applyARM = true;
199  FillDaughters();
200 }
201 
203 {}
204 
205 void G4NuclearDecayChannel::FillDaughterNucleus(G4int index, G4int A, G4int Z,
206  G4double theDaughterExcitation)
207 {
208  // Determine if the proposed daughter nucleus has a sensible A, Z and
209  // excitation energy.
210  if (A < 1 || Z < 0 || theDaughterExcitation < 0.0) {
212  ed << "Inappropriate values of daughter A, Z or excitation: "
213  << A << " , " << Z << " , " << theDaughterExcitation*MeV << " MeV "
214  << G4endl;
215  G4Exception("G4NuclearDecayChannel::FillDaughterNucleus()", "HAD_RDM_006",
216  FatalException, ed);
217  }
218 
219  // Save A and Z to local variables. Find the GROUND STATE of the daughter
220  // nucleus and save this, as an ion, in the array of daughters.
221  daughterA = A;
222  daughterZ = Z;
223  if (Z == 1 && A == 1) {
225  } else if (Z == 0 && A == 1) {
227  } else {
228  G4IonTable *theIonTable =
230  // GetIon with only Z and A arguments returns ground state
231  daughterNucleus = theIonTable->GetIon(daughterZ, daughterA);
232  }
233 
234  daughterExcitation = theDaughterExcitation;
236 }
237 
239 {
240  G4double deltaM;
241  if (decayMode == 1) { // beta- decay
242  deltaM = CLHEP::electron_mass_c2;
243  } else if (decayMode == 2) { // beta+ decay
244  deltaM = 2.*CLHEP::electron_mass_c2;
245  } else if (decayMode < 6 && decayMode > 2) { // EC
246  deltaM = -CLHEP::electron_mass_c2;
247  } else { // all others
248  deltaM = 0.0;
249  }
250 
251  // Mass available for decay in rest frame of parent after correcting for
252  // the appropriate number of electron masses and reserving the daughter
253  // excitation energy to be applied later
254  G4double massOfParent = G4MT_parent->GetPDGMass(); // PDG mass includes excitation energy
255  SetParentMass(massOfParent - deltaM - daughterExcitation);
256 
257  // Define a product vector.
258  G4DecayProducts* products = 0;
259 
260  // Depending upon the number of daughters, select the appropriate decay
261  // kinematics scheme.
262  switch (numberOfDaughters) {
263  case 0:
264  {
266  ed << " No daughters defined " << G4endl;
267  G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_005",
268  JustWarning, ed);
269  }
270  break;
271  case 1:
272  products = OneBodyDecayIt();
273  break;
274  case 2:
275  products = TwoBodyDecayIt();
276  break;
277  case 3:
278  products = BetaDecayIt();
279  break;
280  default:
281  {
283  ed << " More than 3 daughters in decay: N = " << numberOfDaughters
284  << G4endl;
285  G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_007",
286  FatalException, ed);
287  }
288  }
289  if (products == 0) {
291  ed << " Parent nucleus " << *parent_name << " was not decayed " << G4endl;
292  G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_008",
293  JustWarning, ed);
294  DumpInfo();
295  } else {
296  // If the decay is to an excited state of the daughter nuclide, the photon
297  // evaporation process must be applied.
298 
299  // Need to hold the shell idex after ICM
300  G4int shellIndex = -1;
301  if (daughterExcitation > 0.0) {
302  // Pop the daughter nucleus off the product vector to get its 4-momentum
303  dynamicDaughter = products->PopProducts();
304  G4LorentzVector daughterMomentum = dynamicDaughter->Get4Momentum();
305  if (dynamicDaughter) delete dynamicDaughter;
306 
307  // Using daughter nucleus, set up a G4Fragment for photon evaporatation
308  if (decayMode == 0) {
309  G4double exe = ((const G4Ions*)(G4MT_parent))->GetExcitationEnergy();
310  daughterMomentum.setE(daughterMomentum.e() + exe);
311  }
312  G4Fragment nucleus(daughterA, daughterZ, daughterMomentum);
313  G4PhotonEvaporation* deexcitation = new G4PhotonEvaporation;
314  deexcitation->SetVerboseLevel(GetVerboseLevel());
315 
316  // switch on/off internal electron conversion
317  deexcitation->SetICM(applyICM);
318 
319  // In IT mode, we need to force the transition
320  if (decayMode == 0) {
321  deexcitation->RDMForced(true);
322  // at this point, de-excitation will occur even if IT state is long-lived (>1ns)
323  // Why does it need to be forced?
324  } else {
325  // Not forced, but decay will still happen if lifetime < 1 ns,
326  // otherwise no gamma decay is performed
327  deexcitation->RDMForced(false);
328  }
329 
330  //////////////////////////////////////////////////////////////////////////
331  // //
332  // Apply photon evaporation if Isomeric Transition is indicated. //
333  // Use G4PhotonEvaporation::BreakUp() which does only one transition. //
334  // This allows IC to be done. //
335  // //
336  //////////////////////////////////////////////////////////////////////////
337 
338  G4IonTable* theIonTable =
340  G4ParticleDefinition* daughterIon = 0;
341 
342  if (decayMode != 0) {
343  daughterIon = theIonTable->GetIon(daughterZ, daughterA, daughterExcitation);
344  } else {
345  // The fragment vector from photon evaporation contains the list of
346  // evaporated gammas, some of which may have been replaced by conversion
347  // electrons. The last element is the residual nucleus.
348  G4FragmentVector* gammas = deexcitation->BreakUp(nucleus); // Evaporate only one photon
349  G4int nGammas = gammas->size() - 1;
350  G4double finalDaughterExcitation =
351  gammas->operator[](nGammas)->GetExcitationEnergy();
352  // Go through each gamma/e- and add it to the decay product. The
353  // angular distribution of the gammas is isotropic, and the residual
354  // nucleus is assumed not to have suffered any recoil as a result of
355  // this de-excitation.
356  G4double Egamma = 0.0;
357  for (G4int ig = 0; ig < nGammas; ig++) {
358  G4DynamicParticle* theGammaRay =
359  new G4DynamicParticle(gammas->operator[](ig)->GetParticleDefinition(),
360  gammas->operator[](ig)->GetMomentum());
361  theGammaRay -> SetProperTime(gammas->operator[](ig)->GetCreationTime());
362  products->PushProducts (theGammaRay);
363  Egamma = (gammas->operator[](ig)->GetMomentum()).e();
364  }
365  if (finalDaughterExcitation <= 1.0*keV) finalDaughterExcitation = 0;
366 
367  // Get new ion with excitation energy reduced by emitted gamma energy
368  daughterIon =
369  theIonTable->GetIon(daughterZ, daughterA, finalDaughterExcitation);
370  daughterMomentum.setE(daughterMomentum.e() - Egamma);
371 
372  // Delete/reset variables associated with the gammas.
373  while (!gammas->empty() ) {
374  delete *(gammas->end()-1);
375  gammas->pop_back();
376  }
377  delete gammas;
378 
379  } // end if decayMode == 0
380 
381  G4ThreeVector const daughterMomentum1(static_cast<const G4LorentzVector> (daughterMomentum));
382  dynamicDaughter = new G4DynamicParticle(daughterIon, daughterMomentum1);
383  products->PushProducts(dynamicDaughter);
384 
385  // retrieve the ICM shell index
386  shellIndex = deexcitation->GetVacantShellNumber();
387 
388  delete deexcitation;
389  } // if daughter excitation > 0
390 
391  // Now take care of the EC products which have to go through the ARM
392  G4int eShell = -1;
393  if (decayMode == 3 || decayMode == 4 || decayMode == 5) {
394  switch (decayMode)
395  {
396  case KshellEC:
397  {
398  eShell = 0; // --> 0 from 1 (f.lei 30/4/2008)
399  }
400  break;
401  case LshellEC:
402  {
403  eShell = G4int(G4UniformRand()*3)+1;
404  }
405  break;
406  case MshellEC:
407  {
408  // limit the shell index to 6 as specified by the ARM (F.Lei 06/05/2010)
409  // eShell = G4int(G4UniformRand()*5)+4;
410  eShell = G4int(G4UniformRand()*3)+4;
411  }
412  break;
413  case RDM_ERROR:
414  default:
415  G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_009",
416  FatalException, "Incorrect decay mode selection");
417  }
418  } else {
419  // For other cases eShell comes from shellIndex resulting from the photo decay
420  // modeled by G4PhotonEvaporation* de-excitation (see above)
421  eShell = shellIndex;
422  }
423 
424  // now apply ARM if it is requested and there is a vaccancy
425  if (applyARM && eShell != -1) {
426  G4int aZ = daughterZ;
427 
428  // V.Ivanchenko migration to new interface to atomic deexcitation
429  // no check on index of G4MaterialCutsCouple, simplified
430  // check on secondary energy Esec < 0.1 keV
431  G4VAtomDeexcitation* atomDeex =
433  if (atomDeex) {
434  if(atomDeex->IsFluoActive() && aZ > 5 && aZ < 100) { // only applies to 5< Z <100
435  if (eShell >= G4AtomicShells::GetNumberOfShells(aZ)){
436  eShell = G4AtomicShells::GetNumberOfShells(aZ)-1;
437  }
439  const G4AtomicShell* shell = atomDeex->GetAtomicShell(aZ, as);
440  std::vector<G4DynamicParticle*> armProducts;
441  const G4double deexLimit = 0.1*keV;
442  atomDeex->GenerateParticles(&armProducts, shell, aZ, deexLimit, deexLimit);
443  size_t narm = armProducts.size();
444  if (narm > 0) {
445  // L.Desorgher: need to initialize dynamicDaughter in some decay
446  // cases (for example Hg194)
447  dynamicDaughter = products->PopProducts();
449  for (size_t i = 0; i<narm; ++i) {
450  G4DynamicParticle* dp = armProducts[i];
451  G4LorentzVector lv = dp->Get4Momentum().boost(bst);
452  dp->Set4Momentum(lv);
453  products->PushProducts(dp);
454  }
455  products->PushProducts(dynamicDaughter);
456  }
457  }
458  }
459  }
460  } // Parent nucleus decayed
461  /*
462 
463  if (atomDeex && aZ > 5 && aZ < 100) { // only applies to 5< Z <100
464  // Retrieve the corresponding identifier and binding energy of the selected shell
465  const G4AtomicTransitionManager* transitionManager = G4AtomicTransitionManager::Instance();
466 
467  //check that eShell is smaller than the max number of shells
468  //this to avoid a warning from transitionManager(otherwise the output is the same)
469  //Correction to Bug 1662. L Desorgher (04/02/2011)
470  if (eShell >= transitionManager->NumberOfShells(aZ)){
471  eShell=transitionManager->NumberOfShells(aZ)-1;
472  }
473 
474  const G4AtomicShell* shell = transitionManager->Shell(aZ, eShell);
475  G4double bindingEnergy = shell->BindingEnergy();
476  G4int shellId = shell->ShellId();
477 
478  G4AtomicDeexcitation* atomDeex = new G4AtomicDeexcitation();
479  //the default is no Auger electron generation.
480  // Switch it on/off here!
481  atomDeex->ActivateAugerElectronProduction(true);
482  std::vector<G4DynamicParticle*>* armProducts = atomDeex->GenerateParticles(aZ,shellId);
483 
484  // pop up the daughter before insertion;
485  // f.lei (30/04/2008) check if the total kinetic energy is less than
486  // the shell binding energy; if true add the difference to the daughter to conserve the energy
487  dynamicDaughter = products->PopProducts();
488  G4double tARMEnergy = 0.0;
489  for (size_t i = 0; i < armProducts->size(); i++) {
490  products->PushProducts ((*armProducts)[i]);
491  tARMEnergy += (*armProducts)[i]->GetKineticEnergy();
492  }
493  if ((bindingEnergy - tARMEnergy) > 0.1*keV){
494  G4double dEnergy = dynamicDaughter->GetKineticEnergy() + (bindingEnergy - tARMEnergy);
495  dynamicDaughter->SetKineticEnergy(dEnergy);
496  }
497  products->PushProducts(dynamicDaughter);
498 
499 #ifdef G4VERBOSE
500  if (GetVerboseLevel()>0)
501  {
502  G4cout <<"G4NuclearDecayChannel::Selected shell number for ARM = " <<shellId <<G4endl;
503  G4cout <<"G4NuclearDecayChannel::ARM products = " <<armProducts->size()<<G4endl;
504  G4cout <<" The binding energy = " << bindingEnergy << G4endl;
505  G4cout <<" Total ARM particle kinetic energy = " << tARMEnergy << G4endl;
506  }
507 #endif
508 
509  delete armProducts;
510  delete atomDeex;
511  }
512  }
513  */
514 
515  return products;
516 }
517 
518 G4DecayProducts* G4NuclearDecayChannel::BetaDecayIt()
519 {
520  G4double pmass = GetParentMass();
521 
522  G4double daughtermass[3];
523  for (G4int index = 0; index < 3; index++) {
524  daughtermass[index] = G4MT_daughters[index]->GetPDGMass();
525  }
526 
527  // Add excitation energy so daughter can be decayed later
528  daughtermass[1] += daughterExcitation;
529 
530  // Create parent G4DynamicParticle at rest and create products
531  G4ParticleMomentum dummy;
532  G4DynamicParticle parentParticle(G4MT_parent, dummy, 0.0);
533  G4DecayProducts* products = new G4DecayProducts(parentParticle);
534 
535  // faster method as suggested by Dirk Kruecker of FZ-Julich
536  G4double daughtermomentum[3];
537  G4double daughterenergy[3];
538  // Use the histogram distribution to generate the beta energy
539  daughterenergy[0] = Qtransition*RandomEnergy->shoot(G4Random::getTheEngine());
540  daughtermomentum[0] = std::sqrt(daughterenergy[0]*(daughterenergy[0] + 2.*daughtermass[0]) );
541 
542  // neutrino energy distribution is flat within the kinematical limits
543  G4double rd = 2.*G4UniformRand() - 1.;
544  // limits
545  G4double Mme = daughtermass[1] + Qtransition;
546  G4double K = 0.5 - daughtermass[1]*daughtermass[1]/(2*Mme*Mme-4*pmass*daughterenergy[0]);
547 
548  daughterenergy[2] = K * (Mme - daughterenergy[0] + rd*daughtermomentum[0]);
549  daughtermomentum[2] = daughterenergy[2];
550 
551  // the recoil nucleus
552  daughterenergy[1] = Qtransition - daughterenergy[0] - daughterenergy[2];
553  G4double recoilmomentumsquared =
554  daughterenergy[1]*(daughterenergy[1] + 2.0*daughtermass[1]);
555  if (recoilmomentumsquared < 0.0) recoilmomentumsquared = 0.0;
556  daughtermomentum[1] = std::sqrt(recoilmomentumsquared);
557 
558  // output message
559  if (GetVerboseLevel()>1) {
560  G4cout << " G4NuclearDecayChannel::BetaDecayIt() " << G4endl;
561  G4cout <<" e- momentum: " <<daughtermomentum[0]/GeV <<" [GeV/c]" <<G4endl;
562  G4cout <<" daughter momentum: " <<daughtermomentum[1]/GeV <<" [GeV/c]" <<G4endl;
563  G4cout <<" nu momentum: " <<daughtermomentum[2]/GeV <<" [GeV/c]" <<G4endl;
564  G4cout <<" e- energy: " << daughtermass[0] + daughterenergy[0] << G4endl;
565  G4cout <<" daughter energy: " << daughtermass[1] + daughterenergy[1] << G4endl;
566  G4cout <<" nu energy: " << daughtermass[2] + daughterenergy[2] << G4endl;
567  G4cout <<" total of daughter energies: " << daughtermass[0] + daughtermass[1] +
568  daughtermass[2] + daughterenergy[0] + daughterenergy[1] + daughterenergy[2]
569  << G4endl;
570  }
571  //create daughter G4DynamicParticle
572  G4double costheta, sintheta, phi, sinphi, cosphi;
573  G4double costhetan, sinthetan, phin, sinphin, cosphin;
574  costheta = 2.*G4UniformRand()-1.0;
575  sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
576  phi = twopi*G4UniformRand()*rad;
577  sinphi = std::sin(phi);
578  cosphi = std::cos(phi);
579  G4ParticleMomentum direction0(sintheta*cosphi,sintheta*sinphi,costheta);
580  G4DynamicParticle * daughterparticle
581  = new G4DynamicParticle( G4MT_daughters[0], direction0*daughtermomentum[0]);
582  products->PushProducts(daughterparticle);
583 
584  costhetan = (daughtermomentum[1]*daughtermomentum[1]-
585  daughtermomentum[2]*daughtermomentum[2]-
586  daughtermomentum[0]*daughtermomentum[0])/
587  (2.0*daughtermomentum[2]*daughtermomentum[0]);
588 
589  if (costhetan > 1.) costhetan = 1.;
590  if (costhetan < -1.) costhetan = -1.;
591  sinthetan = std::sqrt((1.0-costhetan)*(1.0+costhetan));
592  phin = twopi*G4UniformRand()*rad;
593  sinphin = std::sin(phin);
594  cosphin = std::cos(phin);
595  G4ParticleMomentum direction2;
596  direction2.setX(sinthetan*cosphin*costheta*cosphi -
597  sinthetan*sinphin*sinphi + costhetan*sintheta*cosphi);
598  direction2.setY(sinthetan*cosphin*costheta*sinphi +
599  sinthetan*sinphin*cosphi + costhetan*sintheta*sinphi);
600  direction2.setZ(-sinthetan*cosphin*sintheta + costhetan*costheta);
601  daughterparticle = new G4DynamicParticle(G4MT_daughters[2],
602  direction2*(daughtermomentum[2]/direction2.mag()));
603  products->PushProducts(daughterparticle);
604 
605  daughterparticle =
607  (direction0*daughtermomentum[0] +
608  direction2*(daughtermomentum[2]/direction2.mag()))*(-1.0));
609  products->PushProducts(daughterparticle);
610 
611  if (GetVerboseLevel()>1) {
612  G4cout << "G4NuclearDecayChannel::BetaDecayIt ";
613  G4cout << " create decay products in rest frame " <<G4endl;
614  products->DumpInfo();
615  }
616  return products;
617 }
Hep3Vector boostVector() const
void SetBR(G4double value)
static G4LossTableManager * Instance()
G4bool IsFluoActive() const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4int PushProducts(G4DynamicParticle *aParticle)
G4ParticleDefinition * G4MT_parent
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:449
G4ParticleDefinition ** G4MT_daughters
G4double G4MT_parent_mass
const G4RadioactiveDecayMode decayMode
#define G4ThreadLocal
Definition: tls.hh:52
int G4int
Definition: G4Types.hh:78
void setY(double)
void setZ(double)
void setX(double)
static G4Proton * Definition()
Definition: G4Proton.cc:49
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4IonTable * GetIonTable() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
Definition: G4Ions.hh:51
bool G4bool
Definition: G4Types.hh:79
void DumpInfo() const
#define G4RandGeneral
Definition: Randomize.hh:84
void SetNumberOfDaughters(G4int value)
HepLorentzVector & boost(double, double, double)
G4RadioactiveDecayMode
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
virtual G4FragmentVector * BreakUp(const G4Fragment &nucleus)
G4String * parent_name
int nanosecond
Definition: hepunit.py:82
G4LorentzVector Get4Momentum() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void Set4Momentum(const G4LorentzVector &momentum)
G4int GetVerboseLevel() const
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
G4NuclearDecayChannel(const G4RadioactiveDecayMode &theMode, G4int Verbose)
void SetParent(const G4ParticleDefinition *particle_type)
G4DynamicParticle * PopProducts()
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
static G4Neutron * Definition()
Definition: G4Neutron.cc:54
G4ParticleDefinition * daughterNucleus
#define G4endl
Definition: G4ios.hh:61
static G4ThreadLocal G4DynamicParticle * dynamicDaughter
G4VAtomDeexcitation * AtomDeexcitation()
double G4double
Definition: G4Types.hh:76
static const G4double levelTolerance
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
double mag() const
void SetParentMass(const G4double aParentMass)
G4AtomicShellEnumerator
G4DecayProducts * DecayIt(G4double)
void SetVerboseLevel(G4int verbose)
static const G4double pTolerance
static G4int GetNumberOfShells(G4int Z)