Isis Developer Reference
Isis::NumericalApproximation Class Reference

NumericalApproximation provides various numerical analysis methods of interpolation, extrapolation and approximation of a tabulated set of x, y data. More...

#include <NumericalApproximation.h>

Inheritance diagram for Isis::NumericalApproximation:
Inheritance graph
Collaboration diagram for Isis::NumericalApproximation:
Collaboration graph

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.
 
 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.
 
 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.
 
 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.
 
NumericalApproximationoperator= (const NumericalApproximation &numApMeth)
 NumericalApproximation assigmment operator sets this object "equal to" another.
 
virtual ~NumericalApproximation ()
 Destructor deallocates memory being used.
 
string Name () const
 Get name of interpolating function assigned to object.
 
InterpType InterpolationType ()
 Returns the enumerated type of interpolation chosen.
 
int MinPoints ()
 Minimum number of points required by interpolating function.
 
int MinPoints (NumericalApproximation::InterpType itype)
 Minimum number of points required by interpolating function.
 
void AddData (const double x, const double y)
 Add a datapoint to the set.
 
void AddData (unsigned int n, double *x, double *y)
 Add set of points to the data set using pointer arrays.
 
void AddData (const vector< double > &x, const vector< double > &y)
 Add set of points to the data set using vectors.
 
void SetCubicClampedEndptDeriv (const double yp1, const double ypn)
 Sets the values for the first derivatives of the endpoints of the data set.
 
void AddCubicHermiteDeriv (unsigned int n, double *fprimeOfx)
 Adds values for the first derivatives of the data points.
 
void AddCubicHermiteDeriv (const vector< double > &fprimeOfx)
 Adds values for the first derivatives of the data points.
 
void AddCubicHermiteDeriv (const double fprimeOfx)
 Adds value of a first derivative of a data point.
 
double DomainMinimum ()
 Input data domain minimum value.
 
double DomainMaximum ()
 Input data domain maximum value.
 
bool Contains (double x)
 Returns whether the passed value is an element of the set of x-values in the data set.
 
unsigned int Size ()
 Returns the number of the coordinates added to the data set.
 
double Evaluate (const double a, const ExtrapType &etype=ThrowError)
 Calculates interpolated or extrapolated value of tabulated data set for given domain value.
 
vector< double > Evaluate (const vector< double > &a, const ExtrapType &etype=ThrowError)
 Calculates interpolated value for given set of domain values.
 
vector< double > PolynomialNevilleErrorEstimate ()
 Retrieves the error estimate for the Neville's polynomial interpolation type.
 
vector< double > CubicClampedSecondDerivatives ()
 Retrieves the second derivatives of the data set.
 
double EvaluateCubicHermiteFirstDeriv (const double a)
 Approximates the first derivative of the data set function evaluated at the given domain value for CubicHermite interpolation type.
 
double EvaluateCubicHermiteSecDeriv (const double a)
 Approximates the second derivative of the data set function evaluated at the given domain value for CubicHermite interpolation type.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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).
 
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).
 
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).
 
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).
 
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).
 
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).
 
void Reset ()
 Resets the state of the object.
 
void Reset (NumericalApproximation::InterpType itype)
 Resets the state of the object and resets interpolation type.
 
void SetInterpType (NumericalApproximation::InterpType itype)
 Sets interpolation type.
 

Protected Types

typedef const gsl_interp_type * InterpFunctor
 GSL Interpolation specs.
 
typedef map< InterpType, InterpFunctorFunctorList
 Set up a std::map of GSL interpolator functors. List of function types.
 
typedef FunctorList::const_iterator FunctorConstIter
 GSL Iterator.
 

Protected Member Functions

void Init (NumericalApproximation::InterpType itype)
 Initializes the object upon instantiation.
 
bool GslInterpType (NumericalApproximation::InterpType itype) const
 Returns whether an interpolation type is adapted from the GSL library.
 
void GslAllocation (unsigned int npoints)
 Allocates GSL interpolation functions.
 
void GslDeallocation ()
 Deallocate GSL interpolator resources, if used.
 
InterpFunctor GslFunctor (NumericalApproximation::InterpType itype) const
 Search for a GSL interpolation function.
 
void GslIntegrityCheck (int gsl_status, const char *src, int line)
 Checks the status of the GSL interpolation operations.
 
void ValidateDataSet ()
 Validates the data set before computing interpolation.
 
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.
 
bool GslComputed () const
 Returns whether a GSL interpolation computation of the data set has been performed.
 
void ComputeGsl ()
 Computes the GSL interpolation for a set of (x,y) data points.
 
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.
 
double ValueToExtrapolate (const double a, const ExtrapType &etype)
 Returns the domain value at which to evaluate.
 
double EvaluateCubicNeighborhood (const double a)
 Performs cubic spline interpolation for a neighborhood about a.
 
vector< double > EvaluateCubicNeighborhood (const vector< double > &a, const ExtrapType &etype)
 Performs cubic spline interpolations for neighborhoods about each value of a.
 
double EvaluateCubicClamped (const double a)
 Performs cubic spline interpolation with clamped boundary conditions, if possible.
 
double EvaluateCubicHermite (const double a)
 Performs interpolation using the Hermite cubic polynomial.
 
double EvaluatePolynomialNeville (const double a)
 Performs polynomial interpolation using Neville's algorithm.
 
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.
 
int FindIntervalLowerIndex (const double a)
 Find the index of the x-value in the data set that is just below the input value, a.
 
void ReportException (IException::ErrorType type, const string &method, const string &message, const char *filesrc, int lineno) const
 Generalized error report generator.
 

Protected Attributes

InterpType p_itype
 Interpolation type.
 
vector< double > p_x
 List of X values.
 
vector< double > p_y
 List of Y values.
 
bool p_dataValidated
 Flag variable to determine whether ValidateDataSet() has been called.
 
gsl_interp_accel * p_acc
 Lookup accelorator.
 
gsl_spline * p_interp
 Currently active interpolator.
 
bool p_clampedEndptsSet
 Flag variable to determine whether SetCubicClampedEndptDeriv() has been called after all data was added for CubicClamped interpolation.
 
bool p_clampedComputed
 Flag variable to determine whether ComputeCubicClamped() has been called.
 
double p_clampedDerivFirstPt
 First derivative of first x-value, p_x[0]. This is only used for the CubicClamped interpolation type.
 
double p_clampedDerivLastPt
 First derivative of last x-value, p_x[n-1]. This is only used for the CubicClamped interpolation type.
 
vector< double > p_clampedSecondDerivs
 List of second derivatives evaluated at p_x values. This is only used for the CubicClamped interpolation type.
 
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.
 
