G4LEProtonInelastic.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // Hadronic Process: Low Energy Proton Inelastic Process
00027 // J.L. Chuma, TRIUMF, 19-Nov-1996
00028  
00029 #include <iostream>
00030 
00031 #include "G4LEProtonInelastic.hh"
00032 #include "G4PhysicalConstants.hh"
00033 #include "G4SystemOfUnits.hh"
00034 #include "Randomize.hh"
00035 
00036 G4LEProtonInelastic::G4LEProtonInelastic(const G4String& name)
00037  :G4InelasticInteraction(name)
00038 {
00039   SetMinEnergy(0.0);
00040   SetMaxEnergy(55.*GeV);
00041   G4cout << "WARNING: model G4LEProtonInelastic is being deprecated and will\n"
00042          << "disappear in Geant4 version 10.0"  << G4endl;  
00043 }
00044 
00045 
00046 void G4LEProtonInelastic::ModelDescription(std::ostream& outFile) const
00047 {
00048   outFile << "G4LEProtonInelastic is one of the Low Energy Parameterized\n"
00049           << "(LEP) models used to implement inelastic proton scattering\n"
00050           << "from nuclei.  It is a re-engineered version of the GHEISHA\n"
00051           << "code of H. Fesefeldt.  It divides the initial collision\n"
00052           << "products into backward- and forward-going clusters which are\n"
00053           << "then decayed into final state hadrons.  The model does not\n"
00054           << "conserve energy on an event-by-event basis.  It may be\n"
00055           << "applied to protons with initial energies between 0 and 25\n"
00056           << "GeV.\n";
00057 }
00058 
00059 
00060 G4HadFinalState*
00061 G4LEProtonInelastic::ApplyYourself(const G4HadProjectile& aTrack,
00062                                    G4Nucleus& targetNucleus)
00063 {
00064   theParticleChange.Clear();
00065   const G4HadProjectile *originalIncident = &aTrack;
00066   if (originalIncident->GetKineticEnergy()<= 0.1*MeV) {
00067     theParticleChange.SetStatusChange(isAlive);
00068     theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy());
00069     theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit()); 
00070     return &theParticleChange;      
00071   }
00072 
00073   // Create the target particle
00074 
00075   G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
00076   if (verboseLevel > 1) {
00077     const G4Material *targetMaterial = aTrack.GetMaterial();
00078     G4cout << "G4LEProtonInelastic::ApplyYourself called" << G4endl;
00079     G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
00080     G4cout << "target material = " << targetMaterial->GetName() << ", ";
00081     G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
00082            << G4endl;
00083   }
00084 
00085   if (originalIncident->GetKineticEnergy()/GeV < 0.01+2.*G4UniformRand()/9.) {
00086     SlowProton(originalIncident, targetNucleus);
00087     if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
00088     delete originalTarget;
00089     return &theParticleChange;
00090   }
00091 
00092   // Fermi motion and evaporation
00093   // As of Geant3, the Fermi energy calculation had not been Done
00094 
00095   G4double ek = originalIncident->GetKineticEnergy()/MeV;
00096   G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
00097   G4ReactionProduct modifiedOriginal;
00098   modifiedOriginal = *originalIncident;
00099     
00100   G4double tkin = targetNucleus.Cinema( ek );
00101   ek += tkin;
00102   modifiedOriginal.SetKineticEnergy( ek*MeV );
00103   G4double et = ek + amas;
00104   G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
00105   G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
00106   if (pp > 0.0) {
00107     G4ThreeVector momentum = modifiedOriginal.GetMomentum();
00108     modifiedOriginal.SetMomentum( momentum * (p/pp) );
00109   }
00110 
00111   // calculate black track energies
00112 
00113   tkin = targetNucleus.EvaporationEffects( ek );
00114   ek -= tkin;
00115   modifiedOriginal.SetKineticEnergy( ek*MeV );
00116   et = ek + amas;
00117   p = std::sqrt( std::abs((et-amas)*(et+amas)) );
00118   pp = modifiedOriginal.GetMomentum().mag()/MeV;
00119   if (pp > 0.0) {
00120     G4ThreeVector momentum = modifiedOriginal.GetMomentum();
00121     modifiedOriginal.SetMomentum( momentum * (p/pp) );
00122   }
00123   const G4double cutOff = 0.1;
00124   if (modifiedOriginal.GetKineticEnergy()/MeV <= cutOff) {
00125     SlowProton( originalIncident, targetNucleus);
00126     if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
00127     delete originalTarget;
00128     return &theParticleChange;
00129   }
00130   G4ReactionProduct currentParticle = modifiedOriginal;
00131   G4ReactionProduct targetParticle;
00132   targetParticle = *originalTarget;
00133   currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
00134   targetParticle.SetSide( -1 );  // target always goes in backward hemisphere
00135   G4bool incidentHasChanged = false;
00136   G4bool targetHasChanged = false;
00137   G4bool quasiElastic = false;
00138   G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec;  // vec will contain the sec. particles
00139   G4int vecLen = 0;
00140   vec.Initialize( 0 );
00141     
00142   Cascade(vec, vecLen,
00143           originalIncident, currentParticle, targetParticle,
00144           incidentHasChanged, targetHasChanged, quasiElastic);
00145     
00146   CalculateMomenta(vec, vecLen,
00147                    originalIncident, originalTarget, modifiedOriginal,
00148                    targetNucleus, currentParticle, targetParticle,
00149                    incidentHasChanged, targetHasChanged, quasiElastic);
00150     
00151   SetUpChange(vec, vecLen,
00152               currentParticle, targetParticle,
00153               incidentHasChanged);
00154 
00155   if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);    
00156   delete originalTarget;
00157   return &theParticleChange;
00158 }
00159 
00160 
00161 void
00162 G4LEProtonInelastic::SlowProton(const G4HadProjectile* originalIncident,
00163                                 G4Nucleus& targetNucleus)
00164 {
00165     const G4double A = targetNucleus.GetA_asInt();    // atomic weight
00166     const G4double Z = targetNucleus.GetZ_asInt();    // atomic number
00167     //
00168     // calculate Q-value of reactions
00169     //
00170     G4double theAtomicMass = targetNucleus.AtomicMass( A, Z );
00171     G4double massVec[9];
00172     massVec[0] = targetNucleus.AtomicMass( A+1.0, Z+1.0 );
00173     massVec[1] = 0.;
00174     if (A > Z+1.0)
00175        massVec[1] = targetNucleus.AtomicMass( A    , Z+1.0 );
00176     massVec[2] = theAtomicMass;
00177     massVec[3] = 0.;
00178     if (A > 1.0 && A-1.0 > Z) 
00179        massVec[3] = targetNucleus.AtomicMass( A-1.0, Z );
00180     massVec[4] = 0.;
00181     if (A > 2.0 && A-2.0 > Z) 
00182        massVec[4] = targetNucleus.AtomicMass( A-2.0, Z     );
00183     massVec[5] = 0.;
00184     if (A > 3.0 && Z > 1.0 && A-3.0 > Z-1.0) 
00185        massVec[5] = targetNucleus.AtomicMass( A-3.0, Z-1.0 );
00186     massVec[6] = 0.;
00187     if (A > 1.0 && A-1.0 > Z+1.0) 
00188        massVec[6] = targetNucleus.AtomicMass( A-1.0, Z+1.0 );
00189     massVec[7] = massVec[3];
00190     massVec[8] = 0.;
00191     if (A > 1.0 && Z > 1.0) 
00192        massVec[8] = targetNucleus.AtomicMass( A-1.0, Z-1.0 );
00193     
00194     G4FastVector<G4ReactionProduct,4> vec;  // vec will contain the secondary particles
00195     G4int vecLen = 0;
00196     vec.Initialize( 0 );
00197     
00198     theReactionDynamics.NuclearReaction( vec, vecLen, originalIncident,
00199                                          targetNucleus, theAtomicMass, massVec );
00200     
00201     theParticleChange.SetStatusChange( stopAndKill );
00202     theParticleChange.SetEnergyChange( 0.0 );
00203     
00204     G4DynamicParticle *pd;
00205     for( G4int i=0; i<vecLen; ++i )
00206     {
00207       pd = new G4DynamicParticle();
00208       pd->SetDefinition( vec[i]->GetDefinition() );
00209       pd->SetMomentum( vec[i]->GetMomentum() );
00210       theParticleChange.AddSecondary( pd );
00211       delete vec[i];
00212     }
00213   }
00214  
00215  void
00216   G4LEProtonInelastic::Cascade(
00217    G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
00218    G4int &vecLen,
00219    const G4HadProjectile *originalIncident,
00220    G4ReactionProduct &currentParticle,
00221    G4ReactionProduct &targetParticle,
00222    G4bool &incidentHasChanged,
00223    G4bool &targetHasChanged,
00224    G4bool &quasiElastic )
00225   {
00226     // derived from original FORTRAN code CASP by H. Fesefeldt (13-Sep-1987)
00227     //
00228     // Proton undergoes interaction with nucleon within a nucleus.  Check if it is
00229     // energetically possible to produce pions/kaons.  In not, assume nuclear excitation
00230     // occurs and input particle is degraded in energy. No other particles are produced.
00231     // If reaction is possible, find the correct number of pions/protons/neutrons
00232     // produced using an interpolation to multiplicity data.  Replace some pions or
00233     // protons/neutrons by kaons or strange baryons according to the average
00234     // multiplicity per Inelastic reaction.
00235     //
00236     // the center of mass energy is based on the initial energy, before
00237     // Fermi motion and evaporation effects are taken into account
00238     //
00239     const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
00240     const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
00241     const G4double targetMass = targetParticle.GetMass()/MeV;
00242     G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
00243                                         targetMass*targetMass +
00244                                         2.0*targetMass*etOriginal );
00245     G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
00246     if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV )
00247     {      // not energetically possible to produce pion(s)
00248       quasiElastic = true;
00249       return;
00250     }
00251     static G4bool first = true;
00252     const G4int numMul = 1200;
00253     const G4int numSec = 60;
00254     static G4double protmul[numMul], protnorm[numSec]; // proton constants
00255     static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
00256     // npos = number of pi+, nneg = number of pi-, nzero = number of pi0
00257     G4int counter, nt=0, npos=0, nneg=0, nzero=0;
00258     const G4double c = 1.25;    
00259     const G4double b[] = { 0.70, 0.35 };
00260     if( first )      // compute normalization constants, this will only be Done once
00261     {
00262       first = false;
00263       G4int i;
00264       for( i=0; i<numMul; ++i )protmul[i] = 0.0;
00265       for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
00266       counter = -1;
00267       for( npos=0; npos<(numSec/3); ++npos )
00268       {
00269         for( nneg=std::max(0,npos-2); nneg<=npos; ++nneg )
00270         {
00271           for( nzero=0; nzero<numSec/3; ++nzero )
00272           {
00273             if( ++counter < numMul )
00274             {
00275               nt = npos+nneg+nzero;
00276               if( nt>0 && nt<=numSec )
00277               {
00278                 protmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[0],c) /
00279                   ( theReactionDynamics.Factorial(2-npos+nneg)*
00280                     theReactionDynamics.Factorial(npos-nneg) );
00281                 protnorm[nt-1] += protmul[counter];
00282               }
00283             }
00284           }
00285         }
00286       }
00287       for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
00288       for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
00289       counter = -1;
00290       for( npos=0; npos<numSec/3; ++npos )
00291       {
00292         for( nneg=std::max(0,npos-1); nneg<=(npos+1); ++nneg )
00293         {
00294           for( nzero=0; nzero<numSec/3; ++nzero )
00295           {
00296             if( ++counter < numMul )
00297             {
00298               nt = npos+nneg+nzero;
00299               if( nt>0 && nt<=numSec )
00300               {
00301                 neutmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[1],c) /
00302                   ( theReactionDynamics.Factorial(1-npos+nneg)*
00303                     theReactionDynamics.Factorial(1+npos-nneg) );
00304                 neutnorm[nt-1] += neutmul[counter];
00305               }
00306             }
00307           }
00308         }
00309       }
00310       for( i=0; i<numSec; ++i )
00311       {
00312         if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
00313         if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
00314       }
00315     }   // end of initialization
00316     
00317     const G4double expxu =  82.;          // upper bound for arg. of exp
00318     const G4double expxl = -expxu;        // lower bound for arg. of exp
00319     G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
00320     G4ParticleDefinition *aProton = G4Proton::Proton();
00321     G4int ieab = static_cast<G4int>(availableEnergy*5.0/GeV);
00322     const G4double supp[] = {0.,0.4,0.55,0.65,0.75,0.82,0.86,0.90,0.94,0.98};
00323     G4double test, w0, wp, wt, wm;
00324     if( (availableEnergy < 2.0*GeV) && (G4UniformRand() >= supp[ieab]) )
00325     {
00326       // suppress high multiplicity events at low momentum
00327       // only one pion will be produced
00328       
00329       npos = nneg = nzero = 0;
00330       if( targetParticle.GetDefinition() == aProton )
00331       {
00332         test = std::exp( std::min( expxu, std::max(
00333          expxl, -(1.0+b[0])*(1.0+b[0])/(2.0*c*c) ) ) );
00334         w0 = test/2.0;
00335         wp = test;
00336         if( G4UniformRand() < w0/(w0+wp) )
00337           nzero = 1;
00338         else
00339           npos = 1;
00340       }
00341       else // target is a neutron
00342       {
00343         test = std::exp( std::min( expxu, std::max(
00344          expxl, -(1.0+b[1])*(1.0+b[1])/(2.0*c*c) ) ) );
00345         w0 = test;
00346         wp = test/2.0;        
00347         test = std::exp( std::min( expxu, std::max(
00348          expxl, -(-1.0+b[1])*(-1.0+b[1])/(2.0*c*c) ) ) );
00349         wm = test/2.0;
00350         wt = w0+wp+wm;
00351         wp += w0;
00352         G4double ran = G4UniformRand();
00353         if( ran < w0/wt )
00354           nzero = 1;
00355         else if( ran < wp/wt )
00356           npos = 1;
00357         else
00358           nneg = 1;
00359       }
00360     }
00361     else // (availableEnergy >= 2.0*GeV) || (random number < supp[ieab])
00362     {
00363       G4double n, anpn;
00364       GetNormalizationConstant( availableEnergy, n, anpn );
00365       G4double ran = G4UniformRand();
00366       G4double dum, excs = 0.0;
00367       if( targetParticle.GetDefinition() == aProton )
00368       {
00369         counter = -1;
00370         for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
00371         {
00372           for( nneg=std::max(0,npos-2); nneg<=npos && ran>=excs; ++nneg )
00373           {
00374             for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
00375             {
00376               if( ++counter < numMul )
00377               {
00378                 nt = npos+nneg+nzero;
00379                 if( nt>0 && nt<=numSec )
00380                 {
00381                   test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
00382                   dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
00383                   if( std::fabs(dum) < 1.0 )
00384                   {
00385                     if( test >= 1.0e-10 )excs += dum*test;
00386                   } else {
00387                     excs += dum*test;
00388                   }
00389                 }
00390               }
00391             }
00392           }
00393         }
00394         if( ran >= excs )  // 3 previous loops continued to the end
00395         {
00396           quasiElastic = true;
00397           return;
00398         }
00399       }
00400       else  // target must be a neutron
00401       {
00402         counter = -1;
00403         for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
00404         {
00405           for( nneg=std::max(0,npos-1); nneg<=(npos+1) && ran>=excs; ++nneg )
00406           {
00407             for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
00408             {
00409               if( ++counter < numMul )
00410               {
00411                 nt = npos+nneg+nzero;
00412                 if( nt>0 && nt<=numSec )
00413                 {
00414                   test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
00415                   dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
00416                   if( std::fabs(dum) < 1.0 )
00417                   {
00418                     if( test >= 1.0e-10 )excs += dum*test;
00419                   } else {
00420                     excs += dum*test;
00421                   }
00422                 }
00423               }
00424             }
00425           }
00426         }
00427         if( ran >= excs )  // 3 previous loops continued to the end
00428         {
00429           quasiElastic = true;
00430           return;
00431         }
00432       }
00433       npos--; nneg--; nzero--;
00434     }
00435     if( targetParticle.GetDefinition() == aProton )
00436     {
00437       switch( npos-nneg )
00438       {
00439        case 1:
00440          if( G4UniformRand() < 0.5 )
00441          {
00442            targetParticle.SetDefinitionAndUpdateE( aNeutron );
00443            targetHasChanged = true;
00444          } else {
00445            currentParticle.SetDefinitionAndUpdateE( aNeutron );
00446            incidentHasChanged = true;
00447          }
00448          break;
00449        case 2:
00450          currentParticle.SetDefinitionAndUpdateE( aNeutron );
00451          targetParticle.SetDefinitionAndUpdateE( aNeutron );
00452          incidentHasChanged = true;
00453          targetHasChanged = true;
00454          break;
00455        default:
00456          break;
00457       }
00458     }
00459     else  // target is a neutron
00460     {
00461       switch( npos-nneg )
00462       {
00463        case 0:
00464          if( G4UniformRand() < 0.333333 )
00465          {
00466            currentParticle.SetDefinitionAndUpdateE( aNeutron );
00467            targetParticle.SetDefinitionAndUpdateE( aProton );
00468            incidentHasChanged = true;
00469            targetHasChanged = true;
00470          }
00471          break;
00472        case 1:
00473          currentParticle.SetDefinitionAndUpdateE( aNeutron );
00474          incidentHasChanged = true;
00475          break;
00476        default:
00477          targetParticle.SetDefinitionAndUpdateE( aProton );
00478          targetHasChanged = true;
00479          break;
00480       }
00481     }
00482     SetUpPions( npos, nneg, nzero, vec, vecLen );
00483     return;
00484   }
00485 
00486  /* end of file */
00487  

Generated on Mon May 27 17:48:46 2013 for Geant4 by  doxygen 1.4.7