/[ascend]/trunk/base/generic/solver/ida.c
ViewVC logotype

Contents of /trunk/base/generic/solver/ida.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 930 - (show annotations) (download) (as text)
Wed Nov 22 13:09:09 2006 UTC (17 years, 9 months ago) by johnpye
File MIME type: text/x-csrc
File size: 17918 byte(s)
solve.py raises ImportError if 'browser' object is not available (right approach?)
Freeing some variables in dsgsat2.a4c.
Returing Py_None from extpy routine in 'browser' object not defined (eg during non-GUI unit testing)
Error reporting from extpy import handler (ongoing)
Timeout in versioncheck (when server unavailable)
A little more tinkering with IDA.
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 <compiler/instance_enum.h>
43 #include "var.h"
44 #include "rel.h"
45 #include "discrete.h"
46 #include "conditional.h"
47 #include "logrel.h"
48 #include "bnd.h"
49 #include "linsol.h"
50 #include "linsolqr.h"
51 #include "slv_common.h"
52 #include "slv_client.h"
53 #include "relman.h"
54
55 /* SUNDIALS includes */
56 #ifdef ASC_WITH_IDA
57 # include <sundials/sundials_config.h>
58 # include <ida/ida.h>
59 # include <nvector/nvector_serial.h>
60 # include <ida/ida_spgmr.h>
61 # ifndef IDA_SUCCESS
62 # error "Failed to include SUNDIALS IDA header file"
63 # endif
64 #endif
65
66 /*
67 for the benefit of build tools that didn't sniff the SUNDIALS version, we
68 assume version 2.2.x (and thence possible errors).
69 */
70 #ifndef SUNDIALS_VERSION_MINOR
71 # ifdef __GNUC__
72 # warning "GUESSING SUNDIALS VERSION 2.2"
73 # endif
74 # define SUNDIALS_VERSION_MINOR 2
75 #endif
76 #ifndef SUNDIALS_VERSION_MAJOR
77 # define SUNDIALS_VERSION_MAJOR 2
78 #endif
79
80 /* check that we've got what we expect now */
81 #ifndef ASC_IDA_H
82 # error "Failed to include ASCEND IDA header file"
83 #endif
84
85 /**
86 Struct containing any stuff that IDA needs that doesn't fit into the
87 common IntegratorSystem struct.
88 */
89 typedef struct{
90 struct rel_relation **rellist; /**< NULL terminated list of rels */
91 struct var_variable **varlist; /**< NULL terminated list of vars. ONLY USED FOR DEBUGGING -- get rid of it! */
92 int nrels;
93 int safeeval; /**< whether to pass the 'safe' flag to relman_eval */
94 } IntegratorIdaData;
95
96 /*-------------------------------------------------------------
97 FORWARD DECLS
98 */
99 /* residual function forward declaration */
100 int integrator_ida_fex(realtype tt, N_Vector yy, N_Vector yp, N_Vector rr, void *res_data);
101
102 int integrator_ida_jvex(realtype tt, N_Vector yy, N_Vector yp, N_Vector rr
103 , N_Vector v, N_Vector Jv, realtype c_j
104 , void *jac_data, N_Vector tmp1, N_Vector tmp2
105 );
106
107 /* error handler forward declaration */
108 void integrator_ida_error(int error_code
109 , const char *module, const char *function
110 , char *msg, void *eh_data
111 );
112
113 /*-------------------------------------------------------------
114 SETUP/TEARDOWN ROUTINES
115 */
116 void integrator_ida_create(IntegratorSystem *blsys){
117 CONSOLE_DEBUG("ALLOCATING IDA ENGINE DATA");
118 IntegratorIdaData *enginedata;
119 enginedata = ASC_NEW(IntegratorIdaData);
120 enginedata->rellist = NULL;
121 enginedata->varlist = NULL;
122 enginedata->safeeval = 1;
123 blsys->enginedata = (void *)enginedata;
124 }
125
126 void integrator_ida_free(void *enginedata){
127 CONSOLE_DEBUG("DELETING IDA ENGINE DATA");
128 IntegratorIdaData *d = (IntegratorIdaData *)enginedata;
129 /* note, we don't own the rellist, so don't need to free it */
130 ASC_FREE(d);
131 }
132
133 IntegratorIdaData *integrator_ida_enginedata(IntegratorSystem *blsys){
134 IntegratorIdaData *d;
135 assert(blsys!=NULL);
136 assert(blsys->enginedata!=NULL);
137 assert(blsys->engine==INTEG_IDA);
138 d = ((IntegratorIdaData *)(blsys->enginedata));
139 assert(d->safeeval = 1);
140 return d;
141 }
142
143 /*-------------------------------------------------------------
144 MAIN IDA SOLVER ROUTINE, see IDA manual, sec 5.4, p. 27 ff.
145 */
146
147 /* return 1 on success */
148 int integrator_ida_solve(
149 IntegratorSystem *blsys
150 , unsigned long start_index
151 , unsigned long finish_index
152 ){
153 void *ida_mem;
154 int size, flag, t_index;
155 realtype t0, reltol, t, tret, tout1;
156 N_Vector y0, yp0, abstol, ypret, yret;
157 IntegratorIdaData *enginedata;
158
159 CONSOLE_DEBUG("STARTING IDA...");
160
161 enginedata = integrator_ida_enginedata(blsys);
162 CONSOLE_DEBUG("safeeval = %d",enginedata->safeeval);
163
164 /* store reference to list of relations (in enginedata) */
165 enginedata->nrels = slv_get_num_solvers_rels(blsys->system);
166 enginedata->rellist = slv_get_solvers_rel_list(blsys->system);
167 enginedata->varlist = slv_get_solvers_var_list(blsys->system);
168 CONSOLE_DEBUG("Number of relations: %d",enginedata->nrels);
169 CONSOLE_DEBUG("Number of dependent vars: %ld",blsys->n_y);
170 size = blsys->n_y;
171
172 if(enginedata->nrels!=size){
173 ERROR_REPORTER_HERE(ASC_USER_ERROR,"Integration problem is not square (%d rels, %d vars)", enginedata->nrels, size);
174 return 0; /* failure */
175 }
176
177 CONSOLE_DEBUG("RETRIEVING t0");
178
179 /* retrieve initial values from the system */
180 t0 = samplelist_get(blsys->samples,start_index);
181
182 CONSOLE_DEBUG("RETRIEVING y0");
183
184 y0 = N_VNew_Serial(size);
185 integrator_get_y(blsys,NV_DATA_S(y0));
186
187 CONSOLE_DEBUG("RETRIEVING yp0");
188
189 yp0 = N_VNew_Serial(size);
190 integrator_get_ydot(blsys,NV_DATA_S(yp0));
191
192 N_VPrint_Serial(yp0);
193 CONSOLE_DEBUG("yp0 is at %p",&yp0);
194
195 /* create IDA object */
196 ida_mem = IDACreate();
197
198 /* retrieve the absolute tolerance values for each variable */
199 abstol = N_VNew_Serial(size);
200 N_VConst(0.1,abstol); /** @TODO fill in the abstol values from the model */
201 reltol = 0.001;
202
203 /* allocate internal memory */
204 flag = IDAMalloc(ida_mem, &integrator_ida_fex, t0, y0, yp0, IDA_SV, reltol, abstol);
205 if(flag==IDA_MEM_NULL){
206 ERROR_REPORTER_HERE(ASC_PROG_ERR,"ida_mem is NULL");
207 return 0;
208 }else if(flag==IDA_MEM_FAIL){
209 ERROR_REPORTER_HERE(ASC_PROG_ERR,"Unable to allocate memory (IDAMalloc)");
210 return 0;
211 }else if(flag==IDA_ILL_INPUT){
212 ERROR_REPORTER_HERE(ASC_PROG_ERR,"Invalid input to IDAMalloc");
213 return 0;
214 }/* else success */
215
216 /* set optional inputs... */
217 IDASetErrHandlerFn(ida_mem, &integrator_ida_error, (void *)blsys);
218 IDASetRdata(ida_mem, (void *)blsys);
219 IDASetMaxStep(ida_mem, integrator_get_maxstep(blsys));
220 IDASetInitStep(ida_mem, integrator_get_stepzero(blsys));
221 IDASetMaxNumSteps(ida_mem, integrator_get_maxsubsteps(blsys));
222 /* there's no capability for setting *minimum* step size in IDA */
223
224 CONSOLE_DEBUG("ASSIGNING LINEAR SOLVER");
225
226 /* attach linear solver module, using the default value of maxl */
227 flag = IDASpgmr(ida_mem, 0);
228 if(flag==IDASPILS_MEM_NULL){
229 ERROR_REPORTER_HERE(ASC_PROG_ERR,"ida_mem is NULL");
230 return 0;
231 }else if(flag==IDASPILS_MEM_FAIL){
232 ERROR_REPORTER_HERE(ASC_PROG_ERR,"Unable to allocate memory (IDASpgmr)");
233 return 0;
234 }/* else success */
235
236 /* assign the J*v function */
237 #if 0
238 flag = IDASpilsSetJacTimesVecFn(ida_mem, &integrator_ida_jvex, (void *)blsys);
239 if(flag==IDASPILS_MEM_NULL){
240 ERROR_REPORTER_HERE(ASC_PROG_ERR,"ida_mem is NULL");
241 return 0;
242 }else if(flag==IDASPILS_LMEM_NULL){
243 ERROR_REPORTER_HERE(ASC_PROG_ERR,"IDASPILS linear solver has not been initialized");
244 return 0;
245 }/* else success */
246 #endif
247
248 /* set linear solver optional inputs...
249
250 ...nothing here at the moment...
251
252 */
253
254 /* correct initial values, given derivatives */
255 blsys->currentstep=0;
256 t_index=start_index+1;
257 tout1 = samplelist_get(blsys->samples, t_index);
258
259 #if SUNDIALS_VERSION_MAJOR==2 && SUNDIALS_VERSION_MINOR==3
260 flag = IDACalcIC(ida_mem, IDA_Y_INIT, tout1);
261 #else
262 flag = IDACalcIC(ida_mem, t0, y0, yp0, IDA_Y_INIT, tout1);
263 #endif
264
265 if(flag!=IDA_SUCCESS){
266 ERROR_REPORTER_HERE(ASC_PROG_ERR,"Unable to solve initial values (IDACalcIC)");
267 return 0;
268 }/* else success */
269
270 CONSOLE_DEBUG("INITIAL CONDITIONS SOLVED :-)");
271
272 /* optionally, specify ROO-FINDING PROBLEM */
273
274 /* -- set up the IntegratorReporter */
275 integrator_output_init(blsys);
276
277 /* -- store the initial values of all the stuff */
278 integrator_output_write(blsys);
279 integrator_output_write_obs(blsys);
280
281 /* specify where the returned values should be stored */
282 yret = y0;
283 ypret = yp0;
284
285 /* advance solution in time, return values as yret and derivatives as ypret */
286 blsys->currentstep=1;
287 for(t_index=start_index+1;t_index <= finish_index;++t_index, ++blsys->currentstep){
288 t = samplelist_get(blsys->samples, t_index);
289
290 CONSOLE_DEBUG("SOLVING UP TO t = %f", t);
291
292 flag = IDASolve(ida_mem, t, &tret, yret, ypret, IDA_NORMAL);
293
294 /* pass the values of everything back to the compiler */
295 integrator_set_t(blsys, (double)tret);
296 integrator_set_y(blsys, NV_DATA_S(yret));
297 integrator_set_ydot(blsys, NV_DATA_S(ypret));
298
299 if(flag<0){
300 ERROR_REPORTER_HERE(ASC_PROG_ERR,"Failed to solve t = %f (IDASolve), error %d", t, flag);
301 break;
302 }
303
304 /* -- do something so that blsys knows the values of tret, yret and ypret */
305
306 /* -- store the current values of all the stuff */
307 integrator_output_write(blsys);
308 integrator_output_write_obs(blsys);
309 }
310
311 /* -- close the IntegratorReporter */
312 integrator_output_close(blsys);
313
314 if(flag < 0){
315 ERROR_REPORTER_HERE(ASC_PROG_ERR,"Solving aborted while attempting t = %f", t);
316 return 0;
317 }
318
319 /* get optional outputs */
320
321 /* free solution memory */
322 N_VDestroy_Serial(yret);
323 N_VDestroy_Serial(ypret);
324
325 /* free solver memory */
326 IDAFree(ida_mem);
327
328 /* all done */
329 return 1;
330 }
331
332 /*--------------------------------------------------
333 RESIDUALS AND JACOBIAN
334 */
335 /**
336 Function to evaluate system residuals, in the form required for IDA.
337
338 Given tt, yy and yp, we need to evaluate and return rr.
339
340 @param tt current value of indep variable (time)
341 @param yy current values of dependent variable vector
342 @param yp current values of derivatives of dependent variables
343 @param rr the output residual vector (is we're returning data to)
344 @param res_data pointer to our stuff (blsys in this case).
345
346 @return 0 on success, positive on recoverable error, and
347 negative on unrecoverable error.
348 */
349 int integrator_ida_fex(realtype tt, N_Vector yy, N_Vector yp, N_Vector rr, void *res_data){
350 IntegratorSystem *blsys;
351 IntegratorIdaData *enginedata;
352 int i, calc_ok, is_error;
353 struct rel_relation** relptr;
354 double resid;
355 char *relname;
356
357 blsys = (IntegratorSystem *)res_data;
358 enginedata = integrator_ida_enginedata(blsys);
359
360 /* fprintf(stderr,"\n\n"); */
361 /* CONSOLE_DEBUG("ABOUT TO EVALUTE RESIDUALS..."); */
362
363 /* pass the values of everything back to the compiler */
364 integrator_set_t(blsys, (double)tt);
365 integrator_set_y(blsys, NV_DATA_S(yy));
366 integrator_set_ydot(blsys, NV_DATA_S(yp));
367
368 /* revaluate the system residuals using the new data */
369 is_error = 0;
370 relptr = enginedata->rellist;
371
372 /* CONSOLE_DEBUG("IDA requests residuals of length %lu",NV_LENGTH_S(rr)); */
373 if(NV_LENGTH_S(rr)!=enginedata->nrels){
374 ERROR_REPORTER_HERE(ASC_PROG_ERR,"Invalid residuals nrels!=length(rr)");
375 return -1; /* unrecoverable */
376 }
377
378 Asc_SignalHandlerPush(SIGFPE,SIG_IGN);
379 if (setjmp(g_fpe_env)==0) {
380 for(i=0, relptr = enginedata->rellist;
381 i< enginedata->nrels && relptr != NULL;
382 ++i, ++relptr
383 ){
384 resid = relman_eval(*relptr, &calc_ok, enginedata->safeeval);
385
386 relname = rel_make_name(blsys->system, *relptr);
387 /* if(calc_ok){
388 CONSOLE_DEBUG("residual[%d:\"%s\"] = %f",i,relname,resid);
389 }else{
390 CONSOLE_DEBUG("residual[%d:\"%s\"] = %f (ERROR)",i,relname,resid);
391 }*/
392 ASC_FREE(relname);
393
394 NV_Ith_S(rr,i) = resid;
395 if(!calc_ok){
396 /* presumable some output already made? */
397 is_error = 1;
398 }
399 }
400 }else{
401 CONSOLE_DEBUG("FLOATING POINT ERROR WITH i=%d",i);
402 }
403 Asc_SignalHandlerPop(SIGFPE,SIG_IGN);
404
405 if(is_error)CONSOLE_DEBUG("SOME ERRORS FOUND IN EVALUATION");
406 return is_error;
407 }
408
409 /**
410 Function to evaluate the product J*v, in the form required for IDA (see IDASpilsSetJacTimesVecFn)
411
412 Given tt, yy, yp, rr and v, we need to evaluate and return Jv.
413
414 @param tt current value of the independent variable (time, t)
415 @param yy current value of the dependent variable vector, y(t).
416 @param yp current value of y'(t).
417 @param rr current value of the residual vector F(t, y, y').
418 @param v the vector by which the Jacobian must be multiplied to the right.
419 @param Jv the output vector computed
420 @param c_j the scalar in the system Jacobian, proportional to the inverse of the step size ($ \alpha$ in Eq. (3.5) ).
421 @param jac_data pointer to our stuff (blsys in this case, passed into IDA via IDASp*SetJacTimesVecFn.)
422 @param tmp1 @see tmp2
423 @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.
424 @return 0 on success
425 */
426 int integrator_ida_jvex(realtype tt, N_Vector yy, N_Vector yp, N_Vector rr
427 , N_Vector v, N_Vector Jv, realtype c_j
428 , void *jac_data, N_Vector tmp1, N_Vector tmp2
429 ){
430 IntegratorSystem *blsys;
431 IntegratorIdaData *enginedata;
432 int i, j, is_error=0;
433 struct rel_relation** relptr;
434 char *relname, *varname;
435 int status;
436 double Jv_i;
437 int var_yindex;
438
439 int *variables;
440 double *derivatives;
441 var_filter_t filter;
442 int count;
443
444 fprintf(stderr,"\n--------------\n");
445 CONSOLE_DEBUG("EVALUTING JACOBIAN...");
446
447 blsys = (IntegratorSystem *)jac_data;
448 enginedata = integrator_ida_enginedata(blsys);
449
450 /* pass the values of everything back to the compiler */
451 integrator_set_t(blsys, (double)tt);
452 integrator_set_y(blsys, NV_DATA_S(yy));
453 integrator_set_ydot(blsys, NV_DATA_S(yp));
454 /* no real use for residuals (rr) here, I don't think? */
455
456 /* allocate space for returns from relman_diff2: we *should* be able to use 'tmp1' and 'tmp2' here... */
457 variables = ASC_NEW_ARRAY(int, NV_LENGTH_S(yy) * 2);
458 derivatives = ASC_NEW_ARRAY(double, NV_LENGTH_S(yy) * 2);
459
460 /* evaluate the derivatives... */
461 /* J = dG_dy = dF_dy + alpha * dF_dyp */
462
463 filter.matchbits = VAR_SVAR;
464 filter.matchvalue = VAR_SVAR;
465
466 CONSOLE_DEBUG("PRINTING VALUES OF 'v' VECTOR (length %ld)",NV_LENGTH_S(v));
467 for(i=0; i<NV_LENGTH_S(v); ++i){
468 CONSOLE_DEBUG("v[%d] = %f",i,NV_Ith_S(v,i));
469 }
470
471 Asc_SignalHandlerPush(SIGFPE,SIG_IGN);
472 if (setjmp(g_fpe_env)==0) {
473 for(i=0, relptr = enginedata->rellist;
474 i< enginedata->nrels && relptr != NULL;
475 ++i, ++relptr
476 ){
477 fprintf(stderr,"\n");
478 relname = rel_make_name(blsys->system, *relptr);
479 CONSOLE_DEBUG("RELATION %d '%s'",i,relname);
480 ASC_FREE(relname);
481
482 /* get derivatives for this particular relation */
483 status = relman_diff2(*relptr, &filter, derivatives, variables, &count, enginedata->safeeval);
484 CONSOLE_DEBUG("Got derivatives against %d matching variables", count);
485
486 for(j=0;j<count;++j){
487 varname = var_make_name(blsys->system, enginedata->varlist[variables[j]]);
488 CONSOLE_DEBUG("derivatives[%d] = %f (variable %d, '%s')",j,derivatives[j],variables[j],varname);
489 ASC_FREE(varname);
490 }
491
492 if(!status){
493 CONSOLE_DEBUG("Derivatives for relation %d OK",i);
494 }else{
495 CONSOLE_DEBUG("ERROR calculating derivatives for relation %d",i);
496 break;
497 }
498
499 /*
500 Now we have the derivatives wrt each alg/diff variable in the
501 present equation. variables[] points into the varlist. need
502 a mapping from the varlist to the y and ydot lists.
503 */
504
505 Jv_i = 0;
506 for(j=0; j < count; ++j){
507 /* CONSOLE_DEBUG("j = %d, variables[j] = %d, n_y = %ld", j, variables[j], blsys->n_y); */
508 varname = var_make_name(blsys->system, enginedata->varlist[variables[j]]);
509 if(varname){
510 CONSOLE_DEBUG("Variable %d '%s' derivative = %f", variables[j],varname,derivatives[j]);
511 ASC_FREE(varname);
512 }else{
513 CONSOLE_DEBUG("Variable %d (UNKNOWN!): derivative = %f",variables[j],derivatives[j]);
514 }
515
516 var_yindex = blsys->y_id[variables[j]];
517 CONSOLE_DEBUG("j = %d: variables[j] = %d, y_id = %d",j,variables[j],var_yindex);
518
519 if(var_yindex >= 0){
520 CONSOLE_DEBUG("j = %d: algebraic, deriv[j] = %f, v[%d] = %f",j,derivatives[j], var_yindex, NV_Ith_S(v,var_yindex));
521 Jv_i += derivatives[j] * NV_Ith_S(v,var_yindex);
522 }else{
523 var_yindex = -var_yindex-1;
524 CONSOLE_DEBUG("j = %d: differential, deriv[j] = %f, v[%d] = %f",j,derivatives[j], var_yindex, NV_Ith_S(v,var_yindex));
525 Jv_i += derivatives[j] * NV_Ith_S(v,var_yindex) / c_j;
526 }
527 }
528
529 NV_Ith_S(Jv,i) = Jv_i;
530 CONSOLE_DEBUG("(J*v)[%d] = %f", i, Jv_i);
531
532 if(status){
533 /* presumably some error_reporter will already have been made*/
534 is_error = 1;
535 }
536 }
537 }else{
538 CONSOLE_DEBUG("FLOATING POINT ERROR WITH i=%d",i);
539 }
540 Asc_SignalHandlerPop(SIGFPE,SIG_IGN);
541
542 if(is_error)CONSOLE_DEBUG("SOME ERRORS FOUND IN EVALUATION");
543
544
545
546 return is_error;
547 }
548
549 /*----------------------------------------------
550 ERROR REPORTING
551 */
552 /**
553 Error message reporter function to be passed to IDA. All error messages
554 will trigger a call to this function, so we should find everything
555 appearing on the console (in the case of Tcl/Tk) or in the errors/warnings
556 panel (in the case of PyGTK).
557 */
558 void integrator_ida_error(int error_code
559 , const char *module, const char *function
560 , char *msg, void *eh_data
561 ){
562 IntegratorSystem *blsys;
563 error_severity_t sev;
564
565 /* cast back the IntegratorSystem, just in case we need it */
566 blsys = (IntegratorSystem *)eh_data;
567
568 /* severity depends on the sign of the error_code value */
569 if(error_code <= 0){
570 sev = ASC_PROG_ERR;
571 }else{
572 sev = ASC_PROG_WARNING;
573 }
574
575 /* use our all-purpose error reporting to get stuff back to the GUI */
576 error_reporter(sev,module,0,function,"%s (error %d)",msg,error_code);
577 }

john.pye@anu.edu.au
ViewVC Help
Powered by ViewVC 1.1.22