G4CascadeData.hh

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id$
00027 //
00028 // 20100507  M. Kelsey -- Use template arguments to dimension const-refs
00029 //              to arrays,for use in passing to functions as dimensioned.
00030 //              Add two additional optional(!) template args for piN/NN.
00031 //              Add new data member "sum" to separate summed xsec values
00032 //              from measured inclusive (tot) cross-sections.  Add two
00033 //              ctors to pass inclusive xsec array as input (for piN/NN).
00034 // 20100611  M. Kelsey -- Work around Intel ICC compiler warning about
00035 //              index[] subscripts out of range.  Dimension to full [9].
00036 // 20100803  M. Kelsey -- Add printing function for debugging, split
00037 //              implementation code to .icc file.  Add name argument.
00038 // 20110718  M. Kelsey -- Add inelastic cross-section sum to deal with
00039 //              suppressing elastic scattering off free nucleons (hydrogen)
00040 // 20110719  M. Kelsey -- Add ctor argument for two-body initial state
00041 // 20110725  M. Kelsey -- Save initial state as data member
00042 // 20110923  M. Kelsey -- Add optional ostream& argument to print() fns
00043 
00044 #ifndef G4_CASCADE_DATA_HH
00045 #define G4_CASCADE_DATA_HH
00046 
00047 #include "globals.hh"
00048 #include "G4CascadeSampler.hh"          /* To get number of energy bins */
00049 #include "G4String.hh"
00050 
00051 
00052 template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8=0,int N9=0>
00053 struct G4CascadeData
00054 {
00055   // NOTE: Need access to N2 by value to initialize index array
00056   enum { N02=N2, N23=N2+N3, N24=N23+N4, N25=N24+N5, N26=N25+N6, N27=N26+N7,
00057          N28=N27+N8, N29=N28+N9 };
00058 
00059   enum { N8D=N8?N8:1, N9D=N9?N9:1 };    // SPECIAL: Can't dimension arrays [0]
00060 
00061   enum { NM=N9?8:N8?7:6, NXS=N29 };     // Multiplicity and cross-section bins
00062 
00063   G4int index[9];                       // Start and stop indices to xsec's
00064   G4double multiplicities[NM][NE];      // Multiplicity distributions
00065 
00066   const G4int (&x2bfs)[N2][2];          // Initialized from file-scope inputs
00067   const G4int (&x3bfs)[N3][3];
00068   const G4int (&x4bfs)[N4][4];
00069   const G4int (&x5bfs)[N5][5];
00070   const G4int (&x6bfs)[N6][6];
00071   const G4int (&x7bfs)[N7][7];
00072   const G4int (&x8bfs)[N8D][8];         // These may not be used if mult==7
00073   const G4int (&x9bfs)[N9D][9];
00074   const G4double (&crossSections)[NXS][NE];
00075 
00076   G4double sum[NE];                     // Summed cross-sections, computed
00077   const G4double (&tot)[NE];            // Inclusive cross-sections (from input)
00078 
00079   G4double inelastic[NE];               // Sum of only inelastic channels
00080 
00081   static const G4int empty8bfs[1][8];   // For multiplicity==7 case
00082   static const G4int empty9bfs[1][9];
00083 
00084   const G4String name;                  // For diagnostic purposes
00085   const G4int initialState;             // For registration in lookup table
00086 
00087   G4int maxMultiplicity() const { return NM+1; }  // Used by G4CascadeFunctions
00088 
00089   // Dump multiplicty tables to specified stream
00090   void print(std::ostream& os=G4cout) const;
00091   void print(G4int mult, std::ostream& os) const;
00092   void printXsec(const G4double (&xsec)[NE], std::ostream& os) const;
00093 
00094   // Constructor for kaon/hyperon channels, with multiplicity <= 7
00095   G4CascadeData(const G4int (&the2bfs)[N2][2], const G4int (&the3bfs)[N3][3],
00096                 const G4int (&the4bfs)[N4][4], const G4int (&the5bfs)[N5][5],
00097                 const G4int (&the6bfs)[N6][6], const G4int (&the7bfs)[N7][7],
00098                 const G4double (&xsec)[NXS][NE], G4int ini,
00099                 const G4String& aName="G4CascadeData")
00100     : x2bfs(the2bfs), x3bfs(the3bfs), x4bfs(the4bfs), x5bfs(the5bfs),
00101       x6bfs(the6bfs), x7bfs(the7bfs), x8bfs(empty8bfs), x9bfs(empty9bfs),
00102       crossSections(xsec), tot(sum), name(aName), initialState(ini) {
00103     initialize();
00104   }
00105 
00106   // Constructor for kaon/hyperon channels, with multiplicity <= 7 and inclusive
00107   G4CascadeData(const G4int (&the2bfs)[N2][2], const G4int (&the3bfs)[N3][3],
00108                 const G4int (&the4bfs)[N4][4], const G4int (&the5bfs)[N5][5],
00109                 const G4int (&the6bfs)[N6][6], const G4int (&the7bfs)[N7][7],
00110                 const G4double (&xsec)[NXS][NE], const G4double (&theTot)[NE],
00111                 G4int ini, const G4String& aName="G4CascadeData")
00112     : x2bfs(the2bfs), x3bfs(the3bfs), x4bfs(the4bfs), x5bfs(the5bfs),
00113       x6bfs(the6bfs), x7bfs(the7bfs), x8bfs(empty8bfs), x9bfs(empty9bfs),
00114       crossSections(xsec), tot(theTot), name(aName), initialState(ini) {
00115     initialize();
00116   }
00117 
00118   // Constructor for pion/nucleon channels, with multiplicity > 7
00119   G4CascadeData(const G4int (&the2bfs)[N2][2], const G4int (&the3bfs)[N3][3],
00120                 const G4int (&the4bfs)[N4][4], const G4int (&the5bfs)[N5][5],
00121                 const G4int (&the6bfs)[N6][6], const G4int (&the7bfs)[N7][7],
00122                 const G4int (&the8bfs)[N8D][8], const G4int (&the9bfs)[N9D][9],
00123                 const G4double (&xsec)[NXS][NE], G4int ini,
00124                 const G4String& aName="G4CascadeData")
00125     : x2bfs(the2bfs), x3bfs(the3bfs), x4bfs(the4bfs), x5bfs(the5bfs),
00126       x6bfs(the6bfs), x7bfs(the7bfs), x8bfs(the8bfs), x9bfs(the9bfs),
00127       crossSections(xsec), tot(sum), name(aName), initialState(ini) {
00128     initialize();
00129   }
00130 
00131   // Constructor for pion/nucleon channels, with multiplicity > 7 and inclusive
00132   G4CascadeData(const G4int (&the2bfs)[N2][2], const G4int (&the3bfs)[N3][3],
00133                 const G4int (&the4bfs)[N4][4], const G4int (&the5bfs)[N5][5],
00134                 const G4int (&the6bfs)[N6][6], const G4int (&the7bfs)[N7][7],
00135                 const G4int (&the8bfs)[N8D][8], const G4int (&the9bfs)[N9D][9],
00136                 const G4double (&xsec)[NXS][NE], const G4double (&theTot)[NE],
00137                 G4int ini, const G4String& aName="G4CascadeData")
00138     : x2bfs(the2bfs), x3bfs(the3bfs), x4bfs(the4bfs), x5bfs(the5bfs),
00139       x6bfs(the6bfs), x7bfs(the7bfs), x8bfs(the8bfs), x9bfs(the9bfs),
00140       crossSections(xsec), tot(theTot), name(aName), initialState(ini) {
00141     initialize();
00142   }
00143 
00144   void initialize();                    // Fill summed arrays from input
00145 };
00146 
00147 // Dummy arrays for use when optional template arguments are skipped
00148 template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
00149 const G4int G4CascadeData<NE,N2,N3,N4,N5,N6,N7,N8,N9>::empty8bfs[1][8] = {{0}};
00150 
00151 template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
00152 const G4int G4CascadeData<NE,N2,N3,N4,N5,N6,N7,N8,N9>::empty9bfs[1][9] = {{0}};
00153 
00154 // GCC and other compilers require template implementations here
00155 #include "G4CascadeData.icc"
00156 
00157 #endif

Generated on Mon May 27 17:47:49 2013 for Geant4 by  doxygen 1.4.7