Geant4  10.02.p01
G4FermiPhaseSpaceDecay.hh
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 // $Id: G4ExcitationHandler.hh,v 1.13 2010-11-17 16:20:31 vnivanch Exp $
27 //
28 // Hadronic Process: Phase space decay for the Fermi BreakUp model
29 // by V. Lara
30 //
31 // Modifications:
32 // 01.04.2011 General cleanup by V.Ivanchenko:
33 // - IsotropicVector is inlined
34 // - Momentum computation return zero or positive value
35 // - DumpProblem method is added providing more information
36 // - Reduced usage of exotic std functions
37 //
38 
39 #ifndef G4FermiPhaseSpaceDecay_hh
40 #define G4FermiPhaseSpaceDecay_hh 1
41 
42 #include "G4LorentzVector.hh"
43 #include "G4ThreeVector.hh"
44 #include "Randomize.hh"
45 #include "G4Pow.hh"
46 
47 #include <vector>
48 #include <CLHEP/Units/PhysicalConstants.h>
49 #include <CLHEP/Random/RandomEngine.h>
50 
52 {
53 public:
54 
57 
58  inline std::vector<G4LorentzVector*>*
59  Decay(G4double parent_mass,
60  const std::vector<G4double>& fragment_masses) const;
61 
62 private:
63 
64  inline G4double PtwoBody(G4double E, G4double P1, G4double P2) const;
65 
66  inline G4ThreeVector IsotropicVector(G4double Magnitude,
67  CLHEP::HepRandomEngine*) const;
68 
69  inline G4double BetaKopylov(G4int, CLHEP::HepRandomEngine*) const;
70 
71  std::vector<G4LorentzVector*> *
72  KopylovNBodyDecay(G4double, const std::vector<G4double>&) const;
73 
74  void DumpProblem(G4double E, G4double P1, G4double P2, G4double P) const;
75 
80 
82 };
83 
84 inline G4double
86 {
87  G4double res = 0.0;
88  G4double P = (E+P1+P2)*(E+P1-P2)*(E-P1+P2)*(E-P1-P2)/(4.0*E*E);
89  if (P>0.0) { res = std::sqrt(P); }
90  else { DumpProblem(E,P1,P2,P); }
91  return res;
92 }
93 
94 inline std::vector<G4LorentzVector*>*
96  const std::vector<G4double>& fragment_masses) const
97 {
98  return KopylovNBodyDecay(parent_mass, fragment_masses);
99 }
100 
101 inline G4double
103  CLHEP::HepRandomEngine* rndmEngine) const
104 {
105  G4int N = 3*K - 5;
106  G4double xN = G4double(N);
107  G4double F;
108  // VI variant
109  G4double Fmax = std::sqrt(g4pow->powN(xN/(xN + 1),N)/(xN + 1));
110  G4double chi;
111  do {
112  chi = rndmEngine->flat();
113  F = std::sqrt(g4pow->powN(chi,N)*(1-chi));
114  // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
115  } while ( Fmax*rndmEngine->flat() > F);
116  return chi;
117 }
118 
119 inline G4ThreeVector
121  CLHEP::HepRandomEngine* rndmEngine) const
122  // Samples a isotropic random vectorwith a magnitud given by Magnitude.
123  // By default Magnitude = 1.0
124 {
125  G4double CosTheta = 2.0*rndmEngine->flat() - 1.0;
126  G4double SinTheta = std::sqrt((1. - CosTheta)*(1. + CosTheta));
127  G4double Phi = CLHEP::twopi*rndmEngine->flat();
128  G4ThreeVector Vector(Magnitude*std::cos(Phi)*SinTheta,
129  Magnitude*std::sin(Phi)*SinTheta,
130  Magnitude*CosTheta);
131  return Vector;
132 }
133 
134 #endif
const G4FermiPhaseSpaceDecay & operator=(const G4FermiPhaseSpaceDecay &)
G4double powN(G4double x, G4int n) const
Definition: G4Pow.cc:128
static const G4double * P1[nN]
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Pow.hh:56
G4bool operator!=(const G4FermiPhaseSpaceDecay &)
int G4int
Definition: G4Types.hh:78
static double P[]
G4double PtwoBody(G4double E, G4double P1, G4double P2) const
void DumpProblem(G4double E, G4double P1, G4double P2, G4double P) const
G4bool operator==(const G4FermiPhaseSpaceDecay &)
bool G4bool
Definition: G4Types.hh:79
static const double twopi
Definition: G4SIunits.hh:75
G4double BetaKopylov(G4int, CLHEP::HepRandomEngine *) const
std::vector< G4LorentzVector * > * Decay(G4double parent_mass, const std::vector< G4double > &fragment_masses) const
G4ThreeVector IsotropicVector(G4double Magnitude, CLHEP::HepRandomEngine *) const
static const G4double * P2[nN]
double G4double
Definition: G4Types.hh:76
std::vector< G4LorentzVector * > * KopylovNBodyDecay(G4double, const std::vector< G4double > &) const