Geant4  10.01.p02
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 //
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 #include "G4BetaFermiFunction.hh"
76 #include "G4PhotonEvaporation.hh"
77 
78 #include "G4VAtomDeexcitation.hh"
79 #include "G4AtomicShells.hh"
80 #include "G4LossTableManager.hh"
81 
82 //Model const parameters
85 //const G4bool G4NuclearDecayChannel:: FermiOn = true;
86 
87 //This is a kind of "cache"
89 
90 // Constructor for one decay product (the nucleus)
93  G4int Verbose,
94  const G4ParticleDefinition* theParentNucleus,
95  const G4double theBR,
96  const G4double theQtransition,
97  const G4int A, const G4int Z,
98  const G4double theDaughterExcitation)
99  :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode),
100  Qtransition(theQtransition), RandomEnergy(0)
101 {
102 #ifdef G4VERBOSE
103  if (GetVerboseLevel() > 1) {
104  G4cout << "G4NuclearDecayChannel constructor for " << G4int(theMode)
105  << G4endl;
106  }
107 #endif
108  SetParent(theParentNucleus);
109  FillParent();
110  G4MT_parent_mass = theParentNucleus->GetPDGMass();
111  SetBR(theBR);
113  FillDaughterNucleus(0, A, Z, theDaughterExcitation);
115  applyICM = true;
116  applyARM = true;
117  FillDaughters();
118 }
119 
120 // Constructor for a daughter nucleus and one other particle.
121 //
124  G4int Verbose,
125  const G4ParticleDefinition* theParentNucleus,
126  const G4double theBR,
127  const G4double theQtransition,
128  const G4int A, const G4int Z,
129  const G4double theDaughterExcitation,
130  const G4String theDaughterName1)
131  :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode),
132  Qtransition(theQtransition), RandomEnergy(0)
133 {
134 #ifdef G4VERBOSE
135  if (GetVerboseLevel() > 1) {
136  G4cout << "G4NuclearDecayChannel constructor for " << G4int(theMode)
137  << G4endl;
138  }
139 #endif
140  SetParent(theParentNucleus);
141  FillParent();
142  G4MT_parent_mass = theParentNucleus->GetPDGMass();
143  SetBR(theBR);
145  SetDaughter(0, theDaughterName1);
146  FillDaughterNucleus(1, A, Z, theDaughterExcitation);
148  applyICM = true;
149  applyARM = true;
150  FillDaughters();
151 }
152 
153 // Constructor for a daughter nucleus and two other particles
154 //
157  G4int Verbose,
158  const G4ParticleDefinition *theParentNucleus,
159  const G4double theBR,
160  G4double /* theFFN */,
161  G4bool /* betaS */,
162  G4RandGeneral* randBeta,
163  const G4double theQtransition,
164  const G4int A, const G4int Z,
165  const G4double theDaughterExcitation,
166  const G4String theDaughterName1,
167  const G4String theDaughterName2)
168  :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode),
169  Qtransition(theQtransition)
170 {
171 #ifdef G4VERBOSE
172  if (GetVerboseLevel() > 1) {
173  G4cout << "G4NuclearDecayChannel constructor for " << G4int(theMode)
174  << G4endl;
175  }
176 #endif
177  SetParent(theParentNucleus);
178  FillParent();
179  G4MT_parent_mass = theParentNucleus->GetPDGMass();
180  SetBR (theBR);
182  SetDaughter(0, theDaughterName1);
183  SetDaughter(2, theDaughterName2);
184  FillDaughterNucleus(1, A, Z, theDaughterExcitation);
185  RandomEnergy = randBeta;
186 
188  applyICM = true;
189  applyARM = true;
190  FillDaughters();
191 }
192 
194 {}
195 
197  const G4double theDaughterExcitation)
198 {
199  // Determine if the proposed daughter nucleus has a sensible A, Z and
200  // excitation energy.
201  if (A < 1 || Z < 0 || theDaughterExcitation < 0.0) {
203  ed << "Inappropriate values of daughter A, Z or excitation: "
204  << A << " , " << Z << " , " << theDaughterExcitation*MeV << " MeV "
205  << G4endl;
206  G4Exception("G4NuclearDecayChannel::FillDaughterNucleus()", "HAD_RDM_006",
207  FatalException, ed);
208  }
209 
210  // Save A and Z to local variables. Find the GROUND STATE of the daughter
211  // nucleus and save this, as an ion, in the array of daughters.
212  daughterA = A;
213  daughterZ = Z;
214  if (Z == 1 && A == 1) {
216  } else if (Z == 0 && A == 1) {
218  } else {
219  G4IonTable *theIonTable =
221  // GetIon with only Z and A arguments returns ground state
222  daughterNucleus = theIonTable->GetIon(daughterZ, daughterA);
223  }
224 
225  daughterExcitation = theDaughterExcitation;
227 }
228 
230 {
231  G4double deltaM;
232  if (decayMode == 1) { // beta- decay
233  deltaM = CLHEP::electron_mass_c2;
234  } else if (decayMode == 2) { // beta+ decay
235  deltaM = 2.*CLHEP::electron_mass_c2;
236  } else if (decayMode < 6 && decayMode > 2) { // EC
237  deltaM = -CLHEP::electron_mass_c2;
238  } else { // all others
239  deltaM = 0.0;
240  }
241 
242  // Mass available for decay in rest frame of parent after correcting for
243  // the appropriate number of electron masses and reserving the daughter
244  // excitation energy to be applied later
245  G4double massOfParent = G4MT_parent->GetPDGMass(); // PDG mass includes excitation energy
246  SetParentMass(massOfParent - deltaM - daughterExcitation);
247 
248  // Define a product vector.
249  G4DecayProducts* products = 0;
250 
251  // Depending upon the number of daughters, select the appropriate decay
252  // kinematics scheme.
253  switch (numberOfDaughters) {
254  case 0:
255  {
257  ed << " No daughters defined " << G4endl;
258  G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_005",
259  JustWarning, ed);
260  }
261  break;
262  case 1:
263  products = OneBodyDecayIt();
264  break;
265  case 2:
266  products = TwoBodyDecayIt();
267  break;
268  case 3:
269  products = BetaDecayIt();
270  break;
271  default:
272  {
274  ed << " More than 3 daughters in decay: N = " << numberOfDaughters
275  << G4endl;
276  G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_007",
277  FatalException, ed);
278  }
279  }
280 
281  if (products == 0) {
283  ed << " Parent nucleus " << *parent_name << " was not decayed " << G4endl;
284  G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_008",
285  JustWarning, ed);
286  DumpInfo();
287  } else {
288  // If the decay is to an excited state of the daughter nuclide, the photon
289  // evaporation process must be applied.
290 
291  // Need to hold the shell idex after ICM
292  G4int shellIndex = -1;
293  if (daughterExcitation > 0.0) {
294  // Pop the daughter nucleus off the product vector to get its 4-momentum
295  dynamicDaughter = products->PopProducts();
296  G4LorentzVector daughterMomentum = dynamicDaughter->Get4Momentum();
297  if (dynamicDaughter) delete dynamicDaughter;
298 
299  // Using daughter nucleus, set up a G4Fragment for photon evaporatation
300  if (decayMode == 0) {
301  G4double exe = ((const G4Ions*)(G4MT_parent))->GetExcitationEnergy();
302  daughterMomentum.setE(daughterMomentum.e() + exe);
303  }
304  G4Fragment nucleus(daughterA, daughterZ, daughterMomentum);
305  G4PhotonEvaporation* deexcitation = new G4PhotonEvaporation;
306  deexcitation->SetVerboseLevel(GetVerboseLevel());
307 
308  // switch on/off internal electron conversion
309  deexcitation->SetICM(applyICM);
310 
311  // In IT mode, we need to force the transition
312  if (decayMode == 0) {
313  deexcitation->RDMForced(true);
314  // at this point, de-excitation will occur even if IT state is long-lived (>1ns)
315  // Why does it need to be forced?
316  } else {
317  // Not forced, but decay will still happen if lifetime < 1 ns,
318  // otherwise no gamma decay is performed
319  deexcitation->RDMForced(false);
320  }
321 
323  // //
324  // Apply photon evaporation if Isomeric Transition is indicated. //
325  // Use G4PhotonEvaporation::BreakUp() which does only one transition. //
326  // This allows IC to be done. //
327  // //
329 
330  G4IonTable* theIonTable =
332  G4ParticleDefinition* daughterIon = 0;
333 
334  if (decayMode != 0) {
335  daughterIon = theIonTable->GetIon(daughterZ, daughterA, daughterExcitation);
336  } else {
337  // The fragment vector from photon evaporation contains the list of
338  // evaporated gammas, some of which may have been replaced by conversion
339  // electrons. The last element is the residual nucleus.
340  // Note: try photoEvapProducts as a name instead of gammas
341  G4FragmentVector* gammas = deexcitation->BreakUp(nucleus);
342  G4int nFrags = G4int(gammas->size());
343  G4double eOrGammaEnergy = 0.0;
344 
345  if (nFrags < 1) {
347  ed << nFrags << " No fragments produced by photon evaporation. " << G4endl;
348  G4Exception("G4NuclearDecayChannel::DecayIt()","HAD_RDM_012",
349  FatalException, ed);
350  } else if (nFrags > 1) {
351  // Add gamma/e- to the decay product. The angular distribution of this
352  // particle is assumed to be isotropic
353  G4Fragment* eOrGamma;
354  G4DynamicParticle* eOrGammaDyn;
355  for (G4int i = 0; i < nFrags - 1; i++) {
356  eOrGamma = gammas->operator[](i);
357  eOrGammaDyn = new G4DynamicParticle(eOrGamma->GetParticleDefinition(),
358  eOrGamma->GetMomentum() );
359  eOrGammaDyn->SetProperTime(eOrGamma->GetCreationTime() );
360  products->PushProducts(eOrGammaDyn);
361  eOrGammaEnergy += eOrGamma->GetMomentum().e();
362  }
363  }
364 
365  G4double finalDaughterExcitation =
366  gammas->operator[](nFrags-1)->GetExcitationEnergy();
367  if (finalDaughterExcitation <= 1.0*keV) finalDaughterExcitation = 0;
368 
369  // Get new ion with excitation energy reduced by emitted gamma energy
370  daughterIon =
371  theIonTable->GetIon(daughterZ, daughterA, finalDaughterExcitation);
372  daughterMomentum.setE(daughterMomentum.e() - eOrGammaEnergy);
373 
374  // Delete/reset variables associated with the gammas.
375  while (!gammas->empty() ) {
376  delete *(gammas->end()-1);
377  gammas->pop_back();
378  }
379  delete gammas;
380  } // end if decayMode == 0
381 
382  G4ThreeVector const daughterMomentum1(static_cast<const G4LorentzVector> (daughterMomentum));
383  dynamicDaughter = new G4DynamicParticle(daughterIon, daughterMomentum1);
384  products->PushProducts(dynamicDaughter);
385 
386  // retrieve the ICM shell index
387  shellIndex = deexcitation->GetVacantShellNumber();
388 
389  delete deexcitation;
390  } // if daughter excitation > 0
391 
392  // Now take care of the EC products which have to go through the ARM
393  G4int eShell = -1;
394  if (decayMode == 3 || decayMode == 4 || decayMode == 5) {
395  switch (decayMode)
396  {
397  case KshellEC:
398  {
399  eShell = 0; // --> 0 from 1 (f.lei 30/4/2008)
400  }
401  break;
402  case LshellEC:
403  {
404  eShell = G4int(G4UniformRand()*3)+1;
405  }
406  break;
407  case MshellEC:
408  {
409  // limit the shell index to 6 as specified by the ARM (F.Lei 06/05/2010)
410  // eShell = G4int(G4UniformRand()*5)+4;
411  eShell = G4int(G4UniformRand()*3)+4;
412  }
413  break;
414  case RDM_ERROR:
415  default:
416  G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_009",
417  FatalException, "Incorrect decay mode selection");
418  }
419  } else {
420  // For other cases eShell comes from shellIndex resulting from the photo decay
421  // modeled by G4PhotonEvaporation* de-excitation (see above)
422  eShell = shellIndex;
423  }
424 
425  // now apply ARM if it is requested and there is a vaccancy
426  if (applyARM && eShell != -1) {
427  G4int aZ = daughterZ;
428 
429  // V.Ivanchenko migration to new interface to atomic deexcitation
430  // no check on index of G4MaterialCutsCouple, simplified
431  // check on secondary energy Esec < 0.1 keV
432  G4VAtomDeexcitation* atomDeex =
434  if (atomDeex) {
435  if(atomDeex->IsFluoActive() && aZ > 5 && aZ < 100) { // only applies to 5< Z <100
436  if (eShell >= G4AtomicShells::GetNumberOfShells(aZ)){
437  eShell = G4AtomicShells::GetNumberOfShells(aZ)-1;
438  }
440  const G4AtomicShell* shell = atomDeex->GetAtomicShell(aZ, as);
441  std::vector<G4DynamicParticle*> armProducts;
442  const G4double deexLimit = 0.1*keV;
443  atomDeex->GenerateParticles(&armProducts, shell, aZ, deexLimit, deexLimit);
444  size_t narm = armProducts.size();
445  if (narm > 0) {
446  // L.Desorgher: need to initialize dynamicDaughter in some decay
447  // cases (for example Hg194)
448  dynamicDaughter = products->PopProducts();
449  G4ThreeVector bst = dynamicDaughter->Get4Momentum().boostVector();
450  for (size_t i = 0; i<narm; ++i) {
451  G4DynamicParticle* dp = armProducts[i];
452  G4LorentzVector lv = dp->Get4Momentum().boost(bst);
453  dp->Set4Momentum(lv);
454  products->PushProducts(dp);
455  }
456  products->PushProducts(dynamicDaughter);
457  }
458  }
459  }
460  }
461  } // Parent nucleus decayed
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 
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  // 0 = electron, 1 = daughter, 2 = neutrino
540  daughterenergy[0] = Qtransition*RandomEnergy->shoot(G4Random::getTheEngine());
541  daughtermomentum[0] = std::sqrt(daughterenergy[0]*(daughterenergy[0] + 2.*daughtermass[0]) );
542 
543  // neutrino energy distribution is flat within the kinematical limits
544  G4double rd = 2.*G4UniformRand() - 1.;
545  // limits
546  G4double Mme = daughtermass[1] + Qtransition;
547  G4double K = 0.5 - daughtermass[1]*daughtermass[1]/(2*Mme*Mme-4*pmass*daughterenergy[0]);
548 
549  daughterenergy[2] = K * (Mme - daughterenergy[0] + rd*daughtermomentum[0]);
550  daughtermomentum[2] = daughterenergy[2];
551 
552  // the recoil nucleus
553  daughterenergy[1] = Qtransition - daughterenergy[0] - daughterenergy[2];
554  G4double recoilmomentumsquared =
555  daughterenergy[1]*(daughterenergy[1] + 2.0*daughtermass[1]);
556  if (recoilmomentumsquared < 0.0) recoilmomentumsquared = 0.0;
557  daughtermomentum[1] = std::sqrt(recoilmomentumsquared);
558 
559  // output message
560  if (GetVerboseLevel()>1) {
561  G4cout << " G4NuclearDecayChannel::BetaDecayIt() " << G4endl;
562  G4cout <<" e- momentum: " <<daughtermomentum[0]/GeV <<" [GeV/c]" <<G4endl;
563  G4cout <<" daughter momentum: " <<daughtermomentum[1]/GeV <<" [GeV/c]" <<G4endl;
564  G4cout <<" nu momentum: " <<daughtermomentum[2]/GeV <<" [GeV/c]" <<G4endl;
565  G4cout <<" e- energy: " << daughtermass[0] + daughterenergy[0] << G4endl;
566  G4cout <<" daughter energy: " << daughtermass[1] + daughterenergy[1] << G4endl;
567  G4cout <<" nu energy: " << daughtermass[2] + daughterenergy[2] << G4endl;
568  G4cout <<" total of daughter energies: " << daughtermass[0] + daughtermass[1] +
569  daughtermass[2] + daughterenergy[0] + daughterenergy[1] + daughterenergy[2]
570  << G4endl;
571  }
572  //create daughter G4DynamicParticle
573  G4double costheta, sintheta, phi, sinphi, cosphi;
574  G4double costhetan, sinthetan, phin, sinphin, cosphin;
575  costheta = 2.*G4UniformRand()-1.0;
576  sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
577  phi = twopi*G4UniformRand()*rad;
578  sinphi = std::sin(phi);
579  cosphi = std::cos(phi);
580 // electron chosen isotropically
581  G4ParticleMomentum direction0(sintheta*cosphi,sintheta*sinphi,costheta);
582  G4DynamicParticle * daughterparticle
583  = new G4DynamicParticle( G4MT_daughters[0], direction0*daughtermomentum[0]);
584  products->PushProducts(daughterparticle);
585  // cos of angle between electron and neutrino
586  costhetan = (daughtermomentum[1]*daughtermomentum[1]-
587  daughtermomentum[2]*daughtermomentum[2]-
588  daughtermomentum[0]*daughtermomentum[0])/
589  (2.0*daughtermomentum[2]*daughtermomentum[0]);
590 
591  if (costhetan > 1.) costhetan = 1.;
592  if (costhetan < -1.) costhetan = -1.;
593  sinthetan = std::sqrt((1.0-costhetan)*(1.0+costhetan));
594  phin = twopi*G4UniformRand()*rad;
595  sinphin = std::sin(phin);
596  cosphin = std::cos(phin);
597  G4ParticleMomentum direction2;
598  direction2.setX(sinthetan*cosphin*costheta*cosphi -
599  sinthetan*sinphin*sinphi + costhetan*sintheta*cosphi);
600  direction2.setY(sinthetan*cosphin*costheta*sinphi +
601  sinthetan*sinphin*cosphi + costhetan*sintheta*sinphi);
602  direction2.setZ(-sinthetan*cosphin*sintheta + costhetan*costheta);
603  daughterparticle = new G4DynamicParticle(G4MT_daughters[2],
604  direction2*(daughtermomentum[2]/direction2.mag()));
605  products->PushProducts(daughterparticle);
606  // daughter nucleus p = - (p_e + p_nu )
607  daughterparticle =
609  (direction0*daughtermomentum[0] +
610  direction2*(daughtermomentum[2]/direction2.mag()))*(-1.0));
611  products->PushProducts(daughterparticle);
612 
613  if (GetVerboseLevel()>1) {
614  G4cout << "G4NuclearDecayChannel::BetaDecayIt ";
615  G4cout << " create decay products in rest frame " <<G4endl;
616  products->DumpInfo();
617  }
618  return products;
619 }
static const double MeV
Definition: G4SIunits.hh:193
void SetBR(G4double value)
static G4LossTableManager * Instance()
G4bool IsFluoActive() const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
CLHEP::Hep3Vector G4ThreeVector
G4int PushProducts(G4DynamicParticle *aParticle)
G4ParticleDefinition * G4MT_parent
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:463
G4ParticleDefinition ** G4MT_daughters
const G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Fragment.hh:403
G4double G4MT_parent_mass
const G4RadioactiveDecayMode decayMode
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
static G4Proton * Definition()
Definition: G4Proton.cc:49
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4IonTable * GetIonTable() const
#define G4UniformRand()
Definition: Randomize.hh:93
G4GLOB_DLL std::ostream G4cout
Definition: G4Ions.hh:51
G4double GetCreationTime() const
Definition: G4Fragment.hh:413
bool G4bool
Definition: G4Types.hh:79
void DumpInfo() const
#define G4RandGeneral
Definition: Randomize.hh:90
void SetNumberOfDaughters(G4int value)
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:276
G4RadioactiveDecayMode
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
virtual G4FragmentVector * BreakUp(const G4Fragment &nucleus)
G4String * parent_name
static const double GeV
Definition: G4SIunits.hh:196
static const G4double A[nN]
G4LorentzVector Get4Momentum() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const double rad
Definition: G4SIunits.hh:130
void Set4Momentum(const G4LorentzVector &momentum)
static const double nanosecond
Definition: G4SIunits.hh:139
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)
G4DecayProducts * BetaDecayIt()
static G4Neutron * Definition()
Definition: G4Neutron.cc:54
G4ParticleDefinition * daughterNucleus
#define G4endl
Definition: G4ios.hh:61
static G4ThreadLocal G4DynamicParticle * dynamicDaughter
G4VAtomDeexcitation * AtomDeexcitation()
static const double keV
Definition: G4SIunits.hh:195
double G4double
Definition: G4Types.hh:76
static const G4double levelTolerance
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
G4ThreeVector G4ParticleMomentum
void SetProperTime(G4double)
void SetParentMass(const G4double aParentMass)
G4AtomicShellEnumerator
G4DecayProducts * DecayIt(G4double)
void SetVerboseLevel(G4int verbose)
static const G4double pTolerance
CLHEP::HepLorentzVector G4LorentzVector
void FillDaughterNucleus(G4int index, G4int A, G4int Z, const G4double theDaughterExcitation)
static G4int GetNumberOfShells(G4int Z)