3984
|
1 /* |
|
2 |
|
3 Copyright (C) 2002 John W. Eaton |
|
4 |
|
5 This file is part of Octave. |
|
6 |
|
7 Octave is free software; you can redistribute it and/or modify it |
|
8 under the terms of the GNU General Public License as published by the |
|
9 Free Software Foundation; either version 2, or (at your option) any |
|
10 later version. |
|
11 |
|
12 Octave is distributed in the hope that it will be useful, but WITHOUT |
|
13 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or |
|
14 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License |
|
15 for more details. |
|
16 |
|
17 You should have received a copy of the GNU General Public License |
|
18 along with Octave; see the file COPYING. If not, write to the Free |
|
19 Software Foundation, 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. |
|
20 |
|
21 */ |
|
22 |
4192
|
23 #if defined (__GNUG__) && defined (USE_PRAGMA_INTERFACE_IMPLEMENTATION) |
3984
|
24 #pragma implementation |
|
25 #endif |
|
26 |
|
27 #ifdef HAVE_CONFIG_H |
|
28 #include <config.h> |
|
29 #endif |
|
30 |
|
31 #include <cfloat> |
|
32 #include <cmath> |
|
33 |
|
34 // For instantiating the Array<Matrix> object. |
|
35 #include "Array.h" |
|
36 #include "Array.cc" |
|
37 |
|
38 #include "ODESSA.h" |
|
39 #include "f77-fcn.h" |
|
40 #include "lo-error.h" |
4051
|
41 #include "lo-sstream.h" |
4153
|
42 #include "quit.h" |
3984
|
43 |
|
44 typedef int (*odessa_fcn_ptr) (int*, const double&, double*, |
3991
|
45 double*, double*); |
3984
|
46 |
|
47 typedef int (*odessa_jac_ptr) (int*, const double&, double*, |
3991
|
48 double*, const int&, const int&, |
|
49 double*, const int&); |
3984
|
50 |
|
51 typedef int (*odessa_dfdp_ptr) (int*, const double&, double*, |
3991
|
52 double*, double*, const int&); |
3984
|
53 |
|
54 |
|
55 extern "C" |
|
56 int F77_FUNC (odessa, ODESSA) (odessa_fcn_ptr, odessa_dfdp_ptr, int*, |
|
57 double*, double*, double&, double&, |
|
58 int&, double&, const double*, int&, |
|
59 int&, int*, double*, int&, int*, int&, |
|
60 odessa_jac_ptr, int&); |
|
61 |
|
62 template class Array<Matrix>; |
|
63 |
|
64 static ODESFunc::ODES_fsub user_fsub; |
|
65 static ODESFunc::ODES_bsub user_bsub; |
|
66 static ODESFunc::ODES_jsub user_jsub; |
|
67 |
|
68 |
|
69 static int |
|
70 odessa_f (int* neq, const double& t, double *state, |
|
71 double *par, double *fval) |
|
72 { |
4180
|
73 BEGIN_INTERRUPT_WITH_EXCEPTIONS; |
|
74 |
3984
|
75 int n = neq[0]; |
|
76 int n_par = neq[1]; |
|
77 |
|
78 // Load the state and parameter arrays as Octave objects |
|
79 |
|
80 ColumnVector tmp_state (n); |
|
81 for (int i = 0; i < n; i++) |
|
82 tmp_state(i) = state[i]; |
|
83 |
|
84 ColumnVector tmp_param (n_par); |
|
85 for (int i = 0; i < n_par; i++) |
|
86 tmp_param(i) = par[i]; |
|
87 |
4141
|
88 ColumnVector tmp_fval = user_fsub (tmp_state, t, tmp_param); |
3984
|
89 |
4152
|
90 if (tmp_fval.length () == 0) |
4153
|
91 octave_jump_to_enclosing_context (); |
4152
|
92 else |
3984
|
93 { |
4152
|
94 for (int i = 0; i < n; i++) |
|
95 fval[i] = tmp_fval(i); |
3984
|
96 } |
|
97 |
4180
|
98 END_INTERRUPT_WITH_EXCEPTIONS; |
|
99 |
3984
|
100 return 0; |
|
101 } |
|
102 |
|
103 static int |
|
104 odessa_j (int* neq, const double& t, double *state, |
4152
|
105 double *par, const int& ml, const int& mu, double *pd, |
|
106 const int& nrowpd) |
3984
|
107 { |
4180
|
108 BEGIN_INTERRUPT_WITH_EXCEPTIONS; |
|
109 |
3984
|
110 int n = neq[0]; |
|
111 int n_par = neq[1]; |
|
112 |
|
113 // Load the state and parameter arrays as Octave objects |
|
114 ColumnVector tmp_state (n); |
|
115 for (int i = 0; i < n; i++) |
|
116 tmp_state(i) = state[i]; |
|
117 |
|
118 ColumnVector tmp_param (n_par); |
|
119 for (int i = 0; i < n_par; i++) |
|
120 tmp_param(i) = par[i]; |
|
121 |
4141
|
122 Matrix tmp_fval = user_jsub (tmp_state, t, tmp_param); |
3984
|
123 |
4152
|
124 if (tmp_fval.length () == 0) |
4153
|
125 octave_jump_to_enclosing_context (); |
4152
|
126 else |
|
127 { |
|
128 for (int j = 0; j < n; j++) |
|
129 for (int i = 0; i < nrowpd; i++) |
|
130 pd[nrowpd*j+i] = tmp_fval(i,j); |
|
131 } |
3984
|
132 |
4180
|
133 END_INTERRUPT_WITH_EXCEPTIONS; |
|
134 |
3984
|
135 return 0; |
|
136 } |
|
137 |
|
138 static int |
|
139 odessa_b (int* neq, const double& t, double *state, |
|
140 double *par, double *dfdp, const int& jpar) |
|
141 |
|
142 { |
4180
|
143 BEGIN_INTERRUPT_WITH_EXCEPTIONS; |
|
144 |
3984
|
145 int n = neq[0]; |
|
146 int n_par = neq[1]; |
|
147 |
|
148 // Load the state and parameter arrays as Octave objects |
|
149 ColumnVector tmp_state (n); |
|
150 for (int i = 0; i < n; i++) |
|
151 tmp_state(i) = state[i]; |
|
152 |
|
153 ColumnVector tmp_param (n_par); |
|
154 for (int i = 0; i < n_par; i++) |
|
155 tmp_param(i) = par[i]; |
|
156 |
4141
|
157 ColumnVector tmp_fval = user_bsub (tmp_state, t, tmp_param, jpar); |
3984
|
158 |
4152
|
159 if (tmp_fval.length () == 0) |
4153
|
160 octave_jump_to_enclosing_context (); |
4152
|
161 else |
|
162 { |
|
163 for (int i = 0; i < n; i++) |
|
164 dfdp[i] = tmp_fval(i); |
|
165 } |
3984
|
166 |
4180
|
167 END_INTERRUPT_WITH_EXCEPTIONS; |
|
168 |
3984
|
169 return 0; |
|
170 } |
|
171 |
|
172 ODESSA::ODESSA (void) : ODES (), ODESSA_options () |
|
173 { |
|
174 neq.resize(2); |
|
175 n = size (); |
|
176 |
|
177 iopt.resize(4); |
|
178 |
|
179 itask = 1; |
|
180 iopt(0) = 0; |
|
181 isopt = 0; |
|
182 iopt(1) = isopt; |
|
183 npar = 0; |
|
184 neq(0) = n; |
|
185 neq(1) = npar; |
|
186 |
|
187 sanity_checked = false; |
|
188 } |
|
189 |
|
190 ODESSA::ODESSA (const ColumnVector& state, double time, ODESFunc& f) |
|
191 : ODES (state, time, f), ODESSA_options () |
|
192 { |
|
193 neq.resize(2); |
|
194 n = size (); |
|
195 |
|
196 iopt.resize(4); |
|
197 itask = 1; |
|
198 iopt(0) = 0; |
|
199 isopt = 0; |
|
200 iopt(1) = isopt; |
|
201 |
|
202 sanity_checked = false; |
|
203 |
|
204 npar = 0; |
|
205 neq(0) = n; |
|
206 neq(1) = npar; |
|
207 |
|
208 y.resize (n, 1, 0.0); |
|
209 } |
|
210 |
|
211 ODESSA::ODESSA (const ColumnVector& state, const ColumnVector& theta, |
|
212 const Matrix& sensitivity_guess, double time, ODESFunc& f) |
|
213 : ODES (state, theta, time, f) |
|
214 { |
|
215 initialized = false; |
|
216 |
|
217 neq.resize(2); |
|
218 n = state.length(); |
|
219 npar = theta.length(); |
|
220 |
|
221 neq(0) = n; |
|
222 neq(1) = npar; |
|
223 |
|
224 sx0 = sensitivity_guess; |
|
225 par.resize (npar); |
|
226 |
|
227 for (int i = 0; i < npar; i++) |
|
228 { |
|
229 par(i) = theta(i); |
|
230 } |
|
231 |
|
232 sanity_checked = false; |
|
233 |
|
234 npar = theta.length (); |
|
235 |
|
236 iopt.resize(4); |
|
237 itask = 1; |
|
238 iopt(0) = 0; |
|
239 isopt = 1; |
|
240 iopt(1) = isopt; |
|
241 |
|
242 y.resize (n, npar+1, 0.0); |
|
243 } |
|
244 |
|
245 void |
|
246 ODESSA::integrate (double tout) |
|
247 { |
|
248 ODESSA_result retval; |
|
249 if (! initialized) |
|
250 { |
|
251 |
|
252 for (int i = 0; i < n; i++) |
|
253 y(i,0) = x(i); |
|
254 |
|
255 if (npar > 0) |
|
256 { |
|
257 for (int j = 0; j < npar; j++) |
|
258 for (int i = 0; i < n; i++) |
|
259 y(i,j+1) = sx0(i,j); |
|
260 } |
|
261 |
|
262 integration_error = false; |
|
263 |
|
264 user_fsub = ODESFunc::fsub_function (); |
|
265 user_bsub = ODESFunc::bsub_function (); |
|
266 user_jsub = ODESFunc::jsub_function (); |
|
267 |
|
268 int idf; |
|
269 |
|
270 if (user_bsub) |
|
271 idf = 1; |
|
272 else |
|
273 idf = 0; |
|
274 |
|
275 iopt(2) = idf; |
|
276 |
|
277 |
|
278 if (restart) |
|
279 { |
|
280 restart = false; |
|
281 istate = 1; |
|
282 } |
|
283 |
|
284 if (integration_method () == "stiff") |
|
285 { |
|
286 if (user_jsub) |
|
287 method_flag = 21; |
|
288 else |
|
289 method_flag = 22; |
|
290 if (isopt) |
|
291 { |
|
292 liw = 21 + n + npar; |
|
293 lrw = 22 + 8*(npar+1)*n + n*n + n; |
|
294 } |
|
295 else |
|
296 { |
|
297 liw = 20 + n; |
|
298 lrw = 22 + 9*n + n*n; |
|
299 } |
|
300 } |
|
301 else |
|
302 { |
|
303 if (isopt) |
|
304 { |
|
305 if (user_jsub) |
|
306 method_flag = 11; |
|
307 else |
|
308 method_flag = 12; |
|
309 liw = 21 + n + npar; |
|
310 lrw = 22 + 15*(npar+1)*n + n*n + n; |
|
311 } |
|
312 else |
|
313 { |
|
314 method_flag = 10; |
|
315 liw = 20 + n; |
|
316 lrw = 22 + 16 * n; |
|
317 } |
|
318 } |
|
319 |
|
320 if (iwork.length () != liw) |
|
321 { |
|
322 iwork.resize (liw); |
|
323 |
|
324 for (int i = 4; i < 10; i++) |
|
325 iwork.elem (i) = 0; |
|
326 } |
|
327 |
|
328 if (rwork.length () != lrw) |
|
329 { |
|
330 rwork.resize (lrw); |
|
331 |
|
332 for (int i = 4; i < 10; i++) |
|
333 rwork.elem (i) = 0.0; |
|
334 } |
|
335 |
|
336 initialized = true; |
|
337 } |
|
338 integration_error = false; |
|
339 |
|
340 // NOTE: this won't work if LSODE passes copies of the state vector. |
|
341 // In that case we have to create a temporary vector object |
|
342 // and copy. |
|
343 |
|
344 |
|
345 if (! sanity_checked) |
|
346 { |
4141
|
347 ColumnVector fval = user_fsub (x, t, theta); |
3984
|
348 |
|
349 if (fval.length () != x.length ()) |
|
350 { |
|
351 (*current_liboctave_error_handler) |
|
352 ("odessa: inconsistent sizes for state and residual vectors"); |
|
353 |
|
354 integration_error = true; |
|
355 return; |
|
356 } |
|
357 |
|
358 sanity_checked = true; |
|
359 } |
|
360 |
|
361 if (stop_time_set) |
|
362 { |
|
363 itask = 4; |
|
364 rwork.elem (0) = stop_time; |
|
365 iopt(0) = 1; |
|
366 } |
|
367 else |
|
368 { |
|
369 itask = 1; |
|
370 } |
|
371 |
|
372 double rel_tol = relative_tolerance (); |
|
373 |
|
374 Array<double> abs_tol = absolute_tolerance (); //note; this should |
|
375 // be a matrix, not a vector |
|
376 |
|
377 int abs_tol_len = abs_tol.length (); |
|
378 |
|
379 int itol; |
|
380 |
|
381 if (abs_tol_len == 1) |
|
382 itol = 1; |
|
383 else if (abs_tol_len == n) |
|
384 itol = 2; |
|
385 else |
|
386 { |
|
387 (*current_liboctave_error_handler) |
3996
|
388 ("odessa: inconsistent sizes for state and absolute tolerance vectors"); |
3984
|
389 |
|
390 integration_error = 1; |
|
391 return; |
|
392 } |
|
393 |
|
394 if (initial_step_size () >= 0.0) |
|
395 { |
|
396 rwork.elem (4) = initial_step_size (); |
|
397 iopt(0) = 1; |
|
398 } |
|
399 |
|
400 if (maximum_step_size () >= 0.0) |
|
401 { |
|
402 rwork.elem (5) = maximum_step_size (); |
|
403 iopt(0) = 1; |
|
404 } |
|
405 |
|
406 if (minimum_step_size () >= 0.0) |
|
407 { |
|
408 rwork.elem (6) = minimum_step_size (); |
|
409 iopt(0) = 1; |
|
410 } |
|
411 |
|
412 if (step_limit () > 0) |
|
413 { |
|
414 iwork.elem (5) = step_limit (); |
|
415 iopt(0) = 1; |
|
416 } |
|
417 |
|
418 py = y.fortran_vec (); |
|
419 piwork = iwork.fortran_vec (); |
|
420 prwork = rwork.fortran_vec (); |
|
421 ppar = par.fortran_vec (); |
|
422 piopt = iopt.fortran_vec (); |
|
423 pneq = neq.fortran_vec (); |
|
424 |
|
425 const double *pabs_tol = abs_tol.fortran_vec (); |
|
426 |
3997
|
427 F77_XFCN (odessa, ODESSA, (odessa_f, odessa_b, pneq, py, ppar, t, |
|
428 tout, itol, rel_tol, pabs_tol, itask, |
|
429 istate, piopt, prwork, lrw, piwork, liw, |
|
430 odessa_j, method_flag)); |
3984
|
431 |
|
432 if (f77_exception_encountered) |
|
433 { |
|
434 integration_error = true; |
|
435 (*current_liboctave_error_handler) ("unrecoverable error in odessa"); |
|
436 } |
|
437 else |
|
438 { |
|
439 switch (istate) |
|
440 { |
3996
|
441 case 1: // prior to initial integration step. |
|
442 case 2: // odessa was successful. |
|
443 t = tout; |
|
444 break; |
|
445 |
|
446 case -1: // excess work done on this call (perhaps wrong mf). |
|
447 case -2: // excess accuracy requested (tolerances too small). |
|
448 case -3: // illegal input detected (see printed message). |
|
449 case -4: // repeated error test failures (check all inputs). |
3984
|
450 case -5: // repeated convergence failures (perhaps bad jacobian |
|
451 // supplied or wrong choice of mf or tolerances). |
3996
|
452 case -6: // error weight became zero during problem. (solution |
|
453 // component i vanished, and atol or atol(i) = 0.) |
|
454 case -13: // Return requested in user-supplied function. |
3984
|
455 integration_error = 1; |
|
456 break; |
|
457 |
3996
|
458 default: |
|
459 integration_error = 1; |
3984
|
460 (*current_liboctave_error_handler) |
3996
|
461 ("unrecognized value of istate (= %d) returned from odessa", |
|
462 istate); |
3984
|
463 break; |
|
464 } |
|
465 } |
|
466 } |
|
467 |
|
468 std::string |
|
469 ODESSA::error_message (void) const |
|
470 { |
|
471 std::string retval; |
|
472 |
4051
|
473 OSSTREAM buf; |
|
474 buf << t << OSSTREAM_ENDS; |
|
475 std::string t_curr = OSSTREAM_STR (buf); |
|
476 OSSTREAM_FREEZE (buf); |
4043
|
477 |
3984
|
478 switch (istate) |
|
479 { |
|
480 case 1: |
3996
|
481 retval = "prior to initial integration step"; |
3984
|
482 break; |
|
483 |
|
484 case 2: |
|
485 retval = "successful exit"; |
|
486 break; |
|
487 |
|
488 case 3: |
|
489 retval = "prior to continuation call with modified parameters"; |
|
490 break; |
|
491 |
3996
|
492 case -1: |
4043
|
493 retval = std::string ("excess work on this call (t = ") |
|
494 + t_curr + "; perhaps wrong integration method)"; |
3996
|
495 break; |
|
496 |
|
497 case -2: |
|
498 retval = "excess accuracy requested (tolerances too small)"; |
|
499 break; |
|
500 |
|
501 case -3: |
|
502 retval = "invalid input detected (see printed message)"; |
|
503 break; |
|
504 |
|
505 case -4: |
4043
|
506 retval = std::string ("repeated error test failures (t = ") |
|
507 + t_curr + "check all inputs)"; |
3996
|
508 break; |
|
509 |
|
510 case -5: |
4043
|
511 retval = std::string ("repeated convergence failures (t = ") |
|
512 + t_curr |
|
513 + "perhaps bad jacobian supplied or wrong choice of integration method or tolerances)"; |
3996
|
514 break; |
|
515 |
|
516 case -6: |
4043
|
517 retval = std::string ("error weight became zero during problem. (t = ") |
|
518 + t_curr |
|
519 + "; solution component i vanished, and atol or atol(i) == 0)"; |
3996
|
520 break; |
|
521 |
|
522 case -13: |
4043
|
523 retval = "return requested in user-supplied function (t = " |
|
524 + t_curr + ")"; |
3996
|
525 break; |
|
526 |
3984
|
527 default: |
|
528 retval = "unknown error state"; |
|
529 break; |
|
530 } |
|
531 |
|
532 return retval; |
|
533 } |
|
534 |
|
535 |
|
536 ODESSA_result |
|
537 ODESSA::integrate (const ColumnVector& tout) |
|
538 { |
|
539 ODESSA_result retval; |
|
540 |
|
541 Matrix x_out; |
|
542 |
|
543 Array<Matrix> x_s_out; |
|
544 |
|
545 int n_out = tout.capacity (); |
|
546 |
|
547 if (n_out > 0 && n > 0) |
|
548 { |
|
549 x_out.resize (n_out, n); |
|
550 |
|
551 x_s_out.resize (npar, Matrix (n_out, n, 0.0)); |
|
552 |
|
553 for (int j = 0; j < n_out; j++) |
|
554 { |
|
555 integrate (tout(j)); |
|
556 |
|
557 if (integration_error) |
|
558 { |
|
559 retval = ODESSA_result (x_out, x_s_out); |
|
560 return retval; |
|
561 } |
|
562 |
|
563 for (int i = 0; i < n; i++) |
|
564 { |
|
565 x_out(j,i) = y(i,0); |
|
566 |
|
567 for (int k = 0; k < npar; k++) |
|
568 { |
|
569 x_s_out(k)(j,i) = y(i,k+1); |
|
570 } |
|
571 } |
|
572 } |
|
573 } |
|
574 |
|
575 retval = ODESSA_result (x_out, x_s_out); |
|
576 |
|
577 return retval; |
|
578 } |
|
579 |
|
580 ODESSA_result |
|
581 ODESSA::integrate (const ColumnVector& tout, const ColumnVector& tcrit) |
|
582 { |
|
583 ODESSA_result retval; |
|
584 |
|
585 Matrix x_out; |
|
586 |
|
587 Array<Matrix> x_s_out; |
|
588 |
|
589 int n_out = tout.capacity (); |
|
590 |
|
591 if (n_out > 0 && n > 0) |
|
592 { |
|
593 x_out.resize (n_out, n); |
|
594 |
|
595 x_s_out.resize (npar, Matrix (n_out, n, 0.0)); |
|
596 |
|
597 int n_crit = tcrit.capacity (); |
|
598 |
|
599 if (n_crit > 0) |
|
600 { |
|
601 int i_crit = 0; |
|
602 int i_out = 0; |
|
603 double next_crit = tcrit(0); |
|
604 double next_out; |
|
605 while (i_out < n_out) |
|
606 { |
|
607 bool do_restart = false; |
|
608 |
|
609 next_out = tout(i_out); |
|
610 if (i_crit < n_crit) |
|
611 next_crit = tcrit(i_crit); |
|
612 |
|
613 int save_output; |
|
614 double t_out; |
|
615 |
|
616 if (next_crit == next_out) |
|
617 { |
|
618 set_stop_time (next_crit); |
|
619 t_out = next_out; |
|
620 save_output = true; |
|
621 i_out++; |
|
622 i_crit++; |
|
623 do_restart = true; |
|
624 } |
|
625 else if (next_crit < next_out) |
|
626 { |
|
627 if (i_crit < n_crit) |
|
628 { |
|
629 set_stop_time (next_crit); |
|
630 t_out = next_crit; |
|
631 save_output = false; |
|
632 i_crit++; |
|
633 do_restart = true; |
|
634 } |
|
635 else |
|
636 { |
|
637 clear_stop_time (); |
|
638 t_out = next_out; |
|
639 save_output = true; |
|
640 i_out++; |
|
641 } |
|
642 } |
|
643 else |
|
644 { |
|
645 set_stop_time (next_crit); |
|
646 t_out = next_out; |
|
647 save_output = true; |
|
648 i_out++; |
|
649 } |
|
650 integrate (t_out); |
|
651 if (integration_error) |
|
652 { |
|
653 retval = ODESSA_result (x_out, x_s_out); |
|
654 return retval; |
|
655 } |
|
656 if (save_output) |
|
657 { |
|
658 for (int i = 0; i < n; i++) |
|
659 { |
|
660 x_out(i_out-1,i) = y(i,0); |
|
661 |
|
662 for (int k = 0; k < npar; k++) |
|
663 { |
|
664 x_s_out(k)(i_out-1,i) = y(i,k+1); |
|
665 } |
|
666 } |
|
667 } |
|
668 |
|
669 if (do_restart) |
|
670 force_restart (); |
|
671 } |
|
672 |
|
673 retval = ODESSA_result (x_out, x_s_out); |
|
674 } |
|
675 else |
|
676 { |
|
677 retval = integrate (tout); |
|
678 |
|
679 if (integration_error) |
|
680 return retval; |
|
681 } |
|
682 } |
|
683 |
|
684 return retval; |
|
685 } |
|
686 |
|
687 /* |
|
688 ;;; Local Variables: *** |
|
689 ;;; mode: C++ *** |
|
690 ;;; End: *** |
|
691 */ |