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

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