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

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