Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:58:00

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 //
0027 // 20100803  M. Kelsey -- Move implementations to this .icc file.  Use name
0028 //      string to report output
0029 // 20110718  M. Kelsey -- Add inelastic cross-section sum to deal with
0030 //      suppressing elastic scattering off free nucleons (hydrogen)
0031 // 20110719  M. Kelsey -- Add ctor argument for two-body initial state
0032 // 20110725  M. Kelsey -- Save initial state as data member
0033 // 20110728  M. Kelsey -- Fix Coverity #22954, recursive #include.
0034 // 20110923  M. Kelsey -- Add optional ostream& argument to print() fns,
0035 //      drop all "inline" keywords
0036 // 20120608  M. Kelsey -- Fix variable-name "shadowing" compiler warnings.
0037 // 20130627  M. Kelsey -- Use new function to print particle name strings.
0038 
0039 #ifndef G4_CASCADE_DATA_ICC
0040 #define G4_CASCADE_DATA_ICC
0041 
0042 #include "G4InuclParticleNames.hh"
0043 #include <iostream>
0044 #include <iomanip>
0045 
0046 
0047 // Fill cumulative cross-section arrays from input data
0048 template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
0049 void G4CascadeData<NE,N2,N3,N4,N5,N6,N7,N8,N9>::initialize() {
0050   // Initialize index offsets for cross-section array (can't do globally)
0051   index[0] = 0;   index[1] = N02; index[2] = N23; index[3] = N24;
0052   index[4] = N25; index[5] = N26; index[6] = N27; index[7] = N28;
0053   index[8] = N29;
0054 
0055   // Initialize multiplicity array
0056   for (G4int im = 0; im < NM; im++) {
0057     G4int start = index[im];
0058     G4int stop = index[im+1];
0059     for (G4int k = 0; k < NE; k++) {
0060       multiplicities[im][k] = 0.0;
0061       for (G4int i = start; i < stop; i++) {
0062     multiplicities[im][k] += crossSections[i][k];
0063       }
0064     }
0065   }
0066 
0067   // Initialize total cross section arrays
0068   for (G4int k = 0; k < NE; k++) {
0069     sum[k] = 0.0;
0070     for (G4int im = 0; im < NM; im++) {
0071       sum[k] += multiplicities[im][k];
0072     }
0073   }
0074 
0075   // Identify elastic scattering channel and subtract from inclusive
0076   G4int i2b = 0;
0077   for (i2b=index[0]; i2b<index[1]; i2b++) {
0078     if (x2bfs[i2b][0]*x2bfs[i2b][1] == initialState) break; // Found it
0079   }
0080 
0081   for (G4int k = 0; k < NE; k++) {
0082     if (i2b<index[1]) inelastic[k] = tot[k] - crossSections[i2b][k];
0083     else inelastic[k] = tot[k];     // FIXME:  No elastic channel in table!
0084   }
0085 }
0086 
0087 
0088 // Dump complete final state tables, all multiplicities
0089 template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
0090 void G4CascadeData<NE,N2,N3,N4,N5,N6,N7,N8,N9>::print(std::ostream& os) const {
0091   os << "\n " << name << " Total cross section:" << G4endl;
0092   printXsec(tot, os);
0093   os << "\n Summed cross section:" << G4endl;
0094   printXsec(sum, os);
0095   os << "\n Inelastic cross section:" << G4endl;
0096   printXsec(inelastic, os);
0097   os << "\n Individual channel cross sections" << G4endl;
0098   
0099   for (int im=2; im<NM+2; im++) print(im, os);
0100   return;
0101 }
0102 
0103 // Dump tables for specified multiplicity
0104 template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
0105 void G4CascadeData<NE,N2,N3,N4,N5,N6,N7,N8,N9>::
0106 print(G4int mult, std::ostream& os) const {
0107   if (mult < 0) {       // Old interface used mult == -1 for all
0108     print(os);
0109     return;
0110   }
0111 
0112   G4int im = mult-2;        // Convert multiplicity to array index
0113 
0114   G4int start = index[im];
0115   G4int stop = index[im+1];
0116   os << "\n Mulitplicity " << mult << " (indices " << start << " to "
0117      << stop-1 << ") summed cross section:" << G4endl;
0118 
0119   printXsec(multiplicities[im], os);
0120   
0121   for (G4int i=start; i<stop; i++) {
0122     G4int ichan=i-start;
0123     os << "\n final state x" << mult << "bfs[" << ichan << "] : ";
0124     for (G4int fsi=0; fsi<mult; fsi++) {
0125       switch (mult) {
0126       case 2: os << " " << G4InuclParticleNames::nameShort(x2bfs[ichan][fsi]);
0127     break;
0128       case 3: os << " " << G4InuclParticleNames::nameShort(x3bfs[ichan][fsi]);
0129     break;
0130       case 4: os << " " << G4InuclParticleNames::nameShort(x4bfs[ichan][fsi]);
0131     break;
0132       case 5: os << " " << G4InuclParticleNames::nameShort(x5bfs[ichan][fsi]);
0133     break;
0134       case 6: os << " " << G4InuclParticleNames::nameShort(x6bfs[ichan][fsi]);
0135     break;
0136       case 7: os << " " << G4InuclParticleNames::nameShort(x7bfs[ichan][fsi]);
0137     break;
0138       case 8: os << " " << G4InuclParticleNames::nameShort(x8bfs[ichan][fsi]);
0139     break;
0140       case 9: os << " " << G4InuclParticleNames::nameShort(x9bfs[ichan][fsi]);
0141     break;
0142       default: ;
0143       }
0144     }
0145     os << " -- cross section [" << i << "]:" << G4endl;
0146     printXsec(crossSections[i], os);
0147   }
0148 }
0149 
0150 // Dump individual cross-section table, two lines of 12 values
0151 template <int NE,int N2,int N3,int N4,int N5,int N6,int N7,int N8,int N9>
0152 void G4CascadeData<NE,N2,N3,N4,N5,N6,N7,N8,N9>::
0153 printXsec(const G4double (&xsec)[NE], std::ostream& os) const {
0154   for (G4int k=0; k<NE; k++) {
0155     os << " " << std::setw(6) << xsec[k];
0156     if ((k+1)%10 == 0) os << G4endl;
0157   }
0158   os << G4endl;
0159 }
0160 
0161 #endif  /* G4_CASCADE_DATA_ICC */