G4CascadeData.icc

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 // 20100803  M. Kelsey -- Move implementations to this .icc file.  Use name
00029 //              string to report output
00030 // 20110718  M. Kelsey -- Add inelastic cross-section sum to deal with
00031 //              suppressing elastic scattering off free nucleons (hydrogen)
00032 // 20110719  M. Kelsey -- Add ctor argument for two-body initial state
00033 // 20110725  M. Kelsey -- Save initial state as data member
00034 // 20110728  M. Kelsey -- Fix Coverity #22954, recursive #include.
00035 // 20110923  M. Kelsey -- Add optional ostream& argument to print() fns,
00036 //              drop all "inline" keywords
00037 // 20120608  M. Kelsey -- Fix variable-name "shadowing" compiler warnings.
00038 
00039 #ifndef G4_CASCADE_DATA_ICC
00040 #define G4_CASCADE_DATA_ICC
00041 
00042 #include <iostream>
00043 #include <iomanip>
00044 
00045 
00046 // Fill cumulative cross-section arrays from input data
00047 template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
00048 void G4CascadeData<NE,N2,N3,N4,N5,N6,N7,N8,N9>::initialize() {
00049   // Initialize index offsets for cross-section array (can't do globally)
00050   index[0] = 0;   index[1] = N02; index[2] = N23; index[3] = N24;
00051   index[4] = N25; index[5] = N26; index[6] = N27; index[7] = N28;
00052   index[8] = N29;
00053 
00054   // Initialize multiplicity array
00055   for (G4int im = 0; im < NM; im++) {
00056     G4int start = index[im];
00057     G4int stop = index[im+1];
00058     for (G4int k = 0; k < NE; k++) {
00059       multiplicities[im][k] = 0.0;
00060       for (G4int i = start; i < stop; i++) {
00061         multiplicities[im][k] += crossSections[i][k];
00062       }
00063     }
00064   }
00065 
00066   // Initialize total cross section arrays
00067   for (G4int k = 0; k < NE; k++) {
00068     sum[k] = 0.0;
00069     for (G4int im = 0; im < NM; im++) {
00070       sum[k] += multiplicities[im][k];
00071     }
00072   }
00073 
00074   // Identify elastic scattering channel and subtract from inclusive
00075   G4int i2b = 0;
00076   for (i2b=index[0]; i2b<index[1]; i2b++) {
00077     if (x2bfs[i2b][0]*x2bfs[i2b][1] == initialState) break;     // Found it
00078   }
00079 
00080   for (G4int k = 0; k < NE; k++) {
00081     if (i2b<index[1]) inelastic[k] = tot[k] - crossSections[i2b][k];
00082     else inelastic[k] = tot[k];         // FIXME:  No elastic channel in table!
00083   }
00084 }
00085 
00086 
00087 // Dump complete final state tables, all multiplicities
00088 template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
00089 void G4CascadeData<NE,N2,N3,N4,N5,N6,N7,N8,N9>::print(std::ostream& os) const {
00090   os << "\n " << name << " Total cross section:" << G4endl;
00091   printXsec(tot, os);
00092   os << "\n Summed cross section:" << G4endl;
00093   printXsec(sum, os);
00094   os << "\n Inelastic cross section:" << G4endl;
00095   printXsec(inelastic, os);
00096   os << "\n Individual channel cross sections" << G4endl;
00097   
00098   for (int im=2; im<NM+2; im++) print(im, os);
00099   return;
00100 }
00101 
00102 // Dump tables for specified multiplicity
00103 template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
00104 void G4CascadeData<NE,N2,N3,N4,N5,N6,N7,N8,N9>::
00105 print(G4int mult, std::ostream& os) const {
00106   if (mult < 0) {               // Old interface used mult == -1 for all
00107     print(os);
00108     return;
00109   }
00110 
00111   G4int im = mult-2;            // Convert multiplicity to array index
00112 
00113   G4int start = index[im];
00114   G4int stop = index[im+1];
00115   os << "\n Mulitplicity " << mult << " (indices " << start << " to "
00116          << stop-1 << ") summed cross section:" << G4endl;
00117 
00118   printXsec(multiplicities[im], os);
00119   
00120   for (G4int i=start; i<stop; i++) {
00121     G4int ichan=i-start;
00122     os << "\n final state x" << mult << "bfs[" << ichan << "] : ";
00123     for (G4int fsi=0; fsi<mult; fsi++) {
00124       switch (mult) {
00125       case 2: os << " " << x2bfs[ichan][fsi]; break;
00126       case 3: os << " " << x3bfs[ichan][fsi]; break;
00127       case 4: os << " " << x4bfs[ichan][fsi]; break;
00128       case 5: os << " " << x5bfs[ichan][fsi]; break;
00129       case 6: os << " " << x6bfs[ichan][fsi]; break;
00130       case 7: os << " " << x7bfs[ichan][fsi]; break;
00131       case 8: os << " " << x8bfs[ichan][fsi]; break;
00132       case 9: os << " " << x9bfs[ichan][fsi]; break;
00133       default: ;
00134       }
00135     }
00136     os << " -- cross section [" << i << "]:" << G4endl;
00137     printXsec(crossSections[i], os);
00138   }
00139 }
00140 
00141 // Dump individual cross-section table, two lines of 12 values
00142 template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
00143 void G4CascadeData<NE,N2,N3,N4,N5,N6,N7,N8,N9>::
00144 printXsec(const G4double (&xsec)[NE], std::ostream& os) const {
00145   for (G4int k=0; k<NE; k++) {
00146     os << std::setw(6) << xsec[k];              // Use sign-gap as separator
00147     if ((k+1)%12 == 0) os << G4endl;
00148   }
00149   os << G4endl;
00150 }
00151 
00152 #endif  /* G4_CASCADE_DATA_ICC */

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