Geant4  10.02
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 "G4PhotonEvaporation.hh"
76 
77 #include "G4VAtomDeexcitation.hh"
78 #include "G4AtomicShells.hh"
79 #include "G4LossTableManager.hh"
80 
81 //Model const parameters
84 //const G4bool G4NuclearDecayChannel:: FermiOn = true;
85 
86 //This is a kind of "cache"
88 
89 // Constructor for one decay product (the nucleus)
92  G4int Verbose,
93  const G4ParticleDefinition* theParentNucleus,
94  const G4double theBR,
95  const G4double theQtransition,
96  const G4int A, const G4int Z,
97  const G4double theDaughterExcitation)
98  :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode),
99  Qtransition(theQtransition), RandomEnergy(0)
100 {
101 #ifdef G4VERBOSE
102  if (GetVerboseLevel() > 1) {
103  G4cout << "G4NuclearDecayChannel constructor for " << G4int(theMode)
104  << G4endl;
105  }
106 #endif
107  SetParent(theParentNucleus);
108  FillParent();
109  G4MT_parent_mass = theParentNucleus->GetPDGMass();
110  SetBR(theBR);
112  FillDaughterNucleus(0, A, Z, theDaughterExcitation);
114  applyICM = true;
115  applyARM = true;
116  FillDaughters();
117 }
118 
119 // Constructor for a daughter nucleus and one other particle.
120 //
123  G4int Verbose,
124  const G4ParticleDefinition* theParentNucleus,
125  const G4double theBR,
126  const G4double theQtransition,
127  const G4int A, const G4int Z,
128  const G4double theDaughterExcitation,
129  const G4String theDaughterName1)
130  :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode),
131  Qtransition(theQtransition), RandomEnergy(0)
132 {
133 #ifdef G4VERBOSE
134  if (GetVerboseLevel() > 1) {
135  G4cout << "G4NuclearDecayChannel constructor for " << G4int(theMode)
136  << G4endl;
137  }
138 #endif
139  SetParent(theParentNucleus);
140  FillParent();
141  G4MT_parent_mass = theParentNucleus->GetPDGMass();
142  SetBR(theBR);
144  SetDaughter(0, theDaughterName1);
145  FillDaughterNucleus(1, A, Z, theDaughterExcitation);
147  applyICM = true;
148  applyARM = true;
149  FillDaughters();
150 }
151 
152 // Constructor for a daughter nucleus and two other particles
153 //
156  G4int Verbose,
157  const G4ParticleDefinition *theParentNucleus,
158  const G4double theBR,
159  G4double /* theFFN */,
160  G4bool /* betaS */,
161  G4RandGeneral* randBeta,
162  const G4double theQtransition,
163  const G4int A, const G4int Z,
164  const G4double theDaughterExcitation,
165  const G4String theDaughterName1,
166  const G4String theDaughterName2)
167  :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode),
168  Qtransition(theQtransition)
169 {
170 #ifdef G4VERBOSE
171  if (GetVerboseLevel() > 1) {
172  G4cout << "G4NuclearDecayChannel constructor for " << G4int(theMode)
173  << G4endl;
174  }
175 #endif
176  SetParent(theParentNucleus);
177  FillParent();
178  G4MT_parent_mass = theParentNucleus->GetPDGMass();
179  SetBR (theBR);
181  SetDaughter(0, theDaughterName1);
182  SetDaughter(2, theDaughterName2);
183  FillDaughterNucleus(1, A, Z, theDaughterExcitation);
184  RandomEnergy = randBeta;
185 
187  applyICM = true;
188  applyARM = true;
189  FillDaughters();
190 }
191 
193 {}
194 
196  const G4double theDaughterExcitation)
197 {
198  // Determine if the proposed daughter nucleus has a sensible A, Z and
199  // excitation energy.
200  if (A < 1 || Z < 0 || theDaughterExcitation < 0.0) {
202  ed << "Inappropriate values of daughter A, Z or excitation: "
203  << A << " , " << Z << " , " << theDaughterExcitation*MeV << " MeV "
204  << G4endl;
205  G4Exception("G4NuclearDecayChannel::FillDaughterNucleus()", "HAD_RDM_006",
206  FatalException, ed);
207  }
208 
209  // Save A and Z to local variables. Find the GROUND STATE of the daughter
210  // nucleus and save this, as an ion, in the array of daughters.
211  daughterA = A;
212  daughterZ = Z;
213  if (Z == 1 && A == 1) {
215  } else if (Z == 0 && A == 1) {
217  } else {
218  G4IonTable *theIonTable =
220  // GetIon with only Z and A arguments returns ground state
221  daughterNucleus = theIonTable->GetIon(daughterZ, daughterA);
222  }
223 
224  daughterExcitation = theDaughterExcitation;
226 }
227 
229 {
230  G4double deltaM;
231  if (decayMode == 1) { // beta- decay
232  deltaM = CLHEP::electron_mass_c2;
233  } else if (decayMode == 2) { // beta+ decay
234  deltaM = 2.*CLHEP::electron_mass_c2;
235  } else if (decayMode < 6 && decayMode > 2) { // EC
236  deltaM = -CLHEP::electron_mass_c2;
237  } else { // all others
238  deltaM = 0.0;
239  }
240 
241  // Mass available for decay in rest frame of parent after correcting for
242  // the appropriate number of electron masses and reserving the daughter
243  // excitation energy to be applied later
244  G4double massOfParent = G4MT_parent->GetPDGMass(); // PDG mass includes excitation energy
245  SetParentMass(massOfParent - deltaM - daughterExcitation);
246 
247  // Define a product vector.
248  G4DecayProducts* products = 0;
249 
250  // Depending upon the number of daughters, select the appropriate decay
251  // kinematics scheme.
252  switch (numberOfDaughters) {
253  case 0:
254  {
256  ed << " No daughters defined " << G4endl;
257  G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_005",
258  JustWarning, ed);
259  }
260  break;
261  case 1:
262  products = OneBodyDecayIt();
263  break;
264  case 2:
265  products = TwoBodyDecayIt();
266  break;
267  case 3:
268  products = BetaDecayIt();
269  break;
270  default:
271  {
273  ed << " More than 3 daughters in decay: N = " << numberOfDaughters
274  << G4endl;
275  G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_007",
276  FatalException, ed);
277  }
278  }
279 
280  if (products == 0) {
282  ed << " Parent nucleus " << *parent_name << " was not decayed " << G4endl;
283  G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_008",
284  JustWarning, ed);
285  DumpInfo();
286  } else {
287  // If the decay is to an excited state of the daughter nuclide, the photon
288  // evaporation process must be applied.
289 
290  // Need to hold the shell idex after ICM
291  G4int shellIndex = -1;
292  if (daughterExcitation > 0.0) {
293  // Pop the daughter nucleus off the product vector to get its 4-momentum
294  dynamicDaughter = products->PopProducts();
295  G4LorentzVector daughterMomentum = dynamicDaughter->Get4Momentum();
296  if (dynamicDaughter) delete dynamicDaughter;
297 
298  // Using daughter nucleus, set up a G4Fragment for photon evaporatation
299  if (decayMode == 0) {
300  G4double exe = ((const G4Ions*)(G4MT_parent))->GetExcitationEnergy();
301  daughterMomentum.setE(daughterMomentum.e() + exe);
302  }
303  G4Fragment nucleus(daughterA, daughterZ, daughterMomentum);
304 
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 
337  } else {
338  // The fragment vector from photon evaporation contains the list of
339  // evaporated gammas, some of which may have been replaced by conversion
340  // electrons. The last element is the residual nucleus.
341  // Note: try photoEvapProducts as a name instead of gammas
342  G4FragmentVector* gammas = deexcitation->BreakUp(nucleus);
343  G4int nFrags = G4int(gammas->size());
344  G4double eOrGammaEnergy = 0.0;
345 
346  if (nFrags < 1) {
348  ed << nFrags << " No fragments produced by photon evaporation. " << G4endl;
349  G4Exception("G4NuclearDecayChannel::DecayIt()","HAD_RDM_012",
350  FatalException, ed);
351  } else if (nFrags > 1) {
352  // Add gamma/e- to the decay product. The angular distribution of this
353  // particle is assumed to be isotropic
354  G4Fragment* eOrGamma;
355  G4DynamicParticle* eOrGammaDyn;
356  for (G4int i = 0; i < nFrags - 1; i++) {
357  eOrGamma = gammas->operator[](i);
358  eOrGammaDyn = new G4DynamicParticle(eOrGamma->GetParticleDefinition(),
359  eOrGamma->GetMomentum() );
360  eOrGammaDyn->SetProperTime(eOrGamma->GetCreationTime() );
361  products->PushProducts(eOrGammaDyn);
362  eOrGammaEnergy += eOrGamma->GetMomentum().e();
363  }
364  }
365 
366  G4double finalDaughterExcitation =
367  gammas->operator[](nFrags-1)->GetExcitationEnergy();
368  if (finalDaughterExcitation <= 1.0*keV) finalDaughterExcitation = 0;
369 
370  // Get new ion with excitation energy reduced by emitted gamma energy
371  daughterIon =
372  theIonTable->GetIon(daughterZ, daughterA, finalDaughterExcitation);
373  daughterMomentum.setE(daughterMomentum.e() - eOrGammaEnergy);
374 
375  // Delete/reset variables associated with the gammas.
376  while (!gammas->empty() ) { /* Loop checking, 01.09.2015, D.Wright */
377  delete *(gammas->end()-1);
378  gammas->pop_back();
379  }
380  delete gammas;
381  } // end if decayMode == 0
382 
383  G4ThreeVector const daughterMomentum1(static_cast<const G4LorentzVector> (daughterMomentum));
384  dynamicDaughter = new G4DynamicParticle(daughterIon, daughterMomentum1);
385  products->PushProducts(dynamicDaughter);
386 
387  // retrieve the ICM shell index
388  shellIndex = deexcitation->GetVacantShellNumber();
389 
390  delete deexcitation;
391  } // if daughter excitation > 0
392 
393  // Now take care of the EC products which have to go through the ARM
394  G4int eShell = -1;
395  if (decayMode == 3 || decayMode == 4 || decayMode == 5) {
396  switch (decayMode)
397  {
398  case KshellEC:
399  {
400  eShell = 0; // --> 0 from 1 (f.lei 30/4/2008)
401  }
402  break;
403  case LshellEC:
404  {
405  eShell = G4int(G4UniformRand()*3)+1;
406  }
407  break;
408  case MshellEC:
409  {
410  // limit the shell index to 6 as specified by the ARM (F.Lei 06/05/2010)
411  // eShell = G4int(G4UniformRand()*5)+4;
412  eShell = G4int(G4UniformRand()*3)+4;
413  }
414  break;
415  case RDM_ERROR:
416  default:
417  G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_009",
418  FatalException, "Incorrect decay mode selection");
419  }
420  } else {
421  // For other cases eShell comes from shellIndex resulting from the photo decay
422  // modeled by G4PhotonEvaporation* de-excitation (see above)
423  eShell = shellIndex;
424  }
425 
426  // now apply ARM if it is requested and there is a vaccancy
427  if (applyARM && eShell != -1) {
428  G4int aZ = daughterZ;
429 
430  // V.Ivanchenko migration to new interface to atomic deexcitation
431  // no check on index of G4MaterialCutsCouple, simplified
432  // check on secondary energy Esec < 0.1 keV
433  G4VAtomDeexcitation* atomDeex =
435  if (atomDeex) {
436  if(atomDeex->IsFluoActive() && aZ > 5 && aZ < 100) { // only applies to 5< Z <100
437  if (eShell >= G4AtomicShells::GetNumberOfShells(aZ)){
438  eShell = G4AtomicShells::GetNumberOfShells(aZ)-1;
439  }
441  const G4AtomicShell* shell = atomDeex->GetAtomicShell(aZ, as);
442  std::vector<G4DynamicParticle*> armProducts;
443 
444  // VI, SI
445  // Allows fixing of Bugzilla 1727
446  //const G4double deexLimit = 0.1*keV;
447  G4double deexLimit = 0.1*keV;
448  if (G4EmParameters::Instance()->DeexcitationIgnoreCut()) deexLimit =0.;
449  //
450 
451  atomDeex->GenerateParticles(&armProducts, shell, aZ, deexLimit, deexLimit);
452  size_t narm = armProducts.size();
453  if (narm > 0) {
454  // L.Desorgher: need to initialize dynamicDaughter in some decay
455  // cases (for example Hg194)
456  dynamicDaughter = products->PopProducts();
457  G4ThreeVector bst = dynamicDaughter->Get4Momentum().boostVector();
458  for (size_t i = 0; i<narm; ++i) {
459  G4DynamicParticle* dp = armProducts[i];
460  G4LorentzVector lv = dp->Get4Momentum().boost(bst);
461  dp->Set4Momentum(lv);
462  products->PushProducts(dp);
463  }
464  products->PushProducts(dynamicDaughter);
465  }
466  }
467  }
468  }
469  } // Parent nucleus decayed
470  /*
471  if (atomDeex && aZ > 5 && aZ < 100) { // only applies to 5< Z <100
472  // Retrieve the corresponding identifier and binding energy of the selected shell
473  const G4AtomicTransitionManager* transitionManager = G4AtomicTransitionManager::Instance();
474 
475  //check that eShell is smaller than the max number of shells
476  //this to avoid a warning from transitionManager(otherwise the output is the same)
477  //Correction to Bug 1662. L Desorgher (04/02/2011)
478  if (eShell >= transitionManager->NumberOfShells(aZ)){
479  eShell=transitionManager->NumberOfShells(aZ)-1;
480  }
481 
482  const G4AtomicShell* shell = transitionManager->Shell(aZ, eShell);
483  G4double bindingEnergy = shell->BindingEnergy();
484  G4int shellId = shell->ShellId();
485 
486  G4AtomicDeexcitation* atomDeex = new G4AtomicDeexcitation();
487  //the default is no Auger electron generation.
488  // Switch it on/off here!
489  atomDeex->ActivateAugerElectronProduction(true);
490  std::vector<G4DynamicParticle*>* armProducts = atomDeex->GenerateParticles(aZ,shellId);
491 
492  // pop up the daughter before insertion;
493  // f.lei (30/04/2008) check if the total kinetic energy is less than
494  // the shell binding energy; if true add the difference to the daughter to conserve the energy
495  dynamicDaughter = products->PopProducts();
496  G4double tARMEnergy = 0.0;
497  for (size_t i = 0; i < armProducts->size(); i++) {
498  products->PushProducts ((*armProducts)[i]);
499  tARMEnergy += (*armProducts)[i]->GetKineticEnergy();
500  }
501  if ((bindingEnergy - tARMEnergy) > 0.1*keV){
502  G4double dEnergy = dynamicDaughter->GetKineticEnergy() + (bindingEnergy - tARMEnergy);
503  dynamicDaughter->SetKineticEnergy(dEnergy);
504  }
505  products->PushProducts(dynamicDaughter);
506 
507 #ifdef G4VERBOSE
508  if (GetVerboseLevel()>0)
509  {
510  G4cout <<"G4NuclearDecayChannel::Selected shell number for ARM = " <<shellId <<G4endl;
511  G4cout <<"G4NuclearDecayChannel::ARM products = " <<armProducts->size()<<G4endl;
512  G4cout <<" The binding energy = " << bindingEnergy << G4endl;
513  G4cout <<" Total ARM particle kinetic energy = " << tARMEnergy << G4endl;
514  }
515 #endif
516 
517  delete armProducts;
518  delete atomDeex;
519  }
520  }
521  */
522 
523  return products;
524 }
525 
527 {
528  G4double pmass = GetParentMass();
529 
530  G4double daughtermass[3];
531  for (G4int index = 0; index < 3; index++) {
532  daughtermass[index] = G4MT_daughters[index]->GetPDGMass();
533  }
534 
535  // Add excitation energy so daughter can be decayed later
536  daughtermass[1] += daughterExcitation;
537 
538  // Create parent G4DynamicParticle at rest and create products
539  G4ParticleMomentum dummy;
540  G4DynamicParticle parentParticle(G4MT_parent, dummy, 0.0);
541  G4DecayProducts* products = new G4DecayProducts(parentParticle);
542 
543  // faster method as suggested by Dirk Kruecker of FZ-Julich
544  G4double daughtermomentum[3];
545  G4double daughterenergy[3];
546  // Use the histogram distribution to generate the beta energy
547  // 0 = electron, 1 = daughter, 2 = neutrino
548  daughterenergy[0] = Qtransition*RandomEnergy->shoot(G4Random::getTheEngine());
549  daughtermomentum[0] = std::sqrt(daughterenergy[0]*(daughterenergy[0] + 2.*daughtermass[0]) );
550 
551  // neutrino energy distribution is flat within the kinematical limits
552  G4double rd = 2.*G4UniformRand() - 1.;
553  // limits
554  G4double Mme = daughtermass[1] + Qtransition;
555  G4double K = 0.5 - daughtermass[1]*daughtermass[1]/(2*Mme*Mme-4*pmass*daughterenergy[0]);
556 
557  daughterenergy[2] = K * (Mme - daughterenergy[0] + rd*daughtermomentum[0]);
558  daughtermomentum[2] = daughterenergy[2];
559 
560  // the recoil nucleus
561  daughterenergy[1] = Qtransition - daughterenergy[0] - daughterenergy[2];
562  G4double recoilmomentumsquared =
563  daughterenergy[1]*(daughterenergy[1] + 2.0*daughtermass[1]);
564  if (recoilmomentumsquared < 0.0) recoilmomentumsquared = 0.0;
565  daughtermomentum[1] = std::sqrt(recoilmomentumsquared);
566 
567  // output message
568  if (GetVerboseLevel()>1) {
569  G4cout << " G4NuclearDecayChannel::BetaDecayIt() " << G4endl;
570  G4cout <<" e- momentum: " <<daughtermomentum[0]/GeV <<" [GeV/c]" <<G4endl;
571  G4cout <<" daughter momentum: " <<daughtermomentum[1]/GeV <<" [GeV/c]" <<G4endl;
572  G4cout <<" nu momentum: " <<daughtermomentum[2]/GeV <<" [GeV/c]" <<G4endl;
573  G4cout <<" e- energy: " << daughtermass[0] + daughterenergy[0] << G4endl;
574  G4cout <<" daughter energy: " << daughtermass[1] + daughterenergy[1] << G4endl;
575  G4cout <<" nu energy: " << daughtermass[2] + daughterenergy[2] << G4endl;
576  G4cout <<" total of daughter energies: " << daughtermass[0] + daughtermass[1] +
577  daughtermass[2] + daughterenergy[0] + daughterenergy[1] + daughterenergy[2]
578  << G4endl;
579  }
580  //create daughter G4DynamicParticle
581  G4double costheta, sintheta, phi, sinphi, cosphi;
582  G4double costhetan, sinthetan, phin, sinphin, cosphin;
583  costheta = 2.*G4UniformRand()-1.0;
584  sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
585  phi = twopi*G4UniformRand()*rad;
586  sinphi = std::sin(phi);
587  cosphi = std::cos(phi);
588 // electron chosen isotropically
589  G4ParticleMomentum direction0(sintheta*cosphi,sintheta*sinphi,costheta);
590  G4DynamicParticle * daughterparticle
591  = new G4DynamicParticle( G4MT_daughters[0], direction0*daughtermomentum[0]);
592  products->PushProducts(daughterparticle);
593  // cos of angle between electron and neutrino
594  costhetan = (daughtermomentum[1]*daughtermomentum[1]-
595  daughtermomentum[2]*daughtermomentum[2]-
596  daughtermomentum[0]*daughtermomentum[0])/
597  (2.0*daughtermomentum[2]*daughtermomentum[0]);
598 
599  if (costhetan > 1.) costhetan = 1.;
600  if (costhetan < -1.) costhetan = -1.;
601  sinthetan = std::sqrt((1.0-costhetan)*(1.0+costhetan));
602  phin = twopi*G4UniformRand()*rad;
603  sinphin = std::sin(phin);
604  cosphin = std::cos(phin);
605  G4ParticleMomentum direction2;
606  direction2.setX(sinthetan*cosphin*costheta*cosphi -
607  sinthetan*sinphin*sinphi + costhetan*sintheta*cosphi);
608  direction2.setY(sinthetan*cosphin*costheta*sinphi +
609  sinthetan*sinphin*cosphi + costhetan*sintheta*sinphi);
610  direction2.setZ(-sinthetan*cosphin*sintheta + costhetan*costheta);
611  daughterparticle = new G4DynamicParticle(G4MT_daughters[2],
612  direction2*(daughtermomentum[2]/direction2.mag()));
613  products->PushProducts(daughterparticle);
614  // daughter nucleus p = - (p_e + p_nu )
615  daughterparticle =
617  (direction0*daughtermomentum[0] +
618  direction2*(daughtermomentum[2]/direction2.mag()))*(-1.0));
619  products->PushProducts(daughterparticle);
620 
621  if (GetVerboseLevel()>1) {
622  G4cout << "G4NuclearDecayChannel::BetaDecayIt ";
623  G4cout << " create decay products in rest frame " <<G4endl;
624  products->DumpInfo();
625  }
626  return products;
627 }
virtual void SetICM(G4bool)
static const double MeV
Definition: G4SIunits.hh:211
void SetBR(G4double value)
static G4LossTableManager * Instance()
G4bool IsFluoActive() const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
CLHEP::Hep3Vector G4ThreeVector
virtual void RDMForced(G4bool)
G4int PushProducts(G4DynamicParticle *aParticle)
G4ParticleDefinition * G4MT_parent
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:491
G4ParticleDefinition ** G4MT_daughters
const G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Fragment.hh:395
G4double G4MT_parent_mass
const G4RadioactiveDecayMode decayMode
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
virtual G4FragmentVector * BreakUp(const G4Fragment &theNucleus)
static G4Proton * Definition()
Definition: G4Proton.cc:49
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4IonTable * GetIonTable() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
Definition: G4Ions.hh:51
G4double GetCreationTime() const
Definition: G4Fragment.hh:405
bool G4bool
Definition: G4Types.hh:79
void DumpInfo() const
#define G4RandGeneral
Definition: Randomize.hh:94
void SetNumberOfDaughters(G4int value)
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:284
G4RadioactiveDecayMode
static const double twopi
Definition: G4SIunits.hh:75
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:63
G4String * parent_name
static const double GeV
Definition: G4SIunits.hh:214
G4LorentzVector Get4Momentum() const
G4int GetVacantShellNumber() 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:148
void Set4Momentum(const G4LorentzVector &momentum)
static const double nanosecond
Definition: G4SIunits.hh:157
G4int GetVerboseLevel() const
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
G4NuclearDecayChannel(const G4RadioactiveDecayMode &theMode, G4int Verbose)
void SetParent(const G4ParticleDefinition *particle_type)
static G4EmParameters * Instance()
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:213
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)