Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ChannelingFastSimInterpolation.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26
27#ifndef G4ChannelingFastSimInterpolation_h
28#define G4ChannelingFastSimInterpolation_h
29
30#include "G4PhysicsVector.hh"
31#include "G4Physics2DVector.hh"
32#include "G4ThreeVector.hh"
34
35/** \file G4ChannelingFastSimInterpolation.hh
36* \brief Definition of the G4ChannelingFastSimInterpolation class
37* The class includes spline interpolation coefficients for the important functions
38* needed in Channeling FastSimulation model, i.e. electric fields, nuclear and
39* electron densities and minimum energy of ionization. All the functions have
40* transverse coordinates as arguments (x for crystal planes, x and y for crystal axes).
41* All the functions are calculated in the co-rotating reference system (along crystal
42* planes/axes).
43*/
44
46
47public:
49 G4double dy0,
50 G4int nPointsx0,
51 G4int nPointsy0,
52 G4int iModel0);
54
55 ///Get Spline Function
57
58 ///Set spline coefficients
60 G4double BI0,
61 G4double CI0,
62 G4double DI0,
63 G4int i);
64 void SetCoefficients2D(G4double AI3D0,
65 G4double BI3D0,
66 G4double CI3D0,
67 G4int i,
68 G4int j,
69 G4int k);
70
71private:
72 G4double Spline1D(G4double xx);
73 G4double Spline2D(G4double xx, G4double yy);// cubic spline of 2-variable function
74
75 G4double fDx=0., fDy=0.; //channel width and height
76 G4double fStepi=0., fStepj=0.; //interpolation steps in x and y, respectively
77 G4double fStepi2=0.; //=fStepi*fStepi
78 G4int nPointsx=0, nPointsy=0; //number of interpolation nodes in x and y, respectively
79
80 std::vector <G4double> fAI;
81 std::vector <G4double> fBI;
82 std::vector <G4double> fCI;
83 std::vector <G4double> fDI;
84
85 std::vector<std::vector<G4double>> fAI3D;
86 std::vector<std::vector<G4double>> fBI3D;
87 std::vector<std::vector<G4double>> fCI3D;
88 std::vector<std::vector<G4double>> fAI3D3;
89 std::vector<std::vector<G4double>> fBI3D3;
90 std::vector<std::vector<G4double>> fCI3D3;
91
92 G4int iModel=1;
93
94};
95
96#endif
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
void SetCoefficients1D(G4double AI0, G4double BI0, G4double CI0, G4double DI0, G4int i)
Set spline coefficients.
void SetCoefficients2D(G4double AI3D0, G4double BI3D0, G4double CI3D0, G4int i, G4int j, G4int k)
G4ChannelingFastSimInterpolation(G4double dx0, G4double dy0, G4int nPointsx0, G4int nPointsy0, G4int iModel0)
G4double GetIF(G4double xx, G4double yy)
Get Spline Function.