vector< double > p_fprimeOfx
 List of first derivatives corresponding to each x value in the data set (i.e. each value in p_x)
 

Static Protected Attributes

static FunctorList p_interpFunctors
 Maintains list of interpolator options.
 

Detailed Description

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 Isis 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.

Differentiation methods require the parameters a (domain value at which the derivative will be evaluated), n (number of points used in the difference formula), and h (step size of points). Note: a is denoted x0 in the above formulas.
Numerical Differentiation Type Difference Formulas Methods Available Description
Backward Difference
  • 2-point backward difference.

    \[ f\prime(x_0) \approx \frac{1}{h}[f(x_0) - f(x_0 - h)] \]

  • 3-point backward difference.

    \[ f\prime(x_0) \approx \frac{1}{2h}[3f(x_0) - 4f(x_0 - h) + f(x_0 - 2h)] \]

  • 3-point backward second difference.

    \[ f\prime\prime(x_0) \approx \frac{1}{h^2}[f(x_0) - 2f(x_0 - h) + f(x_0 - 2h)] \]

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
  • 2-point forward difference.

    \[ f\prime(x_0) \approx \frac{1}{h}[f(x_0 + h) - f(x_0)] \]

  • 3-point forward difference.

    \[ f\prime(x_0) \approx \frac{1}{2h}[-f(x_0 + 2h) + 4f(x_0 + h) - 3f(x_0)] \]

  • 3-point forward second difference.

    \[ f\prime\prime(x_0) \approx \frac{1}{h^2}[f(x_0 + 2h) - 2f(x_0 + h) + f(x_0)] \]

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
  • 3-point center difference.

    \[ f\prime(x_0) \approx \frac{1}{2h}[f(x_0 + h) - f(x_0 - h)] \]

  • 5-point center difference.

    \[ f\prime(x_0) \approx \frac{1}{12h}[-f(x_0 + 2h) + 8f(x_0 +h) - 8f(x_0 - h) + f(x_0 - 2h)] \]

  • 3-point center second difference.

    \[ f\prime\prime(x_0) \approx \frac{1}{h^2}[f(x_0 + h) - 2f(x_0) + f(x_0 - h)] \]

  • 5-point center second difference.

    \[ f\prime\prime(x_0) \approx \frac{1}{12h^2}[-f(x_0 + 2h) + 16f(x_0 +h) - 30f(x_0) + 16f(x_0 - h) - f(x_0 - 2h)] \]

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
  • (i - (n - 1)/2)h and evaluating them at 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.

Integration methods require the parameters a, b (interval of domain over which to integrate).
Numerical Integration Type Integration Formulas Methods Available Description
Trapezoidal Rule
  • 2-point Newton-Cotes trapezoidal rule:

    \[ \int_{a}^b f(x)dx \approx \frac{h}{2}[f(a) + f(b)] \]

    where h = b - a.
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
  • 3-point Newton-Cotes, Simpson's Rule

    \[ \int_{a}^b f(x)dx \approx \frac{h}{3}[f(a) + 4f(a+h) + f(a+2h)] \]

    where h = (b - a)/2.
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
  • 4-point Newton-Cotes, Simpson's 3/8 Rule

    \[ \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.
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
  • 5-point Newton-Cotes, Boole's Rule

    \[ \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.
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
  • Extended closed trapezoidal rule

    \[ \int_{a}^b f(x)dx \approx h[\frac12 f(x_0) + f(x_1) + ... + f(x_{n-2} + \frac12 f(x_{n-1})] \]

    where h = (b - a)/4, x0 = a, and xn-1 = b.
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 for Integration
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:

inline double f(double x) { return (x + cos ((x * x))); }
for (int x = 0; x < 200 ; x++) {
spline.addPoint(x, f(x));
}
spline.SetCubicClampedEndptDeriv(0,0);
double yinterp = spline.Evaluate(50.7);
double yextrap =
double deriv = spline.CenterFirstDifference(10);
double integ =
spline.RombergsMethod(spline.DomainMinimum(),spline.DomainMaximum());
NumericalApproximation provides various numerical analysis methods of interpolation,...
Definition NumericalApproximation.h:726
@ CubicClamped
Cubic Spline interpolation with clamped boundary conditions.
Definition NumericalApproximation.h:735
@ Extrapolate
Evaluate() attempts to extrapolate if a is outside of the domain. This is only valid for NumericalApp...
Definition NumericalApproximation.h:814

To compute the same data set using the Akima spline, just use the following:

double yinterp = spline.Evaluate(50.7);
double yextrap =
double deriv = spline.CenterFirstDifference(10);
double integ =
spline.RombergsMethod(spline.DomainMinimum(),spline.DomainMaximum());
@ Akima
Non-rounded Akima Spline interpolation with natural boundary conditions.
Definition NumericalApproximation.h:739
@ NearestEndpoint
Evaluate() returns the y-value of the nearest endpoint if a is outside of the domain.
Definition NumericalApproximation.h:817

Caveats

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.

Author
2008-11-05 Janet Barrett, Kris Becker, K Teal Thompson, Jeannie Walldren

Member Typedef Documentation

◆ FunctorConstIter

typedef FunctorList::const_iterator Isis::NumericalApproximation::FunctorConstIter
protected

GSL Iterator.

◆ FunctorList

Set up a std::map of GSL interpolator functors. List of function types.

◆ InterpFunctor

typedef const gsl_interp_type* Isis::NumericalApproximation::InterpFunctor
protected

GSL Interpolation specs.

Member Enumeration Documentation

◆ ExtrapType

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.

◆ InterpType

This enum defines the types of interpolation supported in this class.

Enumerator
Linear 

Linear interpolation.

Polynomial 

Polynomial interpolation.

PolynomialNeville 

Polynomial interpolation using Neville's algorithm.

CubicNatural 

Cubic Spline interpolation with natural boundary conditions.

CubicClamped 

Cubic Spline interpolation with clamped boundary conditions.

CubicNatPeriodic 

Cubic Spline interpolation with periodic boundary conditions.

CubicNeighborhood 

Cubic Spline interpolation using 4-pt Neighborhoods with natural boundary conditions.

CubicHermite 

Cubic Spline interpolation using the Hermite cubic polynomial.

Akima 

Non-rounded Akima Spline interpolation with natural boundary conditions.

AkimaPeriodic 

Non-rounded Akima Spline interpolation with periodic boundary conditions.

Constructor & Destructor Documentation

◆ NumericalApproximation() [1/4]

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.

Parameters
itypeNumericalApproximation::InterpType enum value to be assigned to this object (Default: CubicNatural)
Exceptions
Isis::IException::Programmer"Unable to construct NumericalApproximation object"

References _FILEINFO_, Init(), and p_dataValidated.

◆ NumericalApproximation() [2/4]

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.

Parameters
nNumber of data points to be used.
xArray of domain values to add to data set
yArray of range values corresponding to each domain value in x
itypeNumericalApproximation::InterpType enum value to be assigned to this object (Default: CubicNatural)
Exceptions
Isis::IException::Programmer"Unable to construct NumericalApproximation object"

References _FILEINFO_, AddData(), Init(), and ValidateDataSet().

◆ NumericalApproximation() [3/4]

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.

Parameters
xVector of domain values to add to data set
yVector of range values corresponding to each domain value in x
itypeNumericalApproximation::InterpType enum value to be assigned to this object (Default: CubicNatural)
Exceptions
Isis::IException::Programmer"Unable to construct NumericalApproximation object"

References _FILEINFO_, AddData(), Init(), and ValidateDataSet().

◆ NumericalApproximation() [4/4]

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.

Parameters
oldObjectNumericalApproximation object to gleen data from for creation of this object
Exceptions
Isis::IException::Programmer"Unable to construct NumericalApproximation object"

References _FILEINFO_, Init(), p_clampedComputed, p_clampedDerivFirstPt, p_clampedDerivLastPt, p_clampedEndptsSet, p_clampedSecondDerivs, p_dataValidated, p_fprimeOfx, p_polyNevError, p_x, and p_y.

◆ ~NumericalApproximation()

Isis::NumericalApproximation::~NumericalApproximation ( )
virtual

Destructor deallocates memory being used.

Explicit operations are required to free resources of the interpolation

References GslDeallocation(), GslInterpType(), and p_itype.

Member Function Documentation

◆ AddCubicHermiteDeriv() [1/3]

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.

Parameters
fprimeOfxDerivative value to be added.
Exceptions
Isis::IException::Programmer"Method only used for cspline-Hermite interpolation type"

References _FILEINFO_, CubicHermite, Name(), p_fprimeOfx, p_itype, Isis::IException::Programmer, and ReportException().

◆ AddCubicHermiteDeriv() [2/3]

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.

Parameters
fprimeOfxVector of derivative values to be added.
Exceptions
Isis::IException::Programmer"Method only used for cspline-Hermite interpolation type"

References _FILEINFO_, CubicHermite, Name(), p_fprimeOfx, p_itype, Isis::IException::Programmer, and ReportException().

◆ AddCubicHermiteDeriv() [3/3]

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.

Parameters
nNumber of derivative values to be added.
fprimeOfxArray of derivative values to be added.
Exceptions
Isis::IException::Programmer"Method only used for cspline-Hermite interpolation type"

References _FILEINFO_, CubicHermite, Name(), p_fprimeOfx, p_itype, Isis::IException::Programmer, and ReportException().

◆ AddData() [1/3]

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.

Parameters
xDomain value to add to data set
yRange value corresponding to domain value x

References p_acc, p_clampedComputed, p_clampedDerivFirstPt, p_clampedDerivLastPt, p_clampedEndptsSet, p_clampedSecondDerivs, p_dataValidated, p_interp, p_polyNevError, p_x, and p_y.

Referenced by Isis::AtmosModel::GenerateAhTable(), Isis::AtmosModel::GenerateHahgTables(), Isis::LunarLambertEmpirical::LunarLambertEmpirical(), Isis::MinnaertEmpirical::MinnaertEmpirical(), NumericalApproximation(), and NumericalApproximation().

◆ AddData() [2/3]

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.

Parameters
xVector of domain values to add to data set
yVector of range values corresponding to each domain value in x
Exceptions
Isis::IException::Programmer"Invalid arguments. The number of x-values does not equal the number of y-values."

References _FILEINFO_, p_acc, p_clampedComputed, p_clampedDerivFirstPt, p_clampedDerivLastPt, p_clampedEndptsSet, p_clampedSecondDerivs, p_dataValidated, p_interp, p_polyNevError, p_x, p_y, Isis::IException::Programmer, and ReportException().

◆ AddData() [3/3]

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.

Parameters
nNumber of data points to be used.
xArray of domain values to add to data set
yArray of range values corresponding to each domain value in x

References p_acc, p_clampedComputed, p_clampedDerivFirstPt, p_clampedDerivLastPt, p_clampedEndptsSet, p_clampedSecondDerivs, p_dataValidated, p_interp, p_polyNevError, p_x, and p_y.

◆ BackwardFirstDifference()

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:

  • 2-point backward difference.

    \[ f\prime(a) \approx \frac{1}{h}[f(a) - f(a - h)] \]

  • 3-point backward difference.

    \[ f\prime(a) \approx \frac{1}{2h}[3f(a) - 4f(a - h) + f(a - 2h)] \]

Parameters
aDomain value at which first deriviative is evaluated.
nThe number of points used in the formula
hDistance between nearest points in the formula
Returns
double First derivative approximation for the given domain value
Exceptions
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"Unable to calculate backward first difference for given (a,n,h)"

References _FILEINFO_, DomainMaximum(), DomainMinimum(), Evaluate(), InsideDomain(), p_dataValidated, Isis::IException::Programmer, ReportException(), and ValidateDataSet().

◆ BackwardSecondDifference()

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:

  • 3-point backward second difference.

    \[ f\prime\prime(a) \approx \frac{1}{h^2}[f(a) - 2f(a - h) + f(a - 2h)] \]

Parameters
aDomain value at which second deriviative is evaluated.
nThe number of points used in the formula.
hDistance between nearest points in the formula.
Returns
double Second derivative approximation for the given domain value
Exceptions
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_, DomainMaximum(), DomainMinimum(), Evaluate(), InsideDomain(), p_dataValidated, Isis::IException::Programmer, ReportException(), and ValidateDataSet().

◆ BoolesRule()

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.

Parameters
aLower bound of interval.
bUpper bound of interval.
Returns
double Approximate integral of data set function from a to b.
Exceptions
Isis::IException::Programmer"Unable to calculate the integral on the interval (a,b) using Boole's rule"

References _FILEINFO_, and EvaluateForIntegration().

◆ CenterFirstDifference()

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:

  • 3-point center difference.

    \[ f\prime(a) \approx \frac{1}{2h}[f(a + h) - f(a - h)] \]

  • 5-point center difference.

    \[ f\prime(a) \approx \frac{1}{12h}[-f(a + 2h) + 8f(a +h) - 8f(a - h) + f(a - 2h)] \]

Parameters
aDomain value at which first deriviative is evaluated.
nThe number of points used in the formula.
hDistance between nearest points in the formula.
Returns
double First derivative approximation for the given domain value
Exceptions
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_, DomainMaximum(), DomainMinimum(), Evaluate(), InsideDomain(), p_dataValidated, Isis::IException::Programmer, ReportException(), and ValidateDataSet().

◆ CenterSecondDifference()

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:

  • 3-point center second difference.

    \[ f\prime\prime(a) \approx \frac{1}{h^2}[f(a + h) - 2f(a) + f(a - h)] \]

  • 5-point center second difference.

    \[ f\prime\prime(a) \approx \frac{1}{12h^2}[-f(a + 2h) + 16f(a +h) - 30f(a) + 16f(a - h) - f(a - 2h)] \]

Parameters
aDomain value at which second deriviative is evaluated.
nThe number of points used in the formula.
hDistance between nearest points in the formula.
Returns
double Second derivative approximation for the given domain value
Exceptions
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_, DomainMaximum(), DomainMinimum(), Evaluate(), InsideDomain(), p_dataValidated, Isis::IException::Programmer, ReportException(), and ValidateDataSet().

◆ ComputeCubicClamped()

void Isis::NumericalApproximation::ComputeCubicClamped ( )
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.

Exceptions
Issis::IException::Programmer"Must use SetCubicClampedEndptDeriv() before computing cubic spline with clamped boundary"
Issis::IException::Programmer"Unable to compute cubic clamped interpolation"
See also
SetCubicClampedEndptDeriv()
EvaluateCubicClamped()

References _FILEINFO_, p_clampedComputed, p_clampedDerivFirstPt, p_clampedDerivLastPt, p_clampedEndptsSet, p_clampedSecondDerivs, p_dataValidated, p_x, p_y, Isis::IException::Programmer, ReportException(), Size(), and ValidateDataSet().

Referenced by CubicClampedSecondDerivatives(), and Evaluate().

◆ ComputeGsl()

void Isis::NumericalApproximation::ComputeGsl ( )
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.

Exceptions
Isis::IException::Programmer"Unable to compute GSL interpolation"

References _FILEINFO_, GslAllocation(), GslComputed(), GslIntegrityCheck(), p_dataValidated, p_interp, p_x, p_y, Size(), and ValidateDataSet().

Referenced by DomainMaximum(), DomainMinimum(), Evaluate(), GslFirstDerivative(), GslIntegral(), and GslSecondDerivative().

◆ Contains()

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.

Parameters
xValue to search for in the data set.
Returns
bool Whether the passed value is contained in the x-values of the data set.

References p_x.

◆ CubicClampedSecondDerivatives()

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.

Returns
vector <double> List of second derivatives for each p_x value in the data set.
Exceptions
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."
See also
ComputeCubicClamped()

References _FILEINFO_, ComputeCubicClamped(), CubicClamped, Name(), p_clampedComputed, p_clampedSecondDerivs, p_itype, Isis::IException::Programmer, and ReportException().

◆ DomainMaximum()

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.

Returns
double Maximum domain value
Exceptions
Isis::IException::Programmer"Unable to calculate the domain minimum for the data set."

References _FILEINFO_, ComputeGsl(), GslComputed(), GslInterpType(), p_dataValidated, p_interp, p_itype, p_x, and ValidateDataSet().

Referenced by BackwardFirstDifference(), BackwardSecondDifference(), CenterFirstDifference(), CenterSecondDifference(), EvaluateForIntegration(), ForwardFirstDifference(), ForwardSecondDifference(), GslFirstDerivative(), GslIntegral(), GslSecondDerivative(), InsideDomain(), and ValueToExtrapolate().

◆ DomainMinimum()

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.

Returns
double Minimum domain value
Exceptions
Isis::IException::Programmer"Unable to calculate the domain maximum for the data set."

References _FILEINFO_, ComputeGsl(), GslComputed(), GslInterpType(), p_dataValidated, p_interp, p_itype, p_x, and ValidateDataSet().

Referenced by BackwardFirstDifference(), BackwardSecondDifference(), CenterFirstDifference(), CenterSecondDifference(), EvaluateForIntegration(), FindIntervalLowerIndex(), ForwardFirstDifference(), ForwardSecondDifference(), GslFirstDerivative(), GslIntegral(), GslSecondDerivative(), InsideDomain(), and ValueToExtrapolate().

◆ Evaluate() [1/2]

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.

Parameters
aDomain value from which to interpolate a corresponding range value
etypeNumericalApproximation::ExtrapType enum value indicates how to evaluate if a falls outside of the domain. (Default: ThrowError)
Returns
double Range value interpolated or extrapolated for the given domain value
Exceptions
Isis::IException::Programmer"Unable to evaluate the function at the point a"

References _FILEINFO_, ComputeCubicClamped(), ComputeGsl(), CubicClamped, CubicHermite, CubicNeighborhood, EvaluateCubicClamped(), EvaluateCubicHermite(), EvaluateCubicNeighborhood(), EvaluatePolynomialNeville(), GslComputed(), GslIntegrityCheck(), InsideDomain(), p_acc, p_clampedComputed, p_interp, p_itype, p_polyNevError, PolynomialNeville, and ValueToExtrapolate().

Referenced by Isis::HapkeAtm1::AtmosModelAlgorithm(), Isis::HapkeAtm2::AtmosModelAlgorithm(), BackwardFirstDifference(), BackwardSecondDifference(), CenterFirstDifference(), CenterSecondDifference(), Evaluate(), EvaluateForIntegration(), ForwardFirstDifference(), ForwardSecondDifference(), Isis::LunarLambertEmpirical::PhotoModelAlgorithm(), Isis::MinnaertEmpirical::PhotoModelAlgorithm(), and RefineExtendedTrap().

◆ Evaluate() [2/2]

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.

Parameters
aVector of domain values from which to interpolate a vector of corresponding range values
etypeNumericalApproximation::ExtrapType enum value indicates how to evaluate any value of a that falls outside of the domain. (Default: ThrowError)
Returns
vector <double> Set of interpolated range values corresponding to the elements of the vector of domain values given
Exceptions
Isis::IException::Programmer"Unable to evaluate the function at the vector of points, a."

References _FILEINFO_, CubicNeighborhood, Evaluate(), EvaluateCubicNeighborhood(), EvaluatePolynomialNeville(), InsideDomain(), p_itype, p_polyNevError, PolynomialNeville, and ValueToExtrapolate().

◆ EvaluateCubicClamped()

double Isis::NumericalApproximation::EvaluateCubicClamped ( const double a)
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.

Parameters
aDomain value from which to interpolate a corresponding range value
Returns
double Value returned from interpolation or extrapolation.
See also
Evaluate()
ComputeCubicClamped()
SetCubicClampedEndptDeriv()

References p_clampedSecondDerivs, p_x, p_y, and Size().

Referenced by Evaluate().

◆ EvaluateCubicHermite()

double Isis::NumericalApproximation::EvaluateCubicHermite ( const double a)
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.

Parameters
aDomain value from which to interpolate a corresponding range value
Returns
double Value returned from interpolation or extrapolation
Exceptions
Isis::IException::User"Invalid arguments. The size of the first derivative vector does not match the number of (x,y) data points."
See also
http://mathworld.wolfram.com/HermitesInterpolatingPolynomial.html
http://en.wikipedia.org/wiki/Lagrange_Polynomial
Evaluate()
AddCubicHermiteDeriv()

References _FILEINFO_, FindIntervalLowerIndex(), p_fprimeOfx, p_x, p_y, ReportException(), Size(), and Isis::IException::User.

Referenced by Evaluate().

◆ EvaluateCubicHermiteFirstDeriv()

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_, CubicHermite, FindIntervalLowerIndex(), Name(), p_fprimeOfx, p_itype, p_x, p_y, ReportException(), Size(), and Isis::IException::User.

◆ EvaluateCubicHermiteSecDeriv()

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_, CubicHermite, FindIntervalLowerIndex(), Name(), p_fprimeOfx, p_itype, p_x, p_y, ReportException(), Size(), and Isis::IException::User.

◆ EvaluateCubicNeighborhood() [1/2]

double Isis::NumericalApproximation::EvaluateCubicNeighborhood ( const double a)
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.

Parameters
aDomain value from which to interpolate a corresponding range value.
Returns
double Result of interpolated value of f(a).
Exceptions
Isis::IException::Programmer"Unable to evaluate cubic neighborhood interpolation at a"
See also
Evaluate()

References _FILEINFO_, CubicNatural, NearestEndpoint, p_x, p_y, and Size().

Referenced by Evaluate(), and Evaluate().

◆ EvaluateCubicNeighborhood() [2/2]

vector< double > Isis::NumericalApproximation::EvaluateCubicNeighborhood ( const vector< double > & a,
const ExtrapType & etype )
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.

Parameters
aVector of domain values from which to interpolate a vector of corresponding range values
etypeNumericalApproximation::ExtrapType enum value indicates how to evaluate if a falls outside of the domain. (Default: ThrowError)
Returns
vector <double> Result of interpolated value of f at each value of a.
Exceptions
Isis::IException::Programmer"Unable to evaluate cubic neighborhood interpolation at the values in the vector, a"
See also
Evaluate()

References _FILEINFO_, CubicNatural, InsideDomain(), NearestEndpoint, p_x, p_y, Size(), and ValueToExtrapolate().

◆ EvaluateForIntegration()

vector< double > Isis::NumericalApproximation::EvaluateForIntegration ( const double a,
const double b,
const unsigned int n )
protected

Evaluates data set in order to have enough data points to approximate the function to be integrated.

Parameters
aLower bound of interval.
bUpper bound of interval.
nNumber of points used in Newton-Cotes formula.
Returns
vector <double> Array of y values to be used in numerical integration
Exceptions
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_, DomainMaximum(), DomainMinimum(), Evaluate(), InsideDomain(), Isis::IException::Programmer, ReportException(), and Size().

Referenced by BoolesRule(), Simpsons3PointRule(), Simpsons4PointRule(), and TrapezoidalRule().

◆ EvaluatePolynomialNeville()

double Isis::NumericalApproximation::EvaluatePolynomialNeville ( const double a)
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.

Parameters
aDomain value from which to interpolate a corresponding range value
Returns
double Value returned from interpolation or extrapolation
See also
http://mathworld.wolfram.com/NevillesAlgorithm.html
Evaluate()
PolynomialNevilleErrorEstimate()

References p_polyNevError, p_x, p_y, and Size().

Referenced by Evaluate(), and Evaluate().

◆ FindIntervalLowerIndex()

int Isis::NumericalApproximation::FindIntervalLowerIndex ( const double a)
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.

Parameters
aDomain value around which the interval lies
Returns
int Index of the x-value that is the lower endpoint of the interval of data points that surrounds a. Returns 0 if a is below domain min and Size()-2 if a is above domain max.

References DomainMinimum(), InsideDomain(), p_x, and Size().

Referenced by EvaluateCubicHermite(), EvaluateCubicHermiteFirstDeriv(), and EvaluateCubicHermiteSecDeriv().

◆ ForwardFirstDifference()

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:

  • 2-point forward difference.

    \[ f\prime(a) \approx \frac{1}{h}[f(a + h) - f(a)] \]

  • 3-point forward difference.

    \[ f\prime(a) \approx \frac{1}{2h}[-f(a + 2h) + 4f(a + h) - 3f(a)] \]

Parameters
aDomain value at which first deriviative is evaluated.
nThe number of points used in the formula.
hDistance between nearest points in the formula.
Returns
double First derivative approximation for the given domain value
Exceptions
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_, DomainMaximum(), DomainMinimum(), Evaluate(), InsideDomain(), p_dataValidated, Isis::IException::Programmer, ReportException(), and ValidateDataSet().

◆ ForwardSecondDifference()

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:

  • 3-point forward second difference.

    \[ f\prime\prime(a) \approx \frac{1}{h^2}[f(a + 2h) - 2f(a + h) + f(a)] \]

Parameters
aDomain value at which second deriviative is evaluated.
nThe number of points used in the formula.
hDistance between nearest points in the formula.
Returns
double Second derivative approximation for the given domain value
Exceptions
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_, DomainMaximum(), DomainMinimum(), Evaluate(), InsideDomain(), p_dataValidated, Isis::IException::Programmer, ReportException(), and ValidateDataSet().

◆ GslAllocation()

void Isis::NumericalApproximation::GslAllocation ( unsigned int npoints)
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.

Parameters
npointsNumber of points to allocate for the GSL interpolator

References GslDeallocation(), GslFunctor(), p_acc, p_interp, and p_itype.

Referenced by ComputeGsl().

◆ GslComputed()

bool Isis::NumericalApproximation::GslComputed ( ) const
protected

Returns whether a GSL interpolation computation of the data set has been performed.

This method is only applicable to GSL interpolation types.

Returns
bool True if GSL interpolation has been computed.
Exceptions
Isis::IException::Programmer"Method only valid for GSL interpolation types"
See also
ComputeGsl()

References _FILEINFO_, GslInterpType(), Name(), p_acc, p_interp, p_itype, and Isis::IException::Programmer.

Referenced by ComputeGsl(), DomainMaximum(), DomainMinimum(), Evaluate(), GslFirstDerivative(), GslIntegral(), and GslSecondDerivative().

◆ GslDeallocation()

void Isis::NumericalApproximation::GslDeallocation ( )
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.

References p_acc, and p_interp.

Referenced by GslAllocation(), operator=(), Reset(), and ~NumericalApproximation().

◆ GslFirstDerivative()

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.

Parameters
aDomain value at which first deriviative is evaluated.
Returns
double First derivative approximation for the given domain value, if valid.
Exceptions
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"
Isis::IException::Programmer"Unable to compute the first derivative at a. GSL integrity check failed"

References _FILEINFO_, ComputeGsl(), DomainMaximum(), DomainMinimum(), GslComputed(), GslIntegrityCheck(), GslInterpType(), InsideDomain(), Name(), p_acc, p_interp, p_itype, Isis::IException::Programmer, and ReportException().

◆ GslFunctor()

NumericalApproximation::InterpFunctor Isis::NumericalApproximation::GslFunctor ( NumericalApproximation::InterpType itype) const
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.

Parameters
itypeType of GSL interpolator to find
Returns
NumericalApproximation::InterpFunctor Pointer to the GSL spline interpolator construct.
Exceptions
Isis::IException::Programmer"Invalid argument. Unable to find GSL interpolator"

References _FILEINFO_, p_interpFunctors, Isis::IException::Programmer, and ReportException().

Referenced by GslAllocation(), MinPoints(), Name(), and SetInterpType().

◆ GslIntegral()

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.

Parameters
aLower endpoint at which integral is evaluated.
bUpper endpoint at which integral is evaluated.
Returns
double Integral approximation for the given interval, if valid.
Exceptions
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"
Isis::IException::Programmer"Unable to compute the integral on the interval (a,b). GSL integrity check failed."

References _FILEINFO_, ComputeGsl(), DomainMaximum(), DomainMinimum(), GslComputed(), GslIntegrityCheck(), GslInterpType(), InsideDomain(), Name(), p_acc, p_interp, p_itype, Isis::IException::Programmer, and ReportException().

◆ GslIntegrityCheck()

void Isis::NumericalApproximation::GslIntegrityCheck ( int gsl_status,
const char * src,
int line )
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.

Parameters
gsl_statusReturn status of a GSL function call
srcName of the calling source invoking the check. This allows more precise determination of the error.
lineLine of the calling source that invoked the check
Exceptions
Isis::IException::Programmer"GSL error occured"

References Isis::IException::Programmer, and ReportException().

Referenced by ComputeGsl(), Evaluate(), GslFirstDerivative(), GslIntegral(), and GslSecondDerivative().

◆ GslInterpType()

bool Isis::NumericalApproximation::GslInterpType ( NumericalApproximation::InterpType itype) const
protected

Returns whether an interpolation type is adapted from the GSL library.

GSL interpolation types include the following:

  • Linear = 0
  • Polynomial = 1
  • CubicNatural = 3
  • CubicNatPeriodic = 5
  • Akima = 7
  • AkimaPeriodic = 8
Parameters
itypeInterpolation type to be compared to GSL list
Returns
bool True if interpolation type is GSL

References Akima, AkimaPeriodic, CubicNatPeriodic, CubicNatural, Linear, and Polynomial.

Referenced by DomainMaximum(), DomainMinimum(), GslComputed(), GslFirstDerivative(), GslIntegral(), GslSecondDerivative(), MinPoints(), operator=(), RefineExtendedTrap(), Reset(), SetInterpType(), ValueToExtrapolate(), and ~NumericalApproximation().

◆ GslSecondDerivative()

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.

Parameters
aDomain value at which second deriviative is evaluated.
Returns
double Second derivative approximation for the given domain value, if valid.
Exceptions
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_, ComputeGsl(), DomainMaximum(), DomainMinimum(), GslComputed(), GslIntegrityCheck(), GslInterpType(), InsideDomain(), Name(), p_acc, p_interp, p_itype, Isis::IException::Programmer, and ReportException().

◆ Init()

void Isis::NumericalApproximation::Init ( NumericalApproximation::InterpType itype)
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.

Parameters
itypeNumericalApproximation::InterpType enum value to be assigned to this object
Exceptions
Isis::IException::Programmer"Unable to initialize NumericalApproximation object"

References _FILEINFO_, Akima, AkimaPeriodic, CubicNatPeriodic, CubicNatural, Linear, p_acc, p_interp, p_interpFunctors, Polynomial, and SetInterpType().

Referenced by NumericalApproximation(), NumericalApproximation(), NumericalApproximation(), and NumericalApproximation().

◆ InsideDomain()

bool Isis::NumericalApproximation::InsideDomain ( const double a)
protected

Returns whether the passed value is greater than or equal to the domain minimum and less than or equal to the domain maximum.

Parameters
aValue to be verified as valid domain value.
Returns
bool True if passed parameter is within the domain.
Exceptions
Isis::IException::Programmer"Unable to compute domain boundaries"

References _FILEINFO_, DomainMaximum(), and DomainMinimum().

Referenced by BackwardFirstDifference(), BackwardSecondDifference(), CenterFirstDifference(), CenterSecondDifference(), Evaluate(), Evaluate(), EvaluateCubicNeighborhood(), EvaluateForIntegration(), FindIntervalLowerIndex(), ForwardFirstDifference(), ForwardSecondDifference(), GslFirstDerivative(), GslIntegral(), and GslSecondDerivative().

◆ InterpolationType()

InterpType Isis::NumericalApproximation::InterpolationType ( )
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:

  • Linear = 0
  • Polynomial = 1
  • PolynomialNeville = 2
  • CubicNatural = 3
  • CubicClamped = 4
  • CubicNatPeriodic = 5
  • CubicNeighborhood = 6
  • CubicHermite = 7
  • Akima = 8
  • AkimaPeriodic = 9
Returns
NumericalApproximation::InterpType Currently assigned interpolation type

References p_itype.

◆ MinPoints() [1/2]

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.

Returns
int Minimum number of data points required for the interpolation function
Exceptions
Isis::IException::Programmer"Unable to calculate minimum points."

References _FILEINFO_, CubicClamped, CubicHermite, CubicNeighborhood, GslFunctor(), p_itype, and PolynomialNeville.

Referenced by ValidateDataSet().

◆ MinPoints() [2/2]

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.

Parameters
itypeType of interpolation function for which to return minimum points
Returns
int Minimum number of data points required for the interpolation function
Exceptions
Isis::IException::Programmer"Invalid argument. Unknown interpolation type"
Isis::IException::Programmer"Unable to calculate minimum points"

References _FILEINFO_, CubicClamped, CubicHermite, CubicNeighborhood, GslInterpType(), PolynomialNeville, and Isis::IException::Programmer.

◆ Name()

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.

Returns
std::string Name of the interpolation function used/to be used.
Exceptions
Isis::IException::Programmer"Unable to retrieve numerical approximation name",

References _FILEINFO_, CubicClamped, CubicHermite, CubicNatural, CubicNeighborhood, GslFunctor(), p_itype, and PolynomialNeville.

Referenced by AddCubicHermiteDeriv(), AddCubicHermiteDeriv(), AddCubicHermiteDeriv(), CubicClampedSecondDerivatives(), EvaluateCubicHermiteFirstDeriv(), EvaluateCubicHermiteSecDeriv(), GslComputed(), GslFirstDerivative(), GslIntegral(), GslSecondDerivative(), PolynomialNevilleErrorEstimate(), SetCubicClampedEndptDeriv(), ValidateDataSet(), and ValueToExtrapolate().

◆ operator=()

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.

Parameters
oldObjectNumericalApproximation object to copy data from for initialization
Exceptions
Isis::IException::Programmer"Unable to copy NumericalApproximation object"

References _FILEINFO_, GslDeallocation(), GslInterpType(), p_clampedDerivFirstPt, p_clampedDerivLastPt, p_clampedSecondDerivs, p_dataValidated, p_fprimeOfx, p_itype, p_polyNevError, p_x, p_y, and SetInterpType().

◆ PolynomialNevilleErrorEstimate()

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.

Returns
vector <double> Estimation of difference between actual value of f(a) and the polynomial Neville interpolated return value of Evaluate(), where f is the function that defines the given data set.
Exceptions
Isis::IException::Programmer"Method only used for polynomial-Neville's interpolation type"
Isis::IException::Programmer"Error not calculated"

References _FILEINFO_, Name(), p_itype, p_polyNevError, PolynomialNeville, Isis::IException::Programmer, and ReportException().

◆ RefineExtendedTrap()

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:

double result;
for(int j = 1; j <= M+1; j++){
result = RefineExtendedTrap(a,b,result,j);
}
double RefineExtendedTrap(double a, double b, double s, unsigned int n)
Calculates refinements extended trapezoidal rule to approximate the integral of the interpolated data...
Definition NumericalApproximation.cpp:2096

Note: Although this method may be used to approximate an integral, as described above, it is most often used by RombergsMethod() to integrate.

Parameters
aLower limit of integration
bUpper limit of integration
sPrevious value of refinement
nNumber of partitions to use when integrating
Returns
double Integral (refined) approximation of the function on the interval (a, b)
Exceptions
Isis::IException::Programmer"Unable to calculate the integral on the interval (a,b) using the extended trapezoidal rule"
See also
RombergsMethod()

References _FILEINFO_, CubicNeighborhood, Evaluate(), Extrapolate, GslInterpType(), NearestEndpoint, and p_itype.

Referenced by RombergsMethod().

◆ ReportException()

void Isis::NumericalApproximation::ReportException ( IException::ErrorType type,
const string & methodName,
const string & message,
const char * filesrc,
int lineno ) const
protected

Generalized error report generator.

This method is used throughout this class to standardize error reporting as a convenience to its implementor.

Parameters
typeType of Isis::IException
methodNameName of method where exception originated
messageError context string provided by the caller
filesrcName of the file the error occured in.
linenoLine number of the calling source where the error occured.
Exceptions
Isis::IException::errTypeequal to type and error message equal to (methodName + " - " + message)

Referenced by AddCubicHermiteDeriv(), AddCubicHermiteDeriv(), AddCubicHermiteDeriv(), AddData(), BackwardFirstDifference(), BackwardSecondDifference(), CenterFirstDifference(), CenterSecondDifference(), ComputeCubicClamped(), CubicClampedSecondDerivatives(), EvaluateCubicHermite(), EvaluateCubicHermiteFirstDeriv(), EvaluateCubicHermiteSecDeriv(), EvaluateForIntegration(), ForwardFirstDifference(), ForwardSecondDifference(), GslFirstDerivative(), GslFunctor(), GslIntegral(), GslIntegrityCheck(), GslSecondDerivative(), PolynomialNevilleErrorEstimate(), SetCubicClampedEndptDeriv(), SetInterpType(), ValidateDataSet(), and ValueToExtrapolate().

◆ Reset() [1/2]

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.

References GslDeallocation(), GslInterpType(), p_clampedComputed, p_clampedDerivFirstPt, p_clampedDerivLastPt, p_clampedEndptsSet, p_clampedSecondDerivs, p_dataValidated, p_fprimeOfx, p_itype, p_polyNevError, p_x, and p_y.

Referenced by Isis::AtmosModel::GenerateAhTable(), Isis::AtmosModel::GenerateHahgTables(), Isis::LunarLambertEmpirical::LunarLambertEmpirical(), Isis::MinnaertEmpirical::MinnaertEmpirical(), Reset(), Isis::LunarLambertEmpirical::~LunarLambertEmpirical(), and Isis::MinnaertEmpirical::~MinnaertEmpirical().

◆ Reset() [2/2]

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.

Parameters
itypeNumericalApproximation::InterpType enum value to be assigned to this object
Exceptions
Isis::IException::Programmer"Unable to reset interpolation type"

References _FILEINFO_, Reset(), and SetInterpType().

◆ RombergsMethod()

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.

Parameters
aLower limit of integration
bUpper limit of integration
Returns
double Integral approximation of the function on the interval (a, b)
Exceptions
Isis::IException::Programmer"Failed to converge."
Isis::IException::Programmer"Unable to calculate the integral on the interval (a,b) using Romberg's method"
See also
http://mathworld.wolfram.com/RombergIntegration.html
RefineExtendedTrap()

References _FILEINFO_, Extrapolate, PolynomialNeville, Isis::IException::Programmer, and RefineExtendedTrap().

◆ SetCubicClampedEndptDeriv()

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.

Parameters
yp1First derivative of the function evaluated at the domain minimum.
ypnFirst derivative of the function evaluated at the domain maximum.
Exceptions
Isis::IException::Programmer"Method only used for cspline-clamped interpolation type"

References _FILEINFO_, CubicClamped, Name(), p_clampedDerivFirstPt, p_clampedDerivLastPt, p_clampedEndptsSet, p_itype, Isis::IException::Programmer, and ReportException().

Referenced by Isis::AtmosModel::GenerateAhTable(), Isis::AtmosModel::GenerateHahgTables(), Isis::LunarLambertEmpirical::LunarLambertEmpirical(), and Isis::MinnaertEmpirical::MinnaertEmpirical().

◆ SetInterpType()

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.

Parameters
itypeNumericalApproximation::InterpType enum value to be assigned to this object
Exceptions
Isis::IException::Programmer"Invalid argument. Unknown interpolation type"
Isis::IException::Programmer"Unable to set interpolation type"
See also
Reset()

References _FILEINFO_, GslFunctor(), GslInterpType(), p_clampedComputed, p_clampedDerivFirstPt, p_clampedDerivLastPt, p_clampedEndptsSet, p_clampedSecondDerivs, p_dataValidated, p_fprimeOfx, p_itype, p_polyNevError, Isis::IException::Programmer, and ReportException().

Referenced by Isis::AtmosModel::GenerateAhTable(), Isis::AtmosModel::GenerateHahgTables(), Init(), Isis::LunarLambertEmpirical::LunarLambertEmpirical(), Isis::MinnaertEmpirical::MinnaertEmpirical(), operator=(), and Reset().

◆ Simpsons3PointRule()

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.

Parameters
aLower bound of interval.
bUpper bound of interval.
Returns
double Approximate integral of data set function from a to b.
Exceptions
Isis::IException::Programmer"Unable to calculate the integral on the interval (a,b) using Simpson's 3 point rule"

References _FILEINFO_, and EvaluateForIntegration().

◆ Simpsons4PointRule()

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.

Parameters
aLower bound of interval.
bUpper bound of interval.
Returns
double Approximate integral of data set function from a to b.
Exceptions
Isis::IException::Programmer"Unable to calculate the integral on the interval (a,b) using Simpson's 4 point rule"

References _FILEINFO_, and EvaluateForIntegration().

◆ Size()

unsigned int Isis::NumericalApproximation::Size ( )
inline

◆ TrapezoidalRule()

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.

Parameters
aLower bound of interval.
bUpper bound of interval.
Returns
double Approximate integral of data set function from a to b.
Exceptions
Isis::IException::Programmer"Unable to calculate the integral on the interval (a,b) using the trapezoidal rule"

References _FILEINFO_, and EvaluateForIntegration().

◆ ValidateDataSet()

void Isis::NumericalApproximation::ValidateDataSet ( )
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}) \).

