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
|
<!DOCTYPE html>
<html>
<!-- Created by GNU Texinfo 7.1.1, https://www.gnu.org/software/texinfo/ -->
<head>
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
<title>Differential-Algebraic Equations (GNU Octave (version 10.3.0))</title>
<meta name="description" content="Differential-Algebraic Equations (GNU Octave (version 10.3.0))">
<meta name="keywords" content="Differential-Algebraic Equations (GNU Octave (version 10.3.0))">
<meta name="resource-type" content="document">
<meta name="distribution" content="global">
<meta name="Generator" content="makeinfo">
<meta name="viewport" content="width=device-width,initial-scale=1">
<link href="index.html" rel="start" title="Top">
<link href="Concept-Index.html" rel="index" title="Concept Index">
<link href="index.html#SEC_Contents" rel="contents" title="Table of Contents">
<link href="Differential-Equations.html" rel="up" title="Differential Equations">
<link href="Matlab_002dcompatible-solvers.html" rel="next" title="Matlab-compatible solvers">
<link href="Ordinary-Differential-Equations.html" rel="prev" title="Ordinary Differential Equations">
<style type="text/css">
<!--
a.copiable-link {visibility: hidden; text-decoration: none; line-height: 0em}
div.example {margin-left: 3.2em}
span:hover a.copiable-link {visibility: visible}
strong.def-name {font-family: monospace; font-weight: bold; font-size: larger}
-->
</style>
<link rel="stylesheet" type="text/css" href="octave.css">
</head>
<body lang="en">
<div class="section-level-extent" id="Differential_002dAlgebraic-Equations">
<div class="nav-panel">
<p>
Next: <a href="Matlab_002dcompatible-solvers.html" accesskey="n" rel="next">Matlab-compatible solvers</a>, Previous: <a href="Ordinary-Differential-Equations.html" accesskey="p" rel="prev">Ordinary Differential Equations</a>, Up: <a href="Differential-Equations.html" accesskey="u" rel="up">Differential Equations</a> [<a href="index.html#SEC_Contents" title="Table of contents" rel="contents">Contents</a>][<a href="Concept-Index.html" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<h3 class="section" id="Differential_002dAlgebraic-Equations-1"><span>24.2 Differential-Algebraic Equations<a class="copiable-link" href="#Differential_002dAlgebraic-Equations-1"> ¶</a></span></h3>
<p>The function <code class="code">daspk</code> can be used to solve DAEs of the form
</p>
<div class="example">
<pre class="example-preformatted">0 = f (x-dot, x, t), x(t=0) = x_0, x-dot(t=0) = x-dot_0
</pre></div>
<p>where
<em class="math">x-dot</em>
is the derivative of <em class="math">x</em>. The equation is solved using
Petzold’s DAE solver <small class="sc">DASPK</small>.
</p>
<a class="anchor" id="XREFdaspk"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-daspk"><span><code class="def-type">[<var class="var">x</var>, <var class="var">xdot</var>, <var class="var">istate</var>, <var class="var">msg</var>] =</code> <strong class="def-name">daspk</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">x_0</var>, <var class="var">xdot_0</var>, <var class="var">t</var>, <var class="var">t_crit</var>)</code><a class="copiable-link" href="#index-daspk"> ¶</a></span></dt>
<dd><p>Solve a set of differential-algebraic equations.
</p>
<p><code class="code">daspk</code> solves the set of equations
</p>
<div class="example">
<pre class="example-preformatted">0 = f (x, xdot, t)
</pre></div>
<p>with
</p>
<div class="example">
<pre class="example-preformatted">x(t_0) = x_0, xdot(t_0) = xdot_0
</pre></div>
<p>The solution is returned in the matrices <var class="var">x</var> and <var class="var">xdot</var>,
with each row in the result matrices corresponding to one of the
elements in the vector <var class="var">t</var>. The first element of <var class="var">t</var>
should be <em class="math">t_0</em> and correspond to the initial state of the
system <var class="var">x_0</var> and its derivative <var class="var">xdot_0</var>, so that the first
row of the output <var class="var">x</var> is <var class="var">x_0</var> and the first row
of the output <var class="var">xdot</var> is <var class="var">xdot_0</var>.
</p>
<p>The first argument, <var class="var">fcn</var>, is a string, inline, or function handle
that names the function <em class="math">f</em> to call to compute the vector of
residuals for the set of equations. It must have the form
</p>
<div class="example">
<pre class="example-preformatted"><var class="var">res</var> = f (<var class="var">x</var>, <var class="var">xdot</var>, <var class="var">t</var>)
</pre></div>
<p>in which <var class="var">x</var>, <var class="var">xdot</var>, and <var class="var">res</var> are vectors, and <var class="var">t</var> is a
scalar.
</p>
<p>If <var class="var">fcn</var> 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 <em class="math">f</em> described above, and the second element names a
function to compute the modified Jacobian
</p>
<div class="example">
<div class="group"><pre class="example-preformatted"> df df
jac = -- + c ------
dx d xdot
</pre></div></div>
<p>The modified Jacobian function must have the form
</p>
<div class="example">
<div class="group"><pre class="example-preformatted">
<var class="var">jac</var> = j (<var class="var">x</var>, <var class="var">xdot</var>, <var class="var">t</var>, <var class="var">c</var>)
</pre></div></div>
<p>The second and third arguments to <code class="code">daspk</code> 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.
</p>
<p>The set of initial states and derivatives are not strictly required to
be consistent. If they are not consistent, you must use the
<code class="code">daspk_options</code> function to provide additional information so
that <code class="code">daspk</code> can compute a consistent starting point.
</p>
<p>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.
</p>
<p>After a successful computation, the value of <var class="var">istate</var> will be
greater than zero (consistent with the Fortran version of <small class="sc">DASPK</small>).
</p>
<p>If the computation is not successful, the value of <var class="var">istate</var> will be
less than zero and <var class="var">msg</var> will contain additional information.
</p>
<p>You can use the function <code class="code">daspk_options</code> to set optional
parameters for <code class="code">daspk</code>.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFdassl">dassl</a>.
</p></dd></dl>
<a class="anchor" id="XREFdaspk_005foptions"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-daspk_005foptions"><span><strong class="def-name">daspk_options</strong> <code class="def-code-arguments">()</code><a class="copiable-link" href="#index-daspk_005foptions"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-daspk_005foptions-1"><span><code class="def-type">val =</code> <strong class="def-name">daspk_options</strong> <code class="def-code-arguments">(<var class="var">opt</var>)</code><a class="copiable-link" href="#index-daspk_005foptions-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-daspk_005foptions-2"><span><strong class="def-name">daspk_options</strong> <code class="def-code-arguments">(<var class="var">opt</var>, <var class="var">val</var>)</code><a class="copiable-link" href="#index-daspk_005foptions-2"> ¶</a></span></dt>
<dd><p>Query or set options for the function <code class="code">daspk</code>.
</p>
<p>When called with no arguments, the names of all available options and
their current values are displayed.
</p>
<p>Given one argument, return the value of the option <var class="var">opt</var>.
</p>
<p>When called with two arguments, <code class="code">daspk_options</code> sets the option
<var class="var">opt</var> to value <var class="var">val</var>.
</p>
<p>Options include
</p>
<dl class="table">
<dt><code class="code">"absolute tolerance"</code></dt>
<dd><p>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.
</p>
</dd>
<dt><code class="code">"relative tolerance"</code></dt>
<dd><p>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.
</p>
<p>The local error test applied at each integration step is
</p>
<div class="example">
<div class="group"><pre class="example-preformatted"> abs (local error in x(i))
<= rtol(i) * abs (Y(i)) + atol(i)
</pre></div></div>
</dd>
<dt><code class="code">"compute consistent initial condition"</code></dt>
<dd><p>Denoting the differential variables in the state vector by ‘<samp class="samp">Y_d</samp>’
and the algebraic variables by ‘<samp class="samp">Y_a</samp>’, <code class="code">ddaspk</code> can solve
one of two initialization problems:
</p>
<ol class="enumerate">
<li> Given Y_d, calculate Y_a and Y’_d
</li><li> Given Y’, calculate Y.
</li></ol>
<p>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).
</p>
<p>If this option is set to a nonzero value, you must also set the
<code class="code">"algebraic variables"</code> option to declare which variables in the
problem are algebraic.
</p>
</dd>
<dt><code class="code">"use initial condition heuristics"</code></dt>
<dd><p>Set to a nonzero value to use the initial condition heuristics options
described below.
</p>
</dd>
<dt><code class="code">"initial condition heuristics"</code></dt>
<dd><p>A vector of the following parameters that can be used to control the
initial condition calculation.
</p>
<dl class="table">
<dt><code class="code">MXNIT</code></dt>
<dd><p>Maximum number of Newton iterations (default is 5).
</p>
</dd>
<dt><code class="code">MXNJ</code></dt>
<dd><p>Maximum number of Jacobian evaluations (default is 6).
</p>
</dd>
<dt><code class="code">MXNH</code></dt>
<dd><p>Maximum number of values of the artificial stepsize parameter to be
tried if the <code class="code">"compute consistent initial condition"</code> option has
been set to 1 (default is 5).
</p>
<p>Note that the maximum total number of Newton iterations allowed is
<code class="code">MXNIT*MXNJ*MXNH</code> if the <code class="code">"compute consistent initial
condition"</code> option has been set to 1 and <code class="code">MXNIT*MXNJ</code> if it is
set to 2.
</p>
</dd>
<dt><code class="code">LSOFF</code></dt>
<dd><p>Set to a nonzero value to disable the linesearch algorithm (default is
0).
</p>
</dd>
<dt><code class="code">STPTOL</code></dt>
<dd><p>Minimum scaled step in linesearch algorithm (default is eps^(2/3)).
</p>
</dd>
<dt><code class="code">EPINIT</code></dt>
<dd><p>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 class="code">EPINIT*EPCON</code>,
where <code class="code">EPCON</code> = 0.33 is the analogous test constant used in the
time steps. The default is <code class="code">EPINIT</code> = 0.01.
</p></dd>
</dl>
</dd>
<dt><code class="code">"print initial condition info"</code></dt>
<dd><p>Set this option to a nonzero value to display detailed information
about the initial condition calculation (default is 0).
</p>
</dd>
<dt><code class="code">"exclude algebraic variables from error test"</code></dt>
<dd><p>Set to a nonzero value to exclude algebraic variables from the error
test. You must also set the <code class="code">"algebraic variables"</code> option to
declare which variables in the problem are algebraic (default is 0).
</p>
</dd>
<dt><code class="code">"algebraic variables"</code></dt>
<dd><p>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).
</p>
<p>This option is required by the
<code class="code">"compute consistent initial condition"</code> and
<code class="code">"exclude algebraic variables from error test"</code> options.
</p>
</dd>
<dt><code class="code">"enforce inequality constraints"</code></dt>
<dd><p>Set to one of the following values to enforce the inequality
constraints specified by the <code class="code">"inequality constraint types"</code>
option (default is 0).
</p>
<ol class="enumerate">
<li> To have constraint checking only in the initial condition calculation.
</li><li> To enforce constraint checking during the integration.
</li><li> To enforce both options 1 and 2.
</li></ol>
</dd>
<dt><code class="code">"inequality constraint types"</code></dt>
<dd><p>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
</p>
<dl class="table">
<dt>-2</dt>
<dd><p>Less than zero.
</p>
</dd>
<dt>-1</dt>
<dd><p>Less than or equal to zero.
</p>
</dd>
<dt>0</dt>
<dd><p>Not constrained.
</p>
</dd>
<dt>1</dt>
<dd><p>Greater than or equal to zero.
</p>
</dd>
<dt>2</dt>
<dd><p>Greater than zero.
</p></dd>
</dl>
<p>This option only has an effect if the
<code class="code">"enforce inequality constraints"</code> option is nonzero.
</p>
</dd>
<dt><code class="code">"initial step size"</code></dt>
<dd><p>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).
</p>
</dd>
<dt><code class="code">"maximum order"</code></dt>
<dd><p>Restrict the maximum order of the solution method. This option must
be between 1 and 5, inclusive (default is 5).
</p>
</dd>
<dt><code class="code">"maximum step size"</code></dt>
<dd><p>Setting the maximum stepsize will avoid passing over very large
regions (default is not specified).
</p></dd>
</dl>
</dd></dl>
<p>Octave also includes <small class="sc">DASSL</small>, an earlier version of <small class="sc">DASPK</small>,
and <small class="sc">DASRT</small>, which can be used to solve DAEs with constraints
(stopping conditions).
</p>
<a class="anchor" id="XREFdassl"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-dassl"><span><code class="def-type">[<var class="var">x</var>, <var class="var">xdot</var>, <var class="var">istate</var>, <var class="var">msg</var>] =</code> <strong class="def-name">dassl</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">x_0</var>, <var class="var">xdot_0</var>, <var class="var">t</var>, <var class="var">t_crit</var>)</code><a class="copiable-link" href="#index-dassl"> ¶</a></span></dt>
<dd><p>Solve a set of differential-algebraic equations.
</p>
<p><code class="code">dassl</code> solves the set of equations
</p>
<div class="example">
<pre class="example-preformatted">0 = f (x, xdot, t)
</pre></div>
<p>with
</p>
<div class="example">
<pre class="example-preformatted">x(t_0) = x_0, xdot(t_0) = xdot_0
</pre></div>
<p>The solution is returned in the matrices <var class="var">x</var> and <var class="var">xdot</var>,
with each row in the result matrices corresponding to one of the
elements in the vector <var class="var">t</var>. The first element of <var class="var">t</var>
should be <em class="math">t_0</em> and correspond to the initial state of the
system <var class="var">x_0</var> and its derivative <var class="var">xdot_0</var>, so that the first
row of the output <var class="var">x</var> is <var class="var">x_0</var> and the first row
of the output <var class="var">xdot</var> is <var class="var">xdot_0</var>.
</p>
<p>The first argument, <var class="var">fcn</var>, is a string, inline, or function handle
that names the function <em class="math">f</em> to call to compute the vector of
residuals for the set of equations. It must have the form
</p>
<div class="example">
<pre class="example-preformatted"><var class="var">res</var> = f (<var class="var">x</var>, <var class="var">xdot</var>, <var class="var">t</var>)
</pre></div>
<p>in which <var class="var">x</var>, <var class="var">xdot</var>, and <var class="var">res</var> are vectors, and <var class="var">t</var> is a
scalar.
</p>
<p>If <var class="var">fcn</var> 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 <em class="math">f</em> described above, and the second element names a
function to compute the modified Jacobian
</p>
<div class="example">
<div class="group"><pre class="example-preformatted"> df df
jac = -- + c ------
dx d xdot
</pre></div></div>
<p>The modified Jacobian function must have the form
</p>
<div class="example">
<div class="group"><pre class="example-preformatted">
<var class="var">jac</var> = j (<var class="var">x</var>, <var class="var">xdot</var>, <var class="var">t</var>, <var class="var">c</var>)
</pre></div></div>
<p>The second and third arguments to <code class="code">dassl</code> 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.
</p>
<p>The set of initial states and derivatives are not strictly required to
be consistent. In practice, however, <small class="sc">DASSL</small> 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.
</p>
<p>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.
</p>
<p>After a successful computation, the value of <var class="var">istate</var> will be
greater than zero (consistent with the Fortran version of <small class="sc">DASSL</small>).
</p>
<p>If the computation is not successful, the value of <var class="var">istate</var> will be
less than zero and <var class="var">msg</var> will contain additional information.
</p>
<p>You can use the function <code class="code">dassl_options</code> to set optional
parameters for <code class="code">dassl</code>.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFdaspk">daspk</a>, <a class="ref" href="#XREFdasrt">dasrt</a>, <a class="ref" href="Ordinary-Differential-Equations.html#XREFlsode">lsode</a>.
</p></dd></dl>
<a class="anchor" id="XREFdassl_005foptions"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-dassl_005foptions"><span><strong class="def-name">dassl_options</strong> <code class="def-code-arguments">()</code><a class="copiable-link" href="#index-dassl_005foptions"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-dassl_005foptions-1"><span><code class="def-type">val =</code> <strong class="def-name">dassl_options</strong> <code class="def-code-arguments">(<var class="var">opt</var>)</code><a class="copiable-link" href="#index-dassl_005foptions-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-dassl_005foptions-2"><span><strong class="def-name">dassl_options</strong> <code class="def-code-arguments">(<var class="var">opt</var>, <var class="var">val</var>)</code><a class="copiable-link" href="#index-dassl_005foptions-2"> ¶</a></span></dt>
<dd><p>Query or set options for the function <code class="code">dassl</code>.
</p>
<p>When called with no arguments, the names of all available options and
their current values are displayed.
</p>
<p>Given one argument, return the value of the option <var class="var">opt</var>.
</p>
<p>When called with two arguments, <code class="code">dassl_options</code> sets the option
<var class="var">opt</var> to value <var class="var">val</var>.
</p>
<p>Options include
</p>
<dl class="table">
<dt><code class="code">"absolute tolerance"</code></dt>
<dd><p>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.
</p>
</dd>
<dt><code class="code">"relative tolerance"</code></dt>
<dd><p>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.
</p>
<p>The local error test applied at each integration step is
</p>
<div class="example">
<div class="group"><pre class="example-preformatted"> abs (local error in x(i))
<= rtol(i) * abs (Y(i)) + atol(i)
</pre></div></div>
</dd>
<dt><code class="code">"compute consistent initial condition"</code></dt>
<dd><p>If nonzero, <code class="code">dassl</code> 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.
</p>
</dd>
<dt><code class="code">"enforce nonnegativity constraints"</code></dt>
<dd><p>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.
</p>
</dd>
<dt><code class="code">"initial step size"</code></dt>
<dd><p>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.
</p>
</dd>
<dt><code class="code">"maximum order"</code></dt>
<dd><p>Restrict the maximum order of the solution method. This option must
be between 1 and 5, inclusive.
</p>
</dd>
<dt><code class="code">"maximum step size"</code></dt>
<dd><p>Setting the maximum stepsize will avoid passing over very large
regions (default is not specified).
</p>
</dd>
<dt><code class="code">"step limit"</code></dt>
<dd><p>Maximum number of integration steps to attempt on a single call to the
underlying Fortran code.
</p></dd>
</dl>
</dd></dl>
<a class="anchor" id="XREFdasrt"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-dasrt"><span><code class="def-type">[<var class="var">x</var>, <var class="var">xdot</var>, <var class="var">t_out</var>, <var class="var">istat</var>, <var class="var">msg</var>] =</code> <strong class="def-name">dasrt</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">g</var>, <var class="var">x_0</var>, <var class="var">xdot_0</var>, <var class="var">t</var>)</code><a class="copiable-link" href="#index-dasrt"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-dasrt-1"><span><code class="def-type">… =</code> <strong class="def-name">dasrt</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">g</var>, <var class="var">x_0</var>, <var class="var">xdot_0</var>, <var class="var">t</var>, <var class="var">t_crit</var>)</code><a class="copiable-link" href="#index-dasrt-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-dasrt-2"><span><code class="def-type">… =</code> <strong class="def-name">dasrt</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">x_0</var>, <var class="var">xdot_0</var>, <var class="var">t</var>)</code><a class="copiable-link" href="#index-dasrt-2"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-dasrt-3"><span><code class="def-type">… =</code> <strong class="def-name">dasrt</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">x_0</var>, <var class="var">xdot_0</var>, <var class="var">t</var>, <var class="var">t_crit</var>)</code><a class="copiable-link" href="#index-dasrt-3"> ¶</a></span></dt>
<dd><p>Solve a set of differential-algebraic equations.
</p>
<p><code class="code">dasrt</code> solves the set of equations
</p>
<div class="example">
<pre class="example-preformatted">0 = f (x, xdot, t)
</pre></div>
<p>with
</p>
<div class="example">
<pre class="example-preformatted">x(t_0) = x_0, xdot(t_0) = xdot_0
</pre></div>
<p>with functional stopping criteria (root solving).
</p>
<p>The solution is returned in the matrices <var class="var">x</var> and <var class="var">xdot</var>,
with each row in the result matrices corresponding to one of the
elements in the vector <var class="var">t_out</var>. The first element of <var class="var">t</var>
should be <em class="math">t_0</em> and correspond to the initial state of the
system <var class="var">x_0</var> and its derivative <var class="var">xdot_0</var>, so that the first
row of the output <var class="var">x</var> is <var class="var">x_0</var> and the first row
of the output <var class="var">xdot</var> is <var class="var">xdot_0</var>.
</p>
<p>The vector <var class="var">t</var> provides an upper limit on the length of the
integration. If the stopping condition is met, the vector
<var class="var">t_out</var> will be shorter than <var class="var">t</var>, and the final element of
<var class="var">t_out</var> will be the point at which the stopping condition was met,
and may not correspond to any element of the vector <var class="var">t</var>.
</p>
<p>The first argument, <var class="var">fcn</var>, is a string, inline, or function handle
that names the function <em class="math">f</em> to call to compute the vector of
residuals for the set of equations. It must have the form
</p>
<div class="example">
<pre class="example-preformatted"><var class="var">res</var> = f (<var class="var">x</var>, <var class="var">xdot</var>, <var class="var">t</var>)
</pre></div>
<p>in which <var class="var">x</var>, <var class="var">xdot</var>, and <var class="var">res</var> are vectors, and <var class="var">t</var> is a
scalar.
</p>
<p>If <var class="var">fcn</var> 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 <em class="math">f</em> described above, and the second element names a
function to compute the modified Jacobian
</p>
<div class="example">
<div class="group"><pre class="example-preformatted"> df df
jac = -- + c ------
dx d xdot
</pre></div></div>
<p>The modified Jacobian function must have the form
</p>
<div class="example">
<div class="group"><pre class="example-preformatted">
<var class="var">jac</var> = j (<var class="var">x</var>, <var class="var">xdot</var>, <var class="var">t</var>, <var class="var">c</var>)
</pre></div></div>
<p>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
</p>
<div class="example">
<pre class="example-preformatted"><var class="var">g_out</var> = g (<var class="var">x</var>, <var class="var">t</var>)
</pre></div>
<p>and return a vector of the constraint function values.
If the value of any of the constraint functions changes sign, <small class="sc">DASRT</small>
will attempt to stop the integration at the point of the sign change.
</p>
<p>If the name of the constraint function is omitted, <code class="code">dasrt</code> solves
the same problem as <code class="code">daspk</code> or <code class="code">dassl</code>.
</p>
<p>Note that because of numerical errors in the constraint functions
due to round-off and integration error, <small class="sc">DASRT</small> may return false
roots, or return the same root at two or more nearly equal values of
<var class="var">T</var>. If such false roots are suspected, the user should consider
smaller error tolerances or higher precision in the evaluation of the
constraint functions.
</p>
<p>If a root of some constraint function defines the end of the problem,
the input to <small class="sc">DASRT</small> should nevertheless allow integration to a
point slightly past that root, so that <small class="sc">DASRT</small> can locate the root
by interpolation.
</p>
<p>The third and fourth arguments to <code class="code">dasrt</code> 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.
</p>
<p>The set of initial states and derivatives are not strictly required to
be consistent. In practice, however, <small class="sc">DASSL</small> 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.
</p>
<p>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.
</p>
<p>After a successful computation, the value of <var class="var">istate</var> will be
greater than zero (consistent with the Fortran version of <small class="sc">DASSL</small>).
</p>
<p>If the computation is not successful, the value of <var class="var">istate</var> will be
less than zero and <var class="var">msg</var> will contain additional information.
</p>
<p>You can use the function <code class="code">dasrt_options</code> to set optional
parameters for <code class="code">dasrt</code>.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFdasrt_005foptions">dasrt_options</a>, <a class="ref" href="#XREFdaspk">daspk</a>, <a class="ref" href="#XREFdasrt">dasrt</a>, <a class="ref" href="Ordinary-Differential-Equations.html#XREFlsode">lsode</a>.
</p></dd></dl>
<a class="anchor" id="XREFdasrt_005foptions"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-dasrt_005foptions"><span><strong class="def-name">dasrt_options</strong> <code class="def-code-arguments">()</code><a class="copiable-link" href="#index-dasrt_005foptions"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-dasrt_005foptions-1"><span><code class="def-type">val =</code> <strong class="def-name">dasrt_options</strong> <code class="def-code-arguments">(<var class="var">opt</var>)</code><a class="copiable-link" href="#index-dasrt_005foptions-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-dasrt_005foptions-2"><span><strong class="def-name">dasrt_options</strong> <code class="def-code-arguments">(<var class="var">opt</var>, <var class="var">val</var>)</code><a class="copiable-link" href="#index-dasrt_005foptions-2"> ¶</a></span></dt>
<dd><p>Query or set options for the function <code class="code">dasrt</code>.
</p>
<p>When called with no arguments, the names of all available options and
their current values are displayed.
</p>
<p>Given one argument, return the value of the option <var class="var">opt</var>.
</p>
<p>When called with two arguments, <code class="code">dasrt_options</code> sets the option
<var class="var">opt</var> to value <var class="var">val</var>.
</p>
<p>Options include
</p>
<dl class="table">
<dt><code class="code">"absolute tolerance"</code></dt>
<dd><p>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.
</p>
</dd>
<dt><code class="code">"relative tolerance"</code></dt>
<dd><p>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.
</p>
<p>The local error test applied at each integration step is
</p>
<div class="example">
<div class="group"><pre class="example-preformatted"> abs (local error in x(i)) <= ...
rtol(i) * abs (Y(i)) + atol(i)
</pre></div></div>
</dd>
<dt><code class="code">"initial step size"</code></dt>
<dd><p>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.
</p>
</dd>
<dt><code class="code">"maximum order"</code></dt>
<dd><p>Restrict the maximum order of the solution method. This option must
be between 1 and 5, inclusive.
</p>
</dd>
<dt><code class="code">"maximum step size"</code></dt>
<dd><p>Setting the maximum stepsize will avoid passing over very large
regions.
</p>
</dd>
<dt><code class="code">"step limit"</code></dt>
<dd><p>Maximum number of integration steps to attempt on a single call to the
underlying Fortran code.
</p></dd>
</dl>
</dd></dl>
<p>See K. E. Brenan, et al., <cite class="cite">Numerical Solution of Initial-Value
Problems in Differential-Algebraic Equations</cite>, North-Holland (1989),
DOI: <a class="url" href="https://doi.org/10.1137/1.9781611971224">https://doi.org/10.1137/1.9781611971224</a>,
for more information about the implementation of <small class="sc">DASSL</small>.
</p>
</div>
<hr>
<div class="nav-panel">
<p>
Next: <a href="Matlab_002dcompatible-solvers.html">Matlab-compatible solvers</a>, Previous: <a href="Ordinary-Differential-Equations.html">Ordinary Differential Equations</a>, Up: <a href="Differential-Equations.html">Differential Equations</a> [<a href="index.html#SEC_Contents" title="Table of contents" rel="contents">Contents</a>][<a href="Concept-Index.html" title="Index" rel="index">Index</a>]</p>
</div>
</body>
</html>
|