Better fix to Fit
[Math-GSL.git] / Interp.i
blob4e2785095befae28a25691fca45a85489067c590
1 %module "Math::GSL::Interp"
3 %include "typemaps.i"
4 %include "gsl_typemaps.i"
6 %{
7 #include "gsl/gsl_types.h"
8 #include "gsl/gsl_interp.h"
9 %}
10 %include "gsl/gsl_types.h"
11 %include "gsl/gsl_interp.h"
14 %perlcode %{
15 @EXPORT_OK = qw/
16 gsl_interp_accel_alloc
17 gsl_interp_accel_find
18 gsl_interp_accel_reset
19 gsl_interp_accel_free
20 gsl_interp_alloc
21 gsl_interp_init
22 gsl_interp_name
23 gsl_interp_min_size
24 gsl_interp_eval_e
25 gsl_interp_eval
26 gsl_interp_eval_deriv_e
27 gsl_interp_eval_deriv
28 gsl_interp_eval_deriv2_e
29 gsl_interp_eval_deriv2
30 gsl_interp_eval_integ_e
31 gsl_interp_eval_integ
32 gsl_interp_free
33 gsl_interp_bsearch
34 $gsl_interp_linear
35 $gsl_interp_polynomial
36 $gsl_interp_cspline
37 $gsl_interp_cspline_periodic
38 $gsl_interp_akima
39 $gsl_interp_akima_periodic
41 %EXPORT_TAGS = ( all => [ @EXPORT_OK ] );
43 __END__
45 =head1 NAME
47 Math::GSL::Interp - Functions for performing interpolation
49 =head1 SYNOPSIS
51 use Math::GSL::Interp qw/:all/;
53 =head1 DESCRIPTION
55 Here is a list of all the functions included in this module :
57 =over 1
59 =item C<gsl_interp_accel_alloc()> - This function returns a pointer to an accelerator object, which is a kind of iterator for interpolation lookups. It tracks the state of lookups, thus allowing for application of various acceleration strategies.
61 =item C<gsl_interp_accel_find($a, $x_array, $size, $x)> - This function performs a lookup action on the data array $x_array of size $size, using the given accelerator $a. This is how lookups are performed during evaluation of an interpolation. The function returns an index i such that $x_array[i] <= $x < $x_array[i+1].
63 =item C<gsl_interp_accel_reset>
65 =item C<gsl_interp_accel_free($a)> - This function frees the accelerator object $a.
67 =item C<gsl_interp_alloc($T, $alloc)> - This function returns a newly allocated interpolation object of type $T for $size data-points. $T must be one of the constants below.
69 =item C<gsl_interp_init($interp, $xa, $ya, $size)> - This function initializes the interpolation object interp for the data (xa,ya) where xa and ya are arrays of size size. The interpolation object (gsl_interp) does not save the data arrays xa and ya and only stores the static state computed from the data. The xa data array is always assumed to be strictly ordered, with increasing x values; the behavior for other arrangements is not defined.
71 =item C<gsl_interp_name>
73 =item C<gsl_interp_min_size>
75 =item C<gsl_interp_eval_e>
77 =item C<gsl_interp_eval>
79 =item C<gsl_interp_eval_deriv_e>
81 =item C<gsl_interp_eval_deriv>
83 =item C<gsl_interp_eval_deriv2_e>
85 =item C<gsl_interp_eval_deriv2>
87 =item C<gsl_interp_eval_integ_e>
89 =item C<gsl_interp_eval_integ>
91 =item C<gsl_interp_free>
93 =item C<gsl_interp_bsearch>
95 =back
97 This module also includes the following constants :
99 =over 1
101 =item C<$gsl_interp_linear> - Linear interpolation
103 =item C<$gsl_interp_polynomial> - 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.
105 =item C<$gsl_interp_cspline> - 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 first point and last point.
107 =item C<$gsl_interp_cspline_periodic> - 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.
109 =item C<$gsl_interp_akima> - Non-rounded Akima spline with natural boundary conditions. This method uses the non-rounded corner algorithm of Wodicka.
111 =item C<$gsl_interp_akima_periodic> - Non-rounded Akima spline with periodic boundary conditions. This method uses the non-rounded corner algorithm of Wodicka.
113 =back
115 =head1 AUTHORS
117 Jonathan Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>
119 =head1 COPYRIGHT AND LICENSE
121 Copyright (C) 2008 Jonathan Leto and Thierry Moisan
123 This program is free software; you can redistribute it and/or modify it
124 under the same terms as Perl itself.
126 =cut