pchip_module Module

A Fortran package for piecewise cubic Hermite interpolation of data.

Author

  • Fritsch, F. N., (LLNL) -- original author
  • Oct 2019 : Jacob Williams, Extensive refactoring and modernization of the SLATEC code.

Uses

  • module~~pchip_module~~UsesGraph module~pchip_module pchip_module iso_fortran_env iso_fortran_env module~pchip_module->iso_fortran_env

Variables

Type Visibility Attributes Name Initial
real(kind=wp), private, parameter :: zero = 0.0_wp
real(kind=wp), private, parameter :: half = 0.5_wp
real(kind=wp), private, parameter :: one = 1.0_wp
real(kind=wp), private, parameter :: two = 2.0_wp
real(kind=wp), private, parameter :: three = 3.0_wp
real(kind=wp), private, parameter :: four = 4.0_wp
real(kind=wp), private, parameter :: six = 6.0_wp
real(kind=wp), private, parameter :: ten = 10.0_wp
real(kind=wp), private, parameter :: d1mach4 = epsilon(one)

d1mach(4) -- the largest relative spacing


Functions

private pure function dpchsd(s1, s2, h1, h2)

inline function for weighted average of slopes.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: s1
real(kind=wp), intent(in) :: s2
real(kind=wp), intent(in) :: h1
real(kind=wp), intent(in) :: h2

Return Value real(kind=wp)

private function dchfcm(d1, d2, delta) result(ismon)

Check a single cubic for monotonicity

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: d1

derivative value at the end of an interval.

real(kind=wp), intent(in) :: d2

derivative value at the end of an interval.

real(kind=wp), intent(in) :: delta

the data slope over that interval

Return Value integer

indicates the monotonicity of the cubic segment:

  • ISMON = -3 if function is probably decreasing
  • ISMON = -1 if function is strictly decreasing
  • ISMON = 0 if function is constant
  • ISMON = 1 if function is strictly increasing
  • ISMON = 2 if function is non-monotonic
  • ISMON = 3 if function is probably increasing

If ABS(ISMON)=3, the derivative values are too close to the boundary of the monotonicity region to declare monotonicity in the presence of roundoff error.

private pure function dchfie(x1, x2, f1, f2, d1, d2, a, b) result(integral)

Cubic Hermite Function Integral Evaluator

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: x1

endpoints if interval of definition of cubic

real(kind=wp), intent(in) :: x2

endpoints if interval of definition of cubic

real(kind=wp), intent(in) :: f1

function values at the ends of the interval

real(kind=wp), intent(in) :: f2

function values at the ends of the interval

real(kind=wp), intent(in) :: d1

derivative values at the ends of the interval

real(kind=wp), intent(in) :: d2

derivative values at the ends of the interval

real(kind=wp), intent(in) :: a

endpoints of interval of integration

real(kind=wp), intent(in) :: b

endpoints of interval of integration

Return Value real(kind=wp)

private function dpchdf(k, x, s, ierr) result(deriv)

Computes divided differences for dpchce and dpchsp

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: k

is the order of the desired derivative approximation. K must be at least 3 (error return if not).

real(kind=wp), intent(in) :: x(k)

contains the K values of the independent variable. X need not be ordered, but the values MUST be distinct. (Not checked here.)

real(kind=wp), intent(inout) :: s(k)

contains the associated slope values: S(I) = (F(I+1)-F(I))/(X(I+1)-X(I)), I=1(1)K-1. (Note that S need only be of length K-1.) Will be destroyed on output.

integer, intent(out) :: ierr

will be set to -1 if K<2 .

Return Value real(kind=wp)

will be set to the desired derivative approximation if IERR=0 or to zero if IERR=-1.

public function dpchia(n, x, f, d, incfd, skip, a, b, ierr) result(value)

Piecewise Cubic Hermite Integrator, Arbitrary limits

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

number of data points. (Error return if N<2).

real(kind=wp), intent(in) :: x(*)

array of independent variable values. The elements of X must be strictly increasing: X(I-1) < X(I), I = 2(1)N. (Error return if not.)

real(kind=wp), intent(in) :: f(incfd,*)

array of function values. F(1+(I-1)*INCFD) is the value corresponding to X(I).

real(kind=wp), intent(in) :: d(incfd,*)

array of derivative values. D(1+(I-1)*INCFD) is the value corresponding to X(I).

integer, intent(in) :: incfd

increment between successive values in F and D. (Error return if INCFD<1).

logical, intent(inout) :: skip

logical variable which should be set to .TRUE. if the user wishes to skip checks for validity of preceding parameters, or to .FALSE. otherwise. This will save time in case these checks have already been performed (say, in dpchim or dpchic). SKIP will be set to .TRUE. on return with IERR>=0 .

real(kind=wp), intent(in) :: a

the limits of integration. NOTE: There is no requirement that [A,B] be contained in [X(1),X(N)]. However, the resulting integral value will be highly suspect, if not.

real(kind=wp), intent(in) :: b

the limits of integration. NOTE: There is no requirement that [A,B] be contained in [X(1),X(N)]. However, the resulting integral value will be highly suspect, if not.

integer, intent(inout) :: ierr

error flag.

Read more…

Return Value real(kind=wp)

value of the requested integral.

public function dpchid(n, x, f, d, incfd, skip, ia, ib, ierr) result(value)

Piecewise Cubic Hermite Integrator, Data limits

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

number of data points. (Error return if N<2).

real(kind=wp), intent(in) :: x(*)

array of independent variable values. The elements of X must be strictly increasing: X(I-1) < X(I), I = 2(1)N. (Error return if not.)

real(kind=wp), intent(in) :: f(incfd,*)

array of function values. F(1+(I-1)*INCFD) is the value corresponding to X(I).

real(kind=wp), intent(in) :: d(incfd,*)

array of derivative values. D(1+(I-1)*INCFD) is the value corresponding to X(I).

integer, intent(in) :: incfd

increment between successive values in F and D. (Error return if INCFD<1).

logical, intent(inout) :: skip

logical variable which should be set to .TRUE. if the user wishes to skip checks for validity of preceding parameters, or to .FALSE. otherwise. This will save time in case these checks have already been performed (say, in dpchim or dpchic). SKIP will be set to .TRUE. on return with IERR = 0 or -4.

integer, intent(in) :: ia

indices in X-array for the limits of integration. both must be in the range [1,N]. (Error return if not.) No restrictions on their relative values.

integer, intent(in) :: ib

indices in X-array for the limits of integration. both must be in the range [1,N]. (Error return if not.) No restrictions on their relative values.

integer, intent(out) :: ierr

error flag.

Read more…

Return Value real(kind=wp)

value of the requested integral.

private pure function dpchst(arg1, arg2) result(s)

DPCHIP Sign-Testing Routine

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: arg1
real(kind=wp), intent(in) :: arg2

Return Value real(kind=wp)


Subroutines

public subroutine dchfdv(x1, x2, f1, f2, d1, d2, ne, xe, fe, de, next, ierr)

Cubic Hermite Function and Derivative Evaluator

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: x1

initial endpoint of interval of definition of cubic. (Error return if X1==X2.)

real(kind=wp), intent(in) :: x2

final endpoint of interval of definition of cubic. (Error return if X1==X2.)

real(kind=wp), intent(in) :: f1

value of function at X1.

real(kind=wp), intent(in) :: f2

value of function at X2.

real(kind=wp), intent(in) :: d1

value of derivative at X1.

real(kind=wp), intent(in) :: d2

value of derivative at X2.

integer, intent(in) :: ne

number of evaluation points. (Error return if NE<1).

real(kind=wp), intent(in), dimension(*) :: xe

array of points at which the functions are to be evaluated. If any of the XE are outside the interval [X1,X2], a warning error is returned in NEXT.

real(kind=wp), intent(out), dimension(*) :: fe

array of values of the cubic function defined by X1,X2, F1,F2, D1,D2 at the points XE.

real(kind=wp), intent(out), dimension(*) :: de

array of values of the first derivative of the same function at the points XE.

integer, intent(out), dimension(2) :: next

array indicating number of extrapolation points:

Read more…
integer, intent(out) :: ierr

error flag.

Read more…

public subroutine dchfev(x1, x2, f1, f2, d1, d2, ne, xe, fe, next, ierr)

Cubic Hermite Function Evaluator

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: x1

initial endpoint of interval of definition of cubic. (Error return if X1==X2.)

real(kind=wp), intent(in) :: x2

final endpoint of interval of definition of cubic. (Error return if X1==X2.)

real(kind=wp), intent(in) :: f1

value of function at X1.

real(kind=wp), intent(in) :: f2

value of function at X2.

real(kind=wp), intent(in) :: d1

value of derivative at X1.

real(kind=wp), intent(in) :: d2

value of derivative at X2.

integer, intent(in) :: ne

number of evaluation points. (Error return if NE<1).

real(kind=wp), intent(in), dimension(*) :: xe

array of points at which the function is to be evaluated. If any of the XE are outside the interval [X1,X2], a warning error is returned in NEXT.

real(kind=wp), intent(out), dimension(*) :: fe

array of values of the cubic function defined by X1,X2, F1,F2, D1,D2 at the points XE.

integer, intent(out), dimension(2) :: next

integer array indicating number of extrapolation points:

Read more…
integer, intent(out) :: ierr

error flag.

Read more…

public subroutine dpchbs(n, x, f, d, incfd, knotyp, nknots, t, bcoef, ndim, kord, ierr)

Piecewise Cubic Hermite to B-Spline converter.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

the number of data points, N>=2 . (not checked)

real(kind=wp), intent(in) :: x(*)

array of independent variable values. The elements of X must be strictly increasing: X(I-1) < X(I), I = 2(1)N. (not checked) nmax, the dimension of X, must be >=N.

real(kind=wp), intent(in) :: f(incfd,*)

array of dependent variable values. F(1+(I-1)*INCFD) is the value corresponding to X(I). nmax, the second dimension of F, must be >=N.

real(kind=wp), intent(in) :: d(incfd,*)

array of derivative values at the data points. D(1+(I-1)*INCFD) is the value corresponding to X(I). nmax, the second dimension of D, must be >=N.

integer, intent(in) :: incfd

the increment between successive values in F and D. This argument is provided primarily for 2-D applications. It may have the value 1 for one-dimensional applications, in which case F and D may be singly-subscripted arrays.

integer, intent(in) :: knotyp

a flag to control the knot sequence. The knot sequence T is normally computed from X by putting a double knot at each X and setting the end knot pairs according to the value of KNOTYP:

Read more…
integer, intent(inout) :: nknots

the number of knots.

Read more…
real(kind=wp), intent(inout) :: t(*)

array of 2*N+4 knots for the B-representation.

Read more…
real(kind=wp), intent(out) :: bcoef(*)

array of 2*N B-spline coefficients.

integer, intent(out) :: ndim

the dimension of the B-spline space. (Set to 2*N.)

integer, intent(out) :: kord

the order of the B-spline. (Set to 4.)

integer, intent(out) :: ierr

an error flag:

Read more…

private subroutine dpchce(ic, vc, n, x, h, slope, d, incfd, ierr)

Set boundary conditions for dpchic

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: ic(2)

array of length 2 specifying desired boundary conditions:

Read more…
real(kind=wp), intent(in) :: vc(2)

array of length 2 specifying desired boundary values.

Read more…
integer, intent(in) :: n

number of data points. (assumes N>=2)

real(kind=wp), intent(in) :: x(*)

array of independent variable values. (the elements of X are assumed to be strictly increasing.)

real(kind=wp), intent(in) :: h(*)

array of interval lengths.

real(kind=wp), intent(in) :: slope(*)

array of data slopes. If the data are (X(I),Y(I)), I=1(1)N, then these inputs are:

Read more…
real(kind=wp), intent(inout) :: d(incfd,*)

(input) real8 array of derivative values at the data points. The value corresponding to X(I) must be stored in D(1+(I-1)INCFD), I=1(1)N.

