| 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: G4HEProtonInelastic.cc,v 1.16 2010/11/29 05:44:44 dennis Exp $
|
|---|
| 27 | // GEANT4 tag $Name: geant4-09-04-ref-00 $
|
|---|
| 28 | //
|
|---|
| 29 |
|
|---|
| 30 | #include "globals.hh"
|
|---|
| 31 | #include "G4ios.hh"
|
|---|
| 32 |
|
|---|
| 33 | // G4 Process: Gheisha High Energy Collision model.
|
|---|
| 34 | // This includes the high energy cascading model, the two-body-resonance model
|
|---|
| 35 | // and the low energy two-body model. Not included are the low energy stuff
|
|---|
| 36 | // like nuclear reactions, nuclear fission without any cascading and all
|
|---|
| 37 | // processes for particles at rest.
|
|---|
| 38 | // First work done by J.L.Chuma and F.W.Jones, TRIUMF, June 96.
|
|---|
| 39 | // H. Fesefeldt, RWTH-Aachen, 23-October-1996
|
|---|
| 40 | // Last modified: 29-July-1998.
|
|---|
| 41 |
|
|---|
| 42 | #include "G4HEProtonInelastic.hh"
|
|---|
| 43 |
|
|---|
| 44 | G4HadFinalState*
|
|---|
| 45 | G4HEProtonInelastic::ApplyYourself(const G4HadProjectile& aTrack,
|
|---|
| 46 | G4Nucleus& targetNucleus)
|
|---|
| 47 | {
|
|---|
| 48 | G4HEVector* pv = new G4HEVector[MAXPART];
|
|---|
| 49 | const G4HadProjectile* aParticle = &aTrack;
|
|---|
| 50 | const G4double A = targetNucleus.GetN();
|
|---|
| 51 | const G4double Z = targetNucleus.GetZ();
|
|---|
| 52 | G4HEVector incidentParticle(aParticle);
|
|---|
| 53 |
|
|---|
| 54 | G4double atomicNumber = Z;
|
|---|
| 55 | G4double atomicWeight = A;
|
|---|
| 56 |
|
|---|
| 57 | if (verboseLevel > 1)
|
|---|
| 58 | G4cout << "Z , A = " << atomicNumber << " " << atomicWeight << G4endl;
|
|---|
| 59 |
|
|---|
| 60 | G4int incidentCode = incidentParticle.getCode();
|
|---|
| 61 | G4double incidentMass = incidentParticle.getMass();
|
|---|
| 62 | G4double incidentTotalEnergy = incidentParticle.getEnergy();
|
|---|
| 63 | G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
|
|---|
| 64 | G4double incidentKineticEnergy = incidentParticle.getKineticEnergy();
|
|---|
| 65 |
|
|---|
| 66 | if (incidentKineticEnergy < 1.)
|
|---|
| 67 | G4cout << "GHEProtonInelastic: incident energy < 1 GeV" << G4endl;
|
|---|
| 68 |
|
|---|
| 69 | if (verboseLevel > 1) {
|
|---|
| 70 | G4cout << "G4HEProtonInelastic::ApplyYourself" << G4endl;
|
|---|
| 71 | G4cout << "incident particle " << incidentParticle.getName() << " "
|
|---|
| 72 | << "mass " << incidentMass << " "
|
|---|
| 73 | << "kinetic energy " << incidentKineticEnergy
|
|---|
| 74 | << G4endl;
|
|---|
| 75 | G4cout << "target material with (A,Z) = ("
|
|---|
| 76 | << atomicWeight << "," << atomicNumber << ")" << G4endl;
|
|---|
| 77 | }
|
|---|
| 78 |
|
|---|
| 79 | G4double inelasticity = NuclearInelasticity(incidentKineticEnergy,
|
|---|
| 80 | atomicWeight, atomicNumber);
|
|---|
| 81 | if (verboseLevel > 1)
|
|---|
| 82 | G4cout << "nuclear inelasticity = " << inelasticity << G4endl;
|
|---|
| 83 |
|
|---|
| 84 | incidentKineticEnergy -= inelasticity;
|
|---|
| 85 |
|
|---|
| 86 | G4double excitationEnergyGNP = 0.;
|
|---|
| 87 | G4double excitationEnergyDTA = 0.;
|
|---|
| 88 |
|
|---|
| 89 | G4double excitation = NuclearExcitation(incidentKineticEnergy,
|
|---|
| 90 | atomicWeight, atomicNumber,
|
|---|
| 91 | excitationEnergyGNP,
|
|---|
| 92 | excitationEnergyDTA);
|
|---|
| 93 |
|
|---|
| 94 | if (verboseLevel > 1)
|
|---|
| 95 | G4cout << "nuclear excitation = " << excitation << " "
|
|---|
| 96 | << excitationEnergyGNP << " " << excitationEnergyDTA << G4endl;
|
|---|
| 97 |
|
|---|
| 98 | incidentKineticEnergy -= excitation;
|
|---|
| 99 | incidentTotalEnergy = incidentKineticEnergy + incidentMass;
|
|---|
| 100 | incidentTotalMomentum = std::sqrt( (incidentTotalEnergy-incidentMass)
|
|---|
| 101 | *(incidentTotalEnergy+incidentMass));
|
|---|
| 102 |
|
|---|
| 103 | G4HEVector targetParticle;
|
|---|
| 104 | if (G4UniformRand() < atomicNumber/atomicWeight) {
|
|---|
| 105 | targetParticle.setDefinition("Proton");
|
|---|
| 106 | } else {
|
|---|
| 107 | targetParticle.setDefinition("Neutron");
|
|---|
| 108 | }
|
|---|
| 109 |
|
|---|
| 110 | G4double targetMass = targetParticle.getMass();
|
|---|
| 111 | G4double centerOfMassEnergy = std::sqrt(incidentMass*incidentMass
|
|---|
| 112 | + targetMass*targetMass
|
|---|
| 113 | + 2.0*targetMass*incidentTotalEnergy);
|
|---|
| 114 | G4double availableEnergy = centerOfMassEnergy - targetMass - incidentMass;
|
|---|
| 115 |
|
|---|
| 116 | // In the original Gheisha code, the inElastic flag was defined as follows:
|
|---|
| 117 | // G4bool inElastic = InElasticCrossSectionInFirstInt
|
|---|
| 118 | // (availableEnergy, incidentCode, incidentTotalMomentum);
|
|---|
| 119 | // For unknown reasons, it was replaced by the following code in Geant???
|
|---|
| 120 |
|
|---|
| 121 | G4bool inElastic = true;
|
|---|
| 122 | // if (G4UniformRand() < elasticCrossSection/totalCrossSection) inElastic = false;
|
|---|
| 123 |
|
|---|
| 124 | vecLength = 0;
|
|---|
| 125 |
|
|---|
| 126 | if (verboseLevel > 1)
|
|---|
| 127 | G4cout << "ApplyYourself: CallFirstIntInCascade for particle "
|
|---|
| 128 | << incidentCode << G4endl;
|
|---|
| 129 |
|
|---|
| 130 | G4bool successful = false;
|
|---|
| 131 |
|
|---|
| 132 | FirstIntInCasProton(inElastic, availableEnergy, pv, vecLength,
|
|---|
| 133 | incidentParticle, targetParticle, atomicWeight);
|
|---|
| 134 |
|
|---|
| 135 | if (verboseLevel > 1)
|
|---|
| 136 | G4cout << "ApplyYourself::StrangeParticlePairProduction" << G4endl;
|
|---|
| 137 |
|
|---|
| 138 | if ((vecLength > 0) && (availableEnergy > 1.))
|
|---|
| 139 | StrangeParticlePairProduction(availableEnergy, centerOfMassEnergy,
|
|---|
| 140 | pv, vecLength,
|
|---|
| 141 | incidentParticle, targetParticle);
|
|---|
| 142 |
|
|---|
| 143 | HighEnergyCascading(successful, pv, vecLength,
|
|---|
| 144 | excitationEnergyGNP, excitationEnergyDTA,
|
|---|
| 145 | incidentParticle, targetParticle,
|
|---|
| 146 | atomicWeight, atomicNumber);
|
|---|
| 147 | if (!successful)
|
|---|
| 148 | HighEnergyClusterProduction(successful, pv, vecLength,
|
|---|
| 149 | excitationEnergyGNP, excitationEnergyDTA,
|
|---|
| 150 | incidentParticle, targetParticle,
|
|---|
| 151 | atomicWeight, atomicNumber);
|
|---|
| 152 | if (!successful)
|
|---|
| 153 | MediumEnergyCascading(successful, pv, vecLength,
|
|---|
| 154 | excitationEnergyGNP, excitationEnergyDTA,
|
|---|
| 155 | incidentParticle, targetParticle,
|
|---|
| 156 | atomicWeight, atomicNumber);
|
|---|
| 157 |
|
|---|
| 158 | if (!successful)
|
|---|
| 159 | MediumEnergyClusterProduction(successful, pv, vecLength,
|
|---|
| 160 | excitationEnergyGNP, excitationEnergyDTA,
|
|---|
| 161 | incidentParticle, targetParticle,
|
|---|
| 162 | atomicWeight, atomicNumber);
|
|---|
| 163 | if (!successful)
|
|---|
| 164 | QuasiElasticScattering(successful, pv, vecLength,
|
|---|
| 165 | excitationEnergyGNP, excitationEnergyDTA,
|
|---|
| 166 | incidentParticle, targetParticle,
|
|---|
| 167 | atomicWeight, atomicNumber);
|
|---|
| 168 | if (!successful)
|
|---|
| 169 | ElasticScattering(successful, pv, vecLength,
|
|---|
| 170 | incidentParticle,
|
|---|
| 171 | atomicWeight, atomicNumber);
|
|---|
| 172 |
|
|---|
| 173 | if (!successful)
|
|---|
| 174 | G4cout << "GHEInelasticInteraction::ApplyYourself fails to produce final state particles"
|
|---|
| 175 | << G4endl;
|
|---|
| 176 |
|
|---|
| 177 | FillParticleChange(pv, vecLength);
|
|---|
| 178 | delete [] pv;
|
|---|
| 179 | theParticleChange.SetStatusChange(stopAndKill);
|
|---|
| 180 | return &theParticleChange;
|
|---|
| 181 | }
|
|---|
| 182 |
|
|---|
| 183 |
|
|---|
| 184 | void
|
|---|
| 185 | G4HEProtonInelastic::FirstIntInCasProton(G4bool& inElastic,
|
|---|
| 186 | const G4double availableEnergy,
|
|---|
| 187 | G4HEVector pv[],
|
|---|
| 188 | G4int& vecLen,
|
|---|
| 189 | const G4HEVector& incidentParticle,
|
|---|
| 190 | const G4HEVector& targetParticle,
|
|---|
| 191 | const G4double atomicWeight)
|
|---|
| 192 |
|
|---|
| 193 | // Proton undergoes interaction with nucleon within a nucleus. Check if it is
|
|---|
| 194 | // energetically possible to produce pions/kaons. In not, assume nuclear excitation
|
|---|
| 195 | // occurs and input particle is degraded in energy. No other particles are produced.
|
|---|
| 196 | // If reaction is possible, find the correct number of pions/protons/neutrons
|
|---|
| 197 | // produced using an interpolation to multiplicity data. Replace some pions or
|
|---|
| 198 | // protons/neutrons by kaons or strange baryons according to the average
|
|---|
| 199 | // multiplicity per inelastic reaction.
|
|---|
| 200 | {
|
|---|
| 201 | static const G4double expxu = std::log(MAXFLOAT); // upper bound for arg. of exp
|
|---|
| 202 | static const G4double expxl = -expxu; // lower bound for arg. of exp
|
|---|
| 203 |
|
|---|
| 204 | static const G4double protb = 0.7;
|
|---|
| 205 | static const G4double neutb = 0.35;
|
|---|
| 206 | static const G4double c = 1.25;
|
|---|
| 207 |
|
|---|
| 208 | static const G4int numMul = 1200;
|
|---|
| 209 | static const G4int numSec = 60;
|
|---|
| 210 |
|
|---|
| 211 | G4int neutronCode = Neutron.getCode();
|
|---|
| 212 | G4int protonCode = Proton.getCode();
|
|---|
| 213 | G4double pionMass = PionPlus.getMass();
|
|---|
| 214 |
|
|---|
| 215 | G4int targetCode = targetParticle.getCode();
|
|---|
| 216 | G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
|
|---|
| 217 |
|
|---|
| 218 | static G4bool first = true;
|
|---|
| 219 | static G4double protmul[numMul], protnorm[numSec]; // proton constants
|
|---|
| 220 | static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
|
|---|
| 221 |
|
|---|
| 222 | // misc. local variables
|
|---|
| 223 | // np = number of pi+, nm = number of pi-, nz = number of pi0
|
|---|
| 224 |
|
|---|
| 225 | G4int i, counter, nt, np, nm, nz;
|
|---|
| 226 |
|
|---|
| 227 | if (first) {
|
|---|
| 228 | // compute normalization constants, this will only be done once
|
|---|
| 229 | first = false;
|
|---|
| 230 | for (i=0; i<numMul; i++) protmul[i] = 0.0;
|
|---|
| 231 | for (i=0; i<numSec; i++) protnorm[i] = 0.0;
|
|---|
| 232 | counter = -1;
|
|---|
| 233 | for (np=0; np<(numSec/3); np++) {
|
|---|
| 234 | for (nm=Imax(0,np-2); nm<=np; nm++) {
|
|---|
| 235 | for (nz=0; nz<numSec/3; nz++) {
|
|---|
| 236 | if (++counter < numMul) {
|
|---|
| 237 | nt = np+nm+nz;
|
|---|
| 238 | if ( (nt>0) && (nt<=numSec) ) {
|
|---|
| 239 | protmul[counter] = pmltpc(np,nm,nz,nt,protb,c)
|
|---|
| 240 | /(Factorial(2-np+nm)*Factorial(np-nm)) ;
|
|---|
| 241 | protnorm[nt-1] += protmul[counter];
|
|---|
| 242 | }
|
|---|
| 243 | }
|
|---|
| 244 | }
|
|---|
| 245 | }
|
|---|
| 246 | }
|
|---|
| 247 |
|
|---|
| 248 | for( i=0; i<numMul; i++ )neutmul[i] = 0.0;
|
|---|
| 249 | for( i=0; i<numSec; i++ )neutnorm[i] = 0.0;
|
|---|
| 250 | counter = -1;
|
|---|
| 251 | for (np=0; np<numSec/3; np++) {
|
|---|
| 252 | for (nm=Imax(0,np-1); nm<=(np+1); nm++) {
|
|---|
| 253 | for (nz=0; nz<numSec/3; nz++) {
|
|---|
| 254 | if (++counter < numMul) {
|
|---|
| 255 | nt = np+nm+nz;
|
|---|
| 256 | if ( (nt>0) && (nt<=numSec) ) {
|
|---|
| 257 | neutmul[counter] = pmltpc(np,nm,nz,nt,neutb,c)
|
|---|
| 258 | /(Factorial(1-np+nm)*Factorial(1+np-nm));
|
|---|
| 259 | neutnorm[nt-1] += neutmul[counter];
|
|---|
| 260 | }
|
|---|
| 261 | }
|
|---|
| 262 | }
|
|---|
| 263 | }
|
|---|
| 264 | }
|
|---|
| 265 | for (i=0; i<numSec; i++) {
|
|---|
| 266 | if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
|
|---|
| 267 | if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
|
|---|
| 268 | }
|
|---|
| 269 | } // end of initialization
|
|---|
| 270 |
|
|---|
| 271 |
|
|---|
| 272 | // initialize the first two places
|
|---|
| 273 | // the same as beam and target
|
|---|
| 274 | pv[0] = incidentParticle;
|
|---|
| 275 | pv[1] = targetParticle;
|
|---|
| 276 | vecLen = 2;
|
|---|
| 277 |
|
|---|
| 278 | if( !inElastic )
|
|---|
| 279 | { // quasi-elastic scattering, no pions produced
|
|---|
| 280 | if( targetCode == neutronCode )
|
|---|
| 281 | {
|
|---|
| 282 | G4double cech[] = {0.50, 0.45, 0.40, 0.35, 0.30, 0.25, 0.06, 0.04, 0.005, 0.};
|
|---|
| 283 | G4int iplab = G4int( Amin( 9.0, incidentTotalMomentum*2.5 ) );
|
|---|
| 284 | if( G4UniformRand() < cech[iplab]/std::pow(atomicWeight,0.42) )
|
|---|
| 285 | { // charge exchange pi+ n -> pi0 p
|
|---|
| 286 | pv[0] = PionZero;
|
|---|
| 287 | pv[1] = Proton;
|
|---|
| 288 | }
|
|---|
| 289 | }
|
|---|
| 290 | return;
|
|---|
| 291 | }
|
|---|
| 292 | else if (availableEnergy <= pionMass)
|
|---|
| 293 | return;
|
|---|
| 294 |
|
|---|
| 295 | // inelastic scattering
|
|---|
| 296 |
|
|---|
| 297 | np = 0, nm = 0, nz = 0;
|
|---|
| 298 | G4double eab = availableEnergy;
|
|---|
| 299 | G4int ieab = G4int( eab*5.0 );
|
|---|
| 300 |
|
|---|
| 301 | G4double supp[] = {0., 0.4, 0.55, 0.65, 0.75, 0.82, 0.86, 0.90, 0.94, 0.98};
|
|---|
| 302 | if( (ieab <= 9) && (G4UniformRand() >= supp[ieab]) )
|
|---|
| 303 | {
|
|---|
| 304 | // suppress high multiplicity events at low momentum
|
|---|
| 305 | // only one additional pion will be produced
|
|---|
| 306 | G4double w0, wp, wm, wt, ran;
|
|---|
| 307 | if( targetCode == protonCode ) // target is a proton
|
|---|
| 308 | {
|
|---|
| 309 | w0 = - sqr(1.+protb)/(2.*c*c);
|
|---|
| 310 | wp = w0 = std::exp(w0);
|
|---|
| 311 | if( G4UniformRand() < w0/(w0+wp) )
|
|---|
| 312 | { np = 0; nm = 0; nz = 1; }
|
|---|
| 313 | else
|
|---|
| 314 | { np = 1; nm = 0; nz = 0; }
|
|---|
| 315 | }
|
|---|
| 316 | else
|
|---|
| 317 | { // target is a neutron
|
|---|
| 318 | w0 = -sqr(1.+neutb)/(2.*c*c);
|
|---|
| 319 | w0 = std::exp(w0);
|
|---|
| 320 | wp = w0/2.;
|
|---|
| 321 | wm = -sqr(-1.+neutb)/(2.*c*c);
|
|---|
| 322 | wm = std::exp(wm)/2.;
|
|---|
| 323 | wt = w0+wp+wm;
|
|---|
| 324 | wp = w0+wp;
|
|---|
| 325 | ran = G4UniformRand();
|
|---|
| 326 | if( ran < w0/wt)
|
|---|
| 327 | { np = 0; nm = 0; nz = 1; }
|
|---|
| 328 | else if( ran < wp/wt)
|
|---|
| 329 | { np = 1; nm = 0; nz = 0; }
|
|---|
| 330 | else
|
|---|
| 331 | { np = 0; nm = 1; nz = 0; }
|
|---|
| 332 | }
|
|---|
| 333 | }
|
|---|
| 334 | else
|
|---|
| 335 | {
|
|---|
| 336 | // number of total particles vs. centre of mass Energy - 2*proton mass
|
|---|
| 337 |
|
|---|
| 338 | G4double aleab = std::log(availableEnergy);
|
|---|
| 339 | G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
|
|---|
| 340 | + aleab*(0.117712+0.0136912*aleab))) - 2.0;
|
|---|
| 341 |
|
|---|
| 342 | // normalization constant for kno-distribution.
|
|---|
| 343 | // calculate first the sum of all constants, check for numerical problems.
|
|---|
| 344 | G4double test, dum, anpn = 0.0;
|
|---|
| 345 |
|
|---|
| 346 | for (nt=1; nt<=numSec; nt++) {
|
|---|
| 347 | test = std::exp( Amin( expxu, Amax( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
|
|---|
| 348 | dum = pi*nt/(2.0*n*n);
|
|---|
| 349 | if (std::fabs(dum) < 1.0) {
|
|---|
| 350 | if( test >= 1.0e-10 )anpn += dum*test;
|
|---|
| 351 | } else {
|
|---|
| 352 | anpn += dum*test;
|
|---|
| 353 | }
|
|---|
| 354 | }
|
|---|
| 355 |
|
|---|
| 356 | G4double ran = G4UniformRand();
|
|---|
| 357 | G4double excs = 0.0;
|
|---|
| 358 | if( targetCode == protonCode )
|
|---|
| 359 | {
|
|---|
| 360 | counter = -1;
|
|---|
| 361 | for (np=0; np<numSec/3; np++) {
|
|---|
| 362 | for (nm=Imax(0,np-2); nm<=np; nm++) {
|
|---|
| 363 | for (nz=0; nz<numSec/3; nz++) {
|
|---|
| 364 | if (++counter < numMul) {
|
|---|
| 365 | nt = np+nm+nz;
|
|---|
| 366 | if ( (nt>0) && (nt<=numSec) ) {
|
|---|
| 367 | test = std::exp( Amin( expxu, Amax( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
|
|---|
| 368 | dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
|
|---|
| 369 | if (std::fabs(dum) < 1.0) {
|
|---|
| 370 | if( test >= 1.0e-10 )excs += dum*test;
|
|---|
| 371 | } else {
|
|---|
| 372 | excs += dum*test;
|
|---|
| 373 | }
|
|---|
| 374 | if (ran < excs) goto outOfLoop; //------------------>
|
|---|
| 375 | }
|
|---|
| 376 | }
|
|---|
| 377 | }
|
|---|
| 378 | }
|
|---|
| 379 | }
|
|---|
| 380 |
|
|---|
| 381 | // 3 previous loops continued to the end
|
|---|
| 382 | inElastic = false; // quasi-elastic scattering
|
|---|
| 383 | return;
|
|---|
| 384 | }
|
|---|
| 385 | else
|
|---|
| 386 | { // target must be a neutron
|
|---|
| 387 | counter = -1;
|
|---|
| 388 | for (np=0; np<numSec/3; np++) {
|
|---|
| 389 | for (nm=Imax(0,np-1); nm<=(np+1); nm++) {
|
|---|
| 390 | for (nz=0; nz<numSec/3; nz++) {
|
|---|
| 391 | if (++counter < numMul) {
|
|---|
| 392 | nt = np+nm+nz;
|
|---|
| 393 | if ( (nt>=1) && (nt<=numSec) ) {
|
|---|
| 394 | test = std::exp( Amin( expxu, Amax( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
|
|---|
| 395 | dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
|
|---|
| 396 | if (std::fabs(dum) < 1.0) {
|
|---|
| 397 | if( test >= 1.0e-10 )excs += dum*test;
|
|---|
| 398 | } else {
|
|---|
| 399 | excs += dum*test;
|
|---|
| 400 | }
|
|---|
| 401 | if (ran < excs) goto outOfLoop; // ------------->
|
|---|
| 402 | }
|
|---|
| 403 | }
|
|---|
| 404 | }
|
|---|
| 405 | }
|
|---|
| 406 | }
|
|---|
| 407 | // 3 previous loops continued to the end
|
|---|
| 408 | inElastic = false; // quasi-elastic scattering.
|
|---|
| 409 | return;
|
|---|
| 410 | }
|
|---|
| 411 | }
|
|---|
| 412 | outOfLoop: // <-----------------------------------------------
|
|---|
| 413 |
|
|---|
| 414 | if( targetCode == protonCode)
|
|---|
| 415 | {
|
|---|
| 416 | if( np == nm)
|
|---|
| 417 | {
|
|---|
| 418 | }
|
|---|
| 419 | else if (np == (1+nm))
|
|---|
| 420 | {
|
|---|
| 421 | if( G4UniformRand() < 0.5)
|
|---|
| 422 | {
|
|---|
| 423 | pv[1] = Neutron;
|
|---|
| 424 | }
|
|---|
| 425 | else
|
|---|
| 426 | {
|
|---|
| 427 | pv[0] = Neutron;
|
|---|
| 428 | }
|
|---|
| 429 | }
|
|---|
| 430 | else
|
|---|
| 431 | {
|
|---|
| 432 | pv[0] = Neutron;
|
|---|
| 433 | pv[1] = Neutron;
|
|---|
| 434 | }
|
|---|
| 435 | }
|
|---|
| 436 | else
|
|---|
| 437 | {
|
|---|
| 438 | if( np == nm)
|
|---|
| 439 | {
|
|---|
| 440 | if( G4UniformRand() < 0.25)
|
|---|
| 441 | {
|
|---|
| 442 | pv[0] = Neutron;
|
|---|
| 443 | pv[1] = Proton;
|
|---|
| 444 | }
|
|---|
| 445 | else
|
|---|
| 446 | {
|
|---|
| 447 | }
|
|---|
| 448 | }
|
|---|
| 449 | else if ( np == (1+nm))
|
|---|
| 450 | {
|
|---|
| 451 | pv[0] = Neutron;
|
|---|
| 452 | }
|
|---|
| 453 | else
|
|---|
| 454 | {
|
|---|
| 455 | pv[1] = Proton;
|
|---|
| 456 | }
|
|---|
| 457 | }
|
|---|
| 458 |
|
|---|
| 459 |
|
|---|
| 460 | nt = np + nm + nz;
|
|---|
| 461 | while ( nt > 0)
|
|---|
| 462 | {
|
|---|
| 463 | G4double ran = G4UniformRand();
|
|---|
| 464 | if ( ran < (G4double)np/nt)
|
|---|
| 465 | {
|
|---|
| 466 | if( np > 0 )
|
|---|
| 467 | { pv[vecLen++] = PionPlus;
|
|---|
| 468 | np--;
|
|---|
| 469 | }
|
|---|
| 470 | }
|
|---|
| 471 | else if ( ran < (G4double)(np+nm)/nt)
|
|---|
| 472 | {
|
|---|
| 473 | if( nm > 0 )
|
|---|
| 474 | {
|
|---|
| 475 | pv[vecLen++] = PionMinus;
|
|---|
| 476 | nm--;
|
|---|
| 477 | }
|
|---|
| 478 | }
|
|---|
| 479 | else
|
|---|
| 480 | {
|
|---|
| 481 | if( nz > 0 )
|
|---|
| 482 | {
|
|---|
| 483 | pv[vecLen++] = PionZero;
|
|---|
| 484 | nz--;
|
|---|
| 485 | }
|
|---|
| 486 | }
|
|---|
| 487 | nt = np + nm + nz;
|
|---|
| 488 | }
|
|---|
| 489 | if (verboseLevel > 1)
|
|---|
| 490 | {
|
|---|
| 491 | G4cout << "Particles produced: " ;
|
|---|
| 492 | G4cout << pv[0].getName() << " " ;
|
|---|
| 493 | G4cout << pv[1].getName() << " " ;
|
|---|
| 494 | for (i=2; i < vecLen; i++)
|
|---|
| 495 | {
|
|---|
| 496 | G4cout << pv[i].getName() << " " ;
|
|---|
| 497 | }
|
|---|
| 498 | G4cout << G4endl;
|
|---|
| 499 | }
|
|---|
| 500 | return;
|
|---|
| 501 | }
|
|---|
| 502 |
|
|---|
| 503 |
|
|---|
| 504 |
|
|---|
| 505 |
|
|---|
| 506 |
|
|---|
| 507 |
|
|---|
| 508 |
|
|---|
| 509 |
|
|---|
| 510 |
|
|---|