1 |
/* ASCEND modelling environment |
2 |
Copyright (C) 2006 Carnegie Mellon University |
3 |
|
4 |
This program is free software; you can redistribute it and/or modify |
5 |
it under the terms of the GNU General Public License as published by |
6 |
the Free Software Foundation; either version 2, or (at your option) |
7 |
any later version. |
8 |
|
9 |
This program is distributed in the hope that it will be useful, |
10 |
but WITHOUT ANY WARRANTY; without even the implied warranty of |
11 |
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
12 |
GNU General Public License for more details. |
13 |
|
14 |
You should have received a copy of the GNU General Public License |
15 |
along with this program; if not, write to the Free Software |
16 |
Foundation, Inc., 59 Temple Place - Suite 330, |
17 |
Boston, MA 02111-1307, USA. |
18 |
*//** |
19 |
@file |
20 |
Access to the IDA integrator for ASCEND. IDA is a DAE solver that comes |
21 |
as part of the GPL-licensed SUNDIALS solver package from LLNL. |
22 |
@see http://www.llnl.gov/casc/sundials/ |
23 |
*//* |
24 |
by John Pye, May 2006 |
25 |
*/ |
26 |
|
27 |
/* |
28 |
Be careful with the following. This file requires both the 'ida.h' from |
29 |
SUNDIALS as well as the 'ida.h' from ASCEND. Make sure that we're getting |
30 |
both of these; if you get problems check your build tool for the paths being |
31 |
passed to the C preprocessor. |
32 |
*/ |
33 |
|
34 |
/* standard includes */ |
35 |
#include <signal.h> |
36 |
|
37 |
/* ASCEND includes */ |
38 |
#include "ida.h" |
39 |
#include <utilities/error.h> |
40 |
#include <utilities/ascConfig.h> |
41 |
#include <utilities/ascSignal.h> |
42 |
#include <utilities/ascPanic.h> |
43 |
#include <compiler/instance_enum.h> |
44 |
#include "var.h" |
45 |
#include "rel.h" |
46 |
#include "discrete.h" |
47 |
#include "conditional.h" |
48 |
#include "logrel.h" |
49 |
#include "bnd.h" |
50 |
#include "linsol.h" |
51 |
#include "linsolqr.h" |
52 |
#include "slv_common.h" |
53 |
#include "slv_client.h" |
54 |
#include "relman.h" |
55 |
|
56 |
/* SUNDIALS includes */ |
57 |
#ifdef ASC_WITH_IDA |
58 |
# include <sundials/sundials_config.h> |
59 |
# include <sundials/sundials_dense.h> |
60 |
# include <ida/ida.h> |
61 |
# include <nvector/nvector_serial.h> |
62 |
# include <ida/ida_spgmr.h> |
63 |
# include <ida/ida_spbcgs.h> |
64 |
# include <ida/ida_sptfqmr.h> |
65 |
# include <ida/ida_dense.h> |
66 |
# ifndef IDA_SUCCESS |
67 |
# error "Failed to include SUNDIALS IDA header file" |
68 |
# endif |
69 |
#endif |
70 |
|
71 |
/* |
72 |
for cases where we don't have SUNDIALS_VERSION_MINOR defined, guess version 2.2 |
73 |
*/ |
74 |
#ifndef SUNDIALS_VERSION_MINOR |
75 |
# ifdef __GNUC__ |
76 |
# warning "GUESSING SUNDIALS VERSION 2.2" |
77 |
# endif |
78 |
# define SUNDIALS_VERSION_MINOR 2 |
79 |
#endif |
80 |
#ifndef SUNDIALS_VERSION_MAJOR |
81 |
# define SUNDIALS_VERSION_MAJOR 2 |
82 |
#endif |
83 |
|
84 |
/* check that we've got what we expect now */ |
85 |
#ifndef ASC_IDA_H |
86 |
# error "Failed to include ASCEND IDA header file" |
87 |
#endif |
88 |
|
89 |
#define FEX_DEBUG |
90 |
#define JEX_DEBUG |
91 |
|
92 |
const IntegratorInternals integrator_ida_internals = { |
93 |
integrator_ida_create |
94 |
,integrator_ida_params_default |
95 |
,integrator_analyse_dae /* note, this routine is back in integrator.c */ |
96 |
,integrator_ida_solve |
97 |
,integrator_ida_free |
98 |
,INTEG_IDA |
99 |
,"IDA" |
100 |
}; |
101 |
|
102 |
/** |
103 |
Struct containing any stuff that IDA needs that doesn't fit into the |
104 |
common IntegratorSystem struct. |
105 |
*/ |
106 |
typedef struct{ |
107 |
struct rel_relation **rellist; /**< NULL terminated list of rels */ |
108 |
struct var_variable **varlist; /**< NULL terminated list of vars. ONLY USED FOR DEBUGGING -- get rid of it! */ |
109 |
struct bnd_boundary **bndlist; /**< NULL-terminated list of boundaries, for use in the root-finding code */ |
110 |
int nrels; |
111 |
int safeeval; /**< whether to pass the 'safe' flag to relman_eval */ |
112 |
} IntegratorIdaData; |
113 |
|
114 |
/*------------------------------------------------------------- |
115 |
FORWARD DECLS |
116 |
*/ |
117 |
/* residual function forward declaration */ |
118 |
int integrator_ida_fex(realtype tt, N_Vector yy, N_Vector yp, N_Vector rr, void *res_data); |
119 |
|
120 |
int integrator_ida_jvex(realtype tt, N_Vector yy, N_Vector yp, N_Vector rr |
121 |
, N_Vector v, N_Vector Jv, realtype c_j |
122 |
, void *jac_data, N_Vector tmp1, N_Vector tmp2 |
123 |
); |
124 |
|
125 |
/* error handler forward declaration */ |
126 |
void integrator_ida_error(int error_code |
127 |
, const char *module, const char *function |
128 |
, char *msg, void *eh_data |
129 |
); |
130 |
|
131 |
int integrator_ida_djex(long int Neq, realtype tt |
132 |
, N_Vector yy, N_Vector yp, N_Vector rr |
133 |
, realtype c_j, void *jac_data, DenseMat Jac |
134 |
, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3 |
135 |
); |
136 |
|
137 |
/*------------------------------------------------------------- |
138 |
SETUP/TEARDOWN ROUTINES |
139 |
*/ |
140 |
void integrator_ida_create(IntegratorSystem *blsys){ |
141 |
CONSOLE_DEBUG("ALLOCATING IDA ENGINE DATA"); |
142 |
IntegratorIdaData *enginedata; |
143 |
enginedata = ASC_NEW(IntegratorIdaData); |
144 |
enginedata->rellist = NULL; |
145 |
enginedata->varlist = NULL; |
146 |
enginedata->safeeval = 0; |
147 |
blsys->enginedata = (void *)enginedata; |
148 |
integrator_ida_params_default(blsys); |
149 |
} |
150 |
|
151 |
void integrator_ida_free(void *enginedata){ |
152 |
CONSOLE_DEBUG("DELETING IDA ENGINE DATA"); |
153 |
IntegratorIdaData *d = (IntegratorIdaData *)enginedata; |
154 |
/* note, we don't own the rellist, so don't need to free it */ |
155 |
ASC_FREE(d); |
156 |
} |
157 |
|
158 |
IntegratorIdaData *integrator_ida_enginedata(IntegratorSystem *blsys){ |
159 |
IntegratorIdaData *d; |
160 |
assert(blsys!=NULL); |
161 |
assert(blsys->enginedata!=NULL); |
162 |
assert(blsys->engine==INTEG_IDA); |
163 |
d = ((IntegratorIdaData *)(blsys->enginedata)); |
164 |
return d; |
165 |
} |
166 |
|
167 |
/*------------------------------------------------------------- |
168 |
PARAMETERS FOR IDA |
169 |
*/ |
170 |
|
171 |
enum ida_parameters{ |
172 |
IDA_PARAM_LINSOLVER |
173 |
,IDA_PARAM_MAXL |
174 |
,IDA_PARAM_AUTODIFF |
175 |
,IDA_PARAM_CALCIC |
176 |
,IDA_PARAM_SAFEEVAL |
177 |
,IDA_PARAM_RTOL |
178 |
,IDA_PARAM_ATOL |
179 |
,IDA_PARAM_ATOLVECT |
180 |
,IDA_PARAM_GSMODIFIED |
181 |
,IDA_PARAMS_SIZE |
182 |
}; |
183 |
|
184 |
/** |
185 |
Here the full set of parameters is defined, along with upper/lower bounds, |
186 |
etc. The values are stuck into the blsys->params structure. |
187 |
|
188 |
@return 0 on success |
189 |
*/ |
190 |
int integrator_ida_params_default(IntegratorSystem *blsys){ |
191 |
asc_assert(blsys!=NULL); |
192 |
asc_assert(blsys->engine==INTEG_IDA); |
193 |
slv_parameters_t *p; |
194 |
p = &(blsys->params); |
195 |
|
196 |
slv_destroy_parms(p); |
197 |
|
198 |
if(p->parms==NULL){ |
199 |
CONSOLE_DEBUG("params NULL"); |
200 |
p->parms = ASC_NEW_ARRAY(struct slv_parameter, IDA_PARAMS_SIZE); |
201 |
if(p->parms==NULL)return -1; |
202 |
p->dynamic_parms = 1; |
203 |
}else{ |
204 |
CONSOLE_DEBUG("params not NULL"); |
205 |
} |
206 |
|
207 |
/* reset the number of parameters to zero so that we can check it at the end */ |
208 |
p->num_parms = 0; |
209 |
|
210 |
slv_param_bool(p,IDA_PARAM_AUTODIFF |
211 |
,(SlvParameterInitBool){{"autodiff" |
212 |
,"Use auto-diff?",1 |
213 |
,"Use automatic differentiation of expressions (1) or use numerical derivatives (0)" |
214 |
}, TRUE} |
215 |
); |
216 |
|
217 |
slv_param_bool(p,IDA_PARAM_CALCIC |
218 |
,(SlvParameterInitBool){{"calcic" |
219 |
,"Calculate initial conditions?",1 |
220 |
,"Use IDA to calculate initial conditions (1) or else assume that the model will already be solved for this case (0)" |
221 |
}, FALSE} |
222 |
); |
223 |
|
224 |
slv_param_bool(p,IDA_PARAM_SAFEEVAL |
225 |
,(SlvParameterInitBool){{"safeeval" |
226 |
,"Use safe evaluation?",1 |
227 |
,"Use 'safe' function evaluation routines (TRUE) or allow ASCEND to " |
228 |
"throw SIGFPE errors which will then halt integration." |
229 |
}, FALSE} |
230 |
); |
231 |
|
232 |
|
233 |
slv_param_bool(p,IDA_PARAM_ATOLVECT |
234 |
,(SlvParameterInitBool){{"atolvect" |
235 |
,"Use 'ode_atol' values as specified?",1 |
236 |
,"If TRUE, values of 'ode_atol' are taken from your model and used " |
237 |
" in the integration. If FALSE, a scalar absolute tolerance value" |
238 |
" is shared by all variables. See IDA manual, section 5.5.1" |
239 |
}, TRUE } |
240 |
); |
241 |
|
242 |
slv_param_real(p,IDA_PARAM_ATOL |
243 |
,(SlvParameterInitReal){{"atol" |
244 |
,"Scalar absolute error tolerance",1 |
245 |
,"Value of the scalar absolute error tolerance. See also 'atolvect'." |
246 |
" See IDA manual, section 5.5.1" |
247 |
}, 1e-5, DBL_MIN, DBL_MAX } |
248 |
); |
249 |
|
250 |
slv_param_real(p,IDA_PARAM_RTOL |
251 |
,(SlvParameterInitReal){{"rtol" |
252 |
,"Scalar relative error tolerance",1 |
253 |
,"Value of the scalar relative error tolerance. (Note that for IDA," |
254 |
" it's not possible to set per-variable relative tolerances as it is" |
255 |
" with LSODE)." |
256 |
" See IDA manual, section 5.5.1" |
257 |
}, 1e-4, 0, DBL_MAX } |
258 |
); |
259 |
|
260 |
slv_param_char(p,IDA_PARAM_LINSOLVER |
261 |
,(SlvParameterInitChar){{"linsolver" |
262 |
,"Linear solver",1 |
263 |
,"See IDA manual, section 5.5.3." |
264 |
}, "SPGMR"}, (char *[]){"DENSE","BAND","SPGMR","SPBCG","SPTFQMR",NULL} |
265 |
); |
266 |
|
267 |
slv_param_int(p,IDA_PARAM_MAXL |
268 |
,(SlvParameterInitInt){{"maxl" |
269 |
,"Maximum Krylov dimension",0 |
270 |
,"The maximum dimension of Krylov space used by the linear solver" |
271 |
" (for SPGMR, SPBCG, SPTFQMR) with IDA. See IDA manual section 5.5." |
272 |
" The default of 0 results in IDA using its internal default, which" |
273 |
" is currently a value of 5." |
274 |
}, 0, 0, 20 } |
275 |
); |
276 |
|
277 |
slv_param_bool(p,IDA_PARAM_GSMODIFIED |
278 |
,(SlvParameterInitBool){{"gsmodified" |
279 |
,"Use modified Gram-Schmidt orthogonalisation in SPGMR?",2 |
280 |
,"TRUE = GS_MODIFIED, FALSE = GS_CLASSICAL. See IDA manual section" |
281 |
" 5.5.6.6. Only applies when linsolve=SPGMR is selected." |
282 |
}, TRUE} |
283 |
); |
284 |
|
285 |
asc_assert(p->num_parms == IDA_PARAMS_SIZE); |
286 |
|
287 |
CONSOLE_DEBUG("Created %d params", p->num_parms); |
288 |
|
289 |
return 0; |
290 |
} |
291 |
|
292 |
/*------------------------------------------------------------- |
293 |
MAIN IDA SOLVER ROUTINE, see IDA manual, sec 5.4, p. 27 ff. |
294 |
*/ |
295 |
|
296 |
/* return 1 on success */ |
297 |
int integrator_ida_solve( |
298 |
IntegratorSystem *blsys |
299 |
, unsigned long start_index |
300 |
, unsigned long finish_index |
301 |
){ |
302 |
void *ida_mem; |
303 |
int size, flag, t_index; |
304 |
realtype t0, reltol, abstol, t, tret, tout1; |
305 |
N_Vector y0, yp0, abstolvect, ypret, yret; |
306 |
IntegratorIdaData *enginedata; |
307 |
char *linsolver; |
308 |
int maxl; |
309 |
|
310 |
CONSOLE_DEBUG("STARTING IDA..."); |
311 |
|
312 |
enginedata = integrator_ida_enginedata(blsys); |
313 |
|
314 |
enginedata->safeeval = SLV_PARAM_BOOL(&(blsys->params),IDA_PARAM_SAFEEVAL); |
315 |
CONSOLE_DEBUG("safeeval = %d",enginedata->safeeval); |
316 |
|
317 |
/* store reference to list of relations (in enginedata) */ |
318 |
enginedata->nrels = slv_get_num_solvers_rels(blsys->system); |
319 |
enginedata->rellist = slv_get_solvers_rel_list(blsys->system); |
320 |
enginedata->varlist = slv_get_solvers_var_list(blsys->system); |
321 |
enginedata->bndlist = slv_get_solvers_bnd_list(blsys->system); |
322 |
|
323 |
CONSOLE_DEBUG("Number of relations: %d",enginedata->nrels); |
324 |
CONSOLE_DEBUG("Number of dependent vars: %ld",blsys->n_y); |
325 |
size = blsys->n_y; |
326 |
|
327 |
if(enginedata->nrels!=size){ |
328 |
ERROR_REPORTER_HERE(ASC_USER_ERROR,"Integration problem is not square (%d rels, %d vars)", enginedata->nrels, size); |
329 |
return 0; /* failure */ |
330 |
} |
331 |
|
332 |
/* retrieve initial values from the system */ |
333 |
|
334 |
/** @TODO fix this, the starting time != first sample */ |
335 |
t0 = integrator_get_t(blsys); |
336 |
CONSOLE_DEBUG("RETRIEVED t0 = %f",t0); |
337 |
|
338 |
CONSOLE_DEBUG("RETRIEVING y0"); |
339 |
|
340 |
y0 = N_VNew_Serial(size); |
341 |
integrator_get_y(blsys,NV_DATA_S(y0)); |
342 |
|
343 |
CONSOLE_DEBUG("RETRIEVING yp0"); |
344 |
|
345 |
yp0 = N_VNew_Serial(size); |
346 |
integrator_get_ydot(blsys,NV_DATA_S(yp0)); |
347 |
|
348 |
N_VPrint_Serial(yp0); |
349 |
CONSOLE_DEBUG("yp0 is at %p",&yp0); |
350 |
|
351 |
/* create IDA object */ |
352 |
ida_mem = IDACreate(); |
353 |
|
354 |
/* relative error tolerance */ |
355 |
reltol = SLV_PARAM_REAL(&(blsys->params),IDA_PARAM_RTOL); |
356 |
CONSOLE_DEBUG("rtol = %8.2e",reltol); |
357 |
|
358 |
/* allocate internal memory */ |
359 |
if(SLV_PARAM_BOOL(&(blsys->params),IDA_PARAM_ATOLVECT)){ |
360 |
/* vector of absolute tolerances */ |
361 |
CONSOLE_DEBUG("USING VECTOR OF ATOL VALUES"); |
362 |
abstolvect = N_VNew_Serial(size); |
363 |
integrator_get_atol(blsys,NV_DATA_S(abstolvect)); |
364 |
|
365 |
flag = IDAMalloc(ida_mem, &integrator_ida_fex, t0, y0, yp0, IDA_SV, reltol, abstolvect); |
366 |
|
367 |
N_VDestroy_Serial(abstolvect); |
368 |
}else{ |
369 |
/* scalar absolute tolerance (one value for all) */ |
370 |
CONSOLE_DEBUG("USING SCALAR ATOL VALUE = %8.2e",abstol); |
371 |
abstol = SLV_PARAM_REAL(&(blsys->params),IDA_PARAM_ATOL); |
372 |
flag = IDAMalloc(ida_mem, &integrator_ida_fex, t0, y0, yp0, IDA_SS, reltol, &abstol); |
373 |
} |
374 |
|
375 |
if(flag==IDA_MEM_NULL){ |
376 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"ida_mem is NULL"); |
377 |
return 0; |
378 |
}else if(flag==IDA_MEM_FAIL){ |
379 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Unable to allocate memory (IDAMalloc)"); |
380 |
return 0; |
381 |
}else if(flag==IDA_ILL_INPUT){ |
382 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Invalid input to IDAMalloc"); |
383 |
return 0; |
384 |
}/* else success */ |
385 |
|
386 |
/* set optional inputs... */ |
387 |
IDASetErrHandlerFn(ida_mem, &integrator_ida_error, (void *)blsys); |
388 |
IDASetRdata(ida_mem, (void *)blsys); |
389 |
IDASetMaxStep(ida_mem, integrator_get_maxstep(blsys)); |
390 |
IDASetInitStep(ida_mem, integrator_get_stepzero(blsys)); |
391 |
IDASetMaxNumSteps(ida_mem, integrator_get_maxsubsteps(blsys)); |
392 |
if(integrator_get_minstep(blsys)>0){ |
393 |
ERROR_REPORTER_HERE(ASC_PROG_NOTE,"IDA does not support minstep (ignored)\n"); |
394 |
} |
395 |
/* there's no capability for setting *minimum* step size in IDA */ |
396 |
|
397 |
|
398 |
/* attach linear solver module, using the default value of maxl */ |
399 |
linsolver = SLV_PARAM_CHAR(&(blsys->params),IDA_PARAM_LINSOLVER); |
400 |
CONSOLE_DEBUG("ASSIGNING LINEAR SOLVER '%s'",linsolver); |
401 |
if(strcmp(linsolver,"DENSE")==0){ |
402 |
CONSOLE_DEBUG("DENSE DIRECT SOLVER, size = %d",size); |
403 |
flag = IDADense(ida_mem, size); |
404 |
switch(flag){ |
405 |
case IDADENSE_SUCCESS: break; |
406 |
case IDADENSE_MEM_NULL: ERROR_REPORTER_HERE(ASC_PROG_ERR,"ida_mem is NULL"); return 0; |
407 |
case IDADENSE_ILL_INPUT: ERROR_REPORTER_HERE(ASC_PROG_ERR,"IDADENSE is not compatible with current nvector module"); return 0; |
408 |
case IDADENSE_MEM_FAIL: ERROR_REPORTER_HERE(ASC_PROG_ERR,"Memory allocation failed for IDADENSE"); return 0; |
409 |
default: ERROR_REPORTER_HERE(ASC_PROG_ERR,"bad return"); return 0; |
410 |
} |
411 |
|
412 |
if(SLV_PARAM_BOOL(&(blsys->params),IDA_PARAM_AUTODIFF)){ |
413 |
CONSOLE_DEBUG("USING AUTODIFF"); |
414 |
flag = IDADenseSetJacFn(ida_mem, &integrator_ida_djex, (void *)blsys); |
415 |
switch(flag){ |
416 |
case IDADENSE_SUCCESS: break; |
417 |
default: ERROR_REPORTER_HERE(ASC_PROG_ERR,"Failed IDADenseSetJacFn"); return 0; |
418 |
} |
419 |
}else{ |
420 |
CONSOLE_DEBUG("USING NUMERICAL DIFF"); |
421 |
} |
422 |
}else{ |
423 |
/* remaining methods are all SPILS */ |
424 |
CONSOLE_DEBUG("IDA SPILS"); |
425 |
|
426 |
maxl = SLV_PARAM_INT(&(blsys->params),IDA_PARAM_MAXL); |
427 |
CONSOLE_DEBUG("maxl = %d",maxl); |
428 |
|
429 |
if(strcmp(linsolver,"SPGMR")==0){ |
430 |
CONSOLE_DEBUG("IDA SPGMR"); |
431 |
flag = IDASpgmr(ida_mem, maxl); /* 0 means use the default max Krylov dimension of 5 */ |
432 |
}else if(strcmp(linsolver,"SPBCG")==0){ |
433 |
CONSOLE_DEBUG("IDA SPBCG"); |
434 |
flag = IDASpbcg(ida_mem, maxl); |
435 |
}else if(strcmp(linsolver,"SPTFQMR")==0){ |
436 |
CONSOLE_DEBUG("IDA SPTFQMR"); |
437 |
flag = IDASptfqmr(ida_mem,maxl); |
438 |
}else{ |
439 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Unknown IDA linear solver choice '%s'",linsolver); |
440 |
return 0; |
441 |
} |
442 |
|
443 |
if(flag==IDASPILS_MEM_NULL){ |
444 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"ida_mem is NULL"); |
445 |
return 0; |
446 |
}else if(flag==IDASPILS_MEM_FAIL){ |
447 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Unable to allocate memory (IDASpgmr)"); |
448 |
return 0; |
449 |
}/* else success */ |
450 |
|
451 |
/* assign the J*v function */ |
452 |
if(SLV_PARAM_BOOL(&(blsys->params),IDA_PARAM_AUTODIFF)){ |
453 |
CONSOLE_DEBUG("USING AUTODIFF"); |
454 |
flag = IDASpilsSetJacTimesVecFn(ida_mem, &integrator_ida_jvex, (void *)blsys); |
455 |
if(flag==IDASPILS_MEM_NULL){ |
456 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"ida_mem is NULL"); |
457 |
return 0; |
458 |
}else if(flag==IDASPILS_LMEM_NULL){ |
459 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"IDASPILS linear solver has not been initialized"); |
460 |
return 0; |
461 |
}/* else success */ |
462 |
}else{ |
463 |
CONSOLE_DEBUG("USING NUMERICAL DIFF"); |
464 |
} |
465 |
|
466 |
if(strcmp(linsolver,"SPGMR")==0){ |
467 |
/* select Gram-Schmidt orthogonalisation */ |
468 |
if(SLV_PARAM_BOOL(&(blsys->params),IDA_PARAM_GSMODIFIED)){ |
469 |
CONSOLE_DEBUG("USING MODIFIED GS"); |
470 |
flag = IDASpilsSetGSType(ida_mem,MODIFIED_GS); |
471 |
if(flag!=IDASPILS_SUCCESS){ |
472 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Failed to set GS_MODIFIED"); |
473 |
return 0; |
474 |
} |
475 |
}else{ |
476 |
CONSOLE_DEBUG("USING CLASSICAL GS"); |
477 |
flag = IDASpilsSetGSType(ida_mem,CLASSICAL_GS); |
478 |
if(flag!=IDASPILS_SUCCESS){ |
479 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Failed to set GS_MODIFIED"); |
480 |
return 0; |
481 |
} |
482 |
} |
483 |
} |
484 |
} |
485 |
|
486 |
/* set linear solver optional inputs... |
487 |
...nothing here at the moment... |
488 |
*/ |
489 |
|
490 |
/* calculate initial conditions */ |
491 |
if(SLV_PARAM_BOOL(&(blsys->params),IDA_PARAM_CALCIC)){ |
492 |
CONSOLE_DEBUG("Solving initial conditions (given derivatives)"); |
493 |
|
494 |
blsys->currentstep=0; |
495 |
t_index=start_index; |
496 |
tout1 = samplelist_get(blsys->samples, t_index); |
497 |
|
498 |
CONSOLE_DEBUG("SOLVING INITIAL CONDITIONS IDACalcIC (tout1 = %f)", tout1); |
499 |
|
500 |
/* correct initial values, given derivatives */ |
501 |
# if SUNDIALS_VERSION_MAJOR==2 && SUNDIALS_VERSION_MINOR==3 |
502 |
/* note the new API from version 2.3 and onwards */ |
503 |
flag = IDACalcIC(ida_mem, IDA_Y_INIT, tout1); |
504 |
# else |
505 |
flag = IDACalcIC(ida_mem, t0, y0, yp0, IDA_Y_INIT, tout1); |
506 |
# endif |
507 |
|
508 |
if(flag!=IDA_SUCCESS){ |
509 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Unable to solve initial values (IDACalcIC)"); |
510 |
return 0; |
511 |
}/* else success */ |
512 |
|
513 |
CONSOLE_DEBUG("INITIAL CONDITIONS SOLVED :-)"); |
514 |
}else{ |
515 |
CONSOLE_DEBUG("Not solving initial conditions (see IDA parameter 'calcic')"); |
516 |
} |
517 |
|
518 |
/* optionally, specify ROO-FINDING PROBLEM */ |
519 |
|
520 |
/* -- set up the IntegratorReporter */ |
521 |
integrator_output_init(blsys); |
522 |
|
523 |
/* -- store the initial values of all the stuff */ |
524 |
integrator_output_write(blsys); |
525 |
integrator_output_write_obs(blsys); |
526 |
|
527 |
/* specify where the returned values should be stored */ |
528 |
yret = y0; |
529 |
ypret = yp0; |
530 |
|
531 |
/* advance solution in time, return values as yret and derivatives as ypret */ |
532 |
blsys->currentstep=1; |
533 |
for(t_index=start_index+1;t_index <= finish_index;++t_index, ++blsys->currentstep){ |
534 |
t = samplelist_get(blsys->samples, t_index); |
535 |
t0 = integrator_get_t(blsys); |
536 |
asc_assert(t > t0); |
537 |
|
538 |
CONSOLE_DEBUG("Integratoring from t0 = %f to t = %f", t0, t); |
539 |
|
540 |
flag = IDASolve(ida_mem, t, &tret, yret, ypret, IDA_NORMAL); |
541 |
|
542 |
/* pass the values of everything back to the compiler */ |
543 |
integrator_set_t(blsys, (double)tret); |
544 |
integrator_set_y(blsys, NV_DATA_S(yret)); |
545 |
integrator_set_ydot(blsys, NV_DATA_S(ypret)); |
546 |
|
547 |
if(flag<0){ |
548 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Failed to solve t = %f (IDASolve), error %d", t, flag); |
549 |
break; |
550 |
} |
551 |
|
552 |
/* -- do something so that blsys knows the values of tret, yret and ypret */ |
553 |
|
554 |
/* -- store the current values of all the stuff */ |
555 |
integrator_output_write(blsys); |
556 |
integrator_output_write_obs(blsys); |
557 |
} |
558 |
|
559 |
/* -- close the IntegratorReporter */ |
560 |
integrator_output_close(blsys); |
561 |
|
562 |
if(flag < 0){ |
563 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Solving aborted while attempting t = %f", t); |
564 |
return 0; |
565 |
} |
566 |
|
567 |
/* get optional outputs */ |
568 |
|
569 |
/* free solution memory */ |
570 |
N_VDestroy_Serial(yret); |
571 |
N_VDestroy_Serial(ypret); |
572 |
|
573 |
/* free solver memory */ |
574 |
IDAFree(ida_mem); |
575 |
|
576 |
/* all done */ |
577 |
return 1; |
578 |
} |
579 |
|
580 |
/*-------------------------------------------------- |
581 |
RESIDUALS AND JACOBIAN |
582 |
*/ |
583 |
/** |
584 |
Function to evaluate system residuals, in the form required for IDA. |
585 |
|
586 |
Given tt, yy and yp, we need to evaluate and return rr. |
587 |
|
588 |
@param tt current value of indep variable (time) |
589 |
@param yy current values of dependent variable vector |
590 |
@param yp current values of derivatives of dependent variables |
591 |
@param rr the output residual vector (is we're returning data to) |
592 |
@param res_data pointer to our stuff (blsys in this case). |
593 |
|
594 |
@return 0 on success, positive on recoverable error, and |
595 |
negative on unrecoverable error. |
596 |
*/ |
597 |
int integrator_ida_fex(realtype tt, N_Vector yy, N_Vector yp, N_Vector rr, void *res_data){ |
598 |
IntegratorSystem *blsys; |
599 |
IntegratorIdaData *enginedata; |
600 |
int i, calc_ok, is_error; |
601 |
struct rel_relation** relptr; |
602 |
double resid; |
603 |
char *relname; |
604 |
#ifdef FEX_DEBUG |
605 |
char *varname; |
606 |
#endif |
607 |
|
608 |
blsys = (IntegratorSystem *)res_data; |
609 |
enginedata = integrator_ida_enginedata(blsys); |
610 |
|
611 |
#ifdef FEX_DEBUG |
612 |
/* fprintf(stderr,"\n\n"); */ |
613 |
CONSOLE_DEBUG("EVALUTE RESIDUALS..."); |
614 |
#endif |
615 |
|
616 |
/* pass the values of everything back to the compiler */ |
617 |
integrator_set_t(blsys, (double)tt); |
618 |
integrator_set_y(blsys, NV_DATA_S(yy)); |
619 |
integrator_set_ydot(blsys, NV_DATA_S(yp)); |
620 |
|
621 |
if(NV_LENGTH_S(rr)!=enginedata->nrels){ |
622 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Invalid residuals nrels!=length(rr)"); |
623 |
return -1; /* unrecoverable */ |
624 |
} |
625 |
|
626 |
/** |
627 |
@TODO does this function (fex) do bounds checking already? |
628 |
*/ |
629 |
|
630 |
/* evaluate each residual in the rellist */ |
631 |
is_error = 0; |
632 |
relptr = enginedata->rellist; |
633 |
|
634 |
if(enginedata->safeeval){ |
635 |
Asc_SignalHandlerPush(SIGFPE,SIG_IGN); |
636 |
}else{ |
637 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"SETTING TO CATCH SIGFPE..."); |
638 |
Asc_SignalHandlerPushDefault(SIGFPE); |
639 |
} |
640 |
|
641 |
if (setjmp(g_fpe_env)==0) { |
642 |
for(i=0, relptr = enginedata->rellist; |
643 |
i< enginedata->nrels && relptr != NULL; |
644 |
++i, ++relptr |
645 |
){ |
646 |
resid = relman_eval(*relptr, &calc_ok, enginedata->safeeval); |
647 |
|
648 |
NV_Ith_S(rr,i) = resid; |
649 |
if(!calc_ok){ |
650 |
relname = rel_make_name(blsys->system, *relptr); |
651 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Calculation error in rel '%s'",relname); |
652 |
ASC_FREE(relname); |
653 |
/* presumable some output already made? */ |
654 |
is_error = 1; |
655 |
}/*else{ |
656 |
CONSOLE_DEBUG("Calc OK"); |
657 |
}*/ |
658 |
} |
659 |
}else{ |
660 |
relname = rel_make_name(blsys->system, *relptr); |
661 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Floating point error (SIGFPE) in rel '%s'",relname); |
662 |
ASC_FREE(relname); |
663 |
is_error = 1; |
664 |
} |
665 |
|
666 |
if(enginedata->safeeval){ |
667 |
Asc_SignalHandlerPop(SIGFPE,SIG_IGN); |
668 |
}else{ |
669 |
Asc_SignalHandlerPopDefault(SIGFPE); |
670 |
} |
671 |
|
672 |
#ifdef FEX_DEBUG |
673 |
/* output residuals to console */ |
674 |
CONSOLE_DEBUG("RESIDUAL OUTPUT"); |
675 |
fprintf(stderr,"index\t%20s\t%20s\t%s\n","y","ydot","resid"); |
676 |
for(i=0; i<blsys->n_y; ++i){ |
677 |
varname = var_make_name(blsys->system,blsys->y[i]); |
678 |
fprintf(stderr,"%d\t%10s=%10f\t",i,varname,NV_Ith_S(yy,i)); |
679 |
if(blsys->ydot[i]){ |
680 |
varname = var_make_name(blsys->system,blsys->ydot[i]); |
681 |
fprintf(stderr,"%10s=%10f\t",varname,NV_Ith_S(yp,i)); |
682 |
}else{ |
683 |
fprintf(stderr,"diff(%4s)=%10f\t",varname,NV_Ith_S(yp,i)); |
684 |
} |
685 |
ASC_FREE(varname); |
686 |
relname = rel_make_name(blsys->system,enginedata->rellist[i]); |
687 |
fprintf(stderr,"'%s'=%f\n",relname,NV_Ith_S(rr,i)); |
688 |
} |
689 |
#endif |
690 |
|
691 |
if(is_error){ |
692 |
return 1; |
693 |
} |
694 |
|
695 |
#ifdef FEX_DEBUG |
696 |
CONSOLE_DEBUG("RESIDUAL OK"); |
697 |
#endif |
698 |
return 0; |
699 |
} |
700 |
|
701 |
/** |
702 |
Dense Jacobian evaluation. Only suitable for small problems! |
703 |
*/ |
704 |
int integrator_ida_djex(long int Neq, realtype tt |
705 |
, N_Vector yy, N_Vector yp, N_Vector rr |
706 |
, realtype c_j, void *jac_data, DenseMat Jac |
707 |
, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3 |
708 |
){ |
709 |
IntegratorSystem *blsys; |
710 |
IntegratorIdaData *enginedata; |
711 |
char *relname; |
712 |
#ifdef JEX_DEBUG |
713 |
char *varname; |
714 |
#endif |
715 |
int status; |
716 |
struct rel_relation **relptr; |
717 |
int i; |
718 |
var_filter_t filter = {VAR_SVAR, VAR_SVAR}; |
719 |
double *derivatives; |
720 |
int *variables; |
721 |
int count, j; |
722 |
long var_yindex; |
723 |
|
724 |
blsys = (IntegratorSystem *)jac_data; |
725 |
enginedata = integrator_ida_enginedata(blsys); |
726 |
|
727 |
/* allocate space for returns from relman_diff2: we *should* be able to use 'tmp1' and 'tmp2' here... */ |
728 |
variables = ASC_NEW_ARRAY(int, NV_LENGTH_S(yy) * 2); |
729 |
derivatives = ASC_NEW_ARRAY(double, NV_LENGTH_S(yy) * 2); |
730 |
|
731 |
/* pass the values of everything back to the compiler */ |
732 |
integrator_set_t(blsys, (double)tt); |
733 |
integrator_set_y(blsys, NV_DATA_S(yy)); |
734 |
integrator_set_ydot(blsys, NV_DATA_S(yp)); |
735 |
|
736 |
#ifdef JEX_DEBUG |
737 |
/* print vars */ |
738 |
for(i=0; i < blsys->n_y; ++i){ |
739 |
varname = var_make_name(blsys->system, blsys->y[i]); |
740 |
CONSOLE_DEBUG("%s = %f = %f",varname,NV_Ith_S(yy,i),var_value(blsys->y[i])); |
741 |
ASC_FREE(varname); |
742 |
} |
743 |
|
744 |
/* print derivatives */ |
745 |
for(i=0; i < blsys->n_y; ++i){ |
746 |
if(blsys->ydot[i]){ |
747 |
varname = var_make_name(blsys->system, blsys->ydot[i]); |
748 |
CONSOLE_DEBUG("%s = %f =%f",varname,NV_Ith_S(yp,i),var_value(blsys->ydot[i])); |
749 |
ASC_FREE(varname); |
750 |
}else{ |
751 |
varname = var_make_name(blsys->system, blsys->y[i]); |
752 |
CONSOLE_DEBUG("diff(%s) = %f",varname,NV_Ith_S(yp,i)); |
753 |
ASC_FREE(varname); |
754 |
} |
755 |
} |
756 |
|
757 |
/* print step size */ |
758 |
CONSOLE_DEBUG("<c_j> = %f",c_j); |
759 |
#endif |
760 |
|
761 |
/* build up the dense jacobian matrix... */ |
762 |
status = 0; |
763 |
for(i=0, relptr = enginedata->rellist; |
764 |
i< enginedata->nrels && relptr != NULL; |
765 |
++i, ++relptr |
766 |
){ |
767 |
/* get derivatives for this particular relation */ |
768 |
status = relman_diff2(*relptr, &filter, derivatives, variables, &count, enginedata->safeeval); |
769 |
|
770 |
if(status){ |
771 |
relname = rel_make_name(blsys->system, *relptr); |
772 |
CONSOLE_DEBUG("ERROR calculating derivatives for relation '%s'",relname); |
773 |
ASC_FREE(relname); |
774 |
break; |
775 |
} |
776 |
|
777 |
/* output what's going on here ... */ |
778 |
#ifdef JEX_DEBUG |
779 |
relname = rel_make_name(blsys->system, *relptr); |
780 |
CONSOLE_DEBUG("RELATION %d '%s'",i,relname); |
781 |
fprintf(stderr,"%d: '%s': ",i,relname); |
782 |
ASC_FREE(relname); |
783 |
for(j=0;j<count;++j){ |
784 |
varname = var_make_name(blsys->system, enginedata->varlist[variables[j]]); |
785 |
var_yindex = blsys->y_id[variables[j]]; |
786 |
if(var_yindex >=0){ |
787 |
fprintf(stderr," var[%d]='%s'=y[%ld]",variables[j],varname,var_yindex); |
788 |
}else{ |
789 |
fprintf(stderr," var[%d]='%s'=ydot[%ld]",variables[j],varname,-var_yindex-1); |
790 |
} |
791 |
ASC_FREE(varname); |
792 |
} |
793 |
fprintf(stderr,"\n"); |
794 |
#endif |
795 |
/* insert values into the Jacobian row in appropriate spots (can assume Jac starts with zeros -- IDA manual) */ |
796 |
for(j=0; j < count; ++j){ |
797 |
var_yindex = blsys->y_id[variables[j]]; |
798 |
ASC_ASSERT_RANGE(var_yindex, -Jac->N, Jac->N); |
799 |
if(var_yindex >= 0){ |
800 |
asc_assert(blsys->y[var_yindex]==enginedata->varlist[variables[j]]); |
801 |
DENSE_ELEM(Jac,i,var_yindex) += derivatives[j]; |
802 |
}else{ |
803 |
asc_assert(blsys->ydot[-var_yindex-1]==enginedata->varlist[variables[j]]); |
804 |
DENSE_ELEM(Jac,i,-var_yindex-1) += derivatives[j] * c_j; |
805 |
} |
806 |
} |
807 |
} |
808 |
|
809 |
#ifdef JEX_DEBUG |
810 |
CONSOLE_DEBUG("PRINTING JAC"); |
811 |
fprintf(stderr,"\t"); |
812 |
for(j=0; j < blsys->n_y; ++j){ |
813 |
if(j)fprintf(stderr,"\t"); |
814 |
varname = var_make_name(blsys->system,blsys->y[j]); |
815 |
fprintf(stderr,"%11s",varname); |
816 |
ASC_FREE(varname); |
817 |
} |
818 |
fprintf(stderr,"\n"); |
819 |
for(i=0; i < enginedata->nrels; ++i){ |
820 |
relname = rel_make_name(blsys->system, enginedata->rellist[i]); |
821 |
fprintf(stderr,"%s\t",relname); |
822 |
ASC_FREE(relname); |
823 |
|
824 |
for(j=0; j < blsys->n_y; ++j){ |
825 |
if(j)fprintf(stderr,"\t"); |
826 |
fprintf(stderr,"%11.2e",DENSE_ELEM(Jac,i,j)); |
827 |
} |
828 |
fprintf(stderr,"\n"); |
829 |
} |
830 |
#endif |
831 |
|
832 |
if(status){ |
833 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"There were derivative evaluation errors in the dense jacobian"); |
834 |
return 1; |
835 |
} |
836 |
|
837 |
#ifdef FEX_DEBUG |
838 |
CONSOLE_DEBUG("DJEX RETURNING 0"); |
839 |
#endif |
840 |
return 0; |
841 |
} |
842 |
|
843 |
/** |
844 |
Function to evaluate the product J*v, in the form required for IDA (see IDASpilsSetJacTimesVecFn) |
845 |
|
846 |
Given tt, yy, yp, rr and v, we need to evaluate and return Jv. |
847 |
|
848 |
@param tt current value of the independent variable (time, t) |
849 |
@param yy current value of the dependent variable vector, y(t). |
850 |
@param yp current value of y'(t). |
851 |
@param rr current value of the residual vector F(t, y, y'). |
852 |
@param v the vector by which the Jacobian must be multiplied to the right. |
853 |
@param Jv the output vector computed |
854 |
@param c_j the scalar in the system Jacobian, proportional to the inverse of the step size ($ \alpha$ in Eq. (3.5) ). |
855 |
@param jac_data pointer to our stuff (blsys in this case, passed into IDA via IDASp*SetJacTimesVecFn.) |
856 |
@param tmp1 @see tmp2 |
857 |
@param tmp2 (as well as tmp1) pointers to memory allocated for variables of type N_Vector for use here as temporary storage or work space. |
858 |
@return 0 on success |
859 |
*/ |
860 |
int integrator_ida_jvex(realtype tt, N_Vector yy, N_Vector yp, N_Vector rr |
861 |
, N_Vector v, N_Vector Jv, realtype c_j |
862 |
, void *jac_data, N_Vector tmp1, N_Vector tmp2 |
863 |
){ |
864 |
IntegratorSystem *blsys; |
865 |
IntegratorIdaData *enginedata; |
866 |
int i, j, is_error=0; |
867 |
struct rel_relation** relptr; |
868 |
char *relname; |
869 |
int status; |
870 |
double Jv_i; |
871 |
long var_yindex; |
872 |
|
873 |
int *variables; |
874 |
double *derivatives; |
875 |
var_filter_t filter; |
876 |
int count; |
877 |
|
878 |
#ifdef JEX_DEBUG |
879 |
CONSOLE_DEBUG("EVALUATING JACOBIAN..."); |
880 |
#endif |
881 |
|
882 |
blsys = (IntegratorSystem *)jac_data; |
883 |
enginedata = integrator_ida_enginedata(blsys); |
884 |
|
885 |
/* pass the values of everything back to the compiler */ |
886 |
integrator_set_t(blsys, (double)tt); |
887 |
integrator_set_y(blsys, NV_DATA_S(yy)); |
888 |
integrator_set_ydot(blsys, NV_DATA_S(yp)); |
889 |
/* no real use for residuals (rr) here, I don't think? */ |
890 |
|
891 |
/* allocate space for returns from relman_diff2: we *should* be able to use 'tmp1' and 'tmp2' here... */ |
892 |
variables = ASC_NEW_ARRAY(int, NV_LENGTH_S(yy) * 2); |
893 |
derivatives = ASC_NEW_ARRAY(double, NV_LENGTH_S(yy) * 2); |
894 |
|
895 |
/* evaluate the derivatives... */ |
896 |
/* J = dG_dy = dF_dy + alpha * dF_dyp */ |
897 |
|
898 |
filter.matchbits = VAR_SVAR; |
899 |
filter.matchvalue = VAR_SVAR; |
900 |
|
901 |
Asc_SignalHandlerPushDefault(SIGFPE); |
902 |
if (setjmp(g_fpe_env)==0) { |
903 |
for(i=0, relptr = enginedata->rellist; |
904 |
i< enginedata->nrels && relptr != NULL; |
905 |
++i, ++relptr |
906 |
){ |
907 |
/* get derivatives for this particular relation */ |
908 |
status = relman_diff2(*relptr, &filter, derivatives, variables, &count, enginedata->safeeval); |
909 |
/* CONSOLE_DEBUG("Got derivatives against %d matching variables", count); */ |
910 |
|
911 |
if(status){ |
912 |
relname = rel_make_name(blsys->system, *relptr); |
913 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Calculation error in rel '%s'",relname); |
914 |
ASC_FREE(relname); |
915 |
is_error = 1; |
916 |
break; |
917 |
} |
918 |
|
919 |
/* |
920 |
Now we have the derivatives wrt each alg/diff variable in the |
921 |
present equation. variables[] points into the varlist. need |
922 |
a mapping from the varlist to the y and ydot lists. |
923 |
*/ |
924 |
|
925 |
Jv_i = 0; |
926 |
for(j=0; j < count; ++j){ |
927 |
/* CONSOLE_DEBUG("j = %d, variables[j] = %d, n_y = %ld", j, variables[j], blsys->n_y); |
928 |
varname = var_make_name(blsys->system, enginedata->varlist[variables[j]]); |
929 |
if(varname){ |
930 |
CONSOLE_DEBUG("Variable %d '%s' derivative = %f", variables[j],varname,derivatives[j]); |
931 |
ASC_FREE(varname); |
932 |
}else{ |
933 |
CONSOLE_DEBUG("Variable %d (UNKNOWN!): derivative = %f",variables[j],derivatives[j]); |
934 |
} |
935 |
*/ |
936 |
|
937 |
var_yindex = blsys->y_id[variables[j]]; |
938 |
CONSOLE_DEBUG("j = %d: variables[j] = %d, y_id = %ld",j,variables[j],var_yindex); |
939 |
ASC_ASSERT_RANGE(-var_yindex-1, -NV_LENGTH_S(v),NV_LENGTH_S(v)); |
940 |
|
941 |
if(var_yindex >= 0){ |
942 |
#ifdef JEX_DEBUG |
943 |
asc_assert(blsys->y[var_yindex]==enginedata->varlist[variables[j]]); |
944 |
fprintf(stderr,"Jv[%d] += %f (dF[%d]/dy[%ld] = %f, v[%ld] = %f)\n", i |
945 |
, derivatives[j] * NV_Ith_S(v,var_yindex) |
946 |
, i, var_yindex, derivatives[j] |
947 |
, var_yindex, NV_Ith_S(v,var_yindex) |
948 |
); |
949 |
#endif |
950 |
Jv_i += derivatives[j] * NV_Ith_S(v,var_yindex); |
951 |
}else{ |
952 |
#ifdef JEX_DEBUG |
953 |
fprintf(stderr,"Jv[%d] += %f (dF[%d]/dydot[%ld] = %f, v[%ld] = %f)\n", i |
954 |
, derivatives[j] * NV_Ith_S(v,-var_yindex-1) |
955 |
, i, -var_yindex-1, derivatives[j] |
956 |
, -var_yindex-1, NV_Ith_S(v,-var_yindex-1) |
957 |
); |
958 |
#endif |
959 |
asc_assert(blsys->ydot[-var_yindex-1]==enginedata->varlist[variables[j]]); |
960 |
Jv_i += derivatives[j] * NV_Ith_S(v,-var_yindex-1) * c_j; |
961 |
} |
962 |
} |
963 |
|
964 |
NV_Ith_S(Jv,i) = Jv_i; |
965 |
#ifdef JEX_DEBUG |
966 |
relname = rel_make_name(blsys->system, *relptr); |
967 |
CONSOLE_DEBUG("'%s': Jv[%d] = %f", relname, i, NV_Ith_S(Jv,i)); |
968 |
ASC_FREE(relname); |
969 |
return 1; |
970 |
#endif |
971 |
} |
972 |
}else{ |
973 |
relname = rel_make_name(blsys->system, *relptr); |
974 |
ERROR_REPORTER_HERE(ASC_PROG_ERR,"Floating point error (SIGFPE) in rel '%s'",relname); |
975 |
ASC_FREE(relname); |
976 |
is_error = 1; |
977 |
} |
978 |
Asc_SignalHandlerPopDefault(SIGFPE); |
979 |
|
980 |
if(is_error){ |
981 |
CONSOLE_DEBUG("SOME ERRORS FOUND IN EVALUATION"); |
982 |
return 1; |
983 |
} |
984 |
return 0; |
985 |
} |
986 |
|
987 |
/*---------------------------------------------- |
988 |
ERROR REPORTING |
989 |
*/ |
990 |
/** |
991 |
Error message reporter function to be passed to IDA. All error messages |
992 |
will trigger a call to this function, so we should find everything |
993 |
appearing on the console (in the case of Tcl/Tk) or in the errors/warnings |
994 |
panel (in the case of PyGTK). |
995 |
*/ |
996 |
void integrator_ida_error(int error_code |
997 |
, const char *module, const char *function |
998 |
, char *msg, void *eh_data |
999 |
){ |
1000 |
IntegratorSystem *blsys; |
1001 |
error_severity_t sev; |
1002 |
|
1003 |
/* cast back the IntegratorSystem, just in case we need it */ |
1004 |
blsys = (IntegratorSystem *)eh_data; |
1005 |
|
1006 |
/* severity depends on the sign of the error_code value */ |
1007 |
if(error_code <= 0){ |
1008 |
sev = ASC_PROG_ERR; |
1009 |
}else{ |
1010 |
sev = ASC_PROG_WARNING; |
1011 |
} |
1012 |
|
1013 |
/* use our all-purpose error reporting to get stuff back to the GUI */ |
1014 |
error_reporter(sev,module,0,function,"%s (error %d)",msg,error_code); |
1015 |
} |