File: diffeq.texi

package info (click to toggle)
octave 4.0.3-3
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 94,200 kB
  • ctags: 52,925
  • sloc: cpp: 316,850; ansic: 43,469; fortran: 23,670; sh: 13,805; yacc: 8,204; objc: 7,939; lex: 3,631; java: 2,127; makefile: 1,746; perl: 1,022; awk: 988
file content (992 lines) | stat: -rw-r--r-- 31,103 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
@c DO NOT EDIT!  Generated automatically by munge-texi.pl.

@c Copyright (C) 1996-2015 John W. Eaton
@c
@c This file is part of Octave.
@c
@c Octave is free software; you can redistribute it and/or modify it
@c under the terms of the GNU General Public License as published by the
@c Free Software Foundation; either version 3 of the License, or (at
@c your option) any later version.
@c
@c Octave is distributed in the hope that it will be useful, but WITHOUT
@c ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
@c FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
@c for more details.
@c
@c You should have received a copy of the GNU General Public License
@c along with Octave; see the file COPYING.  If not, see
@c <http://www.gnu.org/licenses/>.

@node Differential Equations
@chapter Differential Equations

Octave has built-in functions for solving ordinary differential equations,
and differential-algebraic equations.
All solvers are based on reliable ODE routines written in Fortran.

@menu
* Ordinary Differential Equations::
* Differential-Algebraic Equations::
@end menu

@cindex differential equations
@cindex ODE
@cindex DAE

@node Ordinary Differential Equations
@section Ordinary Differential Equations

The function @code{lsode} can be used to solve ODEs of the form
@tex
$$
 {dx\over dt} = f (x, t)
$$
@end tex
@ifnottex

@example
@group
dx
-- = f (x, t)
dt
@end group
@end example

@end ifnottex