Read more…
integer, intent(in) :: incfd

increment between successive values in D. This argument is provided primarily for 2-D applications.

integer, intent(out) :: ierr

error flag.

Read more…

private subroutine dpchci(n, h, slope, d, incfd)

Set interior derivatives for DPCHIC

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

number of data points. If N=2, simply does linear interpolation.

real(kind=wp), intent(in) :: h(*)

array of interval lengths.

real(kind=wp), intent(in) :: slope(*)

array of data slopes. If the data are (X(I),Y(I)), I=1(1)N, then these inputs are:

Read more…
real(kind=wp), intent(out) :: d(incfd,*)

array of derivative values at data points. If the data are monotonic, these values will determine a a monotone cubic Hermite function. The value corresponding to X(I) is stored in D(1+(I-1)*INCFD), I=1(1)N. No other entries in D are changed.

integer, intent(in) :: incfd

increment between successive values in D. This argument is provided primarily for 2-D applications.

public subroutine dpchcm(n, x, f, d, incfd, skip, ismon, ierr)

Check a cubic Hermite function for monotonicity

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

the number of data points. (Error return if N<2).

real(kind=wp), intent(in) :: x(n)

array of independent variable values. The elements of X must be strictly increasing: X(I-1) < X(I), I = 2(1)N. (Error return if not.)

real(kind=wp), intent(in) :: f(incfd,n)

array of function values. F(1+(I-1)*INCFD) is the value corresponding to X(I).

real(kind=wp), intent(in) :: d(incfd,n)

array of derivative values. D(1+(I-1)*INCFD) is is the value corresponding to X(I).

integer, intent(in) :: incfd

the increment between successive values in F and D. (Error return if INCFD<1).

logical, intent(inout) :: skip

logical variable which should be set to .TRUE. if the user wishes to skip checks for validity of preceding parameters, or to .FALSE. otherwise. This will save time in case these checks have already been performed. SKIP will be set to .TRUE. on normal return.

integer, intent(out) :: ismon(n)

array indicating on which intervals the PCH function defined by N, X, F, D is monotonic. For data interval [X(I),X(I+1)]:

Read more…
integer, intent(out) :: ierr

error flag.

Read more…

private subroutine dpchcs(switch, n, h, slope, d, incfd, ierr)

Monotonicity Switch Derivative Setter

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: switch

indicates the amount of control desired over local excursions from data.

integer, intent(in) :: n

number of data points. (assumes N>2).

real(kind=wp), intent(in) :: h(*)

array of interval lengths.

real(kind=wp), intent(in) :: slope(*)

array of data slopes. If the data are (X(I),Y(I)), I=1(1)N, then these inputs are:

Read more…
real(kind=wp), intent(inout) :: d(incfd,*)

(input) array of derivative values at the data points, as determined by dpchci.

Read more…
integer, intent(in) :: incfd

increment between successive values in D. This argument is provided primarily for 2-D applications.

integer, intent(out) :: ierr

error flag. should be zero. If negative, trouble in dpchsw. (should never happen.)

public subroutine dpchfd(n, x, f, d, incfd, skip, ne, xe, fe, de, ierr)

Piecewise Cubic Hermite Function and Derivative evaluator

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

number of data points. (Error return if N<2).

real(kind=wp), intent(in) :: x(*)

array of independent variable values. The elements of X must be strictly increasing: X(I-1) < X(I), I = 2(1)N. (Error return if not.)

real(kind=wp), intent(in) :: f(incfd,*)

array of function values. F(1+(I-1)*INCFD) is the value corresponding to X(I).

real(kind=wp), intent(in) :: d(incfd,*)

array of derivative values. D(1+(I-1)*INCFD) is the value corresponding to X(I).

integer, intent(in) :: incfd

increment between successive values in F and D. (Error return if INCFD<1).

logical, intent(inout) :: skip

logical variable which should be set to .TRUE. if the user wishes to skip checks for validity of preceding parameters, or to .FALSE. otherwise. This will save time in case these checks have already been performed (say, in dpchim or dpchic). SKIP will be set to .TRUE. on normal return.

