summaryrefslogtreecommitdiff
path: root/gsl-1.9/ode-initval/rk2simp.c
blob: 61248acc431351c973c4bd9a6ee4930f83bc80f6 (plain)
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
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
/* ode-initval/rk2simp.c
 * 
 * Copyright (C) 2004 Tuomo Keskitalo
 * 
 * This program 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 2 of the License, or (at
 * your option) any later version.
 * 
 * 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.
 * 
 * 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., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
 */

/* Runge-Kutta 2, Gaussian implicit. Also known as implicit midpoint rule.

   Non-linear equations solved by linearization, LU-decomposition
   and matrix inversion. For reference, see eg.

   Ascher, U.M., Petzold, L.R., Computer methods for ordinary
   differential and differential-algebraic equations, SIAM,
   Philadelphia, 1998.
 */

#include <config.h>
#include <stdlib.h>
#include <string.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_errno.h>
#include <gsl/gsl_odeiv.h>
#include <gsl/gsl_linalg.h>

#include "odeiv_util.h"

typedef struct
{
  double *Y1;
  double *y0;
  double *y0_orig;
  double *ytmp;
  double *dfdy;                 /* Jacobian */
  double *dfdt;                 /* time derivatives, not used */
  double *y_onestep;
  gsl_permutation *p;
}
rk2simp_state_t;

static void *
rk2simp_alloc (size_t dim)
{
  rk2simp_state_t *state =
    (rk2simp_state_t *) malloc (sizeof (rk2simp_state_t));

  if (state == 0)
    {
      GSL_ERROR_NULL ("failed to allocate space for rk2simp_state",
                      GSL_ENOMEM);
    }

  state->Y1 = (double *) malloc (dim * sizeof (double));

  if (state->Y1 == 0)
    {
      free (state);
      GSL_ERROR_NULL ("failed to allocate space for Y1", GSL_ENOMEM);
    }

  state->y0 = (double *) malloc (dim * sizeof (double));

  if (state->y0 == 0)
    {
      free (state->Y1);
      free (state);
      GSL_ERROR_NULL ("failed to allocate space for y0", GSL_ENOMEM);
    }

  state->y0_orig = (double *) malloc (dim * sizeof (double));

  if (state->y0_orig == 0)
    {
      free (state->Y1);
      free (state->y0);
      free (state);
      GSL_ERROR_NULL ("failed to allocate space for y0_orig", GSL_ENOMEM);
    }

  state->ytmp = (double *) malloc (dim * sizeof (double));

  if (state->ytmp == 0)
    {
      free (state->Y1);
      free (state->y0);
      free (state->y0_orig);
      free (state);
      GSL_ERROR_NULL ("failed to allocate space for ytmp", GSL_ENOMEM);
    }

  state->dfdy = (double *) malloc (dim * dim * sizeof (double));

  if (state->dfdy == 0)
    {
      free (state->Y1);
      free (state->y0);
      free (state->y0_orig);
      free (state->ytmp);
      free (state);
      GSL_ERROR_NULL ("failed to allocate space for dfdy", GSL_ENOMEM);
    }

  state->dfdt = (double *) malloc (dim * sizeof (double));

  if (state->dfdt == 0)
    {
      free (state->Y1);
      free (state->y0);
      free (state->y0_orig);
      free (state->ytmp);
      free (state->dfdy);
      free (state);
      GSL_ERROR_NULL ("failed to allocate space for dfdt", GSL_ENOMEM);
    }

  state->y_onestep = (double *) malloc (dim * sizeof (double));

  if (state->y_onestep == 0)
    {
      free (state->Y1);
      free (state->y0);
      free (state->y0_orig);
      free (state->ytmp);
      free (state->dfdy);
      free (state->dfdt);
      free (state);
      GSL_ERROR_NULL ("failed to allocate space for y_onestep", GSL_ENOMEM);
    }

  state->p = gsl_permutation_alloc (dim);

  if (state->p == 0)
    {
      free (state->Y1);
      free (state->y0);
      free (state->y0_orig);
      free (state->ytmp);
      free (state->dfdy);
      free (state->dfdt);
      free (state);
      GSL_ERROR_NULL ("failed to allocate space for p", GSL_ENOMEM);
    }

  return state;
}


