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