integer, intent(in) :: ne

number of evaluation points. (Error return if NE<1).

real(kind=wp), intent(in) :: xe(*)

array of points at which the functions are to be evaluated. NOTES:

Read more…
real(kind=wp), intent(out) :: fe(*)

array of values of the cubic Hermite function defined by N, X, F, D at the points XE.

real(kind=wp), intent(out) :: de(*)

array of values of the first derivative of the same function at the points XE.

integer, intent(out) :: ierr

error flag.

Read more…

public subroutine dpchfe(n, x, f, d, incfd, skip, ne, xe, fe, ierr)

Piecewise Cubic Hermite Function Evaluator

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

number of data points. (Error return if N<2).

real(kind=wp), intent(in) :: x(*)

array of independent variable values. The elements of X must be strictly increasing: X(I-1) < X(I), I = 2(1)N. (Error return if not.)

real(kind=wp), intent(in) :: f(incfd,*)

array of function values. F(1+(I-1)*INCFD) is the value corresponding to X(I).

real(kind=wp), intent(in) :: d(incfd,*)

array of derivative values. D(1+(I-1)*INCFD) is the value corresponding to X(I).

integer, intent(in) :: incfd

increment between successive values in F and D. (Error return if INCFD<1).

logical, intent(inout) :: skip

logical variable which should be set to .TRUE. if the user wishes to skip checks for validity of preceding parameters, or to .FALSE. otherwise. This will save time in case these checks have already been performed (say, in dpchim or dpchic). SKIP will be set to .TRUE. on normal return.

integer, intent(in) :: ne

number of evaluation points. (Error return if NE<1).

real(kind=wp), intent(in) :: xe(*)

array of points at which the function is to be evaluated.

Read more…
real(kind=wp), intent(out) :: fe(*)

array of values of the cubic Hermite function defined by N, X, F, D at the points XE.

integer, intent(out) :: ierr

error flag.

Read more…

public subroutine dpchic(ic, vc, switch, n, x, f, d, incfd, wk, nwk, ierr)

Piecewise Cubic Hermite Interpolation Coefficients.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: ic(2)

array of length 2 specifying desired boundary conditions:

Read more…
real(kind=wp), intent(in) :: vc(2)

array of length 2 specifying desired boundary values, as indicated above.

Read more…
real(kind=wp), intent(in) :: switch

indicates desired treatment of points where direction of monotonicity switches: Set SWITCH to zero if interpolant is required to be mono- tonic in each interval, regardless of monotonicity of data.

Read more…
integer, intent(in) :: n

number of data points. (Error return if N<2).

real(kind=wp), intent(in) :: x(*)

array of independent variable values. The elements of X must be strictly increasing: X(I-1) < X(I), I = 2(1)N. (Error return if not.)

real(kind=wp), intent(in) :: f(incfd,*)

array of dependent variable values to be interpolated. F(1+(I-1)*INCFD) is value corresponding to X(I).

real(kind=wp), intent(out) :: d(incfd,*)

array of derivative values at the data points. These values will determine a monotone cubic Hermite function on each subinterval on which the data are monotonic, except possibly adjacent to switches in monotonicity. The value corresponding to X(I) is stored in D(1+(I-1)*INCFD), I=1(1)N. No other entries in D are changed.

integer, intent(in) :: incfd

increment between successive values in F and D. This argument is provided primarily for 2-D applications. (Error return if INCFD<1).

real(kind=wp), intent(inout) :: wk(nwk)

array of working storage. The user may wish to know that the returned values are:

Read more…
integer, intent(in) :: nwk

length of work array. (Error return if NWK<2*(N-1)).

integer, intent(out) :: ierr

error flag.

Read more…

public subroutine dpchim(n, x, f, d, incfd, ierr)

Piecewise Cubic Hermite Interpolation to Monotone data.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

number of data points. (Error return if N<2). If N=2, simply does linear interpolation.

real(kind=wp), intent(in) :: x(*)

