2// ********************************************************************
3// * License and Disclaimer *
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. *
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. *
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// ********************************************************************
27// 20100803 M. Kelsey -- Move implementations to this .icc file. Use name
28// string to report output
29// 20110718 M. Kelsey -- Add inelastic cross-section sum to deal with
30// suppressing elastic scattering off free nucleons (hydrogen)
31// 20110719 M. Kelsey -- Add ctor argument for two-body initial state
32// 20110725 M. Kelsey -- Save initial state as data member
33// 20110728 M. Kelsey -- Fix Coverity #22954, recursive #include.
34// 20110923 M. Kelsey -- Add optional ostream& argument to print() fns,
35// drop all "inline" keywords
36// 20120608 M. Kelsey -- Fix variable-name "shadowing" compiler warnings.
37// 20130627 M. Kelsey -- Use new function to print particle name strings.
39#ifndef G4_CASCADE_DATA_ICC
40#define G4_CASCADE_DATA_ICC
42#include "G4InuclParticleNames.hh"
47// Fill cumulative cross-section arrays from input data
48template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
49void G4CascadeData<NE,N2,N3,N4,N5,N6,N7,N8,N9>::initialize() {
50 // Initialize index offsets for cross-section array (can't do globally)
51 index[0] = 0; index[1] = N02; index[2] = N23; index[3] = N24;
52 index[4] = N25; index[5] = N26; index[6] = N27; index[7] = N28;
55 // Initialize multiplicity array
56 for (G4int im = 0; im < NM; im++) {
57 G4int start = index[im];
58 G4int stop = index[im+1];
59 for (G4int k = 0; k < NE; k++) {
60 multiplicities[im][k] = 0.0;
61 for (G4int i = start; i < stop; i++) {
62 multiplicities[im][k] += crossSections[i][k];
67 // Initialize total cross section arrays
68 for (G4int k = 0; k < NE; k++) {
70 for (G4int im = 0; im < NM; im++) {
71 sum[k] += multiplicities[im][k];
75 // Identify elastic scattering channel and subtract from inclusive
77 for (i2b=index[0]; i2b<index[1]; i2b++) {
78 if (x2bfs[i2b][0]*x2bfs[i2b][1] == initialState) break; // Found it
81 for (G4int k = 0; k < NE; k++) {
82 if (i2b<index[1]) inelastic[k] = tot[k] - crossSections[i2b][k];
83 else inelastic[k] = tot[k]; // FIXME: No elastic channel in table!
88// Dump complete final state tables, all multiplicities
89template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
90void G4CascadeData<NE,N2,N3,N4,N5,N6,N7,N8,N9>::print(std::ostream& os) const {
91 os << "\n " << name << " Total cross section:" << G4endl;
93 os << "\n Summed cross section:" << G4endl;
95 os << "\n Inelastic cross section:" << G4endl;
96 printXsec(inelastic, os);
97 os << "\n Individual channel cross sections" << G4endl;
99 for (int im=2; im<NM+2; im++) print(im, os);
103// Dump tables for specified multiplicity
104template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
105void G4CascadeData<NE,N2,N3,N4,N5,N6,N7,N8,N9>::
106print(G4int mult, std::ostream& os) const {
107 if (mult < 0) { // Old interface used mult == -1 for all
112 G4int im = mult-2; // Convert multiplicity to array index
114 G4int start = index[im];
115 G4int stop = index[im+1];
116 os << "\n Mulitplicity " << mult << " (indices " << start << " to "
117 << stop-1 << ") summed cross section:" << G4endl;
119 printXsec(multiplicities[im], os);
121 for (G4int i=start; i<stop; i++) {
123 os << "\n final state x" << mult << "bfs[" << ichan << "] : ";
124 for (G4int fsi=0; fsi<mult; fsi++) {
126 case 2: os << " " << G4InuclParticleNames::nameShort(x2bfs[ichan][fsi]);
128 case 3: os << " " << G4InuclParticleNames::nameShort(x3bfs[ichan][fsi]);
130 case 4: os << " " << G4InuclParticleNames::nameShort(x4bfs[ichan][fsi]);
132 case 5: os << " " << G4InuclParticleNames::nameShort(x5bfs[ichan][fsi]);
134 case 6: os << " " << G4InuclParticleNames::nameShort(x6bfs[ichan][fsi]);
136 case 7: os << " " << G4InuclParticleNames::nameShort(x7bfs[ichan][fsi]);
138 case 8: os << " " << G4InuclParticleNames::nameShort(x8bfs[ichan][fsi]);
140 case 9: os << " " << G4InuclParticleNames::nameShort(x9bfs[ichan][fsi]);
145 os << " -- cross section [" << i << "]:" << G4endl;
146 printXsec(crossSections[i], os);
150// Dump individual cross-section table, two lines of 12 values
151template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
152void G4CascadeData<NE,N2,N3,N4,N5,N6,N7,N8,N9>::
153printXsec(const G4double (&xsec)[NE], std::ostream& os) const {
154 for (G4int k=0; k<NE; k++) {
155 os << " " << std::setw(6) << xsec[k];
156 if ((k+1)%10 == 0) os << G4endl;
161#endif /* G4_CASCADE_DATA_ICC */