Geant4-11
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// G4NeutronBetaDecayChannel class implementation
27//
28// Author: H.Kurashige, 18 September 2001
29// --------------------------------------------------------------------
30
33#include "G4SystemOfUnits.hh"
34#include "G4DecayProducts.hh"
35#include "G4VDecayChannel.hh"
37#include "Randomize.hh"
38#include "G4RotationMatrix.hh"
39#include "G4LorentzVector.hh"
40#include "G4LorentzRotation.hh"
41
44{
45}
46
48G4NeutronBetaDecayChannel(const G4String& theParentName,
49 G4double theBR)
50 : G4VDecayChannel("Neutron Decay")
51{
52 // set names for daughter particles
53 if (theParentName == "neutron")
54 {
55 SetBR(theBR);
56 SetParent("neutron");
58 SetDaughter(0, "e-");
59 SetDaughter(1, "anti_nu_e");
60 SetDaughter(2, "proton");
61 }
62 else if (theParentName == "anti_neutron")
63 {
64 SetBR(theBR);
65 SetParent("anti_neutron");
67 SetDaughter(0, "e+");
68 SetDaughter(1, "nu_e");
69 SetDaughter(2, "anti_proton");
70 }
71 else
72 {
73#ifdef G4VERBOSE
74 if (GetVerboseLevel()>0)
75 {
76 G4cout << "G4NeutronBetaDecayChannel:: constructor :";
77 G4cout << " parent particle is not neutron but ";
78 G4cout << theParentName << G4endl;
79 }
80#endif
81 }
82}
83
85{
86}
87
90 : G4VDecayChannel(right)
91{
92}
93
96{
97 if (this != &right)
98 {
101 rbranch = right.rbranch;
102
103 // copy parent name
104 delete parent_name;
105 parent_name = new G4String(*right.parent_name);
106
107 // clear daughters_name array
109
110 // recreate array
112 if ( numberOfDaughters >0 )
113 {
115 // copy daughters name
116 for (G4int index=0; index<numberOfDaughters; ++index)
117 {
118 daughters_name[index] = new G4String(*right.daughters_name[index]);
119 }
120 }
121 }
122 return *this;
123}
124
126{
127 // This class describes free neutron beta decay kinematics.
128 // This version neglects neutron/electron polarization
129 // without Coulomb effect
130
131#ifdef G4VERBOSE
132 if (GetVerboseLevel()>1) G4cout << "G4NeutronBetaDecayChannel::DecayIt ";
133#endif
134
137
138 // parent mass
139 G4double parentmass = G4MT_parent->GetPDGMass();
140
141 // daughters'mass
142 G4double daughtermass[3];
143 G4double sumofdaughtermass = 0.0;
144 for (G4int index=0; index<3; ++index)
145 {
146 daughtermass[index] = G4MT_daughters[index]->GetPDGMass();
147 sumofdaughtermass += daughtermass[index];
148 }
149 G4double xmax = parentmass-sumofdaughtermass;
150
151 // create parent G4DynamicParticle at rest
152 G4ThreeVector dummy;
153 G4DynamicParticle* parentparticle
154 = new G4DynamicParticle( G4MT_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 const std::size_t MAX_LOOP=10000;
171 for (std::size_t loop_counter=0; loop_counter<MAX_LOOP; ++loop_counter)
172 {
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 //
184 G4double costheta = 2.*G4UniformRand()-1.0;
185 G4double theta = std::acos(costheta)*rad;
188 rm.rotateY(theta);
189 rm.rotateZ(phi);
190
191 // daughter 0 (electron) in Z direction
192 daughtermomentum[0] = p;
193 G4ThreeVector direction0(0.0, 0.0, 1.0);
194 direction0 = rm * direction0;
195 G4DynamicParticle* daughterparticle0
196 = new G4DynamicParticle(G4MT_daughters[0], direction0*daughtermomentum[0]);
197 products->PushProducts(daughterparticle0);
198
199 // daughter 1 (nutrino) in XZ plane
200 G4double eNu; // Enu
201 eNu = (parentmass-daughtermass[2])
202 * (parentmass+daughtermass[2])+(dm*dm)-2.*parentmass*(x+dm);
203 eNu /= 2.*(parentmass+p*w-(x+dm));
204 G4double cosn = w;
206 G4double sinn = std::sqrt((1.0-cosn)*(1.0+cosn));
207
208 G4ThreeVector direction1(sinn*std::cos(phin), sinn*std::sin(phin), cosn);
209 direction1 = rm * direction1;
210 G4DynamicParticle* daughterparticle1
211 = new G4DynamicParticle( G4MT_daughters[1], direction1*eNu);
212 products->PushProducts(daughterparticle1);
213
214 // daughter 2 (proton) at REST
215 G4double eP; // Eproton
216 eP = parentmass-eNu-(x+dm)-daughtermass[2];
217 G4double pPx = -eNu*sinn;
218 G4double pPz = -p-eNu*cosn;
219 G4double pP = std::sqrt(eP*(eP+2.*daughtermass[2]));
220 G4ThreeVector direction2(pPx/pP*std::cos(phin),
221 pPx/pP*std::sin(phin), pPz/pP);
222 direction2 = rm * direction2;
223 G4DynamicParticle* daughterparticle2
224 = new G4DynamicParticle( G4MT_daughters[2], direction2*pP);
225 products->PushProducts(daughterparticle2);
226
227 // output message
228#ifdef G4VERBOSE
229 if (GetVerboseLevel()>1)
230 {
231 G4cout << "G4NeutronBetaDecayChannel::DecayIt ";
232 G4cout << " create decay products in rest frame " <<G4endl;
233 products->DumpInfo();
234 }
235#endif
236 return products;
237}
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double rad
Definition: G4SIunits.hh:129
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:87
HepRotation & rotateY(double delta)
Definition: Rotation.cc:74
void DumpInfo() const
G4int PushProducts(G4DynamicParticle *aParticle)
G4NeutronBetaDecayChannel & operator=(const G4NeutronBetaDecayChannel &)
virtual G4DecayProducts * DecayIt(G4double)
G4ParticleDefinition ** G4MT_daughters
G4String * parent_name
void SetBR(G4double value)
G4String ** daughters_name
G4int GetVerboseLevel() const
void SetNumberOfDaughters(G4int value)
G4ParticleDefinition * G4MT_parent
void CheckAndFillDaughters()
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
G4String kinematics_name
void SetParent(const G4ParticleDefinition *particle_type)