Exceptions
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_, CubicNatPeriodic, MinPoints(), Name(), p_dataValidated, p_itype, p_x, p_y, PolynomialNeville, Isis::IException::Programmer, ReportException(), and Size().

Referenced by BackwardFirstDifference(), BackwardSecondDifference(), CenterFirstDifference(), CenterSecondDifference(), ComputeCubicClamped(), ComputeGsl(), DomainMaximum(), DomainMinimum(), ForwardFirstDifference(), ForwardSecondDifference(), NumericalApproximation(), and NumericalApproximation().

◆ ValueToExtrapolate()

double Isis::NumericalApproximation::ValueToExtrapolate ( const double a,
const ExtrapType & etype )
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.

Parameters
aValue passed into Evaluate() that falls outside of domain.
etypeNumericalApproximation::ExtrapType enum value indicates how to evaluate if a falls outside of the domain. (Default: ThrowError)
Returns
double Value returned to Evaluate() to be extrapolated.
Exceptions
Isis::IException::Programmer"Invalid argument. Value entered is outside of domain @throws Isis::IException::Programmer "Invalid argument. Cannot extrapolate for interpolation type" (GSL or CubicNeighborhood)
See also
Evaluate()

References _FILEINFO_, CubicNeighborhood, DomainMaximum(), DomainMinimum(), GslInterpType(), Name(), NearestEndpoint, p_itype, Isis::IException::Programmer, ReportException(), and ThrowError.

