2928
|
1 /* |
|
2 |
|
3 Copyright (C) 1996, 1997 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 |
7016
|
9 Free Software Foundation; either version 3 of the License, or (at your |
|
10 option) any later version. |
2928
|
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/>. |
2928
|
20 |
|
21 */ |
|
22 |
|
23 #ifdef HAVE_CONFIG_H |
|
24 #include <config.h> |
|
25 #endif |
|
26 |
|
27 #include <string> |
|
28 |
3567
|
29 #include <iomanip> |
3523
|
30 #include <iostream> |
2928
|
31 |
|
32 #include "DASSL.h" |
|
33 |
|
34 #include "defun-dld.h" |
|
35 #include "error.h" |
|
36 #include "gripes.h" |
|
37 #include "oct-obj.h" |
2968
|
38 #include "ov-fcn.h" |
5729
|
39 #include "ov-cell.h" |
2928
|
40 #include "pager.h" |
3243
|
41 #include "unwind-prot.h" |
2928
|
42 #include "utils.h" |
|
43 #include "variables.h" |
|
44 |
3998
|
45 #include "DASSL-opts.cc" |
|
46 |
2928
|
47 // Global pointer for user defined function required by dassl. |
2968
|
48 static octave_function *dassl_fcn; |
2928
|
49 |
3991
|
50 // Global pointer for optional user defined jacobian function. |
|
51 static octave_function *dassl_jac; |
|
52 |
4140
|
53 // Have we warned about imaginary values returned from user function? |
|
54 static bool warned_fcn_imaginary = false; |
|
55 static bool warned_jac_imaginary = false; |
|
56 |
3243
|
57 // Is this a recursive call? |
|
58 static int call_depth = 0; |
|
59 |
2928
|
60 ColumnVector |
3849
|
61 dassl_user_function (const ColumnVector& x, const ColumnVector& xdot, |
5275
|
62 double t, octave_idx_type& ires) |
2928
|
63 { |
|
64 ColumnVector retval; |
|
65 |
4628
|
66 assert (x.capacity () == xdot.capacity ()); |
2928
|
67 |
|
68 octave_value_list args; |
4628
|
69 |
2928
|
70 args(2) = t; |
4628
|
71 args(1) = xdot; |
|
72 args(0) = x; |
2928
|
73 |
|
74 if (dassl_fcn) |
|
75 { |
3544
|
76 octave_value_list tmp = dassl_fcn->do_multi_index_op (1, args); |
2928
|
77 |
|
78 if (error_state) |
|
79 { |
|
80 gripe_user_supplied_eval ("dassl"); |
|
81 return retval; |
|
82 } |
|
83 |
3849
|
84 int tlen = tmp.length (); |
|
85 if (tlen > 0 && tmp(0).is_defined ()) |
2928
|
86 { |
4140
|
87 if (! warned_fcn_imaginary && tmp(0).is_complex_type ()) |
|
88 { |
|
89 warning ("dassl: ignoring imaginary part returned from user-supplied function"); |
|
90 warned_fcn_imaginary = true; |
|
91 } |
|
92 |
3419
|
93 retval = ColumnVector (tmp(0).vector_value ()); |
2928
|
94 |
3849
|
95 if (tlen > 1) |
|
96 ires = tmp(1).int_value (); |
|
97 |
2928
|
98 if (error_state || retval.length () == 0) |
|
99 gripe_user_supplied_eval ("dassl"); |
|
100 } |
|
101 else |
|
102 gripe_user_supplied_eval ("dassl"); |
|
103 } |
|
104 |
|
105 return retval; |
|
106 } |
|
107 |
3991
|
108 Matrix |
|
109 dassl_user_jacobian (const ColumnVector& x, const ColumnVector& xdot, |
|
110 double t, double cj) |
|
111 { |
|
112 Matrix retval; |
|
113 |
4628
|
114 assert (x.capacity () == xdot.capacity ()); |
3991
|
115 |
|
116 octave_value_list args; |
|
117 |
|
118 args(3) = cj; |
|
119 args(2) = t; |
4628
|
120 args(1) = xdot; |
|
121 args(0) = x; |
3991
|
122 |
|
123 if (dassl_jac) |
|
124 { |
|
125 octave_value_list tmp = dassl_jac->do_multi_index_op (1, args); |
|
126 |
|
127 if (error_state) |
|
128 { |
|
129 gripe_user_supplied_eval ("dassl"); |
|
130 return retval; |
|
131 } |
|
132 |
|
133 int tlen = tmp.length (); |
|
134 if (tlen > 0 && tmp(0).is_defined ()) |
|
135 { |
4140
|
136 if (! warned_jac_imaginary && tmp(0).is_complex_type ()) |
|
137 { |
|
138 warning ("dassl: ignoring imaginary part returned from user-supplied jacobian function"); |
|
139 warned_jac_imaginary = true; |
|
140 } |
|
141 |
3991
|
142 retval = tmp(0).matrix_value (); |
|
143 |
|
144 if (error_state || retval.length () == 0) |
|
145 gripe_user_supplied_eval ("dassl"); |
|
146 } |
|
147 else |
|
148 gripe_user_supplied_eval ("dassl"); |
|
149 } |
|
150 |
|
151 return retval; |
|
152 } |
|
153 |
3323
|
154 #define DASSL_ABORT() \ |
|
155 do \ |
|
156 { \ |
|
157 unwind_protect::run_frame ("Fdassl"); \ |
|
158 return retval; \ |
|
159 } \ |
|
160 while (0) |
|
161 |
|
162 #define DASSL_ABORT1(msg) \ |
|
163 do \ |
|
164 { \ |
3747
|
165 ::error ("dassl: " msg); \ |
3323
|
166 DASSL_ABORT (); \ |
|
167 } \ |
|
168 while (0) |
|
169 |
|
170 #define DASSL_ABORT2(fmt, arg) \ |
|
171 do \ |
|
172 { \ |
3747
|
173 ::error ("dassl: " fmt, arg); \ |
3323
|
174 DASSL_ABORT (); \ |
|
175 } \ |
|
176 while (0) |
|
177 |
3991
|
178 DEFUN_DLD (dassl, args, nargout, |
3373
|
179 "-*- texinfo -*-\n\ |
4115
|
180 @deftypefn {Loadable Function} {[@var{x}, @var{xdot}, @var{istate}, @var{msg}] =} dassl (@var{fcn}, @var{x_0}, @var{xdot_0}, @var{t}, @var{t_crit})\n\ |
|
181 Solve the set of differential-algebraic equations\n\ |
4912
|
182 @iftex\n\ |
4115
|
183 @tex\n\ |
4852
|
184 $$ 0 = f (x, \\dot{x}, t) $$\n\ |
4115
|
185 with\n\ |
|
186 $$ x(t_0) = x_0, \\dot{x}(t_0) = \\dot{x}_0 $$\n\ |
|
187 @end tex\n\ |
4912
|
188 @end iftex\n\ |
4115
|
189 @ifinfo\n\ |
|
190 \n\ |
|
191 @example\n\ |
4698
|
192 0 = f (x, xdot, t)\n\ |
4115
|
193 @end example\n\ |
|
194 \n\ |
4912
|
195 @noindent\n\ |
4115
|
196 with\n\ |
|
197 \n\ |
|
198 @example\n\ |
|
199 x(t_0) = x_0, xdot(t_0) = xdot_0\n\ |
|
200 @end example\n\ |
|
201 \n\ |
|
202 @end ifinfo\n\ |
|
203 The solution is returned in the matrices @var{x} and @var{xdot},\n\ |
|
204 with each row in the result matrices corresponding to one of the\n\ |
3373
|
205 elements in the vector @var{t}. The first element of @var{t}\n\ |
4115
|
206 should be @math{t_0} and correspond to the initial state of the\n\ |
|
207 system @var{x_0} and its derivative @var{xdot_0}, so that the first\n\ |
|
208 row of the output @var{x} is @var{x_0} and the first row\n\ |
|
209 of the output @var{xdot} is @var{xdot_0}.\n\ |
3373
|
210 \n\ |
5729
|
211 The first argument, @var{fcn}, is a string or a two element cell array\n\ |
|
212 of strings, inline or function handle, that names the function, to call\n\ |
|
213 to compute the vector of residuals for the set of equations. It must\n\ |
|
214 have the form\n\ |
3373
|
215 \n\ |
|
216 @example\n\ |
|
217 @var{res} = f (@var{x}, @var{xdot}, @var{t})\n\ |
|
218 @end example\n\ |
2928
|
219 \n\ |
3373
|
220 @noindent\n\ |
4115
|
221 in which @var{x}, @var{xdot}, and @var{res} are vectors, and @var{t} is a\n\ |
3373
|
222 scalar.\n\ |
|
223 \n\ |
4115
|
224 If @var{fcn} is a two-element string array, the first element names\n\ |
|
225 the function @math{f} described above, and the second element names\n\ |
4117
|
226 a function to compute the modified Jacobian\n\ |
4115
|
227 \n\ |
4912
|
228 @iftex\n\ |
4115
|
229 @tex\n\ |
|
230 $$\n\ |
|
231 J = {\\partial f \\over \\partial x}\n\ |
|
232 + c {\\partial f \\over \\partial \\dot{x}}\n\ |
|
233 $$\n\ |
|
234 @end tex\n\ |
4912
|
235 @end iftex\n\ |
|
236 @ifinfo\n\ |
4888
|
237 @example\n\ |
4115
|
238 df df\n\ |
|
239 jac = -- + c ------\n\ |
|
240 dx d xdot\n\ |
|
241 @end example\n\ |
|
242 @end ifinfo\n\ |
|
243 \n\ |
|
244 The modified Jacobian function must have the form\n\ |
|
245 \n\ |
|
246 @example\n\ |
|
247 \n\ |
|
248 @var{jac} = j (@var{x}, @var{xdot}, @var{t}, @var{c})\n\ |
|
249 \n\ |
|
250 @end example\n\ |
|
251 \n\ |
3373
|
252 The second and third arguments to @code{dassl} specify the initial\n\ |
|
253 condition of the states and their derivatives, and the fourth argument\n\ |
4115
|
254 specifies a vector of output times at which the solution is desired,\n\ |
3373
|
255 including the time corresponding to the initial condition.\n\ |
2928
|
256 \n\ |
3373
|
257 The set of initial states and derivatives are not strictly required to\n\ |
|
258 be consistent. In practice, however, @sc{Dassl} is not very good at\n\ |
|
259 determining a consistent set for you, so it is best if you ensure that\n\ |
|
260 the initial values result in the function evaluating to zero.\n\ |
2928
|
261 \n\ |
3373
|
262 The fifth argument is optional, and may be used to specify a set of\n\ |
|
263 times that the DAE solver should not integrate past. It is useful for\n\ |
|
264 avoiding difficulties with singularities and points where there is a\n\ |
|
265 discontinuity in the derivative.\n\ |
3964
|
266 \n\ |
4115
|
267 After a successful computation, the value of @var{istate} will be\n\ |
|
268 greater than zero (consistent with the Fortran version of @sc{Dassl}).\n\ |
|
269 \n\ |
|
270 If the computation is not successful, the value of @var{istate} will be\n\ |
|
271 less than zero and @var{msg} will contain additional information.\n\ |
|
272 \n\ |
3964
|
273 You can use the function @code{dassl_options} to set optional\n\ |
|
274 parameters for @code{dassl}.\n\ |
5694
|
275 @seealso{daspk, dasrt, lsode}\n\ |
5642
|
276 @end deftypefn") |
2928
|
277 { |
|
278 octave_value_list retval; |
|
279 |
4140
|
280 warned_fcn_imaginary = false; |
|
281 warned_jac_imaginary = false; |
|
282 |
3243
|
283 unwind_protect::begin_frame ("Fdassl"); |
2928
|
284 |
3243
|
285 unwind_protect_int (call_depth); |
|
286 call_depth++; |
2928
|
287 |
3243
|
288 if (call_depth > 1) |
3323
|
289 DASSL_ABORT1 ("invalid recursive call"); |
2928
|
290 |
3243
|
291 int nargin = args.length (); |
|
292 |
3991
|
293 if (nargin > 3 && nargin < 6 && nargout < 5) |
2928
|
294 { |
5729
|
295 std::string fcn_name, fname, jac_name, jname; |
3991
|
296 dassl_fcn = 0; |
|
297 dassl_jac = 0; |
|
298 |
|
299 octave_value f_arg = args(0); |
|
300 |
5729
|
301 if (f_arg.is_cell ()) |
|
302 { |
|
303 Cell c = f_arg.cell_value (); |
|
304 if (c.length() == 1) |
|
305 f_arg = c(0); |
|
306 else if (c.length() == 2) |
|
307 { |
|
308 if (c(0).is_function_handle () || c(0).is_inline_function ()) |
|
309 dassl_fcn = c(0).function_value (); |
|
310 else |
|
311 { |
|
312 fcn_name = unique_symbol_name ("__dassl_fcn__"); |
|
313 fname = "function y = "; |
|
314 fname.append (fcn_name); |
|
315 fname.append (" (x, xdot, t) y = "); |
|
316 dassl_fcn = extract_function |
|
317 (c(0), "dassl", fcn_name, fname, "; endfunction"); |
|
318 } |
|
319 |
|
320 if (dassl_fcn) |
|
321 { |
|
322 if (c(1).is_function_handle () || c(1).is_inline_function ()) |
|
323 dassl_jac = c(1).function_value (); |
|
324 else |
|
325 { |
|
326 jac_name = unique_symbol_name ("__dassl_jac__"); |
|
327 jname = "function jac = "; |
|
328 jname.append(jac_name); |
|
329 jname.append (" (x, xdot, t, cj) jac = "); |
|
330 dassl_jac = extract_function |
|
331 (c(1), "dassl", jac_name, jname, "; endfunction"); |
3991
|
332 |
5729
|
333 if (!dassl_jac) |
|
334 { |
|
335 if (fcn_name.length()) |
|
336 clear_function (fcn_name); |
|
337 dassl_fcn = 0; |
|
338 } |
|
339 } |
|
340 } |
|
341 } |
|
342 else |
|
343 DASSL_ABORT1 ("incorrect number of elements in cell array"); |
|
344 } |
3243
|
345 |
5729
|
346 if (!dassl_fcn && ! f_arg.is_cell()) |
|
347 { |
|
348 if (f_arg.is_function_handle () || f_arg.is_inline_function ()) |
|
349 dassl_fcn = f_arg.function_value (); |
|
350 else |
|
351 { |
|
352 switch (f_arg.rows ()) |
|
353 { |
|
354 case 1: |
|
355 do |
|
356 { |
|
357 fcn_name = unique_symbol_name ("__dassl_fcn__"); |
|
358 fname = "function y = "; |
|
359 fname.append (fcn_name); |
|
360 fname.append (" (x, xdot, t) y = "); |
|
361 dassl_fcn = extract_function |
|
362 (f_arg, "dassl", fcn_name, fname, "; endfunction"); |
|
363 } |
|
364 while (0); |
|
365 break; |
3991
|
366 |
5729
|
367 case 2: |
3991
|
368 { |
5729
|
369 string_vector tmp = f_arg.all_strings (); |
|
370 |
|
371 if (! error_state) |
|
372 { |
|
373 fcn_name = unique_symbol_name ("__dassl_fcn__"); |
|
374 fname = "function y = "; |
|
375 fname.append (fcn_name); |
|
376 fname.append (" (x, xdot, t) y = "); |
|
377 dassl_fcn = extract_function |
|
378 (tmp(0), "dassl", fcn_name, fname, "; endfunction"); |
3991
|
379 |
5729
|
380 if (dassl_fcn) |
|
381 { |
|
382 jac_name = unique_symbol_name ("__dassl_jac__"); |
|
383 jname = "function jac = "; |
|
384 jname.append(jac_name); |
|
385 jname.append (" (x, xdot, t, cj) jac = "); |
|
386 dassl_jac = extract_function |
|
387 (tmp(1), "dassl", jac_name, jname, |
|
388 "; endfunction"); |
|
389 |
|
390 if (!dassl_jac) |
|
391 { |
|
392 if (fcn_name.length()) |
|
393 clear_function (fcn_name); |
|
394 dassl_fcn = 0; |
|
395 } |
|
396 } |
|
397 } |
3991
|
398 } |
5729
|
399 } |
|
400 } |
3991
|
401 } |
|
402 |
|
403 if (error_state || ! dassl_fcn) |
3323
|
404 DASSL_ABORT (); |
3243
|
405 |
3419
|
406 ColumnVector state = ColumnVector (args(1).vector_value ()); |
2928
|
407 |
|
408 if (error_state) |
3323
|
409 DASSL_ABORT1 ("expecting state vector as second argument"); |
3243
|
410 |
3419
|
411 ColumnVector deriv (args(2).vector_value ()); |
3243
|
412 |
|
413 if (error_state) |
3323
|
414 DASSL_ABORT1 ("expecting derivative vector as third argument"); |
3243
|
415 |
3419
|
416 ColumnVector out_times (args(3).vector_value ()); |
3243
|
417 |
|
418 if (error_state) |
3323
|
419 DASSL_ABORT1 ("expecting output time vector as fourth argument"); |
2928
|
420 |
3243
|
421 ColumnVector crit_times; |
|
422 int crit_times_set = 0; |
|
423 if (nargin > 4) |
|
424 { |
3419
|
425 crit_times = ColumnVector (args(4).vector_value ()); |
2928
|
426 |
3243
|
427 if (error_state) |
3323
|
428 DASSL_ABORT1 ("expecting critical time vector as fifth argument"); |
2928
|
429 |
3243
|
430 crit_times_set = 1; |
|
431 } |
2928
|
432 |
3243
|
433 if (state.capacity () != deriv.capacity ()) |
3323
|
434 DASSL_ABORT1 ("x and xdot must have the same size"); |
3243
|
435 |
|
436 double tzero = out_times (0); |
2928
|
437 |
3243
|
438 DAEFunc func (dassl_user_function); |
3991
|
439 if (dassl_jac) |
|
440 func.set_jacobian_function (dassl_user_jacobian); |
|
441 |
3243
|
442 DASSL dae (state, deriv, tzero, func); |
3991
|
443 |
4122
|
444 dae.set_options (dassl_opts); |
2928
|
445 |
3243
|
446 Matrix output; |
|
447 Matrix deriv_output; |
|
448 |
|
449 if (crit_times_set) |
|
450 output = dae.integrate (out_times, deriv_output, crit_times); |
|
451 else |
|
452 output = dae.integrate (out_times, deriv_output); |
2928
|
453 |
5729
|
454 if (fcn_name.length()) |
|
455 clear_function (fcn_name); |
|
456 if (jac_name.length()) |
|
457 clear_function (jac_name); |
|
458 |
3243
|
459 if (! error_state) |
|
460 { |
3997
|
461 std::string msg = dae.error_message (); |
|
462 |
|
463 retval(3) = msg; |
|
464 retval(2) = static_cast<double> (dae.integration_state ()); |
|
465 |
|
466 if (dae.integration_ok ()) |
|
467 { |
|
468 retval(1) = deriv_output; |
|
469 retval(0) = output; |
|
470 } |
|
471 else |
|
472 { |
|
473 retval(1) = Matrix (); |
|
474 retval(0) = Matrix (); |
|
475 |
|
476 if (nargout < 3) |
|
477 error ("dassl: %s", msg.c_str ()); |
|
478 } |
3243
|
479 } |
2928
|
480 } |
3243
|
481 else |
5823
|
482 print_usage (); |
3243
|
483 |
|
484 unwind_protect::run_frame ("Fdassl"); |
2928
|
485 |
|
486 return retval; |
|
487 } |
|
488 |
|
489 /* |
|
490 ;;; Local Variables: *** |
|
491 ;;; mode: C++ *** |
|
492 ;;; End: *** |
|
493 */ |