|
||||
File indexing completed on 2025-01-18 09:58:46
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 // Author: Dennis Wright (SLAC) 0027 // Date: 28 January 2013 0028 // 0029 // Description: Templated base class for numerically integrated two-body 0030 // final state angular distributions in Bertini-style cascade 0031 // 0032 // 20130219 M. Kelsey: Rewrite with C-arrays, using template args for 0033 // dimensions (c.f. G4CascadeSampler) 0034 // 20130620 Address Coverity #51342; initialize angDist[] buffer in ctor 0035 0036 #ifndef G4NumIntTwoBodyAngDst_h 0037 #define G4NumIntTwoBodyAngDst_h 1 0038 0039 #include "G4VTwoBodyAngDst.hh" 0040 #include <algorithm> 0041 0042 template <G4int NKEBINS, G4int NANGLES> 0043 class G4NumIntTwoBodyAngDst : public G4VTwoBodyAngDst { 0044 public: 0045 enum { nDists=NKEBINS, nAngles=NANGLES }; // For use in function arguments 0046 0047 G4NumIntTwoBodyAngDst(const G4String& name, 0048 const G4double (&kebins)[nDists], 0049 const G4double (&angles)[nAngles], 0050 const G4double (&dists)[nDists][nAngles], 0051 const G4double highKEscale, G4int verbose = 0) 0052 : G4VTwoBodyAngDst(name, verbose), tcoeff(highKEscale), 0053 labKE(kebins), cosBins(angles), angDists(dists) { 0054 std::fill(angDist, angDist+nAngles, 0.); // Initialize working buffer 0055 } 0056 0057 virtual ~G4NumIntTwoBodyAngDst() {;} 0058 0059 virtual G4double GetCosTheta(const G4double& ekin, const G4double& pcm) const; 0060 0061 protected: 0062 G4double tcoeff; // Fall-off of exponential for high energy ang. dist. 0063 0064 // Kinetic energies at which angular distributions are taken 0065 const G4double (&labKE)[nDists]; 0066 0067 // Bins of the angular distributions in cos(theta) 0068 const G4double (&cosBins)[nAngles]; 0069 0070 // table of numerical normalized indefinite integrals of 0071 // angular distributions vs. KE and angle 0072 const G4double (&angDists)[nDists][nAngles]; 0073 0074 // Compute interpolated angular distribution between energy bins 0075 void Interpolate(const G4double& ekin) const; 0076 mutable G4double angDist[nAngles]; // Reusable buffer 0077 }; 0078 0079 // Implementations must be included for templated classes 0080 #include "G4NumIntTwoBodyAngDst.icc" 0081 0082 #endif /* G4NumIntTwoBodyAngDst_h */
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |