3
|
1 /* |
|
2 |
2847
|
3 Copyright (C) 1996, 1997 John W. Eaton |
3
|
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 |
7016
|
9 Free Software Foundation; either version 3 of the License, or (at your |
|
10 option) any later version. |
3
|
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 |
7016
|
18 along with Octave; see the file COPYING. If not, see |
|
19 <http://www.gnu.org/licenses/>. |
3
|
20 |
|
21 */ |
|
22 |
238
|
23 #ifdef HAVE_CONFIG_H |
1192
|
24 #include <config.h> |
3
|
25 #endif |
|
26 |
1367
|
27 #include <cfloat> |
|
28 #include <cmath> |
|
29 |
5765
|
30 #include <sstream> |
|
31 |
1842
|
32 #include "LSODE.h" |
1847
|
33 #include "f77-fcn.h" |
227
|
34 #include "lo-error.h" |
4180
|
35 #include "quit.h" |
3
|
36 |
5275
|
37 typedef octave_idx_type (*lsode_fcn_ptr) (const octave_idx_type&, const double&, double*, |
|
38 double*, octave_idx_type&); |
3507
|
39 |
5275
|
40 typedef octave_idx_type (*lsode_jac_ptr) (const octave_idx_type&, const double&, double*, |
|
41 const octave_idx_type&, const octave_idx_type&, double*, const |
|
42 octave_idx_type&); |
3507
|
43 |
3
|
44 extern "C" |
4552
|
45 { |
|
46 F77_RET_T |
5275
|
47 F77_FUNC (dlsode, DLSODE) (lsode_fcn_ptr, octave_idx_type&, double*, double&, |
|
48 double&, octave_idx_type&, double&, const double*, octave_idx_type&, |
|
49 octave_idx_type&, octave_idx_type&, double*, octave_idx_type&, octave_idx_type*, octave_idx_type&, |
|
50 lsode_jac_ptr, octave_idx_type&); |
4552
|
51 } |
3
|
52 |
532
|
53 static ODEFunc::ODERHSFunc user_fun; |
|
54 static ODEFunc::ODEJacFunc user_jac; |
3
|
55 static ColumnVector *tmp_x; |
|
56 |
5275
|
57 static octave_idx_type |
|
58 lsode_f (const octave_idx_type& neq, const double& time, double *, |
|
59 double *deriv, octave_idx_type& ierr) |
3
|
60 { |
4180
|
61 BEGIN_INTERRUPT_WITH_EXCEPTIONS; |
|
62 |
2343
|
63 ColumnVector tmp_deriv; |
3
|
64 |
1360
|
65 // NOTE: this won't work if LSODE passes copies of the state vector. |
|
66 // In that case we have to create a temporary vector object |
|
67 // and copy. |
|
68 |
1251
|
69 tmp_deriv = (*user_fun) (*tmp_x, time); |
3
|
70 |
258
|
71 if (tmp_deriv.length () == 0) |
1251
|
72 ierr = -1; |
258
|
73 else |
|
74 { |
5275
|
75 for (octave_idx_type i = 0; i < neq; i++) |
258
|
76 deriv [i] = tmp_deriv.elem (i); |
|
77 } |
3
|
78 |
4180
|
79 END_INTERRUPT_WITH_EXCEPTIONS; |
|
80 |
3
|
81 return 0; |
|
82 } |
|
83 |
5275
|
84 static octave_idx_type |
|
85 lsode_j (const octave_idx_type& neq, const double& time, double *, |
|
86 const octave_idx_type&, const octave_idx_type&, double *pd, const octave_idx_type& nrowpd) |
3
|
87 { |
4180
|
88 BEGIN_INTERRUPT_WITH_EXCEPTIONS; |
|
89 |
1251
|
90 Matrix tmp_jac (neq, neq); |
3
|
91 |
1360
|
92 // NOTE: this won't work if LSODE passes copies of the state vector. |
|
93 // In that case we have to create a temporary vector object |
|
94 // and copy. |
|
95 |
1251
|
96 tmp_jac = (*user_jac) (*tmp_x, time); |
3
|
97 |
5275
|
98 for (octave_idx_type j = 0; j < neq; j++) |
|
99 for (octave_idx_type i = 0; i < neq; i++) |
1251
|
100 pd [nrowpd * j + i] = tmp_jac (i, j); |
3
|
101 |
4180
|
102 END_INTERRUPT_WITH_EXCEPTIONS; |
|
103 |
3
|
104 return 0; |
|
105 } |
|
106 |
|
107 ColumnVector |
1842
|
108 LSODE::do_integrate (double tout) |
3
|
109 { |
1945
|
110 ColumnVector retval; |
|
111 |
5275
|
112 static octave_idx_type nn = 0; |
4049
|
113 |
|
114 if (! initialized || restart || ODEFunc::reset || LSODE_options::reset) |
1945
|
115 { |
4049
|
116 integration_error = false; |
1945
|
117 |
4049
|
118 initialized = true; |
|
119 |
|
120 istate = 1; |
|
121 |
5275
|
122 octave_idx_type n = size (); |
4049
|
123 |
|
124 nn = n; |
3955
|
125 |
5275
|
126 octave_idx_type max_maxord = 0; |
4231
|
127 |
4049
|
128 if (integration_method () == "stiff") |
|
129 { |
4231
|
130 max_maxord = 5; |
|
131 |
4049
|
132 if (jac) |
|
133 method_flag = 21; |
|
134 else |
|
135 method_flag = 22; |
3955
|
136 |
5553
|
137 liw = 20 + n; |
|
138 lrw = 22 + n * (9 + n); |
4049
|
139 } |
|
140 else |
|
141 { |
4231
|
142 max_maxord = 12; |
|
143 |
4049
|
144 method_flag = 10; |
3955
|
145 |
4049
|
146 liw = 20; |
|
147 lrw = 22 + 16 * n; |
|
148 } |
|
149 |
4231
|
150 maxord = maximum_order (); |
|
151 |
5552
|
152 iwork.resize (liw); |
|
153 |
|
154 for (octave_idx_type i = 4; i < 9; i++) |
|
155 iwork(i) = 0; |
|
156 |
|
157 rwork.resize (lrw); |
|
158 |
|
159 for (octave_idx_type i = 4; i < 9; i++) |
|
160 rwork(i) = 0; |
|
161 |
4231
|
162 if (maxord >= 0) |
|
163 { |
|
164 if (maxord > 0 && maxord <= max_maxord) |
|
165 { |
|
166 iwork(4) = maxord; |
|
167 iopt = 1; |
|
168 } |
|
169 else |
|
170 { |
|
171 (*current_liboctave_error_handler) |
|
172 ("lsode: invalid value for maximum order"); |
|
173 integration_error = true; |
|
174 return retval; |
|
175 } |
|
176 } |
|
177 |
4049
|
178 if (stop_time_set) |
|
179 { |
|
180 itask = 4; |
|
181 rwork(0) = stop_time; |
|
182 iopt = 1; |
|
183 } |
|
184 else |
|
185 { |
|
186 itask = 1; |
|
187 } |
258
|
188 |
4049
|
189 px = x.fortran_vec (); |
3
|
190 |
4049
|
191 piwork = iwork.fortran_vec (); |
|
192 prwork = rwork.fortran_vec (); |
|
193 |
|
194 restart = false; |
|
195 |
|
196 // ODEFunc |
3
|
197 |
4049
|
198 // NOTE: this won't work if LSODE passes copies of the state vector. |
|
199 // In that case we have to create a temporary vector object |
|
200 // and copy. |
3
|
201 |
4049
|
202 tmp_x = &x; |
|
203 |
|
204 user_fun = function (); |
|
205 user_jac = jacobian_function (); |
|
206 |
2343
|
207 ColumnVector xdot = (*user_fun) (x, t); |
|
208 |
|
209 if (x.length () != xdot.length ()) |
|
210 { |
|
211 (*current_liboctave_error_handler) |
|
212 ("lsode: inconsistent sizes for state and derivative vectors"); |
|
213 |
3995
|
214 integration_error = true; |
2343
|
215 return retval; |
|
216 } |
|
217 |
4049
|
218 ODEFunc::reset = false; |
|
219 |
|
220 // LSODE_options |
|
221 |
|
222 rel_tol = relative_tolerance (); |
|
223 abs_tol = absolute_tolerance (); |
|
224 |
5275
|
225 octave_idx_type abs_tol_len = abs_tol.length (); |
4049
|
226 |
|
227 if (abs_tol_len == 1) |
|
228 itol = 1; |
|
229 else if (abs_tol_len == n) |
|
230 itol = 2; |
|
231 else |
|
232 { |
|
233 (*current_liboctave_error_handler) |
|
234 ("lsode: inconsistent sizes for state and absolute tolerance vectors"); |
|
235 |
|
236 integration_error = true; |
|
237 return retval; |
|
238 } |
2343
|
239 |
4049
|
240 double iss = initial_step_size (); |
|
241 if (iss >= 0.0) |
|
242 { |
|
243 rwork(4) = iss; |
|
244 iopt = 1; |
|
245 } |
|
246 |
|
247 double maxss = maximum_step_size (); |
|
248 if (maxss >= 0.0) |
|
249 { |
|
250 rwork(5) = maxss; |
|
251 iopt = 1; |
|
252 } |
|
253 |
|
254 double minss = minimum_step_size (); |
|
255 if (minss >= 0.0) |
|
256 { |
|
257 rwork(6) = minss; |
|
258 iopt = 1; |
|
259 } |
|
260 |
5275
|
261 octave_idx_type sl = step_limit (); |
4049
|
262 if (sl > 0) |
|
263 { |
|
264 iwork(5) = sl; |
|
265 iopt = 1; |
|
266 } |
|
267 |
|
268 pabs_tol = abs_tol.fortran_vec (); |
|
269 |
|
270 LSODE_options::reset = false; |
3
|
271 } |
|
272 |
4583
|
273 F77_XFCN (dlsode, DLSODE, (lsode_f, nn, px, t, tout, itol, rel_tol, |
|
274 pabs_tol, itask, istate, iopt, prwork, lrw, |
|
275 piwork, liw, lsode_j, method_flag)); |
3
|
276 |
1945
|
277 if (f77_exception_encountered) |
3178
|
278 { |
3995
|
279 integration_error = true; |
3178
|
280 (*current_liboctave_error_handler) ("unrecoverable error in lsode"); |
|
281 } |
1945
|
282 else |
3
|
283 { |
1945
|
284 switch (istate) |
|
285 { |
3996
|
286 case 1: // prior to initial integration step. |
|
287 case 2: // lsode was successful. |
1945
|
288 retval = x; |
|
289 t = tout; |
|
290 break; |
|
291 |
3996
|
292 case -1: // excess work done on this call (perhaps wrong mf). |
|
293 case -2: // excess accuracy requested (tolerances too small). |
|
294 case -3: // illegal input detected (see printed message). |
|
295 case -4: // repeated error test failures (check all inputs). |
|
296 case -5: // repeated convergence failures (perhaps bad jacobian |
|
297 // supplied or wrong choice of mf or tolerances). |
|
298 case -6: // error weight became zero during problem. (solution |
|
299 // component i vanished, and atol or atol(i) = 0.) |
|
300 case -13: // return requested in user-supplied function. |
|
301 integration_error = true; |
|
302 break; |
|
303 |
|
304 default: |
|
305 integration_error = true; |
227
|
306 (*current_liboctave_error_handler) |
3996
|
307 ("unrecognized value of istate (= %d) returned from lsode", |
|
308 istate); |
1945
|
309 break; |
3
|
310 } |
|
311 } |
|
312 |
1945
|
313 return retval; |
3
|
314 } |
|
315 |
3959
|
316 std::string |
|
317 LSODE::error_message (void) const |
|
318 { |
|
319 std::string retval; |
|
320 |
5765
|
321 std::ostringstream buf; |
|
322 buf << t; |
|
323 std::string t_curr = buf.str (); |
4042
|
324 |
3959
|
325 switch (istate) |
|
326 { |
|
327 case 1: |
3996
|
328 retval = "prior to initial integration step"; |
3959
|
329 break; |
|
330 |
|
331 case 2: |
|
332 retval = "successful exit"; |
|
333 break; |
|
334 |
|
335 case 3: |
|
336 retval = "prior to continuation call with modified parameters"; |
|
337 break; |
|
338 |
3996
|
339 case -1: |
4042
|
340 retval = std::string ("excess work on this call (t = ") |
4043
|
341 + t_curr + "; perhaps wrong integration method)"; |
3996
|
342 break; |
|
343 |
|
344 case -2: |
|
345 retval = "excess accuracy requested (tolerances too small)"; |
|
346 break; |
|
347 |
|
348 case -3: |
|
349 retval = "invalid input detected (see printed message)"; |
|
350 break; |
|
351 |
|
352 case -4: |
4042
|
353 retval = std::string ("repeated error test failures (t = ") |
4043
|
354 + t_curr + "check all inputs)"; |
3996
|
355 break; |
|
356 |
|
357 case -5: |
4042
|
358 retval = std::string ("repeated convergence failures (t = ") |
4043
|
359 + t_curr |
4042
|
360 + "perhaps bad jacobian supplied or wrong choice of integration method or tolerances)"; |
3996
|
361 break; |
|
362 |
|
363 case -6: |
4042
|
364 retval = std::string ("error weight became zero during problem. (t = ") |
4043
|
365 + t_curr |
4042
|
366 + "; solution component i vanished, and atol or atol(i) == 0)"; |
3996
|
367 break; |
|
368 |
|
369 case -13: |
4042
|
370 retval = "return requested in user-supplied function (t = " |
4043
|
371 + t_curr + ")"; |
3996
|
372 break; |
|
373 |
3959
|
374 default: |
|
375 retval = "unknown error state"; |
|
376 break; |
|
377 } |
|
378 |
|
379 return retval; |
|
380 } |
|
381 |
3
|
382 Matrix |
1842
|
383 LSODE::do_integrate (const ColumnVector& tout) |
3
|
384 { |
|
385 Matrix retval; |
4049
|
386 |
5275
|
387 octave_idx_type n_out = tout.capacity (); |
|
388 octave_idx_type n = size (); |
3
|
389 |
|
390 if (n_out > 0 && n > 0) |
|
391 { |
|
392 retval.resize (n_out, n); |
|
393 |
5275
|
394 for (octave_idx_type i = 0; i < n; i++) |
3
|
395 retval.elem (0, i) = x.elem (i); |
|
396 |
5275
|
397 for (octave_idx_type j = 1; j < n_out; j++) |
3
|
398 { |
1842
|
399 ColumnVector x_next = do_integrate (tout.elem (j)); |
258
|
400 |
|
401 if (integration_error) |
|
402 return retval; |
|
403 |
5275
|
404 for (octave_idx_type i = 0; i < n; i++) |
3
|
405 retval.elem (j, i) = x_next.elem (i); |
|
406 } |
|
407 } |
|
408 |
|
409 return retval; |
|
410 } |
|
411 |
|
412 Matrix |
3511
|
413 LSODE::do_integrate (const ColumnVector& tout, const ColumnVector& tcrit) |
3
|
414 { |
|
415 Matrix retval; |
4049
|
416 |
5275
|
417 octave_idx_type n_out = tout.capacity (); |
|
418 octave_idx_type n = size (); |
3
|
419 |
|
420 if (n_out > 0 && n > 0) |
|
421 { |
|
422 retval.resize (n_out, n); |
|
423 |
5275
|
424 for (octave_idx_type i = 0; i < n; i++) |
3
|
425 retval.elem (0, i) = x.elem (i); |
|
426 |
5275
|
427 octave_idx_type n_crit = tcrit.capacity (); |
3
|
428 |
|
429 if (n_crit > 0) |
|
430 { |
5275
|
431 octave_idx_type i_crit = 0; |
|
432 octave_idx_type i_out = 1; |
3
|
433 double next_crit = tcrit.elem (0); |
|
434 double next_out; |
|
435 while (i_out < n_out) |
|
436 { |
3995
|
437 bool do_restart = false; |
3
|
438 |
|
439 next_out = tout.elem (i_out); |
|
440 if (i_crit < n_crit) |
|
441 next_crit = tcrit.elem (i_crit); |
|
442 |
5275
|
443 octave_idx_type save_output; |
3
|
444 double t_out; |
|
445 |
|
446 if (next_crit == next_out) |
|
447 { |
|
448 set_stop_time (next_crit); |
|
449 t_out = next_out; |
|
450 save_output = 1; |
|
451 i_out++; |
|
452 i_crit++; |
3995
|
453 do_restart = true; |
3
|
454 } |
|
455 else if (next_crit < next_out) |
|
456 { |
|
457 if (i_crit < n_crit) |
|
458 { |
|
459 set_stop_time (next_crit); |
|
460 t_out = next_crit; |
|
461 save_output = 0; |
|
462 i_crit++; |
3995
|
463 do_restart = true; |
3
|
464 } |
|
465 else |
|
466 { |
|
467 clear_stop_time (); |
|
468 t_out = next_out; |
|
469 save_output = 1; |
|
470 i_out++; |
|
471 } |
|
472 } |
|
473 else |
|
474 { |
|
475 set_stop_time (next_crit); |
|
476 t_out = next_out; |
|
477 save_output = 1; |
|
478 i_out++; |
|
479 } |
|
480 |
1842
|
481 ColumnVector x_next = do_integrate (t_out); |
3
|
482 |
258
|
483 if (integration_error) |
|
484 return retval; |
|
485 |
3
|
486 if (save_output) |
|
487 { |
5275
|
488 for (octave_idx_type i = 0; i < n; i++) |
3
|
489 retval.elem (i_out-1, i) = x_next.elem (i); |
|
490 } |
|
491 |
|
492 if (do_restart) |
|
493 force_restart (); |
|
494 } |
|
495 } |
|
496 else |
258
|
497 { |
1842
|
498 retval = do_integrate (tout); |
258
|
499 |
|
500 if (integration_error) |
|
501 return retval; |
|
502 } |
3
|
503 } |
|
504 |
|
505 return retval; |
|
506 } |
|
507 |
|
508 /* |
|
509 ;;; Local Variables: *** |
|
510 ;;; mode: C++ *** |
|
511 ;;; End: *** |
|
512 */ |