@noindent
using @nospell{Hindmarsh's} ODE solver @sc{lsode}.



@c lsode libinterp/corefcn/lsode.cc
@anchor{XREFlsode}
@deftypefn  {Built-in Function} {[@var{x}, @var{istate}, @var{msg}] =} lsode (@var{fcn}, @var{x_0}, @var{t})
@deftypefnx {Built-in Function} {[@var{x}, @var{istate}, @var{msg}] =} lsode (@var{fcn}, @var{x_0}, @var{t}, @var{t_crit})
Ordinary Differential Equation (ODE) solver.

The set of differential equations to solve is
@tex
$$ {dx \over dt} = f (x, t) $$
with
$$ x(t_0) = x_0 $$
@end tex
@ifnottex

@example
@group
dx
-- = f (x, t)
dt
@end group
@end example

@noindent
with

@example
x(t_0) = x_0
@end example

@end ifnottex
The solution is returned in the matrix @var{x}, with each row
corresponding to an element of the vector @var{t}.  The first element
of @var{t} should be @math{t_0} and should correspond to the initial
state of the system @var{x_0}, so that the first row of the output
is @var{x_0}.

The first argument, @var{fcn}, is a string, inline, or function handle
that names the function @math{f} to call to compute the vector of right
hand sides for the set of equations.  The function must have the form

@example
@var{xdot} = f (@var{x}, @var{t})
@end example

@noindent
in which @var{xdot} and @var{x} are vectors and @var{t} is a scalar.

If @var{fcn} is a two-element string array or a two-element cell array
of strings, inline functions, or function handles, the first element names
the function @math{f} described above, and the second element names a
function to compute the Jacobian of @math{f}.  The Jacobian function
must have the form

@example
@var{jac} = j (@var{x}, @var{t})
@end example

@noindent
in which @var{jac} is the matrix of partial derivatives
@tex
$$ J = {\partial f_i \over \partial x_j} = \left[\matrix{
{\partial f_1 \over \partial x_1}
  & {\partial f_1 \over \partial x_2}
  & \cdots
  & {\partial f_1 \over \partial x_N} \cr
{\partial f_2 \over \partial x_1}
  & {\partial f_2 \over \partial x_2}
  & \cdots
  & {\partial f_2 \over \partial x_N} \cr
 \vdots & \vdots & \ddots & \vdots \cr
{\partial f_3 \over \partial x_1}
  & {\partial f_3 \over \partial x_2}
  & \cdots
  & {\partial f_3 \over \partial x_N} \cr}\right]$$
@end tex
@ifnottex

@example
@group
             | df_1  df_1       df_1 |
             | ----  ----  ...  ---- |
             | dx_1  dx_2       dx_N |
             |                       |
             | df_2  df_2       df_2 |
             | ----  ----  ...  ---- |
      df_i   | dx_1  dx_2       dx_N |
jac = ---- = |                       |
      dx_j   |  .    .     .    .    |
             |  .    .      .   .    |
             |  .    .       .  .    |
             |                       |
             | df_N  df_N       df_N |
             | ----  ----  ...  ---- |
             | dx_1  dx_2       dx_N |
@end group
@end example

@end ifnottex

The second and third arguments specify the initial state of the system,
@math{x_0}, and the initial value of the independent variable @math{t_0}.

The fourth argument is optional, and may be used to specify a set of
times that the ODE solver should not integrate past.  It is useful for
avoiding difficulties with singularities and points where there is a
discontinuity in the derivative.

After a successful computation, the value of @var{istate} will be 2
(consistent with the Fortran version of @sc{lsode}).

If the computation is not successful, @var{istate} will be something
other than 2 and @var{msg} will contain additional information.

You can use the function @code{lsode_options} to set optional
parameters for @code{lsode}.
@seealso{@ref{XREFdaspk,,daspk}, @ref{XREFdassl,,dassl}, @ref{XREFdasrt,,dasrt}}
@end deftypefn


@c lsode_options libinterp/corefcn/lsode.cc
@anchor{XREFlsode_options}
@deftypefn  {Built-in Function} {} lsode_options ()
@deftypefnx {Built-in Function} {val =} lsode_options (@var{opt})
@deftypefnx {Built-in Function} {} lsode_options (@var{opt}, @var{val})
Query or set options for the function @code{lsode}.

When called with no arguments, the names of all available options and
their current values are displayed.

Given one argument, return the value of the option @var{opt}.

When called with two arguments, @code{lsode_options} sets the option
@var{opt} to value @var{val}.

Options include

@table @code
@item @qcode{"absolute tolerance"}
Absolute tolerance.  May be either vector or scalar.  If a vector, it
must match the dimension of the state vector.

@item @qcode{"relative tolerance"}
Relative tolerance parameter.  Unlike the absolute tolerance, this
parameter may only be a scalar.

The local error test applied at each integration step is

@example
@group
  abs (local error in x(i)) <= ...
      rtol * abs (y(i)) + atol(i)
@end group
@end example

@item @qcode{"integration method"}
A string specifying the method of integration to use to solve the ODE
system.  Valid values are

@table @asis
@item  @qcode{"adams"}
@itemx @qcode{"non-stiff"}
No Jacobian used (even if it is available).

@item  @qcode{"bdf"}
@itemx @qcode{"stiff"}
Use stiff backward differentiation formula (BDF) method.  If a
function to compute the Jacobian is not supplied, @code{lsode} will
compute a finite difference approximation of the Jacobian matrix.
@end table

@item @qcode{"initial step size"}
The step size to be attempted on the first step (default is determined
automatically).

@item @qcode{"maximum order"}
Restrict the maximum order of the solution method.  If using the Adams
method, this option must be between 1 and 12.  Otherwise, it must be
between 1 and 5, inclusive.

@item @qcode{"maximum step size"}
Setting the maximum stepsize will avoid passing over very large
regions  (default is not specified).

@item @qcode{"minimum step size"}
The minimum absolute step size allowed (default is 0).

@item @qcode{"step limit"}
Maximum number of steps allowed (default is 100000).
@end table
@end deftypefn


Here is an example of solving a set of three differential equations using
@code{lsode}.  Given the function

@cindex oregonator

@example
@group
## oregonator differential equation
function xdot = f (x, t)

  xdot = zeros (3,1);

  xdot(1) = 77.27 * (x(2) - x(1)*x(2) + x(1) \
            - 8.375e-06*x(1)^2);
  xdot(2) = (x(3) - x(1)*x(2) - x(2)) / 77.27;
  xdot(3) = 0.161*(x(1) - x(3));

endfunction
@end group
@end example

@noindent
and the initial condition @code{x0 = [ 4; 1.1; 4 ]}, the set of
equations can be integrated using the command

@example
@group
t = linspace (0, 500, 1000);

y = lsode ("f", x0, t);
@end group
@end example

If you try this, you will see that the value of the result changes
dramatically between @var{t} = 0 and 5, and again around @var{t} = 305.
A more efficient set of output points might be

@example
@group
t = [0, logspace(-1, log10(303), 150), \
        logspace(log10(304), log10(500), 150)];
@end group
@end example

See @nospell{Alan C. Hindmarsh},
@cite{ODEPACK, A Systematized Collection of ODE Solvers},
in Scientific Computing, @nospell{R. S. Stepleman}, editor, (1983)
for more information about the inner workings of @code{lsode}.

An m-file for the differential equation used above is included with the
Octave distribution in the examples directory under the name
@file{oregonator.m}.

@node Differential-Algebraic Equations
@section Differential-Algebraic Equations

The function @code{daspk} can be used to solve DAEs of the form
@tex
$$
 0 = f (\dot{x}, x, t), \qquad x(t=0) = x_0, \dot{x}(t=0) = \dot{x}_0
$$
@end tex
@ifnottex

@example
0 = f (x-dot, x, t),    x(t=0) = x_0, x-dot(t=0) = x-dot_0
@end example

@end ifnottex

@noindent
where
@tex
$\dot{x} = {dx \over dt}$
@end tex
@ifnottex
@math{x-dot}
@end ifnottex
is the derivative of @math{x}.  The equation is solved using
@nospell{Petzold's} DAE solver @sc{daspk}.

@c daspk libinterp/corefcn/daspk.cc
@anchor{XREFdaspk}
@deftypefn {Built-in Function} {[@var{x}, @var{xdot}, @var{istate}, @var{msg}] =} daspk (@var{fcn}, @var{x_0}, @var{xdot_0}, @var{t}, @var{t_crit})
Solve the set of differential-algebraic equations
@tex
$$ 0 = f (x, \dot{x}, t) $$
with
$$ x(t_0) = x_0, \dot{x}(t_0) = \dot{x}_0 $$
@end tex
@ifnottex

@example
0 = f (x, xdot, t)
@end example

@noindent
with

@example
x(t_0) = x_0, xdot(t_0) = xdot_0
@end example

@end ifnottex
The solution is returned in the matrices @var{x} and @var{xdot},
with each row in the result matrices corresponding to one of the
elements in the vector @var{t}.  The first element of @var{t}
should be @math{t_0} and correspond to the initial state of the
system @var{x_0} and its derivative @var{xdot_0}, so that the first
row of the output @var{x} is @var{x_0} and the first row
of the output @var{xdot} is @var{xdot_0}.

The first argument, @var{fcn}, is a string, inline, or function handle
that names the function @math{f} to call to compute the vector of
residuals for the set of equations.  It must have the form

@example
@var{res} = f (@var{x}, @var{xdot}, @var{t})
@end example

@noindent
in which @var{x}, @var{xdot}, and @var{res} are vectors, and @var{t} is a
scalar.

If @var{fcn} is a two-element string array or a two-element cell array
of strings, inline functions, or function handles, the first element names
the function @math{f} described above, and the second element names a
function to compute the modified Jacobian
@tex
$$
J = {\partial f \over \partial x}
  + c {\partial f \over \partial \dot{x}}
$$
@end tex
@ifnottex

@example
@group
      df       df
jac = -- + c ------
      dx     d xdot
@end group
@end example

@end ifnottex

The modified Jacobian function must have the form

@example
@group

@var{jac} = j (@var{x}, @var{xdot}, @var{t}, @var{c})

@end group
@end example

The second and third arguments to @code{daspk} specify the initial
condition of the states and their derivatives, and the fourth argument
specifies a vector of output times at which the solution is desired,
including the time corresponding to the initial condition.

The set of initial states and derivatives are not strictly required to
be consistent.  If they are not consistent, you must use the
@code{daspk_options} function to provide additional information so
that @code{daspk} can compute a consistent starting point.

The fifth argument is optional, and may be used to specify a set of
times that the DAE solver should not integrate past.  It is useful for
avoiding difficulties with singularities and points where there is a
discontinuity in the derivative.

After a successful computation, the value of @var{istate} will be
greater than zero (consistent with the Fortran version of @sc{daspk}).

If the computation is not successful, the value of @var{istate} will be
less than zero and @var{msg} will contain additional information.

You can use the function @code{daspk_options} to set optional
parameters for @code{daspk}.
@seealso{@ref{XREFdassl,,dassl}}
@end deftypefn


@c daspk_options libinterp/corefcn/daspk.cc
@anchor{XREFdaspk_options}
@deftypefn  {Built-in Function} {} daspk_options ()
@deftypefnx {Built-in Function} {val =} daspk_options (@var{opt})
@deftypefnx {Built-in Function} {} daspk_options (@var{opt}, @var{val})
Query or set options for the function @code{daspk}.

When called with no arguments, the names of all available options and
their current values are displayed.

Given one argument, return the value of the option @var{opt}.

When called with two arguments, @code{daspk_options} sets the option
@var{opt} to value @var{val}.

Options include

@table @code
@item @qcode{"absolute tolerance"}
Absolute tolerance.  May be either vector or scalar.  If a vector, it
must match the dimension of the state vector, and the relative
tolerance must also be a vector of the same length.

@item @qcode{"relative tolerance"}
Relative tolerance.  May be either vector or scalar.  If a vector, it
must match the dimension of the state vector, and the absolute
tolerance must also be a vector of the same length.

The local error test applied at each integration step is

@example
@group
  abs (local error in x(i))
       <= rtol(i) * abs (Y(i)) + atol(i)
@end group
@end example

@item @qcode{"compute consistent initial condition"}
Denoting the differential variables in the state vector by @samp{Y_d}
and the algebraic variables by @samp{Y_a}, @code{ddaspk} can solve
one of two initialization problems:

@enumerate
@item Given Y_d, calculate Y_a and Y'_d

@item Given Y', calculate Y.
@end enumerate

In either case, initial values for the given components are input, and
initial guesses for the unknown components must also be provided as
input.  Set this option to 1 to solve the first problem, or 2 to solve
the second (the default is 0, so you must provide a set of
initial conditions that are consistent).

If this option is set to a nonzero value, you must also set the
@qcode{"algebraic variables"} option to declare which variables in the
problem are algebraic.

@item @qcode{"use initial condition heuristics"}
Set to a nonzero value to use the initial condition heuristics options
described below.

@item @qcode{"initial condition heuristics"}
A vector of the following parameters that can be used to control the
initial condition calculation.

@table @code
@item MXNIT
Maximum number of Newton iterations (default is 5).

@item MXNJ
Maximum number of Jacobian evaluations (default is 6).

@item MXNH
Maximum number of values of the artificial stepsize parameter to be
tried if the @qcode{"compute consistent initial condition"} option has
been set to 1 (default is 5).

Note that the maximum total number of Newton iterations allowed is
@code{MXNIT*MXNJ*MXNH} if the @qcode{"compute consistent initial
condition"} option has been set to 1 and @code{MXNIT*MXNJ} if it is
set to 2.

@item LSOFF
Set to a nonzero value to disable the linesearch algorithm (default is
0).

@item STPTOL
Minimum scaled step in linesearch algorithm (default is eps^(2/3)).

@item EPINIT
Swing factor in the Newton iteration convergence test.  The test is
applied to the residual vector, premultiplied by the approximate
Jacobian.  For convergence, the weighted RMS norm of this vector
(scaled by the error weights) must be less than @code{EPINIT*EPCON},
where @code{EPCON} = 0.33 is the analogous test constant used in the
time steps.  The default is @code{EPINIT} = 0.01.
@end table

@item @qcode{"print initial condition info"}
Set this option to a nonzero value to display detailed information
about the initial condition calculation (default is 0).

@item @qcode{"exclude algebraic variables from error test"}
Set to a nonzero value to exclude algebraic variables from the error
test.  You must also set the @qcode{"algebraic variables"} option to
declare which variables in the problem are algebraic (default is 0).

@item @qcode{"algebraic variables"}
A vector of the same length as the state vector.  A nonzero element
indicates that the corresponding element of the state vector is an
algebraic variable (i.e., its derivative does not appear explicitly
in the equation set).

This option is required by the
@qcode{"compute consistent initial condition"} and
@qcode{"exclude algebraic variables from error test"} options.

@item @qcode{"enforce inequality constraints"}
Set to one of the following values to enforce the inequality
constraints specified by the @qcode{"inequality constraint types"}
option (default is 0).

@enumerate
@item To have constraint checking only in the initial condition calculation.

@item To enforce constraint checking during the integration.

@item To enforce both options 1 and 2.
@end enumerate

@item @qcode{"inequality constraint types"}
A vector of the same length as the state specifying the type of
inequality constraint.  Each element of the vector corresponds to an
element of the state and should be assigned one of the following
codes

@table @asis
@item -2
Less than zero.

@item -1
Less than or equal to zero.

@item 0
Not constrained.

@item 1
Greater than or equal to zero.

@item 2
Greater than zero.
@end table

This option only has an effect if the
@qcode{"enforce inequality constraints"} option is nonzero.

@item @qcode{"initial step size"}
Differential-algebraic problems may occasionally suffer from severe
scaling difficulties on the first step.  If you know a great deal
about the scaling of your problem, you can help to alleviate this
problem by specifying an initial stepsize (default is computed
automatically).

@item @qcode{"maximum order"}
Restrict the maximum order of the solution method.  This option must
be between 1 and 5, inclusive (default is 5).

@item @qcode{"maximum step size"}
Setting the maximum stepsize will avoid passing over very large
regions (default is not specified).
@end table
@end deftypefn


Octave also includes @sc{dassl}, an earlier version of @sc{daspk},
and @sc{dasrt}, which can be used to solve DAEs with constraints
(stopping conditions).

@c dassl libinterp/corefcn/dassl.cc
@anchor{XREFdassl}
@deftypefn {Built-in Function} {[@var{x}, @var{xdot}, @var{istate}, @var{msg}] =} dassl (@var{fcn}, @var{x_0}, @var{xdot_0}, @var{t}, @var{t_crit})
Solve the set of differential-algebraic equations
@tex
$$ 0 = f (x, \dot{x}, t) $$
with
$$ x(t_0) = x_0, \dot{x}(t_0) = \dot{x}_0 $$
@end tex
@ifnottex

@example
0 = f (x, xdot, t)
@end example

@noindent
with

@example
x(t_0) = x_0, xdot(t_0) = xdot_0
@end example

@end ifnottex
The solution is returned in the matrices @var{x} and @var{xdot},
with each row in the result matrices corresponding to one of the
elements in the vector @var{t}.  The first element of @var{t}
should be @math{t_0} and correspond to the initial state of the
system @var{x_0} and its derivative @var{xdot_0}, so that the first
row of the output @var{x} is @var{x_0} and the first row
of the output @var{xdot} is @var{xdot_0}.

The first argument, @var{fcn}, is a string, inline, or function handle
that names the function @math{f} to call to compute the vector of
residuals for the set of equations.  It must have the form

@example
@var{res} = f (@var{x}, @var{xdot}, @var{t})
@end example

@noindent
in which @var{x}, @var{xdot}, and @var{res} are vectors, and @var{t} is a
scalar.

If @var{fcn} is a two-element string array or a two-element cell array
of strings, inline functions, or function handles, the first element names
the function @math{f} described above, and the second element names a
function to compute the modified Jacobian

@tex
$$
J = {\partial f \over \partial x}
  + c {\partial f \over \partial \dot{x}}
$$
@end tex
@ifnottex

@example
@group
      df       df
jac = -- + c ------
      dx     d xdot
@end group
@end example

@end ifnottex

The modified Jacobian function must have the form

@example
@group

@var{jac} = j (@var{x}, @var{xdot}, @var{t}, @var{c})

@end group
@end example

The second and third arguments to @code{dassl} specify the initial
condition of the states and their derivatives, and the fourth argument
specifies a vector of output times at which the solution is desired,
including the time corresponding to the initial condition.

The set of initial states and derivatives are not strictly required to
be consistent.  In practice, however, @sc{dassl} is not very good at
determining a consistent set for you, so it is best if you ensure that
the initial values result in the function evaluating to zero.

The fifth argument is optional, and may be used to specify a set of
times that the DAE solver should not integrate past.  It is useful for
avoiding difficulties with singularities and points where there is a
discontinuity in the derivative.

After a successful computation, the value of @var{istate} will be
greater than zero (consistent with the Fortran version of @sc{dassl}).

If the computation is not successful, the value of @var{istate} will be
less than zero and @var{msg} will contain additional information.

You can use the function @code{dassl_options} to set optional
parameters for @code{dassl}.
@seealso{@ref{XREFdaspk,,daspk}, @ref{XREFdasrt,,dasrt}, @ref{XREFlsode,,lsode}}
@end deftypefn


@c dassl_options libinterp/corefcn/dassl.cc
@anchor{XREFdassl_options}
@deftypefn  {Built-in Function} {} dassl_options ()
@deftypefnx {Built-in Function} {val =} dassl_options (@var{opt})
@deftypefnx {Built-in Function} {} dassl_options (@var{opt}, @var{val})
Query or set options for the function @code{dassl}.

When called with no arguments, the names of all available options and
their current values are displayed.

Given one argument, return the value of the option @var{opt}.

When called with two arguments, @code{dassl_options} sets the option
@var{opt} to value @var{val}.

Options include

@table @code
@item @qcode{"absolute tolerance"}
Absolute tolerance.  May be either vector or scalar.  If a vector, it
must match the dimension of the state vector, and the relative
tolerance must also be a vector of the same length.

@item @qcode{"relative tolerance"}
Relative tolerance.  May be either vector or scalar.  If a vector, it
must match the dimension of the state vector, and the absolute
tolerance must also be a vector of the same length.

The local error test applied at each integration step is

@example
@group
  abs (local error in x(i))
       <= rtol(i) * abs (Y(i)) + atol(i)
@end group
@end example

@item @qcode{"compute consistent initial condition"}
If nonzero, @code{dassl} will attempt to compute a consistent set of initial
conditions.  This is generally not reliable, so it is best to provide
a consistent set and leave this option set to zero.

@item @qcode{"enforce nonnegativity constraints"}
If you know that the solutions to your equations will always be
non-negative, it may help to set this parameter to a nonzero
value.  However, it is probably best to try leaving this option set to
zero first, and only setting it to a nonzero value if that doesn't
work very well.

@item @qcode{"initial step size"}
Differential-algebraic problems may occasionally suffer from severe
scaling difficulties on the first step.  If you know a great deal
about the scaling of your problem, you can help to alleviate this
problem by specifying an initial stepsize.

@item @qcode{"maximum order"}
Restrict the maximum order of the solution method.  This option must
be between 1 and 5, inclusive.

@item @qcode{"maximum step size"}
Setting the maximum stepsize will avoid passing over very large
regions  (default is not specified).

@item @qcode{"step limit"}
Maximum number of integration steps to attempt on a single call to the
underlying Fortran code.
@end table
@end deftypefn


@c dasrt libinterp/corefcn/dasrt.cc
@anchor{XREFdasrt}
@deftypefn  {Built-in Function} {[@var{x}, @var{xdot}, @var{t_out}, @var{istat}, @var{msg}] =} dasrt (@var{fcn}, [], @var{x_0}, @var{xdot_0}, @var{t})
@deftypefnx {Built-in Function} {@dots{} =} dasrt (@var{fcn}, @var{g}, @var{x_0}, @var{xdot_0}, @var{t})
@deftypefnx {Built-in Function} {@dots{} =} dasrt (@var{fcn}, [], @var{x_0}, @var{xdot_0}, @var{t}, @var{t_crit})
@deftypefnx {Built-in Function} {@dots{} =} dasrt (@var{fcn}, @var{g}, @var{x_0}, @var{xdot_0}, @var{t}, @var{t_crit})
Solve the set of differential-algebraic equations
@tex
$$ 0 = f (x, \dot{x}, t) $$
with
$$ x(t_0) = x_0, \dot{x}(t_0) = \dot{x}_0 $$
@end tex
@ifnottex

@example
0 = f (x, xdot, t)
@end example

@noindent
with

@example
x(t_0) = x_0, xdot(t_0) = xdot_0
@end example

@end ifnottex
with functional stopping criteria (root solving).

The solution is returned in the matrices @var{x} and @var{xdot},
with each row in the result matrices corresponding to one of the
elements in the vector @var{t_out}.  The first element of @var{t}
should be @math{t_0} and correspond to the initial state of the
system @var{x_0} and its derivative @var{xdot_0}, so that the first
row of the output @var{x} is @var{x_0} and the first row
of the output @var{xdot} is @var{xdot_0}.

The vector @var{t} provides an upper limit on the length of the
integration.  If the stopping condition is met, the vector
@var{t_out} will be shorter than @var{t}, and the final element of
@var{t_out} will be the point at which the stopping condition was met,
and may not correspond to any element of the vector @var{t}.

The first argument, @var{fcn}, is a string, inline, or function handle
that names the function @math{f} to call to compute the vector of
residuals for the set of equations.  It must have the form

@example
@var{res} = f (@var{x}, @var{xdot}, @var{t})
@end example

@noindent
in which @var{x}, @var{xdot}, and @var{res} are vectors, and @var{t} is a
scalar.

If @var{fcn} is a two-element string array or a two-element cell array
of strings, inline functions, or function handles, the first element names
the function @math{f} described above, and the second element names a
function to compute the modified Jacobian

@tex
$$
J = {\partial f \over \partial x}
  + c {\partial f \over \partial \dot{x}}
$$
@end tex
@ifnottex

@example
@group
      df       df
jac = -- + c ------
      dx     d xdot
@end group
@end example

@end ifnottex

The modified Jacobian function must have the form

@example
@group

@var{jac} = j (@var{x}, @var{xdot}, @var{t}, @var{c})

@end group
@end example

The optional second argument names a function that defines the
constraint functions whose roots are desired during the integration.
This function must have the form

@example
@var{g_out} = g (@var{x}, @var{t})
@end example

@noindent
and return a vector of the constraint function values.
If the value of any of the constraint functions changes sign, @sc{dasrt}
will attempt to stop the integration at the point of the sign change.

If the name of the constraint function is omitted, @code{dasrt} solves
the same problem as @code{daspk} or @code{dassl}.

Note that because of numerical errors in the constraint functions
due to round-off and integration error, @sc{dasrt} may return false
roots, or return the same root at two or more nearly equal values of
@var{T}.  If such false roots are suspected, the user should consider
smaller error tolerances or higher precision in the evaluation of the
constraint functions.

If a root of some constraint function defines the end of the problem,
the input to @sc{dasrt} should nevertheless allow integration to a
point slightly past that root, so that @sc{dasrt} can locate the root
by interpolation.

The third and fourth arguments to @code{dasrt} specify the initial
condition of the states and their derivatives, and the fourth argument
specifies a vector of output times at which the solution is desired,
including the time corresponding to the initial condition.

The set of initial states and derivatives are not strictly required to
be consistent.  In practice, however, @sc{dassl} is not very good at
determining a consistent set for you, so it is best if you ensure that
the initial values result in the function evaluating to zero.

The sixth argument is optional, and may be used to specify a set of
times that the DAE solver should not integrate past.  It is useful for
avoiding difficulties with singularities and points where there is a
discontinuity in the derivative.

After a successful computation, the value of @var{istate} will be
greater than zero (consistent with the Fortran version of @sc{dassl}).

If the computation is not successful, the value of @var{istate} will be
less than zero and @var{msg} will contain additional information.

You can use the function @code{dasrt_options} to set optional
parameters for @code{dasrt}.
@seealso{@ref{XREFdasrt_options,,dasrt_options}, @ref{XREFdaspk,,daspk}, @ref{XREFdasrt,,dasrt}, @ref{XREFlsode,,lsode}}
@end deftypefn


@c dasrt_options libinterp/corefcn/dasrt.cc
@anchor{XREFdasrt_options}
@deftypefn  {Built-in Function} {} dasrt_options ()
@deftypefnx {Built-in Function} {val =} dasrt_options (@var{opt})
@deftypefnx {Built-in Function} {} dasrt_options (@var{opt}, @var{val})
Query or set options for the function @code{dasrt}.

When called with no arguments, the names of all available options and
their current values are displayed.

Given one argument, return the value of the option @var{opt}.

When called with two arguments, @code{dasrt_options} sets the option
@var{opt} to value @var{val}.

Options include

@table @code
@item @qcode{"absolute tolerance"}
Absolute tolerance.  May be either vector or scalar.  If a vector, it
must match the dimension of the state vector, and the relative
tolerance must also be a vector of the same length.

@item @qcode{"relative tolerance"}
Relative tolerance.  May be either vector or scalar.  If a vector, it
must match the dimension of the state vector, and the absolute
tolerance must also be a vector of the same length.

The local error test applied at each integration step is

@example
@group
  abs (local error in x(i)) <= ...
      rtol(i) * abs (Y(i)) + atol(i)
@end group
@end example

@item @qcode{"initial step size"}
Differential-algebraic problems may occasionally suffer from severe
scaling difficulties on the first step.  If you know a great deal
about the scaling of your problem, you can help to alleviate this
problem by specifying an initial stepsize.

@item @qcode{"maximum order"}
Restrict the maximum order of the solution method.  This option must
be between 1 and 5, inclusive.

@item @qcode{"maximum step size"}
Setting the maximum stepsize will avoid passing over very large
regions.

@item @qcode{"step limit"}
Maximum number of integration steps to attempt on a single call to the
underlying Fortran code.
@end table
@end deftypefn


See @nospell{K. E. Brenan}, et al., @cite{Numerical Solution of Initial-Value
Problems in Differential-Algebraic Equations}, North-Holland (1989) for
more information about the implementation of @sc{dassl}.