1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175
|
<!DOCTYPE html PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN" "http://www.w3.org/TR/html4/loose.dtd">
<html>
<!-- Copyright (C) 1996, 1997, 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010, 2011, 2012, 2013, 2014, 2015, 2016 The GSL Team.
Permission is granted to copy, distribute and/or modify this document
under the terms of the GNU Free Documentation License, Version 1.3 or
any later version published by the Free Software Foundation; with the
Invariant Sections being "GNU General Public License" and "Free Software
Needs Free Documentation", the Front-Cover text being "A GNU Manual",
and with the Back-Cover Text being (a) (see below). A copy of the
license is included in the section entitled "GNU Free Documentation
License".
(a) The Back-Cover Text is: "You have the freedom to copy and modify this
GNU Manual." -->
<!-- Created by GNU Texinfo 5.1, http://www.gnu.org/software/texinfo/ -->
<head>
<title>GNU Scientific Library – Reference Manual: 1D Interpolation Types</title>
<meta name="description" content="GNU Scientific Library – Reference Manual: 1D Interpolation Types">
<meta name="keywords" content="GNU Scientific Library – Reference Manual: 1D Interpolation Types">
<meta name="resource-type" content="document">
<meta name="distribution" content="global">
<meta name="Generator" content="makeinfo">
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
<link href="index.html#Top" rel="start" title="Top">
<link href="Function-Index.html#Function-Index" rel="index" title="Function Index">
<link href="Interpolation.html#Interpolation" rel="up" title="Interpolation">
<link href="1D-Index-Look_002dup-and-Acceleration.html#g_t1D-Index-Look_002dup-and-Acceleration" rel="next" title="1D Index Look-up and Acceleration">
<link href="1D-Interpolation-Functions.html#g_t1D-Interpolation-Functions" rel="previous" title="1D Interpolation Functions">
<style type="text/css">
<!--
a.summary-letter {text-decoration: none}
blockquote.smallquotation {font-size: smaller}
div.display {margin-left: 3.2em}
div.example {margin-left: 3.2em}
div.indentedblock {margin-left: 3.2em}
div.lisp {margin-left: 3.2em}
div.smalldisplay {margin-left: 3.2em}
div.smallexample {margin-left: 3.2em}
div.smallindentedblock {margin-left: 3.2em; font-size: smaller}
div.smalllisp {margin-left: 3.2em}
kbd {font-style:oblique}
pre.display {font-family: inherit}
pre.format {font-family: inherit}
pre.menu-comment {font-family: serif}
pre.menu-preformatted {font-family: serif}
pre.smalldisplay {font-family: inherit; font-size: smaller}
pre.smallexample {font-size: smaller}
pre.smallformat {font-family: inherit; font-size: smaller}
pre.smalllisp {font-size: smaller}
span.nocodebreak {white-space:nowrap}
span.nolinebreak {white-space:nowrap}
span.roman {font-family:serif; font-weight:normal}
span.sansserif {font-family:sans-serif; font-weight:normal}
ul.no-bullet {list-style: none}
-->
</style>
</head>
<body lang="en" bgcolor="#FFFFFF" text="#000000" link="#0000FF" vlink="#800080" alink="#FF0000">
<a name="g_t1D-Interpolation-Types"></a>
<div class="header">
<p>
Next: <a href="1D-Index-Look_002dup-and-Acceleration.html#g_t1D-Index-Look_002dup-and-Acceleration" accesskey="n" rel="next">1D Index Look-up and Acceleration</a>, Previous: <a href="1D-Interpolation-Functions.html#g_t1D-Interpolation-Functions" accesskey="p" rel="previous">1D Interpolation Functions</a>, Up: <a href="Interpolation.html#Interpolation" accesskey="u" rel="up">Interpolation</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<a name="g_t1D-Interpolation-Types-1"></a>
<h3 class="section">28.3 1D Interpolation Types</h3>
<a name="index-gsl_005finterp_005ftype"></a>
<p>The interpolation library provides the following interpolation types:
</p>
<dl>
<dt><a name="index-gsl_005finterp_005flinear"></a>Interpolation Type: <strong>gsl_interp_linear</strong></dt>
<dd><a name="index-linear-interpolation"></a>
<p>Linear interpolation. This interpolation method does not require any
additional memory.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005finterp_005fpolynomial"></a>Interpolation Type: <strong>gsl_interp_polynomial</strong></dt>
<dd><a name="index-polynomial-interpolation"></a>
<p>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.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005finterp_005fcspline"></a>Interpolation Type: <strong>gsl_interp_cspline</strong></dt>
<dd><a name="index-cubic-splines"></a>
<p>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.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005finterp_005fcspline_005fperiodic"></a>Interpolation Type: <strong>gsl_interp_cspline_periodic</strong></dt>
<dd><p>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.
</p>
</dd></dl>
<dl>
<dt><a name="index-gsl_005finterp_005fakima"></a>Interpolation Type: <strong>gsl_interp_akima</strong></dt>
<dd><a name="index-Akima-splines"></a>
<p>Non-rounded Akima spline with natural boundary conditions. This method
uses the non-rounded corner algorithm of Wodicka.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005finterp_005fakima_005fperiodic"></a>Interpolation Type: <strong>gsl_interp_akima_periodic</strong></dt>
<dd><p>Non-rounded Akima spline with periodic boundary conditions. This method
uses the non-rounded corner algorithm of Wodicka.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005finterp_005fsteffen"></a>Interpolation Type: <strong>gsl_interp_steffen</strong></dt>
<dd><p>Steffen’s method guarantees the monotonicity of the interpolating function
between the given data points. Therefore, minima and maxima can only occur
exactly at the data points, and there can never be spurious oscillations
between data points. The interpolated function is piecewise cubic
in each interval. The resulting curve and its first derivative
are guaranteed to be continuous, but the second derivative may be
discontinuous.
</p></dd></dl>
<p>The following related functions are available:
</p>
<dl>
<dt><a name="index-gsl_005finterp_005fname"></a>Function: <em>const char *</em> <strong>gsl_interp_name</strong> <em>(const gsl_interp * <var>interp</var>)</em></dt>
<dd><p>This function returns the name of the interpolation type used by <var>interp</var>.
For example,
</p>
<div class="example">
<pre class="example">printf ("interp uses '%s' interpolation.\n",
gsl_interp_name (interp));
</pre></div>
<p>would print something like,
</p>
<div class="example">
<pre class="example">interp uses 'cspline' interpolation.
</pre></div>
</dd></dl>
<dl>
<dt><a name="index-gsl_005finterp_005fmin_005fsize"></a>Function: <em>unsigned int</em> <strong>gsl_interp_min_size</strong> <em>(const gsl_interp * <var>interp</var>)</em></dt>
<dt><a name="index-gsl_005finterp_005ftype_005fmin_005fsize"></a>Function: <em>unsigned int</em> <strong>gsl_interp_type_min_size</strong> <em>(const gsl_interp_type * <var>T</var>)</em></dt>
<dd><p>These functions return the minimum number of points required by the
interpolation object <var>interp</var> or interpolation type <var>T</var>. For
example, Akima spline interpolation requires a minimum of 5 points.
</p></dd></dl>
<hr>
<div class="header">
<p>
Next: <a href="1D-Index-Look_002dup-and-Acceleration.html#g_t1D-Index-Look_002dup-and-Acceleration" accesskey="n" rel="next">1D Index Look-up and Acceleration</a>, Previous: <a href="1D-Interpolation-Functions.html#g_t1D-Interpolation-Functions" accesskey="p" rel="previous">1D Interpolation Functions</a>, Up: <a href="Interpolation.html#Interpolation" accesskey="u" rel="up">Interpolation</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
</body>
</html>
|