source: trunk/source/processes/hadronic/cross_sections/include/G4UPiNuclearCrossSection.hh@ 1344

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

update ti head

File size: 5.4 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// Calculation of the total, elastic and inelastic cross-sections
27// based on Barashenkov parametrisations of pion data
28//
29// 16.08.06 V.Ivanchenko - first implementation
30// 22.01.07 V.Ivanchenko - add cross section interfaces with Z and A
31// 05.03.07 V.Ivanchenko - add IfZAApplicable
32//
33
34
35#ifndef G4UPiNuclearCrossSection_h
36#define G4UPiNuclearCrossSection_h
37
38#include "G4VCrossSectionDataSet.hh"
39#include "G4DynamicParticle.hh"
40#include "G4DataVector.hh"
41#include "G4HadTmpUtil.hh"
42#include "globals.hh"
43#include <vector>
44
45class G4PhysicsTable;
46
47class G4UPiNuclearCrossSection : public G4VCrossSectionDataSet
48{
49public:
50
51 G4UPiNuclearCrossSection();
52
53 virtual ~G4UPiNuclearCrossSection();
54
55 virtual
56 G4bool IsApplicable(const G4DynamicParticle* aParticle,
57 const G4Element* anElement);
58
59 virtual
60 G4bool IsIsoApplicable(const G4DynamicParticle* aParticle,
61 G4int Z, G4int A);
62
63 virtual
64 G4double GetCrossSection(const G4DynamicParticle* aParticle,
65 const G4Element* anElement, G4double T=0.);
66
67 virtual
68 G4double GetZandACrossSection(const G4DynamicParticle* aParticle,
69 G4int Z, G4int A, G4double T=0.);
70
71 G4double GetElasticCrossSection(const G4DynamicParticle* aParticle,
72 const G4Element* anElement);
73
74 G4double GetElasticCrossSection(const G4DynamicParticle* aParticle,
75 G4int Z, G4int A);
76
77 G4double GetInelasticCrossSection(const G4DynamicParticle* aParticle,
78 const G4Element* anElement);
79
80 G4double GetInelasticCrossSection(const G4DynamicParticle* aParticle,
81 G4int Z, G4int A);
82
83 void BuildPhysicsTable(const G4ParticleDefinition&);
84
85 void DumpPhysicsTable(const G4ParticleDefinition&);
86
87private:
88
89 void Initialise();
90
91 void AddDataSet(const G4String& p, const G4double* tot,
92 const G4double* in, const G4double* e, G4int n);
93
94 G4double Interpolate(G4int Z, G4int A, G4double ekin,
95 G4PhysicsTable*);
96
97 G4int NZ;
98 std::vector<G4int> theZ;
99 G4DataVector theA;
100 G4PhysicsTable* piPlusElastic;
101 G4PhysicsTable* piPlusInelastic;
102 G4PhysicsTable* piMinusElastic;
103 G4PhysicsTable* piMinusInelastic;
104
105 G4double aPower;
106 G4double elow;
107 G4double elowest;
108 G4double APower[93];
109
110 const G4ParticleDefinition* piPlus;
111 const G4ParticleDefinition* piMinus;
112};
113
114inline G4bool G4UPiNuclearCrossSection::IsApplicable(
115 const G4DynamicParticle* part,
116 const G4Element* elm)
117{
118 G4int Z = G4lrint(elm->GetZ());
119 G4int A = G4lrint(elm->GetN());
120 return IsIsoApplicable(part, Z, A);
121}
122
123inline G4bool
124G4UPiNuclearCrossSection::IsIsoApplicable(const G4DynamicParticle* part,
125 G4int Z, G4int)
126{
127 return ((part->GetDefinition() == piMinus ||
128 part->GetDefinition() == piPlus) &&
129 Z > 1);
130}
131
132inline G4double G4UPiNuclearCrossSection::GetCrossSection(
133 const G4DynamicParticle* dp,
134 const G4Element* elm, G4double)
135{
136 G4int Z = G4lrint(elm->GetZ());
137 G4int A = G4lrint(elm->GetN());
138 return GetInelasticCrossSection(dp, Z, A);
139}
140
141inline G4double
142G4UPiNuclearCrossSection::GetZandACrossSection(const G4DynamicParticle* dp,
143 G4int Z, G4int A, G4double)
144{
145 return GetInelasticCrossSection(dp, Z, A);
146}
147
148inline G4double G4UPiNuclearCrossSection::GetInelasticCrossSection(
149 const G4DynamicParticle* dp,
150 const G4Element* elm)
151{
152 G4int Z = G4lrint(elm->GetZ());
153 G4int A = G4lrint(elm->GetN());
154 return GetInelasticCrossSection(dp, Z, A);
155}
156
157inline G4double G4UPiNuclearCrossSection::GetElasticCrossSection(
158 const G4DynamicParticle* dp,
159 const G4Element* elm)
160{
161 G4int Z = G4lrint(elm->GetZ());
162 G4int A = G4lrint(elm->GetN());
163 return GetElasticCrossSection(dp, Z, A);
164}
165
166#endif
Note: See TracBrowser for help on using the repository browser.