Christian Soeller > PDL-2.4.1 > PDL::Slatec

Download:
PDL-2.4.1.tar.gz

Annotate this POD

Related Modules

Math::MatrixReal
Inline::C
HTML::Template
Math::Pari
Text::Template
GD::Graph
more...
By perlmonks.org

CPAN RT

Open  0
View/Report Bugs
Source   Latest Release: PDL-2.007_04

NAME ^

PDL::Slatec - PDL interface to the slatec numerical programming library

SYNOPSIS ^

 use PDL::Slatec;

 ($ndeg, $r, $ierr, $a) = polyfit($x, $y, $w, $maxdeg, $eps);

DESCRIPTION ^

This module serves the dual purpose of providing an interface to parts of the slatec library and showing how to interface PDL to an external library. Using this library requires a fortran compiler; the source for the routines is provided for convenience.

Currently available are routines to: manipulate matrices; calculate FFT's; fit data using polynomials; and interpolate/integrate data using piecewise cubic Hermite interpolation.

Piecewise cubic Hermite interpolation (PCHIP)

PCHIP is the slatec package of routines to perform piecewise cubic Hermite interpolation of data. It features software to produce a monotone and "visually pleasing" interpolant to monotone data. According to Fritsch & Carlson ("Monotone piecewise cubic interpolation", SIAM Journal on Numerical Analysis 17, 2 (April 1980), pp. 238-246), such an interpolant may be more reasonable than a cubic spline if the data contains both "steep" and "flat" sections. Interpolation of cumulative probability distribution functions is another application. These routines are cryptically named (blame FORTRAN), beginning with 'ch', and accept either float or double piddles.

Most of the routines require an integer parameter called check; if set to 0, then no checks on the validity of the input data are made, otherwise these checks are made. The value of check can be set to 0 if a routine such as chim has already been successfully called.

EOD # ' un-confuse emacs

# if define chbs, then add something like the following to point 3: # # or use chbs to convert a PCH function into B-representation # for use with the B-spline routines of slatec # (although no interface to them currently exist). #

# add function definitions after finishing the first pp_addpm(), since this # adds a '=head1 FUNCTIONS' line at the end of the text

pp_addpm(<<'END'); =head2 eigsys

Eigenvalues and eigenvectors of a real positive definite symmetric matrix.

 ($eigvals,$eigvecs) = eigsys($mat)

Note: this function should be extended to calculate only eigenvalues if called in scalar context!

matinv

Inverse of a square matrix

 ($inv) = matinv($mat)

polyfit

Convenience wrapper routine about the polfit slatec function. Separates supplied arguments and return values.

