|
||||
File indexing completed on 2025-01-18 09:58:01
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 #ifndef G4ChannelingFastSimInterpolation_h 0028 #define G4ChannelingFastSimInterpolation_h 0029 0030 #include "G4PhysicsVector.hh" 0031 #include "G4Physics2DVector.hh" 0032 #include "G4ThreeVector.hh" 0033 #include "G4PhysicsLinearVector.hh" 0034 0035 /** \file G4ChannelingFastSimInterpolation.hh 0036 * \brief Definition of the G4ChannelingFastSimInterpolation class 0037 * The class includes spline interpolation coefficients for the important functions 0038 * needed in Channeling FastSimulation model, i.e. electric fields, nuclear and 0039 * electron densities and minimum energy of ionization. All the functions have 0040 * transverse coordinates as arguments (x for crystal planes, x and y for crystal axes). 0041 * All the functions are calculated in the co-rotating reference system (along crystal 0042 * planes/axes). 0043 */ 0044 0045 class G4ChannelingFastSimInterpolation { 0046 0047 public: 0048 G4ChannelingFastSimInterpolation(G4double dx0, 0049 G4double dy0, 0050 G4int nPointsx0, 0051 G4int nPointsy0, 0052 G4int iModel0); 0053 ~G4ChannelingFastSimInterpolation() = default; 0054 0055 ///Get Spline Function 0056 G4double GetIF(G4double xx, G4double yy); 0057 0058 ///Set spline coefficients 0059 void SetCoefficients1D(G4double AI0, 0060 G4double BI0, 0061 G4double CI0, 0062 G4double DI0, 0063 G4int i); 0064 void SetCoefficients2D(G4double AI3D0, 0065 G4double BI3D0, 0066 G4double CI3D0, 0067 G4int i, 0068 G4int j, 0069 G4int k); 0070 0071 private: 0072 G4double Spline1D(G4double xx); 0073 G4double Spline2D(G4double xx, G4double yy);// cubic spline of 2-variable function 0074 0075 G4double fDx=0., fDy=0.; //channel width and height 0076 G4double fStepi=0., fStepj=0.; //interpolation steps in x and y, respectively 0077 G4double fStepi2=0.; //=fStepi*fStepi 0078 G4int nPointsx=0, nPointsy=0; //number of interpolation nodes in x and y, respectively 0079 0080 std::vector <G4double> fAI; 0081 std::vector <G4double> fBI; 0082 std::vector <G4double> fCI; 0083 std::vector <G4double> fDI; 0084 0085 std::vector<std::vector<G4double>> fAI3D; 0086 std::vector<std::vector<G4double>> fBI3D; 0087 std::vector<std::vector<G4double>> fCI3D; 0088 std::vector<std::vector<G4double>> fAI3D3; 0089 std::vector<std::vector<G4double>> fBI3D3; 0090 std::vector<std::vector<G4double>> fCI3D3; 0091 0092 G4int iModel=1; 0093 0094 }; 0095 0096 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |