Geant4  10.01
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 <vector>
43 #include <CLHEP/Units/PhysicalConstants.h>
44 
45 #include "G4LorentzVector.hh"
46 #include "G4ThreeVector.hh"
47 #include "Randomize.hh"
48 #include "G4Pow.hh"
49 
51 {
52 public:
53 
56 
57  inline std::vector<G4LorentzVector*> *
58  Decay(const G4double, const std::vector<G4double>&) const;
59 
60 private:
61 
62  inline G4double PtwoBody(G4double E, G4double P1, G4double P2) const;
63 
64  inline G4ThreeVector IsotropicVector(const G4double Magnitude = 1.0) const;
65 
66  inline G4double BetaKopylov(G4int) const;
67 
68  std::vector<G4LorentzVector*> *
69  TwoBodyDecay(G4double, const std::vector<G4double>&) const;
70 
71  std::vector<G4LorentzVector*> *
72  NBodyDecay(G4double, const std::vector<G4double>&) const;
73 
74  std::vector<G4LorentzVector*> *
75  KopylovNBodyDecay(G4double, const std::vector<G4double>&) const;
76 
77  void DumpProblem(G4double E, G4double P1, G4double P2, G4double P) const;
78 
83 
85 };
86 
87 inline G4double
89 {
90  G4double res = 0.0;
91  G4double P = (E+P1+P2)*(E+P1-P2)*(E-P1+P2)*(E-P1-P2)/(4.0*E*E);
92  if (P>0.0) { res = std::sqrt(P); }
93  else { DumpProblem(E,P1,P2,P); }
94  return res;
95 }
96 
97 inline std::vector<G4LorentzVector*> * G4FermiPhaseSpaceDecay::
98 Decay(G4double parent_mass_parameter, const std::vector<G4double>& fragment_masses) const
99 {
100  return KopylovNBodyDecay(parent_mass_parameter,fragment_masses);
101 }
102 
104 {
105  G4int N = 3*K - 5;
106  G4double xN = G4double(N);
107  G4double F;
108  //G4double Fmax=std::pow((3.*K-5.)/(3.*K-4.),(3.*K-5.)/2.)*std::sqrt(1-((3.*K-5.)/(3.*K-4.)));
109  // VI variant
110  G4double Fmax = std::sqrt(g4pow->powN(xN/(xN + 1),N)/(xN + 1));
111  G4double chi;
112  do {
113  chi = G4UniformRand();
114  F = std::sqrt(g4pow->powN(chi,N)*(1-chi));
115  } while ( Fmax*G4UniformRand() > F);
116  return chi;
117 }
118 
119 inline G4ThreeVector
121  // Samples a isotropic random vectorwith a magnitud given by Magnitude.
122  // By default Magnitude = 1.0
123 {
124  G4double CosTheta = 2.0*G4UniformRand() - 1.0;
125  G4double SinTheta = std::sqrt((1. - CosTheta)*(1. + CosTheta));
126  G4double Phi = CLHEP::twopi*G4UniformRand();
127  G4ThreeVector Vector(Magnitude*std::cos(Phi)*SinTheta,
128  Magnitude*std::sin(Phi)*SinTheta,
129  Magnitude*CosTheta);
130  return Vector;
131 }
132 
133 #endif
const G4FermiPhaseSpaceDecay & operator=(const G4FermiPhaseSpaceDecay &)
std::vector< G4LorentzVector * > * NBodyDecay(G4double, const std::vector< G4double > &) const
G4double powN(G4double x, G4int n) const
Definition: G4Pow.cc:128
static const G4double * P1[nN]
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Pow.hh:56
std::vector< G4LorentzVector * > * TwoBodyDecay(G4double, const std::vector< G4double > &) const
G4bool operator!=(const G4FermiPhaseSpaceDecay &)
int G4int
Definition: G4Types.hh:78
std::vector< G4LorentzVector * > * Decay(const G4double, const std::vector< G4double > &) const
G4double PtwoBody(G4double E, G4double P1, G4double P2) const
void DumpProblem(G4double E, G4double P1, G4double P2, G4double P) const
G4double BetaKopylov(G4int) const
#define G4UniformRand()
Definition: Randomize.hh:95
G4bool operator==(const G4FermiPhaseSpaceDecay &)
bool G4bool
Definition: G4Types.hh:79
G4ThreeVector IsotropicVector(const G4double Magnitude=1.0) const
static const G4double * P2[nN]
double G4double
Definition: G4Types.hh:76
std::vector< G4LorentzVector * > * KopylovNBodyDecay(G4double, const std::vector< G4double > &) const