Referenced by Evaluate(), Evaluate(), and EvaluateCubicNeighborhood().

Member Data Documentation

◆ p_acc

gsl_interp_accel* Isis::NumericalApproximation::p_acc
protected

◆ p_clampedComputed

bool Isis::NumericalApproximation::p_clampedComputed
protected

◆ p_clampedDerivFirstPt

double Isis::NumericalApproximation::p_clampedDerivFirstPt
protected

First derivative of first x-value, p_x[0]. This is only used for the CubicClamped interpolation type.

Referenced by AddData(), AddData(), AddData(), ComputeCubicClamped(), NumericalApproximation(), operator=(), Reset(), SetCubicClampedEndptDeriv(), and SetInterpType().

◆ p_clampedDerivLastPt

double Isis::NumericalApproximation::p_clampedDerivLastPt
protected

First derivative of last x-value, p_x[n-1]. This is only used for the CubicClamped interpolation type.

Referenced by AddData(), AddData(), AddData(), ComputeCubicClamped(), NumericalApproximation(), operator=(), Reset(), SetCubicClampedEndptDeriv(), and SetInterpType().

◆ p_clampedEndptsSet

bool Isis::NumericalApproximation::p_clampedEndptsSet
protected

Flag variable to determine whether SetCubicClampedEndptDeriv() has been called after all data was added for CubicClamped interpolation.

