| 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 | //
|
|---|
| 28 | // Hadronic Process: Low Energy Neutron Inelastic Process
|
|---|
| 29 | // J.L. Chuma, TRIUMF, 04-Feb-1997
|
|---|
| 30 |
|
|---|
| 31 | #include "G4LENeutronInelastic.hh"
|
|---|
| 32 | #include "Randomize.hh"
|
|---|
| 33 | #include "G4Electron.hh"
|
|---|
| 34 | // #include "DumpFrame.hh"
|
|---|
| 35 |
|
|---|
| 36 | G4HadFinalState *
|
|---|
| 37 | G4LENeutronInelastic::ApplyYourself( const G4HadProjectile &aTrack,
|
|---|
| 38 | G4Nucleus &targetNucleus )
|
|---|
| 39 | {
|
|---|
| 40 | theParticleChange.Clear();
|
|---|
| 41 | const G4HadProjectile *originalIncident = &aTrack;
|
|---|
| 42 | //
|
|---|
| 43 | // create the target particle
|
|---|
| 44 | //
|
|---|
| 45 | G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
|
|---|
| 46 |
|
|---|
| 47 | if( verboseLevel > 1 )
|
|---|
| 48 | {
|
|---|
| 49 | const G4Material *targetMaterial = aTrack.GetMaterial();
|
|---|
| 50 | G4cout << "G4LENeutronInelastic::ApplyYourself called" << G4endl;
|
|---|
| 51 | G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
|
|---|
| 52 | G4cout << "target material = " << targetMaterial->GetName() << ", ";
|
|---|
| 53 | G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
|
|---|
| 54 | << G4endl;
|
|---|
| 55 | }
|
|---|
| 56 | /* not true, for example for Fe56, etc..
|
|---|
| 57 | if( originalIncident->GetKineticEnergy()/MeV < 0.000001 )
|
|---|
| 58 | throw G4HadronicException(__FILE__, __LINE__, "G4LENeutronInelastic: should be capture process!");
|
|---|
| 59 | if( originalIncident->Get4Momentum().vect().mag()/MeV < 0.000001 )
|
|---|
| 60 | throw G4HadronicException(__FILE__, __LINE__, "G4LENeutronInelastic: should be capture process!");
|
|---|
| 61 | */
|
|---|
| 62 |
|
|---|
| 63 | G4ReactionProduct modifiedOriginal;
|
|---|
| 64 | modifiedOriginal = *originalIncident;
|
|---|
| 65 | G4ReactionProduct targetParticle;
|
|---|
| 66 | targetParticle = *originalTarget;
|
|---|
| 67 | if( originalIncident->GetKineticEnergy()/GeV < 0.01 + 2.*G4UniformRand()/9. )
|
|---|
| 68 | {
|
|---|
| 69 | SlowNeutron( originalIncident, modifiedOriginal, targetParticle, targetNucleus );
|
|---|
| 70 | delete originalTarget;
|
|---|
| 71 | return &theParticleChange;
|
|---|
| 72 | }
|
|---|
| 73 | //
|
|---|
| 74 | // Fermi motion and evaporation
|
|---|
| 75 | // As of Geant3, the Fermi energy calculation had not been Done
|
|---|
| 76 | //
|
|---|
| 77 | G4double ek = originalIncident->GetKineticEnergy()/MeV;
|
|---|
| 78 | G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
|
|---|
| 79 |
|
|---|
| 80 | G4double tkin = targetNucleus.Cinema( ek );
|
|---|
| 81 | ek += tkin;
|
|---|
| 82 | modifiedOriginal.SetKineticEnergy( ek*MeV );
|
|---|
| 83 | G4double et = ek + amas;
|
|---|
| 84 | G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
|
|---|
| 85 | G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
|
|---|
| 86 | if( pp > 0.0 )
|
|---|
| 87 | {
|
|---|
| 88 | G4ThreeVector momentum = modifiedOriginal.GetMomentum();
|
|---|
| 89 | modifiedOriginal.SetMomentum( momentum * (p/pp) );
|
|---|
| 90 | }
|
|---|
| 91 | //
|
|---|
| 92 | // calculate black track energies
|
|---|
| 93 | //
|
|---|
| 94 | tkin = targetNucleus.EvaporationEffects( ek );
|
|---|
| 95 | ek -= tkin;
|
|---|
| 96 | modifiedOriginal.SetKineticEnergy( ek*MeV );
|
|---|
| 97 | et = ek + amas;
|
|---|
| 98 | p = std::sqrt( std::abs((et-amas)*(et+amas)) );
|
|---|
| 99 | pp = modifiedOriginal.GetMomentum().mag()/MeV;
|
|---|
| 100 | if( pp > 0.0 )
|
|---|
| 101 | {
|
|---|
| 102 | G4ThreeVector momentum = modifiedOriginal.GetMomentum();
|
|---|
| 103 | modifiedOriginal.SetMomentum( momentum * (p/pp) );
|
|---|
| 104 | }
|
|---|
| 105 | const G4double cutOff = 0.1;
|
|---|
| 106 | if( modifiedOriginal.GetKineticEnergy()/MeV <= cutOff )
|
|---|
| 107 | {
|
|---|
| 108 | SlowNeutron( originalIncident, modifiedOriginal, targetParticle, targetNucleus );
|
|---|
| 109 | delete originalTarget;
|
|---|
| 110 | return &theParticleChange;
|
|---|
| 111 | }
|
|---|
| 112 | G4ReactionProduct currentParticle = modifiedOriginal;
|
|---|
| 113 | currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
|
|---|
| 114 | targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
|
|---|
| 115 | G4bool incidentHasChanged = false;
|
|---|
| 116 | G4bool targetHasChanged = false;
|
|---|
| 117 | G4bool quasiElastic = false;
|
|---|
| 118 | G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
|
|---|
| 119 | G4int vecLen = 0;
|
|---|
| 120 | vec.Initialize( 0 );
|
|---|
| 121 |
|
|---|
| 122 | Cascade( vec, vecLen,
|
|---|
| 123 | originalIncident, currentParticle, targetParticle,
|
|---|
| 124 | incidentHasChanged, targetHasChanged, quasiElastic );
|
|---|
| 125 |
|
|---|
| 126 | CalculateMomenta( vec, vecLen,
|
|---|
| 127 | originalIncident, originalTarget, modifiedOriginal,
|
|---|
| 128 | targetNucleus, currentParticle, targetParticle,
|
|---|
| 129 | incidentHasChanged, targetHasChanged, quasiElastic );
|
|---|
| 130 |
|
|---|
| 131 | SetUpChange( vec, vecLen,
|
|---|
| 132 | currentParticle, targetParticle,
|
|---|
| 133 | incidentHasChanged );
|
|---|
| 134 |
|
|---|
| 135 | delete originalTarget;
|
|---|
| 136 | return &theParticleChange;
|
|---|
| 137 | }
|
|---|
| 138 |
|
|---|
| 139 | void
|
|---|
| 140 | G4LENeutronInelastic::SlowNeutron(
|
|---|
| 141 | const G4HadProjectile *originalIncident,
|
|---|
| 142 | G4ReactionProduct &modifiedOriginal,
|
|---|
| 143 | G4ReactionProduct &targetParticle,
|
|---|
| 144 | G4Nucleus &targetNucleus )
|
|---|
| 145 | {
|
|---|
| 146 | const G4double A = targetNucleus.GetN(); // atomic weight
|
|---|
| 147 | const G4double Z = targetNucleus.GetZ(); // atomic number
|
|---|
| 148 |
|
|---|
| 149 | G4double currentKinetic = modifiedOriginal.GetKineticEnergy()/MeV;
|
|---|
| 150 | G4double currentMass = modifiedOriginal.GetMass()/MeV;
|
|---|
| 151 | if( A < 1.5 ) // Hydrogen
|
|---|
| 152 | {
|
|---|
| 153 | //
|
|---|
| 154 | // very simple simulation of scattering angle and energy
|
|---|
| 155 | // nonrelativistic approximation with isotropic angular
|
|---|
| 156 | // distribution in the cms system
|
|---|
| 157 | //
|
|---|
| 158 | G4double cost1, eka = 0.0;
|
|---|
| 159 | while (eka <= 0.0)
|
|---|
| 160 | {
|
|---|
| 161 | cost1 = -1.0 + 2.0*G4UniformRand();
|
|---|
| 162 | eka = 1.0 + 2.0*cost1*A + A*A;
|
|---|
| 163 | }
|
|---|
| 164 | G4double cost = std::min( 1.0, std::max( -1.0, (A*cost1+1.0)/std::sqrt(eka) ) );
|
|---|
| 165 | eka /= (1.0+A)*(1.0+A);
|
|---|
| 166 | G4double ek = currentKinetic*MeV/GeV;
|
|---|
| 167 | G4double amas = currentMass*MeV/GeV;
|
|---|
| 168 | ek *= eka;
|
|---|
| 169 | G4double en = ek + amas;
|
|---|
| 170 | G4double p = std::sqrt(std::abs(en*en-amas*amas));
|
|---|
| 171 | G4double sint = std::sqrt(std::abs(1.0-cost*cost));
|
|---|
| 172 | G4double phi = G4UniformRand()*twopi;
|
|---|
| 173 | G4double px = sint*std::sin(phi);
|
|---|
| 174 | G4double py = sint*std::cos(phi);
|
|---|
| 175 | G4double pz = cost;
|
|---|
| 176 | targetParticle.SetMomentum( px*GeV, py*GeV, pz*GeV );
|
|---|
| 177 | G4double pxO = originalIncident->Get4Momentum().x()/GeV;
|
|---|
| 178 | G4double pyO = originalIncident->Get4Momentum().y()/GeV;
|
|---|
| 179 | G4double pzO = originalIncident->Get4Momentum().z()/GeV;
|
|---|
| 180 | G4double ptO = pxO*pxO + pyO+pyO;
|
|---|
| 181 | if( ptO > 0.0 )
|
|---|
| 182 | {
|
|---|
| 183 | G4double pO = std::sqrt(pxO*pxO+pyO*pyO+pzO*pzO);
|
|---|
| 184 | cost = pzO/pO;
|
|---|
| 185 | sint = 0.5*(std::sqrt(std::abs((1.0-cost)*(1.0+cost)))+std::sqrt(ptO)/pO);
|
|---|
| 186 | G4double ph = pi/2.0;
|
|---|
| 187 | if( pyO < 0.0 )ph = ph*1.5;
|
|---|
| 188 | if( std::abs(pxO) > 0.000001 )ph = std::atan2(pyO,pxO);
|
|---|
| 189 | G4double cosp = std::cos(ph);
|
|---|
| 190 | G4double sinp = std::sin(ph);
|
|---|
| 191 | px = cost*cosp*px - sinp*py+sint*cosp*pz;
|
|---|
| 192 | py = cost*sinp*px + cosp*py+sint*sinp*pz;
|
|---|
| 193 | pz = -sint*px + cost*pz;
|
|---|
| 194 | }
|
|---|
| 195 | else
|
|---|
| 196 | {
|
|---|
| 197 | if( pz < 0.0 )pz *= -1.0;
|
|---|
| 198 | }
|
|---|
| 199 | G4double pu = std::sqrt(px*px+py*py+pz*pz);
|
|---|
| 200 | modifiedOriginal.SetMomentum( targetParticle.GetMomentum() * (p/pu) );
|
|---|
| 201 | modifiedOriginal.SetKineticEnergy( ek*GeV );
|
|---|
| 202 |
|
|---|
| 203 | targetParticle.SetMomentum(
|
|---|
| 204 | originalIncident->Get4Momentum().vect() - modifiedOriginal.GetMomentum() );
|
|---|
| 205 | G4double pp = targetParticle.GetMomentum().mag();
|
|---|
| 206 | G4double tarmas = targetParticle.GetMass();
|
|---|
| 207 | targetParticle.SetTotalEnergy( std::sqrt( pp*pp + tarmas*tarmas ) );
|
|---|
| 208 |
|
|---|
| 209 | theParticleChange.SetEnergyChange( modifiedOriginal.GetKineticEnergy() );
|
|---|
| 210 | G4DynamicParticle *pd = new G4DynamicParticle;
|
|---|
| 211 | pd->SetDefinition( targetParticle.GetDefinition() );
|
|---|
| 212 | pd->SetMomentum( targetParticle.GetMomentum() );
|
|---|
| 213 | theParticleChange.AddSecondary( pd );
|
|---|
| 214 | return;
|
|---|
| 215 | }
|
|---|
| 216 | G4FastVector<G4ReactionProduct,4> vec; // vec will contain the secondary particles
|
|---|
| 217 | G4int vecLen = 0;
|
|---|
| 218 | vec.Initialize( 0 );
|
|---|
| 219 |
|
|---|
| 220 | G4double theAtomicMass = targetNucleus.AtomicMass( A, Z );
|
|---|
| 221 | G4double massVec[9];
|
|---|
| 222 | massVec[0] = targetNucleus.AtomicMass( A+1.0, Z );
|
|---|
| 223 | massVec[1] = theAtomicMass;
|
|---|
| 224 | massVec[2] = 0.;
|
|---|
| 225 | if (Z > 1.0)
|
|---|
| 226 | massVec[2] = targetNucleus.AtomicMass( A , Z-1.0 );
|
|---|
| 227 | massVec[3] = 0.;
|
|---|
| 228 | if (Z > 1.0 && A > 1.0)
|
|---|
| 229 | massVec[3] = targetNucleus.AtomicMass( A-1.0, Z-1.0 );
|
|---|
| 230 | massVec[4] = 0.;
|
|---|
| 231 | if (Z > 1.0 && A > 2.0 && A-2.0 > Z-1.0)
|
|---|
| 232 | massVec[4] = targetNucleus.AtomicMass( A-2.0, Z-1.0 );
|
|---|
| 233 | massVec[5] = 0.;
|
|---|
| 234 | if (Z > 2.0 && A > 3.0 && A-3.0 > Z-2.0)
|
|---|
| 235 | massVec[5] = targetNucleus.AtomicMass( A-3.0, Z-2.0 );
|
|---|
| 236 | massVec[6] = 0.;
|
|---|
| 237 | if (A > 1.0 && A-1.0 > Z)
|
|---|
| 238 | massVec[6] = targetNucleus.AtomicMass( A-1.0, Z );
|
|---|
| 239 | massVec[7] = massVec[3];
|
|---|
| 240 | massVec[8] = 0.;
|
|---|
| 241 | if (Z > 2.0 && A > 1.0)
|
|---|
| 242 | massVec[8] = targetNucleus.AtomicMass( A-1.0, Z-2.0 );
|
|---|
| 243 |
|
|---|
| 244 | theReactionDynamics.NuclearReaction( vec, vecLen, originalIncident,
|
|---|
| 245 | targetNucleus, theAtomicMass, massVec );
|
|---|
| 246 |
|
|---|
| 247 | theParticleChange.SetStatusChange( stopAndKill );
|
|---|
| 248 | theParticleChange.SetEnergyChange( 0.0 );
|
|---|
| 249 |
|
|---|
| 250 | G4DynamicParticle * pd;
|
|---|
| 251 | for( G4int i=0; i<vecLen; ++i )
|
|---|
| 252 | {
|
|---|
| 253 | pd = new G4DynamicParticle();
|
|---|
| 254 | pd->SetDefinition( vec[i]->GetDefinition() );
|
|---|
| 255 | pd->SetMomentum( vec[i]->GetMomentum() );
|
|---|
| 256 | theParticleChange.AddSecondary( pd );
|
|---|
| 257 | delete vec[i];
|
|---|
| 258 | }
|
|---|
| 259 | }
|
|---|
| 260 |
|
|---|
| 261 | void
|
|---|
| 262 | G4LENeutronInelastic::Cascade(
|
|---|
| 263 | G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
|
|---|
| 264 | G4int& vecLen,
|
|---|
| 265 | const G4HadProjectile *originalIncident,
|
|---|
| 266 | G4ReactionProduct ¤tParticle,
|
|---|
| 267 | G4ReactionProduct &targetParticle,
|
|---|
| 268 | G4bool &incidentHasChanged,
|
|---|
| 269 | G4bool &targetHasChanged,
|
|---|
| 270 | G4bool &quasiElastic )
|
|---|
| 271 | {
|
|---|
| 272 | // derived from original FORTRAN code CASN by H. Fesefeldt (13-Sep-1987)
|
|---|
| 273 | //
|
|---|
| 274 | // Neutron undergoes interaction with nucleon within a nucleus. Check if it is
|
|---|
| 275 | // energetically possible to produce pions/kaons. In not, assume nuclear excitation
|
|---|
| 276 | // occurs and input particle is degraded in energy. No other particles are produced.
|
|---|
| 277 | // If reaction is possible, find the correct number of pions/protons/neutrons
|
|---|
| 278 | // produced using an interpolation to multiplicity data. Replace some pions or
|
|---|
| 279 | // protons/neutrons by kaons or strange baryons according to the average
|
|---|
| 280 | // multiplicity per Inelastic reaction.
|
|---|
| 281 | //
|
|---|
| 282 | const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
|
|---|
| 283 | const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
|
|---|
| 284 | const G4double targetMass = targetParticle.GetMass()/MeV;
|
|---|
| 285 | G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
|
|---|
| 286 | targetMass*targetMass +
|
|---|
| 287 | 2.0*targetMass*etOriginal );
|
|---|
| 288 | G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
|
|---|
| 289 | if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV )
|
|---|
| 290 | {
|
|---|
| 291 | quasiElastic = true;
|
|---|
| 292 | return;
|
|---|
| 293 | }
|
|---|
| 294 | static G4bool first = true;
|
|---|
| 295 | const G4int numMul = 1200;
|
|---|
| 296 | const G4int numSec = 60;
|
|---|
| 297 | static G4double protmul[numMul], protnorm[numSec]; // proton constants
|
|---|
| 298 | static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
|
|---|
| 299 | // np = number of pi+, nm = number of pi-, nz = number of pi0
|
|---|
| 300 | G4int counter, nt=0, np=0, nm=0, nz=0;
|
|---|
| 301 | const G4double c = 1.25;
|
|---|
| 302 | const G4double b[] = { 0.35, 0.0 };
|
|---|
| 303 | if( first ) // compute normalization constants, this will only be Done once
|
|---|
| 304 | {
|
|---|
| 305 | first = false;
|
|---|
| 306 | G4int i;
|
|---|
| 307 | for( i=0; i<numMul; ++i )protmul[i] = 0.0;
|
|---|
| 308 | for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
|
|---|
| 309 | counter = -1;
|
|---|
| 310 | for( np=0; np<numSec/3; ++np )
|
|---|
| 311 | {
|
|---|
| 312 | for( nm=std::max(0,np-1); nm<=(np+1); ++nm )
|
|---|
| 313 | {
|
|---|
| 314 | for( nz=0; nz<numSec/3; ++nz )
|
|---|
| 315 | {
|
|---|
| 316 | if( ++counter < numMul )
|
|---|
| 317 | {
|
|---|
| 318 | nt = np+nm+nz;
|
|---|
| 319 | if( nt > 0 )
|
|---|
| 320 | {
|
|---|
| 321 | protmul[counter] = Pmltpc(np,nm,nz,nt,b[0],c) /
|
|---|
| 322 | ( theReactionDynamics.Factorial(1-np+nm)*
|
|---|
| 323 | theReactionDynamics.Factorial(1+np-nm) );
|
|---|
| 324 | protnorm[nt-1] += protmul[counter];
|
|---|
| 325 | }
|
|---|
| 326 | }
|
|---|
| 327 | }
|
|---|
| 328 | }
|
|---|
| 329 | }
|
|---|
| 330 | for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
|
|---|
| 331 | for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
|
|---|
| 332 | counter = -1;
|
|---|
| 333 | for( np=0; np<(numSec/3); ++np )
|
|---|
| 334 | {
|
|---|
| 335 | for( nm=np; nm<=(np+2); ++nm )
|
|---|
| 336 | {
|
|---|
| 337 | for( nz=0; nz<numSec/3; ++nz )
|
|---|
| 338 | {
|
|---|
| 339 | if( ++counter < numMul )
|
|---|
| 340 | {
|
|---|
| 341 | nt = np+nm+nz;
|
|---|
| 342 | if( (nt>0) && (nt<=numSec) )
|
|---|
| 343 | {
|
|---|
| 344 | neutmul[counter] = Pmltpc(np,nm,nz,nt,b[1],c) /
|
|---|
| 345 | ( theReactionDynamics.Factorial(nm-np)*
|
|---|
| 346 | theReactionDynamics.Factorial(2-nm+np) );
|
|---|
| 347 | neutnorm[nt-1] += neutmul[counter];
|
|---|
| 348 | }
|
|---|
| 349 | }
|
|---|
| 350 | }
|
|---|
| 351 | }
|
|---|
| 352 | }
|
|---|
| 353 | for( i=0; i<numSec; ++i )
|
|---|
| 354 | {
|
|---|
| 355 | if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
|
|---|
| 356 | if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
|
|---|
| 357 | }
|
|---|
| 358 | } // end of initialization
|
|---|
| 359 |
|
|---|
| 360 | const G4double expxu = 82.; // upper bound for arg. of exp
|
|---|
| 361 | const G4double expxl = -expxu; // lower bound for arg. of exp
|
|---|
| 362 | G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
|
|---|
| 363 | G4ParticleDefinition *aProton = G4Proton::Proton();
|
|---|
| 364 | G4int ieab = static_cast<G4int>(availableEnergy*5.0/GeV);
|
|---|
| 365 | const G4double supp[] = {0.,0.4,0.55,0.65,0.75,0.82,0.86,0.90,0.94,0.98};
|
|---|
| 366 | G4double test, w0, wp, wt, wm;
|
|---|
| 367 | if( (availableEnergy < 2.0*GeV) && (G4UniformRand() >= supp[ieab]) )
|
|---|
| 368 | {
|
|---|
| 369 | // suppress high multiplicity events at low momentum
|
|---|
| 370 | // only one pion will be produced
|
|---|
| 371 |
|
|---|
| 372 | nm = np = nz = 0;
|
|---|
| 373 | if( targetParticle.GetDefinition() == aNeutron )
|
|---|
| 374 | {
|
|---|
| 375 | test = std::exp( std::min( expxu, std::max( expxl, -(1.0+b[1])*(1.0+b[1])/(2.0*c*c) ) ) );
|
|---|
| 376 | w0 = test/2.0;
|
|---|
| 377 | wm = test;
|
|---|
| 378 | if( G4UniformRand() < w0/(w0+wm) )
|
|---|
| 379 | nz = 1;
|
|---|
| 380 | else
|
|---|
| 381 | nm = 1;
|
|---|
| 382 | } else { // target is a proton
|
|---|
| 383 | test = std::exp( std::min( expxu, std::max( expxl, -(1.0+b[0])*(1.0+b[0])/(2.0*c*c) ) ) );
|
|---|
| 384 | w0 = test;
|
|---|
| 385 | wp = test/2.0;
|
|---|
| 386 | test = std::exp( std::min( expxu, std::max( expxl, -(-1.0+b[0])*(-1.0+b[0])/(2.0*c*c) ) ) );
|
|---|
| 387 | wm = test/2.0;
|
|---|
| 388 | wt = w0+wp+wm;
|
|---|
| 389 | wp += w0;
|
|---|
| 390 | G4double ran = G4UniformRand();
|
|---|
| 391 | if( ran < w0/wt )
|
|---|
| 392 | nz = 1;
|
|---|
| 393 | else if( ran < wp/wt )
|
|---|
| 394 | np = 1;
|
|---|
| 395 | else
|
|---|
| 396 | nm = 1;
|
|---|
| 397 | }
|
|---|
| 398 | } else { // (availableEnergy >= 2.0*GeV) || (random number < supp[ieab])
|
|---|
| 399 | G4double n, anpn;
|
|---|
| 400 | GetNormalizationConstant( availableEnergy, n, anpn );
|
|---|
| 401 | G4double ran = G4UniformRand();
|
|---|
| 402 | G4double dum, excs = 0.0;
|
|---|
| 403 | if( targetParticle.GetDefinition() == aProton )
|
|---|
| 404 | {
|
|---|
| 405 | counter = -1;
|
|---|
| 406 | for( np=0; np<numSec/3 && ran>=excs; ++np )
|
|---|
| 407 | {
|
|---|
| 408 | for( nm=std::max(0,np-1); nm<=(np+1) && ran>=excs; ++nm )
|
|---|
| 409 | {
|
|---|
| 410 | for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
|
|---|
| 411 | {
|
|---|
| 412 | if( ++counter < numMul )
|
|---|
| 413 | {
|
|---|
| 414 | nt = np+nm+nz;
|
|---|
| 415 | if( nt > 0 )
|
|---|
| 416 | {
|
|---|
| 417 | test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
|
|---|
| 418 | dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
|
|---|
| 419 | if( std::fabs(dum) < 1.0 ) {
|
|---|
| 420 | if( test >= 1.0e-10 )excs += dum*test;
|
|---|
| 421 | } else {
|
|---|
| 422 | excs += dum*test;
|
|---|
| 423 | }
|
|---|
| 424 | }
|
|---|
| 425 | }
|
|---|
| 426 | }
|
|---|
| 427 | }
|
|---|
| 428 | }
|
|---|
| 429 | if( ran >= excs ) // 3 previous loops continued to the end
|
|---|
| 430 | {
|
|---|
| 431 | quasiElastic = true;
|
|---|
| 432 | return;
|
|---|
| 433 | }
|
|---|
| 434 | np--; nm--; nz--;
|
|---|
| 435 | } else { // target must be a neutron
|
|---|
| 436 | counter = -1;
|
|---|
| 437 | for( np=0; np<numSec/3 && ran>=excs; ++np )
|
|---|
| 438 | {
|
|---|
| 439 | for( nm=np; nm<=(np+2) && ran>=excs; ++nm )
|
|---|
| 440 | {
|
|---|
| 441 | for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
|
|---|
| 442 | {
|
|---|
| 443 | if( ++counter < numMul )
|
|---|
| 444 | {
|
|---|
| 445 | nt = np+nm+nz;
|
|---|
| 446 | if( (nt>=1) && (nt<=numSec) )
|
|---|
| 447 | {
|
|---|
| 448 | test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
|
|---|
| 449 | dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
|
|---|
| 450 | if( std::fabs(dum) < 1.0 ) {
|
|---|
| 451 | if( test >= 1.0e-10 )excs += dum*test;
|
|---|
| 452 | } else {
|
|---|
| 453 | excs += dum*test;
|
|---|
| 454 | }
|
|---|
| 455 | }
|
|---|
| 456 | }
|
|---|
| 457 | }
|
|---|
| 458 | }
|
|---|
| 459 | }
|
|---|
| 460 | if( ran >= excs ) // 3 previous loops continued to the end
|
|---|
| 461 | {
|
|---|
| 462 | quasiElastic = true;
|
|---|
| 463 | return;
|
|---|
| 464 | }
|
|---|
| 465 | np--; nm--; nz--;
|
|---|
| 466 | }
|
|---|
| 467 | }
|
|---|
| 468 | if( targetParticle.GetDefinition() == aProton )
|
|---|
| 469 | {
|
|---|
| 470 | switch( np-nm )
|
|---|
| 471 | {
|
|---|
| 472 | case 0:
|
|---|
| 473 | if( G4UniformRand() < 0.33 )
|
|---|
| 474 | {
|
|---|
| 475 | currentParticle.SetDefinitionAndUpdateE( aProton );
|
|---|
| 476 | targetParticle.SetDefinitionAndUpdateE( aNeutron );
|
|---|
| 477 | incidentHasChanged = true;
|
|---|
| 478 | targetHasChanged = true;
|
|---|
| 479 | }
|
|---|
| 480 | break;
|
|---|
| 481 | case 1:
|
|---|
| 482 | targetParticle.SetDefinitionAndUpdateE( aNeutron );
|
|---|
| 483 | targetHasChanged = true;
|
|---|
| 484 | break;
|
|---|
| 485 | default:
|
|---|
| 486 | currentParticle.SetDefinitionAndUpdateE( aProton );
|
|---|
| 487 | incidentHasChanged = true;
|
|---|
| 488 | break;
|
|---|
| 489 | }
|
|---|
| 490 | } else { // target must be a neutron
|
|---|
| 491 | switch( np-nm )
|
|---|
| 492 | {
|
|---|
| 493 | case -1: // changed from +1 by JLC, 7Jul97
|
|---|
| 494 | if( G4UniformRand() < 0.5 )
|
|---|
| 495 | {
|
|---|
| 496 | currentParticle.SetDefinitionAndUpdateE( aProton );
|
|---|
| 497 | incidentHasChanged = true;
|
|---|
| 498 | } else {
|
|---|
| 499 | targetParticle.SetDefinitionAndUpdateE( aProton );
|
|---|
| 500 | targetHasChanged = true;
|
|---|
| 501 | }
|
|---|
| 502 | break;
|
|---|
| 503 | case 0:
|
|---|
| 504 | break;
|
|---|
| 505 | default:
|
|---|
| 506 | currentParticle.SetDefinitionAndUpdateE( aProton );
|
|---|
| 507 | targetParticle.SetDefinitionAndUpdateE( aProton );
|
|---|
| 508 | incidentHasChanged = true;
|
|---|
| 509 | targetHasChanged = true;
|
|---|
| 510 | break;
|
|---|
| 511 | }
|
|---|
| 512 | }
|
|---|
| 513 | SetUpPions( np, nm, nz, vec, vecLen );
|
|---|
| 514 | // DEBUG --> DumpFrames::DumpFrame(vec, vecLen);
|
|---|
| 515 | return;
|
|---|
| 516 | }
|
|---|
| 517 |
|
|---|
| 518 | /* end of file */
|
|---|
| 519 |
|
|---|