source: trunk/source/processes/hadronic/cross_sections/include/G4NucleonNuclearCrossSection.hh @ 1315

Last change on this file since 1315 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 5.8 KB
Line 
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// author: Vladimir.Grichine@cern.ch
27//
28// Implements data from: Barashenkov V.S., Nucleon-Nucleus Cross Section,
29// Preprint JINR P2-89-770, p. 12, Dubna 1989 (scanned version from KEK)
30// Based on G. Folger version of G4PiNuclearCrossSection class
31//
32// Modified:
33// 05.03.07 V.Ivanchenko - add IfZAApplicable, remove "debug"
34// 06.03.07 V.Ivanchenko - add GetElasticCrossSection for combined dataset
35//
36
37#ifndef G4NucleonNuclearCrossSection_h
38#define G4NucleonNuclearCrossSection_h
39
40#include "G4VCrossSectionDataSet.hh"
41#include "G4ParticleDefinition.hh"
42
43#include "globals.hh"
44#include "G4PiData.hh"
45#include "G4HadTmpUtil.hh"
46
47class G4NucleonNuclearCrossSection : public G4VCrossSectionDataSet
48{
49public:
50 
51  G4NucleonNuclearCrossSection();
52  virtual ~G4NucleonNuclearCrossSection();
53
54  virtual G4bool IsApplicable(const G4DynamicParticle* aParticle, const G4Element* anElement);
55
56  virtual G4bool IsZAApplicable(const G4DynamicParticle* aParticle, G4double Z, G4double A);
57
58  G4double GetCrossSection(const G4DynamicParticle* aParticle, 
59                           const G4Element* anElement,
60                           G4double T=0.);
61
62  G4double GetIsoZACrossSection(const G4DynamicParticle* aParticle, 
63                                G4double ZZ, G4double AA, G4double T=0. );
64
65  G4double GetElasticCrossSection(const G4DynamicParticle* aParticle, 
66                                  G4double ZZ, G4double AA);
67
68  G4double GetTotalXsc()  { return fTotalXsc;   };
69  G4double GetElasticXsc(){ return fElasticXsc; };
70
71
72  void BuildPhysicsTable(const G4ParticleDefinition&) {}
73  void DumpPhysicsTable(const G4ParticleDefinition&) {}
74 
75private:
76
77  G4double Interpolate(G4int Z1, G4int Z2, G4int Z, G4double x1, G4double x2);
78
79// add Hydrogen from PDG group.
80
81static const G4double e1[44];
82
83static const G4double he_m_t[44];
84static const G4double he_m_in[44];
85static const G4double he_p_in[44];
86
87static const G4double be_m_t[44];
88static const G4double be_m_in[44];
89static const G4double be_p_in[44];
90
91static const G4double c_m_t[44];
92static const G4double c_m_in[44];
93static const G4double c_p_in[44];
94
95
96static const G4double e2[44];
97
98static const G4double n_m_t[44];
99static const G4double n_m_in[44];
100static const G4double n_p_in[44];
101
102static const G4double o_m_t[44];
103static const G4double o_m_in[44];
104static const G4double o_p_in[44];
105
106static const G4double na_m_t[44];
107static const G4double na_m_in[44];
108static const G4double na_p_in[44];
109
110
111static const G4double e3[45];
112
113// static const G4double e3_1[31];
114
115static const G4double al_m_t[45];
116static const G4double al_m_in[45];
117static const G4double al_p_in[45];
118
119static const G4double si_m_t[45];
120static const G4double si_m_in[45];
121static const G4double si_p_in[45];
122
123static const G4double ca_m_t[45];
124static const G4double ca_m_in[45];
125static const G4double ca_p_in[45];
126
127
128static const G4double e4[47];
129
130static const G4double fe_m_t[47];
131static const G4double fe_m_in[47];
132static const G4double fe_p_in[47];
133
134static const G4double cu_m_t[47];
135static const G4double cu_m_in[47];
136static const G4double cu_p_in[47];
137
138static const G4double mo_m_t[47];
139static const G4double mo_m_in[47];
140static const G4double mo_p_in[47];
141
142
143static const G4double e5[48];
144
145static const G4double cd_m_t[48];
146static const G4double cd_m_in[48];
147static const G4double cd_p_in[48];
148
149static const G4double sn_m_t[48];
150static const G4double sn_m_in[48];
151static const G4double sn_p_in[48];
152
153static const G4double w_m_t[48];
154static const G4double w_m_in[48];
155static const G4double w_p_in[48];
156
157
158
159
160static const G4double e6[46];
161
162
163// static const G4double e7[46];
164
165static const G4double pb_m_t[46];
166static const G4double pb_m_in[46];
167static const G4double pb_p_in[46];
168
169static const G4double u_m_t[46];
170static const G4double u_m_in[46];
171static const G4double u_p_in[46];
172
173// vectors for treatment
174
175std::vector< G4int >     theZ;
176std::vector< G4PiData* > thePipData;
177std::vector< G4PiData* > thePimData;
178
179  // cross sections
180
181  G4double fTotalXsc;
182  G4double fElasticXsc;
183
184  // particles
185  const G4ParticleDefinition* theProton;
186  const G4ParticleDefinition* theNeutron;
187
188};
189
190inline
191G4double G4NucleonNuclearCrossSection::GetElasticCrossSection(
192         const G4DynamicParticle* dp, G4double ZZ, G4double AA)
193{
194  GetIsoZACrossSection(dp, ZZ, AA);
195  return fElasticXsc;
196}
197
198#endif
Note: See TracBrowser for help on using the repository browser.