1 |
johnpye |
669 |
/* ASCEND modelling environment |
2 |
|
|
Copyright 1997, Carnegie Mellon University |
3 |
|
|
Copyright (C) 2006 Carnegie Mellon University |
4 |
|
|
|
5 |
|
|
This program is free software; you can redistribute it and/or modify |
6 |
|
|
it under the terms of the GNU General Public License as published by |
7 |
|
|
the Free Software Foundation; either version 2, or (at your option) |
8 |
|
|
any later version. |
9 |
|
|
|
10 |
|
|
This program is distributed in the hope that it will be useful, |
11 |
|
|
but WITHOUT ANY WARRANTY; without even the implied warranty of |
12 |
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
13 |
|
|
GNU General Public License for more details. |
14 |
|
|
|
15 |
|
|
You should have received a copy of the GNU General Public License |
16 |
|
|
along with this program; if not, write to the Free Software |
17 |
|
|
Foundation, Inc., 59 Temple Place - Suite 330, |
18 |
|
|
Boston, MA 02111-1307, USA. |
19 |
|
|
*//** |
20 |
|
|
@file |
21 |
|
|
LSODE integrator. |
22 |
|
|
|
23 |
|
|
(old implementation notes:) |
24 |
|
|
|
25 |
|
|
As fortran io is unreliably portable (vc5+digital fortran) |
26 |
|
|
we have converted xerrwv to xascwv provided here. |
27 |
johnpye |
741 |
|
28 |
johnpye |
669 |
The lsode interface variable t is actually an array of |
29 |
|
|
2 doubles rather than just 1. The first is the the one |
30 |
|
|
used by lsode. The second is used by LSODE_FEX to tell |
31 |
|
|
what the last time it was called was. This is so the |
32 |
|
|
C driver can tell if it needs to resolve d to compute |
33 |
|
|
observation variables. If x[0]==x[1] we save ourselves |
34 |
|
|
a solve. |
35 |
|
|
|
36 |
johnpye |
741 |
@NOTE The above doesn't work since lsode doesn't use the same t internally |
37 |
johnpye |
669 |
that we hand it. |
38 |
|
|
|
39 |
|
|
*//* |
40 |
|
|
by Kirk Abbott and Ben Allan |
41 |
|
|
Created: 1/94 |
42 |
|
|
Version: $Revision: 1.29 $ |
43 |
|
|
Version control file: $RCSfile: Lsode.c,v $ |
44 |
|
|
Date last modified: $Date: 2000/01/25 02:26:31 $ |
45 |
|
|
Last modified by: $Author: ballan $ |
46 |
|
|
*/ |
47 |
|
|
|
48 |
|
|
#ifndef NO_SIGNAL_TRAPS |
49 |
|
|
#include <signal.h> |
50 |
|
|
#include <setjmp.h> |
51 |
|
|
#endif /* NO_SIGNAL_TRAPS */ |
52 |
|
|
|
53 |
|
|
#include <utilities/ascConfig.h> |
54 |
|
|
#include <utilities/error.h> |
55 |
|
|
#include <compiler/instance_enum.h> |
56 |
|
|
#include <utilities/ascSignal.h> |
57 |
|
|
#include <utilities/ascMalloc.h> |
58 |
johnpye |
888 |
#include <utilities/ascPanic.h> |
59 |
johnpye |
669 |
|
60 |
|
|
#include "slv_types.h" |
61 |
|
|
#include "mtx.h" |
62 |
|
|
#include "rel.h" |
63 |
|
|
#include "var.h" |
64 |
|
|
#include "discrete.h" |
65 |
|
|
#include "conditional.h" |
66 |
|
|
#include "bnd.h" |
67 |
|
|
#include "logrel.h" |
68 |
|
|
#include "slv_common.h" |
69 |
|
|
#include "linsol.h" |
70 |
|
|
#include "linsolqr.h" |
71 |
|
|
#include "slv_client.h" |
72 |
johnpye |
1055 |
#include "slv_stdcalls.h" |
73 |
johnpye |
669 |
|
74 |
|
|
#include "integrator.h" |
75 |
|
|
#include "lsode.h" |
76 |
|
|
|
77 |
johnpye |
977 |
const IntegratorInternals integrator_lsode_internals = { |
78 |
|
|
integrator_lsode_create |
79 |
|
|
,integrator_lsode_params_default |
80 |
|
|
,integrator_analyse_ode /* note, this routine is back in integrator.c */ |
81 |
|
|
,integrator_lsode_solve |
82 |
|
|
,integrator_lsode_free |
83 |
|
|
,INTEG_LSODE |
84 |
|
|
,"LSODE" |
85 |
|
|
}; |
86 |
|
|
|
87 |
johnpye |
669 |
/* |
88 |
|
|
#include "Sensitivity.h" |
89 |
|
|
*//* see the packages dir */ |
90 |
|
|
|
91 |
|
|
/* |
92 |
|
|
* NOUNDERBARS --> FORTRAN compiler naming convention for subroutine |
93 |
|
|
* is wierd. WIN32/CRAY is treated as special case |
94 |
|
|
*/ |
95 |
|
|
#ifdef APOLLO |
96 |
|
|
#define NOUNDERBARS TRUE |
97 |
|
|
#endif |
98 |
|
|
#ifdef _HPUX_SOURCE |
99 |
|
|
#define NOUNDERBARS TRUE |
100 |
|
|
#endif |
101 |
|
|
/* AIX xlf will not suffix an underbar on a symbol |
102 |
|
|
* unless xlf is given the ``-qextname'' option |
103 |
|
|
*/ |
104 |
|
|
#ifdef _AIX |
105 |
|
|
#define NOUNDERBARS TRUE |
106 |
|
|
#endif |
107 |
|
|
|
108 |
|
|
#ifdef NOUNDERBARS |
109 |
|
|
#define LSODE lsode |
110 |
|
|
#define LSODE_JEX jex |
111 |
|
|
#define LSODE_FEX fex |
112 |
|
|
#define GETCOMMON get_lsode_common |
113 |
|
|
#define XASCWV xascwv |
114 |
|
|
#else |
115 |
|
|
/* sun, __alpha, __sgi, ... */ |
116 |
|
|
#define LSODE lsode_ |
117 |
|
|
#define LSODE_JEX jex_ |
118 |
|
|
#define LSODE_FEX fex_ |
119 |
|
|
#define GETCOMMON get_lsode_common_ |
120 |
|
|
#define XASCWV xascwv_ |
121 |
|
|
#endif |
122 |
|
|
|
123 |
|
|
#if defined(CRAY) || (defined(__WIN32__) && !defined(__MINGW32_VERSION)) |
124 |
|
|
#undef LSODE |
125 |
|
|
#undef LSODE_JEX |
126 |
|
|
#undef LSODE_FEX |
127 |
|
|
#undef GETCOMMON |
128 |
|
|
#undef XASCWV |
129 |
|
|
#define XASCWV XASCWV |
130 |
|
|
#define LSODE LSODE |
131 |
|
|
#define LSODE_JEX JEX |
132 |
|
|
#define LSODE_FEX FEX |
133 |
|
|
#define GETCOMMON GET_LSODE_COMMON |
134 |
|
|
#endif |
135 |
|
|
|
136 |
|
|
#define DOTIME FALSE |
137 |
|
|
|
138 |
|
|
/* definitions of lsode supported children of atoms, etc */ |
139 |
|
|
/********************************************************************/ |
140 |
|
|
/* solver_var children expected for state variables */ |
141 |
|
|
static symchar *g_symbols[2]; |
142 |
|
|
#define STATERTOL g_symbols[0] |
143 |
|
|
#define STATEATOL g_symbols[1] |
144 |
|
|
static |
145 |
|
|
void InitTolNames(void) |
146 |
|
|
{ |
147 |
|
|
STATERTOL = AddSymbol("ode_rtol"); |
148 |
|
|
STATEATOL = AddSymbol("ode_atol"); |
149 |
|
|
} |
150 |
|
|
|
151 |
|
|
/** |
152 |
johnpye |
741 |
Because LSODE doesn't seem to make an allowance for 'client data' we |
153 |
johnpye |
669 |
have to store this as a 'local global' and fish it out when we're in the |
154 |
|
|
callbacks. |
155 |
|
|
*/ |
156 |
|
|
IntegratorSystem *l_lsode_blsys; |
157 |
|
|
|
158 |
|
|
enum Lsode_enum { |
159 |
|
|
lsode_none, /* true on first call */ |
160 |
|
|
lsode_function, lsode_derivative, /* functions or gradients done */ |
161 |
|
|
lsode_sparse, lsode_dense, /* what type of backend should we */ |
162 |
|
|
lsode_band, /* use for the integrator */ |
163 |
|
|
lsode_ok, lsode_nok /* bad return from func or grad */ |
164 |
|
|
}; |
165 |
|
|
|
166 |
|
|
static struct Lsode_Data { |
167 |
|
|
enum Lsode_enum lastcall; /* type of last call; func or grad */ |
168 |
|
|
enum Lsode_enum status; /* solve status */ |
169 |
|
|
int partitioned; /* partioned func evals or not */ |
170 |
|
|
} lsodesys = {lsode_none, lsode_ok, 1}; |
171 |
|
|
|
172 |
|
|
|
173 |
|
|
/*-------------------------- |
174 |
|
|
Data space for use by LSODE |
175 |
|
|
*/ |
176 |
johnpye |
1050 |
typedef struct IntegratorLsodeDataStruct{ |
177 |
ben.allan |
704 |
long n_eqns; /**< dimension of state vector */ |
178 |
|
|
int *input_indices; /**< vector of state vars indexes */ |
179 |
|
|
int *output_indices; /**< vector of derivative var indexes */ |
180 |
|
|
struct var_variable **y_vars; /**< NULL terminated list of states vars */ |
181 |
|
|
struct var_variable **ydot_vars; /**< NULL terminated list of derivative vars*/ |
182 |
|
|
struct rel_relation **rlist; /**< NULL terminated list of relevant rels |
183 |
|
|
to be differentiated */ |
184 |
|
|
double **dydx_dx; /**< change in derivatives wrt states |
185 |
|
|
I prefer to call this: d(ydot)/dy */ |
186 |
johnpye |
669 |
} IntegratorLsodeData; |
187 |
|
|
|
188 |
|
|
|
189 |
|
|
/*---------------------------- |
190 |
|
|
Function types that LSODE wants to use |
191 |
|
|
*/ |
192 |
|
|
|
193 |
|
|
/** |
194 |
|
|
Type of function used to evaluate derivative system. |
195 |
|
|
*/ |
196 |
johnpye |
888 |
typedef void LsodeEvalFn(int *, double *, double *, double *); |
197 |
johnpye |
669 |
|
198 |
|
|
/** |
199 |
|
|
Type of function used to evaluate jacobian system. |
200 |
|
|
*/ |
201 |
johnpye |
888 |
typedef void LsodeJacobianFn(int *, double *, double *, int *, int *, double *, int *); |
202 |
johnpye |
669 |
|
203 |
|
|
/*---------------------------- |
204 |
|
|
forward declarations |
205 |
|
|
*/ |
206 |
|
|
|
207 |
|
|
int integrator_lsode_setup_diffs(IntegratorSystem *blsys); |
208 |
johnpye |
888 |
static double **lsode_densematrix_create(int nrows, int ncols); |
209 |
|
|
static void lsode_densematrix_destroy(double **matrix,int nrows); |
210 |
johnpye |
669 |
|
211 |
|
|
/** |
212 |
|
|
void LSODE(&fex, &neq, y, &x, &xend, &itol, reltol, abtol, &itask, |
213 |
|
|
&istate, &iopt ,rwork, &lrw, iwork, &liw, &jex, &mf); |
214 |
|
|
|
215 |
|
|
This is a prototype for the *fortran* LSODE function. |
216 |
|
|
|
217 |
|
|
No 'extern' here, so we want linker to complain if no static linkage. |
218 |
|
|
*/ |
219 |
|
|
void LSODE(LsodeEvalFn*,int *neq ,double *y ,double *x |
220 |
|
|
,double *xend |
221 |
|
|
,int *itol ,double *reltol ,double *abtol |
222 |
|
|
,int *itask ,int *istate ,int *iopt |
223 |
|
|
,double *rwork ,int *lrw |
224 |
|
|
,int *iwork ,int *liw |
225 |
|
|
,LsodeJacobianFn *jex ,int *mf |
226 |
|
|
); |
227 |
|
|
|
228 |
|
|
/*------------------------------------------------------ |
229 |
|
|
Memory allocation/free |
230 |
|
|
*/ |
231 |
|
|
|
232 |
|
|
void integrator_lsode_create(IntegratorSystem *blsys){ |
233 |
|
|
IntegratorLsodeData *d; |
234 |
|
|
d = ASC_NEW_CLEAR(IntegratorLsodeData); |
235 |
|
|
d->n_eqns=0; |
236 |
|
|
d->input_indices=NULL; |
237 |
|
|
d->output_indices=NULL; |
238 |
|
|
d->y_vars=NULL; |
239 |
|
|
d->ydot_vars=NULL; |
240 |
|
|
d->rlist=NULL; |
241 |
|
|
d->dydx_dx=NULL; |
242 |
|
|
blsys->enginedata=(void*)d; |
243 |
johnpye |
962 |
integrator_lsode_params_default(blsys); |
244 |
|
|
|
245 |
johnpye |
669 |
} |
246 |
|
|
|
247 |
|
|
/** |
248 |
johnpye |
888 |
Cleanup the data struct that belongs to LSODE |
249 |
johnpye |
669 |
*/ |
250 |
|
|
void integrator_lsode_free(void *enginedata){ |
251 |
|
|
IntegratorLsodeData d; |
252 |
|
|
d = *((IntegratorLsodeData *)enginedata); |
253 |
|
|
|
254 |
|
|
if(d.input_indices)ASC_FREE(d.input_indices); |
255 |
|
|
d.input_indices = NULL; |
256 |
|
|
|
257 |
|
|
if(d.output_indices)ASC_FREE(d.output_indices); |
258 |
|
|
d.output_indices = NULL; |
259 |
|
|
|
260 |
|
|
if(d.y_vars)ASC_FREE(d.y_vars); |
261 |
|
|
d.y_vars = NULL; |
262 |
|
|
|
263 |
|
|
if(d.ydot_vars)ASC_FREE(d.ydot_vars); |
264 |
|
|
d.ydot_vars = NULL; |
265 |
|
|
|
266 |
|
|
if(d.rlist)ASC_FREE(d.rlist); |
267 |
|
|
d.rlist = NULL; |
268 |
|
|
|
269 |
johnpye |
888 |
if(d.dydx_dx)lsode_densematrix_destroy(d.dydx_dx, d.n_eqns); |
270 |
johnpye |
669 |
d.dydx_dx = NULL; |
271 |
|
|
|
272 |
|
|
d.n_eqns = 0L; |
273 |
|
|
} |
274 |
|
|
|
275 |
johnpye |
942 |
/*------------------------------------------------------------------------------ |
276 |
|
|
PARAMETERS |
277 |
|
|
*/ |
278 |
|
|
|
279 |
|
|
enum ida_parameters{ |
280 |
|
|
LSODE_PARAM_TIMING |
281 |
johnpye |
963 |
,LSODE_PARAM_RTOLVECT |
282 |
|
|
,LSODE_PARAM_RTOL |
283 |
|
|
,LSODE_PARAM_ATOLVECT |
284 |
|
|
,LSODE_PARAM_ATOL |
285 |
johnpye |
942 |
,LSODE_PARAMS_SIZE |
286 |
|
|
}; |
287 |
|
|
|
288 |
|
|
/** |
289 |
|
|
Here the full set of parameters is defined, along with upper/lower bounds, |
290 |
|
|
etc. The values are stuck into the blsys->params structure. |
291 |
|
|
|
292 |
|
|
@return 0 on success |
293 |
|
|
*/ |
294 |
|
|
int integrator_lsode_params_default(IntegratorSystem *blsys){ |
295 |
|
|
|
296 |
|
|
asc_assert(blsys!=NULL); |
297 |
johnpye |
962 |
asc_assert(blsys->engine==INTEG_LSODE); |
298 |
johnpye |
942 |
slv_parameters_t *p; |
299 |
|
|
p = &(blsys->params); |
300 |
|
|
|
301 |
|
|
slv_destroy_parms(p); |
302 |
|
|
|
303 |
|
|
if(p->parms==NULL){ |
304 |
|
|
CONSOLE_DEBUG("params NULL"); |
305 |
|
|
p->parms = ASC_NEW_ARRAY(struct slv_parameter, LSODE_PARAMS_SIZE); |
306 |
|
|
if(p->parms==NULL)return -1; |
307 |
|
|
p->dynamic_parms = 1; |
308 |
|
|
}else{ |
309 |
|
|
asc_assert(p->num_parms == LSODE_PARAMS_SIZE); |
310 |
|
|
CONSOLE_DEBUG("reusing parm memory"); |
311 |
|
|
} |
312 |
|
|
|
313 |
|
|
/* reset the number of parameters to zero so that we can check it at the end */ |
314 |
|
|
p->num_parms = 0; |
315 |
|
|
|
316 |
|
|
slv_param_bool(p,LSODE_PARAM_TIMING |
317 |
|
|
,(SlvParameterInitBool){{"timing" |
318 |
|
|
,"Output timing statistics?",1,NULL |
319 |
|
|
}, TRUE} |
320 |
|
|
); |
321 |
|
|
|
322 |
johnpye |
963 |
slv_param_bool(p,LSODE_PARAM_ATOLVECT |
323 |
|
|
,(SlvParameterInitBool){{"atolvect" |
324 |
|
|
,"Use 'ode_atol' values as specified for each var?",1 |
325 |
|
|
,"If TRUE, values of 'ode_atol' are taken from your model and used " |
326 |
|
|
" in the integration. If FALSE, a scalar absolute tolerance (atol)" |
327 |
|
|
" is shared by all variables." |
328 |
|
|
}, TRUE } |
329 |
|
|
); |
330 |
|
|
|
331 |
|
|
slv_param_real(p,LSODE_PARAM_ATOL |
332 |
|
|
,(SlvParameterInitReal){{"atol" |
333 |
johnpye |
962 |
,"Scalar absolute error tolerance",1 |
334 |
|
|
,"Default value of the scalar absolute error tolerance (for cases" |
335 |
|
|
" where not specified in oda_atol var property. See 'lsode.f' for" |
336 |
|
|
" details" |
337 |
|
|
}, 1e-6, DBL_MIN, DBL_MAX } |
338 |
|
|
); |
339 |
|
|
|
340 |
johnpye |
963 |
slv_param_bool(p,LSODE_PARAM_RTOLVECT |
341 |
|
|
,(SlvParameterInitBool){{"rtolvect" |
342 |
|
|
,"Use 'ode_rtol' values as specified for each var?",1 |
343 |
|
|
,"If TRUE, values of 'ode_atol' are taken from your model and used " |
344 |
|
|
" in the integration. If FALSE, a scalar absolute tolerance (rtol)" |
345 |
|
|
" is shared by all variables." |
346 |
|
|
}, TRUE } |
347 |
|
|
); |
348 |
|
|
|
349 |
|
|
slv_param_real(p,LSODE_PARAM_RTOL |
350 |
|
|
,(SlvParameterInitReal){{"rtol" |
351 |
johnpye |
962 |
,"Scalar relative error tolerance",1 |
352 |
|
|
,"Default value of the scalar relative error tolerance (for cases" |
353 |
|
|
" where not specified in oda_rtol var property. See 'lsode.f' for" |
354 |
|
|
" details" |
355 |
|
|
}, 1e-6, DBL_MIN, DBL_MAX } |
356 |
|
|
); |
357 |
|
|
|
358 |
johnpye |
942 |
asc_assert(p->num_parms == LSODE_PARAMS_SIZE); |
359 |
|
|
|
360 |
|
|
CONSOLE_DEBUG("Created %d params", p->num_parms); |
361 |
|
|
|
362 |
|
|
return 0; |
363 |
|
|
} |
364 |
|
|
|
365 |
johnpye |
669 |
/*--------------------------------------------------------- |
366 |
|
|
Couple of matrix methods...? |
367 |
|
|
*/ |
368 |
|
|
|
369 |
johnpye |
888 |
static double **lsode_densematrix_create(int nrows, int ncols){ |
370 |
johnpye |
669 |
int c; |
371 |
|
|
double **result; |
372 |
|
|
assert(nrows>0); |
373 |
|
|
assert(ncols>0); |
374 |
|
|
result = ASC_NEW_ARRAY(double *, nrows); |
375 |
|
|
for (c=0;c<nrows;c++) { |
376 |
|
|
result[c] = ASC_NEW_ARRAY_CLEAR(double, ncols); |
377 |
|
|
} |
378 |
|
|
return result; |
379 |
|
|
} |
380 |
|
|
|
381 |
johnpye |
888 |
static void lsode_densematrix_destroy(double **matrix,int nrows){ |
382 |
johnpye |
669 |
int c; |
383 |
|
|
if (matrix) { |
384 |
|
|
for (c=0;c<nrows;c++) { |
385 |
|
|
if (matrix[c]) { |
386 |
|
|
ascfree((char *)matrix[c]); |
387 |
|
|
} |
388 |
|
|
} |
389 |
|
|
ascfree((char *)matrix); |
390 |
|
|
} |
391 |
ben.allan |
704 |
} |
392 |
johnpye |
669 |
|
393 |
|
|
/*------------------------------------------------------------------------------ |
394 |
|
|
PROBLEM ANALYSIS |
395 |
|
|
*/ |
396 |
|
|
|
397 |
ben.allan |
704 |
/** |
398 |
|
|
@TODO needs work. Assumes struct Instance* and struct var_variable* |
399 |
|
|
are synonymous, which demonstrates the need for a method to take |
400 |
|
|
an instance and ask the solvers for its global or local index |
401 |
|
|
if var and inst are decoupled. |
402 |
|
|
*/ |
403 |
|
|
int integrator_lsode_setup_diffs(IntegratorSystem *blsys) { |
404 |
johnpye |
908 |
/* long n_eqns; */ |
405 |
ben.allan |
704 |
unsigned long nch,i; |
406 |
|
|
|
407 |
|
|
struct var_variable **vp; |
408 |
johnpye |
669 |
int *ip; |
409 |
|
|
|
410 |
|
|
IntegratorLsodeData *enginedata; |
411 |
ben.allan |
704 |
enginedata = (IntegratorLsodeData *)blsys->enginedata; |
412 |
johnpye |
669 |
assert(enginedata!=NULL); |
413 |
|
|
|
414 |
|
|
assert(enginedata->n_eqns==blsys->n_y); |
415 |
ben.allan |
704 |
|
416 |
johnpye |
669 |
/* |
417 |
johnpye |
741 |
Put the |
418 |
ben.allan |
704 |
Let us now process what we consider *inputs* to the problem as |
419 |
|
|
far as ASCEND is concerned; i.e. the state vars or the y_vars's |
420 |
|
|
if you prefer. |
421 |
|
|
*/ |
422 |
johnpye |
669 |
nch = enginedata->n_eqns; |
423 |
|
|
|
424 |
ben.allan |
704 |
|
425 |
|
|
vp = enginedata->y_vars; |
426 |
|
|
ip = enginedata->input_indices; |
427 |
|
|
for (i=0;i<nch;i++) { |
428 |
|
|
*vp = (struct var_variable *)blsys->y[i]; |
429 |
|
|
*ip = var_sindex(*vp); |
430 |
|
|
vp++; |
431 |
|
|
ip++; |
432 |
|
|
} |
433 |
|
|
*vp = NULL; /* terminate */ |
434 |
|
|
|
435 |
|
|
/* |
436 |
|
|
Let us now go for the outputs, ie the derivative terms. |
437 |
|
|
*/ |
438 |
|
|
vp = enginedata->ydot_vars; |
439 |
|
|
ip = enginedata->output_indices; |
440 |
|
|
for (i=0;i<nch;i++) { |
441 |
|
|
*vp = (struct var_variable *)blsys->ydot[i]; |
442 |
|
|
*ip = var_sindex(*vp); |
443 |
|
|
vp++; /* dont assume that a var is synonymous with */ |
444 |
|
|
ip++; /* an Instance; that might/will change soon */ |
445 |
|
|
} |
446 |
|
|
*vp = NULL; /* terminate */ |
447 |
|
|
|
448 |
|
|
return 0; |
449 |
johnpye |
669 |
} |
450 |
|
|
|
451 |
|
|
/** |
452 |
johnpye |
888 |
allocates, fills, and returns the atol vector based on LSODE |
453 |
johnpye |
669 |
|
454 |
johnpye |
963 |
State variables missing child ode_rtol will be defaulted to ATOL |
455 |
johnpye |
669 |
*/ |
456 |
johnpye |
888 |
static double *lsode_get_atol( IntegratorSystem *blsys) { |
457 |
johnpye |
669 |
|
458 |
|
|
struct Instance *tol; |
459 |
|
|
double *atoli; |
460 |
|
|
int i,len; |
461 |
johnpye |
963 |
double atol; |
462 |
johnpye |
669 |
|
463 |
|
|
len = blsys->n_y; |
464 |
johnpye |
963 |
atoli = ASC_NEW_ARRAY(double, blsys->n_y); /* changed, this was n_y+1 before, dunnowi -- JP */ |
465 |
johnpye |
669 |
if (atoli == NULL) { |
466 |
|
|
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Insufficient memory"); |
467 |
|
|
return atoli; |
468 |
|
|
} |
469 |
johnpye |
963 |
|
470 |
|
|
if(!SLV_PARAM_BOOL(&(blsys->params),LSODE_PARAM_ATOLVECT)){ |
471 |
|
|
atol = SLV_PARAM_REAL(&(blsys->params),LSODE_PARAM_ATOL); |
472 |
|
|
CONSOLE_DEBUG("Using ATOL = %f for all vars", atol); |
473 |
|
|
for(i=0; i<len; ++i){ |
474 |
|
|
atoli[i] = atol; |
475 |
|
|
} |
476 |
|
|
}else{ |
477 |
|
|
InitTolNames(); |
478 |
|
|
for (i=0; i<len; i++) { |
479 |
|
|
|
480 |
|
|
tol = ChildByChar(var_instance(blsys->y[i]),STATEATOL); |
481 |
|
|
if (tol == NULL || !AtomAssigned(tol) ) { |
482 |
|
|
atoli[i] = SLV_PARAM_REAL(&(blsys->params),LSODE_PARAM_ATOL); |
483 |
|
|
ERROR_REPORTER_HERE(ASC_PROG_WARNING,"Assuming atol = %3g" |
484 |
|
|
"for ode_atol child undefined for state variable %ld." |
485 |
|
|
,atoli[i], blsys->y_id[i] |
486 |
|
|
); |
487 |
|
|
} else { |
488 |
|
|
atoli[i] = RealAtomValue(tol); |
489 |
|
|
CONSOLE_DEBUG("Using atol %3g for state variable %ld.",atoli[i], blsys->y_id[i]); |
490 |
|
|
} |
491 |
johnpye |
669 |
} |
492 |
|
|
} |
493 |
|
|
return atoli; |
494 |
|
|
} |
495 |
|
|
|
496 |
|
|
/** |
497 |
johnpye |
888 |
Allocates, fills, and returns the rtol vector based on LSODE |
498 |
johnpye |
669 |
|
499 |
johnpye |
963 |
State variables missing child ode_rtol will be defaulted to RTOL |
500 |
johnpye |
669 |
*/ |
501 |
johnpye |
888 |
static double *lsode_get_rtol( IntegratorSystem *blsys) { |
502 |
johnpye |
669 |
|
503 |
|
|
struct Instance *tol; |
504 |
johnpye |
963 |
double rtol, *rtoli; |
505 |
johnpye |
669 |
int i,len; |
506 |
|
|
|
507 |
|
|
len = blsys->n_y; |
508 |
|
|
rtoli = ASC_NEW_ARRAY(double, blsys->n_y+1); |
509 |
|
|
if (rtoli == NULL) { |
510 |
|
|
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Insufficient memory"); |
511 |
|
|
return rtoli; |
512 |
|
|
} |
513 |
johnpye |
963 |
if(!SLV_PARAM_BOOL(&(blsys->params),LSODE_PARAM_RTOLVECT)){ |
514 |
|
|
rtol = SLV_PARAM_REAL(&(blsys->params),LSODE_PARAM_RTOL); |
515 |
|
|
CONSOLE_DEBUG("Using RTOL = %f for all vars", rtol); |
516 |
|
|
for(i=0; i<len; ++i){ |
517 |
|
|
rtoli[i] = rtol; |
518 |
|
|
} |
519 |
|
|
}else{ |
520 |
|
|
InitTolNames(); |
521 |
|
|
for (i=0; i<len; i++) { |
522 |
|
|
tol = ChildByChar(var_instance(blsys->y[i]),STATERTOL); |
523 |
|
|
if (tol == NULL || !AtomAssigned(tol) ) { |
524 |
|
|
rtoli[i] = SLV_PARAM_REAL(&(blsys->params),LSODE_PARAM_RTOL); |
525 |
johnpye |
669 |
|
526 |
johnpye |
963 |
ERROR_REPORTER_HERE(ASC_PROG_WARNING,"Assuming rtol = %3g" |
527 |
|
|
"for ode_rtol child undefined for state variable %ld." |
528 |
|
|
,rtoli[i], blsys->y_id[i] |
529 |
|
|
); |
530 |
|
|
|
531 |
|
|
} else { |
532 |
|
|
rtoli[i] = RealAtomValue(tol); |
533 |
|
|
} |
534 |
johnpye |
669 |
} |
535 |
|
|
} |
536 |
johnpye |
963 |
rtoli[len] = SLV_PARAM_REAL(&(blsys->params),LSODE_PARAM_RTOL); |
537 |
johnpye |
669 |
return rtoli; |
538 |
|
|
} |
539 |
|
|
|
540 |
|
|
/* |
541 |
|
|
Write out a a status message based on the istate parameter. |
542 |
|
|
*/ |
543 |
johnpye |
888 |
static void lsode_write_istate( int istate) { |
544 |
johnpye |
669 |
switch (istate) { |
545 |
|
|
case -1: |
546 |
|
|
FPRINTF(ASCERR,"Excess steps taken on this call (perhaps wrong MF)."); |
547 |
|
|
break; |
548 |
|
|
case -2: |
549 |
|
|
FPRINTF(ASCERR,"Excess accuracy requested (tolerances too small)."); |
550 |
|
|
break; |
551 |
|
|
case -3: |
552 |
johnpye |
888 |
FPRINTF(ASCERR,"Illegal input detected (see console)."); |
553 |
johnpye |
669 |
break; |
554 |
|
|
case -4: |
555 |
johnpye |
888 |
FPRINTF(ASCERR,"Repeated error test failures (check all inputs)."); |
556 |
johnpye |
669 |
break; |
557 |
|
|
case -5: |
558 |
johnpye |
888 |
FPRINTF(ASCERR,"Repeated convergence failures (perhaps bad Jacobian supplied, or wrong choice of MF or tolerances)."); |
559 |
johnpye |
669 |
break; |
560 |
|
|
case -6: |
561 |
johnpye |
888 |
FPRINTF(ASCERR,"Error weight became zero during problem (solution component i vanished, and atol or atol(i) = 0)."); |
562 |
johnpye |
669 |
break; |
563 |
|
|
case -7: |
564 |
johnpye |
888 |
FPRINTF(ASCERR,"Interrupted? User cancelled operation?"); |
565 |
johnpye |
669 |
break; |
566 |
|
|
default: |
567 |
johnpye |
888 |
FPRINTF(ASCERR,"Unknown 'istate' error code %d from LSODE.",istate); |
568 |
johnpye |
669 |
break; |
569 |
|
|
} |
570 |
|
|
} |
571 |
|
|
|
572 |
|
|
/** |
573 |
|
|
Free memory allocated for the LSODE, but first check. |
574 |
|
|
*/ |
575 |
johnpye |
888 |
static void lsode_free_mem(double *y, double *reltol, double *abtol, double *rwork, |
576 |
johnpye |
669 |
int *iwork, double *obs, double *dydx) |
577 |
|
|
{ |
578 |
|
|
if (y != NULL) { |
579 |
|
|
ascfree((double *)y); |
580 |
|
|
} |
581 |
|
|
if (reltol != NULL) { |
582 |
|
|
ascfree((double *)reltol); |
583 |
|
|
} |
584 |
|
|
if (abtol != NULL) { |
585 |
|
|
ascfree((double *)abtol); |
586 |
|
|
} |
587 |
|
|
if (rwork != NULL) { |
588 |
|
|
ascfree((double *)rwork); |
589 |
|
|
} |
590 |
|
|
if (iwork != NULL) { |
591 |
|
|
ascfree((int *)iwork); |
592 |
|
|
} |
593 |
|
|
if (obs != NULL) { |
594 |
|
|
ascfree((double *)obs); |
595 |
|
|
} |
596 |
|
|
if (dydx != NULL) { |
597 |
|
|
ascfree((double *)dydx); |
598 |
|
|
} |
599 |
|
|
} |
600 |
|
|
|
601 |
|
|
/* |
602 |
|
|
********************************************************************* |
603 |
|
|
* This code is provided for the benefit of a temporary |
604 |
|
|
* fix for the derivative problem in Lsode. |
605 |
|
|
* The proper permanent fix for lsode is to dump it in favor of |
606 |
|
|
* cvode or dassl. |
607 |
johnpye |
888 |
* Extended 7/95 baa to deal with linsolqr and lsode. |
608 |
johnpye |
669 |
* It is assumed the system has been solved at the current point. |
609 |
|
|
********************************************************************* |
610 |
|
|
*/ |
611 |
johnpye |
888 |
int lsode_derivatives(slv_system_t sys, double **dy_dx, |
612 |
johnpye |
669 |
int *inputs_ndx_list, int ninputs, |
613 |
|
|
int *outputs_ndx_list, int noutputs) |
614 |
|
|
{ |
615 |
|
|
static int n_calls = 0; |
616 |
|
|
linsolqr_system_t lqr_sys; /* stuff for the linear system & matrix */ |
617 |
|
|
mtx_matrix_t mtx; |
618 |
|
|
int32 capacity; |
619 |
|
|
real64 *scratch_vector = NULL; |
620 |
|
|
int result=0; |
621 |
|
|
|
622 |
|
|
(void)NumberFreeVars(NULL); /* used to re-init the system */ |
623 |
|
|
(void)NumberIncludedRels(NULL); /* used to re-init the system */ |
624 |
|
|
if (!sys) { |
625 |
|
|
FPRINTF(stderr,"The solve system does not exist !\n"); |
626 |
|
|
return 1; |
627 |
|
|
} |
628 |
|
|
|
629 |
|
|
result = Compute_J(sys); |
630 |
|
|
if (result) { |
631 |
|
|
FPRINTF(stderr,"Early termination due to failure in calc Jacobian\n"); |
632 |
|
|
return 1; |
633 |
|
|
} |
634 |
|
|
|
635 |
|
|
lqr_sys = slv_get_linsolqr_sys(sys); /* get the linear system */ |
636 |
|
|
if (lqr_sys==NULL) { |
637 |
|
|
FPRINTF(stderr,"Early termination due to missing linsolqr system.\n"); |
638 |
|
|
return 1; |
639 |
|
|
} |
640 |
|
|
mtx = slv_get_sys_mtx(sys); /* get the matrix */ |
641 |
|
|
if (mtx==NULL) { |
642 |
|
|
FPRINTF(stderr,"Early termination due to missing mtx in linsolqr.\n"); |
643 |
|
|
return 1; |
644 |
|
|
} |
645 |
|
|
capacity = mtx_capacity(mtx); |
646 |
|
|
scratch_vector = ASC_NEW_ARRAY_CLEAR(real64,capacity); |
647 |
|
|
linsolqr_add_rhs(lqr_sys,scratch_vector,FALSE); |
648 |
|
|
|
649 |
|
|
result = LUFactorJacobian(sys); |
650 |
|
|
if (result) { |
651 |
|
|
FPRINTF(stderr,"Early termination due to failure in LUFactorJacobian\n"); |
652 |
|
|
goto error; |
653 |
|
|
} |
654 |
|
|
result = Compute_dy_dx_smart(sys, scratch_vector, dy_dx, |
655 |
|
|
inputs_ndx_list, ninputs, |
656 |
|
|
outputs_ndx_list, noutputs); |
657 |
|
|
|
658 |
|
|
linsolqr_remove_rhs(lqr_sys,scratch_vector); |
659 |
|
|
if (result) { |
660 |
|
|
FPRINTF(stderr,"Early termination due to failure in Compute_dy_dx\n"); |
661 |
|
|
goto error; |
662 |
|
|
} |
663 |
|
|
|
664 |
|
|
error: |
665 |
|
|
n_calls++; |
666 |
|
|
if (scratch_vector) { |
667 |
|
|
ascfree((char *)scratch_vector); |
668 |
|
|
} |
669 |
|
|
return result; |
670 |
|
|
} |
671 |
|
|
|
672 |
|
|
/** |
673 |
|
|
The current way that we are getting the derivatives (if the problem |
674 |
|
|
was solved partitioned) messes up the slv_system so that we *have* |
675 |
|
|
to do a *presolve* rather than a simply a *resolve* before doing |
676 |
|
|
function calls. This code below attempts to handle these cases. |
677 |
|
|
*/ |
678 |
|
|
static void LSODE_FEX( int *n_eq ,double *t ,double *y ,double *ydot) |
679 |
|
|
{ |
680 |
|
|
slv_status_t status; |
681 |
|
|
|
682 |
|
|
/* slv_parameters_t parameters; pity lsode doesn't allow error returns */ |
683 |
johnpye |
908 |
/* int i; */ |
684 |
johnpye |
1049 |
unsigned long res; |
685 |
johnpye |
669 |
|
686 |
|
|
#if DOTIME |
687 |
|
|
double time1,time2; |
688 |
|
|
#endif |
689 |
|
|
|
690 |
johnpye |
938 |
/* CONSOLE_DEBUG("Calling for a function evaluation"); */ |
691 |
johnpye |
903 |
|
692 |
johnpye |
669 |
#if DOTIME |
693 |
|
|
CONSOLE_DEBUG("Calling for a function evaluation"); |
694 |
|
|
time1 = tm_cpu_time(); |
695 |
|
|
#endif |
696 |
|
|
|
697 |
|
|
/* |
698 |
johnpye |
888 |
t[1]=t[0]; can't do this. lsode calls us with a different t than the t we sent in. |
699 |
johnpye |
669 |
*/ |
700 |
|
|
integrator_set_t(l_lsode_blsys, t[0]); |
701 |
|
|
integrator_set_y(l_lsode_blsys, y); |
702 |
|
|
|
703 |
|
|
#if DOTIME |
704 |
|
|
time2 = tm_cpu_time(); |
705 |
|
|
#endif |
706 |
|
|
|
707 |
|
|
switch(lsodesys.lastcall) { |
708 |
|
|
case lsode_none: /* first call */ |
709 |
johnpye |
903 |
CONSOLE_DEBUG("FIRST CALL..."); |
710 |
|
|
|
711 |
johnpye |
669 |
case lsode_derivative: |
712 |
|
|
if (lsodesys.partitioned) { |
713 |
johnpye |
938 |
/* CONSOLE_DEBUG("PRE-SOLVE"); */ |
714 |
johnpye |
669 |
slv_presolve(l_lsode_blsys->system); |
715 |
|
|
} else { |
716 |
johnpye |
938 |
/** @TODO this doesn't ever seem to be called */ |
717 |
johnpye |
903 |
CONSOLE_DEBUG("RE-SOLVE"); |
718 |
johnpye |
669 |
slv_resolve(l_lsode_blsys->system); |
719 |
|
|
} |
720 |
|
|
break; |
721 |
|
|
default: |
722 |
|
|
case lsode_function: |
723 |
|
|
slv_resolve(l_lsode_blsys->system); |
724 |
|
|
break; |
725 |
|
|
} |
726 |
johnpye |
903 |
|
727 |
johnpye |
669 |
slv_solve(l_lsode_blsys->system); |
728 |
|
|
slv_get_status(l_lsode_blsys->system, &status); |
729 |
johnpye |
1055 |
if(slv_check_bounds(l_lsode_blsys->system,0,-1,"")){ |
730 |
|
|
lsodesys.status = lsode_nok; |
731 |
|
|
} |
732 |
|
|
|
733 |
johnpye |
888 |
/* pass the solver status to the integrator */ |
734 |
johnpye |
1049 |
res = integrator_checkstatus(status); |
735 |
johnpye |
669 |
|
736 |
|
|
#if DOTIME |
737 |
|
|
time2 = tm_cpu_time() - time2; |
738 |
|
|
#endif |
739 |
|
|
|
740 |
johnpye |
1049 |
if(res){ |
741 |
|
|
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Failed to solve for derivatives (%d)",res); |
742 |
johnpye |
1055 |
#if 0 |
743 |
johnpye |
669 |
ERROR_REPORTER_START_HERE(ASC_PROG_ERR); |
744 |
|
|
FPRINTF(ASCERR,"Unable to compute the vector of derivatives with the following values for the state variables:\n"); |
745 |
|
|
for (i = 0; i< *n_eq; i++) { |
746 |
|
|
FPRINTF(ASCERR,"y[%4d] = %f\n",i, y[i]); |
747 |
|
|
} |
748 |
|
|
error_reporter_end_flush(); |
749 |
johnpye |
1055 |
#endif |
750 |
johnpye |
669 |
lsodesys.status = lsode_nok; |
751 |
johnpye |
1049 |
}else{ |
752 |
johnpye |
669 |
lsodesys.status = lsode_ok; |
753 |
|
|
} |
754 |
|
|
integrator_get_ydot(l_lsode_blsys, ydot); |
755 |
|
|
|
756 |
|
|
lsodesys.lastcall = lsode_function; |
757 |
|
|
#if DOTIME |
758 |
|
|
time1 = tm_cpu_time() - time1; |
759 |
|
|
CONSOLE_DEBUG("Function evalulation has been completed in time %g. True function call time = %g",time1,time2); |
760 |
|
|
#endif |
761 |
|
|
} |
762 |
|
|
|
763 |
|
|
/** |
764 |
|
|
Evaluate the jacobian |
765 |
|
|
*/ |
766 |
|
|
static void LSODE_JEX(int *neq ,double *t, double *y, |
767 |
|
|
int *ml ,int *mu ,double *pd, int *nrpd) |
768 |
|
|
{ |
769 |
|
|
int nok = 0; |
770 |
|
|
int i,j; |
771 |
|
|
|
772 |
|
|
IntegratorLsodeData enginedata=*((IntegratorLsodeData *)l_lsode_blsys->enginedata); |
773 |
|
|
|
774 |
|
|
UNUSED_PARAMETER(t); |
775 |
|
|
UNUSED_PARAMETER(y); |
776 |
|
|
UNUSED_PARAMETER(ml); |
777 |
|
|
UNUSED_PARAMETER(mu); |
778 |
|
|
|
779 |
johnpye |
938 |
/* CONSOLE_DEBUG("Calling for a gradient evaluation"); */ |
780 |
johnpye |
669 |
#if DOTIME |
781 |
|
|
double time1; |
782 |
|
|
|
783 |
|
|
CONSOLE_DEBUG("Calling for a gradient evaluation"); |
784 |
|
|
time1 = tm_cpu_time(); |
785 |
|
|
#endif |
786 |
|
|
/* |
787 |
|
|
* Make the real call. |
788 |
|
|
*/ |
789 |
johnpye |
888 |
nok = lsode_derivatives(l_lsode_blsys->system |
790 |
johnpye |
669 |
, enginedata.dydx_dx |
791 |
|
|
, enginedata.input_indices |
792 |
|
|
, *neq |
793 |
|
|
, enginedata.output_indices |
794 |
|
|
, *nrpd |
795 |
|
|
); |
796 |
|
|
|
797 |
|
|
if (nok) { |
798 |
|
|
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Error in computing the derivatives for the system. Failing..."); |
799 |
|
|
lsodesys.status = lsode_nok; |
800 |
|
|
lsodesys.lastcall = lsode_derivative; |
801 |
|
|
return; |
802 |
|
|
} else { |
803 |
|
|
lsodesys.status = lsode_ok; |
804 |
|
|
lsodesys.lastcall = lsode_derivative; |
805 |
|
|
} |
806 |
|
|
/* |
807 |
johnpye |
903 |
Map data from C based matrix to Fortan matrix. |
808 |
|
|
We will send in a column major ordering vector for pd. |
809 |
|
|
*/ |
810 |
johnpye |
888 |
for (j=0;j<*neq;j++) { /* loop through columnns */ |
811 |
|
|
for (i=0;i<*nrpd;i++){ /* loop through rows */ |
812 |
|
|
/* CONSOLE_DEBUG("JAC[r=%d,c=%d]=%f",i,j,enginedata.dydx_dx[i][j]); */ |
813 |
johnpye |
669 |
*pd++ = enginedata.dydx_dx[i][j]; |
814 |
|
|
} |
815 |
|
|
} |
816 |
|
|
|
817 |
|
|
#if DOTIME |
818 |
|
|
time1 = tm_cpu_time() - time1; |
819 |
|
|
CONSOLE_DEBUG("Time to do gradient evaluation %g",time1); |
820 |
|
|
#endif |
821 |
|
|
|
822 |
|
|
return; |
823 |
|
|
} |
824 |
|
|
|
825 |
|
|
/** |
826 |
|
|
The public function: here we do the actual integration, I guess. |
827 |
johnpye |
966 |
|
828 |
johnpye |
1049 |
Return 0 on success |
829 |
johnpye |
669 |
*/ |
830 |
|
|
int integrator_lsode_solve(IntegratorSystem *blsys |
831 |
|
|
, unsigned long start_index, unsigned long finish_index |
832 |
|
|
){ |
833 |
|
|
slv_status_t status; |
834 |
|
|
slv_parameters_t params; |
835 |
johnpye |
741 |
IntegratorLsodeData *d; |
836 |
johnpye |
669 |
|
837 |
|
|
double x[2]; |
838 |
|
|
double xend,xprev; |
839 |
|
|
unsigned long nsamples, neq; |
840 |
|
|
long nobs; |
841 |
|
|
int itol, itask, mf, lrw, liw; |
842 |
|
|
unsigned long index; |
843 |
|
|
int istate, iopt; |
844 |
|
|
double * rwork; |
845 |
|
|
int * iwork; |
846 |
|
|
double *y, *abtol, *reltol, *obs, *dydx; |
847 |
|
|
int my_neq; |
848 |
|
|
FILE *y_out =NULL; |
849 |
|
|
FILE *obs_out =NULL; |
850 |
johnpye |
902 |
int reporterstatus; |
851 |
johnpye |
669 |
|
852 |
|
|
/* store the local variable so that we can get at stuff from inside LSODE_FEX. */ |
853 |
|
|
l_lsode_blsys = blsys; |
854 |
|
|
|
855 |
|
|
d = (IntegratorLsodeData *)(blsys->enginedata); |
856 |
|
|
|
857 |
ben.allan |
704 |
/* the numer of equations must be equal to blsys->n_y, the number of states */ |
858 |
|
|
d->n_eqns = blsys->n_y; |
859 |
|
|
assert(d->n_eqns>0); |
860 |
johnpye |
669 |
|
861 |
ben.allan |
704 |
d->input_indices = ASC_NEW_ARRAY_CLEAR(int, d->n_eqns); |
862 |
|
|
d->output_indices = ASC_NEW_ARRAY_CLEAR(int, d->n_eqns); |
863 |
johnpye |
888 |
d->dydx_dx = lsode_densematrix_create(d->n_eqns,d->n_eqns); |
864 |
johnpye |
669 |
|
865 |
|
|
d->y_vars = ASC_NEW_ARRAY(struct var_variable *,d->n_eqns+1); |
866 |
|
|
d->ydot_vars = ASC_NEW_ARRAY(struct var_variable *, d->n_eqns+1); |
867 |
|
|
|
868 |
|
|
integrator_lsode_setup_diffs(blsys); |
869 |
|
|
|
870 |
ben.allan |
704 |
/* this is a lie, but we will keep it. |
871 |
|
|
We handle any linsol/linsolqr based solver. */ |
872 |
|
|
if (strcmp(slv_solver_name(slv_get_selected_solver(blsys->system)),"QRSlv") != 0) { |
873 |
|
|
ERROR_REPORTER_NOLINE(ASC_USER_ERROR,"QRSlv must be selected before integration."); |
874 |
johnpye |
1049 |
return 1; |
875 |
johnpye |
669 |
} |
876 |
|
|
|
877 |
|
|
slv_get_status(l_lsode_blsys->system, &status); |
878 |
|
|
|
879 |
|
|
if (status.struct_singular) { |
880 |
|
|
ERROR_REPORTER_HERE(ASC_USER_ERROR,"Integration will not be performed. The system is structurally singular."); |
881 |
|
|
lsodesys.status = lsode_nok; |
882 |
johnpye |
1049 |
return 2; |
883 |
johnpye |
669 |
} |
884 |
|
|
|
885 |
|
|
#if defined(STATIC_LSOD) || defined (DYNAMIC_LSOD) |
886 |
|
|
|
887 |
|
|
/* here we assume integrators.c is in charge of dynamic loading */ |
888 |
|
|
|
889 |
|
|
slv_get_parameters(blsys->system,¶ms); |
890 |
|
|
lsodesys.partitioned = 1; |
891 |
|
|
|
892 |
|
|
nsamples = integrator_getnsamples(blsys); |
893 |
|
|
if (nsamples <2) { |
894 |
|
|
ERROR_REPORTER_HERE(ASC_USER_ERROR,"Integration will not be performed. The system has no end sample time defined."); |
895 |
|
|
lsodesys.status = lsode_nok; |
896 |
johnpye |
1049 |
return 3; |
897 |
johnpye |
669 |
} |
898 |
|
|
neq = blsys->n_y; |
899 |
|
|
nobs = blsys->n_obs; |
900 |
|
|
|
901 |
johnpye |
966 |
/* samplelist_debug(blsys->samples); */ |
902 |
|
|
|
903 |
|
|
/* x[0] = integrator_get_t(blsys); */ |
904 |
|
|
x[0] = integrator_getsample(blsys, 0); |
905 |
johnpye |
669 |
x[1] = x[0]-1; /* make sure we don't start with wierd x[1] */ |
906 |
|
|
lrw = 22 + 9*neq + neq*neq; |
907 |
|
|
rwork = ASC_NEW_ARRAY_CLEAR(double, lrw+1); |
908 |
|
|
liw = 20 + neq; |
909 |
|
|
iwork = ASC_NEW_ARRAY_CLEAR(int, liw+1); |
910 |
|
|
y = integrator_get_y(blsys, NULL); |
911 |
johnpye |
888 |
reltol = lsode_get_rtol(blsys); |
912 |
|
|
abtol = lsode_get_atol(blsys); |
913 |
johnpye |
669 |
obs = integrator_get_observations(blsys, NULL); |
914 |
|
|
dydx = ASC_NEW_ARRAY_CLEAR(double, neq+1); |
915 |
|
|
if (!y || !obs || !abtol || !reltol || !rwork || !iwork || !dydx) { |
916 |
johnpye |
888 |
lsode_free_mem(y,reltol,abtol,rwork,iwork,obs,dydx); |
917 |
|
|
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Insufficient memory for lsode."); |
918 |
johnpye |
669 |
lsodesys.status = lsode_nok; |
919 |
johnpye |
1049 |
return 4; |
920 |
johnpye |
669 |
} |
921 |
|
|
|
922 |
|
|
/* |
923 |
|
|
Prepare args and call lsode. |
924 |
|
|
*/ |
925 |
|
|
itol = 4; |
926 |
|
|
itask = 1; |
927 |
|
|
istate = 1; |
928 |
|
|
iopt = 1; |
929 |
|
|
rwork[4] = integrator_get_stepzero(blsys); |
930 |
|
|
rwork[5] = integrator_get_maxstep(blsys); |
931 |
|
|
rwork[6] = integrator_get_minstep(blsys); |
932 |
|
|
iwork[5] = integrator_get_maxsubsteps(blsys); |
933 |
johnpye |
888 |
mf = 21; /* 21 = BDF with exact jacobian. 22 = BDF with finite diff Jacobian */ |
934 |
johnpye |
669 |
|
935 |
johnpye |
966 |
if(x[0] > integrator_getsample(blsys, 2)){ |
936 |
|
|
ERROR_REPORTER_HERE(ASC_USER_ERROR,"Invalid initialisation time: exceeds second timestep value"); |
937 |
johnpye |
1049 |
return 5; |
938 |
johnpye |
966 |
} |
939 |
|
|
|
940 |
johnpye |
669 |
/* put the values from derivative system into the record */ |
941 |
|
|
integrator_setsample(blsys, start_index, x[0]); |
942 |
|
|
|
943 |
|
|
integrator_output_init(blsys); |
944 |
|
|
|
945 |
johnpye |
979 |
/* -- store the initial values of all the stuff */ |
946 |
|
|
integrator_output_write(blsys); |
947 |
|
|
integrator_output_write_obs(blsys); |
948 |
|
|
|
949 |
johnpye |
669 |
my_neq = (int)neq; |
950 |
|
|
|
951 |
|
|
/* |
952 |
|
|
First time entering lsode, x is input. After that, |
953 |
|
|
lsode uses x as output (y output is y(x)). To drive |
954 |
|
|
the loop ahead in time, all we need to do is keep upping |
955 |
|
|
xend. |
956 |
|
|
*/ |
957 |
johnpye |
966 |
|
958 |
johnpye |
669 |
blsys->currentstep = 0; |
959 |
|
|
for (index = start_index; index < finish_index; index++, blsys->currentstep++) { |
960 |
|
|
xend = integrator_getsample(blsys, index+1); |
961 |
|
|
xprev = x[0]; |
962 |
johnpye |
966 |
asc_assert(xend > xprev); |
963 |
|
|
/* CONSOLE_DEBUG("LSODE call #%lu: x = [%f,%f]", index,xprev,xend); */ |
964 |
johnpye |
669 |
|
965 |
|
|
# ifndef NO_SIGNAL_TRAPS |
966 |
johnpye |
997 |
if (SETJMP(g_fpe_env)==0) { |
967 |
johnpye |
669 |
# endif /* NO_SIGNAL_TRAPS */ |
968 |
|
|
|
969 |
johnpye |
895 |
/* CONSOLE_DEBUG("Calling LSODE with end-time = %f",xend); */ |
970 |
|
|
/* |
971 |
|
|
switch(mf){ |
972 |
johnpye |
888 |
case 10: |
973 |
|
|
CONSOLE_DEBUG("Non-stiff (Adams) method; no Jacobian will be used"); break; |
974 |
|
|
case 21: |
975 |
|
|
CONSOLE_DEBUG("Stiff (BDF) method, user-supplied full Jacobian"); break; |
976 |
|
|
case 22: |
977 |
|
|
CONSOLE_DEBUG("Stiff (BDF) method, internally generated full Jacobian"); break; |
978 |
|
|
case 24: |
979 |
|
|
CONSOLE_DEBUG("Stiff (BDF) method, user-supplied banded jacobian"); break; |
980 |
|
|
case 25: |
981 |
|
|
CONSOLE_DEBUG("Stiff (BDF) method, internally generated banded jacobian"); break; |
982 |
|
|
default: |
983 |
|
|
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Invalid method id %d for LSODE",mf); |
984 |
johnpye |
895 |
return 0; * failure * |
985 |
johnpye |
888 |
} |
986 |
johnpye |
895 |
*/ |
987 |
johnpye |
669 |
|
988 |
|
|
LSODE(&(LSODE_FEX), &my_neq, y, x, &xend, |
989 |
|
|
&itol, reltol, abtol, &itask, &istate, |
990 |
|
|
&iopt ,rwork, &lrw, iwork, &liw, &(LSODE_JEX), &mf); |
991 |
|
|
|
992 |
|
|
|
993 |
|
|
# ifndef NO_SIGNAL_TRAPS |
994 |
|
|
} else { |
995 |
|
|
FPRINTF(stderr, |
996 |
|
|
"Integration terminated due to float error in LSODE call.\n"); |
997 |
johnpye |
888 |
lsode_free_mem(y,reltol,abtol,rwork,iwork,obs,dydx); |
998 |
johnpye |
669 |
lsodesys.status = lsode_ok; /* clean up before we go */ |
999 |
|
|
lsodesys.lastcall = lsode_none; |
1000 |
|
|
if (y_out!=NULL) { |
1001 |
|
|
fclose(y_out); |
1002 |
|
|
} |
1003 |
|
|
if (obs_out!=NULL) { |
1004 |
|
|
fclose(obs_out); |
1005 |
|
|
} |
1006 |
johnpye |
1049 |
return 6; |
1007 |
johnpye |
669 |
} |
1008 |
|
|
# endif /* NO_SIGNAL_TRAPS */ |
1009 |
|
|
|
1010 |
|
|
/* CONSOLE_DEBUG("AFTER %lu LSODE CALL\n", index); */ |
1011 |
|
|
/* this check is better done in fex,jex, but lsode takes no status */ |
1012 |
|
|
if (Solv_C_CheckHalt()) { |
1013 |
|
|
if (istate >= 0) { |
1014 |
|
|
istate=-7; |
1015 |
|
|
} |
1016 |
|
|
} |
1017 |
|
|
|
1018 |
|
|
if (istate < 0 ) { |
1019 |
|
|
/* some kind of error occurred... */ |
1020 |
|
|
ERROR_REPORTER_START_HERE(ASC_PROG_ERR); |
1021 |
johnpye |
888 |
lsode_write_istate(istate); |
1022 |
|
|
FPRINTF(ASCERR, "\nFurthest point reached was t = %g.\n",x[0]); |
1023 |
johnpye |
669 |
error_reporter_end_flush(); |
1024 |
|
|
|
1025 |
johnpye |
888 |
lsode_free_mem(y,reltol,abtol,rwork,iwork,obs,dydx); |
1026 |
johnpye |
669 |
integrator_output_close(blsys); |
1027 |
johnpye |
1049 |
return 7; |
1028 |
johnpye |
669 |
} |
1029 |
|
|
|
1030 |
|
|
if (lsodesys.status==lsode_nok) { |
1031 |
|
|
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Integration terminated due to an error in derivative computations."); |
1032 |
johnpye |
888 |
lsode_free_mem(y,reltol,abtol,rwork,iwork,obs,dydx); |
1033 |
johnpye |
669 |
lsodesys.status = lsode_ok; /* clean up before we go */ |
1034 |
|
|
lsodesys.lastcall = lsode_none; |
1035 |
|
|
integrator_output_close(blsys); |
1036 |
johnpye |
1049 |
return 8; |
1037 |
johnpye |
669 |
} |
1038 |
|
|
|
1039 |
|
|
integrator_setsample(blsys, index+1, x[0]); |
1040 |
|
|
/* record when lsode actually came back */ |
1041 |
|
|
integrator_set_t(blsys, x[0]); |
1042 |
|
|
integrator_set_y(blsys, y); |
1043 |
|
|
/* put x,y in d in case lsode got x,y by interpolation, as it does */ |
1044 |
|
|
|
1045 |
johnpye |
902 |
reporterstatus = integrator_output_write(blsys); |
1046 |
johnpye |
669 |
|
1047 |
johnpye |
902 |
if(reporterstatus==0){ |
1048 |
|
|
ERROR_REPORTER_HERE(ASC_USER_ERROR,"Integration cancelled"); |
1049 |
|
|
lsode_free_mem(y,reltol,abtol,rwork,iwork,obs,dydx); |
1050 |
|
|
lsodesys.status = lsode_ok; |
1051 |
|
|
lsodesys.lastcall = lsode_none; |
1052 |
|
|
integrator_output_close(blsys); |
1053 |
johnpye |
1049 |
return 9; |
1054 |
johnpye |
902 |
} |
1055 |
|
|
|
1056 |
johnpye |
669 |
if (nobs > 0) { |
1057 |
|
|
# ifndef NO_SIGNAL_TRAPS |
1058 |
johnpye |
997 |
if (SETJMP(g_fpe_env)==0) { |
1059 |
johnpye |
669 |
# endif /* NO_SIGNAL_TRAPS */ |
1060 |
|
|
|
1061 |
|
|
/* solve for obs since d isn't necessarily already |
1062 |
|
|
computed there though lsode's x and y may be. |
1063 |
|
|
Note that since lsode usually steps beyond xend |
1064 |
|
|
x1 usually wouldn't be x0 precisely if the x1/x0 |
1065 |
|
|
scheme worked, which it doesn't anyway. */ |
1066 |
|
|
|
1067 |
|
|
LSODE_FEX(&my_neq, x, y, dydx); |
1068 |
|
|
|
1069 |
|
|
/* calculate observations, if any, at returned x and y. */ |
1070 |
|
|
obs = integrator_get_observations(blsys, obs); |
1071 |
|
|
|
1072 |
|
|
integrator_output_write_obs(blsys); |
1073 |
|
|
|
1074 |
|
|
# ifndef NO_SIGNAL_TRAPS |
1075 |
|
|
} else { |
1076 |
|
|
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Integration terminated due to float error in LSODE FEX call."); |
1077 |
johnpye |
888 |
lsode_free_mem(y,reltol,abtol,rwork,iwork,obs,dydx); |
1078 |
johnpye |
669 |
lsodesys.status = lsode_ok; /* clean up before we go */ |
1079 |
|
|
lsodesys.lastcall = lsode_none; |
1080 |
|
|
integrator_output_close(blsys); |
1081 |
johnpye |
1049 |
return 10; |
1082 |
johnpye |
669 |
} |
1083 |
|
|
# endif /* NO_SIGNAL_TRAPS */ |
1084 |
|
|
} |
1085 |
johnpye |
855 |
/* CONSOLE_DEBUG("Integration completed from %3g to %3g.",xprev,x[0]); */ |
1086 |
johnpye |
669 |
} |
1087 |
|
|
|
1088 |
johnpye |
888 |
CONSOLE_DEBUG("..."); |
1089 |
johnpye |
669 |
CONSOLE_DEBUG("Number of steps taken: %1d.", iwork[10]); |
1090 |
|
|
CONSOLE_DEBUG("Number of function evaluations: %1d.", iwork[11]); |
1091 |
|
|
CONSOLE_DEBUG("Number of Jacobian evaluations: %1d.", iwork[12]); |
1092 |
johnpye |
888 |
CONSOLE_DEBUG("..."); |
1093 |
johnpye |
669 |
|
1094 |
|
|
|
1095 |
johnpye |
888 |
lsode_free_mem(y,reltol,abtol,rwork,iwork,obs,dydx); |
1096 |
johnpye |
669 |
|
1097 |
|
|
/* |
1098 |
|
|
* return the system to its original state. |
1099 |
|
|
*/ |
1100 |
|
|
|
1101 |
|
|
lsodesys.status = lsode_ok; |
1102 |
|
|
lsodesys.lastcall = lsode_none; |
1103 |
|
|
|
1104 |
|
|
integrator_output_close(blsys); |
1105 |
|
|
|
1106 |
johnpye |
888 |
CONSOLE_DEBUG("--- LSODE done ---"); |
1107 |
johnpye |
1049 |
return 0; /* success */ |
1108 |
johnpye |
669 |
|
1109 |
|
|
#else /* STATIC_LSOD || DYNAMIC_LSOD */ |
1110 |
|
|
|
1111 |
|
|
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Integration will not be performed. LSODE binary not available."); |
1112 |
|
|
lsodesys.status = lsode_nok; |
1113 |
johnpye |
1049 |
return 11; |
1114 |
johnpye |
669 |
|
1115 |
|
|
#endif |
1116 |
|
|
} |
1117 |
|
|
|
1118 |
|
|
/** |
1119 |
johnpye |
888 |
Function XASCWV is an error reporting function replacing the XERRWV |
1120 |
|
|
routine in lsode.f. The call signature is the same with the original Fortran |
1121 |
|
|
function. |
1122 |
johnpye |
669 |
|
1123 |
johnpye |
888 |
@see the comments for 'xerrwv' from lsode.f, with which XASCWV is compatible... |
1124 |
johnpye |
669 |
|
1125 |
johnpye |
888 |
@param msg = the message (hollerith literal or integer array). |
1126 |
|
|
@param nmes = the length of msg (number of characters). |
1127 |
|
|
@param nerr = the error number (not used). |
1128 |
|
|
@param level = the error level.. |
1129 |
johnpye |
669 |
0 or 1 means recoverable (control returns to caller). |
1130 |
|
|
2 means fatal (run is aborted--see note below). |
1131 |
johnpye |
888 |
@param ni = number of integers (0, 1, or 2) to be printed with message. |
1132 |
|
|
@param i1,i2 = integers to be printed, depending on ni. |
1133 |
|
|
@param nr = number of reals (0, 1, or 2) to be printed with message. |
1134 |
|
|
@param r1,r2 = reals to be printed, depending on nr. |
1135 |
|
|
*/ |
1136 |
|
|
void XASCWV( char *msg, /* pointer to start of message */ |
1137 |
|
|
int *nmes, /* the length of msg (number of characters) */ |
1138 |
|
|
int *nerr, /* the error number (not used). */ |
1139 |
johnpye |
908 |
int *level, |
1140 |
johnpye |
669 |
int *ni, |
1141 |
|
|
int *i1, |
1142 |
|
|
int *i2, |
1143 |
|
|
int *nr, |
1144 |
|
|
double *r1, |
1145 |
|
|
double *r2 |
1146 |
johnpye |
888 |
){ |
1147 |
johnpye |
966 |
static double r1last; |
1148 |
|
|
|
1149 |
|
|
asc_assert(*level!=2); // LSODE doesn't give level 2 in our version. |
1150 |
|
|
|
1151 |
johnpye |
888 |
switch(*nerr){ |
1152 |
johnpye |
966 |
case 52: |
1153 |
|
|
if(*nr==2){ |
1154 |
|
|
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Illegal t = %f, not in range (t - hu,t) = (%f,%f)", r1last, *r1, *r2); |
1155 |
|
|
return; |
1156 |
|
|
}else if(*nr==1){ |
1157 |
|
|
r1last = *r1; |
1158 |
|
|
return; |
1159 |
|
|
} break; |
1160 |
johnpye |
888 |
case 204: |
1161 |
johnpye |
1049 |
if(*nr==0 && *ni==0)return; |
1162 |
johnpye |
966 |
if(*nr==2){ |
1163 |
|
|
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Error test failed repeatedly or with abs(h)=hmin.\nt=%f and step size h=%f",*r1,*r2); |
1164 |
|
|
return; |
1165 |
|
|
} break; |
1166 |
johnpye |
888 |
case 205: |
1167 |
johnpye |
966 |
if(*nr==2){ |
1168 |
|
|
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Corrector convergence test failed repeatedly or with abs(h)=hmin.\nt=%f and step size h=%f",*r1,*r2); |
1169 |
|
|
return; |
1170 |
|
|
} break; |
1171 |
|
|
case 27: |
1172 |
|
|
if(*nr==1 && *ni==1){ |
1173 |
|
|
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Trouble with INTDY: itask = %d, tout = %f", *i1, *r1); |
1174 |
|
|
return; |
1175 |
|
|
} break; |
1176 |
|
|
} |
1177 |
johnpye |
908 |
|
1178 |
johnpye |
966 |
ERROR_REPORTER_START_NOLINE(ASC_PROG_ERR); |
1179 |
johnpye |
669 |
|
1180 |
johnpye |
966 |
/* note that %.*s means that a string length (integer) and string pointer are being required */ |
1181 |
|
|
FPRINTF(stderr,"LSODE error: (%d) %.*s",*nerr,*nmes,msg); |
1182 |
|
|
if (*ni == 1) { |
1183 |
|
|
FPRINTF(stderr,"\nwhere i1 = %d",*i1); |
1184 |
johnpye |
888 |
} |
1185 |
johnpye |
966 |
if (*ni == 2) { |
1186 |
|
|
FPRINTF(stderr,"\nwhere i1 = %d, i2 = %d",*i1,*i2); |
1187 |
johnpye |
888 |
} |
1188 |
johnpye |
966 |
if (*nr == 1) { |
1189 |
|
|
FPRINTF(stderr,"\nwhere r1 = %.13g", *r1); |
1190 |
|
|
} |
1191 |
|
|
if (*nr == 2) { |
1192 |
|
|
FPRINTF(stderr,"\nwhere r1 = %.13g, r2 = %.13g", *r1,*r2); |
1193 |
|
|
} |
1194 |
johnpye |
888 |
error_reporter_end_flush(); |
1195 |
johnpye |
669 |
} |