G4IonsShenCrossSection.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 // 18-Sep-2003 First version is written by T. Koi
00027 // 12-Nov-2003 Add energy check at lower side T. Koi
00028 // 15-Nov-2006 Above 10GeV/n Cross Section become constant T. Koi (SLAC/SCCS)
00029 // 23-Dec-2006 Isotope dependence adde by D. Wright
00030 // 14-Mar-2011 Moved constructor, destructor and virtual methods to source by V.Ivanchenko
00031 // 19-Aug-2011 V.Ivanchenko move to new design and make x-section per element
00032 //
00033 
00034 #include "G4IonsShenCrossSection.hh"
00035 #include "G4PhysicalConstants.hh"
00036 #include "G4SystemOfUnits.hh"
00037 #include "G4DynamicParticle.hh"
00038 #include "G4NucleiProperties.hh"
00039 #include "G4HadTmpUtil.hh"
00040 #include "G4NistManager.hh"
00041 #include "G4Pow.hh"
00042 
00043 G4IonsShenCrossSection::G4IonsShenCrossSection()
00044   : G4VCrossSectionDataSet("IonsShen"),
00045     upperLimit( 10*GeV ), lowerLimit( 10*MeV ), r0 ( 1.1 )
00046 {}
00047 
00048 G4IonsShenCrossSection::~G4IonsShenCrossSection()
00049 {}
00050 
00051 void
00052 G4IonsShenCrossSection::CrossSectionDescription(std::ostream& outFile) const
00053 {
00054   outFile << "G4IonsShenCrossSection calculates the total reaction cross\n"
00055           << "section for nucleus-nucleus scattering using the Shen\n"
00056           << "parameterization.  It is valid for projectiles and targets of\n"
00057           << "all Z, and projectile energies up to 1 TeV/n.  Above 10 GeV/n"
00058           << "the cross section is constant.  Below 10 MeV/n zero cross\n"
00059           << "is returned.\n";
00060 }
00061    
00062 G4bool G4IonsShenCrossSection::IsElementApplicable(const G4DynamicParticle* aDP, 
00063                                                    G4int, const G4Material*)
00064 {
00065   return (1 <= aDP->GetDefinition()->GetBaryonNumber());
00066 }
00067 
00068 G4double 
00069 G4IonsShenCrossSection::GetElementCrossSection(const G4DynamicParticle* aParticle, 
00070                                                G4int Z, 
00071                                                const G4Material*)
00072 {
00073   G4int A = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z));
00074   return GetIsoCrossSection(aParticle, Z, A);
00075 }
00076 
00077 G4double G4IonsShenCrossSection::GetIsoCrossSection(const G4DynamicParticle* aParticle, 
00078                                                     G4int Zt, G4int At,  
00079                                                     const G4Isotope*,
00080                                                     const G4Element*,
00081                                                     const G4Material*)
00082 
00083 {
00084    G4double xsection = 0.0;
00085 
00086    G4int Ap = aParticle->GetDefinition()->GetBaryonNumber();
00087    G4int Zp = G4lrint(aParticle->GetDefinition()->GetPDGCharge()/eplus); 
00088    G4double ke_per_N = aParticle->GetKineticEnergy() / Ap; 
00089    if ( ke_per_N > upperLimit ) { ke_per_N = upperLimit; }
00090 
00091    // Apply energy check, if less than lower limit then 0 value is returned
00092    //if (  ke_per_N < lowerLimit ) { return xsection; }
00093 
00094    G4Pow* g4pow = G4Pow::GetInstance();
00095   
00096    G4double cubicrAt = g4pow->Z13(At);
00097    G4double cubicrAp = g4pow->Z13(Ap);
00098  
00099    G4double Rt = 1.12 * cubicrAt - 0.94 * ( 1.0 / cubicrAt );
00100    G4double Rp = 1.12 * cubicrAp - 0.94 * ( 1.0 / cubicrAp );
00101 
00102    G4double r = Rt + Rp + 3.2;   // in fm
00103    G4double b = 1.0;   // in MeV/fm
00104    G4double targ_mass = G4NucleiProperties::GetNuclearMass(At, Zt);
00105 
00106    G4double proj_mass = aParticle->GetMass();
00107    G4double proj_momentum = aParticle->GetMomentum().mag();
00108 
00109    G4double Ecm = calEcmValue (proj_mass, targ_mass, proj_momentum); 
00110 
00111    G4double B = 1.44 * Zt * Zp / r - b * Rt * Rp / ( Rt + Rp ); 
00112    if(Ecm <= B) { return xsection; }
00113 
00114    G4double c = calCeValue ( ke_per_N / MeV  );  
00115 
00116    G4double R1 = r0 * (cubicrAt + cubicrAp + 1.85*cubicrAt*cubicrAp/(cubicrAt + cubicrAp) - c); 
00117 
00118    G4double R2 = 1.0 * ( At - 2 * Zt ) * Zp / ( Ap * At );
00119 
00120 
00121    G4double R3 = (0.176 / g4pow->A13(Ecm)) * cubicrAt * cubicrAp /(cubicrAt + cubicrAp);
00122 
00123    G4double R = R1 + R2 + R3;
00124 
00125    xsection = 10 * pi * R * R * ( 1 - B / Ecm );   
00126    xsection = xsection * millibarn;   // mulitply xsection by millibarn
00127 
00128    return xsection; 
00129 }
00130 
00131 G4double
00132 G4IonsShenCrossSection::calEcmValue(const G4double mp, const G4double mt,
00133                                     const G4double Plab)
00134 {
00135    G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
00136    G4double Ecm = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
00137    G4double Pcm = Plab * mt / Ecm;
00138    G4double KEcm = std::sqrt ( Pcm * Pcm + mp * mp ) - mp;
00139    return KEcm;
00140 }
00141 
00142 
00143 G4double G4IonsShenCrossSection::calCeValue(const G4double ke)
00144 {
00145   // Calculate c value 
00146   // This value is indepenent from projectile and target particle 
00147   // ke is projectile kinetic energy per nucleon in the Lab system 
00148   // with MeV unit 
00149   // fitting function is made by T. Koi 
00150   // There are no data below 30 MeV/n in Kox et al., 
00151 
00152    G4double Ce; 
00153    G4double log10_ke = std::log10 ( ke );   
00154    if (log10_ke > 1.5) 
00155    {
00156      Ce = -10.0/std::pow(G4double(log10_ke), G4double(5)) + 2.0;
00157    }
00158    else
00159    {
00160      Ce = (-10.0/std::pow(G4double(1.5), G4double(5) ) + 2.0) /
00161          std::pow(G4double(1.5) , G4double(3)) * std::pow(G4double(log10_ke), G4double(3));
00162    }
00163    return Ce;
00164 }
00165 

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