Geant4  10.02.p01
G4NeutronBetaDecayChannel.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 // $Id: G4NeutronBetaDecayChannel.cc 91896 2015-08-10 09:54:06Z gcosmo $
28 // GEANT4 tag $Name: geant4-09-04-ref-00 $
29 //
30 //
31 // ------------------------------------------------------------
32 // GEANT 4 class header file
33 //
34 // History: first implementation, based on object model of
35 // 18 Sep 2001 H.Kurashige
36 //---
37 // Fix energy of proton and neutrino May 2011 H.Kurashige
38 // Fix direction of proton and neutrino Nov 2013 H.Kurashige
39 // ------------------------------------------------------------
40 
41 #include "G4ParticleDefinition.hh"
42 #include "G4PhysicalConstants.hh"
43 #include "G4SystemOfUnits.hh"
44 #include "G4DecayProducts.hh"
45 #include "G4VDecayChannel.hh"
47 #include "Randomize.hh"
48 #include "G4RotationMatrix.hh"
49 #include "G4LorentzVector.hh"
50 #include "G4LorentzRotation.hh"
51 
53  :G4VDecayChannel(),
54  aENuCorr(-0.102)
55 {
56 }
57 
59  const G4String& theParentName,
60  G4double theBR)
61  :G4VDecayChannel("Neutron Decay"),
62  aENuCorr(-0.102)
63 {
64  // set names for daughter particles
65  if (theParentName == "neutron") {
66  SetBR(theBR);
67  SetParent("neutron");
69  SetDaughter(0, "e-");
70  SetDaughter(1, "anti_nu_e");
71  SetDaughter(2, "proton");
72  } else if (theParentName == "anti_neutron") {
73  SetBR(theBR);
74  SetParent("anti_neutron");
76  SetDaughter(0, "e+");
77  SetDaughter(1, "nu_e");
78  SetDaughter(2, "anti_proton");
79  } else {
80 #ifdef G4VERBOSE
81  if (GetVerboseLevel()>0) {
82  G4cout << "G4NeutronBetaDecayChannel:: constructor :";
83  G4cout << " parent particle is not neutron but ";
84  G4cout << theParentName << G4endl;
85  }
86 #endif
87  }
88 }
89 
91 {
92 }
93 
95  : G4VDecayChannel(right),
96  aENuCorr(-0.102)
97 {
98 }
99 
100 
102 {
103  if (this != &right) {
105  verboseLevel = right.verboseLevel;
106  rbranch = right.rbranch;
107 
108  // copy parent name
109  parent_name = new G4String(*right.parent_name);
110 
111  // clear daughters_name array
113 
114  // recreate array
116  if ( numberOfDaughters >0 ) {
119  //copy daughters name
120  for (G4int index=0; index < numberOfDaughters; index++) {
121  daughters_name[index] = new G4String(*right.daughters_name[index]);
122  }
123  }
124  }
125  return *this;
126 }
127 
129 {
130  // This class describes free neutron beta decay kinemtics.
131  // This version neglects neutron/electron polarization
132  // without Coulomb effect
133 
134 #ifdef G4VERBOSE
135  if (GetVerboseLevel()>1) G4cout << "G4NeutronBetaDecayChannel::DecayIt ";
136 #endif
137 
138  if (G4MT_parent == 0) FillParent();
139  if (G4MT_daughters == 0) FillDaughters();
140 
141  // parent mass
142  G4double parentmass = G4MT_parent->GetPDGMass();
143 
144  //daughters'mass
145  G4double daughtermass[3];
146  G4double sumofdaughtermass = 0.0;
147  for (G4int index=0; index<3; index++){
148  daughtermass[index] = G4MT_daughters[index]->GetPDGMass();
149  sumofdaughtermass += daughtermass[index];
150  }
151  G4double xmax = parentmass-sumofdaughtermass;
152 
153  //create parent G4DynamicParticle at rest
154  G4ThreeVector dummy;
155  G4DynamicParticle * parentparticle = new G4DynamicParticle( G4MT_parent, dummy, 0.0);
156 
157  //create G4Decayproducts
158  G4DecayProducts *products = new G4DecayProducts(*parentparticle);
159  delete parentparticle;
160 
161  // calculate daughter momentum
162  G4double daughtermomentum[3];
163 
164  // calcurate electron energy
165  G4double x; // Ee
166  G4double p; // Pe
167  G4double dm = daughtermass[0]; //Me
168  G4double w; // cosine of e-nu angle
169  G4double r;
170  G4double r0;
171  const size_t MAX_LOOP=10000;
172  for (size_t loop_counter=0; loop_counter <MAX_LOOP; ++loop_counter){
173  x = xmax*G4UniformRand();
174  p = std::sqrt(x*(x+2.0*dm));
175  w = 1.0-2.0*G4UniformRand();
176  r = p*(x+dm)*(xmax-x)*(xmax-x)*(1.0+aENuCorr*p/(x+dm)*w);
177  r0 = G4UniformRand()*(xmax+dm)*(xmax+dm)*xmax*xmax*(1.0+aENuCorr);
178  if (r > r0) break;
179  }
180 
181  //create daughter G4DynamicParticle
182  // rotation materix to lab frame
183  G4double costheta = 2.*G4UniformRand()-1.0;
184  G4double theta = std::acos(costheta)*rad;
186  G4RotationMatrix rm;
187  rm.rotateY(theta);
188  rm.rotateZ(phi);
189 
190  // daughter 0 (electron) in Z direction
191  daughtermomentum[0] = p;
192  G4ThreeVector direction0(0.0, 0.0, 1.0);
193  direction0 = rm * direction0;
194  G4DynamicParticle * daughterparticle0
195  = new G4DynamicParticle( G4MT_daughters[0], direction0*daughtermomentum[0]);
196  products->PushProducts(daughterparticle0);
197 
198  // daughter 1 (nutrino) in XZ plane
199  G4double eNu; // Enu
200  eNu = (parentmass-daughtermass[2])*(parentmass+daughtermass[2])+(dm*dm)-2.*parentmass*(x+dm);
201  eNu /= 2.*(parentmass+p*w-(x+dm));
202  G4double cosn = w;
203  G4double phin = twopi*G4UniformRand()*rad;
204  G4double sinn = std::sqrt((1.0-cosn)*(1.0+cosn));
205 
206  G4ThreeVector direction1(sinn*std::cos(phin), sinn*std::sin(phin), cosn);
207  direction1 = rm * direction1;
208  G4DynamicParticle * daughterparticle1
209  = new G4DynamicParticle( G4MT_daughters[1], direction1*eNu);
210  products->PushProducts(daughterparticle1);
211 
212  // daughter 2 (proton) at REST
213  G4double eP; // Eproton
214  eP = parentmass-eNu-(x+dm)-daughtermass[2];
215  G4double pPx = -eNu*sinn;
216  G4double pPz = -p-eNu*cosn;
217  G4double pP = std::sqrt(eP*(eP+2.*daughtermass[2]));
218  G4ThreeVector direction2(pPx/pP*std::cos(phin), pPx/pP*std::sin(phin), pPz/pP);
219  G4DynamicParticle * daughterparticle2
220  = new G4DynamicParticle( G4MT_daughters[2], direction2);
221  products->PushProducts(daughterparticle2);
222 
223 
224  // output message
225 #ifdef G4VERBOSE
226  if (GetVerboseLevel()>1) {
227  G4cout << "G4NeutronBetaDecayChannel::DecayIt ";
228  G4cout << " create decay products in rest frame " <<G4endl;
229  products->DumpInfo();
230  }
231 #endif
232  return products;
233 }
234 
235 
236 
237 
238 
239 
void SetBR(G4double value)
G4NeutronBetaDecayChannel & operator=(const G4NeutronBetaDecayChannel &)
CLHEP::Hep3Vector G4ThreeVector
CLHEP::HepRotation G4RotationMatrix
G4int PushProducts(G4DynamicParticle *aParticle)
G4ParticleDefinition * G4MT_parent
G4ParticleDefinition ** G4MT_daughters
const G4double w[NPOINTSGL]
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4String kinematics_name
void DumpInfo() const
void SetNumberOfDaughters(G4int value)
static const double twopi
Definition: G4SIunits.hh:75
G4String * parent_name
static const double rad
Definition: G4SIunits.hh:148
G4int GetVerboseLevel() const
G4double GetPDGMass() const
const G4double x[NPOINTSGL]
void SetParent(const G4ParticleDefinition *particle_type)
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4String ** daughters_name
const G4double r0
virtual G4DecayProducts * DecayIt(G4double)