PCHIP
Piecewise Cubic Hermite Interpolant Package
PCHIP
is a FORTRAN77 library which
can construct a piecewise cubic Hermite interpolant to data, and
carry out various related operations, by Fred Fritsch.
Licensing:
The computer code and data files made available on this
web page are distributed under
the GNU LGPL license.
Languages:
PCHIP is available in
a FORTRAN77 version and
a FORTRAN90 version.
Related Data and Programs:
BARYCENTRIC_INTERP_1D,
a FORTRAN77 library which
defines and evaluates the barycentric Lagrange polynomial p(x)
which interpolates a set of data, so that p(x(i)) = y(i).
The barycentric approach means that very high degree polynomials can
safely be used.
CHEBYSHEV_INTERP_1D,
a FORTRAN77 library which
determines the combination of Chebyshev polynomials which
interpolates a set of data, so that p(x(i)) = y(i).
INTERPOLATION,
a dataset directory which
contains datasets to be interpolated.
LAGRANGE_INTERP_1D,
a FORTRAN77 library which
defines and evaluates the Lagrange polynomial p(x)
which interpolates a set of data, so that p(x(i)) = y(i).
NEAREST_INTERP_1D,
a FORTRAN77 library which
interpolates a set of data using a piecewise constant interpolant
defined by the nearest neighbor criterion,
creating graphics files for processing by gnuplot.
NMS,
a FORTRAN77 library which
includes a wide variety of numerical software, including
solvers for linear systems of equations,
a Piecewise Cubic Hermite Interpolation Package (PCHIP),
numerical quadrature, linear least squares data fitting,
the solution of nonlinear equations,
ordinary differential equations (ODE's),
optimization and nonlinear least squares, simulation and random numbers,
trigonometric approximation and Fast Fourier Transforms (FFT).
PWL_INTERP_1D,
a FORTRAN77 library which
interpolates a set of data using a piecewise linear function,
creating graphics files for processing by gnuplot.
TEST_INTERP,
a FORTRAN77 library which
defines a number of test problems for interpolation,
provided as a set of (x,y) data.
TEST_INTERP_1D,
a FORTRAN77 library which
defines test problems for interpolation of data y(x),
depending on a 1D argument.
Reference:
-
Carl deBoor,
A Practical Guide to Splines,
Springer, 2001,
ISBN: 0387953663,
LC: QA1.A647.v27.
-
Carl deBoor,
Package for calculating with B-splines,
SIAM Journal on Numerical Analysis,
Volume 14, Number 3, June 1977, pages 441-472.
-
Fred Fritsch,
Representations for Parametric Cubic Splines,
Computer Aided Geometric Design,
Volume 6, 1989, pages 79-82.
-
Fred Fritsch,
Piecewise Cubic Hermite Interpolation Package, Final Specifications,
Lawrence Livermore National Laboratory,
Computer Documentation UCID-30194, August 1982.
-
Fred Fritsch, Judy Butland,
A Method for Constructing Local Monotone Piecewise
Cubic Interpolants,
SIAM Journal on Scientific and Statistical Computing,
Volume 5, Number 2, June 1984, pages 300-304.
-
Fred Fritsch, Ralph Carlson,
Monotone Piecewise Cubic Interpolation,
SIAM Journal on Numerical Analysis,
Volume 17, Number 2, April 1980, pages 238-246.
-
David Kahaner, Cleve Moler, Steven Nash,
Numerical Methods and Software,
Prentice Hall, 1989,
ISBN: 0-13-627258-4,
LC: TA345.K34.
Source Code:
-
pchip.f, the source code.
-
pchip.sh,
BASH commands to compile the source code.
Examples and Tests:
TEST02 generates some files that can be input to GNUPLOT
to illustrate the problem of interpolating monotonic data with
a piecewise cubic Hermite spline.
TEST03 generates some files that can be input to GNUPLOT
to illustrate the problem of interpolating monotonic data with
a cubic spline.
List of Routines:
-
BVALU evaluates a B-spline or its derivatives.
-
CHFCM checks a single cubic for monotonicity.
-
CHFDV evaluates a cubic polynomial and derivative in Hermite form.
-
CHFEV evaluates a cubic polynomial in Hermite form.
-
CHFIE evaluates the integral of a single cubic for PCHIA.
-
COMP compares actual and expected error flag values.
-
D1MACH returns double precision real machine-dependent constants.
-
DBVALU evaluates a B-spline or its derivatives.
-
DCHFCM checks a single cubic for monotonicity.
-
DCHFDV evaluates a cubic Hermite polynomial and derivative at many points.
-
DCHFEV evaluates a cubic Hermite polynomial at many points.
-
DCHFIE evaluates the integral of a single cubic for DPCHIA.
-
DEVCHK tests the evaluation accuracy of DCHFDV and DCHFEV for DPCHQ1.
-
DEVERK tests error returns from DPCHIP evaluators for DPCHQ1.
-
DEVPCK tests usage of increment argument in DPCHFD and DPCHFE for DPCHQ1.
-
DFDTRU computes exact function values for DEVCHK.
-
DINTRV determines the interval containing, or nearest, a point.
-
DPCHBS: Piecewise cubic Hermite to B-Spline converter.
-
DPCHCE sets boundary conditions for DPCHIC.
-
DPCHCI sets interior derivatives for DPCHIC.
-
DPCHCM checks a cubic Hermite function for monotonicity.
-
DPCHCS adjusts derivative values for DPCHIC.
-
DPCHDF computes divided differences for DPCHCE and DPCHSP.
-
DPCHEV: value and derivative of spline or cubic Hermite at many points.
-
DPCHEZ sets up a spline or cubic Hermite interpolant.
-
DPCHFD evaluates a piecewise cubic Hermite and derivative at many points.
-
DPCHFE evaluates a piecewise cubic Hermite function at many points.
-
DPCHIA evaluates the definite integral of a piecewise cubic Hermite function.
-
DPCHIC sets derivatives for a monotone piecewise cubic Hermite interpolant.
-
DPCHID: integral of a piecewise cubic Hermite function between data points.
-
DPCHIM sets derivatives for a monotone piecewise cubic Hermite interpolant.
-
DPCHKT computes the B-spline knot sequence for DPCHBS.
-
DPCHQ1: Test the PCHIP evaluators DCHFDV, DCHFEV, DPCHFD, DPCHFE.
-
DPCHQ2 tests the PCHIP integrators DPCHIA and DPCHID.
-
DPCHQ3 tests the PCHIP interpolators DPCHIC, DPCHIM, DPCHSP.
-
DPCHQ4 tests the PCHIP monotonicity checker DPCHCM.
-
DPCHQ5 tests the PCH to B-spline conversion routine DPCHBS.
-
DPCHQA: definite integral of spline or piecewise cubic Hermite interpolant.
-
DPCHSP: derivatives for Hermite representation of cubic spline interpolant.
-
DPCHST: carry out a sign test.
-
DPCHSW limits excursion from data for DPCHCS.
-
EVCHK tests evaluation accuracy of CHFDV and DHFEV for PCHQK1.
-
EVERCK tests error returns from PCHIP evaluators for PCHQK1.
-
EVPCCK tests usage of increment argument in PCHFD and PCHFE for PCHQK1.
-
FDTRUE computes exact function values for EVCHCK.
-
FDUMP can be used to print a symbolic dump.
-
GET_UNIT returns a free FORTRAN unit number.
-
I1MACH returns integer machine dependent constants.
-
INTRV seeks the interval containing or nearest to a given point.
-
J4SAVE saves or recalls global variables needed by the error handler.
-
PCHBS: Piecewise cubic Hermite to B-Spline converter.
-
PCHCE sets boundary conditions for PCHIC.
-
PCHCI sets interior derivatives for PCHIC.
-
PCHCM checks a cubic Hermite function for monotonicity.
-
PCHCS adjusts derivative values for PCHIC.
-
PCHDF computes divided differences for PCHCE and PCHSP.
-
PCHEV: value and derivative of spline or cubic Hermite at many points.
-
PCHEZ sets up a spline or cubic Hermite interpolant.
-
PCHFD evaluates a piecewise cubic Hermite and derivative at many points.
-
PCHFE evaluates a piecewise cubic Hermite function at many points.
-
PCHIA evaluates the definite integral of a piecewise cubic Hermite function.
-
PCHIC sets derivatives for a monotone piecewise cubic Hermite interpolant.
-
PCHID: integral of a piecewise cubic Hermite function between data points.
-
PCHIM sets derivatives for a monotone piecewise cubic Hermite interpolant.
-
PCHKT computes the B-spline knot sequence for PCHBS.
-
PCHQA: definite integral of spline or piecewise cubic Hermite interpolant.
-
PCHQK1 tests the PCHIP evaluators CHFDV, CHFEV, PCHFD and PCHFE.
-
PCHQK2 tests the PCHIP integrators PCHIA and PCHID.
-
PCHQK3 tests the PCHIP interpolators PCHIC, PCHIM, PCHSP.
-
PCHQK4 tests the PCHIP monotonicity checker PCHCM.
-
PCHQK5 tests the PCH to B-spline conversion routine PCHBS.
-
PCHSP: derivatives for Hermite representation of cubic spline interpolant.
-
PCHST: carry out a sign test.
-
PCHSW limits excursion from data for PCHCS.
-
R1MACH returns single precision real machine constants.
-
RAND generates a uniformly distributed random number.
-
TIMESTAMP prints out the current YMDHMS date as a timestamp.
-
XERCNT allows user control over handling of errors.
-
XERDMP prints the error tables and clears them.
-
XERHLT halts program execution after a fatal error.
-
XERMSG processes error messages.
-
XERPRN prints error messages from XERMSG.
-
XERSVE records that an error has occurred.
-
XGETF returns the current value of the error control flag.
-
XGETUA returns error unit numbers.
-
XSETF sets the error control flag.
You can go up one level to
the FORTRAN77 source codes.
Last revised on 06 April 2015.