Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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// G4DataInterpolation
27//
28// Class description:
29//
30// The class consists of some methods for data interpolations and
31// extrapolations. The methods based mainly on recommendations given in the
32// book: An introduction to NUMERICAL METHODS IN C++, B.H. Flowers,
33// Claredon Press, Oxford, 1995.
34
35// Author: V.Grichine, 03.04.1997
36// --------------------------------------------------------------------
37#ifndef G4DATAINTERPOLATION_HH
38#define G4DATAINTERPOLATION_HH 1
39
40#include "globals.hh"
41
43{
44 public:
45 G4DataInterpolation(G4double pX[], G4double pY[], G4int number);
46 // Constructor for initializing data members.
47
48 G4DataInterpolation(G4double pX[], G4double pY[], G4int number,
49 G4double pFirstDerStart, G4double pFirstDerFinish);
50 // Constructor for cubic spline interpolation. It creates fSecond Deivative
51 // array as well as fArgument and fFunction.
52
54 // Destructor deletes dynamically created arrays for data members: fArgument,
55 // fFunction and fSecondDerivative, all have dimension of fNumber.
56
59 // Copy constructor and assignement operator not allowed.
60
62 // This function returns the value P(pX), where P(x) is polynom of fNumber-1
63 // degree such that P(fArgument[i]) = fFunction[i], for i = 0, ..., fNumber-1.
64
65 void PolIntCoefficient(G4double cof[]) const;
66 // Given arrays fArgument[0,..,fNumber-1] and fFunction[0,..,fNumber-1], this
67 // function calculates an array of coefficients.
68 // The coefficients don't provide usually (fNumber>10) better accuracy for
69 // polynom interpolation, as compared with PolynomInterpolation() function.
70 // They could be used instead for derivate calculations and some other
71 // applications.
72
74 // The function returns diagonal rational function (Bulirsch and Stoer
75 // algorithm 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
77 // with polynomial interpolation.
78
80 // Cubic spline interpolation in point pX for function given by the table:
81 // fArgument, fFunction. The constructor, which creates fSecondDerivative,
82 // must be called before. The function works optimal, if sequential calls
83 // are in random values of pX.
84
85 G4double FastCubicSpline(G4double pX, G4int index) const;
86 // Return cubic spline interpolation in the point pX which is located between
87 // fArgument[index] and fArgument[index+1]. It is usually called in sequence
88 // of known from external analysis values of index.
89
91 // Given argument pX, returns index k, so that pX bracketed by fArgument[k]
92 // and fArgument[k+1].
93
94 void CorrelatedSearch(G4double pX, G4int& index) const;
95 // Given a value pX, returns a value 'index' such that pX is between
96 // fArgument[index] and fArgument[index+1]. fArgument MUST BE MONOTONIC,
97 // either increasing or decreasing. If index = -1 or fNumber, this indicates
98 // that pX is out of range. The value index on input is taken as the initial
99 // approximation for index on output.
100
101 private:
102 // pointers to data table to be interpolated for y[i] and x[i] respectively
103 G4double* fArgument = nullptr;
104 G4double* fFunction = nullptr;
105
106 G4double* fSecondDerivative = nullptr;
107
108 G4int fNumber = 0; // the corresponding table size
109};
110
111#endif
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4double RationalPolInterpolation(G4double pX, G4double &deltaY) const
void CorrelatedSearch(G4double pX, G4int &index) const
G4double FastCubicSpline(G4double pX, G4int index) const
G4DataInterpolation(const G4DataInterpolation &)=delete
G4double PolynomInterpolation(G4double pX, G4double &deltaY) const
G4int LocateArgument(G4double pX) const
G4DataInterpolation & operator=(const G4DataInterpolation &)=delete
G4double CubicSplineInterpolation(G4double pX) const
void PolIntCoefficient(G4double cof[]) const