source: trunk/source/processes/hadronic/cross_sections/include/G4PiNuclearCrossSection.hh@ 1347

Last change on this file since 1347 was 1340, checked in by garnier, 15 years ago

update ti head

File size: 6.1 KB
RevLine 
[819]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#ifndef G4PiNuclearCrossSection_h
27#define G4PiNuclearCrossSection_h
28
29#include "G4VCrossSectionDataSet.hh"
30
31#include "globals.hh"
32#include "G4PionMinus.hh"
33#include "G4PionPlus.hh"
34#include "G4PiData.hh"
35#include "G4HadTmpUtil.hh"
36
37class G4PiNuclearCrossSection : public G4VCrossSectionDataSet
38{
39 public:
40
41 G4PiNuclearCrossSection();
42 virtual ~G4PiNuclearCrossSection();
43
[1340]44 G4bool IsApplicable(const G4DynamicParticle* aParticle,
45 const G4Element* anElement)
[819]46 {
47 G4bool result = false;
48 if(aParticle->GetDefinition() == G4PionMinus::PionMinus()) result=true;
49 if(aParticle->GetDefinition() == G4PionPlus::PionPlus()) result=true;
50 if(G4lrint(anElement->GetZ()) == 1) result = false;
[1228]51 if(aParticle->GetKineticEnergy() > 99.9*TeV) result=false;
[819]52 return result;
53 }
54
[1340]55 G4bool
56 IsIsoApplicable(const G4DynamicParticle* particle, G4int ZZ, G4int /*AA*/)
[819]57 {
58 G4bool result = false;
59 if(particle->GetDefinition() == G4PionMinus::PionMinus()) result=true;
60 if(particle->GetDefinition() == G4PionPlus::PionPlus()) result=true;
[1340]61 if(ZZ == 1) result = false;
[1228]62 if(particle->GetKineticEnergy() > 99.9*TeV) result=false;
[819]63 return result;
64 }
65
[1340]66
[819]67 G4double GetCrossSection(const G4DynamicParticle* particle,
68 const G4Element* element,
69 G4double temperature)
70 {
[1340]71 G4int Z = G4lrint(element->GetZ());
72 G4int A = G4lrint(element->GetN());
73 return GetZandACrossSection(particle, Z, A, temperature);
[819]74 }
75
[1340]76 G4double GetZandACrossSection(const G4DynamicParticle* aParticle,
77 G4int ZZ, G4int AA,
[819]78 G4double /*aTemperature*/);
79
[1340]80 G4double GetTotalXsc() {return fTotalXsc;};
81 G4double GetElasticXsc() {return fElasticXsc;};
[819]82
83
84 void BuildPhysicsTable(const G4ParticleDefinition&) {}
85 void DumpPhysicsTable(const G4ParticleDefinition&) {}
86
87 private:
88 G4double Interpolate(G4int Z1, G4int Z2, G4int Z, G4double x1, G4double x2);
89
90// add Hydrogen from PDG group.
91
92static const G4double e1[38];
93static const G4double he_t[38];
94static const G4double he_in[38];
95static const G4double be_m_t[38];
96static const G4double be_m_in[38];
97static const G4double be_p_t[24];
98static const G4double be_p_in[24];
99static const G4double e2[39];
100static const G4double c_m_t[39];
101static const G4double c_m_in[39];
102static const G4double c_p_t[24];
103static const G4double c_p_in[24];
104static const G4double n_m_t[39];
105static const G4double n_m_in[39];
106static const G4double n_p_t[27];
107static const G4double n_p_in[27];
108static const G4double e3[31];
109static const G4double o_m_t[31];
110static const G4double o_m_in[31];
111static const G4double o_p_t[20];
112static const G4double o_p_in[20];
113static const G4double na_m_t[31];
114static const G4double na_m_in[31];
115static const G4double na_p_t[22];
116static const G4double na_p_in[22];
117static const G4double e3_1[31];
118static const G4double al_m_t[31];
119static const G4double al_m_in[31];
120static const G4double al_p_t[21];
121static const G4double al_p_in[21];
122static const G4double ca_m_t[31];
123static const G4double ca_m_in[31];
124static const G4double ca_p_t[23];
125static const G4double ca_p_in[23];
126
127static const G4double e4[32];
128static const G4double fe_m_t[32];
129static const G4double fe_m_in[32];
130static const G4double fe_p_t[25];
131static const G4double fe_p_in[25];
132static const G4double cu_m_t[32];
133static const G4double cu_m_in[32];
134static const G4double cu_p_t[25];
135static const G4double cu_p_in[25];
136static const G4double e5[34];
137static const G4double mo_m_t[34];
138static const G4double mo_m_in[34];
139static const G4double mo_p_t[27];
140static const G4double mo_p_in[27];
141static const G4double cd_m_t[34];
142static const G4double cd_m_in[34];
143static const G4double cd_p_t[28];
144static const G4double cd_p_in[28];
145static const G4double e6[35];
146static const G4double sn_m_t[35];
147static const G4double sn_m_in[35];
148static const G4double sn_p_t[29];
149static const G4double sn_p_in[29];
150static const G4double w_m_t[35];
151static const G4double w_m_in[35];
152static const G4double w_p_t[30];
153static const G4double w_p_in[30];
154static const G4double e7[35];
155static const G4double pb_m_t[35];
156static const G4double pb_m_in[35];
157static const G4double pb_p_t[30];
158static const G4double pb_p_in[30];
159static const G4double u_m_t[35];
160static const G4double u_m_in[35];
161static const G4double u_p_t[30];
162static const G4double u_p_in[30];
163
164std::vector<G4int> theZ;
165std::vector<G4PiData *> thePipData;
166std::vector<G4PiData *> thePimData;
167
168 // cross sections
169
170 G4double fTotalXsc;
171 G4double fElasticXsc;
172
173
174};
175
176#endif
Note: See TracBrowser for help on using the repository browser.