Fit discrete data in a least squares sense by polynomials in one variable. Handles threading correctly--one can pass in a 2D PDL (as $y) and it will pass back a 2D PDL, the rows of which are the polynomial regression results (in $r corresponding to the rows of $y.

 ($ndeg, $r, $ierr, $a) = polyfit($x, $y, $w, $maxdeg, $eps);

 where on input:

 C<x> and C<y> are the values to fit to a polynomial.
 C<w> are weighting factors
 C<maxdeg> is the maximum degree of polynomial to use and 
 C<eps> is the required degree of fit.

 and on output:

 C<ndeg> is the degree of polynomial actually used
 C<r> is the values of the fitted polynomial 
 C<ierr> is a return status code, and
 C<a> is some working array or other
 C<eps> is modified to contain the rms error of the fit.

This version of polyfit handles bad values correctly. It strips them out of the $x variable and creates an appropriate $y variable containing indices of the non-bad members of $x before calling the Slatec routine polfit.

polycoef

Convenience wrapper routine around the pcoef slatec function. Separates supplied arguments and return values.

Convert the polyfit/polfit coefficients to Taylor series form.

 $tc = polycoef($l, $c, $a);

polyvalue

Convenience wrapper routine around the pvalue slatec function. Separates supplied arguments and return values.

For multiple input x positions, a corresponding y position is calculated.

The derivatives PDL is one dimensional (of size nder) if a single x position is supplied, two dimensional if more than one x position is supplied.

Use the coefficients generated by polyfit (or polfit) to evaluate the polynomial fit of degree l, along with the first nder of its derivatives, at a specified point.

 ($yfit, $yp) = polyvalue($l, $nder, $x, $a);

detslatec

compute the determinant of an invertible matrix

  $mat = zeroes(5,5); $mat->diagonal(0,1) .= 1; # unity matrix
  $det = detslatec $mat;

Usage:

  $determinant = detslatec $matrix;
  Signature: detslatec(mat(n,m); [o] det())

detslatec computes the determinant of an invertible matrix and barfs if the matrix argument provided is non-invertible. The matrix threads as usual.

This routine was previously known as det which clashes now with det:PDL::MatrixOps/det which is provided by PDL::MatrixOps. For the moment PDL::Slatec will also load PDL::MatrixOps thereby making sure that older scripts work.

# switch has become mflag, since `switch' is a reserved word in # C. # # can not say (nwk=2*n) --- the rhs has to equal a number # -> could Basic/Gen/PP/Dims.pm be hacked to allow this? # # I didn't have much success with preceeding wk by [t] # defslatec( 'chic', {S => 'pchic', D => 'dpchic'}, 'IntFlag ic (two=2); Mat vc (two=2); Mat mflag (); Dim n; Mat x (n); Mat f (n); Mat [o] d (n); Incfd dummy; Mat wk (nwk); Dim nwk; IntFlag [o] ierr (); ', 'Calculate the derivatives at the given points ($x,$f, where $x is strictly increasing). Control over the boundary conditions is given by the $ic and $vc piddles, and the value of $mflag determines the treatment of points where monotoncity switches direction. A simpler, more restricted, interface is available using chim.

The first and second elements of $ic determine the boundary conditions at the start and end of the data respectively. If the value is 0, then the default condition, as used by chim, is adopted. If greater than zero, no adjustment for monotonicity is made, otherwise if less than zero the derivative will be adjusted. The allowed magnitudes for ic(0) are:

The values for ic(1) are the same as above, except that the first-derivative value is stored in vc(1) for cases 1 and 2. The values of $vc need only be set if options 1 or 2 are chosen for $ic.

Set $mflag = 0 if interpolant is required to be monotonic in each interval, regardless of the data. This causes $d to be set to 0 at all switch points. Set $mflag to be non-zero to use a formula based on the 3-point difference formula at switch points. If $mflag > 0, then the interpolant at swich points is forced to not deviate from the data by more than $mflag*dfloc, where dfloc is the maximum of the change of $f on this interval and its two immediate neighbours. If $mflag < 0, no such control is to be imposed.

The piddle $wk is only needed for work space. However, I could not get it to work as a temporary variable, so you must supply it; it is a 1D piddle with 2*n elements.

Error status returned by $ierr:

# as above, have made wk an actual piddle, rather than a [t] defslatec( 'chsp', {S => 'pchsp', D => 'dpchsp'}, 'IntFlag ic (two=2); Mat vc (two=2); Dim n; Mat x (n); Mat f (n); Mat [o] d (n); Incfd dummy; Mat wk (nwk); Dim nwk; IntFlag [o] ierr (); ', 'Calculate the derivatives, using cubic spline interpolation, at the given points ($x,$f), with the specified boundary conditions. Control over the boundary conditions is given by the $ic and $vc piddles. The resulting values - $x,$f,$d - can be used in all the functions which expect a cubic Hermite function.

The first and second elements of $ic determine the boundary conditions at the start and end of the data respectively. The allowed values for ic(0) are:

The values for ic(1) are the same as above, except that the first-derivative value is stored in vc(1) for cases 1 and 2. The values of $vc need only be set if options 1 or 2 are chosen for $ic.

The piddle $wk is only needed for work space. However, I could not get it to work as a temporary variable, so you must supply it; it is a 1D piddle with 2*n elements.

Error status returned by $ierr:

defslatec( 'chfd', {S => 'pchfd', D => 'dpchfd'}, 'Dim n; Mat x (n); Mat f (n); Mat d (n); Incfd dummy; CheckFlag check (); Dim ne; Mat xe (ne); Mat [o] fe (ne); Mat [o] de (ne); IntFlag [o] ierr (); ', 'Given a piecewise cubic Hermite function - such as from chim - evaluate the function ($fe) and derivative ($de) at a set of points ($xe). If function values alone are required, use chfe. Set check to 0 to skip checks on the input data.

Error status returned by $ierr:

defslatec( 'chfe', {S => 'pchfe', D => 'dpchfe'}, 'Dim n; Mat x (n); Mat f (n); Mat d (n); Incfd dummy; CheckFlag check (); Dim ne; Mat xe (ne); Mat [o] fe (ne); IntFlag [o] ierr (); ', 'Given a piecewise cubic Hermite function - such as from chim - evaluate the function ($fe) at a set of points ($xe). If derivative values are also required, use chfd. Set check to 0 to skip checks on the input data.

Error status returned by $ierr:

defslatec( 'chia', {S => 'pchia', D => 'dpchia'}, 'Dim n; Mat x (n); Mat f (n); Mat d (n); Incfd dummy; CheckFlag check (); Mat a (); Mat b (); FuncRet [o] ans (); IntFlag [o] ierr (); ', 'Evaluate the definite integral of a a piecewise cubic Hermite function over an arbitrary interval, given by [$a,$b]. See chid if the integration limits are data points. Set check to 0 to skip checks on the input data.

The values of $a and $b do not have to lie within $x, although the resulting integral value will be highly suspect if they are not.

Error status returned by $ierr:

defslatec( 'chid', {S => 'pchid', D => 'dpchid'}, 'Dim n; Mat x (n); Mat f (n); Mat d (n); Incfd dummy; CheckFlag check (); FortranIndex ia (); FortranIndex ib (); FuncRet [o] ans (); IntFlag [o] ierr (); ', 'Evaluate the definite integral of a a piecewise cubic Hermite function between x($ia) and x($ib).

See chia for integration between arbitrary limits.

Although using a fortran routine, the values of $ia and $ib are zero offset. Set check to 0 to skip checks on the input data.

Error status returned by $ierr:

defslatec( 'chcm', {S => 'pchcm', D => 'dpchcm'}, 'Dim n; Mat x (n); Mat f (n); Mat d (n); Incfd dummy; CheckFlag check (); IntFlag [o] ismon (n); IntFlag [o] ierr (); ', 'The outout piddle $ismon indicates over which intervals the function is monotonic. Set check to 0 to skip checks on the input data.

For the data interval [x(i),x(i+1)], the values of ismon(i) can be:

If abs(ismon(i)) == 3, the derivative values are near the boundary of the monotonicity region. A small increase produces non-monotonicity, whereas a decrease produces strict monotonicity.

The above applies to i = 0 .. nelem($x)-1. The last element of $ismon indicates whether the entire function is monotonic over $x.

Error status returned by $ierr:

# XXX tsize = 2*n+4 # bsize = 2*n # # ndim gets set to 2*n # # Changed by routine: # nknots # t defslatec( 'chbs', {S => 'pchbs', D => 'dpchbs'}, 'Dim n; Mat x (n); Mat f (n); Mat d (n); Incfd dummy; IntFlag knotyp (); IntFlag nknots (); Mat t (tsize); Mat [o] bcoef (bsize); IntFlag [o] ndim (); IntFlag [o] kord (); IntFlag [o] ierr (); ', 'The resulting B-spline representation of the data (i.e. nknots, t, bcoeff, ndim, and kord) can be evaluated by bvalu (which is currently not available).

Array sizes: tsize = 2*n + 4, bsize = 2*n, and ndim = 2*n.

knotyp is a flag which controls 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 (where m = ndim = 2*n):

nknots is the number of knots and may be changed by the routine. If knotyp >= 0, nknots will be set to ndim+4, otherwise it is an input variable, and an error will occur if its value is not equal to ndim+4.

t is the array of 2*n+4 knots for the B-representation and may be changed by the routine. If knotyp >= 0, t will be changed so that the interior double knots are equal to the x-values and the boundary knots set as indicated above, otherwise it is assumed that t was set by a previous call (no check is made to verify that the data forms a legitimate knot sequence).

Error status returned by $ierr:

AUTHOR ^

Copyright (C) 1997 Tuomas J. Lukka. Copyright (C) 2000 Tim Jenness, Doug Burke. All rights reserved. There is no warranty. You are allowed to redistribute this software / documentation under certain conditions. For details, see the file COPYING in the PDL distribution. If this file is separated from the PDL distribution, the copyright notice should be included in the file.

syntax highlighting: