Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FermiConfigurationList.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: G4VFermiBreakUp.cc,v 1.5 2006-06-29 20:13:13 gunter Exp $
27// GEANT4 tag $Name: not supported by cvs2svn $
28//
29// Hadronic Process: Nuclear De-excitations
30// by V. Lara (Nov 1998)
31//
32// Modifications:
33// 01.04.2011 General cleanup by V.Ivanchenko - more clean usage of static
34// 23.04.2011 V.Ivanchenko: make this class to be responsible for
35// selection of decay channel and decay
36
37#include <set>
38
42#include "Randomize.hh"
43#include "G4Pow.hh"
44
45const G4double G4FermiConfigurationList::Kappa = 6.0;
46const G4double G4FermiConfigurationList::r0 = 1.3*CLHEP::fermi;
47
49{
51 g4pow = G4Pow::GetInstance();
52 Coef = 0.6*(CLHEP::elm_coupling/r0)/g4pow->Z13(1+G4int(Kappa));
53 ConstCoeff = g4pow->powN(r0/hbarc,3)*Kappa*std::sqrt(2.0/pi)/3.0;
54
55 // 16 is the max number
56 nmax = 50;
57 NormalizedWeights.resize(nmax,0.0);
58}
59
61{}
62
64G4FermiConfigurationList::CoulombBarrier(
65 const std::vector<const G4VFermiFragment*>& conf)
66{
67 // Calculates Coulomb Barrier (MeV) for given channel with K fragments.
68 G4int SumA = 0;
69 G4int SumZ = 0;
70 G4double CoulombEnergy = 0.;
71 size_t nn = conf.size();
72 for (size_t i=0; i<nn; ++i) {
73 G4int z = conf[i]->GetZ();
74 G4int a = conf[i]->GetA();
75 CoulombEnergy += G4double(z*z)/g4pow->Z13(a);
76 SumA += a;
77 SumZ += z;
78 }
79 CoulombEnergy -= SumZ*SumZ/g4pow->Z13(SumA);
80 return -Coef * CoulombEnergy;
81}
82
84G4FermiConfigurationList::DecayProbability(G4int A, G4double TotalE,
86 // Decay probability for a given channel with K fragments
87{
88 // A: Atomic Weight
89 // TotalE: Total energy of nucleus
90
91 G4double KineticEnergy = TotalE; // MeV
92 G4double ProdMass = 1.0;
93 G4double SumMass = 0.0;
94 G4double S_n = 1.0;
95 std::set<G4int> combSet;
96 std::multiset<G4int> combmSet;
97
98 const std::vector<const G4VFermiFragment*> flist =
99 conf->GetFragmentList();
100
101 // Number of fragments
102 G4int K = flist.size();
103
104 for (G4int i=0; i<K; ++i) {
105 G4int a = flist[i]->GetA();
106 combSet.insert(a);
107 combmSet.insert(a);
108 G4double mass = flist[i]->GetFragmentMass();
109 ProdMass *= mass;
110 SumMass += mass;
111 // Spin factor S_n
112 S_n *= flist[i]->GetPolarization();
113 KineticEnergy -= mass + flist[i]->GetExcitationEnergy();
114 }
115
116 // Check that there is enough energy to produce K fragments
117 KineticEnergy -= CoulombBarrier(flist);
118 if (KineticEnergy <= 0.0) { return 0.0; }
119
120 G4double MassFactor = ProdMass/SumMass;
121 MassFactor *= std::sqrt(MassFactor);
122
123 // This is the constant (doesn't depend on nucleus) part
124 G4double Coeff = g4pow->powN(ConstCoeff*A, K-1);
125
126 //JMQ 111009 Bug fixed: gamma function for odd K was wrong by a factor 2
127 // new implementation explicitely according to standard properties of Gamma function
128 // Calculation of 1/Gamma(3(k-1)/2)
129 G4double Gamma = 1.0;
130 // G4double arg = 3.0*(K-1)/2.0 - 1.0;
131 // while (arg > 1.1)
132 // {
133 // Gamma *= arg;
134 // arg--;
135 // }
136 // if ((K-1)%2 == 1) Gamma *= std::sqrt(pi);
137
138 if ((K-1)%2 != 1)
139
140 {
141 G4double arg = 3.0*(K-1)/2.0 - 1.0;
142 while (arg > 1.1)
143 {
144 Gamma *= arg;
145 arg--;
146 }
147 }
148 else {
149 G4double n= 3.0*K/2.0-2.0;
150 G4double arg2=2*n-1;
151 while (arg2>1.1)
152 {
153 Gamma *= arg2;
154 arg2 -= 2;
155 }
156 Gamma *= std::sqrt(pi)/g4pow->powZ(2,n);
157 }
158 // end of new implementation of Gamma function
159
160
161 // Permutation Factor G_n
162 G4double G_n = 1.0;
163 for (std::set<G4int>::iterator itr = combSet.begin();
164 itr != combSet.end(); ++itr)
165 {
166 for (G4int ni = combmSet.count(*itr); ni > 1; ni--) { G_n *= ni; }
167 }
168
169 G4double Weight = Coeff * MassFactor * (S_n / G_n) / Gamma;
170 Weight *= std::sqrt(g4pow->powN(KineticEnergy,3*(K-1)))/KineticEnergy;
171
172 return Weight;
173}
174
177{
178 // Calculate Momenta of K fragments
179 G4double M = theNucleus.GetMomentum().m();
180 const std::vector<const G4VFermiFragment*>* conf =
181 SelectConfiguration(theNucleus.GetZ_asInt(),
182 theNucleus.GetA_asInt(), M);
183
184
185 G4FragmentVector* theResult = new G4FragmentVector();
186 size_t nn = conf->size();
187 if(1 >= nn) {
188 theResult->push_back(new G4Fragment(theNucleus));
189 delete conf;
190 return theResult;
191 }
192
193 G4ThreeVector boostVector = theNucleus.GetMomentum().boostVector();
194 std::vector<G4double> mr;
195 mr.reserve(nn);
196 for(size_t i=0; i<nn; ++i) {
197 mr.push_back( (*conf)[i]->GetTotalEnergy() );
198 }
199 std::vector<G4LorentzVector*>* mom = thePhaseSpace.Decay(M,mr);
200 if(!mom) {
201 delete conf;
202 return theResult;
203 }
204
205 size_t nmom = mom->size();
206
207 // Go back to the Lab Frame
208 if(0 < nmom) {
209 for (size_t j=0; j<nmom; ++j) {
210 G4LorentzVector* FourMomentum = (*mom)[j];
211
212 // Lorentz boost
213 FourMomentum->boost(boostVector);
214
215 G4FragmentVector* fragment = (*conf)[j]->GetFragment(*FourMomentum);
216
217 size_t nfrag = fragment->size();
218 for (size_t k=0; k<nfrag; ++k) { theResult->push_back((*fragment)[k]); }
219 delete fragment;
220 delete (*mom)[j];
221 }
222 }
223
224 delete mom;
225 delete conf;
226 return theResult;
227}
228
229const std::vector<const G4VFermiFragment*>*
230G4FermiConfigurationList::SelectConfiguration(G4int Z, G4int A, G4double mass)
231{
232 std::vector<const G4VFermiFragment*>* res =
233 new std::vector<const G4VFermiFragment*>;
234 const std::vector<G4FermiConfiguration*>* conflist =
235 thePool->GetConfigurationList(Z, A, mass);
236 if(!conflist) { return res; }
237 size_t nn = conflist->size();
238 if(0 < nn) {
239 size_t idx = 0;
240 if(1 < nn) {
241 if(nn > nmax) {
242 nmax = nn;
243 NormalizedWeights.resize(nmax,0.0);
244 }
245 G4double prob = 0.0;
246 for(size_t i=0; i<nn; ++i) {
247 prob += DecayProbability(A, mass, (*conflist)[i]);
248 NormalizedWeights[i] = prob;
249 }
250 prob *= G4UniformRand();
251 for(idx=0; idx<nn; ++idx) {
252 if(NormalizedWeights[idx] >= prob) { break; }
253 }
254 }
255 const std::vector<const G4VFermiFragment*> flist =
256 (*conflist)[idx]->GetFragmentList();
257 size_t nf = flist.size();
258 for(size_t i=0; i<nf; ++i) { res->push_back(flist[i]); }
259 //G4cout << "FermiBreakUp: " << nn << " conf; selected idx= "
260 // << idx << " Nprod= " << nf << G4endl;
261 }
262 delete conflist;
263 return res;
264}
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
G4FragmentVector * GetFragments(const G4Fragment &theNucleus)
const std::vector< const G4VFermiFragment * > & GetFragmentList()
const std::vector< G4FermiConfiguration * > * GetConfigurationList(G4int Z, G4int A, G4double mass)
static G4FermiFragmentsPool * Instance()
std::vector< G4LorentzVector * > * Decay(const G4double, const std::vector< G4double > &) const
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:251
G4int GetZ_asInt() const
Definition: G4Fragment.hh:223
G4int GetA_asInt() const
Definition: G4Fragment.hh:218
static G4Pow * GetInstance()
Definition: G4Pow.cc:50
G4double powN(G4double x, G4int n)
Definition: G4Pow.cc:98
G4double Z13(G4int Z)
Definition: G4Pow.hh:110
G4double powZ(G4int Z, G4double y)
Definition: G4Pow.hh:180