3912
|
1 /* |
|
2 |
|
3 Copyright (C) 1996, 1997, 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 |
|
23 #ifdef HAVE_CONFIG_H |
|
24 #include <config.h> |
|
25 #endif |
|
26 |
|
27 #include <string> |
|
28 |
|
29 #include <iomanip> |
|
30 #include <iostream> |
|
31 |
|
32 #include "DASPK.h" |
|
33 |
|
34 #include "defun-dld.h" |
|
35 #include "error.h" |
|
36 #include "gripes.h" |
|
37 #include "oct-obj.h" |
|
38 #include "ov-fcn.h" |
|
39 #include "pager.h" |
|
40 #include "unwind-prot.h" |
|
41 #include "utils.h" |
|
42 #include "variables.h" |
|
43 |
|
44 // Global pointer for user defined function required by daspk. |
|
45 static octave_function *daspk_fcn; |
|
46 |
|
47 static DASPK_options daspk_opts; |
|
48 |
|
49 // Is this a recursive call? |
|
50 static int call_depth = 0; |
|
51 |
|
52 ColumnVector |
|
53 daspk_user_function (const ColumnVector& x, const ColumnVector& xdot, |
|
54 double t, int& ires) |
|
55 { |
|
56 ColumnVector retval; |
|
57 |
|
58 int nstates = x.capacity (); |
|
59 |
|
60 assert (nstates == xdot.capacity ()); |
|
61 |
|
62 octave_value_list args; |
|
63 args(2) = t; |
|
64 |
|
65 if (nstates > 1) |
|
66 { |
|
67 Matrix m1 (nstates, 1); |
|
68 Matrix m2 (nstates, 1); |
|
69 for (int i = 0; i < nstates; i++) |
|
70 { |
|
71 m1 (i, 0) = x (i); |
|
72 m2 (i, 0) = xdot (i); |
|
73 } |
|
74 octave_value state (m1); |
|
75 octave_value deriv (m2); |
|
76 args(1) = deriv; |
|
77 args(0) = state; |
|
78 } |
|
79 else |
|
80 { |
|
81 double d1 = x (0); |
|
82 double d2 = xdot (0); |
|
83 octave_value state (d1); |
|
84 octave_value deriv (d2); |
|
85 args(1) = deriv; |
|
86 args(0) = state; |
|
87 } |
|
88 |
|
89 if (daspk_fcn) |
|
90 { |
|
91 octave_value_list tmp = daspk_fcn->do_multi_index_op (1, args); |
|
92 |
|
93 if (error_state) |
|
94 { |
|
95 gripe_user_supplied_eval ("daspk"); |
|
96 return retval; |
|
97 } |
|
98 |
|
99 int tlen = tmp.length (); |
|
100 if (tlen > 0 && tmp(0).is_defined ()) |
|
101 { |
|
102 retval = ColumnVector (tmp(0).vector_value ()); |
|
103 |
|
104 if (tlen > 1) |
|
105 ires = tmp(1).int_value (); |
|
106 |
|
107 if (error_state || retval.length () == 0) |
|
108 gripe_user_supplied_eval ("daspk"); |
|
109 } |
|
110 else |
|
111 gripe_user_supplied_eval ("daspk"); |
|
112 } |
|
113 |
|
114 return retval; |
|
115 } |
|
116 |
|
117 #define DASPK_ABORT() \ |
|
118 do \ |
|
119 { \ |
|
120 unwind_protect::run_frame ("Fdaspk"); \ |
|
121 return retval; \ |
|
122 } \ |
|
123 while (0) |
|
124 |
|
125 #define DASPK_ABORT1(msg) \ |
|
126 do \ |
|
127 { \ |
|
128 ::error ("daspk: " msg); \ |
|
129 DASPK_ABORT (); \ |
|
130 } \ |
|
131 while (0) |
|
132 |
|
133 #define DASPK_ABORT2(fmt, arg) \ |
|
134 do \ |
|
135 { \ |
|
136 ::error ("daspk: " fmt, arg); \ |
|
137 DASPK_ABORT (); \ |
|
138 } \ |
|
139 while (0) |
|
140 |
|
141 DEFUN_DLD (daspk, args, , |
|
142 "-*- texinfo -*-\n\ |
|
143 @deftypefn {Loadable Function} {[@var{x}, @var{xdot}] =} daspk (@var{fcn}, @var{x0}, @var{xdot0}, @var{t}, @var{t_crit})\n\ |
|
144 Return a matrix of states and their first derivatives with respect to\n\ |
|
145 @var{t}. Each row in the result matrices correspond to one of the\n\ |
|
146 elements in the vector @var{t}. The first element of @var{t}\n\ |
|
147 corresponds to the initial state @var{x0} and derivative @var{xdot0}, so\n\ |
|
148 that the first row of the output @var{x} is @var{x0} and the first row\n\ |
|
149 of the output @var{xdot} is @var{xdot0}.\n\ |
|
150 \n\ |
|
151 The first argument, @var{fcn}, is a string that names the function to\n\ |
|
152 call to compute the vector of residuals for the set of equations.\n\ |
|
153 It must have the form\n\ |
|
154 \n\ |
|
155 @example\n\ |
|
156 @var{res} = f (@var{x}, @var{xdot}, @var{t})\n\ |
|
157 @end example\n\ |
|
158 \n\ |
|
159 @noindent\n\ |
|
160 where @var{x}, @var{xdot}, and @var{res} are vectors, and @var{t} is a\n\ |
|
161 scalar.\n\ |
|
162 \n\ |
|
163 The second and third arguments to @code{daspk} specify the initial\n\ |
|
164 condition of the states and their derivatives, and the fourth argument\n\ |
|
165 specifies a vector of output times at which the solution is desired, \n\ |
|
166 including the time corresponding to the initial condition.\n\ |
|
167 \n\ |
|
168 The set of initial states and derivatives are not strictly required to\n\ |
|
169 be consistent. In practice, however, @sc{Dassl} is not very good at\n\ |
|
170 determining a consistent set for you, so it is best if you ensure that\n\ |
|
171 the initial values result in the function evaluating to zero.\n\ |
|
172 \n\ |
|
173 The fifth argument is optional, and may be used to specify a set of\n\ |
|
174 times that the DAE solver should not integrate past. It is useful for\n\ |
|
175 avoiding difficulties with singularities and points where there is a\n\ |
|
176 discontinuity in the derivative.\n\ |
|
177 @end deftypefn") |
|
178 { |
|
179 octave_value_list retval; |
|
180 |
|
181 unwind_protect::begin_frame ("Fdaspk"); |
|
182 |
|
183 unwind_protect_int (call_depth); |
|
184 call_depth++; |
|
185 |
|
186 if (call_depth > 1) |
|
187 DASPK_ABORT1 ("invalid recursive call"); |
|
188 |
|
189 int nargin = args.length (); |
|
190 |
|
191 if (nargin > 3 && nargin < 6) |
|
192 { |
|
193 daspk_fcn = extract_function |
|
194 (args(0), "daspk", "__daspk_fcn__", |
|
195 "function res = __daspk_fcn__ (x, xdot, t) res = ", |
|
196 "; endfunction"); |
|
197 |
|
198 if (! daspk_fcn) |
|
199 DASPK_ABORT (); |
|
200 |
|
201 ColumnVector state = ColumnVector (args(1).vector_value ()); |
|
202 |
|
203 if (error_state) |
|
204 DASPK_ABORT1 ("expecting state vector as second argument"); |
|
205 |
|
206 ColumnVector deriv (args(2).vector_value ()); |
|
207 |
|
208 if (error_state) |
|
209 DASPK_ABORT1 ("expecting derivative vector as third argument"); |
|
210 |
|
211 ColumnVector out_times (args(3).vector_value ()); |
|
212 |
|
213 if (error_state) |
|
214 DASPK_ABORT1 ("expecting output time vector as fourth argument"); |
|
215 |
|
216 ColumnVector crit_times; |
|
217 int crit_times_set = 0; |
|
218 if (nargin > 4) |
|
219 { |
|
220 crit_times = ColumnVector (args(4).vector_value ()); |
|
221 |
|
222 if (error_state) |
|
223 DASPK_ABORT1 ("expecting critical time vector as fifth argument"); |
|
224 |
|
225 crit_times_set = 1; |
|
226 } |
|
227 |
|
228 if (state.capacity () != deriv.capacity ()) |
|
229 DASPK_ABORT1 ("x and xdot must have the same size"); |
|
230 |
|
231 double tzero = out_times (0); |
|
232 |
|
233 DAEFunc func (daspk_user_function); |
|
234 DASPK dae (state, deriv, tzero, func); |
|
235 dae.copy (daspk_opts); |
|
236 |
|
237 Matrix output; |
|
238 Matrix deriv_output; |
|
239 |
|
240 if (crit_times_set) |
|
241 output = dae.integrate (out_times, deriv_output, crit_times); |
|
242 else |
|
243 output = dae.integrate (out_times, deriv_output); |
|
244 |
|
245 if (! error_state) |
|
246 { |
|
247 retval.resize (2); |
|
248 |
|
249 retval(0) = output; |
|
250 retval(1) = deriv_output; |
|
251 } |
|
252 } |
|
253 else |
|
254 print_usage ("daspk"); |
|
255 |
|
256 unwind_protect::run_frame ("Fdaspk"); |
|
257 |
|
258 return retval; |
|
259 } |
|
260 |
|
261 typedef void (DASPK_options::*d_set_opt_mf) (double); |
|
262 typedef double (DASPK_options::*d_get_opt_mf) (void); |
|
263 |
|
264 #define MAX_TOKENS 3 |
|
265 |
|
266 struct DASPK_OPTIONS |
|
267 { |
|
268 const char *keyword; |
|
269 const char *kw_tok[MAX_TOKENS + 1]; |
|
270 int min_len[MAX_TOKENS + 1]; |
|
271 int min_toks_to_match; |
|
272 d_set_opt_mf d_set_fcn; |
|
273 d_get_opt_mf d_get_fcn; |
|
274 }; |
|
275 |
|
276 static DASPK_OPTIONS daspk_option_table [] = |
|
277 { |
|
278 { "absolute tolerance", |
|
279 { "absolute", "tolerance", 0, 0, }, |
|
280 { 1, 0, 0, 0, }, 1, |
|
281 &DASPK_options::set_absolute_tolerance, |
|
282 &DASPK_options::absolute_tolerance, }, |
|
283 |
|
284 { "initial step size", |
|
285 { "initial", "step", "size", 0, }, |
|
286 { 1, 0, 0, 0, }, 1, |
|
287 &DASPK_options::set_initial_step_size, |
|
288 &DASPK_options::initial_step_size, }, |
|
289 |
|
290 { "maximum step size", |
|
291 { "maximum", "step", "size", 0, }, |
|
292 { 2, 0, 0, 0, }, 1, |
|
293 &DASPK_options::set_maximum_step_size, |
|
294 &DASPK_options::maximum_step_size, }, |
|
295 |
|
296 { "relative tolerance", |
|
297 { "relative", "tolerance", 0, 0, }, |
|
298 { 1, 0, 0, 0, }, 1, |
|
299 &DASPK_options::set_relative_tolerance, |
|
300 &DASPK_options::relative_tolerance, }, |
|
301 |
|
302 { 0, |
|
303 { 0, 0, 0, 0, }, |
|
304 { 0, 0, 0, 0, }, 0, |
|
305 0, 0, }, |
|
306 }; |
|
307 |
|
308 static void |
|
309 print_daspk_option_list (std::ostream& os) |
|
310 { |
|
311 print_usage ("daspk_options", 1); |
|
312 |
|
313 os << "\n" |
|
314 << "Options for daspk include:\n\n" |
|
315 << " keyword value\n" |
|
316 << " ------- -----\n\n"; |
|
317 |
|
318 DASPK_OPTIONS *list = daspk_option_table; |
|
319 |
|
320 const char *keyword; |
|
321 while ((keyword = list->keyword) != 0) |
|
322 { |
|
323 os << " " |
|
324 << std::setiosflags (std::ios::left) << std::setw (40) |
|
325 << keyword |
|
326 << std::resetiosflags (std::ios::left) |
|
327 << " "; |
|
328 |
|
329 double val = (daspk_opts.*list->d_get_fcn) (); |
|
330 if (val < 0.0) |
|
331 os << "computed automatically"; |
|
332 else |
|
333 os << val; |
|
334 |
|
335 os << "\n"; |
|
336 list++; |
|
337 } |
|
338 |
|
339 os << "\n"; |
|
340 } |
|
341 |
|
342 static void |
|
343 set_daspk_option (const std::string& keyword, double val) |
|
344 { |
|
345 DASPK_OPTIONS *list = daspk_option_table; |
|
346 |
|
347 while (list->keyword != 0) |
|
348 { |
|
349 if (keyword_almost_match (list->kw_tok, list->min_len, keyword, |
|
350 list->min_toks_to_match, MAX_TOKENS)) |
|
351 { |
|
352 (daspk_opts.*list->d_set_fcn) (val); |
|
353 |
|
354 return; |
|
355 } |
|
356 list++; |
|
357 } |
|
358 |
|
359 warning ("daspk_options: no match for `%s'", keyword.c_str ()); |
|
360 } |
|
361 |
|
362 static octave_value_list |
|
363 show_daspk_option (const std::string& keyword) |
|
364 { |
|
365 octave_value retval; |
|
366 |
|
367 DASPK_OPTIONS *list = daspk_option_table; |
|
368 |
|
369 while (list->keyword != 0) |
|
370 { |
|
371 if (keyword_almost_match (list->kw_tok, list->min_len, keyword, |
|
372 list->min_toks_to_match, MAX_TOKENS)) |
|
373 { |
|
374 double val = (daspk_opts.*list->d_get_fcn) (); |
|
375 if (val < 0.0) |
|
376 retval = "computed automatically"; |
|
377 else |
|
378 retval = val; |
|
379 |
|
380 return retval; |
|
381 } |
|
382 list++; |
|
383 } |
|
384 |
|
385 warning ("daspk_options: no match for `%s'", keyword.c_str ()); |
|
386 |
|
387 return retval; |
|
388 } |
|
389 |
|
390 DEFUN_DLD (daspk_options, args, , |
|
391 "-*- texinfo -*-\n\ |
|
392 @deftypefn {Loadable Function} {} daspk_options (@var{opt}, @var{val})\n\ |
|
393 When called with two arguments, this function allows you set options\n\ |
|
394 parameters for the function @code{lsode}. Given one argument,\n\ |
|
395 @code{daspk_options} returns the value of the corresponding option. If\n\ |
|
396 no arguments are supplied, the names of all the available options and\n\ |
|
397 their current values are displayed.\n\ |
|
398 @end deftypefn") |
|
399 { |
|
400 octave_value_list retval; |
|
401 |
|
402 int nargin = args.length (); |
|
403 |
|
404 if (nargin == 0) |
|
405 { |
|
406 print_daspk_option_list (octave_stdout); |
|
407 return retval; |
|
408 } |
|
409 else if (nargin == 1 || nargin == 2) |
|
410 { |
|
411 std::string keyword = args(0).string_value (); |
|
412 |
|
413 if (! error_state) |
|
414 { |
|
415 if (nargin == 1) |
|
416 return show_daspk_option (keyword); |
|
417 else |
|
418 { |
|
419 double val = args(1).double_value (); |
|
420 |
|
421 if (! error_state) |
|
422 { |
|
423 set_daspk_option (keyword, val); |
|
424 return retval; |
|
425 } |
|
426 } |
|
427 } |
|
428 } |
|
429 |
|
430 print_usage ("daspk_options"); |
|
431 |
|
432 return retval; |
|
433 } |
|
434 |
|
435 /* |
|
436 ;;; Local Variables: *** |
|
437 ;;; mode: C++ *** |
|
438 ;;; End: *** |
|
439 */ |