Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 // 27 //--------------------------------------------------------------------- 28 // 29 // Geant4 Class : G4ChipsElasticModel 30 // 31 // Author : V.Ivanchenko 29 June 2009 32 // 33 // Modified: 34 // 13.01.10: M.Kosov: Use G4Q(Pr/Neut)ElasticCS instead of G4QElasticCS 35 // 36 //--------------------------------------------------------------------- 37 // CHIPS model of hadron elastic scattering 38 // 39 40 #include "G4ChipsElasticModel.hh" 41 #include "G4ParticleDefinition.hh" 42 #include <iostream> 43 44 #include "G4CrossSectionDataSetRegistry.hh" 45 46 G4ChipsElasticModel::G4ChipsElasticModel() : G4HadronElastic("hElasticCHIPS") 47 { 48 pxsManager = (G4ChipsProtonElasticXS*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsProtonElasticXS::Default_Name()); 49 nxsManager = (G4ChipsNeutronElasticXS*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsNeutronElasticXS::Default_Name()); 50 PBARxsManager = (G4ChipsAntiBaryonElasticXS*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsAntiBaryonElasticXS::Default_Name()); 51 PIPxsManager = (G4ChipsPionPlusElasticXS*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsPionPlusElasticXS::Default_Name()); 52 PIMxsManager = (G4ChipsPionMinusElasticXS*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsPionMinusElasticXS::Default_Name()); 53 KPxsManager = (G4ChipsKaonPlusElasticXS*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsKaonPlusElasticXS::Default_Name()); 54 KMxsManager = (G4ChipsKaonMinusElasticXS*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsKaonMinusElasticXS::Default_Name()); 55 //Description(); 56 } 57 58 G4ChipsElasticModel::~G4ChipsElasticModel() 59 {} 60 61 void G4ChipsElasticModel::ModelDescription(std::ostream& outFile) const 62 { 63 64 outFile << "The G4ChipsElasticModel model performs hadron-nucleus elastic\n" 65 << "scattering using the parameterized elastic cross sections\n" 66 << "of M. Kossov\n"; 67 68 } 69 70 71 G4double 72 G4ChipsElasticModel::SampleInvariantT(const G4ParticleDefinition* p, 73 G4double plab, G4int Z, G4int A) 74 { 75 G4int N = A - Z; 76 if(Z == 1 && N == 2) { N = 1; } 77 else if(Z == 2 && N == 1) { N = 2; } 78 G4int projPDG = p->GetPDGEncoding(); 79 G4double cs = 0.; 80 if (projPDG==2212) { cs = pxsManager->GetChipsCrossSection(plab,Z,N,projPDG); } 81 else if(projPDG==2112) { cs = nxsManager->GetChipsCrossSection(plab,Z,N,projPDG); } 82 else if(projPDG==-2212){ cs = PBARxsManager->GetChipsCrossSection(plab,Z,N,projPDG); } //Pbar 83 else if(projPDG== 211) { cs = PIPxsManager->GetChipsCrossSection(plab,Z,N,projPDG); } // Pi+ 84 else if(projPDG==-211) { cs = PIMxsManager->GetChipsCrossSection(plab,Z,N,projPDG); } // Pi- 85 else if(projPDG== 321) { cs = KPxsManager->GetChipsCrossSection(plab,Z,N,projPDG); } // K+ 86 else if(projPDG==-321) { cs = KMxsManager->GetChipsCrossSection(plab,Z,N,projPDG); } // K- 87 88 G4double t = 0.0; 89 if(cs > 0.0) 90 { 91 if (projPDG== 2212) { t = pxsManager->GetExchangeT(Z,N,projPDG); } 92 else if(projPDG== 2112) { t = nxsManager->GetExchangeT(Z,N,projPDG); } 93 else if(projPDG==-2212) { t = PBARxsManager->GetExchangeT(Z,N,projPDG); } // Pbar 94 else if(projPDG== 211) { t = PIPxsManager->GetExchangeT(Z,N,projPDG); } // Pi+ 95 else if(projPDG== -211) { t = PIMxsManager->GetExchangeT(Z,N,projPDG); } // Pi- 96 else if(projPDG== 321) { t = KPxsManager->GetExchangeT(Z,N,projPDG); } // K+ 97 else if(projPDG== -321) { t = KMxsManager->GetExchangeT(Z,N,projPDG); } // K- 98 } 99 else { t = G4HadronElastic::SampleInvariantT(p, plab, Z, A); } 100 return t; 101 } 102 103