[819] | 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 | // |
---|
[1315] | 27 | // $Id: G4FermiConfiguration.cc,v 1.13 2010/04/26 11:14:28 vnivanch Exp $ |
---|
| 28 | // GEANT4 tag $Name: geant4-09-04-beta-cand-01 $ |
---|
[819] | 29 | // |
---|
| 30 | // Hadronic Process: Nuclear De-excitations |
---|
| 31 | // by V. Lara (Nov 1998) |
---|
| 32 | // |
---|
[1315] | 33 | // J. M. Quesada (12 October 2009) new implementation of Gamma function in configuration weight |
---|
| 34 | // J. M. Quesada (09 March 2010) Kappa is set to 6. |
---|
[819] | 35 | |
---|
| 36 | #include "G4FermiConfiguration.hh" |
---|
| 37 | #include "G4FermiPhaseSpaceDecay.hh" |
---|
| 38 | #include <set> |
---|
| 39 | |
---|
| 40 | // Kappa = V/V_0 it is used in calculation of Coulomb energy |
---|
| 41 | // Kappa is adimensional |
---|
[1315] | 42 | // JMQ 090310 according to model developer (A. Botvina) no theoretical constraint for kappa below 10 |
---|
| 43 | // kappa values larger than 1 seem to provide better results. 6 is a good choice |
---|
| 44 | // const G4double G4FermiConfiguration::Kappa = 1.0; |
---|
| 45 | const G4double G4FermiConfiguration::Kappa = 6.0; |
---|
[819] | 46 | |
---|
| 47 | // r0 is the nuclear radius |
---|
| 48 | const G4double G4FermiConfiguration::r0 = 1.3*fermi; |
---|
| 49 | |
---|
| 50 | |
---|
| 51 | G4double G4FermiConfiguration::CoulombBarrier(void) |
---|
| 52 | { |
---|
| 53 | // Calculates Coulomb Barrier (MeV) for given channel with K fragments. |
---|
| 54 | const G4double Coef = (3./5.)*(elm_coupling/r0)*std::pow(1./(1.+Kappa), 1./3.); |
---|
| 55 | |
---|
| 56 | G4double SumA = 0; |
---|
| 57 | G4double SumZ = 0; |
---|
| 58 | G4double CoulombEnergy = 0.; |
---|
| 59 | for (std::vector<const G4VFermiFragment*>::iterator i = Configuration.begin(); |
---|
| 60 | i != Configuration.end(); i++) |
---|
| 61 | { |
---|
| 62 | G4double z = static_cast<G4double>((*i)->GetZ()); |
---|
| 63 | G4double a = static_cast<G4double>((*i)->GetA()); |
---|
| 64 | CoulombEnergy += (z*z) / std::pow(a, 1./3.); |
---|
| 65 | SumA += a; |
---|
| 66 | SumZ += z; |
---|
| 67 | } |
---|
| 68 | CoulombEnergy -= SumZ*SumZ/std::pow(SumA, 1./3.); |
---|
| 69 | return -Coef * CoulombEnergy; |
---|
| 70 | } |
---|
| 71 | |
---|
| 72 | |
---|
| 73 | |
---|
| 74 | |
---|
| 75 | G4double G4FermiConfiguration::DecayProbability(const G4int A, const G4double TotalE) |
---|
| 76 | // Decay probability for a given channel with K fragments |
---|
| 77 | { |
---|
| 78 | // A: Atomic Weight |
---|
| 79 | // TotalE: Total energy of nucleus |
---|
| 80 | |
---|
| 81 | |
---|
| 82 | G4double KineticEnergy = TotalE; // MeV |
---|
| 83 | G4double ProdMass = 1.0; |
---|
| 84 | G4double SumMass = 0.0; |
---|
| 85 | G4double S_n = 1.0; |
---|
| 86 | std::set<G4int> combSet; |
---|
| 87 | std::multiset<G4int> combmSet; |
---|
| 88 | |
---|
| 89 | for (std::vector<const G4VFermiFragment*>::iterator i = Configuration.begin(); |
---|
| 90 | i != Configuration.end(); i++) |
---|
| 91 | { |
---|
| 92 | G4int a = (*i)->GetA(); |
---|
| 93 | combSet.insert(a); |
---|
| 94 | combmSet.insert(a); |
---|
| 95 | G4double m = (*i)->GetFragmentMass(); |
---|
| 96 | ProdMass *= m; |
---|
| 97 | SumMass += m; |
---|
| 98 | // Spin factor S_n |
---|
| 99 | S_n *= (*i)->GetPolarization(); |
---|
| 100 | KineticEnergy -= m + (*i)->GetExcitationEnergy(); |
---|
| 101 | } |
---|
| 102 | |
---|
| 103 | // Check that there is enough energy to produce K fragments |
---|
| 104 | if (KineticEnergy <= 0.0) return 0.0; |
---|
| 105 | if ((KineticEnergy -= this->CoulombBarrier()) <= 0.0) return 0.0; |
---|
| 106 | |
---|
| 107 | |
---|
| 108 | G4double MassFactor = ProdMass/SumMass; |
---|
| 109 | MassFactor *= std::sqrt(MassFactor); |
---|
| 110 | |
---|
| 111 | // Number of fragments |
---|
| 112 | G4int K = Configuration.size(); |
---|
| 113 | |
---|
| 114 | // This is the constant (doesn't depend on nucleus) part |
---|
| 115 | const G4double ConstCoeff = std::pow(r0/hbarc,3.0)*Kappa*std::sqrt(2.0/pi)/3.0; |
---|
| 116 | G4double Coeff = std::pow(ConstCoeff*A,K-1); |
---|
| 117 | |
---|
[1196] | 118 | //JMQ 111009 Bug fixed: gamma function for odd K was wrong by a factor 2 |
---|
| 119 | // new implementation explicitely according to standard properties of Gamma function |
---|
[819] | 120 | // Calculation of 1/Gamma(3(k-1)/2) |
---|
| 121 | G4double Gamma = 1.0; |
---|
[1196] | 122 | // G4double arg = 3.0*(K-1)/2.0 - 1.0; |
---|
| 123 | // while (arg > 1.1) |
---|
| 124 | // { |
---|
| 125 | // Gamma *= arg; |
---|
| 126 | // arg--; |
---|
| 127 | // } |
---|
| 128 | // if ((K-1)%2 == 1) Gamma *= std::sqrt(pi); |
---|
| 129 | |
---|
| 130 | if ((K-1)%2 != 1) |
---|
| 131 | |
---|
[819] | 132 | { |
---|
[1196] | 133 | G4double arg = 3.0*(K-1)/2.0 - 1.0; |
---|
| 134 | while (arg > 1.1) |
---|
| 135 | { |
---|
| 136 | Gamma *= arg; |
---|
| 137 | arg--; |
---|
| 138 | } |
---|
[819] | 139 | } |
---|
[1196] | 140 | else { |
---|
| 141 | G4double n= 3.0*K/2.0-2.0; |
---|
| 142 | G4double arg2=2*n-1; |
---|
| 143 | while (arg2>1.1) |
---|
| 144 | { |
---|
| 145 | Gamma*=arg2; |
---|
| 146 | arg2-=2; |
---|
| 147 | } |
---|
[1228] | 148 | Gamma=Gamma/std::pow(2.,n)*std::sqrt(pi); |
---|
[1196] | 149 | } |
---|
| 150 | // end of new implementation of Gamma function |
---|
[819] | 151 | |
---|
| 152 | |
---|
| 153 | // Permutation Factor G_n |
---|
| 154 | G4double G_n = 1.0; |
---|
| 155 | for (std::set<G4int>::iterator s = combSet.begin(); s != combSet.end(); ++s) |
---|
| 156 | { |
---|
| 157 | for (G4int ni = combmSet.count(*s); ni > 1; ni--) G_n *= ni; |
---|
| 158 | } |
---|
| 159 | |
---|
| 160 | G4double Weight = Coeff * MassFactor * (S_n / G_n) / Gamma; |
---|
| 161 | Weight *= std::pow(KineticEnergy,3.0*(K-1)/2.0)/KineticEnergy; |
---|
| 162 | |
---|
| 163 | return Weight; |
---|
| 164 | } |
---|
| 165 | |
---|
| 166 | |
---|
| 167 | G4FragmentVector * G4FermiConfiguration::GetFragments(const G4Fragment & theNucleus) |
---|
| 168 | { |
---|
| 169 | |
---|
| 170 | G4FermiPhaseSpaceDecay thePhaseSpace; |
---|
| 171 | |
---|
| 172 | // Calculate Momenta of K fragments |
---|
| 173 | G4double M = theNucleus.GetMomentum().m(); |
---|
| 174 | std::vector<G4double> m; |
---|
| 175 | m.reserve(Configuration.size()); |
---|
| 176 | std::vector<const G4VFermiFragment*>::iterator i; |
---|
| 177 | for (i = Configuration.begin(); i != Configuration.end(); ++i) |
---|
| 178 | { |
---|
| 179 | m.push_back( (*i)->GetTotalEnergy() ); |
---|
| 180 | } |
---|
| 181 | std::vector<G4LorentzVector*>* MomentumComponents = |
---|
| 182 | thePhaseSpace.Decay(M,m); |
---|
| 183 | |
---|
| 184 | G4FragmentVector * theResult = new G4FragmentVector; |
---|
| 185 | |
---|
| 186 | G4ThreeVector boostVector = theNucleus.GetMomentum().boostVector(); |
---|
| 187 | |
---|
| 188 | // Go back to the Lab Frame |
---|
| 189 | for (i = Configuration.begin(); i != Configuration.end(); ++i) |
---|
| 190 | { |
---|
| 191 | #ifdef G4NO_ISO_VECDIST |
---|
| 192 | std::vector<const G4VFermiFragment*>::difference_type n = 0; |
---|
| 193 | std::distance(Configuration.begin(), i, n); |
---|
| 194 | G4LorentzVector FourMomentum(*(MomentumComponents->operator[](n))); |
---|
| 195 | #else |
---|
| 196 | G4LorentzVector FourMomentum(*(MomentumComponents-> |
---|
| 197 | operator[](std::distance(Configuration.begin(),i)))); |
---|
| 198 | #endif |
---|
| 199 | |
---|
| 200 | // Lorentz boost |
---|
| 201 | FourMomentum.boost(boostVector); |
---|
| 202 | |
---|
| 203 | G4FragmentVector * fragment = (*i)->GetFragment(FourMomentum); |
---|
| 204 | |
---|
| 205 | |
---|
| 206 | for (G4FragmentVector::reverse_iterator ri = fragment->rbegin(); |
---|
| 207 | ri != fragment->rend(); ++ri) |
---|
| 208 | { |
---|
| 209 | theResult->push_back(*ri); |
---|
| 210 | } |
---|
| 211 | delete fragment; |
---|
| 212 | } |
---|
| 213 | |
---|
| 214 | if (!MomentumComponents->empty()) |
---|
| 215 | { |
---|
| 216 | std::for_each(MomentumComponents->begin(),MomentumComponents->end(), |
---|
| 217 | DeleteFragment()); |
---|
| 218 | } |
---|
| 219 | |
---|
| 220 | delete MomentumComponents; |
---|
| 221 | |
---|
| 222 | return theResult; |
---|
| 223 | } |
---|
| 224 | |
---|
| 225 | |
---|
| 226 | G4ParticleMomentum G4FermiConfiguration::IsotropicVector(const G4double Magnitude) |
---|
| 227 | // Samples a isotropic random vectorwith a magnitud given by Magnitude. |
---|
| 228 | // By default Magnitude = 1.0 |
---|
| 229 | { |
---|
| 230 | G4double CosTheta = 1.0 - 2.0*G4UniformRand(); |
---|
| 231 | G4double SinTheta = std::sqrt(1.0 - CosTheta*CosTheta); |
---|
| 232 | G4double Phi = twopi*G4UniformRand(); |
---|
| 233 | G4ParticleMomentum Vector(Magnitude*std::cos(Phi)*SinTheta, |
---|
| 234 | Magnitude*std::sin(Phi)*SinTheta, |
---|
| 235 | Magnitude*CosTheta); |
---|
| 236 | return Vector; |
---|
| 237 | } |
---|