array of independent variable values. The elements of X must be strictly increasing: X(I-1) < X(I), I = 2(1)N. (Error return if not.)

real(kind=wp), intent(in) :: f(incfd,*)

array of dependent variable values to be interpolated. F(1+(I-1)*INCFD) is value corresponding to X(I). DPCHIM is designed for monotonic data, but it will work for any F-array. It will force extrema at points where monotonicity switches direction. If some other treatment of switch points is desired, dpchic should be used instead.

real(kind=wp), intent(out) :: d(incfd,*)

array of derivative values at the data points. If the data are monotonic, these values will determine a monotone cubic Hermite function. The value corresponding to X(I) is stored in D(1+(I-1)*INCFD), I=1(1)N. No other entries in D are changed.

integer, intent(in) :: incfd

increment between successive values in F and D. This argument is provided primarily for 2-D applications. (Error return if INCFD<1).

integer, intent(out) :: ierr

error flag.

Read more…

private subroutine dpchkt(n, x, knotyp, t)

Set a knot sequence for the B-spline representation of a PCH function with breakpoints X. All knots will be at least double. Endknots are set as:

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n
real(kind=wp), intent(in) :: x(*)
integer, intent(in) :: knotyp
real(kind=wp), intent(out) :: t(*)

public subroutine dpchsp(ic, vc, n, x, f, d, incfd, wk, nwk, ierr)

Piecewise Cubic Hermite Spline

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: ic(2)

integer array of length 2 specifying desired boundary conditions:

Read more…
real(kind=wp), intent(in) :: vc(2)

array of length 2 specifying desired boundary values, as indicated above.

Read more…
integer, intent(in) :: n

number of data points. (Error return if N<2).

real(kind=wp), intent(in) :: x(*)

array of independent variable values. The elements of X must be strictly increasing: X(I-1) < X(I), I = 2(1)N. (Error return if not.)

real(kind=wp), intent(in) :: f(incfd,*)

array of dependent variable values to be interpolated. F(1+(I-1)*INCFD) is value corresponding to X(I).

real(kind=wp), intent(out) :: d(incfd,*)

array of derivative values at the data points. These values will determine the cubic spline interpolant with the requested boundary conditions. The value corresponding to X(I) is stored in D(1+(I-1)*INCFD), I=1(1)N. No other entries in D are changed.

integer, intent(in) :: incfd

increment between successive values in F and D. This argument is provided primarily for 2-D applications. (Error return if INCFD<1).

real(kind=wp), intent(inout) :: wk(2,*)

array of working storage.

integer, intent(in) :: nwk

length of work array. (Error return if NWK<2*N).

integer, intent(out) :: ierr

error flag.

Read more…

private subroutine dpchsw(dfmax, iextrm, d1, d2, h, slope, ierr)

Switch Excursion Limiter.

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: dfmax

maximum allowed difference between F(IEXTRM) and the cubic determined by derivative values D1,D2. (assumes DFMAX>0.)

integer, intent(in) :: iextrm

index of the extreme data value. (assumes IEXTRM = 1 or 2 . Any value /=1 is treated as 2.)

real(kind=wp), intent(inout) :: d1

derivative values at the ends of the interval. (Assumes D1*D2 <= 0.) (output) may be modified if necessary to meet the restriction imposed by DFMAX.

real(kind=wp), intent(inout) :: d2

derivative values at the ends of the interval. (Assumes D1*D2 <= 0.) (output) may be modified if necessary to meet the restriction imposed by DFMAX.

real(kind=wp), intent(in) :: h

interval length. (Assumes H>0.)

real(kind=wp), intent(in) :: slope

data slope on the interval.

integer, intent(out) :: ierr

error flag. should be zero.

Read more…

private subroutine xermsg(librar, subrou, messg, nerr, level)

Error reporter.

Read more…

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: librar

name of the library

character(len=*), intent(in) :: subrou

name of the routine that detected the error

character(len=*), intent(in) :: messg

text of the error or warning message

integer, intent(in) :: nerr

error code

integer, intent(in) :: level

value in the range 0 to 2 that indicates the level (severity) of the error.