Referenced by AddData(), AddData(), AddData(), ComputeCubicClamped(), NumericalApproximation(), Reset(), SetCubicClampedEndptDeriv(), and SetInterpType().

◆ p_clampedSecondDerivs

vector<double> Isis::NumericalApproximation::p_clampedSecondDerivs
protected

List of second derivatives evaluated at p_x values. This is only used for the CubicClamped interpolation type.

Referenced by AddData(), AddData(), AddData(), ComputeCubicClamped(), CubicClampedSecondDerivatives(), EvaluateCubicClamped(), NumericalApproximation(), operator=(), Reset(), and SetInterpType().

◆ p_dataValidated

◆ p_fprimeOfx

vector<double> Isis::NumericalApproximation::p_fprimeOfx
protected

◆ p_interp

gsl_spline* Isis::NumericalApproximation::p_interp
protected

◆ p_interpFunctors

NumericalApproximation::FunctorList Isis::NumericalApproximation::p_interpFunctors
staticprotected

Maintains list of interpolator options.

Referenced by GslFunctor(), and Init().

◆ p_itype

◆ p_polyNevError

vector<double> Isis::NumericalApproximation::p_polyNevError
protected

Estimate of error for interpolation evaluated at x. This is only used for the PolynomialNeville interpolation type. 91 taken from AtmosModel.

Referenced by AddData(), AddData(), AddData(), Evaluate(), Evaluate(), EvaluatePolynomialNeville(), NumericalApproximation(), operator=(), PolynomialNevilleErrorEstimate(), Reset(), and SetInterpType().

◆ p_x

◆ p_y


The documentation for this class was generated from the following files: