G4ChipsElasticModel.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 //
00030 // Geant4 Class : G4ChipsElasticModel
00031 //
00032 // Author : V.Ivanchenko 29 June 2009 
00033 //  
00034 // Modified:
00035 // 13.01.10: M.Kosov: Use G4Q(Pr/Neut)ElasticCS instead of G4QElasticCS
00036 //
00037 //---------------------------------------------------------------------
00038 // CHIPS model of hadron elastic scattering
00039 //
00040 
00041 #include "G4ChipsElasticModel.hh"
00042 #include "G4ParticleDefinition.hh"
00043 #include <iostream>
00044 
00045 #include "G4CrossSectionDataSetRegistry.hh"
00046 
00047 G4ChipsElasticModel::G4ChipsElasticModel() : G4HadronElastic("hElasticCHIPS")
00048 {
00049     pxsManager    = (G4ChipsProtonElasticXS*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsProtonElasticXS::Default_Name());
00050     nxsManager    = (G4ChipsNeutronElasticXS*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsNeutronElasticXS::Default_Name());
00051     PBARxsManager = (G4ChipsAntiBaryonElasticXS*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsAntiBaryonElasticXS::Default_Name());
00052     PIPxsManager  = (G4ChipsPionPlusElasticXS*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsPionPlusElasticXS::Default_Name());
00053     PIMxsManager  = (G4ChipsPionMinusElasticXS*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsPionMinusElasticXS::Default_Name());
00054     KPxsManager   = (G4ChipsKaonPlusElasticXS*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsKaonPlusElasticXS::Default_Name());
00055     KMxsManager   = (G4ChipsKaonMinusElasticXS*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsKaonMinusElasticXS::Default_Name());
00056   //Description();
00057 }
00058 
00059 G4ChipsElasticModel::~G4ChipsElasticModel()
00060 {}
00061 
00062 void G4ChipsElasticModel::Description() const
00063 {
00064   char* dirName = getenv("G4PhysListDocDir");
00065   if (dirName) {
00066     std::ofstream outFile;
00067     G4String outFileName = GetModelName() + ".html";
00068     G4String pathName = G4String(dirName) + "/" + outFileName;
00069     outFile.open(pathName);
00070     outFile << "<html>\n";
00071     outFile << "<head>\n";
00072 
00073     outFile << "<title>Description of G4ChipsElasticModel</title>\n";
00074     outFile << "</head>\n";
00075     outFile << "<body>\n";
00076 
00077     outFile << "The G4ChipsElasticModel model performs hadron-nucleus elastic\n"
00078             << "scattering using the parameterized elastic cross sections\n"
00079             << "of M. Kossov\n";
00080 
00081     outFile << "</body>\n";
00082     outFile << "</html>\n";
00083     outFile.close();
00084   }
00085 }
00086 
00087 
00088 G4double 
00089 G4ChipsElasticModel::SampleInvariantT(const G4ParticleDefinition* p, 
00090                                  G4double plab, G4int Z, G4int A)
00091 {
00092   G4int N = A - Z;
00093   if(Z == 1 && N == 2)      { N = 1; }
00094   else if(Z == 2 && N == 1) { N = 2; }
00095   G4int projPDG = p->GetPDGEncoding();
00096   G4double cs = 0.;
00097   if     (projPDG==2212) { cs = pxsManager->GetChipsCrossSection(plab,Z,N,projPDG); }
00098   else if(projPDG==2112) { cs = nxsManager->GetChipsCrossSection(plab,Z,N,projPDG); }
00099   else if(projPDG==-2212){ cs = PBARxsManager->GetChipsCrossSection(plab,Z,N,projPDG); } //Pbar
00100   else if(projPDG== 211) { cs = PIPxsManager->GetChipsCrossSection(plab,Z,N,projPDG); } // Pi+
00101   else if(projPDG==-211) { cs = PIMxsManager->GetChipsCrossSection(plab,Z,N,projPDG); } // Pi-
00102   else if(projPDG== 321) { cs = KPxsManager->GetChipsCrossSection(plab,Z,N,projPDG); } // K+
00103   else if(projPDG==-321) { cs = KMxsManager->GetChipsCrossSection(plab,Z,N,projPDG); } // K-
00104 
00105   G4double t = 0.0;
00106   if(cs > 0.0)
00107   {
00108     if     (projPDG== 2212) { t = pxsManager->GetExchangeT(Z,N,projPDG); }
00109     else if(projPDG== 2112) { t = nxsManager->GetExchangeT(Z,N,projPDG); }
00110     else if(projPDG==-2212) { t = PBARxsManager->GetExchangeT(Z,N,projPDG); } // Pbar
00111     else if(projPDG==  211) { t = PIPxsManager->GetExchangeT(Z,N,projPDG); }  // Pi+
00112     else if(projPDG== -211) { t = PIMxsManager->GetExchangeT(Z,N,projPDG); }  // Pi-
00113     else if(projPDG==  321) { t = KPxsManager->GetExchangeT(Z,N,projPDG); }   // K+
00114     else if(projPDG== -321) { t = KMxsManager->GetExchangeT(Z,N,projPDG); }   // K-
00115   }
00116   else  { t = G4HadronElastic::SampleInvariantT(p, plab, Z, A); }
00117   return t;
00118 }
00119 

Generated on Mon May 27 17:47:52 2013 for Geant4 by  doxygen 1.4.7