TMB Documentation  v1.9.11
splines.hpp
Go to the documentation of this file.
1 
5 namespace tmbutils {
6 
7 // Copyright (C) 2013-2015 Kasper Kristensen
8 // License: GPL-2
9 
10 /*
11  * R : A Computer Language for Statistical Data Analysis
12  * Copyright (C) 1995, 1996 Robert Gentleman and Ross Ihaka
13  * Copyright (C) 1998--2001 Robert Gentleman, Ross Ihaka and the
14  * R Development Core Team
15  *
16  * This program is free software; you can redistribute it and/or modify
17  * it under the terms of the GNU General Public License as published by
18  * the Free Software Foundation; either version 2 of the License, or
19  * (at your option) any later version.
20  *
21  * This program is distributed in the hope that it will be useful,
22  * but WITHOUT ANY WARRANTY; without even the implied warranty of
23  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24  * GNU General Public License for more details.
25  *
26  * You should have received a copy of the GNU General Public License
27  * along with this program; if not, a copy is available at
28  * http://www.r-project.org/Licenses/
29  */
30 
31 /* Spline Interpolation
32  * --------------------
33  * C code to perform spline fitting and interpolation.
34  * There is code here for:
35  *
36  * 1. Natural splines.
37  *
38  * 2. Periodic splines
39  *
40  * 3. Splines with end-conditions determined by fitting
41  * cubics in the start and end intervals (Forsythe et al).
42  *
43  *
44  * Computational Techniques
45  * ------------------------
46  * A special LU decomposition for symmetric tridiagonal matrices
47  * is used for all computations, except for periodic splines where
48  * Choleski is more efficient.
49  */
50 
58 template <class Type>
59 class splinefun
60 {
61 private:
62  /* Input to R's "spline_coef" */
63  int method[1];
64  int n[1];
65  Type *x;
66  Type *y;
67  Type *b;
68  Type *c;
69  Type *d;
70  Type *e;
71 
72  /* Not used */
73  //int errno, EDOM;
74 
75  /* Memory management helpers */
76  void clear() {
77  if (*n != 0) {
78  delete [] x;
79  delete [] y;
80  delete [] b;
81  delete [] c;
82  delete [] d;
83  delete [] e;
84  }
85  }
86  void alloc(int n_) {
87  *n = n_;
88  x = new Type[*n];
89  y = new Type[*n];
90  b = new Type[*n];
91  c = new Type[*n];
92  d = new Type[*n];
93  e = new Type[*n];
94  }
95  void copy_from(const splinefun &fun) {
96  *method = fun.method[0];
97  *n = fun.n[0];
98  alloc(*n);
99  for(int i = 0; i < *n; i++) {
100  x[i] = fun.x[i];
101  y[i] = fun.y[i];
102  b[i] = fun.b[i];
103  c[i] = fun.c[i];
104  d[i] = fun.d[i];
105  e[i] = fun.e[i];
106  }
107  }
108 
109 public:
110  /* Construct spline */
111  splinefun() {
112  x = y = b = c = d = e = NULL;
113  *method = *n = 0;
114  };
115  splinefun(const splinefun &fun){
116  copy_from(fun);
117  }
118  splinefun& operator=(const splinefun &x) {
119  if (this != &x) {
120  (*this).clear();
121  (*this).copy_from(x);
122  }
123  return *this;
124  }
125  ~splinefun() {
126  clear();
127  };
138  const vector<Type> &y_,
139  int method_ = 3) {
140  method[0] = method_;
141  n[0] = x_.size();
142  alloc( x_.size() );
143  for(int i=0; i < *n; i++) {
144  x[i] = x_[i];
145  y[i] = y_[i];
146  }
147  spline_coef(method, n, x, y, b, c, d, e);
148  }
149 
151  Type operator()(const Type &x_) {
152  Type u[1];
153  Type v[1];
154  int nu[1];
155  u[0] = x_;
156  nu[0] = 1;
157  spline_eval(method, nu, u, v,
158  n, x, y, b, c, d);
159  return v[0];
160  }
163  vector<Type> y(x.size());
164  for (int i=0; i<x.size(); i++) y[i] = (*this)(x[i]);
165  return y;
166  }
167 
168  /* ------------------------------------------------------------------
169  * The following is copy-pasted from R-package "stats" where "double"
170  * has been replaced by "Type".
171  *-------------------------------------------------------------------
172  */
173 
174  /*
175  * Natural Splines
176  * ---------------
177  * Here the end-conditions are determined by setting the second
178  * derivative of the spline at the end-points to equal to zero.
179  *
180  * There are n-2 unknowns (y[i]'' at x[2], ..., x[n-1]) and n-2
181  * equations to determine them. Either Choleski or Gaussian
182  * elimination could be used.
183  */
184 
185  void natural_spline(int n, Type *x, Type *y, Type *b, Type *c, Type *d)
186  {
187  int nm1, i;
188  Type t;
189 
190  x--; y--; b--; c--; d--;
191 
192  if(n < 2) {
193  //errno = EDOM;
194  return;
195  }
196 
197  if(n < 3) {
198  t = (y[2] - y[1]);
199  b[1] = t / (x[2]-x[1]);
200  b[2] = b[1];
201  c[1] = c[2] = d[1] = d[2] = 0.0;
202  return;
203  }
204 
205  nm1 = n - 1;
206 
207  /* Set up the tridiagonal system */
208  /* b = diagonal, d = offdiagonal, c = right hand side */
209 
210  d[1] = x[2] - x[1];
211  c[2] = (y[2] - y[1])/d[1];
212  for( i=2 ; i<n ; i++) {
213  d[i] = x[i+1] - x[i];
214  b[i] = 2.0 * (d[i-1] + d[i]);
215  c[i+1] = (y[i+1] - y[i])/d[i];
216  c[i] = c[i+1] - c[i];
217  }
218 
219  /* Gaussian elimination */
220 
221  for(i=3 ; i<n ; i++) {
222  t = d[i-1]/b[i-1];
223  b[i] = b[i] - t*d[i-1];
224  c[i] = c[i] - t*c[i-1];
225  }
226 
227  /* Backward substitution */
228 
229  c[nm1] = c[nm1]/b[nm1];
230  for(i=n-2 ; i>1 ; i--)
231  c[i] = (c[i]-d[i]*c[i+1])/b[i];
232 
233  /* End conditions */
234 
235  c[1] = c[n] = 0.0;
236 
237  /* Get cubic coefficients */
238 
239  b[1] = (y[2] - y[1])/d[1] - d[i] * c[2];
240  c[1] = 0.0;
241  d[1] = c[2]/d[1];
242  b[n] = (y[n] - y[nm1])/d[nm1] + d[nm1] * c[nm1];
243  for(i=2 ; i<n ; i++) {
244  b[i] = (y[i+1]-y[i])/d[i] - d[i]*(c[i+1]+2.0*c[i]);
245  d[i] = (c[i+1]-c[i])/d[i];
246  c[i] = 3.0*c[i];
247  }
248  c[n] = 0.0;
249  d[n] = 0.0;
250 
251  return;
252  }
253 
254  /*
255  * Splines a la Forsythe Malcolm and Moler
256  * ---------------------------------------
257  * In this case the end-conditions are determined by fitting
258  * cubic polynomials to the first and last 4 points and matching
259  * the third derivitives of the spline at the end-points to the
260  * third derivatives of these cubics at the end-points.
261  */
262 
263  void fmm_spline(int n, Type *x, Type *y, Type *b, Type *c, Type *d)
264  {
265  int nm1, i;
266  Type t;
267 
268  /* Adjustment for 1-based arrays */
269 
270  x--; y--; b--; c--; d--;
271 
272  if(n < 2) {
273  //errno = EDOM;
274  return;
275  }
276 
277  if(n < 3) {
278  t = (y[2] - y[1]);
279  b[1] = t / (x[2]-x[1]);
280  b[2] = b[1];
281  c[1] = c[2] = d[1] = d[2] = 0.0;
282  return;
283  }
284 
285  nm1 = n - 1;
286 
287  /* Set up tridiagonal system */
288  /* b = diagonal, d = offdiagonal, c = right hand side */
289 
290  d[1] = x[2] - x[1];
291  c[2] = (y[2] - y[1])/d[1];/* = +/- Inf for x[1]=x[2] -- problem? */
292  for(i=2 ; i<n ; i++) {
293  d[i] = x[i+1] - x[i];
294  b[i] = 2.0 * (d[i-1] + d[i]);
295  c[i+1] = (y[i+1] - y[i])/d[i];
296  c[i] = c[i+1] - c[i];
297  }
298 
299  /* End conditions. */
300  /* Third derivatives at x[0] and x[n-1] obtained */
301  /* from divided differences */
302 
303  b[1] = -d[1];
304  b[n] = -d[nm1];
305  c[1] = c[n] = 0.0;
306  if(n > 3) {
307  c[1] = c[3]/(x[4]-x[2]) - c[2]/(x[3]-x[1]);
308  c[n] = c[nm1]/(x[n] - x[n-2]) - c[n-2]/(x[nm1]-x[n-3]);
309  c[1] = c[1]*d[1]*d[1]/(x[4]-x[1]);
310  c[n] = -c[n]*d[nm1]*d[nm1]/(x[n]-x[n-3]);
311  }
312 
313  /* Gaussian elimination */
314 
315  for(i=2 ; i<=n ; i++) {
316  t = d[i-1]/b[i-1];
317  b[i] = b[i] - t*d[i-1];
318  c[i] = c[i] - t*c[i-1];
319  }
320 
321  /* Backward substitution */
322 
323  c[n] = c[n]/b[n];
324  for(i=nm1 ; i>=1 ; i--)
325  c[i] = (c[i]-d[i]*c[i+1])/b[i];
326 
327  /* c[i] is now the sigma[i-1] of the text */
328  /* Compute polynomial coefficients */
329 
330  b[n] = (y[n] - y[n-1])/d[n-1] + d[n-1]*(c[n-1]+ 2.0*c[n]);
331  for(i=1 ; i<=nm1 ; i++) {
332  b[i] = (y[i+1]-y[i])/d[i] - d[i]*(c[i+1]+2.0*c[i]);
333  d[i] = (c[i+1]-c[i])/d[i];
334  c[i] = 3.0*c[i];
335  }
336  c[n] = 3.0*c[n];
337  d[n] = d[nm1];
338  return;
339  }
340 
341 
342  /*
343  * Periodic Spline
344  * ---------------
345  * The end conditions here match spline (and its derivatives)
346  * at x[1] and x[n].
347  *
348  * Note: There is an explicit check that the user has supplied
349  * data with y[1] equal to y[n].
350  */
351 
352  void periodic_spline(int n, Type *x, Type *y,
353  Type *b, Type *c, Type *d, Type *e)
354  {
355  Type s;
356  int i, nm1;
357 
358  /* Adjustment for 1-based arrays */
359 
360  x--; y--; b--; c--; d--; e--;
361 
362  if(n < 2 || y[1] != y[n]) {
363  errno = EDOM;
364  return;
365  }
366 
367  if(n == 2) {
368  b[1] = b[2] = c[1] = c[2] = d[1] = d[2] = 0.0;
369  return;
370  } else if (n == 3) {
371  b[1] = b[2] = b[3] = -(y[1] - y[2])*(x[1] - 2*x[2] + x[3])/(x[3]-x[2])/(x[2]-x[1]);
372  c[1] = -3*(y[1]-y[2])/(x[3]-x[2])/(x[2]-x[1]);
373  c[2] = -c[1];
374  c[3] = c[1];
375  d[1] = -2*c[1]/3/(x[2]-x[1]);
376  d[2] = -d[1]*(x[2]-x[1])/(x[3]-x[2]);
377  d[3] = d[1];
378  return;
379  }
380 
381  /* else --------- n >= 4 --------- */
382  nm1 = n-1;
383 
384  /* Set up the matrix system */
385  /* A = diagonal B = off-diagonal C = rhs */
386 
387 #define A b
388 #define B d
389 #define C c
390 
391  B[1] = x[2] - x[1];
392  B[nm1]= x[n] - x[nm1];
393  A[1] = 2.0 * (B[1] + B[nm1]);
394  C[1] = (y[2] - y[1])/B[1] - (y[n] - y[nm1])/B[nm1];
395 
396  for(i = 2; i < n; i++) {
397  B[i] = x[i+1] - x[i];
398  A[i] = 2.0 * (B[i] + B[i-1]);
399  C[i] = (y[i+1] - y[i])/B[i] - (y[i] - y[i-1])/B[i-1];
400  }
401 
402  /* Choleski decomposition */
403 
404 #define L b
405 #define M d
406 #define E e
407 
408  L[1] = sqrt(A[1]);
409  E[1] = (x[n] - x[nm1])/L[1];
410  s = 0.0;
411  for(i = 1; i <= nm1 - 2; i++) {
412  M[i] = B[i]/L[i];
413  if(i != 1) E[i] = -E[i-1] * M[i-1] / L[i];
414  L[i+1] = sqrt(A[i+1]-M[i]*M[i]);
415  s = s + E[i] * E[i];
416  }
417  M[nm1-1] = (B[nm1-1] - E[nm1-2] * M[nm1-2])/L[nm1-1];
418  L[nm1] = sqrt(A[nm1] - M[nm1-1]*M[nm1-1] - s);
419 
420  /* Forward Elimination */
421 
422 #define Y c
423 #define D c
424 
425  Y[1] = D[1]/L[1];
426  s = 0.0;
427  for(i=2 ; i<=nm1-1 ; i++) {
428  Y[i] = (D[i] - M[i-1]*Y[i-1])/L[i];
429  s = s + E[i-1] * Y[i-1];
430  }
431  Y[nm1] = (D[nm1] - M[nm1-1] * Y[nm1-1] - s) / L[nm1];
432 
433 #define X c
434 
435  X[nm1] = Y[nm1]/L[nm1];
436  X[nm1-1] = (Y[nm1-1] - M[nm1-1] * X[nm1])/L[nm1-1];
437  for(i=nm1-2 ; i>=1 ; i--)
438  X[i] = (Y[i] - M[i] * X[i+1] - E[i] * X[nm1])/L[i];
439 
440  /* Wrap around */
441 
442  X[n] = X[1];
443 
444  /* Compute polynomial coefficients */
445 
446  for(i=1 ; i<=nm1 ; i++) {
447  s = x[i+1] - x[i];
448  b[i] = (y[i+1]-y[i])/s - s*(c[i+1]+2.0*c[i]);
449  d[i] = (c[i+1]-c[i])/s;
450  c[i] = 3.0*c[i];
451  }
452  b[n] = b[1];
453  c[n] = c[1];
454  d[n] = d[1];
455  return;
456  }
457 #undef A
458 #undef B
459 #undef C
460 #undef L
461 #undef M
462 #undef E
463 #undef Y
464 #undef D
465 #undef X
466 
467  void spline_coef(int *method, int *n, Type *x, Type *y,
468  Type *b, Type *c, Type *d, Type *e)
469  {
470  switch(*method) {
471  case 1:
472  periodic_spline(*n, x, y, b, c, d, e); break;
473 
474  case 2:
475  natural_spline(*n, x, y, b, c, d); break;
476 
477  case 3:
478  fmm_spline(*n, x, y, b, c, d); break;
479  }
480  }
481 
482  void spline_eval(int *method, int *nu, Type *u, Type *v,
483  int *n, Type *x, Type *y, Type *b, Type *c, Type *d)
484  {
485  /* Evaluate v[l] := spline(u[l], ...), l = 1,..,nu, i.e. 0:(nu-1)
486  * Nodes x[i], coef (y[i]; b[i],c[i],d[i]); i = 1,..,n , i.e. 0:(*n-1)
487  */
488  const int n_1 = *n - 1;
489  int i, j, k, l;
490  Type ul, dx, tmp;
491 
492  if(*method == 1 && *n > 1) { /* periodic */
493  dx = x[n_1] - x[0];
494  for(l = 0; l < *nu; l++) {
495 
496  /* WARNING - "fmod(AD<double>,AD<double>)" is not defined */
497  v[l] = fmod(asDouble(u[l]-x[0]), asDouble(dx));
498  if(v[l] < 0.0) v[l] += dx;
499  v[l] += x[0];
500  }
501  }
502  else {
503  for(l = 0; l < *nu; l++)
504  v[l] = u[l];
505  }
506 
507  i = 0;
508  for(l = 0; l < *nu; l++) {
509  ul = v[l];
510  if(ul < x[i] || (i < n_1 && x[i+1] < ul)) {
511  /* reset i such that x[i] <= ul <= x[i+1] : */
512  i = 0;
513  j = *n;
514  do {
515  k = (i+j)/2;
516  if(ul < x[k]) j = k;
517  else i = k;
518  }
519  while(j > i+1);
520  }
521  dx = ul - x[i];
522  /* for natural splines extrapolate linearly left */
523  tmp = (*method == 2 && ul < x[0]) ? 0.0 : d[i];
524 
525  v[l] = y[i] + dx*(b[i] + dx*(c[i] + dx*tmp));
526  }
527  }
528 
529 };
530 
531 }
splinefun(const vector< Type > &x_, const vector< Type > &y_, int method_=3)
Construct spline function object.
Definition: splines.hpp:137
Spline Interpolation.
Definition: splines.hpp:59
Vector class used by TMB.
Definition: tmbutils.hpp:18
Type operator()(const Type &x_)
Evaluate spline - scalar argument case.
Definition: splines.hpp:151
Utility functions for TMB (automatically included)
Definition: concat.hpp:5
License: GPL v2