CoastalME (Coastal Modelling Environment)
Simulates the long-term behaviour of complex coastlines
Loading...
Searching...
No Matches
hermite_cubic.cpp
Go to the documentation of this file.
1
10
11/* ===============================================================================================================================
12 This file is part of CoastalME, the Coastal Modelling Environment.
13
14 CoastalME is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.
15
16 This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
17
18 You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
19===============================================================================================================================*/
20#include <iostream>
21using std::cerr;
22using std::endl;
23using std::ios;
24
25#include "hermite_cubic.h"
26
27//===============================================================================================================================
28//
29// Purpose:
30//
31// R8VEC_BRACKET3 finds the interval containing or nearest a given value.
32//
33// Discussion:
34//
35// An R8VEC is a vector of R8's.
36//
37// The routine always returns the index LEFT of the sorted array
38// T with the property that either
39// * T is contained in the interval [ T[LEFT], T[LEFT+1] ], or
40// * T < T[LEFT] = T[0], or
41// * T > T[LEFT+1] = T[N-1].
42//
43// The routine is useful for interpolation problems, where
44// the abscissa must be located within an interval of data
45// abscissas for interpolation, or the "nearest" interval
46// to the (extreme) abscissa must be found so that extrapolation
47// can be carried out.
48//
49// This version of the function has been revised so that the value of
50// LEFT that is returned uses the 0-based indexing natural to C++.
51//
52// Licensing:
53//
54// This code is distributed under the GNU LGPL license.
55//
56// Modified:
57//
58// 30 April 2009
59//
60// Author:
61//
62// John Burkardt
63//
64// Parameters:
65//
66// Input, int N, length of the input array.
67//
68// Input, double T[N], an array that has been sorted into ascending order.
69//
70// Input, double TVAL, a value to be bracketed by entries of T.
71//
72// Input/output, int *LEFT.
73// On input, if 0 <= LEFT <= N-2, LEFT is taken as a suggestion for the
74// interval [ T[LEFT-1] T[LEFT] ] in which TVAL lies. This interval
75// is searched first, followed by the appropriate interval to the left
76// or right. After that, a binary search is used.
77// On output, LEFT is set so that the interval [ T[LEFT], T[LEFT+1] ]
78// is the closest to TVAL; it either contains TVAL, or else TVAL
79// lies outside the interval [ T[0], T[N-1] ].
80//
81//===============================================================================================================================
83//===============================================================================================================================
84void r8vec_bracket3(int const n, double const* t, double const tval, int* left)
85{
86 int high;
87 int low;
88 int mid;
89
90 // Check the input data
91 if (n < 2)
92 {
93 cerr << endl
94 << "R8VEC_BRACKET3 - Fatal error! N must be at least 2." << endl;
95 return;
96 }
97
98 // If *left is not between 0 and n-2, set it to the middle value
99 if ((*left < 0) || (n - 2 < *left))
100 {
101 *left = (n - 1) / 2;
102 }
103
104 // CASE 1: TVAL < T[*LEFT]: search for TVAL in (T[I],T[I+1]), for I = 0 to *LEFT-1.
105 if (tval < t[*left])
106 {
107 if (*left == 0)
108 {
109 return;
110 }
111
112 else if (*left == 1)
113 {
114 *left = 0;
115 return;
116 }
117
118 else if (t[*left - 1] <= tval)
119 {
120 *left = *left - 1;
121 return;
122 }
123
124 else if (tval <= t[1])
125 {
126 *left = 0;
127 return;
128 }
129
130 //
131 // ...Binary search for TVAL in (T[I],T[I+1]), for I = 1 to *LEFT-2.
132 //
133 low = 1;
134 high = *left - 2;
135
136 for (;;)
137 {
138 if (low == high)
139 {
140 *left = low;
141 return;
142 }
143
144 mid = (low + high + 1) / 2;
145
146 if (t[mid] <= tval)
147 {
148 low = mid;
149 }
150
151 else
152 {
153 high = mid - 1;
154 }
155 }
156 }
157
158 //
159 // CASE 2: T[*LEFT+1] < TVAL:
160 // Search for TVAL in (T[I],T[I+1]) for intervals I = *LEFT+1 to N-2.
161 //
162 else if (t[*left + 1] < tval)
163 {
164 if (*left == n - 2)
165 {
166 return;
167 }
168
169 else if (*left == n - 3)
170 {
171 *left = *left + 1;
172 return;
173 }
174
175 else if (tval <= t[*left + 2])
176 {
177 *left = *left + 1;
178 return;
179 }
180
181 else if (t[n - 2] <= tval)
182 {
183 *left = n - 2;
184 return;
185 }
186
187 // ...Binary search for TVAL in (T[I],T[I+1]) for intervals I = *LEFT+2 to N-3.
188 low = *left + 2;
189 high = n - 3;
190
191 for (;;)
192 {
193 if (low == high)
194 {
195 *left = low;
196 return;
197 }
198
199 mid = (low + high + 1) / 2;
200
201 if (t[mid] <= tval)
202 {
203 low = mid;
204 }
205
206 else
207 {
208 high = mid - 1;
209 }
210 }
211 }
212
213 // CASE 3: T[*LEFT] <= TVAL <= T[*LEFT+1]: T is just where the user said it might be.
214 else
215 {
216 }
217
218 return;
219}
220
221//===============================================================================================================================
222//
223// Purpose:
224//
225// HERMITE_CUBIC_VALUE evaluates a Hermite cubic polynomial.
226//
227// Discussion:
228//
229// The input arguments can be vectors.
230//
231// Licensing:
232//
233// This code is distributed under the GNU LGPL license.
234//
235// Modified:
236//
237// 13 February 2011
238//
239// Author:
240//
241// John Burkardt
242//
243// Reference:
244//
245// Fred Fritsch, Ralph Carlson,
246// Monotone Piecewise Cubic Interpolation,
247// SIAM Journal on Numerical Analysis,
248// Volume 17, Number 2, April 1980, pages 238-246.
249//
250// Parameters:
251//
252// Input, double X1, F1, D1, the left endpoint, function value
253// and derivative.
254//
255// Input, double X2, F2, D2, the right endpoint, function value
256// and derivative.
257//
258// Input, int N, the number of evaluation points.
259//
260// Input, double X[N], the points at which the Hermite cubic
261// is to be evaluated.
262//
263// Output, double F[N], D[N], S[N], T[N], the value and first
264// three derivatives of the Hermite cubic at X.
265//
266//===============================================================================================================================
268void hermite_cubic_value(double const x1, double const f1, double const d1, double const x2, double const f2, double const d2, int const n, double const* x, double* f, double* d, double* s, double* t)
269{
270 double c2;
271 double c3;
272 double df;
273 double h;
274 int i;
275
276 h = x2 - x1;
277 df = (f2 - f1) / h;
278
279 c2 = -(2.0 * d1 - 3.0 * df + d2) / h;
280 c3 = (d1 - 2.0 * df + d2) / h / h;
281
282 for (i = 0; i < n; i++)
283 {
284 f[i] = f1 + (x[i] - x1) * (d1 + (x[i] - x1) * (c2 + (x[i] - x1) * c3));
285 d[i] = d1 + (x[i] - x1) * (2.0 * c2 + (x[i] - x1) * 3.0 * c3);
286 s[i] = 2.0 * c2 + (x[i] - x1) * 6.0 * c3;
287 t[i] = 6.0 * c3;
288 }
289
290 return;
291}
292
293//===============================================================================================================================
294//
295// Purpose:
296//
297// HERMITE_CUBIC_SPLINE_VALUE evaluates a Hermite cubic spline.
298//
299// Licensing:
300//
301// This code is distributed under the GNU LGPL license.
302//
303// Modified:
304//
305// 13 February 2011
306//
307// Author:
308//
309// John Burkardt
310//
311// Reference:
312//
313// Fred Fritsch, Ralph Carlson,
314// Monotone Piecewise Cubic Interpolation,
315// SIAM Journal on Numerical Analysis,
316// Volume 17, Number 2, April 1980, pages 238-246.
317//
318// Parameters:
319//
320// Input, int NN, the number of data points.
321//
322// Input, double XN[NN], the coordinates of the data points.
323// The entries in XN must be in strictly ascending order.
324//
325// Input, double FN[NN], the function values.
326//
327// Input, double DN[NN], the derivative values.
328//
329// Input, int N, the number of sample points.
330//
331// Input, double X[N], the coordinates of the sample points.
332//
333// Output, double F[N], the function value at the sample points.
334//
335// Output, double D[N], the derivative value at the sample points.
336//
337// Output, double S[N], the second derivative value at the
338// sample points.
339//
340// Output, double T[N], the third derivative value at the
341// sample points.
342//
343//===============================================================================================================================
345//===============================================================================================================================
346void hermite_cubic_spline_value(int const nn, double* const xn, double* const fn, double* const dn, int const n, double* const x, double* f, double* d, double* s, double* t)
347{
348 int left = n / 2;
349
350 for (int i = 0; i < n; i++)
351 {
352 r8vec_bracket3(nn, xn, x[i], &left);
353
354 hermite_cubic_value(xn[left], fn[left], dn[left], xn[left + 1], fn[left + 1], dn[left + 1], 1, x + i, f + i, d + i, s + i, t + i);
355 }
356
357 return;
358}
void hermite_cubic_value(double const x1, double const f1, double const d1, double const x2, double const f2, double const d2, int const n, double const *x, double *f, double *d, double *s, double *t)
This is part of a C++ library from http://people.sc.fsu.edu/~jburkardt/cpp_src/hermite_cubic/hermite_...
void r8vec_bracket3(int const n, double const *t, double const tval, int *left)
This is part of a C++ library from http://people.sc.fsu.edu/~jburkardt/cpp_src/hermite_cubic/hermite_...
void hermite_cubic_spline_value(int const nn, double *const xn, double *const fn, double *const dn, int const n, double *const x, double *f, double *d, double *s, double *t)
This is part of a C++ library from http://people.sc.fsu.edu/~jburkardt/cpp_src/hermite_cubic/hermite_...
Definitions of some routines from the hermite_cubic library.