Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FermiPhaseSpaceDecay.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// $Id: G4ExcitationHandler.hh,v 1.13 2010-11-17 16:20:31 vnivanch Exp $
27// GEANT4 tag $Name: not supported by cvs2svn $
28//
29// Hadronic Process: Phase space decay for the Fermi BreakUp model
30// by V. Lara
31//
32// Modifications:
33// 01.04.2011 General cleanup by V.Ivanchenko:
34// - IsotropicVector is inlined
35// - Momentum computation return zero or positive value
36// - DumpProblem method is added providing more information
37// - Reduced usage of exotic std functions
38
39#include <numeric>
40
42#include "G4SystemOfUnits.hh"
44
46{
47 g4pow = G4Pow::GetInstance();
48}
49
51{
52}
53
54std::vector<G4LorentzVector*> *
55G4FermiPhaseSpaceDecay::KopylovNBodyDecay(const G4double M,
56 const std::vector<G4double>& mr) const
57 // Calculates momentum for N fragments (Kopylov's method of sampling is used)
58{
59 size_t N = mr.size();
60
61 std::vector<G4LorentzVector*>* P = new std::vector<G4LorentzVector*>(N, 0);
62
63 G4double mtot = std::accumulate( mr.begin(), mr.end(), 0.0);
64 G4double mu = mtot;
65 G4double PFragMagCM = 0.0;
66 G4double Mass = M;
67 G4double T = Mass-mtot;
68 G4LorentzVector PFragCM(0.0,0.0,0.0,0.0);
69 G4LorentzVector PRestCM(0.0,0.0,0.0,0.0);
70 G4LorentzVector PRestLab(0.0,0.0,0.0,Mass);
71
72 for (size_t k = N-1; k>0; --k)
73 {
74 mu -= mr[k];
75 if (k>1) { T *= BetaKopylov(k); }
76 else { T = 0.0; }
77
78 G4double RestMass = mu + T;
79
80 PFragMagCM = PtwoBody(Mass,mr[k],RestMass);
81
82 // Create a unit vector with a random direction isotropically distributed
83 G4ThreeVector RandVector(IsotropicVector(PFragMagCM));
84
85 PFragCM.setVect(RandVector);
86 PFragCM.setE(std::sqrt(PFragMagCM*PFragMagCM + mr[k]*mr[k]));
87
88 PRestCM.setVect(-RandVector);
89 PRestCM.setE(std::sqrt(PFragMagCM*PFragMagCM + RestMass*RestMass));
90
91
92 G4ThreeVector BoostV = PRestLab.boostVector();
93
94 PFragCM.boost(BoostV);
95 PRestCM.boost(BoostV);
96 PRestLab = PRestCM;
97
98 (*P)[k] = new G4LorentzVector(PFragCM);
99
100 Mass = RestMass;
101 }
102
103 (*P)[0] = new G4LorentzVector(PRestLab);
104
105 return P;
106}
107
108
109std::vector<G4LorentzVector*> *
110G4FermiPhaseSpaceDecay::NBodyDecay(G4double M, const std::vector<G4double>& mr) const
111{
112 // Number of fragments
113 size_t N = mr.size();
114 size_t i, j;
115 // Total Daughters Mass
116 G4double mtot = std::accumulate( mr.begin(), mr.end(), 0.0);
117 G4double Emax = M - mtot + mr[0];
118 G4double Emin = 0.0;
119 G4double Wmax = 1.0;
120 for (i = 1; i < N; i++)
121 {
122 Emax += mr[i];
123 Emin += mr[i-1];
124 Wmax *= PtwoBody(Emax, Emin, mr[i]);
125 }
126
127 G4int ntries = 0;
128 G4double weight = 1.0;
129 std::vector<G4double> p(N, 0.0);
130 std::vector<G4double> r(N,0.0);
131 std::vector<G4double> vm(N, 0.0);
132 r[N-1] = 1.0;
133
134 do
135 {
136 // Sample uniform random numbers in increasing order
137 for (i = 1; i < N-1; i++) { r[i] = G4UniformRand(); }
138 std::sort(r.begin(),r.end(), std::less<G4double>());
139
140 // Calculate virtual masses
141 std::partial_sum(mr.begin(), mr.end(), vm.begin());
142 std::transform(r.begin(), r.end(), r.begin(),
143 std::bind2nd(std::multiplies<G4double>(), M-mtot));
144 std::transform(r.begin(), r.end(), vm.begin(), vm.begin(), std::plus<G4double>());
145
146 // Calcualte daughter momenta
147 weight = 1.0;
148 for (j = 0; j < N-1; j++)
149 {
150 p[j] = PtwoBody(vm[j+1],vm[j],mr[j+1]);
151 weight *= p[j];
152 }
153 p[N-1] = PtwoBody(vm[N-2],mr[N-2],mr[N-1]);
154
155
156 if (ntries++ > 1000000)
157 {
158 throw G4HadronicException(__FILE__, __LINE__, "Failed to decay");
159 }
160 }
161 while ( weight < G4UniformRand()*Wmax );
162
163 std::vector<G4LorentzVector*> * P = new std::vector<G4LorentzVector*>(N, 0);
164
165 G4ThreeVector a3P = IsotropicVector(p[0]);
166
167 (*P)[0] = new G4LorentzVector( a3P, std::sqrt(a3P.mag2()+mr[0]*mr[0]) );
168 (*P)[1] = new G4LorentzVector(-a3P, std::sqrt(a3P.mag2()+mr[1]*mr[1]) );
169 for (i = 2; i < N; i++)
170 {
171 a3P = IsotropicVector(p[i-1]);
172 (*P)[i] = new G4LorentzVector(a3P, std::sqrt(a3P.mag2() + mr[i]*mr[i]));
173 G4ThreeVector Beta = -((*P)[i]->boostVector());
174 // boost already created particles
175 for (j = 0; j < i; j++)
176 {
177 (*P)[j]->boost(Beta);
178 }
179 }
180
181 return P;
182}
183
184std::vector<G4LorentzVector*> *
185G4FermiPhaseSpaceDecay::TwoBodyDecay(G4double M,
186 const std::vector<G4double>& mass) const
187{
188 G4double m0 = mass.front();
189 G4double m1 = mass.back();
190 G4double mom = PtwoBody(M,m0,m1);
191 G4ThreeVector p = IsotropicVector(mom);
192
194 P41->setVect(p);
195 P41->setE(std::sqrt(mom*mom + m0*m0));
196
198 P42->setVect(-p);
199 P42->setE(std::sqrt(mom*mom + m1*m1));
200
201 std::vector<G4LorentzVector*> * result = new std::vector<G4LorentzVector*>;
202 result->push_back(P41);
203 result->push_back(P42);
204 return result;
205}
206
207void
208G4FermiPhaseSpaceDecay::DumpProblem(G4double E, G4double P1, G4double P2,
209 G4double P) const
210{
211 G4cout << "G4FermiPhaseSpaceDecay: problem of decay of M(GeV)= " << E/GeV
212 << " on M1(GeV)= " << P1/GeV << " and M2(GeV)= " << P2/GeV
213 << " P(MeV)= " << P/MeV << " < 0" << G4endl;
214 // exception only if the problem is numerically significant
215 if(P < -CLHEP::eV) {
216 throw G4HadronicException(__FILE__, __LINE__,"Error in decay kinematics");
217 }
218}
219
220
CLHEP::HepLorentzVector G4LorentzVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
double mag2() const
void setVect(const Hep3Vector &)
static G4Pow * GetInstance()
Definition: G4Pow.cc:50