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