|
||||
File indexing completed on 2025-01-18 09:58:50
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: Michael Kelsey (SLAC) 0027 // Date: 20 February 2013 0028 // 0029 // Description: Templated base class for exponential parametrization of 0030 // two-body angular distributions in Bertini-style cascade 0031 // 0032 // 20130227 Renamed from "ParamExp" to "ParamExp" to fix misattribution. 0033 0034 #ifndef G4ParamExpTwoBodyAngDst_h 0035 #define G4ParamExpTwoBodyAngDst_h 1 0036 0037 #include "G4VTwoBodyAngDst.hh" 0038 #include "G4CascadeInterpolator.hh" 0039 0040 0041 template <G4int NKEBINS> 0042 class G4ParamExpTwoBodyAngDst : public G4VTwoBodyAngDst { 0043 public: 0044 enum { nKEbins=NKEBINS }; // For use in function arguments 0045 0046 G4ParamExpTwoBodyAngDst(const G4String& name, 0047 const G4double (&kebins)[nKEbins], 0048 const G4double (&pFrac)[nKEbins], 0049 const G4double (&pA)[nKEbins], 0050 const G4double (&pC)[nKEbins], 0051 const G4double (&pCos)[nKEbins], 0052 G4int verbose = 0) 0053 : G4VTwoBodyAngDst(name, verbose), labKE(kebins), angleCut(pFrac), 0054 smallScale(pA), largeScale(pC), cosScale(pCos), interpolator(kebins) {;} 0055 0056 virtual ~G4ParamExpTwoBodyAngDst() {;} 0057 0058 virtual G4double GetCosTheta(const G4double& ekin, const G4double& pcm) const; 0059 0060 protected: 0061 // Kinetic energies at which angular distributions are taken 0062 const G4double (&labKE)[nKEbins]; 0063 0064 // Cut-off for small- vs. large-angle functions at each energy 0065 const G4double (&angleCut)[nKEbins]; 0066 0067 // Scale factors for small- vs. large-angle regions 0068 const G4double (&smallScale)[nKEbins]; 0069 const G4double (&largeScale)[nKEbins]; 0070 0071 // Scale factor for exponential mapping to cos(theta) 0072 const G4double (&cosScale)[nKEbins]; 0073 0074 G4CascadeInterpolator<NKEBINS> interpolator; 0075 }; 0076 0077 // Implementations must be included for templated classes 0078 #include "G4ParamExpTwoBodyAngDst.icc" 0079 0080 #endif /* G4ParamExpTwoBodyAngDst_h */
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |