source: trunk/source/processes/hadronic/models/cascade/cascade/include/G4CascadeData.hh @ 1340

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

update ti head

File size: 6.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// $Id: G4CascadeData.hh,v 1.10 2010/08/04 05:28:24 mkelsey Exp $
27// GEANT4 tag: $Name: hadr-casc-V09-03-85 $
28//
29// 20100507  M. Kelsey -- Use template arguments to dimension const-refs
30//              to arrays,for use in passing to functions as dimensioned.
31//              Add two additional optional(!) template args for piN/NN.
32//              Add new data member "sum" to separate summed xsec values
33//              from measured inclusive (tot) cross-sections.  Add two
34//              ctors to pass inclusive xsec array as input (for piN/NN).
35// 20100611  M. Kelsey -- Work around Intel ICC compiler warning about
36//              index[] subscripts out of range.  Dimension to full [9].
37// 20100803  M. Kelsey -- Add printing function for debugging, split
38//              implementation code to .icc file.  Add name argument.
39
40#ifndef G4_CASCADE_DATA_HH
41#define G4_CASCADE_DATA_HH
42
43#include "globals.hh"
44#include "G4CascadeSampler.hh"          /* To get number of energy bins */
45#include "G4String.hh"
46
47
48template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8=0,int N9=0>
49struct G4CascadeData
50{
51  // NOTE: Need access to N2 by value to initialize index array
52  enum { N02=N2, N23=N2+N3, N24=N23+N4, N25=N24+N5, N26=N25+N6, N27=N26+N7,
53         N28=N27+N8, N29=N28+N9 };
54
55  enum { N8D=N8?N8:1, N9D=N9?N9:1 };    // SPECIAL: Can't dimension arrays [0]
56
57  enum { NM=N9?8:N8?7:6, NXS=N29 };     // Multiplicity and cross-section bins
58
59  G4int index[9];                       // Start and stop indices to xsec's
60  G4double multiplicities[NM][NE];      // Multiplicity distributions
61
62  const G4int (&x2bfs)[N2][2];          // Initialized from file-scope inputs
63  const G4int (&x3bfs)[N3][3];
64  const G4int (&x4bfs)[N4][4];
65  const G4int (&x5bfs)[N5][5];
66  const G4int (&x6bfs)[N6][6];
67  const G4int (&x7bfs)[N7][7];
68  const G4int (&x8bfs)[N8D][8];         // These may not be used if mult==7
69  const G4int (&x9bfs)[N9D][9];
70  const G4double (&crossSections)[NXS][NE];
71
72  G4double sum[NE];                     // Summed cross-sections, computed
73  const G4double (&tot)[NE];            // Inclusive cross-sections (from input)
74
75  static const G4int empty8bfs[1][8];   // For multiplicity==7 case
76  static const G4int empty9bfs[1][9];
77
78  const G4String name;                  // For diagnostic purposes
79
80  G4int maxMultiplicity() const { return NM+1; }  // Used by G4CascadeFunctions
81
82  void print(G4int mult=-1) const;      // Dump multiplicty tables (-1 == all)
83  void printXsec(const G4double (&xsec)[NE]) const;
84
85  // Constructor for kaon/hyperon channels, with multiplicity <= 7
86  G4CascadeData(const G4int (&the2bfs)[N2][2], const G4int (&the3bfs)[N3][3],
87                const G4int (&the4bfs)[N4][4], const G4int (&the5bfs)[N5][5],
88                const G4int (&the6bfs)[N6][6], const G4int (&the7bfs)[N7][7],
89                const G4double (&xsec)[NXS][NE],
90                const G4String& aName="G4CascadeData")
91    : x2bfs(the2bfs), x3bfs(the3bfs), x4bfs(the4bfs), x5bfs(the5bfs),
92      x6bfs(the6bfs), x7bfs(the7bfs), x8bfs(empty8bfs), x9bfs(empty9bfs),
93      crossSections(xsec), tot(sum), name(aName) { initialize(); }
94
95  // Constructor for kaon/hyperon channels, with multiplicity <= 7 and inclusive
96  G4CascadeData(const G4int (&the2bfs)[N2][2], const G4int (&the3bfs)[N3][3],
97                const G4int (&the4bfs)[N4][4], const G4int (&the5bfs)[N5][5],
98                const G4int (&the6bfs)[N6][6], const G4int (&the7bfs)[N7][7],
99                const G4double (&xsec)[NXS][NE], const G4double (&theTot)[NE],
100                const G4String& aName="G4CascadeData")
101    : x2bfs(the2bfs), x3bfs(the3bfs), x4bfs(the4bfs), x5bfs(the5bfs),
102      x6bfs(the6bfs), x7bfs(the7bfs), x8bfs(empty8bfs), x9bfs(empty9bfs),
103      crossSections(xsec), tot(theTot), name(aName) { initialize(); }
104
105  // Constructor for pion/nuleon channels, with multiplicity > 7
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 G4int (&the8bfs)[N8D][8], const G4int (&the9bfs)[N9D][9],
110                const G4double (&xsec)[NXS][NE],
111                const G4String& aName="G4CascadeData")
112    : x2bfs(the2bfs), x3bfs(the3bfs), x4bfs(the4bfs), x5bfs(the5bfs),
113      x6bfs(the6bfs), x7bfs(the7bfs), x8bfs(the8bfs), x9bfs(the9bfs),
114      crossSections(xsec), tot(sum), name(aName) { initialize(); }
115
116  // Constructor for pion/nuleon channels, with multiplicity > 7 and inclusive
117  G4CascadeData(const G4int (&the2bfs)[N2][2], const G4int (&the3bfs)[N3][3],
118                const G4int (&the4bfs)[N4][4], const G4int (&the5bfs)[N5][5],
119                const G4int (&the6bfs)[N6][6], const G4int (&the7bfs)[N7][7],
120                const G4int (&the8bfs)[N8D][8], const G4int (&the9bfs)[N9D][9],
121                const G4double (&xsec)[NXS][NE], const G4double (&theTot)[NE],
122                const G4String& aName="G4CascadeData")
123    : x2bfs(the2bfs), x3bfs(the3bfs), x4bfs(the4bfs), x5bfs(the5bfs),
124      x6bfs(the6bfs), x7bfs(the7bfs), x8bfs(the8bfs), x9bfs(the9bfs),
125      crossSections(xsec), tot(theTot), name(aName) { initialize(); }
126
127  void initialize();                    // Fill summed arrays from input
128};
129
130// Dummy arrays for use when optional template arguments are skipped
131template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
132const G4int G4CascadeData<NE,N2,N3,N4,N5,N6,N7,N8,N9>::empty8bfs[1][8] = {{0}};
133
134template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
135const G4int G4CascadeData<NE,N2,N3,N4,N5,N6,N7,N8,N9>::empty9bfs[1][9] = {{0}};
136
137// GCC and other compilers require template implementations here
138#include "G4CascadeData.icc"
139
140#endif
Note: See TracBrowser for help on using the repository browser.