Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4CHIPSElastic.cc
Go to the documentation of this file.
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 // $Id$
27 //
28 //---------------------------------------------------------------------
29 //
30 // Geant4 Class : G4CHIPSElastic
31 //
32 // Author : V.Ivanchenko 29 June 2009
33 //
34 // Modified:
35 // 13.01.10: M.Kosov: Use G4Q(Pr/Neut)ElasticCS instead of G4QElasticCS
36 //
37 //---------------------------------------------------------------------
38 // CHIPS model of hadron elastic scattering
39 //
40 
41 #include "G4CHIPSElastic.hh"
42 #include "G4VQCrossSection.hh"
43 #include "G4ParticleDefinition.hh"
46 
47 #include "G4QAntiBaryonElasticCrossSection.hh" // Uzhi
48 #include "G4QPionPlusElasticCrossSection.hh" // Uzhi
49 #include "G4QPionMinusElasticCrossSection.hh" // Uzhi
50 #include "G4QKaonPlusElasticCrossSection.hh" // Uzhi
51 #include "G4QKaonMinusElasticCrossSection.hh" // Uzhi
52 #include <iostream>
53 
54 
55 G4VQCrossSection* G4CHIPSElastic::pxsManager = 0;
56 G4VQCrossSection* G4CHIPSElastic::nxsManager = 0;
57 
58 G4VQCrossSection* G4CHIPSElastic::PBARxsManager = 0; // Uzhi
59 G4VQCrossSection* G4CHIPSElastic::PIPxsManager = 0;
60 G4VQCrossSection* G4CHIPSElastic::PIMxsManager = 0;
61 G4VQCrossSection* G4CHIPSElastic::KPxsManager = 0;
62 G4VQCrossSection* G4CHIPSElastic::KMxsManager = 0;
63 
65 {
66  if(!pxsManager)
67  {
70 
71  PBARxsManager = G4QAntiBaryonElasticCrossSection::GetPointer(); // Uzhi
72  PIPxsManager = G4QPionPlusElasticCrossSection::GetPointer(); // Uzhi
73  PIMxsManager = G4QPionMinusElasticCrossSection::GetPointer(); // Uzhi
74  KPxsManager = G4QKaonPlusElasticCrossSection::GetPointer(); // Uzhi
75  KMxsManager = G4QKaonMinusElasticCrossSection::GetPointer(); // Uzhi
76  }
77  //Description();
78 }
79 
81 {}
82 
84 {
85  char* dirName = getenv("G4PhysListDocDir");
86  if (dirName) {
87  std::ofstream outFile;
88  G4String outFileName = GetModelName() + ".html";
89  G4String pathName = G4String(dirName) + "/" + outFileName;
90  outFile.open(pathName);
91  outFile << "<html>\n";
92  outFile << "<head>\n";
93 
94  outFile << "<title>Description of G4CHIPSElastic</title>\n";
95  outFile << "</head>\n";
96  outFile << "<body>\n";
97 
98  outFile << "The G4CHIPSElastic model performs hadron-nucleus elastic\n"
99  << "scattering using the parameterized elastic cross sections\n"
100  << "of M. Kossov\n";
101 
102  outFile << "</body>\n";
103  outFile << "</html>\n";
104  outFile.close();
105  }
106 }
107 
108 
109 G4double
111  G4double plab, G4int Z, G4int A)
112 {
113  G4int N = A - Z;
114  if(Z == 1 && N == 2) { N = 1; }
115  else if(Z == 2 && N == 1) { N = 2; }
116  G4int projPDG = p->GetPDGEncoding();
117  G4double cs = 0.;
118  if (projPDG==2212) { cs = pxsManager->GetCrossSection(false,plab,Z,N,projPDG); }
119  else if(projPDG==2112) { cs = nxsManager->GetCrossSection(false,plab,Z,N,projPDG); }
120  else if(projPDG==-2212){ cs = PBARxsManager->GetCrossSection(false,plab,Z,N,projPDG); } //Pbar
121  else if(projPDG== 211) { cs = PIPxsManager->GetCrossSection(false,plab,Z,N,projPDG); } // Pi+
122  else if(projPDG==-211) { cs = PIMxsManager->GetCrossSection(false,plab,Z,N,projPDG); } // Pi-
123  else if(projPDG== 321) { cs = KPxsManager->GetCrossSection(false,plab,Z,N,projPDG); } // K+
124  else if(projPDG==-321) { cs = KMxsManager->GetCrossSection(false,plab,Z,N,projPDG); } // K-
125 
126  G4double t = 0.0;
127  if(cs > 0.0)
128  {
129  if (projPDG== 2212) { t = pxsManager->GetExchangeT(Z,N,projPDG); }
130  else if(projPDG== 2112) { t = nxsManager->GetExchangeT(Z,N,projPDG); }
131  else if(projPDG==-2212) { t = PBARxsManager->GetExchangeT(Z,N,projPDG); } // Pbar
132  else if(projPDG== 211) { t = PIPxsManager->GetExchangeT(Z,N,projPDG); } // Pi+
133  else if(projPDG== -211) { t = PIMxsManager->GetExchangeT(Z,N,projPDG); } // Pi-
134  else if(projPDG== 321) { t = KPxsManager->GetExchangeT(Z,N,projPDG); } // K+
135  else if(projPDG== -321) { t = KMxsManager->GetExchangeT(Z,N,projPDG); } // K-
136  }
137  else { t = G4HadronElastic::SampleInvariantT(p, plab, Z, A); }
138  return t;
139 }
140