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