static int
rk2simp_step (double *y, rk2simp_state_t * state,
              const double h, const double t,
              const size_t dim, const gsl_odeiv_system * sys)
{
  /* Makes a Runge-Kutta 2nd order semi-implicit advance with step size h.
     y0 is initial values of variables y. 

     The linearized semi-implicit equations to calculate are:

     Y1 = y0 + h/2 * (1 - h/2 * df/dy)^(-1) * f(t + h/2, y0)

     y = y0 + h * f(t + h/2, Y1)
   */

  const double *y0 = state->y0;
  double *Y1 = state->Y1;
  double *ytmp = state->ytmp;

  size_t i;
  int s, ps;

  gsl_matrix_view J = gsl_matrix_view_array (state->dfdy, dim, dim);

  /* First solve Y1. 
     Calculate the inverse matrix (1 - h/2 * df/dy)^-1 
   */

  /* Create matrix to J */

  s = GSL_ODEIV_JA_EVAL (sys, t, y0, state->dfdy, state->dfdt);

  if (s != GSL_SUCCESS)
    {
      return s;
    }

  gsl_matrix_scale (&J.matrix, -h / 2.0);
  gsl_matrix_add_diagonal(&J.matrix, 1.0);

  /* Invert it by LU-decomposition to invmat */

  s += gsl_linalg_LU_decomp (&J.matrix, state->p, &ps);

  if (s != GSL_SUCCESS)
    {
      return GSL_EFAILED;
    }

  /* Evaluate f(t + h/2, y0) */

  s = GSL_ODEIV_FN_EVAL (sys, t + 0.5 * h, y0, ytmp);

  if (s != GSL_SUCCESS)
    {
      return s;
    }

  /* Calculate Y1 = y0 + h/2 * ((1-h/2 * df/dy)^-1) ytmp */

  {
    gsl_vector_const_view y0_view = gsl_vector_const_view_array(y0, dim);
    gsl_vector_view ytmp_view = gsl_vector_view_array(ytmp, dim);
    gsl_vector_view Y1_view = gsl_vector_view_array(Y1, dim);

    s = gsl_linalg_LU_solve (&J.matrix, state->p, 
                             &ytmp_view.vector, &Y1_view.vector);
      
    gsl_vector_scale (&Y1_view.vector, 0.5 * h);
    gsl_vector_add (&Y1_view.vector, &y0_view.vector);
  }

  /* And finally evaluation of f(t + h/2, Y1) and calculation of y */

  s = GSL_ODEIV_FN_EVAL (sys, t + 0.5 * h, Y1, ytmp);

  if (s != GSL_SUCCESS)
    {
      return s;
    }

  for (i = 0; i < dim; i++)
    {
      y[i] = y0[i] + h * ytmp[i];
    }

  return s;
}

static int
rk2simp_apply (void *vstate, size_t dim, double t, double h,
               double y[], double yerr[], const double dydt_in[],
               double dydt_out[], const gsl_odeiv_system * sys)
{
  rk2simp_state_t *state = (rk2simp_state_t *) vstate;

  size_t i;

  double *y0 = state->y0;
  double *y0_orig = state->y0_orig;
  double *y_onestep = state->y_onestep;

  /* Error estimation is done by step doubling procedure */

  DBL_MEMCPY (y0, y, dim);

  /* Save initial values in case of failure */
  DBL_MEMCPY (y0_orig, y, dim);

  /* First traverse h with one step (save to y_onestep) */
  DBL_MEMCPY (y_onestep, y, dim);

  {
    int s = rk2simp_step (y_onestep, state, h, t, dim, sys);

    if (s != GSL_SUCCESS)
      {
        return s;
      }
  }

  /* Then with two steps with half step length (save to y) */

  {
    int s = rk2simp_step (y, state, h / 2.0, t, dim, sys);

    if (s != GSL_SUCCESS)
      {
        /* Restore original y vector */
        DBL_MEMCPY (y, y0_orig, dim);
        return s;
      }
  }

  DBL_MEMCPY (y0, y, dim);

  {
    int s = rk2simp_step (y, state, h / 2.0, t + h / 2.0, dim, sys);

    if (s != GSL_SUCCESS)
      {
        /* Restore original y vector */
        DBL_MEMCPY (y, y0_orig, dim);
        return s;
      }
  }

  /* Derivatives at output */

  if (dydt_out != NULL)
    {
      int s = GSL_ODEIV_FN_EVAL (sys, t + h, y, dydt_out);

      if (s != GSL_SUCCESS)
        {
          /* Restore original y vector */
          DBL_MEMCPY (y, y0_orig, dim);
          return s;
        }
    }

  /* Error estimation */

  for (i = 0; i < dim; i++)
    {
      yerr[i] = 4.0 * (y[i] - y_onestep[i]) / 3.0;
    }

  return GSL_SUCCESS;
}


static int
rk2simp_reset (void *vstate, size_t dim)
{
  rk2simp_state_t *state = (rk2simp_state_t *) vstate;

  DBL_ZERO_MEMSET (state->Y1, dim);
  DBL_ZERO_MEMSET (state->y0, dim);
  DBL_ZERO_MEMSET (state->y0_orig, dim);
  DBL_ZERO_MEMSET (state->ytmp, dim);
  DBL_ZERO_MEMSET (state->dfdt, dim * dim);
  DBL_ZERO_MEMSET (state->dfdt, dim);
  DBL_ZERO_MEMSET (state->y_onestep, dim);

  return GSL_SUCCESS;
}

static unsigned int
rk2simp_order (void *vstate)
{
  rk2simp_state_t *state = (rk2simp_state_t *) vstate;
  state = 0;                    /* prevent warnings about unused parameters */
  return 2;
}

static void
rk2simp_free (void *vstate)
{
  rk2simp_state_t *state = (rk2simp_state_t *) vstate;
  free (state->Y1);
  free (state->y0);
  free (state->y0_orig);
  free (state->ytmp);
  free (state->dfdy);
  free (state->dfdt);
  free (state->y_onestep);
  gsl_permutation_free (state->p);
  free (state);
}

static const gsl_odeiv_step_type rk2simp_type = {
  "rk2simp",                    /* name */
  0,                            /* can use dydt_in? */
  1,                            /* gives exact dydt_out? */
  &rk2simp_alloc,
  &rk2simp_apply,
  &rk2simp_reset,
  &rk2simp_order,
  &rk2simp_free
};

const gsl_odeiv_step_type *gsl_odeiv_step_rk2simp = &rk2simp_type;