[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 | // |
---|
[1340] | 26 | // $Id: G4PreCompoundEmission.cc,v 1.32 2010/09/01 15:11:10 vnivanch Exp $ |
---|
| 27 | // GEANT4 tag $Name: geant4-09-03-ref-09 $ |
---|
[819] | 28 | // |
---|
[1055] | 29 | // ------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // GEANT4 Class file |
---|
| 32 | // |
---|
| 33 | // |
---|
| 34 | // File name: G4PreCompoundEmission |
---|
| 35 | // |
---|
| 36 | // Author: V.Lara |
---|
| 37 | // |
---|
| 38 | // Modified: |
---|
[1315] | 39 | // 15.01.2010 J.M.Quesada added protection against unphysical values of parameter an |
---|
| 40 | // 19.01.2010 V.Ivanchenko simplified computation of parameter an, sample cosTheta |
---|
| 41 | // instead of theta; protect all calls to sqrt |
---|
[1340] | 42 | // 20.08.2010 V.Ivanchenko added G4Pow and G4PreCompoundParameters pointers |
---|
| 43 | // use int Z and A and cleanup |
---|
[1055] | 44 | // |
---|
[819] | 45 | |
---|
| 46 | #include "G4PreCompoundEmission.hh" |
---|
| 47 | #include "G4PreCompoundParameters.hh" |
---|
| 48 | #include "G4PreCompoundEmissionFactory.hh" |
---|
| 49 | #include "G4HETCEmissionFactory.hh" |
---|
[1340] | 50 | #include "G4HadronicException.hh" |
---|
| 51 | #include "G4Pow.hh" |
---|
| 52 | #include "Randomize.hh" |
---|
[819] | 53 | |
---|
| 54 | G4PreCompoundEmission::G4PreCompoundEmission() |
---|
| 55 | { |
---|
| 56 | theFragmentsFactory = new G4PreCompoundEmissionFactory(); |
---|
[1340] | 57 | theFragmentsVector = |
---|
| 58 | new G4PreCompoundFragmentVector(theFragmentsFactory->GetFragmentVector()); |
---|
| 59 | g4pow = G4Pow::GetInstance(); |
---|
| 60 | theParameters = G4PreCompoundParameters::GetAddress(); |
---|
[819] | 61 | } |
---|
| 62 | |
---|
| 63 | G4PreCompoundEmission::~G4PreCompoundEmission() |
---|
| 64 | { |
---|
[1340] | 65 | if (theFragmentsFactory) { delete theFragmentsFactory; } |
---|
| 66 | if (theFragmentsVector) { delete theFragmentsVector; } |
---|
[819] | 67 | } |
---|
| 68 | |
---|
| 69 | void G4PreCompoundEmission::SetDefaultModel() |
---|
| 70 | { |
---|
[1340] | 71 | if (theFragmentsFactory) { delete theFragmentsFactory; } |
---|
[819] | 72 | theFragmentsFactory = new G4PreCompoundEmissionFactory(); |
---|
| 73 | if (theFragmentsVector) |
---|
| 74 | { |
---|
| 75 | theFragmentsVector->SetVector(theFragmentsFactory->GetFragmentVector()); |
---|
| 76 | } |
---|
| 77 | else |
---|
| 78 | { |
---|
[1340] | 79 | theFragmentsVector = |
---|
| 80 | new G4PreCompoundFragmentVector(theFragmentsFactory->GetFragmentVector()); |
---|
[819] | 81 | } |
---|
| 82 | return; |
---|
| 83 | } |
---|
| 84 | |
---|
| 85 | void G4PreCompoundEmission::SetHETCModel() |
---|
| 86 | { |
---|
| 87 | if (theFragmentsFactory) delete theFragmentsFactory; |
---|
| 88 | theFragmentsFactory = new G4HETCEmissionFactory(); |
---|
| 89 | if (theFragmentsVector) |
---|
| 90 | { |
---|
| 91 | theFragmentsVector->SetVector(theFragmentsFactory->GetFragmentVector()); |
---|
| 92 | } |
---|
| 93 | else |
---|
| 94 | { |
---|
[1340] | 95 | theFragmentsVector = |
---|
| 96 | new G4PreCompoundFragmentVector(theFragmentsFactory->GetFragmentVector()); |
---|
[819] | 97 | } |
---|
| 98 | return; |
---|
| 99 | } |
---|
| 100 | |
---|
| 101 | G4ReactionProduct * G4PreCompoundEmission::PerformEmission(G4Fragment & aFragment) |
---|
| 102 | { |
---|
| 103 | // Choose a Fragment for emission |
---|
[1340] | 104 | G4VPreCompoundFragment * thePreFragment = theFragmentsVector->ChooseFragment(); |
---|
| 105 | if (thePreFragment == 0) |
---|
[819] | 106 | { |
---|
[1340] | 107 | G4cout << "G4PreCompoundEmission::PerformEmission : I couldn't choose a fragment\n" |
---|
[819] | 108 | << "while trying to de-excite\n" |
---|
[1340] | 109 | << aFragment << G4endl; |
---|
[819] | 110 | throw G4HadronicException(__FILE__, __LINE__, ""); |
---|
| 111 | } |
---|
| 112 | // Kinetic Energy of emitted fragment |
---|
[1340] | 113 | G4double kinEnergyOfEmittedFragment = thePreFragment->GetKineticEnergy(aFragment); |
---|
| 114 | if(kinEnergyOfEmittedFragment < 0.0) { kinEnergyOfEmittedFragment = 0.0; } |
---|
[819] | 115 | |
---|
| 116 | // Calculate the fragment momentum (three vector) |
---|
[1340] | 117 | AngularDistribution(thePreFragment,aFragment,kinEnergyOfEmittedFragment); |
---|
[819] | 118 | |
---|
| 119 | // Mass of emittef fragment |
---|
[1340] | 120 | G4double EmittedMass = thePreFragment->GetNuclearMass(); |
---|
[819] | 121 | |
---|
| 122 | // Now we can calculate the four momentum |
---|
| 123 | // both options are valid and give the same result but 2nd one is faster |
---|
[1340] | 124 | G4LorentzVector Emitted4Momentum(theFinalMomentum, |
---|
| 125 | EmittedMass + kinEnergyOfEmittedFragment); |
---|
[819] | 126 | |
---|
| 127 | // Perform Lorentz boost |
---|
[1340] | 128 | G4LorentzVector Rest4Momentum = aFragment.GetMomentum(); |
---|
| 129 | Emitted4Momentum.boost(Rest4Momentum.boostVector()); |
---|
[819] | 130 | |
---|
| 131 | // Set emitted fragment momentum |
---|
[1340] | 132 | thePreFragment->SetMomentum(Emitted4Momentum); |
---|
[819] | 133 | |
---|
| 134 | // NOW THE RESIDUAL NUCLEUS |
---|
| 135 | // ------------------------ |
---|
| 136 | |
---|
[1340] | 137 | Rest4Momentum -= Emitted4Momentum; |
---|
[819] | 138 | |
---|
| 139 | // Update nucleus parameters: |
---|
| 140 | // -------------------------- |
---|
[962] | 141 | |
---|
[819] | 142 | // Number of excitons |
---|
| 143 | aFragment.SetNumberOfParticles(aFragment.GetNumberOfParticles()- |
---|
[1340] | 144 | thePreFragment->GetA()); |
---|
[819] | 145 | // Number of charges |
---|
| 146 | aFragment.SetNumberOfCharged(aFragment.GetNumberOfCharged()- |
---|
[1340] | 147 | thePreFragment->GetZ()); |
---|
[819] | 148 | |
---|
[1340] | 149 | // Z and A |
---|
| 150 | aFragment.SetZandA_asInt(thePreFragment->GetRestZ(), |
---|
| 151 | thePreFragment->GetRestA()); |
---|
[819] | 152 | |
---|
[1340] | 153 | // Update nucleus momentum |
---|
| 154 | // A check on consistence of Z, A, and mass will be performed |
---|
| 155 | aFragment.SetMomentum(Rest4Momentum); |
---|
| 156 | |
---|
| 157 | // Create a G4ReactionProduct |
---|
| 158 | G4ReactionProduct * MyRP = thePreFragment->GetReactionProduct(); |
---|
[962] | 159 | |
---|
[1340] | 160 | //G4cout << "G4PreCompoundEmission::Fragment emitted" << G4endl; |
---|
| 161 | //G4cout << thePreFragment << G4endl; |
---|
[819] | 162 | |
---|
| 163 | return MyRP; |
---|
| 164 | } |
---|
| 165 | |
---|
[1340] | 166 | void |
---|
| 167 | G4PreCompoundEmission::AngularDistribution(G4VPreCompoundFragment* thePreFragment, |
---|
[1315] | 168 | const G4Fragment& aFragment, |
---|
[1340] | 169 | G4double ekin) |
---|
[819] | 170 | { |
---|
[1340] | 171 | G4int p = aFragment.GetNumberOfParticles(); |
---|
| 172 | G4int h = aFragment.GetNumberOfHoles(); |
---|
[819] | 173 | G4double U = aFragment.GetExcitationEnergy(); |
---|
[1315] | 174 | |
---|
[819] | 175 | // Emission particle separation energy |
---|
[1340] | 176 | G4double Bemission = thePreFragment->GetBindingEnergy(); |
---|
[1315] | 177 | |
---|
[819] | 178 | // Fermi energy |
---|
[1340] | 179 | G4double Ef = theParameters->GetFermiEnergy(); |
---|
[819] | 180 | |
---|
| 181 | // |
---|
| 182 | // G4EvaporationLevelDensityParameter theLDP; |
---|
| 183 | // G4double g = (6.0/pi2)*aFragment.GetA()* |
---|
[1315] | 184 | |
---|
[1340] | 185 | G4double g = (6.0/pi2)*aFragment.GetA_asInt()*theParameters->GetLevelDensity(); |
---|
[819] | 186 | |
---|
| 187 | // Average exciton energy relative to bottom of nuclear well |
---|
[1340] | 188 | G4double Eav = 2*p*(p+1)/((p+h)*g); |
---|
[819] | 189 | |
---|
| 190 | // Excitation energy relative to the Fermi Level |
---|
| 191 | G4double Uf = std::max(U - (p - h)*Ef , 0.0); |
---|
| 192 | // G4double Uf = U - KineticEnergyOfEmittedFragment - Bemission; |
---|
| 193 | |
---|
[1340] | 194 | G4double w_num = rho(p+1, h, g, Uf, Ef); |
---|
| 195 | G4double w_den = rho(p, h, g, Uf, Ef); |
---|
[819] | 196 | if (w_num > 0.0 && w_den > 0.0) |
---|
| 197 | { |
---|
| 198 | Eav *= (w_num/w_den); |
---|
| 199 | Eav += - Uf/(p+h) + Ef; |
---|
| 200 | } |
---|
| 201 | else |
---|
| 202 | { |
---|
| 203 | Eav = Ef; |
---|
| 204 | } |
---|
| 205 | |
---|
[1315] | 206 | // VI + JMQ 19/01/2010 update computation of the parameter an |
---|
| 207 | // |
---|
| 208 | G4double an = 0.0; |
---|
| 209 | G4double Eeff = ekin + Bemission + Ef; |
---|
| 210 | if(ekin > DBL_MIN && Eeff > DBL_MIN) { |
---|
| 211 | |
---|
| 212 | G4double zeta = std::max(1.0,9.3/std::sqrt(ekin/MeV)); |
---|
| 213 | |
---|
[1340] | 214 | // This should be the projectile energy. If I would know which is |
---|
| 215 | // the projectile (proton, neutron) I could remove the binding energy. |
---|
| 216 | // But, what happens if INC precedes precompound? This approximation |
---|
| 217 | // seems to work well enough |
---|
| 218 | G4double ProjEnergy = aFragment.GetExcitationEnergy(); |
---|
[1315] | 219 | |
---|
[1340] | 220 | an = 3*std::sqrt((ProjEnergy+Ef)*Eeff)/(zeta*Eav); |
---|
| 221 | |
---|
[1315] | 222 | G4int ne = aFragment.GetNumberOfExcitons() - 1; |
---|
| 223 | if ( ne > 1 ) { an /= (G4double)ne; } |
---|
[819] | 224 | |
---|
[1315] | 225 | // protection of exponent |
---|
| 226 | if ( an > 10. ) { an = 10.; } |
---|
[1196] | 227 | } |
---|
| 228 | |
---|
[1315] | 229 | // sample cosine of theta and not theta as in old versions |
---|
| 230 | G4double random = G4UniformRand(); |
---|
| 231 | G4double cost; |
---|
| 232 | |
---|
| 233 | if(an < 0.1) { cost = 1. - 2*random; } |
---|
| 234 | else { |
---|
| 235 | G4double exp2an = std::exp(-2*an); |
---|
| 236 | cost = 1. + std::log(1-random*(1-exp2an))/an; |
---|
| 237 | if(cost > 1.) { cost = 1.; } |
---|
| 238 | else if(cost < -1.) {cost = -1.; } |
---|
| 239 | } |
---|
| 240 | |
---|
[1340] | 241 | G4double phi = CLHEP::twopi*G4UniformRand(); |
---|
[819] | 242 | |
---|
| 243 | // Calculate the momentum magnitude of emitted fragment |
---|
[1340] | 244 | G4double pmag = std::sqrt(ekin*(ekin + 2.0*thePreFragment->GetNuclearMass())); |
---|
[819] | 245 | |
---|
[1315] | 246 | G4double sint = std::sqrt((1.0-cost)*(1.0+cost)); |
---|
[819] | 247 | |
---|
[1340] | 248 | theFinalMomentum.set(pmag*std::cos(phi)*sint,pmag*std::sin(phi)*sint,pmag*cost); |
---|
| 249 | |
---|
[819] | 250 | // theta is the angle wrt the incident direction |
---|
[1340] | 251 | G4ThreeVector theIncidentDirection = aFragment.GetMomentum().vect().unit(); |
---|
| 252 | theFinalMomentum.rotateUz(theIncidentDirection); |
---|
[819] | 253 | } |
---|
| 254 | |
---|
[1340] | 255 | G4double G4PreCompoundEmission::rho(G4int p, G4int h, G4double g, |
---|
| 256 | G4double E, G4double Ef) const |
---|
[1315] | 257 | { |
---|
| 258 | // 25.02.2010 V.Ivanchenko added more protections |
---|
| 259 | G4double Aph = (p*p + h*h + p - 3.0*h)/(4.0*g); |
---|
| 260 | // G4double alpha = (p*p + h*h)/(2.0*g); |
---|
[962] | 261 | |
---|
[1315] | 262 | if ( E - Aph < 0.0) { return 0.0; } |
---|
[962] | 263 | |
---|
[1340] | 264 | G4double logConst = (p+h)*std::log(g) |
---|
| 265 | - g4pow->logfactorial(p+h-1) - g4pow->logfactorial(p) - g4pow->logfactorial(h); |
---|
[962] | 266 | |
---|
[1315] | 267 | // initialise values using j=0 |
---|
[962] | 268 | |
---|
| 269 | G4double t1=1; |
---|
| 270 | G4double t2=1; |
---|
[1340] | 271 | G4double logt3 = (p+h-1) * std::log(E-Aph) + logConst; |
---|
[1315] | 272 | const G4double logmax = 200.; |
---|
| 273 | if(logt3 > logmax) { logt3 = logmax; } |
---|
| 274 | G4double tot = std::exp( logt3 ); |
---|
[962] | 275 | |
---|
[1315] | 276 | // and now sum rest of terms |
---|
| 277 | // 25.02.2010 V.Ivanchenko change while to for loop and cleanup |
---|
| 278 | G4double Eeff = E - Aph; |
---|
| 279 | for(G4int j=1; j<=h; ++j) |
---|
[962] | 280 | { |
---|
[1315] | 281 | Eeff -= Ef; |
---|
| 282 | if(Eeff < 0.0) { break; } |
---|
| 283 | t1 *= -1.; |
---|
| 284 | t2 *= (G4double)(h+1-j)/(G4double)j; |
---|
| 285 | logt3 = (p+h-1) * std::log( Eeff) + logConst; |
---|
| 286 | if(logt3 > logmax) { logt3 = logmax; } |
---|
| 287 | tot += t1*t2*std::exp(logt3); |
---|
[962] | 288 | } |
---|
| 289 | |
---|
| 290 | return tot; |
---|
| 291 | } |
---|