source: trunk/source/processes/hadronic/models/de_excitation/fermi_breakup/src/G4FermiPhaseSpaceDecay.cc @ 1315

Last change on this file since 1315 was 1315, checked in by garnier, 14 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 7.3 KB
Line 
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// Hadronic Process: Nuclear De-excitations
28// by V. Lara
29
30
31#include "G4FermiPhaseSpaceDecay.hh"
32#include "G4HadronicException.hh"
33#include "Randomize.hh"
34
35#include <algorithm>
36#include <numeric>
37#include <functional>
38
39G4FermiPhaseSpaceDecay::G4FermiPhaseSpaceDecay()
40{
41}
42
43G4FermiPhaseSpaceDecay::~G4FermiPhaseSpaceDecay()
44{
45}
46
47std::vector<G4LorentzVector*> *
48G4FermiPhaseSpaceDecay::KopylovNBodyDecay(const G4double M, const std::vector<G4double>& m) const
49  // Calculates momentum for N fragments (Kopylov's method of sampling is used)
50{
51  G4int N = m.size();
52
53  std::vector<G4LorentzVector*>* P = new std::vector<G4LorentzVector*>;
54  P->insert(P->begin(), N, static_cast<G4LorentzVector*>(0));
55
56  G4double mtot = std::accumulate( m.begin(), m.end(), 0.0);
57  G4double mu = mtot;
58  G4double PFragMagCM = 0.0;
59  G4double Mass = M;
60  G4double T = M-mtot;
61  G4LorentzVector PFragCM(0.0,0.0,0.0,0.0);
62  G4LorentzVector PFragLab(0.0,0.0,0.0,0.0);
63  G4LorentzVector PRestCM(0.0,0.0,0.0,0.0);
64  G4LorentzVector PRestLab(0.0,0.0,0.0,Mass);
65
66  for (G4int k = N-1; k > 0; k--)
67    {
68      mu -= m[k];
69      if (k>1) T *= BetaKopylov(k);
70      else T = 0.0;
71     
72      G4double RestMass = mu + T;
73
74      PFragMagCM = PtwoBody(Mass,m[k],RestMass);
75      if (PFragMagCM < 0) 
76        {
77          throw G4HadronicException(__FILE__, __LINE__, "G4FermiPhaseSpaceDecay::KopylovNBodyDecay: Error sampling fragments momenta!!");
78        }
79     
80
81      // Create a unit vector with a random direction isotropically distributed
82      G4ParticleMomentum RandVector(IsotropicVector(PFragMagCM));
83
84      PFragCM.setVect(RandVector);
85      PFragCM.setE(std::sqrt(RandVector.mag2()+m[k]*m[k]));
86
87      PRestCM.setVect(-RandVector);
88      PRestCM.setE(std::sqrt(RandVector.mag2()+RestMass*RestMass));
89
90
91      G4ThreeVector BoostV = PRestLab.boostVector();
92
93      PFragLab = PFragCM;
94      PFragLab.boost(BoostV);
95      PRestLab = PRestCM;
96      PRestLab.boost(BoostV);
97
98      P->operator[](k) = new G4LorentzVector(PFragLab);
99     
100      Mass = RestMass;
101    }
102
103    P->operator[](0) = new G4LorentzVector(PRestLab);
104
105  return P;
106
107}
108
109
110
111std::vector<G4LorentzVector*> *
112G4FermiPhaseSpaceDecay::NBodyDecay(const G4double M,  const std::vector<G4double>& m) const
113{
114  // Number of fragments
115  G4int N = m.size();
116  G4int i, j;
117  // Total Daughters Mass
118  G4double mtot = std::accumulate( m.begin(), m.end(), 0.0);
119  G4double Emax = M - mtot + m[0];
120  G4double Emin = 0.0;
121  G4double Wmax = 1.0;
122  for (i = 1; i < N; i++)
123    {
124      Emax += m[i];
125      Emin += m[i-1];
126      Wmax *= this->PtwoBody(Emax, Emin, m[i]);
127    }
128
129  G4int ntries = 0;
130  G4double weight = 1.0;
131  std::vector<G4double> p(N);
132  do
133    {
134      // Sample uniform random numbers in increasing order
135      std::vector<G4double> r;
136      r.reserve(N);
137      r.push_back(0.0);
138      for (i = 1; i < N-1; i++) r.push_back(G4UniformRand());
139      r.push_back(1.0);
140      std::sort(r.begin(),r.end(), std::less<G4double>());
141     
142      // Calculate virtual masses
143      std::vector<G4double> vm(N);
144      vm[0] = 0.0;
145      std::partial_sum(m.begin(), m.end(), vm.begin());
146      std::transform(r.begin(), r.end(), r.begin(), std::bind2nd(std::multiplies<G4double>(), M-mtot));
147      std::transform(r.begin(), r.end(), vm.begin(), vm.begin(), std::plus<G4double>());
148      r.clear();
149
150      // Calcualte daughter momenta
151      weight = 1.0;
152      for (j = 0; j < N-1; j++)
153        {
154          p[j] = PtwoBody(vm[j+1],vm[j],m[j+1]); 
155          if (p[j] < 0.0)
156            {
157              G4cerr << "G4FermiPhaseSpaceDecay::Decay: Daughter momentum less than zero\n";
158              weight = 0.0;
159              break;
160            }
161          else
162            {
163              weight *= p[j]; 
164            }
165        }
166      p[N-1] = PtwoBody(vm[N-2], m[N-2], m[N-1]);
167
168           
169      if (ntries++ > 1000000)
170        {
171          throw G4HadronicException(__FILE__, __LINE__, "G4FermiPhaseSpaceDecay::Decay: Cannot determine decay kinematics");
172        }
173    }
174  while ( weight < G4UniformRand()*Wmax );
175
176  std::vector<G4LorentzVector*> * P = new std::vector<G4LorentzVector*>; 
177  P->insert(P->begin(),N, static_cast<G4LorentzVector*>(0));
178
179  G4ParticleMomentum a3P = this->IsotropicVector(p[0]);
180
181  P->operator[](0) = new G4LorentzVector( a3P, std::sqrt(a3P.mag2()+m[0]*m[0]) ); 
182  P->operator[](1) = new G4LorentzVector(-a3P, std::sqrt(a3P.mag2()+m[1]*m[1]) );
183  for (i = 2; i < N; i++)
184    {
185      a3P = this->IsotropicVector(p[i-1]);
186      P->operator[](i) = new G4LorentzVector(a3P, std::sqrt(a3P.mag2() + m[i]*m[i]));
187      G4ThreeVector Beta = (-1.0)*P->operator[](i)->boostVector();
188      // boost already created particles
189      for (j = 0; j < i; j++)
190        {
191          P->operator[](j)->boost(Beta);
192        }
193    }
194 
195  return P;
196}
197
198
199
200
201std::vector<G4LorentzVector*> *
202G4FermiPhaseSpaceDecay::
203TwoBodyDecay(const G4double M, const std::vector<G4double>& m) const
204{
205  G4double m0 = m.front();
206  G4double m1 = m.back();
207  G4double psqr = this->PtwoBody(M,m0,m1);
208  G4ParticleMomentum p = this->IsotropicVector(std::sqrt(psqr));
209
210  G4LorentzVector * P41 = new G4LorentzVector;
211  P41->setVect(p);
212  P41->setE(std::sqrt(psqr+m0*m0));
213
214  G4LorentzVector * P42 = new G4LorentzVector;
215  P42->setVect(-p);
216  P42->setE(std::sqrt(psqr+m1*m1));
217
218  std::vector<G4LorentzVector*> * result = new std::vector<G4LorentzVector*>;
219  result->push_back(P41);
220  result->push_back(P42);
221  return result;
222}
223
224
225
226
227G4ParticleMomentum G4FermiPhaseSpaceDecay::IsotropicVector(const G4double Magnitude) const
228  // Samples a isotropic random vectorwith a magnitud given by Magnitude.
229  // By default Magnitude = 1.0
230{
231  G4double CosTheta = 1.0 - 2.0*G4UniformRand();
232  G4double SinTheta = std::sqrt(1.0 - CosTheta*CosTheta);
233  G4double Phi = twopi*G4UniformRand();
234  G4ParticleMomentum Vector(Magnitude*std::cos(Phi)*SinTheta,
235                            Magnitude*std::sin(Phi)*SinTheta,
236                            Magnitude*CosTheta);
237  return Vector;
238}
239
240
Note: See TracBrowser for help on using the repository browser.