Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CascadeData.hh
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//
27// 20100507 M. Kelsey -- Use template arguments to dimension const-refs
28// to arrays,for use in passing to functions as dimensioned.
29// Add two additional optional(!) template args for piN/NN.
30// Add new data member "sum" to separate summed xsec values
31// from measured inclusive (tot) cross-sections. Add two
32// ctors to pass inclusive xsec array as input (for piN/NN).
33// 20100611 M. Kelsey -- Work around Intel ICC compiler warning about
34// index[] subscripts out of range. Dimension to full [9].
35// 20100803 M. Kelsey -- Add printing function for debugging, split
36// implementation code to .icc file. Add name argument.
37// 20110718 M. Kelsey -- Add inelastic cross-section sum to deal with
38// suppressing elastic scattering off free nucleons (hydrogen)
39// 20110719 M. Kelsey -- Add ctor argument for two-body initial state
40// 20110725 M. Kelsey -- Save initial state as data member
41// 20110923 M. Kelsey -- Add optional ostream& argument to print() fns
42
43#ifndef G4_CASCADE_DATA_HH
44#define G4_CASCADE_DATA_HH
45
46#include "globals.hh"
47#include "G4CascadeSampler.hh" /* To get number of energy bins */
48#include "G4String.hh"
49
50
51template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8=0,int N9=0>
53{
54 // NOTE: Need access to N2 by value to initialize index array
55 enum { N02=N2, N23=N2+N3, N24=N23+N4, N25=N24+N5, N26=N25+N6, N27=N26+N7,
56 N28=N27+N8, N29=N28+N9 };
57
58 enum { N8D=N8?N8:1, N9D=N9?N9:1 }; // SPECIAL: Can't dimension arrays [0]
59
60 enum { NM=N9?8:N8?7:6, NXS=N29 }; // Multiplicity and cross-section bins
61
62 G4int index[9]; // Start and stop indices to xsec's
63 G4double multiplicities[NM][NE]; // Multiplicity distributions
64
65 const G4int (&x2bfs)[N2][2]; // Initialized from file-scope inputs
66 const G4int (&x3bfs)[N3][3];
67 const G4int (&x4bfs)[N4][4];
68 const G4int (&x5bfs)[N5][5];
69 const G4int (&x6bfs)[N6][6];
70 const G4int (&x7bfs)[N7][7];
71 const G4int (&x8bfs)[N8D][8]; // These may not be used if mult==7
72 const G4int (&x9bfs)[N9D][9];
74
75 G4double sum[NE]; // Summed cross-sections, computed
76 const G4double (&tot)[NE]; // Inclusive cross-sections (from input)
77
78 G4double inelastic[NE]; // Sum of only inelastic channels
79
80 static const G4int empty8bfs[1][8]; // For multiplicity==7 case
81 static const G4int empty9bfs[1][9];
82
83 const G4String name; // For diagnostic purposes
84 const G4int initialState; // For registration in lookup table
85
86 G4int maxMultiplicity() const { return NM+1; } // Used by G4CascadeFunctions
87
88 // Dump multiplicty tables to specified stream
89 void print(std::ostream& os=G4cout) const;
90 void print(G4int mult, std::ostream& os) const;
91 void printXsec(const G4double (&xsec)[NE], std::ostream& os) const;
92
93 // Constructor for kaon/hyperon channels, with multiplicity <= 7
94 G4CascadeData(const G4int (&the2bfs)[N2][2], const G4int (&the3bfs)[N3][3],
95 const G4int (&the4bfs)[N4][4], const G4int (&the5bfs)[N5][5],
96 const G4int (&the6bfs)[N6][6], const G4int (&the7bfs)[N7][7],
97 const G4double (&xsec)[NXS][NE], G4int ini,
98 const G4String& aName="G4CascadeData")
99 : x2bfs(the2bfs), x3bfs(the3bfs), x4bfs(the4bfs), x5bfs(the5bfs),
100 x6bfs(the6bfs), x7bfs(the7bfs), x8bfs(empty8bfs), x9bfs(empty9bfs),
101 crossSections(xsec), tot(sum), name(aName), initialState(ini) {
102 initialize();
103 }
104
105 // Constructor for kaon/hyperon channels, with multiplicity <= 7 and inclusive
106 G4CascadeData(const G4int (&the2bfs)[N2][2], const G4int (&the3bfs)[N3][3],
107 const G4int (&the4bfs)[N4][4], const G4int (&the5bfs)[N5][5],
108 const G4int (&the6bfs)[N6][6], const G4int (&the7bfs)[N7][7],
109 const G4double (&xsec)[NXS][NE], const G4double (&theTot)[NE],
110 G4int ini, const G4String& aName="G4CascadeData")
111 : x2bfs(the2bfs), x3bfs(the3bfs), x4bfs(the4bfs), x5bfs(the5bfs),
112 x6bfs(the6bfs), x7bfs(the7bfs), x8bfs(empty8bfs), x9bfs(empty9bfs),
113 crossSections(xsec), tot(theTot), name(aName), initialState(ini) {
114 initialize();
115 }
116
117 // Constructor for pion/nucleon channels, with multiplicity > 7
118 G4CascadeData(const G4int (&the2bfs)[N2][2], const G4int (&the3bfs)[N3][3],
119 const G4int (&the4bfs)[N4][4], const G4int (&the5bfs)[N5][5],
120 const G4int (&the6bfs)[N6][6], const G4int (&the7bfs)[N7][7],
121 const G4int (&the8bfs)[N8D][8], const G4int (&the9bfs)[N9D][9],
122 const G4double (&xsec)[NXS][NE], G4int ini,
123 const G4String& aName="G4CascadeData")
124 : x2bfs(the2bfs), x3bfs(the3bfs), x4bfs(the4bfs), x5bfs(the5bfs),
125 x6bfs(the6bfs), x7bfs(the7bfs), x8bfs(the8bfs), x9bfs(the9bfs),
126 crossSections(xsec), tot(sum), name(aName), initialState(ini) {
127 initialize();
128 }
129
130 // Constructor for pion/nucleon channels, with multiplicity > 7 and inclusive
131 G4CascadeData(const G4int (&the2bfs)[N2][2], const G4int (&the3bfs)[N3][3],
132 const G4int (&the4bfs)[N4][4], const G4int (&the5bfs)[N5][5],
133 const G4int (&the6bfs)[N6][6], const G4int (&the7bfs)[N7][7],
134 const G4int (&the8bfs)[N8D][8], const G4int (&the9bfs)[N9D][9],
135 const G4double (&xsec)[NXS][NE], const G4double (&theTot)[NE],
136 G4int ini, const G4String& aName="G4CascadeData")
137 : x2bfs(the2bfs), x3bfs(the3bfs), x4bfs(the4bfs), x5bfs(the5bfs),
138 x6bfs(the6bfs), x7bfs(the7bfs), x8bfs(the8bfs), x9bfs(the9bfs),
139 crossSections(xsec), tot(theTot), name(aName), initialState(ini) {
140 initialize();
141 }
142
143 void initialize(); // Fill summed arrays from input
144};
145
146// Dummy arrays for use when optional template arguments are skipped
147template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
149
150template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
152
153// GCC and other compilers require template implementations here
154#include "G4CascadeData.icc"
155
156#endif
@ NE
Definition: Evaluator.cc:65
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
G4CascadeData(const G4int(&the2bfs)[N2][2], const G4int(&the3bfs)[N3][3], const G4int(&the4bfs)[N4][4], const G4int(&the5bfs)[N5][5], const G4int(&the6bfs)[N6][6], const G4int(&the7bfs)[N7][7], const G4int(&the8bfs)[N8D][8], const G4int(&the9bfs)[N9D][9], const G4double(&xsec)[NXS][NE], const G4double(&theTot)[NE], G4int ini, const G4String &aName="G4CascadeData")
G4double sum[NE]
static const G4int empty8bfs[1][8]
G4CascadeData(const G4int(&the2bfs)[N2][2], const G4int(&the3bfs)[N3][3], const G4int(&the4bfs)[N4][4], const G4int(&the5bfs)[N5][5], const G4int(&the6bfs)[N6][6], const G4int(&the7bfs)[N7][7], const G4double(&xsec)[NXS][NE], G4int ini, const G4String &aName="G4CascadeData")
void print(std::ostream &os=G4cout) const
void initialize()
G4int maxMultiplicity() const
const G4int(& x9bfs)[N9D][9]
G4CascadeData(const G4int(&the2bfs)[N2][2], const G4int(&the3bfs)[N3][3], const G4int(&the4bfs)[N4][4], const G4int(&the5bfs)[N5][5], const G4int(&the6bfs)[N6][6], const G4int(&the7bfs)[N7][7], const G4int(&the8bfs)[N8D][8], const G4int(&the9bfs)[N9D][9], const G4double(&xsec)[NXS][NE], G4int ini, const G4String &aName="G4CascadeData")
G4double inelastic[NE]
const G4int initialState
G4double multiplicities[NM][NE]
G4CascadeData(const G4int(&the2bfs)[N2][2], const G4int(&the3bfs)[N3][3], const G4int(&the4bfs)[N4][4], const G4int(&the5bfs)[N5][5], const G4int(&the6bfs)[N6][6], const G4int(&the7bfs)[N7][7], const G4double(&xsec)[NXS][NE], const G4double(&theTot)[NE], G4int ini, const G4String &aName="G4CascadeData")
const G4double(& crossSections)[NXS][NE]
const G4int(& x5bfs)[N5][5]
const G4int(& x7bfs)[N7][7]
const G4double(& tot)[NE]
void printXsec(const G4double(&xsec)[NE], std::ostream &os) const
const G4int(& x3bfs)[N3][3]
void print(G4int mult, std::ostream &os) const
G4int index[9]
const G4String name
static const G4int empty9bfs[1][9]
const G4int(& x8bfs)[N8D][8]
const G4int(& x6bfs)[N6][6]
const G4int(& x2bfs)[N2][2]
const G4int(& x4bfs)[N4][4]