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

Last change on this file since 962 was 819, checked in by garnier, 16 years ago

import all except CVS

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