Geant4  10.00.p01
G4DataInterpolation.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 // $Id: G4DataInterpolation.hh 67970 2013-03-13 10:10:06Z gcosmo $
28 //
29 // Class description:
30 //
31 // The class consists of some methods for data interpolations and extrapolations.
32 // The methods based mainly on recommendations given in the book : An introduction to
33 // NUMERICAL METHODS IN C++, B.H. Flowers, Claredon Press, Oxford, 1995
34 //
35 // ------------------------------ Data members: ---------------------------------
36 //
37 // fArgument and fFunction - pointers to data table to be interpolated
38 // for y[i] and x[i] respectively
39 // fNumber - the corresponding table size
40 // ......
41 // G4DataInterpolation( G4double pX[], G4double pY[], G4int number )
42 //
43 // Constructor for initializing of fArgument, fFunction and fNumber data members:
44 // ......
45 // G4DataInterpolation( G4double pX[], G4double pY[], G4int number,
46 // G4double pFirstDerStart, G4double pFirstDerFinish )
47 //
48 // Constructor for cubic spline interpolation. It creates the array
49 // fSecondDerivative[0,...fNumber-1] which is used in this interpolation by
50 // the function:
51 // ....
52 // ~G4DataInterpolation()
53 //
54 // Destructor deletes dynamically created arrays for data members: fArgument,
55 // fFunction and fSecondDerivative, all have dimension of fNumber
56 //
57 // ------------------------------ Methods: ----------------------------------------
58 //
59 // G4double PolynomInterpolation(G4double pX, G4double& deltaY ) const
60 //
61 // This function returns the value P(pX), where P(x) is polynom of fNumber-1 degree
62 // such that P(fArgument[i]) = fFunction[i], for i = 0, ..., fNumber-1 .
63 // ........
64 // void PolIntCoefficient( G4double cof[]) const
65 //
66 // Given arrays fArgument[0,..,fNumber-1] and fFunction[0,..,fNumber-1] , this
67 // function calculates an array of coefficients. The coefficients don't provide
68 // usually (fNumber>10) better accuracy for polynom interpolation, as compared with
69 // PolynomInterpolation function. They could be used instead for derivate
70 // calculations and some other applications.
71 // .........
72 // G4double RationalPolInterpolation(G4double pX, G4double& deltaY ) const
73 //
74 // The function returns diagonal rational function (Bulirsch and Stoer algorithm
75 // of Neville type) Pn(x)/Qm(x) where P and Q are polynoms.
76 // Tests showed the method is not stable and hasn't advantage if compared with
77 // polynomial interpolation
78 // ................
79 // G4double CubicSplineInterpolation(G4double pX) const
80 //
81 // Cubic spline interpolation in point pX for function given by the table:
82 // fArgument, fFunction. The constructor, which creates fSecondDerivative, must be
83 // called before. The function works optimal, if sequential calls are in random
84 // values of pX.
85 // ..................
86 // G4double FastCubicSpline(G4double pX, G4int index) const
87 //
88 // Return cubic spline interpolation in the point pX which is located between
89 // fArgument[index] and fArgument[index+1]. It is usually called in sequence of
90 // known from external analysis values of index.
91 // .........
92 // G4int LocateArgument(G4double pX) const
93 //
94 // Given argument pX, returns index k, so that pX bracketed by fArgument[k] and
95 // fArgument[k+1]
96 // ......................
97 // void CorrelatedSearch( G4double pX, G4int& index ) const
98 //
99 // Given a value pX, returns a value 'index' such that pX is between fArgument[index]
100 // and fArgument[index+1]. fArgument MUST BE MONOTONIC, either increasing or
101 // decreasing. If index = -1 or fNumber, this indicates that pX is out of range.
102 // The value index on input is taken as the initial approximation for index on
103 // output.
104 
105 // --------------------------------- History: --------------------------------------
106 //
107 // 3.4.97 V.Grichine (Vladimir.Grichine@cern.ch)
108 //
109 
110 
111 #ifndef G4DATAINTERPOLATION_HH
112 #define G4DATAINTERPOLATION_HH
113 
114 #include "globals.hh"
115 
117 {
118 public:
120  G4double pY[],
121  G4int number );
122 
123 // Constructor for cubic spline interpolation. It creates fSecond Deivative array
124 // as well as fArgument and fFunction
125 
127  G4double pY[],
128  G4int number,
129  G4double pFirstDerStart,
130  G4double pFirstDerFinish ) ;
131 
133 
135  G4double& deltaY ) const ;
136 
137  void PolIntCoefficient( G4double cof[]) const ;
138 
140  G4double& deltaY ) const ;
141 
143 
145  G4int index ) const ;
146 
147  G4int LocateArgument( G4double pX ) const ;
148 
149  void CorrelatedSearch( G4double pX,
150  G4int& index ) const ;
151 
152 private:
153 
156 
157 private:
162 } ;
163 
164 #endif
G4double RationalPolInterpolation(G4double pX, G4double &deltaY) const
G4double FastCubicSpline(G4double pX, G4int index) const
void CorrelatedSearch(G4double pX, G4int &index) const
G4DataInterpolation(G4double pX[], G4double pY[], G4int number)
G4int LocateArgument(G4double pX) const
int G4int
Definition: G4Types.hh:78
G4double CubicSplineInterpolation(G4double pX) const
void PolIntCoefficient(G4double cof[]) const
G4DataInterpolation & operator=(const G4DataInterpolation &)
double G4double
Definition: G4Types.hh:76
G4double PolynomInterpolation(G4double pX, G4double &deltaY) const