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
|
<!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>Matlab-compatible solvers (GNU Octave (version 10.3.0))</title>
<meta name="description" content="Matlab-compatible solvers (GNU Octave (version 10.3.0))">
<meta name="keywords" content="Matlab-compatible solvers (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="Differential_002dAlgebraic-Equations.html" rel="prev" title="Differential-Algebraic 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}
ul.mark-bullet {list-style-type: disc}
-->
</style>
<link rel="stylesheet" type="text/css" href="octave.css">
</head>
<body lang="en">
<div class="section-level-extent" id="Matlab_002dcompatible-solvers">
<div class="nav-panel">
<p>
Previous: <a href="Differential_002dAlgebraic-Equations.html" accesskey="p" rel="prev">Differential-Algebraic 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="Matlab_002dcompatible-solvers-1"><span>24.3 Matlab-compatible solvers<a class="copiable-link" href="#Matlab_002dcompatible-solvers-1"> ¶</a></span></h3>
<p>Octave also provides a set of solvers for initial value problems for ordinary
differential equations (ODEs) that have a <small class="sc">MATLAB</small>-compatible interface.
The options for this class of methods are set using the functions.
</p>
<ul class="itemize mark-bullet">
<li><a class="ref" href="#XREFodeset">odeset</a>
</li><li><a class="ref" href="#XREFodeget">odeget</a>
</li></ul>
<p>Currently implemented solvers are:
</p>
<ul class="itemize mark-bullet">
<li>Runge-Kutta methods
<ul class="itemize mark-bullet">
<li><a class="ref" href="#XREFode45">ode45</a> integrates a system of non-stiff ODEs or
index-1 differential-algebraic equations (DAEs) using the high-order,
variable-step Dormand-Prince method. It requires six function
evaluations per integration step, but may take larger steps on smooth
problems than <code class="code">ode23</code>: potentially offering improved efficiency at
smaller tolerances.
</li><li><a class="ref" href="#XREFode23">ode23</a> integrates a system of non-stiff ODEs or (or
index-1 DAEs). It uses the third-order Bogacki-Shampine method
and adapts the local step size in order to satisfy a user-specified
tolerance. The solver requires three function evaluations per integration
step.
</li><li><a class="ref" href="#XREFode23s">ode23s</a> integrates a system of stiff ODEs (or
index-1 DAEs) using a modified second-order Rosenbrock method.
</li></ul>
</li><li>Linear multistep methods
<ul class="itemize mark-bullet">
<li><a class="ref" href="#XREFode15s">ode15s</a> integrates a system of stiff ODEs (or
index-1 DAEs) using a variable step, variable order method based on
Backward Difference Formulas (BDF).
</li><li><a class="ref" href="#XREFode15i">ode15i</a> integrates a system of fully-implicit ODEs
(or index-1 DAEs) using the same variable step, variable order method as
<code class="code">ode15s</code>. <a class="ref" href="#XREFdecic">decic</a> can be used to compute consistent
initial conditions for <code class="code">ode15i</code>.
</li></ul>
</li></ul>
<p>Detailed information on the solvers are given in L. F. Shampine and
M. W. Reichelt, <cite class="cite">The MATLAB ODE Suite</cite>, SIAM Journal on
Scientific Computing, Vol. 18, 1997, pp. 1–22,
DOI: <a class="url" href="https://doi.org/10.1137/S1064827594276424">https://doi.org/10.1137/S1064827594276424</a>.
</p>
<a class="anchor" id="XREFode45"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-ode45"><span><code class="def-type">[<var class="var">t</var>, <var class="var">y</var>] =</code> <strong class="def-name">ode45</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">trange</var>, <var class="var">init</var>)</code><a class="copiable-link" href="#index-ode45"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode45-1"><span><code class="def-type">[<var class="var">t</var>, <var class="var">y</var>] =</code> <strong class="def-name">ode45</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">trange</var>, <var class="var">init</var>, <var class="var">ode_opt</var>)</code><a class="copiable-link" href="#index-ode45-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode45-2"><span><code class="def-type">[<var class="var">t</var>, <var class="var">y</var>, <var class="var">te</var>, <var class="var">ye</var>, <var class="var">ie</var>] =</code> <strong class="def-name">ode45</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-ode45-2"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode45-3"><span><code class="def-type"><var class="var">solution</var> =</code> <strong class="def-name">ode45</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-ode45-3"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode45-4"><span><strong class="def-name">ode45</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-ode45-4"> ¶</a></span></dt>
<dd>
<p>Solve a set of non-stiff Ordinary Differential Equations (non-stiff ODEs)
with the well known explicit Dormand-Prince method of order 4.
</p>
<p><var class="var">fcn</var> is a function handle, inline function, or string containing the
name of the function that defines the ODE: <code class="code">y' = f(t,y)</code>. The function
must accept two inputs where the first is time <var class="var">t</var> and the second is a
column vector of unknowns <var class="var">y</var>.
</p>
<p><var class="var">trange</var> specifies the time interval over which the ODE will be
evaluated. Typically, it is a two-element vector specifying the initial and
final times (<code class="code">[tinit, tfinal]</code>). If there are more than two elements
then the solution will also be evaluated at these intermediate time
instances.
</p>
<p>By default, <code class="code">ode45</code> uses an adaptive timestep with the
<code class="code">integrate_adaptive</code> algorithm. The tolerance for the timestep
computation may be changed by using the options <code class="code">"RelTol"</code> and
<code class="code">"AbsTol"</code>.
</p>
<p><var class="var">init</var> contains the initial value for the unknowns. If it is a row
vector then the solution <var class="var">y</var> will be a matrix in which each column is
the solution for the corresponding initial value in <var class="var">init</var>.
</p>
<p>The optional fourth argument <var class="var">ode_opt</var> specifies non-default options to
the ODE solver. It is a structure generated by <code class="code">odeset</code>.
</p>
<p>The function typically returns two outputs. Variable <var class="var">t</var> is a
column vector and contains the times where the solution was found. The
output <var class="var">y</var> is a matrix in which each column refers to a different
unknown of the problem and each row corresponds to a time in <var class="var">t</var>.
</p>
<p>The output can also be returned as a structure <var class="var">solution</var> which has a
field <var class="var">x</var> containing a row vector of times where the solution was
evaluated and a field <var class="var">y</var> containing the solution matrix such that each
column corresponds to a time in <var class="var">x</var>. Use
<code class="code">fieldnames (<var class="var">solution</var>)</code><!-- /@w --> to see the other fields and
additional information returned.
</p>
<p>If no output arguments are requested, and no <code class="code">"OutputFcn"</code> is
specified in <var class="var">ode_opt</var>, then the <code class="code">"OutputFcn"</code> is set to
<code class="code">odeplot</code> and the results of the solver are plotted immediately.
</p>
<p>If using the <code class="code">"Events"</code> option then three additional outputs may be
returned. <var class="var">te</var> holds the time when an Event function returned a zero.
<var class="var">ye</var> holds the value of the solution at time <var class="var">te</var>. <var class="var">ie</var>
contains an index indicating which Event function was triggered in the case
of multiple Event functions.
</p>
<p>Example: Solve the Van der Pol equation
</p>
<div class="example">
<div class="group"><pre class="example-preformatted">fvdp = @(<var class="var">t</var>,<var class="var">y</var>) [<var class="var">y</var>(2); (1 - <var class="var">y</var>(1)^2) * <var class="var">y</var>(2) - <var class="var">y</var>(1)];
[<var class="var">t</var>,<var class="var">y</var>] = ode45 (fvdp, [0, 20], [2, 0]);
</pre></div></div>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFodeset">odeset</a>, <a class="ref" href="#XREFodeget">odeget</a>, <a class="ref" href="#XREFode23">ode23</a>, <a class="ref" href="#XREFode15s">ode15s</a>.
</p></dd></dl>
<a class="anchor" id="XREFode23"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-ode23"><span><code class="def-type">[<var class="var">t</var>, <var class="var">y</var>] =</code> <strong class="def-name">ode23</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">trange</var>, <var class="var">init</var>)</code><a class="copiable-link" href="#index-ode23"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode23-1"><span><code class="def-type">[<var class="var">t</var>, <var class="var">y</var>] =</code> <strong class="def-name">ode23</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">trange</var>, <var class="var">init</var>, <var class="var">ode_opt</var>)</code><a class="copiable-link" href="#index-ode23-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode23-2"><span><code class="def-type">[<var class="var">t</var>, <var class="var">y</var>, <var class="var">te</var>, <var class="var">ye</var>, <var class="var">ie</var>] =</code> <strong class="def-name">ode23</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-ode23-2"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode23-3"><span><code class="def-type"><var class="var">solution</var> =</code> <strong class="def-name">ode23</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-ode23-3"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode23-4"><span><strong class="def-name">ode23</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-ode23-4"> ¶</a></span></dt>
<dd>
<p>Solve a set of non-stiff Ordinary Differential Equations (non-stiff ODEs)
with the well known explicit Bogacki-Shampine method of order 3.
</p>
<p><var class="var">fcn</var> is a function handle, inline function, or string containing the
name of the function that defines the ODE: <code class="code">y' = f(t,y)</code>. The function
must accept two inputs where the first is time <var class="var">t</var> and the second is a
column vector of unknowns <var class="var">y</var>.
</p>
<p><var class="var">trange</var> specifies the time interval over which the ODE will be
evaluated. Typically, it is a two-element vector specifying the initial and
final times (<code class="code">[tinit, tfinal]</code>). If there are more than two elements
then the solution will also be evaluated at these intermediate time
instances.
</p>
<p>By default, <code class="code">ode23</code> uses an adaptive timestep with the
<code class="code">integrate_adaptive</code> algorithm. The tolerance for the timestep
computation may be changed by using the options <code class="code">"RelTol"</code> and
<code class="code">"AbsTol"</code>.
</p>
<p><var class="var">init</var> contains the initial value for the unknowns. If it is a row
vector then the solution <var class="var">y</var> will be a matrix in which each column is
the solution for the corresponding initial value in <var class="var">init</var>.
</p>
<p>The optional fourth argument <var class="var">ode_opt</var> specifies non-default options to
the ODE solver. It is a structure generated by <code class="code">odeset</code>.
</p>
<p>The function typically returns two outputs. Variable <var class="var">t</var> is a
column vector and contains the times where the solution was found. The
output <var class="var">y</var> is a matrix in which each column refers to a different
unknown of the problem and each row corresponds to a time in <var class="var">t</var>.
</p>
<p>The output can also be returned as a structure <var class="var">solution</var> which has a
field <var class="var">x</var> containing a row vector of times where the solution was
evaluated and a field <var class="var">y</var> containing the solution matrix such that each
column corresponds to a time in <var class="var">x</var>. Use
<code class="code">fieldnames (<var class="var">solution</var>)</code><!-- /@w --> to see the other fields and
additional information returned.
</p>
<p>If no output arguments are requested, and no <code class="code">"OutputFcn"</code> is
specified in <var class="var">ode_opt</var>, then the <code class="code">"OutputFcn"</code> is set to
<code class="code">odeplot</code> and the results of the solver are plotted immediately.
</p>
<p>If using the <code class="code">"Events"</code> option then three additional outputs may be
returned. <var class="var">te</var> holds the time when an Event function returned a zero.
<var class="var">ye</var> holds the value of the solution at time <var class="var">te</var>. <var class="var">ie</var>
contains an index indicating which Event function was triggered in the case
of multiple Event functions.
</p>
<p>Example: Solve the Van der Pol equation
</p>
<div class="example">
<div class="group"><pre class="example-preformatted">fvdp = @(<var class="var">t</var>,<var class="var">y</var>) [<var class="var">y</var>(2); (1 - <var class="var">y</var>(1)^2) * <var class="var">y</var>(2) - <var class="var">y</var>(1)];
[<var class="var">t</var>,<var class="var">y</var>] = ode23 (fvdp, [0, 20], [2, 0]);
</pre></div></div>
<p>Reference: For the definition of this method see
<a class="url" href="https://en.wikipedia.org/wiki/List_of_Runge%E2%80%93Kutta_methods">https://en.wikipedia.org/wiki/List_of_Runge%E2%80%93Kutta_methods</a>.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFodeset">odeset</a>, <a class="ref" href="#XREFodeget">odeget</a>, <a class="ref" href="#XREFode45">ode45</a>, <a class="ref" href="#XREFode15s">ode15s</a>.
</p></dd></dl>
<a class="anchor" id="XREFode23s"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-ode23s"><span><code class="def-type">[<var class="var">t</var>, <var class="var">y</var>] =</code> <strong class="def-name">ode23s</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">trange</var>, <var class="var">init</var>)</code><a class="copiable-link" href="#index-ode23s"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode23s-1"><span><code class="def-type">[<var class="var">t</var>, <var class="var">y</var>] =</code> <strong class="def-name">ode23s</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">trange</var>, <var class="var">init</var>, <var class="var">ode_opt</var>)</code><a class="copiable-link" href="#index-ode23s-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode23s-2"><span><code class="def-type">[<var class="var">t</var>, <var class="var">y</var>] =</code> <strong class="def-name">ode23s</strong> <code class="def-code-arguments">(…, <var class="var">par1</var>, <var class="var">par2</var>, …)</code><a class="copiable-link" href="#index-ode23s-2"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode23s-3"><span><code class="def-type">[<var class="var">t</var>, <var class="var">y</var>, <var class="var">te</var>, <var class="var">ye</var>, <var class="var">ie</var>] =</code> <strong class="def-name">ode23s</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-ode23s-3"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode23s-4"><span><code class="def-type"><var class="var">solution</var> =</code> <strong class="def-name">ode23s</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-ode23s-4"> ¶</a></span></dt>
<dd>
<p>Solve a set of stiff Ordinary Differential Equations (stiff ODEs) with a
Rosenbrock method of order (2,3).
</p>
<p><var class="var">fcn</var> is a function handle, inline function, or string containing the
name of the function that defines the ODE: <code class="code">M y' = f(t,y)</code>. The
function must accept two inputs where the first is time <var class="var">t</var> and the
second is a column vector of unknowns <var class="var">y</var>. <var class="var">M</var> is a constant mass
matrix, non-singular and possibly sparse. Set the field <code class="code">"Mass"</code> in
<var class="var">odeopts</var> using <var class="var">odeset</var> to specify a mass matrix.
</p>
<p><var class="var">trange</var> specifies the time interval over which the ODE will be
evaluated. Typically, it is a two-element vector specifying the initial
and final times (<code class="code">[tinit, tfinal]</code>). If there are more than two
elements then the solution will also be evaluated at these intermediate
time instances using an interpolation procedure of the same order as the
one of the solver.
</p>
<p>By default, <code class="code">ode23s</code> uses an adaptive timestep with the
<code class="code">integrate_adaptive</code> algorithm. The tolerance for the timestep
computation may be changed by using the options <code class="code">"RelTol"</code> and
<code class="code">"AbsTol"</code>.
</p>
<p><var class="var">init</var> contains the initial value for the unknowns. If it is a row
vector then the solution <var class="var">y</var> will be a matrix in which each column is
the solution for the corresponding initial value in <var class="var">init</var>.
</p>
<p>The optional fourth argument <var class="var">ode_opt</var> specifies non-default options to
the ODE solver. It is a structure generated by <code class="code">odeset</code>.
<code class="code">ode23s</code> will ignore the following options: <code class="code">"BDF"</code>,
<code class="code">"InitialSlope"</code>, <code class="code">"MassSingular"</code>, <code class="code">"MStateDependence"</code>,
<code class="code">"MvPattern"</code>, <code class="code">"MaxOrder"</code>, <code class="code">"Non-negative"</code>.
</p>
<p>The function typically returns two outputs. Variable <var class="var">t</var> is a
column vector and contains the times where the solution was found. The
output <var class="var">y</var> is a matrix in which each column refers to a different
unknown of the problem and each row corresponds to a time in <var class="var">t</var>. If
<var class="var">trange</var> specifies intermediate time steps, only those will be returned.
</p>
<p>The output can also be returned as a structure <var class="var">solution</var> which has a
field <var class="var">x</var> containing a row vector of times where the solution was
evaluated and a field <var class="var">y</var> containing the solution matrix such that each
column corresponds to a time in <var class="var">x</var>. Use
<code class="code">fieldnames (<var class="var">solution</var>)</code><!-- /@w --> to see the other fields and
additional information returned.
</p>
<p>If using the <code class="code">"Events"</code> option then three additional outputs may be
returned. <var class="var">te</var> holds the time when an Event function returned a zero.
<var class="var">ye</var> holds the value of the solution at time <var class="var">te</var>. <var class="var">ie</var>
contains an index indicating which Event function was triggered in the case
of multiple Event functions.
</p>
<p>Example: Solve the stiff Van der Pol equation
</p>
<div class="example">
<div class="group"><pre class="example-preformatted">f = @(<var class="var">t</var>,<var class="var">y</var>) [<var class="var">y</var>(2); 1000*(1 - <var class="var">y</var>(1)^2) * <var class="var">y</var>(2) - <var class="var">y</var>(1)];
opt = odeset ('Mass', [1 0; 0 1], 'MaxStep', 1e-1);
[vt, vy] = ode23s (f, [0 2000], [2 0], opt);
</pre></div></div>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFodeset">odeset</a>, <a class="ref" href="Differential_002dAlgebraic-Equations.html#XREFdaspk">daspk</a>, <a class="ref" href="Differential_002dAlgebraic-Equations.html#XREFdassl">dassl</a>.
</p></dd></dl>
<a class="anchor" id="XREFode15s"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-ode15s"><span><code class="def-type">[<var class="var">t</var>, <var class="var">y</var>] =</code> <strong class="def-name">ode15s</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">trange</var>, <var class="var">y0</var>)</code><a class="copiable-link" href="#index-ode15s"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode15s-1"><span><code class="def-type">[<var class="var">t</var>, <var class="var">y</var>] =</code> <strong class="def-name">ode15s</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">trange</var>, <var class="var">y0</var>, <var class="var">ode_opt</var>)</code><a class="copiable-link" href="#index-ode15s-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode15s-2"><span><code class="def-type">[<var class="var">t</var>, <var class="var">y</var>, <var class="var">te</var>, <var class="var">ye</var>, <var class="var">ie</var>] =</code> <strong class="def-name">ode15s</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-ode15s-2"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode15s-3"><span><code class="def-type"><var class="var">solution</var> =</code> <strong class="def-name">ode15s</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-ode15s-3"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode15s-4"><span><strong class="def-name">ode15s</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-ode15s-4"> ¶</a></span></dt>
<dd><p>Solve a set of stiff Ordinary Differential Equations (ODEs) or stiff
semi-explicit index 1 Differential Algebraic Equations (DAEs).
</p>
<p><code class="code">ode15s</code> uses a variable step, variable order BDF (Backward
Differentiation Formula) method that ranges from order 1 to 5.
</p>
<p><var class="var">fcn</var> is a function handle, inline function, or string containing the
name of the function that defines the ODE: <code class="code">y' = f(t,y)</code>. The function
must accept two inputs where the first is time <var class="var">t</var> and the second is a
column vector of unknowns <var class="var">y</var>.
</p>
<p><var class="var">trange</var> specifies the time interval over which the ODE will be
evaluated. Typically, it is a two-element vector specifying the initial and
final times (<code class="code">[tinit, tfinal]</code>). If there are more than two elements
then the solution will also be evaluated at these intermediate time
instances.
</p>
<p><var class="var">init</var> contains the initial value for the unknowns. If it is a row
vector then the solution <var class="var">y</var> will be a matrix in which each column is
the solution for the corresponding initial value in <var class="var">init</var>.
</p>
<p>The optional fourth argument <var class="var">ode_opt</var> specifies non-default options to
the ODE solver. It is a structure generated by <code class="code">odeset</code>.
</p>
<p>The function typically returns two outputs. Variable <var class="var">t</var> is a
column vector and contains the times where the solution was found. The
output <var class="var">y</var> is a matrix in which each column refers to a different
unknown of the problem and each row corresponds to a time in <var class="var">t</var>.
</p>
<p>The output can also be returned as a structure <var class="var">solution</var> which has a
field <var class="var">x</var> containing a row vector of times where the solution was
evaluated and a field <var class="var">y</var> containing the solution matrix such that each
column corresponds to a time in <var class="var">x</var>. Use
<code class="code">fieldnames (<var class="var">solution</var>)</code><!-- /@w --> to see the other fields and
additional information returned.
</p>
<p>If no output arguments are requested, and no <code class="code">"OutputFcn"</code> is
specified in <var class="var">ode_opt</var>, then the <code class="code">"OutputFcn"</code> is set to
<code class="code">odeplot</code> and the results of the solver are plotted immediately.
</p>
<p>If using the <code class="code">"Events"</code> option then three additional outputs may be
returned. <var class="var">te</var> holds the time when an Event function returned a zero.
<var class="var">ye</var> holds the value of the solution at time <var class="var">te</var>. <var class="var">ie</var>
contains an index indicating which Event function was triggered in the case
of multiple Event functions.
</p>
<p>Example: Solve Robertson’s equations:
</p>
<div class="example smallexample">
<div class="group"><pre class="example-preformatted">function r = robertson_dae (<var class="var">t</var>, <var class="var">y</var>)
r = [ -0.04*<var class="var">y</var>(1) + 1e4*<var class="var">y</var>(2)*<var class="var">y</var>(3)
0.04*<var class="var">y</var>(1) - 1e4*<var class="var">y</var>(2)*<var class="var">y</var>(3) - 3e7*<var class="var">y</var>(2)^2
<var class="var">y</var>(1) + <var class="var">y</var>(2) + <var class="var">y</var>(3) - 1 ];
endfunction
opt = odeset ("Mass", [1 0 0; 0 1 0; 0 0 0], "MStateDependence", "none");
[<var class="var">t</var>,<var class="var">y</var>] = ode15s (@robertson_dae, [0, 1e3], [1; 0; 0], opt);
</pre></div></div>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFdecic">decic</a>, <a class="ref" href="#XREFodeset">odeset</a>, <a class="ref" href="#XREFodeget">odeget</a>, <a class="ref" href="#XREFode23">ode23</a>, <a class="ref" href="#XREFode45">ode45</a>.
</p></dd></dl>
<a class="anchor" id="XREFode15i"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-ode15i"><span><code class="def-type">[<var class="var">t</var>, <var class="var">y</var>] =</code> <strong class="def-name">ode15i</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">trange</var>, <var class="var">y0</var>, <var class="var">yp0</var>)</code><a class="copiable-link" href="#index-ode15i"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode15i-1"><span><code class="def-type">[<var class="var">t</var>, <var class="var">y</var>] =</code> <strong class="def-name">ode15i</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">trange</var>, <var class="var">y0</var>, <var class="var">yp0</var>, <var class="var">ode_opt</var>)</code><a class="copiable-link" href="#index-ode15i-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode15i-2"><span><code class="def-type">[<var class="var">t</var>, <var class="var">y</var>, <var class="var">te</var>, <var class="var">ye</var>, <var class="var">ie</var>] =</code> <strong class="def-name">ode15i</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-ode15i-2"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode15i-3"><span><code class="def-type"><var class="var">solution</var> =</code> <strong class="def-name">ode15i</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-ode15i-3"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-ode15i-4"><span><strong class="def-name">ode15i</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-ode15i-4"> ¶</a></span></dt>
<dd><p>Solve a set of fully-implicit Ordinary Differential Equations (ODEs) or
index 1 Differential Algebraic Equations (DAEs).
</p>
<p><code class="code">ode15i</code> uses a variable step, variable order BDF (Backward
Differentiation Formula) method that ranges from order 1 to 5.
</p>
<p><var class="var">fcn</var> is a function handle, inline function, or string containing the
name of the function that defines the ODE: <code class="code">0 = f(t,y,yp)</code>. The
function must accept three inputs where the first is time <var class="var">t</var>, the
second is the function value <var class="var">y</var> (a column vector), and the third
is the derivative value <var class="var">yp</var> (a column vector).
</p>
<p><var class="var">trange</var> specifies the time interval over which the ODE will be
evaluated. Typically, it is a two-element vector specifying the initial and
final times (<code class="code">[tinit, tfinal]</code>). If there are more than two elements
then the solution will also be evaluated at these intermediate time
instances.
</p>
<p><var class="var">y0</var> and <var class="var">yp0</var> contain the initial values for the unknowns <var class="var">y</var>
and <var class="var">yp</var>. If they are row vectors then the solution <var class="var">y</var> will be a
matrix in which each column is the solution for the corresponding initial
value in <var class="var">y0</var> and <var class="var">yp0</var>.
</p>
<p><var class="var">y0</var> and <var class="var">yp0</var> must be consistent initial conditions, meaning that
<code class="code">f(t,y0,yp0) = 0</code> is satisfied. The function <code class="code">decic</code> may be used
to compute consistent initial conditions given initial guesses.
</p>
<p>The optional fifth argument <var class="var">ode_opt</var> specifies non-default options to
the ODE solver. It is a structure generated by <code class="code">odeset</code>.
</p>
<p>The function typically returns two outputs. Variable <var class="var">t</var> is a
column vector and contains the times where the solution was found. The
output <var class="var">y</var> is a matrix in which each column refers to a different
unknown of the problem and each row corresponds to a time in <var class="var">t</var>.
</p>
<p>The output can also be returned as a structure <var class="var">solution</var> which has a
field <var class="var">x</var> containing a row vector of times where the solution was
evaluated and a field <var class="var">y</var> containing the solution matrix such that each
column corresponds to a time in <var class="var">x</var>. Use
<code class="code">fieldnames (<var class="var">solution</var>)</code><!-- /@w --> to see the other fields and
additional information returned.
</p>
<p>If no output arguments are requested, and no <code class="code">"OutputFcn"</code> is
specified in <var class="var">ode_opt</var>, then the <code class="code">"OutputFcn"</code> is set to
<code class="code">odeplot</code> and the results of the solver are plotted immediately.
</p>
<p>If using the <code class="code">"Events"</code> option then three additional outputs may be
returned. <var class="var">te</var> holds the time when an Event function returned a zero.
<var class="var">ye</var> holds the value of the solution at time <var class="var">te</var>. <var class="var">ie</var>
contains an index indicating which Event function was triggered in the case
of multiple Event functions.
</p>
<p>Example: Solve Robertson’s equations:
</p>
<div class="example smallexample">
<div class="group"><pre class="example-preformatted">function r = robertson_dae (<var class="var">t</var>, <var class="var">y</var>, <var class="var">yp</var>)
r = [ -(<var class="var">yp</var>(1) + 0.04*<var class="var">y</var>(1) - 1e4*<var class="var">y</var>(2)*<var class="var">y</var>(3))
-(<var class="var">yp</var>(2) - 0.04*<var class="var">y</var>(1) + 1e4*<var class="var">y</var>(2)*<var class="var">y</var>(3) + 3e7*<var class="var">y</var>(2)^2)
<var class="var">y</var>(1) + <var class="var">y</var>(2) + <var class="var">y</var>(3) - 1 ];
endfunction
[<var class="var">t</var>,<var class="var">y</var>] = ode15i (@robertson_dae, [0, 1e3], [1; 0; 0], [-1e-4; 1e-4; 0]);
</pre></div></div>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFdecic">decic</a>, <a class="ref" href="#XREFodeset">odeset</a>, <a class="ref" href="#XREFodeget">odeget</a>.
</p></dd></dl>
<a class="anchor" id="XREFdecic"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-decic"><span><code class="def-type">[<var class="var">y0_new</var>, <var class="var">yp0_new</var>] =</code> <strong class="def-name">decic</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">t0</var>, <var class="var">y0</var>, <var class="var">fixed_y0</var>, <var class="var">yp0</var>, <var class="var">fixed_yp0</var>)</code><a class="copiable-link" href="#index-decic"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-decic-1"><span><code class="def-type">[<var class="var">y0_new</var>, <var class="var">yp0_new</var>] =</code> <strong class="def-name">decic</strong> <code class="def-code-arguments">(<var class="var">fcn</var>, <var class="var">t0</var>, <var class="var">y0</var>, <var class="var">fixed_y0</var>, <var class="var">yp0</var>, <var class="var">fixed_yp0</var>, <var class="var">options</var>)</code><a class="copiable-link" href="#index-decic-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-decic-2"><span><code class="def-type">[<var class="var">y0_new</var>, <var class="var">yp0_new</var>, <var class="var">resnorm</var>] =</code> <strong class="def-name">decic</strong> <code class="def-code-arguments">(…)</code><a class="copiable-link" href="#index-decic-2"> ¶</a></span></dt>
<dd>
<p>Compute consistent implicit ODE initial conditions <var class="var">y0_new</var> and
<var class="var">yp0_new</var> given initial guesses <var class="var">y0</var> and <var class="var">yp0</var>.
</p>
<p>A maximum of <code class="code">length (<var class="var">y0</var>)</code> components between <var class="var">fixed_y0</var> and
<var class="var">fixed_yp0</var> may be chosen as fixed values.
</p>
<p><var class="var">fcn</var> is a function handle. The function must accept three inputs where
the first is time <var class="var">t</var>, the second is a column vector of unknowns
<var class="var">y</var>, and the third is a column vector of unknowns <var class="var">yp</var>.
</p>
<p><var class="var">t0</var> is the initial time such that
<code class="code"><var class="var">fcn</var>(<var class="var">t0</var>, <var class="var">y0_new</var>, <var class="var">yp0_new</var>) = 0</code>, specified as a
scalar.
</p>
<p><var class="var">y0</var> is a vector used as the initial guess for <var class="var">y</var>.
</p>
<p><var class="var">fixed_y0</var> is a vector which specifies the components of <var class="var">y0</var> to
hold fixed. Choose a maximum of <code class="code">length (<var class="var">y0</var>)</code> components between
<var class="var">fixed_y0</var> and <var class="var">fixed_yp0</var> as fixed values.
Set <var class="var">fixed_y0</var>(i) component to 1 if you want to fix the value of
<var class="var">y0</var>(i).
Set <var class="var">fixed_y0</var>(i) component to 0 if you want to allow the value of
<var class="var">y0</var>(i) to change.
</p>
<p><var class="var">yp0</var> is a vector used as the initial guess for <var class="var">yp</var>.
</p>
<p><var class="var">fixed_yp0</var> is a vector which specifies the components of <var class="var">yp0</var> to
hold fixed. Choose a maximum of <code class="code">length (<var class="var">yp0</var>)</code> components
between <var class="var">fixed_y0</var> and <var class="var">fixed_yp0</var> as fixed values.
Set <var class="var">fixed_yp0</var>(i) component to 1 if you want to fix the value of
<var class="var">yp0</var>(i).
Set <var class="var">fixed_yp0</var>(i) component to 0 if you want to allow the value of
<var class="var">yp0</var>(i) to change.
</p>
<p>The optional seventh argument <var class="var">options</var> is a structure array. Use
<code class="code">odeset</code> to generate this structure. The relevant options are
<code class="code">RelTol</code> and <code class="code">AbsTol</code> which specify the error thresholds used to
compute the initial conditions.
</p>
<p>The function typically returns two outputs. Variable <var class="var">y0_new</var> is a
column vector and contains the consistent initial value of <var class="var">y</var>. The
output <var class="var">yp0_new</var> is a column vector and contains the consistent initial
value of <var class="var">yp</var>.
</p>
<p>The optional third output <var class="var">resnorm</var> is the norm of the vector of
residuals. If <var class="var">resnorm</var> is small, <code class="code">decic</code> has successfully
computed the initial conditions. If the value of <var class="var">resnorm</var> is large,
use <code class="code">RelTol</code> and <code class="code">AbsTol</code> to adjust it.
</p>
<p>Example: Compute initial conditions for Robertson’s equations:
</p>
<div class="example smallexample">
<div class="group"><pre class="example-preformatted">function r = robertson_dae (<var class="var">t</var>, <var class="var">y</var>, <var class="var">yp</var>)
r = [ -(<var class="var">yp</var>(1) + 0.04*<var class="var">y</var>(1) - 1e4*<var class="var">y</var>(2)*<var class="var">y</var>(3))
-(<var class="var">yp</var>(2) - 0.04*<var class="var">y</var>(1) + 1e4*<var class="var">y</var>(2)*<var class="var">y</var>(3) + 3e7*<var class="var">y</var>(2)^2)
<var class="var">y</var>(1) + <var class="var">y</var>(2) + <var class="var">y</var>(3) - 1 ];
endfunction
</pre></div><pre class="example-preformatted">[<var class="var">y0_new</var>,<var class="var">yp0_new</var>] = decic (@robertson_dae, 0, [1; 0; 0], [1; 1; 0],
[-1e-4; 1; 0], [0; 0; 0]);
</pre></div>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFode15i">ode15i</a>, <a class="ref" href="#XREFodeset">odeset</a>.
</p></dd></dl>
<a class="anchor" id="XREFodeset"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-odeset"><span><code class="def-type"><var class="var">odestruct</var> =</code> <strong class="def-name">odeset</strong> <code class="def-code-arguments">()</code><a class="copiable-link" href="#index-odeset"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-odeset-1"><span><code class="def-type"><var class="var">odestruct</var> =</code> <strong class="def-name">odeset</strong> <code class="def-code-arguments">(<var class="var">"field1"</var>, <var class="var">value1</var>, <var class="var">"field2"</var>, <var class="var">value2</var>, …)</code><a class="copiable-link" href="#index-odeset-1"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-odeset-2"><span><code class="def-type"><var class="var">odestruct</var> =</code> <strong class="def-name">odeset</strong> <code class="def-code-arguments">(<var class="var">oldstruct</var>, <var class="var">"field1"</var>, <var class="var">value1</var>, <var class="var">"field2"</var>, <var class="var">value2</var>, …)</code><a class="copiable-link" href="#index-odeset-2"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-odeset-3"><span><code class="def-type"><var class="var">odestruct</var> =</code> <strong class="def-name">odeset</strong> <code class="def-code-arguments">(<var class="var">oldstruct</var>, <var class="var">newstruct</var>)</code><a class="copiable-link" href="#index-odeset-3"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-odeset-4"><span><strong class="def-name">odeset</strong> <code class="def-code-arguments">()</code><a class="copiable-link" href="#index-odeset-4"> ¶</a></span></dt>
<dd>
<p>Create or modify an ODE options structure.
</p>
<p>When called with no input argument and one output argument, return a new ODE
options structure that contains all possible fields initialized to their
default values. If no output argument is requested, display a list of
the common ODE solver options along with their default value.
</p>
<p>If called with name-value input argument pairs <var class="var">"field1"</var>,
<var class="var">"value1"</var>, <var class="var">"field2"</var>, <var class="var">"value2"</var>, … return a new
ODE options structure with all the most common option fields
initialized, <strong class="strong">and</strong> set the values of the fields <var class="var">"field1"</var>,
<var class="var">"field2"</var>, … to the values <var class="var">value1</var>, <var class="var">value2</var>,
<small class="enddots">...</small>
</p>
<p>If called with an input structure <var class="var">oldstruct</var> then overwrite the
values of the options <var class="var">"field1"</var>, <var class="var">"field2"</var>, … with
new values <var class="var">value1</var>, <var class="var">value2</var>, … and return the
modified structure.
</p>
<p>When called with two input ODE options structures <var class="var">oldstruct</var> and
<var class="var">newstruct</var> overwrite all values from the structure
<var class="var">oldstruct</var> with new values from the structure <var class="var">newstruct</var>.
Empty values in <var class="var">newstruct</var> will not overwrite values in
<var class="var">oldstruct</var>.
</p>
<p>The most commonly used ODE options, which are always assigned a value
by <code class="code">odeset</code>, are the following:
</p>
<dl class="table">
<dt><code class="code">AbsTol</code>: positive scalar | vector, def. <code class="code">1e-6</code></dt>
<dd><p>Absolute error tolerance.
</p>
</dd>
<dt><code class="code">BDF</code>: {<code class="code">"off"</code>} | <code class="code">"on"</code></dt>
<dd><p>Use BDF formulas in implicit multistep methods.
<em class="emph">Note</em>: This option is not yet implemented.
</p>
</dd>
<dt><code class="code">Events</code>: function_handle</dt>
<dd><p>Event function. An event function must have the form
<code class="code">[value, isterminal, direction] = my_events_f (t, y)</code>
</p>
</dd>
<dt><code class="code">InitialSlope</code>: vector</dt>
<dd><p>Consistent initial slope vector for DAE solvers.
</p>
</dd>
<dt><code class="code">InitialStep</code>: positive scalar</dt>
<dd><p>Initial time step size.
</p>
</dd>
<dt><code class="code">Jacobian</code>: matrix | function_handle</dt>
<dd><p>Jacobian matrix, specified as a constant matrix or a function of
time and state.
</p>
</dd>
<dt><code class="code">JConstant</code>: {<code class="code">"off"</code>} | <code class="code">"on"</code></dt>
<dd><p>Specify whether the Jacobian is a constant matrix or depends on the
state.
</p>
</dd>
<dt><code class="code">JPattern</code>: sparse matrix</dt>
<dd><p>If the Jacobian matrix is sparse and non-constant but maintains a
constant sparsity pattern, specify the sparsity pattern.
</p>
</dd>
<dt><code class="code">Mass</code>: matrix | function_handle</dt>
<dd><p>Mass matrix, specified as a constant matrix or a function of
time and state.
</p>
</dd>
<dt><code class="code">MassSingular</code>: {<code class="code">"maybe"</code>} | <code class="code">"yes"</code> | <code class="code">"on"</code></dt>
<dd><p>Specify whether the mass matrix is singular.
</p>
</dd>
<dt><code class="code">MaxOrder</code>: {<code class="code">5</code>} | <code class="code">4</code> | <code class="code">3</code> | <code class="code">2</code> | <code class="code">1</code></dt>
<dd><p>Maximum order of formula.
</p>
</dd>
<dt><code class="code">MaxStep</code>: positive scalar</dt>
<dd><p>Maximum time step value.
</p>
</dd>
<dt><code class="code">MStateDependence</code>: {<code class="code">"weak"</code>} | <code class="code">"none"</code> | <code class="code">"strong"</code></dt>
<dd><p>Specify whether the mass matrix depends on the state or only on time.
</p>
</dd>
<dt><code class="code">MvPattern</code>: sparse matrix</dt>
<dd><p>If the mass matrix is sparse and non-constant but maintains a
constant sparsity pattern, specify the sparsity pattern.
<em class="emph">Note</em>: This option is not yet implemented.
</p>
</dd>
<dt><code class="code">NonNegative</code>: scalar | vector</dt>
<dd><p>Specify elements of the state vector that are expected to remain
non-negative during the simulation.
</p>
</dd>
<dt><code class="code">NormControl</code>: {<code class="code">"off"</code>} | <code class="code">"on"</code></dt>
<dd><p>Control error relative to the 2-norm of the solution, rather than its
absolute value.
</p>
</dd>
<dt><code class="code">OutputFcn</code>: function_handle</dt>
<dd><p>Function to monitor the state during the simulation. For the form of
the function to use see <a class="pxref" href="#XREFodeplot"><code class="code">odeplot</code></a>.
</p>
</dd>
<dt><code class="code">OutputSel</code>: scalar | vector</dt>
<dd><p>Indices of elements of the state vector to be passed to the output
monitoring function.
</p>
</dd>
<dt><code class="code">Refine</code>: positive scalar</dt>
<dd><p>Specify whether output should be returned only at the end of each
time step or also at intermediate time instances. The value should be
a scalar indicating the number of equally spaced time points to use
within each timestep at which to return output.
</p>
</dd>
<dt><code class="code">RelTol</code>: positive scalar</dt>
<dd><p>Relative error tolerance.
</p>
</dd>
<dt><code class="code">Stats</code>: {<code class="code">"off"</code>} | <code class="code">"on"</code></dt>
<dd><p>Print solver statistics after simulation.
</p>
</dd>
<dt><code class="code">Vectorized</code>: {<code class="code">"off"</code>} | <code class="code">"on"</code></dt>
<dd><p>Specify whether <code class="code">odefcn</code> can be passed multiple values of the
state at once.
</p>
</dd>
</dl>
<p>Field names that are not in the above list are also accepted and
added to the result structure.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFodeget">odeget</a>.
</p></dd></dl>
<a class="anchor" id="XREFodeget"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-odeget"><span><code class="def-type"><var class="var">val</var> =</code> <strong class="def-name">odeget</strong> <code class="def-code-arguments">(<var class="var">ode_opt</var>, <var class="var">field</var>)</code><a class="copiable-link" href="#index-odeget"> ¶</a></span></dt>
<dt class="deftypefnx def-cmd-deftypefn" id="index-odeget-1"><span><code class="def-type"><var class="var">val</var> =</code> <strong class="def-name">odeget</strong> <code class="def-code-arguments">(<var class="var">ode_opt</var>, <var class="var">field</var>, <var class="var">default</var>)</code><a class="copiable-link" href="#index-odeget-1"> ¶</a></span></dt>
<dd>
<p>Query the value of the property <var class="var">field</var> in the ODE options structure
<var class="var">ode_opt</var>.
</p>
<p>If called with two input arguments and the first input argument
<var class="var">ode_opt</var> is an ODE option structure and the second input argument
<var class="var">field</var> is a string specifying an option name, then return the option
value <var class="var">val</var> corresponding to <var class="var">field</var> from <var class="var">ode_opt</var>.
</p>
<p>If called with an optional third input argument, and <var class="var">field</var> is
not set in the structure <var class="var">ode_opt</var>, then return the default value
<var class="var">default</var> instead.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFodeset">odeset</a>.
</p></dd></dl>
<a class="anchor" id="XREFodeplot"></a><span style="display:block; margin-top:-4.5ex;"> </span>
<dl class="first-deftypefn">
<dt class="deftypefn" id="index-odeplot"><span><code class="def-type"><var class="var">stop_solve</var> =</code> <strong class="def-name">odeplot</strong> <code class="def-code-arguments">(<var class="var">t</var>, <var class="var">y</var>, <var class="var">flag</var>)</code><a class="copiable-link" href="#index-odeplot"> ¶</a></span></dt>
<dd>
<p>Open a new figure window and plot the solution of an ode problem at each
time step during the integration.
</p>
<p>The types and values of the input parameters <var class="var">t</var> and <var class="var">y</var> depend on
the input <var class="var">flag</var> that is of type string. Valid values of <var class="var">flag</var>
are:
</p>
<dl class="table">
<dt><samp class="option"><code class="code">"init"</code></samp></dt>
<dd><p>The input <var class="var">t</var> must be a column vector of length 2 with the first and
last time step (<code class="code">[<var class="var">tfirst</var> <var class="var">tlast</var>]</code>. The input <var class="var">y</var>
contains the initial conditions for the ode problem (<var class="var">y0</var>).
</p>
</dd>
<dt><samp class="option"><code class="code">""</code></samp></dt>
<dd><p>The input <var class="var">t</var> must be a scalar double or vector specifying the time(s)
for which the solution in input <var class="var">y</var> was calculated.
</p>
</dd>
<dt><samp class="option"><code class="code">"done"</code></samp></dt>
<dd><p>The inputs should be empty, but are ignored if they are present.
</p></dd>
</dl>
<p><code class="code">odeplot</code> always returns false, i.e., don’t stop the ode solver.
</p>
<p>Example: solve an anonymous implementation of the
<code class="code">"Van der Pol"</code> equation and display the results while
solving.
</p>
<div class="example">
<div class="group"><pre class="example-preformatted">fvdp = @(t,y) [y(2); (1 - y(1)^2) * y(2) - y(1)];
opt = odeset ("OutputFcn", @odeplot, "RelTol", 1e-6);
sol = ode45 (fvdp, [0 20], [2 0], opt);
</pre></div></div>
<p>Background Information:
This function is called by an ode solver function if it was specified in
the <code class="code">"OutputFcn"</code> property of an options structure created with
<code class="code">odeset</code>. The ode solver will initially call the function with the
syntax <code class="code">odeplot ([<var class="var">tfirst</var>, <var class="var">tlast</var>], <var class="var">y0</var>, "init")</code>. The
function initializes internal variables, creates a new figure window, and
sets the x limits of the plot. Subsequently, at each time step during the
integration the ode solver calls <code class="code">odeplot (<var class="var">t</var>, <var class="var">y</var>, [])</code>.
At the end of the solution the ode solver calls
<code class="code">odeplot ([], [], "done")</code> so that odeplot can perform any clean-up
actions required.
</p>
<p><strong class="strong">See also:</strong> <a class="ref" href="#XREFodeset">odeset</a>, <a class="ref" href="#XREFodeget">odeget</a>, <a class="ref" href="#XREFode23">ode23</a>, <a class="ref" href="#XREFode45">ode45</a>.
</p></dd></dl>
</div>
<hr>
<div class="nav-panel">
<p>
Previous: <a href="Differential_002dAlgebraic-Equations.html">Differential-Algebraic 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>
|