Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4CHIPSElasticXS.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 file
31//
32//
33// File name: G4CHIPSElasticXS
34//
35// Author Ivantchenko, Geant4, 3-Aug-09
36//
37// Modifications:
38// 31-05-2011 V.Uzhinsky added anti-baryons, Pi+, Pi-, K+, K- cross sections
39// 23-08-2011 V.Ivanchenko migration to new design and cleanup
40//
41
42#include "G4CHIPSElasticXS.hh"
44#include "G4DynamicParticle.hh"
46#include "G4Element.hh"
47#include "G4Proton.hh"
48#include "G4Neutron.hh"
49#include "G4VQCrossSection.hh"
52
55#include "G4QPionPlusElasticCrossSection.hh" // Uzhi
57#include "G4QKaonPlusElasticCrossSection.hh" // Uzhi
58
60 : G4VCrossSectionDataSet("CHIPSElasticXS"),
61 theProton(G4Proton::Proton()),
62 theNeutron(G4Neutron::Neutron()),
63 thEnergy(19*CLHEP::MeV),
64 isInitialized(false)
65{
66 // verboseLevel = 0;
69
70 PBARxsManager = G4QAntiBaryonElasticCrossSection::GetPointer(); // Uzhi
71 PIPxsManager = G4QPionPlusElasticCrossSection::GetPointer(); // Uzhi
72 PIMxsManager = G4QPionMinusElasticCrossSection::GetPointer(); // Uzhi
73 KPxsManager = G4QKaonPlusElasticCrossSection::GetPointer(); // Uzhi
74 KMxsManager = G4QKaonMinusElasticCrossSection::GetPointer(); // Uzhi
75 //Description();
76 theParticle = 0;
77}
78
80{}
81
82
84{
85 char* dirName = getenv("G4PhysListDocDir");
86 if (dirName) {
87 std::ofstream outFile;
88 G4String outFileName = GetName() + ".html";
89 G4String pathName = G4String(dirName) + "/" + outFileName;
90
91 outFile.open(pathName);
92 outFile << "<html>\n";
93 outFile << "<head>\n";
94
95 outFile << "<title>Description of CHIPS Elastic Cross Section</title>\n";
96 outFile << "</head>\n";
97 outFile << "<body>\n";
98
99 outFile << "G4CHIPSElasticXS provides hadron-nuclear elastic scattering\n"
100 << "cross sections for protons and neutrons with incident energies\n"
101 << "between 19 MeV and X GeV. These cross sections represent\n"
102 << "parameterizations developed by M. Kossov. (more detail)\n";
103
104 outFile << "</body>\n";
105 outFile << "</html>\n";
106 outFile.close();
107 }
108}
109
110G4bool
112 G4int Z, G4int /*A*/,
113 const G4Element*, const G4Material*)
114{
115 return (Z <= 2 && dyn->GetKineticEnergy() > thEnergy);
116}
117
120 G4int Z, G4int A,
121 const G4Isotope*, const G4Element*,
122 const G4Material*)
123{
124 G4int N = A - Z;
125 if(Z == 1) {
126 if(N > 1) { N = 1; }
127 } else if(Z == 2) { N = 2; }
128
129 G4double momentum = dyn->GetTotalMomentum();
130 G4int uPDGcode = dyn->GetPDGcode();
131 G4VQCrossSection* CHIPSmanager = 0;
132 G4double cross = 0.0;
133
134 switch(uPDGcode) {
135 case 2212:
136 CHIPSmanager=pCManager;
137 break;
138 case 2112:
139 CHIPSmanager=nCManager;
140 break;
141 case -2212:
142 CHIPSmanager=PBARxsManager;
143 break;
144 case -2112:
145 CHIPSmanager=PBARxsManager;
146 break;
147 case 211:
148 CHIPSmanager=PIPxsManager;
149 break;
150 case -211:
151 CHIPSmanager=PIMxsManager;
152 break;
153 case 321:
154 CHIPSmanager=KPxsManager;
155 break;
156 case -321:
157 CHIPSmanager=KMxsManager;
158 break;
159 case 130:
160 break;
161 case 310:
162 break;
163 case 311:
164 break;
165 case -311:
166 break;
167 default:
168 throw G4HadronicException(__FILE__, __LINE__,
169 "G4CHIPSElasticXS: not applicable for a particle");
170 return cross;
171 }
172 if(CHIPSmanager) {
173 cross = CHIPSmanager->GetCrossSection(false,momentum,Z,N,uPDGcode);
174 } else {
175 cross = 0.5*(KPxsManager->GetCrossSection(false,momentum,Z,N,uPDGcode) +
176 KMxsManager->GetCrossSection(false,momentum,Z,N,uPDGcode));
177 }
178 return cross;
179}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
virtual ~G4CHIPSElasticXS()
virtual void Description() const
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
G4double GetTotalMomentum() const
G4int GetPDGcode() const
static G4VQCrossSection * GetPointer()
const G4String & GetName() const
virtual G4double GetCrossSection(G4bool, G4double, G4int, G4int, G4int pPDG=0)
Definition: DoubConv.h:17