Geant4  10.00.p03
G4KL3DecayChannel.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: G4KL3DecayChannel.cc 67971 2013-03-13 10:13:24Z gcosmo $
28 //
29 //
30 // ------------------------------------------------------------
31 // GEANT 4 class header file
32 //
33 // History: first implementation, based on object model of
34 // 30 May 1997 H.Kurashige
35 // ------------------------------------------------------------
36 
37 #include "G4ParticleDefinition.hh"
38 #include "G4PhysicalConstants.hh"
39 #include "G4SystemOfUnits.hh"
40 #include "G4DecayProducts.hh"
41 #include "G4VDecayChannel.hh"
42 #include "G4KL3DecayChannel.hh"
43 #include "Randomize.hh"
44 #include "G4LorentzVector.hh"
45 #include "G4LorentzRotation.hh"
46 
48  :G4VDecayChannel(),
49  massK(0.0), pLambda(0.0), pXi0(0.0)
50 {
51  daughterM[idPi] = 0.0;
52  daughterM[idLepton] = 0.0;
53  daughterM[idNutrino] = 0.0;
54 }
55 
56 
58  const G4String& theParentName,
59  G4double theBR,
60  const G4String& thePionName,
61  const G4String& theLeptonName,
62  const G4String& theNutrinoName)
63  :G4VDecayChannel("KL3 Decay",theParentName,
64  theBR, 3,
65  thePionName,theLeptonName,theNutrinoName)
66 {
67  static const G4String K_plus("kaon+");
68  static const G4String K_minus("kaon-");
69  static const G4String K_L("kaon0L");
70  static const G4String Mu_plus("mu+");
71  static const G4String Mu_minus("mu-");
72  static const G4String E_plus("e+");
73  static const G4String E_minus("e-");
74 
75  massK = 0.0;
76  daughterM[idPi] = 0.0;
77  daughterM[idLepton] = 0.0;
78  daughterM[idNutrino] = 0.0;
79 
80  // check modes
81  if ( ((theParentName == K_plus)&&(theLeptonName == E_plus)) ||
82  ((theParentName == K_minus)&&(theLeptonName == E_minus)) ) {
83  // K+- (Ke3)
84  pLambda = 0.0286;
85  pXi0 = -0.35;
86  } else if ( ((theParentName == K_plus)&&(theLeptonName == Mu_plus)) ||
87  ((theParentName == K_minus)&&(theLeptonName == Mu_minus)) ) {
88  // K+- (Kmu3)
89  pLambda = 0.033;
90  pXi0 = -0.35;
91  } else if ( (theParentName == K_L) &&
92  ((theLeptonName == E_plus) ||(theLeptonName == E_minus)) ){
93  // K0L (Ke3)
94  pLambda = 0.0300;
95  pXi0 = -0.11;
96  } else if ( (theParentName == K_L) &&
97  ((theLeptonName == Mu_plus) ||(theLeptonName == Mu_minus)) ){
98  // K0L (Kmu3)
99  pLambda = 0.034;
100  pXi0 = -0.11;
101  } else {
102 #ifdef G4VERBOSE
103  if (GetVerboseLevel()>2) {
104  G4cout << "G4KL3DecayChannel:: constructor :";
105  G4cout << "illegal arguments " << G4endl;;
106  DumpInfo();
107  }
108 #endif
109  // set values for K0L (Ke3) temporarily
110  pLambda = 0.0300;
111  pXi0 = -0.11;
112  }
113 }
114 
116 {
117 }
118 
120  G4VDecayChannel(right),
121  massK(right.massK),
122  pLambda(right.pLambda),
123  pXi0(right.pXi0)
124 {
125  daughterM[idPi] = right.daughterM[idPi];
128 }
129 
131 {
132  if (this != &right) {
134  verboseLevel = right.verboseLevel;
135  rbranch = right.rbranch;
136 
137  // copy parent name
138  parent_name = new G4String(*right.parent_name);
139 
140  // clear daughters_name array
142 
143  // recreate array
145  if ( numberOfDaughters >0 ) {
148  //copy daughters name
149  for (G4int index=0; index < numberOfDaughters; index++) {
150  daughters_name[index] = new G4String(*right.daughters_name[index]);
151  }
152  }
153  massK = right.massK;
154  pLambda = right.pLambda;
155  pXi0 = right.pXi0;
156  daughterM[idPi] = right.daughterM[idPi];
159  }
160  return *this;
161 }
162 
163 
165 {
166  // this version neglects muon polarization
167  // assumes the pure V-A coupling
168  // gives incorrect energy spectrum for Nutrinos
169 #ifdef G4VERBOSE
170  if (GetVerboseLevel()>1) G4cout << "G4KL3DecayChannel::DecayIt " << G4endl;
171 #endif
172 
173  // fill parent particle and its mass
174  if (G4MT_parent == 0) {
175  FillParent();
176  }
178 
179  // fill daughter particles and their mass
180  if (G4MT_daughters == 0) {
181  FillDaughters();
182  }
186 
187  // determine momentum/energy of daughters
188  // according to DalitzDensity
189  G4double daughterP[3], daughterE[3];
190  G4double w;
191  G4double r;
192  do {
193  r = G4UniformRand();
194  PhaseSpace(massK, &daughterM[0], &daughterE[0], &daughterP[0]);
195  w = DalitzDensity(daughterE[idPi],daughterE[idLepton],daughterE[idNutrino]);
196  } while ( r > w);
197 
198  // output message
199 #ifdef G4VERBOSE
200  if (GetVerboseLevel()>1) {
201  G4cout << *daughters_name[0] << ":" << daughterP[0]/GeV << "[GeV/c]" <<G4endl;
202  G4cout << *daughters_name[1] << ":" << daughterP[1]/GeV << "[GeV/c]" <<G4endl;
203  G4cout << *daughters_name[2] << ":" << daughterP[2]/GeV << "[GeV/c]" <<G4endl;
204  }
205 #endif
206  //create parent G4DynamicParticle at rest
207  G4ThreeVector* direction = new G4ThreeVector(1.0,0.0,0.0);
208  G4DynamicParticle * parentparticle = new G4DynamicParticle( G4MT_parent, *direction, 0.0);
209  delete direction;
210 
211  //create G4Decayproducts
212  G4DecayProducts *products = new G4DecayProducts(*parentparticle);
213  delete parentparticle;
214 
215  //create daughter G4DynamicParticle
216  G4double costheta, sintheta, phi, sinphi, cosphi;
217  G4double costhetan, sinthetan, phin, sinphin, cosphin;
218 
219  // pion
220  costheta = 2.*G4UniformRand()-1.0;
221  sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
222  phi = twopi*G4UniformRand()*rad;
223  sinphi = std::sin(phi);
224  cosphi = std::cos(phi);
225  direction = new G4ThreeVector(sintheta*cosphi,sintheta*sinphi,costheta);
226  G4ThreeVector momentum0 = (*direction)*daughterP[0];
227  G4DynamicParticle * daughterparticle
228  = new G4DynamicParticle( G4MT_daughters[0], momentum0);
229  products->PushProducts(daughterparticle);
230 
231  // neutrino
232  costhetan = (daughterP[1]*daughterP[1]-daughterP[2]*daughterP[2]-daughterP[0]*daughterP[0])/(2.0*daughterP[2]*daughterP[0]);
233  sinthetan = std::sqrt((1.0-costhetan)*(1.0+costhetan));
234  phin = twopi*G4UniformRand()*rad;
235  sinphin = std::sin(phin);
236  cosphin = std::cos(phin);
237  direction->setX( sinthetan*cosphin*costheta*cosphi - sinthetan*sinphin*sinphi + costhetan*sintheta*cosphi);
238  direction->setY( sinthetan*cosphin*costheta*sinphi + sinthetan*sinphin*cosphi + costhetan*sintheta*sinphi);
239  direction->setZ( -sinthetan*cosphin*sintheta + costhetan*costheta);
240 
241  G4ThreeVector momentum2 = (*direction)*daughterP[2];
242  daughterparticle = new G4DynamicParticle( G4MT_daughters[2], momentum2);
243  products->PushProducts(daughterparticle);
244 
245  //lepton
246  G4ThreeVector momentum1 = (momentum0 + momentum2) * (-1.0);
247  daughterparticle =
248  new G4DynamicParticle( G4MT_daughters[1], momentum1);
249  products->PushProducts(daughterparticle);
250 
251 #ifdef G4VERBOSE
252  if (GetVerboseLevel()>1) {
253  G4cout << "G4KL3DecayChannel::DecayIt ";
254  G4cout << " create decay products in rest frame " <<G4endl;
255  G4cout << " decay products address=" << products << G4endl;
256  products->DumpInfo();
257  }
258 #endif
259  delete direction;
260  return products;
261 }
262 
264  const G4double* M,
265  G4double* E,
266  G4double* P )
267 // algorism of this code is originally written in GDECA3 of GEANT3
268 {
269 
270  //sum of daughters'mass
271  G4double sumofdaughtermass = 0.0;
272  G4int index;
273  for (index=0; index<3; index++){
274  sumofdaughtermass += M[index];
275  }
276 
277  //calculate daughter momentum
278  // Generate two
279  G4double rd1, rd2, rd;
280  G4double momentummax=0.0, momentumsum = 0.0;
282 
283  do {
284  rd1 = G4UniformRand();
285  rd2 = G4UniformRand();
286  if (rd2 > rd1) {
287  rd = rd1;
288  rd1 = rd2;
289  rd2 = rd;
290  }
291  momentummax = 0.0;
292  momentumsum = 0.0;
293  // daughter 0
294  energy = rd2*(parentM - sumofdaughtermass);
295  P[0] = std::sqrt(energy*energy + 2.0*energy*M[0]);
296  E[0] = energy;
297  if ( P[0] >momentummax )momentummax = P[0];
298  momentumsum += P[0];
299  // daughter 1
300  energy = (1.-rd1)*(parentM - sumofdaughtermass);
301  P[1] = std::sqrt(energy*energy + 2.0*energy*M[1]);
302  E[1] = energy;
303  if ( P[1] >momentummax )momentummax = P[1];
304  momentumsum += P[1];
305  // daughter 2
306  energy = (rd1-rd2)*(parentM - sumofdaughtermass);
307  P[2] = std::sqrt(energy*energy + 2.0*energy*M[2]);
308  E[2] = energy;
309  if ( P[2] >momentummax )momentummax = P[2];
310  momentumsum += P[2];
311  } while (momentummax > momentumsum - momentummax );
312 
313 #ifdef G4VERBOSE
314  if (GetVerboseLevel()>2) {
315  G4cout << "G4KL3DecayChannel::PhaseSpace ";
316  G4cout << "Kon mass:" << parentM/GeV << "GeV/c/c" << G4endl;
317  for (index=0; index<3; index++){
318  G4cout << index << " : " << M[index]/GeV << "GeV/c/c ";
319  G4cout << " : " << E[index]/GeV << "GeV ";
320  G4cout << " : " << P[index]/GeV << "GeV/c " << G4endl;
321  }
322  }
323 #endif
324 }
325 
326 
328 {
329  // KL3 decay Dalitz Plot Density
330  // see Chounet et al Phys. Rep. 4, 201
331  // arguments
332  // Epi: kinetic enregy of pion
333  // El: kinetic enregy of lepton (e or mu)
334  // Enu: kinetic energy of nutrino
335  // constants
336  // pLambda : linear energy dependence of f+
337  // pXi0 : = f+(0)/f-
338  // pNorm : normalization factor
339  // variables
340  // Epi: total energy of pion
341  // El: total energy of lepton (e or mu)
342  // Enu: total energy of nutrino
343 
344  // mass of daughters
345  G4double massPi = daughterM[idPi];
346  G4double massL = daughterM[idLepton];
347  G4double massNu = daughterM[idNutrino];
348 
349  // calcurate total energy
350  Epi = Epi + massPi;
351  El = El + massL;
352  Enu = Enu + massNu;
353 
354  G4double Epi_max = (massK*massK+massPi*massPi-massL*massL)/2.0/massK;
355  G4double E = Epi_max - Epi;
356  G4double q2 = massK*massK + massPi*massPi - 2.0*massK*Epi;
357 
358  G4double F = 1.0 + pLambda*q2/massPi/massPi;
359  G4double Fmax = 1.0;
360  if (pLambda >0.0) Fmax = (1.0 + pLambda*(massK*massK/massPi/massPi+1.0));
361 
362  G4double Xi = pXi0*(1.0 + pLambda*q2/massPi/massPi);
363 
364  G4double coeffA = massK*(2.0*El*Enu-massK*E)+massL*massL*(E/4.0-Enu);
365  G4double coeffB = massL*massL*(Enu-E/2.0);
366  G4double coeffC = massL*massL*E/4.0;
367 
368  G4double RhoMax = (Fmax*Fmax)*(massK*massK*massK/8.0);
369 
370  G4double Rho = (F*F)*(coeffA + coeffB*Xi + coeffC*Xi*Xi);
371 
372 #ifdef G4VERBOSE
373  if (GetVerboseLevel()>2) {
374  G4cout << "G4KL3DecayChannel::DalitzDensity " <<G4endl;
375  G4cout << " Pi[" << massPi/GeV <<"GeV/c/c] :" << Epi/GeV << "GeV" <<G4endl;
376  G4cout << " L[" << massL/GeV <<"GeV/c/c] :" << El/GeV << "GeV" <<G4endl;
377  G4cout << " Nu[" << massNu/GeV <<"GeV/c/c] :" << Enu/GeV << "GeV" <<G4endl;
378  G4cout << " F :" << F << " Fmax :" << Fmax << " Xi :" << Xi << G4endl;
379  G4cout << " A :" << coeffA << " B :" << coeffB << " C :"<< coeffC <<G4endl;
380  G4cout << " Rho :" << Rho << " RhoMax :" << RhoMax << G4endl;
381  }
382 #endif
383  return (Rho/RhoMax);
384 }
385 
386 
virtual G4DecayProducts * DecayIt(G4double)
CLHEP::Hep3Vector G4ThreeVector
G4int PushProducts(G4DynamicParticle *aParticle)
G4ParticleDefinition * G4MT_parent
G4ParticleDefinition ** G4MT_daughters
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
G4double DalitzDensity(G4double Epi, G4double El, G4double Enu)
G4String kinematics_name
void DumpInfo() const
G4String * parent_name
static const double GeV
Definition: G4SIunits.hh:196
static const double rad
Definition: G4SIunits.hh:130
G4int GetVerboseLevel() const
G4double GetPDGMass() const
G4double energy(const ThreeVector &p, const G4double m)
G4KL3DecayChannel & operator=(const G4KL3DecayChannel &)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4String ** daughters_name
void PhaseSpace(G4double Mparent, const G4double *Mdaughter, G4double *Edaughter, G4double *Pdaughter)