G4RPGSigmaMinusInelastic.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 // $Id$
00027 //
00028  
00029 #include "G4RPGSigmaMinusInelastic.hh"
00030 #include "G4PhysicalConstants.hh"
00031 #include "G4SystemOfUnits.hh"
00032 #include "Randomize.hh"
00033  
00034 G4HadFinalState*
00035 G4RPGSigmaMinusInelastic::ApplyYourself( const G4HadProjectile &aTrack,
00036                                           G4Nucleus &targetNucleus )
00037 {
00038   const G4HadProjectile *originalIncident = &aTrack;
00039   if (originalIncident->GetKineticEnergy()<= 0.1*MeV) 
00040   {
00041     theParticleChange.SetStatusChange(isAlive);
00042     theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy());
00043     theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit()); 
00044     return &theParticleChange;      
00045   }
00046 
00047   // create the target particle
00048 
00049   G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
00050     
00051   if( verboseLevel > 1 )
00052   {
00053     const G4Material *targetMaterial = aTrack.GetMaterial();
00054     G4cout << "G4RPGSigmaMinusInelastic::ApplyYourself called" << G4endl;
00055     G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
00056     G4cout << "target material = " << targetMaterial->GetName() << ", ";
00057     G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
00058            << G4endl;
00059   }
00060 
00061   // Fermi motion and evaporation
00062   // As of Geant3, the Fermi energy calculation had not been Done
00063 
00064     G4double ek = originalIncident->GetKineticEnergy()/MeV;
00065     G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
00066     G4ReactionProduct modifiedOriginal;
00067     modifiedOriginal = *originalIncident;
00068     
00069     G4double tkin = targetNucleus.Cinema( ek );
00070     ek += tkin;
00071     modifiedOriginal.SetKineticEnergy( ek*MeV );
00072     G4double et = ek + amas;
00073     G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
00074     G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
00075     if( pp > 0.0 )
00076     {
00077       G4ThreeVector momentum = modifiedOriginal.GetMomentum();
00078       modifiedOriginal.SetMomentum( momentum * (p/pp) );
00079     }
00080     //
00081     // calculate black track energies
00082     //
00083     tkin = targetNucleus.EvaporationEffects( ek );
00084     ek -= tkin;
00085     modifiedOriginal.SetKineticEnergy( ek*MeV );
00086     et = ek + amas;
00087     p = std::sqrt( std::abs((et-amas)*(et+amas)) );
00088     pp = modifiedOriginal.GetMomentum().mag()/MeV;
00089     if( pp > 0.0 )
00090     {
00091       G4ThreeVector momentum = modifiedOriginal.GetMomentum();
00092       modifiedOriginal.SetMomentum( momentum * (p/pp) );
00093     }
00094     G4ReactionProduct currentParticle = modifiedOriginal;
00095     G4ReactionProduct targetParticle;
00096     targetParticle = *originalTarget;
00097     currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
00098     targetParticle.SetSide( -1 );  // target always goes in backward hemisphere
00099     G4bool incidentHasChanged = false;
00100     G4bool targetHasChanged = false;
00101     G4bool quasiElastic = false;
00102     G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec;  // vec will contain the secondary particles
00103     G4int vecLen = 0;
00104     vec.Initialize( 0 );
00105     
00106     const G4double cutOff = 0.1;
00107     if( originalIncident->GetKineticEnergy()/MeV > cutOff )
00108       Cascade( vec, vecLen,
00109                originalIncident, currentParticle, targetParticle,
00110                incidentHasChanged, targetHasChanged, quasiElastic );
00111     
00112     CalculateMomenta( vec, vecLen,
00113                       originalIncident, originalTarget, modifiedOriginal,
00114                       targetNucleus, currentParticle, targetParticle,
00115                       incidentHasChanged, targetHasChanged, quasiElastic );
00116     
00117     SetUpChange( vec, vecLen,
00118                  currentParticle, targetParticle,
00119                  incidentHasChanged );
00120     
00121   delete originalTarget;
00122   return &theParticleChange;
00123 }
00124 
00125 
00126 void G4RPGSigmaMinusInelastic::Cascade(
00127    G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
00128    G4int& vecLen,
00129    const G4HadProjectile *originalIncident,
00130    G4ReactionProduct &currentParticle,
00131    G4ReactionProduct &targetParticle,
00132    G4bool &incidentHasChanged,
00133    G4bool &targetHasChanged,
00134    G4bool &quasiElastic )
00135 {
00136   // Derived from H. Fesefeldt's original FORTRAN code CASSM
00137   //
00138   // SigmaMinus undergoes interaction with nucleon within a nucleus.  Check if it is
00139   // energetically possible to produce pions/kaons.  In not, assume nuclear excitation
00140   // occurs and input particle is degraded in energy. No other particles are produced.
00141   // If reaction is possible, find the correct number of pions/protons/neutrons
00142   // produced using an interpolation to multiplicity data.  Replace some pions or
00143   // protons/neutrons by kaons or strange baryons according to the average
00144   // multiplicity per Inelastic reaction.
00145 
00146   const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
00147   const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
00148   const G4double targetMass = targetParticle.GetMass()/MeV;
00149   G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
00150                                       targetMass*targetMass +
00151                                       2.0*targetMass*etOriginal );
00152   G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
00153   if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV )
00154   {
00155     quasiElastic = true;
00156     return;
00157   }
00158     static G4bool first = true;
00159     const G4int numMul = 1200;
00160     const G4int numSec = 60;
00161     static G4double protmul[numMul], protnorm[numSec]; // proton constants
00162     static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
00163     // np = number of pi+, nneg = number of pi-, nz = number of pi0
00164     G4int counter, nt=0, np=0, nneg=0, nz=0;
00165     G4double test;
00166     const G4double c = 1.25;    
00167     const G4double b[] = { 0.70, 0.70 };
00168     if( first )     // compute normalization constants, this will only be Done once
00169     {
00170       first = false;
00171       G4int i;
00172       for( i=0; i<numMul; ++i )protmul[i] = 0.0;
00173       for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
00174       counter = -1;
00175       for( np=0; np<(numSec/3); ++np )
00176       {
00177         for( nneg=std::max(0,np-1); nneg<=(np+1); ++nneg )
00178         {
00179           for( nz=0; nz<numSec/3; ++nz )
00180           {
00181             if( ++counter < numMul )
00182             {
00183               nt = np+nneg+nz;
00184               if( nt>0 && nt<=numSec )
00185               {
00186                 protmul[counter] = Pmltpc(np,nneg,nz,nt,b[0],c);
00187                 protnorm[nt-1] += protmul[counter];
00188               }
00189             }
00190           }
00191         }
00192       }
00193       for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
00194       for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
00195       counter = -1;
00196       for( np=0; np<numSec/3; ++np )
00197       {
00198         for( nneg=np; nneg<=(np+2); ++nneg )
00199         {
00200           for( nz=0; nz<numSec/3; ++nz )
00201           {
00202             if( ++counter < numMul )
00203             {
00204               nt = np+nneg+nz;
00205               if( nt>0 && nt<=numSec )
00206               {
00207                 neutmul[counter] = Pmltpc(np,nneg,nz,nt,b[1],c);
00208                 neutnorm[nt-1] += neutmul[counter];
00209               }
00210             }
00211           }
00212         }
00213       }
00214       for( i=0; i<numSec; ++i )
00215       {
00216         if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
00217         if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
00218       }
00219     }   // end of initialization
00220     
00221     const G4double expxu = 82.;           // upper bound for arg. of exp
00222     const G4double expxl = -expxu;        // lower bound for arg. of exp
00223     G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
00224     G4ParticleDefinition *aProton = G4Proton::Proton();
00225     G4ParticleDefinition *aLambda = G4Lambda::Lambda();
00226     G4ParticleDefinition *aSigmaZero = G4SigmaZero::SigmaZero();
00227     
00228     // energetically possible to produce pion(s)  -->  inelastic scattering
00229     
00230     G4double n, anpn;
00231     GetNormalizationConstant( availableEnergy, n, anpn );
00232     G4double ran = G4UniformRand();
00233     G4double dum, excs = 0.0;
00234     if( targetParticle.GetDefinition() == aProton )
00235     {
00236       counter = -1;
00237       for( np=0; np<numSec/3 && ran>=excs; ++np )
00238       {
00239         for( nneg=std::max(0,np-1); nneg<=(np+1) && ran>=excs; ++nneg )
00240         {
00241           for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
00242           {
00243             if( ++counter < numMul )
00244             {
00245               nt = np+nneg+nz;
00246               if( nt>0 && nt<=numSec )
00247               {
00248                 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
00249                 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
00250                 if( std::fabs(dum) < 1.0 )
00251                 {
00252                   if( test >= 1.0e-10 )excs += dum*test;
00253                 }
00254                 else
00255                   excs += dum*test;
00256               }
00257             }
00258           }
00259         }
00260       }
00261       if( ran >= excs )  // 3 previous loops continued to the end
00262       {
00263         quasiElastic = true;
00264         return;
00265       }
00266       np--; nneg--; nz--;
00267       G4int ncht = std::max( 1, np-nneg+2 );
00268       switch( ncht )
00269       {
00270        case 1:
00271          if( G4UniformRand() < 0.5 )
00272            currentParticle.SetDefinitionAndUpdateE( aLambda );
00273          else
00274            currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
00275          incidentHasChanged = true;
00276          break;
00277        case 2:
00278          if( G4UniformRand() >= 0.5 )
00279          {
00280            if( G4UniformRand() < 0.5 )
00281              currentParticle.SetDefinitionAndUpdateE( aLambda );
00282            else
00283              currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
00284            incidentHasChanged = true;
00285            targetParticle.SetDefinitionAndUpdateE( aNeutron );
00286            targetHasChanged = true;
00287          }             
00288          break;
00289        default:
00290          targetParticle.SetDefinitionAndUpdateE( aNeutron );
00291          targetHasChanged = true;
00292          break;
00293       }
00294     }
00295     else  // target must be a neutron
00296     {
00297       counter = -1;
00298       for( np=0; np<numSec/3 && ran>=excs; ++np )
00299       {
00300         for( nneg=np; nneg<=(np+2) && ran>=excs; ++nneg )
00301         {
00302           for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
00303           {
00304             if( ++counter < numMul )
00305             {
00306               nt = np+nneg+nz;
00307               if( nt>0 && nt<=numSec )
00308               {
00309                 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
00310                 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
00311                 if( std::fabs(dum) < 1.0 )
00312                 {
00313                   if( test >= 1.0e-10 )excs += dum*test;
00314                 }
00315                 else
00316                   excs += dum*test;
00317               }
00318             }
00319           }
00320         }
00321       }
00322       if( ran >= excs )  // 3 previous loops continued to the end
00323       {
00324         quasiElastic = true;
00325         return;
00326       }
00327       np--; nneg--; nz--;
00328       G4int ncht = std::max( 1, np-nneg+3 );
00329       switch( ncht )
00330       {
00331        case 1:
00332          if( G4UniformRand() < 0.5 )
00333            currentParticle.SetDefinitionAndUpdateE( aLambda );
00334          else
00335            currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
00336          incidentHasChanged = true;
00337          targetParticle.SetDefinitionAndUpdateE( aProton );
00338          targetHasChanged = true;
00339          break;
00340        case 2:
00341          if( G4UniformRand() < 0.5 )
00342          {
00343            if( G4UniformRand() < 0.5 )
00344              currentParticle.SetDefinitionAndUpdateE( aLambda );
00345            else
00346              currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
00347            incidentHasChanged = true;
00348          }
00349          else
00350          {
00351            targetParticle.SetDefinitionAndUpdateE( aProton );
00352            targetHasChanged = true;
00353          }
00354          break;
00355        default:
00356          break;
00357       }
00358   }
00359 
00360   SetUpPions(np, nneg, nz, vec, vecLen);
00361   return;
00362 }
00363 
00364  /* end of file */
00365  

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