G4TripathiCrossSection.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 // Implementation of formulas in analogy to NASA technical paper 3621 by 
00027 // Tripathi, et al.
00028 // 
00029 // 26-Dec-2006 Isotope dependence added by D. Wright
00030 // 19-Aug-2011 V.Ivanchenko move to new design and make x-section per element
00031 //
00032 
00033 #include "G4TripathiCrossSection.hh"
00034 #include "G4PhysicalConstants.hh"
00035 #include "G4SystemOfUnits.hh"
00036 #include "G4DynamicParticle.hh"
00037 #include "G4ParticleTable.hh"
00038 #include "G4IonTable.hh"
00039 #include "G4HadTmpUtil.hh"
00040 #include "G4Proton.hh"
00041 #include "G4NistManager.hh"
00042 
00043 G4TripathiCrossSection::G4TripathiCrossSection()
00044  : G4VCrossSectionDataSet("Tripathi")
00045 {}
00046 
00047 G4TripathiCrossSection::~G4TripathiCrossSection() 
00048 {}
00049 
00050 G4bool 
00051 G4TripathiCrossSection::IsElementApplicable(const G4DynamicParticle* aPart, 
00052                                             G4int, const G4Material*)
00053 {
00054   G4bool result = false;
00055   if ( (aPart->GetDefinition()->GetBaryonNumber()>2.5) &&
00056        ( aPart->GetKineticEnergy()/aPart->GetDefinition()->GetBaryonNumber()<1*GeV) ) {
00057     result = true;
00058   }
00059   return result;
00060 }
00061 
00062 G4double G4TripathiCrossSection::
00063 GetElementCrossSection(const G4DynamicParticle* aPart, G4int ZZ,  
00064                        const G4Material*) 
00065 {
00066   G4double result = 0.;
00067   G4double targetAtomicNumber = G4NistManager::Instance()->GetAtomicMassAmu(ZZ);
00068   G4double nTargetProtons = ZZ;
00069   
00070   G4double kineticEnergy = aPart->GetKineticEnergy()/MeV;
00071   G4double nProjProtons = aPart->GetDefinition()->GetPDGCharge();
00072   G4double projectileAtomicNumber = 
00073     aPart->GetDefinition()->GetBaryonNumber();
00074 
00075   const G4double nuleonRadius=1.1E-15;
00076   const G4double myNuleonRadius=1.36E-15;
00077   
00078   // needs target mass
00079   G4double targetMass = 
00080      G4ParticleTable::GetParticleTable()->GetIonTable()
00081          ->GetIonMass(G4lrint(nTargetProtons), G4lrint(targetAtomicNumber));
00082   G4LorentzVector pTarget(0,0,0,targetMass); 
00083   G4LorentzVector pProjectile(aPart->Get4Momentum());
00084   pTarget = pTarget+pProjectile;
00085   G4double E_cm = (pTarget.mag()-targetMass-pProjectile.m())/MeV;
00086   if(E_cm <= DBL_MIN) { return result; }
00087   // done
00088   G4double r_rms_p = 0.6 * myNuleonRadius * 
00089                                    std::pow(projectileAtomicNumber, 1./3.);
00090   G4double r_rms_t = 0.6 * myNuleonRadius * 
00091                                    std::pow(targetAtomicNumber, 1./3.);
00092   
00093   // done
00094   G4double r_p = 1.29*r_rms_p/nuleonRadius ;
00095   G4double r_t = 1.29*r_rms_t/nuleonRadius;
00096   
00097   // done
00098   G4double Radius = r_p + r_t + 
00099            1.2*(std::pow(targetAtomicNumber, 1./3.) + 
00100             std::pow(projectileAtomicNumber, 1./3.))/std::pow(E_cm, 1./3.);
00101 
00102   //done
00103   G4double B = 1.44*nProjProtons*nTargetProtons/Radius;
00104   if(E_cm <= B) return result; 
00105   // done
00106   G4double Energy = kineticEnergy/projectileAtomicNumber;
00107 
00108   // done
00109   //
00110   // Note that this correction to G4TripathiCrossSection is just to accurately
00111   // reflect Tripathi's algorithm.  However, if you're using alpha 
00112   // particles/protons consider using the more accurate 
00113   // G4TripathiLightCrossSection, which Tripathi developed specifically for 
00114   // light systems.
00115   //
00116 
00117   G4double D;
00118   if (nProjProtons==1 && projectileAtomicNumber==1)
00119   {
00120     D = 2.05;
00121   }
00122   else if (nProjProtons==2 && projectileAtomicNumber==4)
00123   {
00124     D = 2.77-(8.0E-3*targetAtomicNumber)+
00125           (1.8E-5*targetAtomicNumber*targetAtomicNumber)
00126                    - 0.8/(1+std::exp((250.-Energy)/75.));
00127   }
00128   else
00129   {
00130   //
00131   // This is the original value used in the G4TripathiCrossSection 
00132   // implementation, and was used for all projectile/target conditions.  
00133   // I'm not touching this, although judging from Tripathi's paper, this is 
00134   // valid for cases where the nucleon density changes little with A.
00135   // 
00136     D = 1.75;
00137   }
00138   // done
00139   G4double C_E = D * (1-std::exp(-Energy/40.)) - 
00140        0.292*std::exp(-Energy/792.)*std::cos(0.229*std::pow(Energy, 0.453));
00141   
00142   // done
00143   G4double S = std::pow(projectileAtomicNumber, 1./3.)*
00144                std::pow(targetAtomicNumber, 1./3.)/
00145                (std::pow(projectileAtomicNumber, 1./3.) + 
00146                std::pow(targetAtomicNumber, 1./3.)); 
00147   
00148   // done
00149   G4double deltaE = 1.85*S + 0.16*S/std::pow(E_cm,1./3.) - C_E +
00150                     0.91*(targetAtomicNumber-2.*nTargetProtons)*nProjProtons/
00151                     (targetAtomicNumber*projectileAtomicNumber);
00152   
00153   // done 
00154   result = pi * nuleonRadius*nuleonRadius * 
00155            std::pow(( std::pow(targetAtomicNumber, 1./3.) + 
00156                  std::pow(projectileAtomicNumber, 1./3.) + deltaE),2.) * 
00157                  (1-B/E_cm);
00158   
00159   if(result < 0.) { result = 0.; }
00160   return result*m2;
00161 
00162 }
00163 

Generated on Mon May 27 17:50:03 2013 for Geant4 by  doxygen 1.4.7