Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4PionRadiativeDecayChannel.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// G4PionRadiativeDecayChannel class implementation
27// GEANT 4 class header file
28//
29// Author: P.Gumplinger, 30 July 2007
30// Reference: M. Blecher, TRIUMF/PIENU Technote
31// "Inclusion of pi->enug in the Monte Carlo"
32// --------------------------------------------------------------------
33
35
37#include "G4SystemOfUnits.hh"
38#include "Randomize.hh"
39#include "G4DecayProducts.hh"
40#include "G4LorentzVector.hh"
41
42namespace
43{
44 const G4double beta = 3.6612e-03;
45 const G4double cib = 1.16141e-03;
46 const G4double csdp = 3.45055e-02;
47 const G4double csdm = 5.14122e-03;
48 const G4double cif = 4.63543e-05;
49 const G4double cig = 1.78928e-05;
50 const G4double xl = 2.*0.1*MeV/139.57*MeV;
51 const G4double yl = ((1.-xl) + std::sqrt((1-xl)*(1-xl)+4*beta*beta))/2.;
52
53 const G4double xu = 1. - (yl - std::sqrt(yl*yl-4.*beta*beta))/2.;
54 const G4double yu = 1. + beta*beta;
55
56 inline G4double D2W(const G4double x,const G4double y)
57 {
58 return cib*(1.-y)*(1.+((1.-x)*(1.-x)))/((x*x)*(x+y-1.)) +
59 csdp*(1.-x)*((x+y-1.)*(x+y-1.)) +
60 csdm*(1.-x)*((1.-y)*(1.-y)) +
61 cif*(x-1.)*(1.-y)/x +
62 cig*(1.-y)*(1.-x+(x*x)/(x+y-1.))/x;
63 }
64
65 const G4double d2wmax = D2W(xl,yl);
66}
67
68// --------------------------------------------------------------------
71{
72}
73
74// --------------------------------------------------------------------
76G4PionRadiativeDecayChannel(const G4String& theParentName,
77 G4double theBR)
78 : G4VDecayChannel("Radiative Pion Decay", 1)
79{
80 // set names for daughter particles
81 if (theParentName == "pi+")
82 {
83 SetBR(theBR);
84 SetParent("pi+");
86 SetDaughter(0, "e+");
87 SetDaughter(1, "gamma");
88 SetDaughter(2, "nu_e");
89 }
90 else if (theParentName == "pi-")
91 {
92 SetBR(theBR);
93 SetParent("pi-");
95 SetDaughter(0, "e-");
96 SetDaughter(1, "gamma");
97 SetDaughter(2, "anti_nu_e");
98 }
99 else
100 {
101#ifdef G4VERBOSE
102 if (GetVerboseLevel()>0)
103 {
104 G4cout << "G4RadiativePionDecayChannel::G4PionRadiativeDecayChannel()"
105 << G4endl;
106 G4cout << "Parent particle is not charged pion: ";
107 G4cout << theParentName << G4endl;
108 }
109#endif
110 }
111}
112
113// --------------------------------------------------------------------
115{
116}
117
118// --------------------------------------------------------------------
121 : G4VDecayChannel(right)
122{
123}
124
127{
128 if (this != &right)
129 {
132 rbranch = right.rbranch;
133
134 // copy parent name
135 parent_name = new G4String(*right.parent_name);
136
137 // clear daughters_name array
139
140 // recreate array
142 if ( numberOfDaughters >0 )
143 {
144 if (daughters_name != nullptr) ClearDaughtersName();
146 //copy daughters name
147 for (G4int index=0; index<numberOfDaughters; ++index)
148 {
149 daughters_name[index] = new G4String(*right.daughters_name[index]);
150 }
151 }
152 }
153 return *this;
154}
155
156// --------------------------------------------------------------------
158{
159
160#ifdef G4VERBOSE
161 if (GetVerboseLevel()>1)
162 G4cout << "G4PionRadiativeDecayChannel::DecayIt ";
163#endif
164
167
168 // parent mass
169 G4double parentmass = G4MT_parent->GetPDGMass();
170
171 G4double EMPI = parentmass;
172
173 // daughters'mass
174 const G4int N_DAUGHTER=3;
175 G4double daughtermass[N_DAUGHTER];
176 //G4double sumofdaughtermass = 0.0;
177 for (G4int index=0; index<N_DAUGHTER; ++index)
178 {
179 daughtermass[index] = G4MT_daughters[index]->GetPDGMass();
180 //sumofdaughtermass += daughtermass[index];
181 }
182
183 G4double EMASS = daughtermass[0];
184
185 // create parent G4DynamicParticle at rest
186 G4ThreeVector dummy;
187 G4DynamicParticle* parentparticle
188 = new G4DynamicParticle( G4MT_parent, dummy, 0.0);
189 // create G4Decayproducts
190 G4DecayProducts *products = new G4DecayProducts(*parentparticle);
191 delete parentparticle;
192
193 G4double x, y;
194
195 const std::size_t MAX_LOOP=1000;
196
197 for (std::size_t loop_counter1=0; loop_counter1<MAX_LOOP; ++loop_counter1)
198 {
199 for (std::size_t loop_counter2=0; loop_counter2<MAX_LOOP; ++loop_counter2)
200 {
201 x = xl + G4UniformRand()*(xu-xl);
202 y = yl + G4UniformRand()*(yu-yl);
203 if (x+y > 1.) break;
204 }
205 G4double d2w = D2W(x,y);
206 if (d2w > G4UniformRand()*d2wmax) break;
207 }
208
209 // Calculate the angle between positron and photon (cosine)
210 //
211 G4double cthetaGE = (y*(x-2.)+2.*(1.-x+beta*beta)) /
212 (x*std::sqrt(y*y-4.*beta*beta));
213
214 G4double G = x * EMPI/2.;
215 G4double E = y * EMPI/2.;
216
217 if (E < EMASS) E = EMASS;
218
219 // calculate daughter momentum
220 G4double daughtermomentum[2];
221
222 daughtermomentum[0] = std::sqrt(E*E - EMASS*EMASS);
223
224 G4double cthetaE = 2.*G4UniformRand()-1.;
225 G4double sthetaE = std::sqrt(1.-cthetaE*cthetaE);
226
227 G4double phiE = twopi*G4UniformRand()*rad;
228 G4double cphiE = std::cos(phiE);
229 G4double sphiE = std::sin(phiE);
230
231 // Coordinates of the decay positron
232 //
233 G4double px = sthetaE*cphiE;
234 G4double py = sthetaE*sphiE;
235 G4double pz = cthetaE;
236
237 G4ThreeVector direction0(px,py,pz);
238
239 G4DynamicParticle * daughterparticle0
240 = new G4DynamicParticle(G4MT_daughters[0], daughtermomentum[0]*direction0);
241
242 products->PushProducts(daughterparticle0);
243
244 daughtermomentum[1] = G;
245
246 G4double sthetaGE = std::sqrt(1.-cthetaGE*cthetaGE);
247
248 G4double phiGE = twopi*G4UniformRand()*rad;
249 G4double cphiGE = std::cos(phiGE);
250 G4double sphiGE = std::sin(phiGE);
251
252 // Coordinates of the decay gamma with respect to the decay positron
253 //
254 px = sthetaGE*cphiGE;
255 py = sthetaGE*sphiGE;
256 pz = cthetaGE;
257
258 G4ThreeVector direction1(px,py,pz);
259
260 direction1.rotateUz(direction0);
261
262 G4DynamicParticle * daughterparticle1
263 = new G4DynamicParticle(G4MT_daughters[1], daughtermomentum[1]*direction1);
264
265 products->PushProducts(daughterparticle1);
266
267 // output message
268#ifdef G4VERBOSE
269 if (GetVerboseLevel()>1)
270 {
271 G4cout << "G4PionRadiativeDecayChannel::DecayIt() -";
272 G4cout << " create decay products in rest frame " << G4endl;
273 products->DumpInfo();
274 }
275#endif
276
277 return products;
278}
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
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
void DumpInfo() const
G4int PushProducts(G4DynamicParticle *aParticle)
G4PionRadiativeDecayChannel & operator=(const G4PionRadiativeDecayChannel &)
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)