Isis 3 Developer Reference
|
NumericalApproximation provides various numerical analysis methods of interpolation, extrapolation and approximation of a tabulated set of x, y data. More...
#include <NumericalApproximation.h>
Public Types | |
enum | InterpType { Linear, Polynomial, PolynomialNeville, CubicNatural, CubicClamped, CubicNatPeriodic, CubicNeighborhood, CubicHermite, Akima, AkimaPeriodic } |
This enum defines the types of interpolation supported in this class. More... | |
enum | ExtrapType { ThrowError, Extrapolate, NearestEndpoint } |
This enum defines the manner in which a value outside of the domain should be handled if passed to the Evaluate() method. More... | |
Public Member Functions | |
NumericalApproximation (const NumericalApproximation::InterpType &itype=CubicNatural) | |
Default constructor creates NumericalApproximation object. More... | |
NumericalApproximation (unsigned int n, double *x, double *y, const NumericalApproximation::InterpType &itype=CubicNatural) | |
Constructs NumericalApproximation object, sets interpolation type, populates the data set with pointer values. More... | |
NumericalApproximation (const vector< double > &x, const vector< double > &y, const NumericalApproximation::InterpType &itype=CubicNatural) | |
Constructs NumericalApproximation object, sets interpolation type, populates the data set with vector values. More... | |
NumericalApproximation (const NumericalApproximation &dint) | |
Copy constructor creates NumericalApproximation object and sets the copies the class variable values and state of input NumericalApproximation object to the new object. More... | |
NumericalApproximation & | operator= (const NumericalApproximation &numApMeth) |
NumericalApproximation assigmment operator sets this object "equal to" another. More... | |
virtual | ~NumericalApproximation () |
Destructor deallocates memory being used. More... | |
string | Name () const |
Get name of interpolating function assigned to object. More... | |
InterpType | InterpolationType () |
Returns the enumerated type of interpolation chosen. More... | |
int | MinPoints () |
Minimum number of points required by interpolating function. More... | |
int | MinPoints (NumericalApproximation::InterpType itype) |
Minimum number of points required by interpolating function. More... | |
void | AddData (const double x, const double y) |
Add a datapoint to the set. More... | |
void | AddData (unsigned int n, double *x, double *y) |
Add set of points to the data set using pointer arrays. More... | |
void | AddData (const vector< double > &x, const vector< double > &y) |
Add set of points to the data set using vectors. More... | |
void | SetCubicClampedEndptDeriv (const double yp1, const double ypn) |
Sets the values for the first derivatives of the endpoints of the data set. More... | |
void | AddCubicHermiteDeriv (unsigned int n, double *fprimeOfx) |
Adds values for the first derivatives of the data points. More... | |
void | AddCubicHermiteDeriv (const vector< double > &fprimeOfx) |
Adds values for the first derivatives of the data points. More... | |
void | AddCubicHermiteDeriv (const double fprimeOfx) |
Adds value of a first derivative of a data point. More... | |
double | DomainMinimum () |
Input data domain minimum value. More... | |
double | DomainMaximum () |
Input data domain maximum value. More... | |
bool | Contains (double x) |
Returns whether the passed value is an element of the set of x-values in the data set. More... | |
unsigned int | Size () |
Returns the number of the coordinates added to the data set. More... | |
double | Evaluate (const double a, const ExtrapType &etype=ThrowError) |
Calculates interpolated or extrapolated value of tabulated data set for given domain value. More... | |
vector< double > | Evaluate (const vector< double > &a, const ExtrapType &etype=ThrowError) |
Calculates interpolated value for given set of domain values. More... | |
vector< double > | PolynomialNevilleErrorEstimate () |
Retrieves the error estimate for the Neville's polynomial interpolation type. More... | |
vector< double > | CubicClampedSecondDerivatives () |
Retrieves the second derivatives of the data set. More... | |
double | EvaluateCubicHermiteFirstDeriv (const double a) |
Approximates the first derivative of the data set function evaluated at the given domain value for CubicHermite interpolation type. More... | |
double | EvaluateCubicHermiteSecDeriv (const double a) |
Approximates the second derivative of the data set function evaluated at the given domain value for CubicHermite interpolation type. More... | |
double | GslFirstDerivative (const double a) |
Approximates the first derivative of the data set function evaluated at the given domain value for GSL supported interpolation types. More... | |
double | BackwardFirstDifference (const double a, const unsigned int n=3, const double h=0.1) |
Uses an n point backward first difference formula to approximate the first derivative evaluated at a given domain value. More... | |
double | ForwardFirstDifference (const double a, const unsigned int n=3, const double h=0.1) |
Uses an n point forward first difference formula to approximate the first derivative evaluated at a given domain value. More... | |
double | CenterFirstDifference (const double a, const unsigned int n=5, const double h=0.1) |
Uses an n point center first difference formula to approximate the first derivative evaluated at a given domain value. More... | |
double | GslSecondDerivative (const double a) |
Approximates the second derivative of the interpolated data set function evaluated at the given domain value for GSL supported interpolation types. More... | |
double | BackwardSecondDifference (const double a, const unsigned int n=3, const double h=0.1) |
Uses an n point backward second difference formula to approximate the second derivative evaluated at a given domain value. More... | |
double | ForwardSecondDifference (const double a, const unsigned int n=3, const double h=0.1) |
Uses an n point forward second difference formula to approximate the second derivative evaluated at a given domain value. More... | |
double | CenterSecondDifference (const double a, const unsigned int n=5, const double h=0.1) |
Uses an n point center second difference formula to approximate the second derivative evaluated at a given domain value. More... | |
double | GslIntegral (const double a, const double b) |
Approximates the integral of the data set function evaluated on the given interval for GSL supported interpolation types. More... | |
double | TrapezoidalRule (const double a, const double b) |
Uses the trapezoidal rule to approximate the integral of the interpolated data set function on the interval (a, b). More... | |
double | Simpsons3PointRule (const double a, const double b) |
Uses Simpson's 3-point rule to approximate the integral of the interpolated data set function on the interval (a,b). More... | |
double | Simpsons4PointRule (const double a, const double b) |
Uses Simpson's 4-point rule to approximate the integral of the interpolated data set function on the interval (a,b). More... | |
double | BoolesRule (const double a, const double b) |
Uses Boole's Rule to approximate the integral of the interpolated data set function on the interval (a,b). More... | |
double | RefineExtendedTrap (double a, double b, double s, unsigned int n) |
Calculates refinements extended trapezoidal rule to approximate the integral of the interpolated data set function on the interval (a,b). More... | |
double | RombergsMethod (double a, double b) |
Uses Romberg's method to approximate the integral of the interpolated data set function on the interval (a,b). More... | |
void | Reset () |
Resets the state of the object. More... | |
void | Reset (NumericalApproximation::InterpType itype) |
Resets the state of the object and resets interpolation type. More... | |
void | SetInterpType (NumericalApproximation::InterpType itype) |
Sets interpolation type. More... | |
Protected Types | |
typedef const gsl_interp_type * | InterpFunctor |
GSL Interpolation specs. More... | |
typedef map< InterpType, InterpFunctor > | FunctorList |
Set up a std::map of GSL interpolator functors. List of function types. More... | |
typedef FunctorList::const_iterator | FunctorConstIter |
GSL Iterator. More... | |
Protected Member Functions | |
void | Init (NumericalApproximation::InterpType itype) |
Initializes the object upon instantiation. More... | |
bool | GslInterpType (NumericalApproximation::InterpType itype) const |
Returns whether an interpolation type is adapted from the GSL library. More... | |
void | GslAllocation (unsigned int npoints) |
Allocates GSL interpolation functions. More... | |
void | GslDeallocation () |
Deallocate GSL interpolator resources, if used. More... | |
InterpFunctor | GslFunctor (NumericalApproximation::InterpType itype) const |
Search for a GSL interpolation function. More... | |
void | GslIntegrityCheck (int gsl_status, const char *src, int line) |
Checks the status of the GSL interpolation operations. More... | |
void | ValidateDataSet () |
Validates the data set before computing interpolation. More... | |
bool | InsideDomain (const double a) |
Returns whether the passed value is greater than or equal to the domain minimum and less than or equal to the domain maximum. More... | |
bool | GslComputed () const |
Returns whether a GSL interpolation computation of the data set has been performed. More... | |
void | ComputeGsl () |
Computes the GSL interpolation for a set of (x,y) data points. More... | |
void | ComputeCubicClamped () |
Computes the cubic clamped interpolation for a set of (x,y) data points, given the first derivatives of the endpoints of the data set. More... | |
double | ValueToExtrapolate (const double a, const ExtrapType &etype) |
Returns the domain value at which to evaluate. More... | |
double | EvaluateCubicNeighborhood (const double a) |
Performs cubic spline interpolation for a neighborhood about a. More... | |
vector< double > | EvaluateCubicNeighborhood (const vector< double > &a, const ExtrapType &etype) |
Performs cubic spline interpolations for neighborhoods about each value of a. More... | |
double | EvaluateCubicClamped (const double a) |
Performs cubic spline interpolation with clamped boundary conditions, if possible. More... | |
double | EvaluateCubicHermite (const double a) |
Performs interpolation using the Hermite cubic polynomial. More... | |
double | EvaluatePolynomialNeville (const double a) |
Performs polynomial interpolation using Neville's algorithm. More... | |
vector< double > | EvaluateForIntegration (const double a, const double b, const unsigned int n) |
Evaluates data set in order to have enough data points to approximate the function to be integrated. More... | |
int | FindIntervalLowerIndex (const double a) |
Find the index of the x-value in the data set that is just below the input value, a. More... | |
void | ReportException (IException::ErrorType type, const string &method, const string &message, const char *filesrc, int lineno) const |
Generalized error report generator. More... | |
Protected Attributes | |
InterpType | p_itype |
Interpolation type. More... | |
vector< double > | p_x |
List of X values. More... | |
vector< double > | p_y |
List of Y values. More... | |
bool | p_dataValidated |
Flag variable to determine whether ValidateDataSet() has been called. More... | |
gsl_interp_accel * | p_acc |
Lookup accelorator. More... | |
gsl_spline * | p_interp |
Currently active interpolator. More... | |
bool | p_clampedEndptsSet |
Flag variable to determine whether SetCubicClampedEndptDeriv() has been called after all data was added for CubicClamped interpolation. More... | |
bool | p_clampedComputed |
Flag variable to determine whether ComputeCubicClamped() has been called. More... | |
double | p_clampedDerivFirstPt |
First derivative of first x-value, p_x[0]. This is only used for the CubicClamped interpolation type. More... | |
double | p_clampedDerivLastPt |
First derivative of last x-value, p_x[n-1]. This is only used for the CubicClamped interpolation type. More... | |
vector< double > | p_clampedSecondDerivs |
List of second derivatives evaluated at p_x values. This is only used for the CubicClamped interpolation type. More... | |
vector< double > | p_polyNevError |
Estimate of error for interpolation evaluated at x. This is only used for the PolynomialNeville interpolation type. 91 taken from AtmosModel. More... | |
vector< double > | p_fprimeOfx |
List of first derivatives corresponding to each x value in the data set (i.e. each value in p_x) More... | |
Static Protected Attributes | |
static FunctorList | p_interpFunctors |
Maintains list of interpolator options. More... | |
NumericalApproximation provides various numerical analysis methods of interpolation, extrapolation and approximation of a tabulated set of x, y data.
This class contains a merged version of the Isis3 classes DataInterp and NumericalMethods. In addition, some methods from AtmosModel were moved to this class, the CubicNeighborhood interpolation type was adapted from the IDL routine called interpol (using the "/spline" keyword), and several differentiation and integration approximation methods were created.
NumericalApproximation contains 9 types of data interpolation routines. These routines act on x, y pairs of data points. The following forms of data interpolation are supported:
Name | Type (enum) | MinPoints | Description |
---|---|---|---|
Linear | Linear | 2 | Linear interpolation approximates a curve by concatinating line segments between known data points. This results in a continuous curve with a discontinuities of the derivative at the known data points. This interpolation type uses the GSL library routines. |
Polynomial | Polynomial | 3 | Computes a polynomial interpolation. This method should only be used for interpolating small numbers of points because polynomial interpolation introduces large oscillations, even for well-behaved datasets. The number of terms in the interpolating polynomial is equal to the number of points. This interpolation type uses the GSL library routines. |
Neville's Algorithm for Polynomial Interpolation | PolynomialNeville | 3 | Polynomial interpolation using Neville's algorithm. This first fits a polynomial of degree 0 through each point. On the second iteration, the polynomial of adjoining indices are combined to fit through pairs of points. This process repeats until a pyramid of approximations is reached. This is based on the Newton form of the interpolating polynomial and the recursion relation of the divided differences. |
Natural Cubic Spline | CubicNatural | 3 | Cubic spline with natural boundary conditions. The resulting curve is piecewise cubic on each interval with matching first and second derivatives at the supplied data points. The second derivative is chosen to be zero at the endpoints, i.e. the first and last points of the data set. This is also refered to as free boundary condtions. This interpolation type uses the GSL library routines. |
Clamped Cubic Spline | CubicClamped | 3 | Cubic Spline interpolation with clamped boundary conditions. The resulting curve is piecewise cubic on each interval. For this type of boundary condition to hold, it is necessary to have either the values of the derivative at the endpoints or an accurate approximation to those values. In general, clamped boundary conditions lead to more accurate approximations since they include more information about the function. |
Periodic Natural Cubic Spline | CubicNatPeriodic | 2 | Cubic spline with periodic boundary conditions. The resulting curve is piecewise cubic on each interval with matching first and second derivatives at the supplied data points. The derivatives at the first and last points are also matched. Note that the last point in the data must have the same y-value as the first point, otherwise the resulting periodic interpolation will have a discontinuity at the boundary. This interpolation type uses the GSL library routines. |
Cubic Spline about 4-point Neighborhood | CubicNeighborhood | 4 | Cubic Spline interpolation using 4-pt Neighborhoods with natural boundary conditions. This interpolation type is adapted from the IDL interpol.pro application using the "/spline" keyword on an irregular grid. This type of cubic spline fits a natural cubic spline to the 4-point neighborhood of known data points surrounding the x value at which we wish to evaluate. For example, suppose {x0, x1, ..., xn} is the array of known domain values in the data set and \( x_i \leq u < x_{i+1} \) for some i such that \( 0 \leq i \leq n \), then to evaluate the y value at u, this method of interpolation evaluates the ordinary cubic spline consisting of the data set {xi-1, xi, xi+1, xi+2} at u. |
Cubic Spline using Hermite cubic polynomial | CubicHermite | 2 | Cubic Spline interpolation using the Hermite cubic polynomial (also called the unique Hermite interpolating fundamental polynomial of degree 3). This method requires the input of the slopes of the tangent lines of the curve (i.e. velocities, or first derivatives) at each of the data points. This ensures that the approximation has the same "shape" as the curve at the known data points. This interpolation type uses piecewise Hermite cubic polynomials for each interval, (x0, x1), using the known data points, (x0,f(x0)) and (x1,f(x1)), and their derivatives, f '(x0) and f '(x1). The Hermite cubic polynomial is defined \[ H_3(x) = f(x_0)h_0(x)+f(x_1)h_1(x) +f\prime(x_0)\hat{h}_0(x)+f\prime(x_1)\hat{h}_1(x) \] where \( h_k(x) = [1-2(x-x_k)L\prime_k(x_k)](L_k(x))^2 \) and \( \hat{h}_k(x) = (x-x_k)(L_k(x))^2 \) for the kth-Lagrange coefficient polynomials of degree n = 1, \( L_0(x) = \frac{x- x_1}{x_0-x_1}\) and \( L_1(x) = \frac{x- x_0}{x_1-x_0}\). |
Akima | Akima | 5 | Non-rounded Akima spline with natural boundary conditions. This method uses non-rounded corner algorithm of Wodicka. This interpolation type uses the GSL library routines. |
Periodic Akima | AkimaPeriodic | 5 | Non-rounded Akima spline with periodic boundary conditions. This method uses non-rounded corner algorithm of Wodicka. This interpolation type uses the GSL library routines. |
Numerical analysis approximation methods for differentiating and integrating unknown functions represented by a data set are implemented by using the interpolations on the data set created with a NumericalApproximation object. The Lagrange polynomial,
\[ L(x) = \sum_{k=0}^{n} \left(f(x_k) \prod_{i=0,i \neq k}^{n} \frac{x- x_i}{x_k-x_i}\right) \]
is used to determine formulas for numerical differentiation and integration.
Numerical Differentiation
The differentiation methods use difference formulas to approximate first and second derivatives at a domain value for a given a set of known (x,y) data points. Once all of the data points are added, a Isis::NumericalApproximation interpolation is used on the dataset to estimate the function, \(f:x \to y\), mapping x values to corresponding y values. Then, the derivative of f evaluated at x0 is approximated. To do so, a uniform step size of h is chosen to determine the distance between xi and xi+1.
Numerical Differentiation Type | Difference Formulas | Methods Available | Description |
---|---|---|---|
Backward Difference |
| Backward difference formulas use a uniform step-size moving in the negative direction from the given x0. These formulas are derived by differentiating the Lagrange polynomials for xi = x0 - ih and evaluating them atx0. | |
Forward Difference |
| Forward difference formulas use a uniform step-size moving in the positive direction from x0. These formulas are derived by differentiating the Lagrange polynomials for xi = x0 + ih and evaluating them atx0. | |
Center Difference |
| Center difference formulas use a uniform step-size moving in both directions from x0 so this point stays centered. These formulas are derived by differentiating the Lagrange polynomials for xi = x0
| |
GSL Differentiation | Unknown | No documentation was found about the algorithms used for these methods. They may only be called when interpolation type is one of the following: Linear, Polynomial, CubicNatural, CubicNatPeriodic, Akima, or AkimaPeriodic. |
Numerical Integration
The integration methods were derived using Newton-Cotes, or quadrature, formulas for approximating integrals given a set of known (x,y) data points. The x values may be irregularly spaced, but must be unique and added to the dataset in ascending order. Once all of the data points are added, a Isis::NumericalApproximation interpolation is used on the dataset to estimate the function, \(f:x \to y\), mapping x values to corresponding y values. Then, the integral of f on the interval from a to b is approximated. To do so, an algorithm for creating a composite formula by applying the original formula to segments that share a boundary point is used. The NumericalApproximation::InterpType chosen for interpolation computation seems to have little affect on the error between the actual integral and the return values of the integration methods. The BoolesRule() method varies the most in error. Although errors are not high for any of the interpolation types, CubicNatural spline interpolations seem to have the smallest error most often with BoolesRule(). For any other numerical integration method, the errors appear to be identical for any NumericalApproximation::InterpType except Polynomial, which usually has a slightly larger error than the other interpolation types. Note: A portion of this algorithm is derived from the IDL function int_tabulated for Boole's Method.
Numerical Integration Type | Integration Formulas | Methods Available | Description |
---|---|---|---|
Trapezoidal Rule |
| The 2-point closed rule, known as the trapezoidal rule, uses straight line segments between known data points to estimate the area under the curve (integral). This Newton-Cotes formula uses a uniform step-size of h = b - a and is derived by integrating the Lagrange polynomials over the closed interval [a, a+h] for xi = x0 + ih. | |
Simpson's Rule, or Simpson's 3-Point Rule |
| The 3-point closed rule, known as Simpson's Rule or Simpson's 3-Point Rule, uses parabolic arcs between known data points to estimate the area under the curve (integral). This Newton-Cotes formula uses a uniform step-size of h = (b - a)/2 and is derived by integrating the Lagrange polynomials over the closed interval [a, a+2h] for xi = x0 + ih. | |
Simpson's 3/8 Rule, or Simpson's 4-Point Rule |
| The 4-point closed rule, known as Simpson's 3/8 Rule or Simpson's 4-Point Rule, uses cubic curves between known data points to estimate the area under the curve (integral). This Newton-Cotes formula uses a uniform step-size of h = (b - a)/3 and is derived by integrating the Lagrange polynomials over the closed interval [a, a+3h] for xi = x0 + ih. | |
Boole's Rule |
| The 5-point closed rule, known as Boole's Rule, uses quartic curves between known data points to estimate the area under the curve (integral). This Newton-Cotes formula uses a uniform step-size of h = (b - a)/4 and is derived by integrating the Lagrange polynomials over the closed interval [a, a+4h] for xi = x0 + ih. | |
Refinements of Extended Trapezoidal Rule |
| The extended (or composite) trapezoidal rule can be used to with a series of refinements to approximate the integral. The first stage of refinement returns the ordinary trapezoidal estimate of the integral. Subsequent stages will improve the accuracy, where, for the nth stage 2n-2 interior points are added. | |
Romberg's Method |
| Romberg's Method is a potent numerical integration tool. It uses the extended trapezoidal rule and Neville's algorithm for polynomial interpolation. | |
GSL Integration | Unknown | No documentation was found about the algorithm used for this methods. They may only be called when interpolation type is one of the following: Linear, Polynomial, CubicNatural, CubicNatPeriodic, Akima, or AkimaPeriodic. |
Below is an example demonstating the use of this class:
To compute the same data set using the Akima spline, just use the following:
When using this class, there are some important details that require consideration. Several interpolation types of this class use the GSL library. In addition, the GSL library default error handling scheme has implications when used in C++ objects. The default behavior is to terminate the application when an error occurs. Options for developers are to turn off error trapping within the GSL, which means users of the library must do their own error checking, or a specific handler can be specified to replace the default mode. Neither option is well suited in object-oriented implementations because they apply globally, to all users of the library. The approach used here is to turn off error handling and require users to handle their own error checking. This is not entirely safe as any user that does not follow this strategy could change this behavior at runtime. Other options should be explored. An additional option is that the default behavior can be altered prior to building GSL but this is also problematic since this behavior cannot be guaranteed universally.
For these types, interpolation is strictly adhered to and extrapolation is not handled well. If the input to be evaluated is outside of the minimum and maximum x values of the original data points, then the y value corresponding to the x value that is nearest to that input is returned. However, for sufficiently close values, the clamped cubic and polynomial Neville's types are fairly accurate in their extrapolation techniques. All differentiation and integration methods throw an error if the value passed is outside of the domain.
|
protected |
GSL Iterator.
|
protected |
Set up a std::map of GSL interpolator functors. List of function types.
|
protected |
GSL Interpolation specs.
This enum defines the manner in which a value outside of the domain should be handled if passed to the Evaluate() method.
Enumerator | |
---|---|
ThrowError | Evaluate() throws an error if a is outside of the domain. |
Extrapolate | Evaluate() attempts to extrapolate if a is outside of the domain. This is only valid for NumericalApproximation::InterpType CubicClamped or PolynomialNeville and the result will be accurate only if sufficiently close to the domain boundary. |
NearestEndpoint | Evaluate() returns the y-value of the nearest endpoint if a is outside of the domain. |
This enum defines the types of interpolation supported in this class.
Isis::NumericalApproximation::NumericalApproximation | ( | const NumericalApproximation::InterpType & | itype = CubicNatural | ) |
Default constructor creates NumericalApproximation object.
Sets the NumericalApproximation::InterpType to the enumerated value of itype. Default is CubicNatural.
itype | NumericalApproximation::InterpType enum value to be assigned to this object (Default: CubicNatural) |
Isis::IException::Programmer | "Unable to construct NumericalApproximation object" |
References _FILEINFO_, and Isis::IException::errorType().
Isis::NumericalApproximation::NumericalApproximation | ( | unsigned int | n, |
double * | x, | ||
double * | y, | ||
const NumericalApproximation::InterpType & | itype = CubicNatural |
||
) |
Constructs NumericalApproximation object, sets interpolation type, populates the data set with pointer values.
This constructor verifies the NumericalApproximation::InterpType passed in as a parameter. If it passes, it assigns this type to the object and x and y are added to the data set and validated. If the interpolation type is CubicClamped, first derivatives of the endpoints are required before Evaluate() may be called. If no value for itype is given, the default is CubicNatural.
n | Number of data points to be used. |
x | Array of domain values to add to data set |
y | Array of range values corresponding to each domain value in x |
itype | NumericalApproximation::InterpType enum value to be assigned to this object (Default: CubicNatural) |
Isis::IException::Programmer | "Unable to construct NumericalApproximation object" |
References _FILEINFO_, and Isis::IException::errorType().
Isis::NumericalApproximation::NumericalApproximation | ( | const vector< double > & | x, |
const vector< double > & | y, | ||
const NumericalApproximation::InterpType & | itype = CubicNatural |
||
) |
Constructs NumericalApproximation object, sets interpolation type, populates the data set with vector values.
This constructor verifies the interpolation type and checks that the x and y vectors have the same number of elements. If so, the object is assigned the NumericalApproximation::InterpType passed in as a parameter and x and y are added to the data set and validated. If the interpolation type is CubicClamped, first derivatives of the endpoints are required before evaluation is allowed. If no value for itype is given, the default is CubicNatural.
x | Vector of domain values to add to data set |
y | Vector of range values corresponding to each domain value in x |
itype | NumericalApproximation::InterpType enum value to be assigned to this object (Default: CubicNatural) |
Isis::IException::Programmer | "Unable to construct NumericalApproximation object" |
References _FILEINFO_, and Isis::IException::errorType().
Isis::NumericalApproximation::NumericalApproximation | ( | const NumericalApproximation & | oldObject | ) |
Copy constructor creates NumericalApproximation object and sets the copies the class variable values and state of input NumericalApproximation object to the new object.
Copying of class objects must be handled explicitly in order for it to be completed properly. The data, interpolator type and other class variables are copied from the supplied object.
oldObject | NumericalApproximation object to gleen data from for creation of this object |
Isis::IException::Programmer | "Unable to construct NumericalApproximation object" |
References _FILEINFO_, Isis::IException::errorType(), p_clampedComputed, p_clampedDerivFirstPt, p_clampedDerivLastPt, p_clampedEndptsSet, p_clampedSecondDerivs, p_dataValidated, p_fprimeOfx, p_itype, p_polyNevError, p_x, and p_y.
|
virtual |
Destructor deallocates memory being used.
Explicit operations are required to free resources of the interpolation
void Isis::NumericalApproximation::AddCubicHermiteDeriv | ( | unsigned int | n, |
double * | fprimeOfx | ||
) |
Adds values for the first derivatives of the data points.
This method can only be called for cubic Hermite splines, i.e. if NumericalApproximation::InterpType is CubicHermite. These values should be entered in the order of the corresponding data points.
n | Number of derivative values to be added. |
fprimeOfx | Array of derivative values to be added. |
Isis::IException::Programmer | "Method only used for cspline-Hermite interpolation type" |
References _FILEINFO_.
Referenced by Isis::SpicePosition::HermiteIndices(), and Isis::SpicePosition::SetEphemerisTimeHermiteCache().
void Isis::NumericalApproximation::AddCubicHermiteDeriv | ( | const vector< double > & | fprimeOfx | ) |
Adds values for the first derivatives of the data points.
This method can only be called for cubic Hermite splines, i.e. if NumericalApproximation::InterpType is CubicHermite. These values should be entered in the order of the corresponding data points.
fprimeOfx | Vector of derivative values to be added. |
Isis::IException::Programmer | "Method only used for cspline-Hermite interpolation type" |
References _FILEINFO_.
void Isis::NumericalApproximation::AddCubicHermiteDeriv | ( | const double | fprimeOfx | ) |
Adds value of a first derivative of a data point.
This method can only be called for cubic Hermite splines, i.e. if NumericalApproximation::InterpType is CubicHermite. Values should be entered in the order of the corresponding data points.
fprimeOfx | Derivative value to be added. |
Isis::IException::Programmer | "Method only used for cspline-Hermite interpolation type" |
References _FILEINFO_.
void Isis::NumericalApproximation::AddData | ( | const double | x, |
const double | y | ||
) |
Add a datapoint to the set.
This method allows the user to add a point to the set of data in preparation for interpolation over an interval.
Note: All data sets must have unique x values. If the interpolation type is not PolynomialNeville, x values must be sorted in ascending order. If CubicNatPeriodic interpolation type is used, the first and last points added to the data set must have the same y value.
x | Domain value to add to data set |
y | Range value corresponding to domain value x |
Referenced by EvaluateCubicNeighborhood(), Isis::SpicePosition::HermiteIndices(), Isis::LunarLambertEmpirical::LunarLambertEmpirical(), Isis::MinnaertEmpirical::MinnaertEmpirical(), Isis::rebin(), Isis::NumericalAtmosApprox::RombergsMethod(), RombergsMethod(), and Isis::SpicePosition::SetEphemerisTimeHermiteCache().
void Isis::NumericalApproximation::AddData | ( | unsigned int | n, |
double * | x, | ||
double * | y | ||
) |
Add set of points to the data set using pointer arrays.
This method allows the user to add a set of points to the set of data in preparation for interpolation over an interval. This method does not overwrite previously added data. Rather, it adds the new data to the end of the existing data set.
Note: Behavior is undefined if the domain is not sorted in ascending or descending order.
n | Number of data points to be used. |
x | Array of domain values to add to data set |
y | Array of range values corresponding to each domain value in x |
void Isis::NumericalApproximation::AddData | ( | const vector< double > & | x, |
const vector< double > & | y | ||
) |
Add set of points to the data set using vectors.
This method allows the user to add a set of points to the set of data in preparation for interpolation over an interval. This method does not overwrite previously added data. Rather, it adds the new data to the end of the existing data set.
Note: Behavior is undefined if the domain is not sorted in ascending or descending order.
x | Vector of domain values to add to data set |
y | Vector of range values corresponding to each domain value in x |
Isis::IException::Programmer | "Invalid arguments. The number of x-values does not equal the number of y-values." |
References _FILEINFO_.
double Isis::NumericalApproximation::BackwardFirstDifference | ( | const double | a, |
const unsigned int | n = 3 , |
||
const double | h = 0.1 |
||
) |
Uses an n point backward first difference formula to approximate the first derivative evaluated at a given domain value.
This method uses backward first difference formulas to return an approximation of the first derivative of the interpolated data set function evaluated at given a valid domain value, a. Backward difference formulas use n points, with the largest x value at a, for numerical differentiation approximation. This method uses one of the following formulas:
\[ f\prime(a) \approx \frac{1}{h}[f(a) - f(a - h)] \]
\[ f\prime(a) \approx \frac{1}{2h}[3f(a) - 4f(a - h) + f(a - 2h)] \]
a | Domain value at which first deriviative is evaluated. |
n | The number of points used in the formula |
h | Distance between nearest points in the formula |
Isis::IException::Programmer | "Invalid argument. Value entered is outside of domain." |
Isis::IException::Programmer | "Formula steps outside of domain." If a-(n-1)h is less than domain min. |
Isis::IException::Programmer | "Invalid argument. There is no n-point backward difference formula in use." @throws Isis::IException::Programmer "Unable to calculate backward first difference for given (a,n,h)" |
References _FILEINFO_, and Isis::IException::errorType().
double Isis::NumericalApproximation::BackwardSecondDifference | ( | const double | a, |
const unsigned int | n = 3 , |
||
const double | h = 0.1 |
||
) |
Uses an n point backward second difference formula to approximate the second derivative evaluated at a given domain value.
This method uses backward second difference formulas to return an approximation of the second derivative of the interpolated data set function evaluated at given a valid domain value, a. Backward second difference formulas use n points, with the largest x value at a, for numerical differentiation approximation. This method uses the following formula:
\[ f\prime\prime(a) \approx \frac{1}{h^2}[f(a) - 2f(a - h) + f(a - 2h)] \]
a | Domain value at which second deriviative is evaluated. |
n | The number of points used in the formula. |
h | Distance between nearest points in the formula. |
Isis::IException::Programmer | "Invalid argument. Value entered is outside of domain." |
Isis::IException::Programmer | "Formula steps outside of domain." If a-(n-1)h is less than domain min. |
Isis::IException::Programmer | "Invalid argument. There is no n-point backward difference formula in use." |
Isis::IException::Programmer |
References _FILEINFO_, and Isis::IException::errorType().
double Isis::NumericalApproximation::BoolesRule | ( | const double | a, |
const double | b | ||
) |
Uses Boole's Rule to approximate the integral of the interpolated data set function on the interval (a,b).
The Boole's Rule for integration approximation uses a 5-point Newton-Cote formula. This rule states:
\[ \int_{a}^b f(x)dx \approx \frac{2h}{45}[7f(a) + 32f(a+h) + 12f(a+2h) + 32f(a+3h) + 7f(a+4h)] \]
where h = (b - a)/4.
This method uses a composite, or extended, Boole's rule formula to approximate the integral.
Note: The method uses an algorithm that is adapted from the IDL function int_tabulated.pro.
a | Lower bound of interval. |
b | Upper bound of interval. |
Isis::IException::Programmer | "Unable to calculate the integral on the interval (a,b) using Boole's rule" |
References _FILEINFO_, and Isis::IException::errorType().
double Isis::NumericalApproximation::CenterFirstDifference | ( | const double | a, |
const unsigned int | n = 5 , |
||
const double | h = 0.1 |
||
) |
Uses an n point center first difference formula to approximate the first derivative evaluated at a given domain value.
This method uses center first difference formulas to return an approximation of the first derivative of the interpolated data set function evaluated at given a valid domain value, a. Center difference formulas use n points, centered at a, for numerical differentiation approximation. This method uses one of the following formulas:
\[ f\prime(a) \approx \frac{1}{2h}[f(a + h) - f(a - h)] \]
\[ f\prime(a) \approx \frac{1}{12h}[-f(a + 2h) + 8f(a +h) - 8f(a - h) + f(a - 2h)] \]
a | Domain value at which first deriviative is evaluated. |
n | The number of points used in the formula. |
h | Distance between nearest points in the formula. |
Isis::IException::Programmer | "Invalid argument. Value entered is outside of domain." |
Isis::IException::Programmer | "Formula steps outside of domain." If a+(n-1)h is greater than domain max or a-(n-1)h is less than domain min. |
Isis::IException::Programmer | "Invalid argument. There is no n-point center difference formula in use." |
Isis::IException::Programmer | "Unable to calculate center first difference for (a,n,h)" |
References _FILEINFO_, and Isis::IException::errorType().
double Isis::NumericalApproximation::CenterSecondDifference | ( | const double | a, |
const unsigned int | n = 5 , |
||
const double | h = 0.1 |
||
) |
Uses an n point center second difference formula to approximate the second derivative evaluated at a given domain value.
This method uses center second difference formulas to return an approximation of the second derivative of the interpolated data set function evaluated at given a valid domain value, a. Center second difference formulas use n points, centered at a, for numerical differentiation approximation. This method uses one of the following formulas:
\[ f\prime\prime(a) \approx \frac{1}{h^2}[f(a + h) - 2f(a) + f(a - h)] \]
\[ f\prime\prime(a) \approx \frac{1}{12h^2}[-f(a + 2h) + 16f(a +h) - 30f(a) + 16f(a - h) - f(a - 2h)] \]
a | Domain value at which second deriviative is evaluated. |
n | The number of points used in the formula. |
h | Distance between nearest points in the formula. |
Isis::IException::Programmer | "Invalid argument. Value entered is outside of domain." |
Isis::IException::Programmer | "Formula steps outside of domain." If a+(n-1)h is greater than domain max or a-(n-1)h is less than domain min. |
Isis::IException::Programmer | "Invalid argument. There is no n-point center difference formula in use." |
Isis::IException::Programmer | "Unable to calculate center second difference for (a,n,h)" |
References _FILEINFO_, and Isis::IException::errorType().
|
protected |
Computes the cubic clamped interpolation for a set of (x,y) data points, given the first derivatives of the endpoints of the data set.
This protected method is called only if the object is assigned a cubic clamped interpolation type and if it has not already been computed for the given data set. It calculates the second derivatives for each p_x value of the known data set and stores these values in p_clampedSecondDerivs so that the EvaluateCubicClamped() method may be called to interpolate the set using clamped boundary conditions, if possible. This method must be called when all the data points have been added to the object for the data set and after SetCubicClampedEndptDeriv() has been called. If the endpoint derivatives are greater than or equal to 1 x 1030, the routine is signaled to set the corresponding boundary condition for a natural spline, with zero second derivative on that boundary.
Issis::IException::Programmer | "Must use SetCubicClampedEndptDeriv() before computing cubic spline with clamped boundary" @throws Issis::IException::Programmer "Unable to compute cubic clamped interpolation" |
References _FILEINFO_, and Isis::IException::errorType().
|
protected |
Computes the GSL interpolation for a set of (x,y) data points.
This protected method is called only if the object is assigned a GSL interpolation type and if it has not already been computed on the given data set. It will compute the interval of interpolated range values over the given domain. A copy of this data is maintained in the object so the data points do not need to be interpolated on each evaluation of a point unless the data set changes.
Isis::IException::Programmer | "Unable to compute GSL interpolation" |
References _FILEINFO_, and Isis::IException::errorType().
bool Isis::NumericalApproximation::Contains | ( | double | x | ) |
Returns whether the passed value is an element of the set of x-values in the data set.
This method uses a binary search of the set of x-values to determine whether the input is contained in this set.
x | Value to search for in the data set. |
vector< double > Isis::NumericalApproximation::CubicClampedSecondDerivatives | ( | ) |
Retrieves the second derivatives of the data set.
This method returns a vector of the same size as the data set. Each component is the second derivative of the corresponding p_x value , as estimated by the ComputeCubicClamped() method.
Isis::IException::Programmer | "Method only used for cspline-clamped interpolation type" |
Isis::IException::Programmer | "Unable to calculate the second derivatives of the data set for a clamped cubic spline." |
References _FILEINFO_, and Isis::IException::errorType().
double Isis::NumericalApproximation::DomainMaximum | ( | ) |
Input data domain maximum value.
This method validates the data set for the assigned interpolation type and returns the largest value in p_x of the data set.
Isis::IException::Programmer | "Unable to calculate the domain minimum for the data set." |
References _FILEINFO_, and Isis::IException::errorType().
double Isis::NumericalApproximation::DomainMinimum | ( | ) |
Input data domain minimum value.
This method validates the data set for the assigned interpolation type and returns the smallest value in p_x of the data set.
Isis::IException::Programmer | "Unable to calculate the domain maximum for the data set." |
References _FILEINFO_, and Isis::IException::errorType().
double Isis::NumericalApproximation::Evaluate | ( | const double | a, |
const ExtrapType & | etype = ThrowError |
||
) |
Calculates interpolated or extrapolated value of tabulated data set for given domain value.
This method returns the approximate value for f(a), where a is an element near the domain and f is the approximated function for the given data set. If the given value, a, falls outside of the domain, then the extrapoltion type is examined to determine the result. CubicNeighborhood and GSL interpolation types can not extrapolate, so the user must choose to throw an error or return f evaluated at the nearest domain boundary. CubicClamped, CubicHermite, and PolynomialNeville interpolation types can extrapolate a value with accuracy only if a is near enough to the domain boundary. Default NumericalApproximation::ExtrapType is ThrowError.
a | Domain value from which to interpolate a corresponding range value |
etype | NumericalApproximation::ExtrapType enum value indicates how to evaluate if a falls outside of the domain. (Default: ThrowError) |
Isis::IException::Programmer | "Unable to evaluate the function at the point a" |
References _FILEINFO_, and Isis::IException::errorType().
Referenced by Isis::HapkeAtm1::AtmosModelAlgorithm(), Isis::HapkeAtm2::AtmosModelAlgorithm(), EvaluateCubicNeighborhood(), Isis::SpicePosition::HermiteIndices(), Isis::LunarLambertEmpirical::PhotoModelAlgorithm(), Isis::MinnaertEmpirical::PhotoModelAlgorithm(), Isis::rebin(), Isis::NumericalAtmosApprox::RombergsMethod(), RombergsMethod(), and Isis::SpicePosition::SetEphemerisTimeHermiteCache().
vector< double > Isis::NumericalApproximation::Evaluate | ( | const vector< double > & | a, |
const ExtrapType & | etype = ThrowError |
||
) |
Calculates interpolated value for given set of domain values.
This method returns the approximate values for f(ai), where ai is an element in the domain and f is the interpolated function for the given data set. If any of the given values in a fall outside of the domain, then the extrapoltion type is examined to determine the result. CubicNeighborhood and GSL interpolation types can not extrapolate, so the user must choose to throw an error or return f evaluated at the nearest domain boundary. CubicClamped and PolynomialNeville interpolation types can extrapolate a value with accuracy only if ai is near enough to the domain boundary.
a | Vector of domain values from which to interpolate a vector of corresponding range values |
etype | NumericalApproximation::ExtrapType enum value indicates how to evaluate any value of a that falls outside of the domain. (Default: ThrowError) |
Isis::IException::Programmer | "Unable to evaluate the function at the vector of points, a." |
References _FILEINFO_, and Isis::IException::errorType().
|
protected |
Performs cubic spline interpolation with clamped boundary conditions, if possible.
This is a protected method called by Evaluate() if the NumericalApproximation::InterpType is CubicClamped and SetCubicClampedEndptDeriv() has been called. It uses the second derivative vector, p_clampedSecondDerivs, to interpolate the value for f(a) using a clamped cubic spline formula,
Note: If the given value, a, falls outside of the domain, then extrapolation is attempted and the result is accurate only if a is near enough to the domain boundary.
a | Domain value from which to interpolate a corresponding range value |
|
protected |
Performs interpolation using the Hermite cubic polynomial.
This is a protected method called by Evaluate() if the NumericalApproximation::InterpType is CubicHermite. It returns an approximate value for f(a) by using the Hermite cubic polynomial, which uses Lagrange coefficient polynomials. The data points and each corresponding first derivative should have been already added.
a | Domain value from which to interpolate a corresponding range value |
Isis::IException::User | "Invalid arguments. The size of the first derivative vector does not match the number of (x,y) data points." |
References _FILEINFO_.
double Isis::NumericalApproximation::EvaluateCubicHermiteFirstDeriv | ( | const double | a | ) |
Approximates the first derivative of the data set function evaluated at the given domain value for CubicHermite interpolation type.
This method returns an approximation of the first derivative evaluated at given a valid domain value, a. It is able to extrapolate for values not far outside of the domain
References _FILEINFO_.
Referenced by Isis::SpicePosition::SetEphemerisTimeHermiteCache().
double Isis::NumericalApproximation::EvaluateCubicHermiteSecDeriv | ( | const double | a | ) |
Approximates the second derivative of the data set function evaluated at the given domain value for CubicHermite interpolation type.
This method returns an approximation of the second derivative evaluated at given a valid domain value, a. It is able to extrapolate for values not far outside of the domain
References _FILEINFO_.
|
protected |
Performs cubic spline interpolation for a neighborhood about a.
This is a protected method called by Evaluate() if the NumericalApproximation::InterpType is CubicNeighborhood. It uses an algorithm that is adapted from the IDL interpol.pro application using the "/spline" keyword on an irregular grid. This type of cubic spline fits a natural cubic spline to the 4-point neighborhood of known data points surrounding a. For example, suppose {x0, x1, ..., xn} is the array of known domain values in the data set and \( x_i \leq a < x_{i+1}\) for some i such that \( 0 \leq i \leq n \), then f(a) is evaluated by interpolating the natural cubic spline consisting of the data set {xi-1, xi, xi+1, xi+2} at a.
Note: If the given value, a, falls outside of the domain, then f evaluated at the nearest domain boundary is returned.
a | Domain value from which to interpolate a corresponding range value. |
Isis::IException::Programmer | "Unable to evaluate cubic neighborhood interpolation at a" |
References _FILEINFO_, AddData(), Isis::IException::errorType(), and Evaluate().
|
protected |
Performs cubic spline interpolations for neighborhoods about each value of a.
This is a protected method called by Evaluate() if the NumericalApproximation::InterpType is CubicNeighborhood. It uses an algorithm that is adapted from the IDL interpol.pro application using the "/spline" keyword on an irregular grid. For each component of a, this method fits a natural cubic spline using the 4-point neighborhood of known data points surrounding that component. For example, suppose {x0, x1, ..., xn} is the array of known domain values in the data set, then for each component of a, ak, in the domain, there is an i such that \( 0 \leq i \leq n \) and \( x_i \leq a_k < x_{i+1}\). Then, f(ak) is evaluated by interpolating the natural cubic spline consisting of the data set {xi-1, xi, xi+1, xi+2} at ak.
Note: If any given value, ai, falls outside of the domain, then f is evaluated at the nearest domain boundary.
a | Vector of domain values from which to interpolate a vector of corresponding range values |
etype | NumericalApproximation::ExtrapType enum value indicates how to evaluate if a falls outside of the domain. (Default: ThrowError) |
Isis::IException::Programmer | "Unable to evaluate cubic neighborhood interpolation at the values in the vector, a" |
References _FILEINFO_, AddData(), Isis::IException::errorType(), Evaluate(), and Reset().
|
protected |
Evaluates data set in order to have enough data points to approximate the function to be integrated.
a | Lower bound of interval. |
b | Upper bound of interval. |
n | Number of points used in Newton-Cotes formula. |
Isis::IException::Programmer | (When a > b) "Invalid interval entered" |
Isis::IException::Programmer | "Invalid arguments. Interval entered is not contained within domain" |
Isis::IException::Programmer | "Unable to evaluate the data set for integration" |
References _FILEINFO_, and Isis::IException::errorType().
|
protected |
Performs polynomial interpolation using Neville's algorithm.
This is a protected method called by Evaluate() if the NumericalApproximation::InterpType is PolynomialNeville. It uses Neville's algorithm for Lagrange polynomials. It returns a value f(a) and sets the error estimate p_polyNevError. After this method is called, the user may access this error estimate by calling PolynomialNevilleErrorEstimate().
Note: If the given value, a, falls outside of the domain, then extrapolation is attempted and the result is accurate only if a is near enough to the domain boundary.
a | Domain value from which to interpolate a corresponding range value |
|
protected |
Find the index of the x-value in the data set that is just below the input value, a.
This method is used by EvaluateCubicHermite(), EvaluateCubFirstDeriv() and EvaluateCubicHermiteSecDeriv() to determine in which interval of x-values a lies. It returns the index of the lower endpoint of the interval. If a is below the domain minimum, the method returns 0 as the lower index. If a is above the domain maximum, it returns the second to last index of the data set, Size()-2, as the lower index.
a | Domain value around which the interval lies |
double Isis::NumericalApproximation::ForwardFirstDifference | ( | const double | a, |
const unsigned int | n = 3 , |
||
const double | h = 0.1 |
||
) |
Uses an n point forward first difference formula to approximate the first derivative evaluated at a given domain value.
This method uses forward first difference formulas to return an approximation of the first derivative of the interpolated data set function evaluated at given a valid domain value, a. Forward difference formulas use n points, with the smallest x value at a, for numerical differentiation approximation. This method uses one of the following formulas:
\[ f\prime(a) \approx \frac{1}{h}[f(a + h) - f(a)] \]
\[ f\prime(a) \approx \frac{1}{2h}[-f(a + 2h) + 4f(a + h) - 3f(a)] \]
a | Domain value at which first deriviative is evaluated. |
n | The number of points used in the formula. |
h | Distance between nearest points in the formula. |
Isis::IException::Programmer | "Invalid argument. Value entered is outside of domain." |
Isis::IException::Programmer | "Formula steps outside of domain." If a+(n-1)h is greater than domain max. |
Isis::IException::Programmer | "Invalid argument. There is no n-point forward difference formula in use." |
Isis::IException::Programmer | "Unable to calculate forward first difference for (a,n,h)" |
References _FILEINFO_, and Isis::IException::errorType().
double Isis::NumericalApproximation::ForwardSecondDifference | ( | const double | a, |
const unsigned int | n = 3 , |
||
const double | h = 0.1 |
||
) |
Uses an n point forward second difference formula to approximate the second derivative evaluated at a given domain value.
This method uses forward second difference formulas to return an approximation of the second derivative of the interpolated data set function evaluated at given a valid domain value, a. Forward second difference formulas use n points, with the smallest x value at a, for numerical differentiation approximation. This method uses the following formula:
\[ f\prime\prime(a) \approx \frac{1}{h^2}[f(a + 2h) - 2f(a + h) + f(a)] \]
a | Domain value at which second deriviative is evaluated. |
n | The number of points used in the formula. |
h | Distance between nearest points in the formula. |
Isis::IException::Programmer | "Invalid argument. Value entered is outside of domain." |
Isis::IException::Programmer | "Formula steps outside of domain." If a+(n-1)h is greater than domain max. |
Isis::IException::Programmer | "Invalid argument. There is no n-point forward difference formula in use." |
Isis::IException::Programmer | "Unable to calculate forward second difference for (a,n,h)" |
References _FILEINFO_, and Isis::IException::errorType().
|
protected |
Allocates GSL interpolation functions.
This method is called within Compute to allocate pointers to a GSL interpolation object and to a GSL accelerator object (for interpolation look-ups). This is only called for the GSL interpolation types. If it is deemed invalid, an exception will be thrown.
npoints | Number of points to allocate for the GSL interpolator |
|
protected |
Returns whether a GSL interpolation computation of the data set has been performed.
This method is only applicable to GSL interpolation types.
Isis::IException::Programmer | "Method only valid for GSL interpolation types" |
References _FILEINFO_.
|
protected |
Deallocate GSL interpolator resources, if used.
If a GSL interpolator function has been allocated, this routine will free its resources and reset internal pointers to reflect this state and provide a mechanism to test its state. Otherwise, this method sets the internal GSL pointers to their default, 0.
double Isis::NumericalApproximation::GslFirstDerivative | ( | const double | a | ) |
Approximates the first derivative of the data set function evaluated at the given domain value for GSL supported interpolation types.
This method returns an approximation of the first derivative evaluated at given a valid domain value, a. It is a wrapper for the GSL subroutine gsl_spline_eval_deriv_e(). If the NumericalApproximation::InterpType is not a GSL type, then this method throws an error. No documentation was found concerning the algorithm used by this method.
a | Domain value at which first deriviative is evaluated. |
Isis::IException::Programmer | "Invalid argument. Value entered is outside of domain." |
Isis::IException::Programmer | "Cannot use this method for interpolation type" (If the interpolation type is not GSL) |
Isis::IException::Programmer | "Unable to compute the first derivative at a using the GSL interpolation" @throws Isis::IException::Programmer "Unable to compute the first derivative at a. GSL integrity check failed" |
References _FILEINFO_, and Isis::IException::errorType().
|
protected |
Search for a GSL interpolation function.
This method searches the supported GSL options table for a given interpolation function as requested by the caller. If it is not found, an exception will be thrown indicating the erroneous request.
itype | Type of GSL interpolator to find |
Isis::IException::Programmer | "Invalid argument. Unable to find GSL interpolator" |
References _FILEINFO_.
Referenced by MinPoints().
double Isis::NumericalApproximation::GslIntegral | ( | const double | a, |
const double | b | ||
) |
Approximates the integral of the data set function evaluated on the given interval for GSL supported interpolation types.
This method returns an approximation of the integral evaluated on a given valid domain interval, (a,b). It is a wrapper for the GSL subroutine gsl_spline_eval_integ_e(). If the NumericalApproximation::InterpType is not a GSL type, then this method throws an error. No documentation was found concerning the algorithm used by this method.
a | Lower endpoint at which integral is evaluated. |
b | Upper endpoint at which integral is evaluated. |
Isis::IException::Programmer | (When a > b) "Invalid interval entered" |
Isis::IException::Programmer | "Invalid arguments. Interval entered is not contained within domain" |
Isis::IException::Programmer | "Cannot use this method for interpolation type" (If the interpolation type is not GSL) |
Isis::IException::Programmer | "Unable to compute the integral on the interval (a,b) using GSL interpolation" @throws Isis::IException::Programmer "Unable to compute the integral on the interval (a,b). GSL integrity check failed." |
References _FILEINFO_, and Isis::IException::errorType().
|
protected |
Checks the status of the GSL interpolation operations.
This method takes a return status from a GSL call and determines if it is completed successfully. This implementation currently allows the GSL_DOMAIN error to propagate through sucessfully as the domain can be checked by the caller if they deem this necessary.
gsl_status | Return status of a GSL function call |
src | Name of the calling source invoking the check. This allows more precise determination of the error. |
line | Line of the calling source that invoked the check |
Isis::IException::Programmer | "GSL error occured" |
|
protected |
Returns whether an interpolation type is adapted from the GSL library.
GSL interpolation types include the following:
itype | Interpolation type to be compared to GSL list |
double Isis::NumericalApproximation::GslSecondDerivative | ( | const double | a | ) |
Approximates the second derivative of the interpolated data set function evaluated at the given domain value for GSL supported interpolation types.
This method returns an approximation of the second derivative evaluated at given a valid domain value, a. It is a wrapper for the GSL subroutine gsl_spline_eval_deriv2_e(). If the NumericalApproximation::InterpType is not a GSL type, then this method throws an error. No documentation was found concerning the algorithm used by this method.
a | Domain value at which second deriviative is evaluated. |
Isis::IException::Programmer | "Invalid argument. Value entered is outside of domain." |
Isis::IException::Programmer | "Cannot use this method for interpolation type" (If the interpolation type is not GSL) |
Isis::IException::Programmer | "Unable to compute the second derivative at a using GSL interpolation" |
Isis::IException::Programmer | "Unable to compute the second derivative at a. GSL integrity check failed" |
References _FILEINFO_, and Isis::IException::errorType().
|
protected |
Initializes the object upon instantiation.
This method sets up the initial state of the object, typically at instantiation. It populates the interpolation function table identifying which options are available to the users of this class.
GSL error handling is turned off - upon repeated instantiation of this object. The GSL default error handling, termination of the application via an abort() is unacceptable. This calls adopts an alternative policy provided by the GSL whereby error checking must be done by the calling environment. This has an unfortunate drawback in that it is not enforceable in an object oriented environment that utilizes the GSL in disjoint classes.
itype | NumericalApproximation::InterpType enum value to be assigned to this object |
Isis::IException::Programmer | "Unable to initialize NumericalApproximation object" |
References _FILEINFO_, and Isis::IException::errorType().
|
protected |
Returns whether the passed value is greater than or equal to the domain minimum and less than or equal to the domain maximum.
a | Value to be verified as valid domain value. |
Isis::IException::Programmer | "Unable to compute domain boundaries" |
References _FILEINFO_, and Isis::IException::errorType().
|
inline |
Returns the enumerated type of interpolation chosen.
This method can be selected after all the points are added in the interpolation by using Compute() method. Note that this prints out as an integer representaion of the enumerated type:
int Isis::NumericalApproximation::MinPoints | ( | ) |
Minimum number of points required by interpolating function.
This method returns the minimum number of points that are required by the interpolating function in order for it to be applied to a data set. It returns the number of the of the currently selected/active interpolation function.
Isis::IException::Programmer | "Unable to calculate minimum points." |
References _FILEINFO_, and Isis::IException::errorType().
int Isis::NumericalApproximation::MinPoints | ( | NumericalApproximation::InterpType | itype | ) |
Minimum number of points required by interpolating function.
This method returns the minimum number of points that are required by the interpolating function in order for it to be applied to a data set. It returns the minimum number for the specifed interpolation function as specified by the caller.
itype | Type of interpolation function for which to return minimum points |
Isis::IException::Programmer | "Invalid argument. Unknown interpolation type" |
Isis::IException::Programmer | "Unable to calculate minimum points" |
References _FILEINFO_, Isis::IException::errorType(), and GslFunctor().
string Isis::NumericalApproximation::Name | ( | ) | const |
Get name of interpolating function assigned to object.
This method returns the name of the interpolation type that is currently assigned to this object as a std::string. This may be called without adding any points. If called before computation, the result reflects the name of the function chosen at instantiation.
Isis::IException::Programmer | "Unable to retrieve numerical approximation name", |
References _FILEINFO_, and Isis::IException::errorType().
NumericalApproximation & Isis::NumericalApproximation::operator= | ( | const NumericalApproximation & | oldObject | ) |
NumericalApproximation assigmment operator sets this object "equal to" another.
Assigning one object to the other requires a deallocation of the current state of this object and reinitialization of data from the right hand object. This may be a bit costly if the number of values is large.
oldObject | NumericalApproximation object to copy data from for initialization |
Isis::IException::Programmer | "Unable to copy NumericalApproximation object" |
References _FILEINFO_, Isis::IException::errorType(), p_clampedDerivFirstPt, p_clampedDerivLastPt, p_clampedSecondDerivs, p_dataValidated, p_fprimeOfx, p_itype, p_polyNevError, p_x, and p_y.
vector< double > Isis::NumericalApproximation::PolynomialNevilleErrorEstimate | ( | ) |
Retrieves the error estimate for the Neville's polynomial interpolation type.
This method must be called after the Evaluate() method has been invoked. If the Evaluate() method is passed a vector, this error will contain an error estimate for each of the values of the passed in vector. Otherwise, it will contain a single value.
Isis::IException::Programmer | "Method only used for polynomial-Neville's interpolation type" |
Isis::IException::Programmer | "Error not calculated" |
References _FILEINFO_.
Referenced by Isis::NumericalAtmosApprox::RombergsMethod(), and RombergsMethod().
double Isis::NumericalApproximation::RefineExtendedTrap | ( | double | a, |
double | b, | ||
double | s, | ||
unsigned int | n | ||
) |
Calculates refinements extended trapezoidal rule to approximate the integral of the interpolated data set function on the interval (a,b).
This method calculates the nth stage of refinement of an extended trapezoidal rule. When called with n = 1, the method returns the non-composite trapezoidal estimate of the integral. Subsequent calls with n = 2, 3, ... (in sequential order) will improve the accuracy by adding 2n-2 additional interior points. This method can be used to integrate by the extended trapeziodal rule if you know the number of steps you want to take. For example, if you want 2M + 1, use the following code:
Note: Although this method may be used to approximate an integral, as described above, it is most often used by RombergsMethod() to integrate.
a | Lower limit of integration |
b | Upper limit of integration |
s | Previous value of refinement |
n | Number of partitions to use when integrating |
Isis::IException::Programmer | "Unable to calculate the integral on the interval (a,b) using the extended trapezoidal rule" |
References _FILEINFO_, and Isis::IException::errorType().
|
protected |
Generalized error report generator.
This method is used throughout this class to standardize error reporting as a convenience to its implementor.
type | Type of Isis::IException |
methodName | Name of method where exception originated |
message | Error context string provided by the caller |
filesrc | Name of the file the error occured in. |
lineno | Line number of the calling source where the error occured. |
Isis::IException::errType | equal to type and error message equal to (methodName + " - " + message) |
void Isis::NumericalApproximation::Reset | ( | ) |
Resets the state of the object.
This method deallocates the internal state of the object and clears the data set and class variables. The object is returned to its original state, so data points must be entered before computing again. This does not clear or reset the interpolation type.
Referenced by Isis::HapkeAtm1::AtmosModelAlgorithm(), Isis::HapkeAtm2::AtmosModelAlgorithm(), EvaluateCubicNeighborhood(), Isis::AtmosModel::GenerateHahgTables(), Isis::AtmosModel::GenerateHahgTablesShadow(), Isis::LunarLambertEmpirical::LunarLambertEmpirical(), Isis::MinnaertEmpirical::MinnaertEmpirical(), Isis::NumericalAtmosApprox::RombergsMethod(), RombergsMethod(), Isis::LunarLambertEmpirical::~LunarLambertEmpirical(), and Isis::MinnaertEmpirical::~MinnaertEmpirical().
void Isis::NumericalApproximation::Reset | ( | NumericalApproximation::InterpType | itype | ) |
Resets the state of the object and resets interpolation type.
This method will clear the data set, reset the validation status of the data to false, reset the interpolation type, clear class variables and deallocate (inactivate) the internal state of the object. The object is returned to its original state, so data points must be entered before computing again.
itype | NumericalApproximation::InterpType enum value to be assigned to this object |
Isis::IException::Programmer | "Unable to reset interpolation type" |
References _FILEINFO_, and Isis::IException::errorType().
double Isis::NumericalApproximation::RombergsMethod | ( | double | a, |
double | b | ||
) |
Uses Romberg's method to approximate the integral of the interpolated data set function on the interval (a,b).
This method returns the integral of the function from a to b using Romberg's method for Numerical Integration of order 2K, where, e.g., K=2 is simpson's rule. This is a generalization of the trapezoidal rule. Romberg Integration uses a series of refinements on the extended (or composite) trapezoidal rule to reduce error terms. This method makes use of Neville's algorithm for polynomial interpolation to extrapolate successive refinements.
a | Lower limit of integration |
b | Upper limit of integration |
Isis::IException::Programmer | "Failed to converge." |
Isis::IException::Programmer | "Unable to calculate the integral on the interval (a,b) using Romberg's method" |
References _FILEINFO_, AddData(), Isis::IException::errorType(), Evaluate(), PolynomialNevilleErrorEstimate(), and Reset().
void Isis::NumericalApproximation::SetCubicClampedEndptDeriv | ( | const double | yp1, |
const double | ypn | ||
) |
Sets the values for the first derivatives of the endpoints of the data set.
This method can only be called for cubic splines with clamped boundary conditions, i.e. if NumericalApproximation::InterpType is CubicClamped.
yp1 | First derivative of the function evaluated at the domain minimum. |
ypn | First derivative of the function evaluated at the domain maximum. |
Isis::IException::Programmer | "Method only used for cspline-clamped interpolation type" |
References _FILEINFO_.
Referenced by Isis::LunarLambertEmpirical::LunarLambertEmpirical(), and Isis::MinnaertEmpirical::MinnaertEmpirical().
void Isis::NumericalApproximation::SetInterpType | ( | NumericalApproximation::InterpType | itype | ) |
Sets interpolation type.
This method will reset the interpolation type to that of the input parameter and deallocate memory. Unlike the Reset() method, this method does NOT discard the data set. However, since interpolations have differing requirements for valid data sets, the data set stored is not yet validated for the new interpolation type. Other class variables are cleared if they are interpolation type dependent.
itype | NumericalApproximation::InterpType enum value to be assigned to this object |
Isis::IException::Programmer | "Invalid argument. Unknown interpolation type" |
Isis::IException::Programmer | "Unable to set interpolation type" |
References _FILEINFO_, and Isis::IException::errorType().
Referenced by Isis::LunarLambertEmpirical::LunarLambertEmpirical(), and Isis::MinnaertEmpirical::MinnaertEmpirical().
double Isis::NumericalApproximation::Simpsons3PointRule | ( | const double | a, |
const double | b | ||
) |
Uses Simpson's 3-point rule to approximate the integral of the interpolated data set function on the interval (a,b).
The Simpson's 3-Point Rule for numerical integration uses a 3-point Newton-Cote formula. This rule states:
\[ \int_{a}^b f(x)dx \approx \frac{h}{3}[f(a) + f(a+h) + 4*f(a+2h)] \]
where h = (b - a)/2. This method uses a composite, or extended, Simpson's rule algorithm to approximate the integral.
a | Lower bound of interval. |
b | Upper bound of interval. |
Isis::IException::Programmer | "Unable to calculate the integral on the interval (a,b) using Simpson's 3 point rule" |
References _FILEINFO_, and Isis::IException::errorType().
double Isis::NumericalApproximation::Simpsons4PointRule | ( | const double | a, |
const double | b | ||
) |
Uses Simpson's 4-point rule to approximate the integral of the interpolated data set function on the interval (a,b).
The Simpson's 4-point Rule for numerical integration uses a 4-point Newton-Cote formula and is sometimes called the Simpson's 3/8 Rule. This rule states:
\[ \int_{a}^b f(x)dx \approx \frac{3h}{8}[f(a) + 3f(a+h) + 3f(a+2h) + f(a+3h)] \]
where h = (b - a)/3. This method uses a composite, or extended, Simpson's 3/8 rule algorithm to approximate the integral.
a | Lower bound of interval. |
b | Upper bound of interval. |
Isis::IException::Programmer | "Unable to calculate the integral on the interval (a,b) using Simpson's 4 point rule" |
References _FILEINFO_, and Isis::IException::errorType().
|
inline |
Returns the number of the coordinates added to the data set.
double Isis::NumericalApproximation::TrapezoidalRule | ( | const double | a, |
const double | b | ||
) |
Uses the trapezoidal rule to approximate the integral of the interpolated data set function on the interval (a, b).
The trapeziod rule for integration approximation uses a 2-point Newton-Cote formula. This rule states:
\[ \int_{a}^b f(x)dx \approx \frac{h}{2}[f(a) + f(b)] \]
where h = b - a. This method uses a composite, or extended, trapeziodal rule algorithm to approximate the integral.
a | Lower bound of interval. |
b | Upper bound of interval. |
Isis::IException::Programmer | "Unable to calculate the integral on the interval (a,b) using the trapezoidal rule" |
References _FILEINFO_, and Isis::IException::errorType().
|
protected |
Validates the data set before computing interpolation.
This method is called from the ComputeCubicClamped() and ComputeGsl() methods to verify that the data set contains the minimum number of required points and that the components of the vector of domain values are unique. For all interpolation types other than polynomial-Neville's, the method verifies that the vector of domain values are also sorted in ascending order. For CubicNatPeriodic interpolation type, this method verifies that the first and last p_y values are equal, i.e. \( f(x_0) = f(x_{n-1}) \).
Isis::IException::Programmer | "Interpolation requires a minimim of data points" |
Isis::IException::Programmer | "Invalid data set, x-values must be unique" |
Isis::IException::Programmer | "Invalid data set, x-values must be in ascending order" (if interpolation type is not polynomial-Neville's) |
Isis::IException::Programmer | "First and last points of the data set must have the same y-value" (if interpolation type is cubic-periodic) |
References _FILEINFO_.
|
protected |
Returns the domain value at which to evaluate.
This protected method is called by Evaluate() if a falls outside of the domain of p_x values in the data set. The return value is determined by the NumericalApproximation::ExtrapType. If it is ThrowError, an error is thrown indicating that a is out of the domain. If it is NearestEndpoint, then the nearest domain boundary value is returned. Otherwise, i.e. Extrapolate, a is returned as long as the NumericalApproximation::InterpType is not GSL or cubic neighborhood.
a | Value passed into Evaluate() that falls outside of domain. |
etype | NumericalApproximation::ExtrapType enum value indicates how to evaluate if a falls outside of the domain. (Default: ThrowError) |
Isis::IException::Programmer | "Invalid argument. Value entered is outside of domain |
Isis::IException::Programmer | "Invalid argument. Cannot extrapolate for interpolation type" (GSL or CubicNeighborhood) |
References _FILEINFO_.
|
protected |
Lookup accelorator.
|
protected |
Flag variable to determine whether ComputeCubicClamped() has been called.
Referenced by NumericalApproximation().
|
protected |
First derivative of first x-value, p_x[0]. This is only used for the CubicClamped interpolation type.
Referenced by NumericalApproximation(), and operator=().
|
protected |
First derivative of last x-value, p_x[n-1]. This is only used for the CubicClamped interpolation type.
Referenced by NumericalApproximation(), and operator=().
|
protected |
Flag variable to determine whether SetCubicClampedEndptDeriv() has been called after all data was added for CubicClamped interpolation.
Referenced by NumericalApproximation().
|
protected |
List of second derivatives evaluated at p_x values. This is only used for the CubicClamped interpolation type.
Referenced by NumericalApproximation(), and operator=().
|
protected |
Flag variable to determine whether ValidateDataSet() has been called.
Referenced by NumericalApproximation(), and operator=().
|
protected |
List of first derivatives corresponding to each x value in the data set (i.e. each value in p_x)
Referenced by NumericalApproximation(), and operator=().
|
protected |
Currently active interpolator.
|
staticprotected |
Maintains list of interpolator options.
|
protected |
Interpolation type.
Referenced by NumericalApproximation(), and operator=().
|
protected |
Estimate of error for interpolation evaluated at x. This is only used for the PolynomialNeville interpolation type. 91 taken from AtmosModel.
Referenced by NumericalApproximation(), and operator=().
|
protected |
List of X values.
Referenced by NumericalApproximation(), and operator=().
|
protected |
List of Y values.
Referenced by NumericalApproximation(), and operator=().