Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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$
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 // ------------------------------------------------------------
39 
40 #include "G4ParticleDefinition.hh"
41 #include "G4PhysicalConstants.hh"
42 #include "G4SystemOfUnits.hh"
43 #include "G4DecayProducts.hh"
44 #include "G4VDecayChannel.hh"
46 #include "Randomize.hh"
47 #include "G4RotationMatrix.hh"
48 #include "G4LorentzVector.hh"
49 #include "G4LorentzRotation.hh"
50 
52  :G4VDecayChannel(),
53  aENuCorr(-0.102)
54 {
55 }
56 
58  const G4String& theParentName,
59  G4double theBR)
60  :G4VDecayChannel("Neutron Decay"),
61  aENuCorr(-0.102)
62 {
63  // set names for daughter particles
64  if (theParentName == "neutron") {
65  SetBR(theBR);
66  SetParent("neutron");
68  SetDaughter(0, "e-");
69  SetDaughter(1, "anti_nu_e");
70  SetDaughter(2, "proton");
71  } else if (theParentName == "anti_neutron") {
72  SetBR(theBR);
73  SetParent("anti_neutron");
75  SetDaughter(0, "e+");
76  SetDaughter(1, "nu_e");
77  SetDaughter(2, "anti_proton");
78  } else {
79 #ifdef G4VERBOSE
80  if (GetVerboseLevel()>0) {
81  G4cout << "G4NeutronBetaDecayChannel:: constructor :";
82  G4cout << " parent particle is not neutron but ";
83  G4cout << theParentName << G4endl;
84  }
85 #endif
86  }
87 }
88 
90 {
91 }
92 
94  : G4VDecayChannel(right),
95  aENuCorr(-0.102)
96 {
97 }
98 
99 
101 {
102  if (this != &right) {
104  verboseLevel = right.verboseLevel;
105  rbranch = right.rbranch;
106 
107  // copy parent name
108  parent_name = new G4String(*right.parent_name);
109 
110  // clear daughters_name array
112 
113  // recreate array
115  if ( numberOfDaughters >0 ) {
118  //copy daughters name
119  for (G4int index=0; index < numberOfDaughters; index++) {
121  }
122  }
123  }
124  return *this;
125 }
126 
128 {
129  // This class describes free neutron beta decay kinemtics.
130  // This version neglects neutron/electron polarization
131  // without Coulomb effect
132 
133 #ifdef G4VERBOSE
134  if (GetVerboseLevel()>1) G4cout << "G4NeutronBetaDecayChannel::DecayIt ";
135 #endif
136 
137  if (parent == 0) FillParent();
138  if (daughters == 0) FillDaughters();
139 
140  // parent mass
141  G4double parentmass = parent->GetPDGMass();
142 
143  //daughters'mass
144  G4double daughtermass[3];
145  G4double sumofdaughtermass = 0.0;
146  for (G4int index=0; index<3; index++){
147  daughtermass[index] = daughters[index]->GetPDGMass();
148  sumofdaughtermass += daughtermass[index];
149  }
150  G4double xmax = parentmass-sumofdaughtermass;
151 
152  //create parent G4DynamicParticle at rest
153  G4ThreeVector dummy;
154  G4DynamicParticle * parentparticle = new G4DynamicParticle( parent, dummy, 0.0);
155 
156  //create G4Decayproducts
157  G4DecayProducts *products = new G4DecayProducts(*parentparticle);
158  delete parentparticle;
159 
160  // calculate daughter momentum
161  G4double daughtermomentum[3];
162 
163  // calcurate electron energy
164  G4double x; // Ee
165  G4double p; // Pe
166  G4double dm = daughtermass[0]; //Me
167  G4double w; // cosine of e-nu angle
168  G4double r;
169  G4double r0;
170  do {
171  x = xmax*G4UniformRand();
172  p = std::sqrt(x*(x+2.0*dm));
173  w = 1.0-2.0*G4UniformRand();
174  r = p*(x+dm)*(xmax-x)*(xmax-x)*(1.0+aENuCorr*p/(x+dm)*w);
175  r0 = G4UniformRand()*(xmax+dm)*(xmax+dm)*xmax*xmax*(1.0+aENuCorr);
176  } while (r < r0);
177 
178 
179  //create daughter G4DynamicParticle
180  // rotation materix to lab frame
181  G4double costheta = 2.*G4UniformRand()-1.0;
182  G4double theta = std::acos(costheta)*rad;
184  G4RotationMatrix rm;
185  rm.rotateY(theta);
186  rm.rotateZ(phi);
187 
188  // daughter 0 (electron) in Z direction
189  daughtermomentum[0] = p;
190  G4ThreeVector direction0(0.0, 0.0, 1.0);
191  direction0 = rm * direction0;
192  G4DynamicParticle * daughterparticle0
193  = new G4DynamicParticle( daughters[0], direction0*daughtermomentum[0]);
194  products->PushProducts(daughterparticle0);
195 
196  // daughter 1 (nutrino) in XZ plane
197  G4double eNu; // Enu
198  eNu = (parentmass-daughtermass[2])*(parentmass+daughtermass[2])+(dm*dm)-2.*parentmass*(x+dm);
199  eNu /= 2.*(parentmass+p*w-(x+dm));
200  G4double cosn = w;
201  G4double sinn = std::sqrt((1.0-cosn)*(1.0+cosn));
202 
203  G4ThreeVector direction1(sinn, 0.0, cosn);
204  direction1 = rm * direction1;
205  G4DynamicParticle * daughterparticle1
206  = new G4DynamicParticle( daughters[1], direction1*eNu);
207  products->PushProducts(daughterparticle1);
208 
209  // daughter 2 (proton) at REST
210  G4double eP; // Eproton
211  eP = parentmass-eNu-(x+dm)-daughtermass[2];
212  G4double pPx = -eNu*sinn;
213  G4double pPz = -p-eNu*cosn;
214  G4double pP = std::sqrt(eP*(eP+2.*daughtermass[2]));
215  G4ThreeVector direction2(pPx/pP, 0.0, pPz/pP);
216  G4DynamicParticle * daughterparticle2
217  = new G4DynamicParticle( daughters[2], direction2);
218  products->PushProducts(daughterparticle2);
219 
220 
221  // output message
222 #ifdef G4VERBOSE
223  if (GetVerboseLevel()>1) {
224  G4cout << "G4NeutronBetaDecayChannel::DecayIt ";
225  G4cout << " create decay products in rest frame " <<G4endl;
226  products->DumpInfo();
227  }
228 #endif
229  return products;
230 }
231 
232 
233 
234 
235 
236