Scilab Website | Contribute with GitLab | Mailing list archives | ATOMS toolboxes
Scilab Online Help
2024.0.0 - Русский


CVode

CVode (short for C-language Variable-coefficients ODE solver) is a numerical solver providing an efficient and stable method to solve Ordinary Differential Equations (ODEs) Initial Value Problems. It uses either BDF or Adams as implicit integration method, and Newton or Functional iterations.

Description

Called by xcos, CVode (short for C-language Variable-coefficients ODE solver) is a numerical solver providing an efficient and stable method to solve Initial Value Problems of the form:

\begin{eqnarray}
                \dot{y} = f(t,y), \hspace{3 mm} y(t_0) = y_0, \hspace{3 mm} y \in R^N
                \end{eqnarray}

Starting with y0 , CVode approximates yn+1 with the formula:

\begin{eqnarray}
                \sum_{i=0}^{K_1} \alpha_{n,i} y_{n-i} + h_n\sum_{i=0}^{K_2} \beta_{n,i} \dot{y}_{n-i} = 0,\hspace{10 mm} (1)
                \end{eqnarray}

with yn the approximation of y(tn) , and hn = tn - tn-1 the step size.

These implicit methods are characterized by their respective order q, which indicates the number of intermediate points required to compute yn+1 .

This is where the difference between BDF and Adams intervenes (Backward Differentiation Formula and Adams-Moulton formula):

If the problem is stiff, the user should select BDF:

  • q, the order of the method, is set between 1 and 5 (automated),
  • K1 = q and K2 = 0.

In the case of nonstiffness, Adams is preferred:

  • q is set between 1 and 12 (automated),
  • K1 = 1 and K2 = q.

The coefficients are fixed, uniquely determined by the method type, its order, the history of the step sizes, and the normalization αn, 0 = -1 .

For either choice and at each step, injecting this integration in (1) yields the nonlinear system:

G(y_n)\equiv y_n-h_n\beta_{n,0}f(t_n,y_n)-a_n=0, \hspace{2 mm} where \hspace{2 mm} a_n\equiv \sum_{i>0} (\alpha_{n,i} y_{n-i} + h_n\beta_{n,i}\dot{y}_{n-i})

This system can be solved by either Functional or Newton iterations, described hereafter.

In both following cases, the initial "predicted" yn(0) is explicitly computed from the history data, by adding derivatives.

  • Functional: this method only involves evaluations of f, it simply computes yn(0) by iterating the formula:

    y_{n(m+1)} = h_n β_{n,0} f(t_n,y_{n(m)}) + a_n where \hspace{2 mm} a_n\equiv \sum_{i>0} (\alpha_{n,i} y_{n-i} + h_n\beta_{n,i}\dot{y}_{n-i})

  • Newton: here, we use an implemented direct dense solver on the linear system:

    M[y_{n(m+1)}-y_{n(m)}]=-G(y_{n(m)}), \hspace{4 mm} M \approx I-\gamma J, \hspace{2 mm} J=\frac{\partial f}{\partial y}, \hspace{2 mm} and \hspace{2 mm} \gamma = h_n\beta_{n,0}

In both situations, CVode uses the history array to control the local error yn(m) - yn(0) and recomputes hn if that error is not satisfying.

The recommended choices are BDF / Newton for stiff problems and Adams / Functional for the nonstiff ones.

The function is called in between activations, because a discrete activation may change the system.

Following the criticality of the event (its effect on the continuous problem), we either relaunch the solver with different start and final times as if nothing happened, or, if the system has been modified, we need to "cold-restart" the problem by reinitializing it anew and relaunching the solver.

Averagely, CVode accepts tolerances up to 10-16. Beyond that, it returns a Too much accuracy requested error.

Examples

The integral block returns its continuous state, we can evaluate it with BDF / Newton by running the example:

// Import the diagram and set the ending time
loadScicos();
loadXcosLibs();
importXcosDiagram("SCI/modules/xcos/examples/solvers/ODE_Example.zcos");
scs_m.props.tf = 5000;

// Select the solver BDF / Newton
scs_m.props.tol(6) = 1;

// Start the timer, launch the simulation and display time
tic();
try xcos_simulate(scs_m, 4); catch disp(lasterror()); end
t = toc();
disp("Time for BDF / Newton:", t);

The Scilab console displays:

Time for BDF / Newton:
 13.438
            

Now, in the following script, we compare the time difference between the methods by running the example with the four solvers in turn: Open the script

Results:

Time for BDF / Newton:
 18.894

Time for BDF / Functional:
 18.382

Time for Adams / Newton:
 10.368

Time for Adams / Functional:
 9.815
            

The results show that for a simple nonstiff continuous problem, Adams / Functional is fastest.

See also

  • LSodar — LSodar (short for Livermore Solver for Ordinary Differential equations, with Automatic method switching for stiff and nonstiff problems, and with Root-finding) is a numerical solver providing an efficient and stable method to solve Ordinary Differential Equations (ODEs) Initial Value Problems.
  • Runge-Kutta 4(5) — Runge-Kutta is a numerical solver providing an efficient explicit method to solve Ordinary Differential Equations (ODEs) Initial Value Problems.
  • Dormand-Prince 4(5) — Dormand-Prince is a numerical solver providing an efficient explicit method to solve Ordinary Differential Equations (ODEs) Initial Value Problems.
  • Implicit Runge-Kutta 4(5) — Numerical solver providing an efficient and stable implicit method to solve Ordinary Differential Equations (ODEs) Initial Value Problems.
  • Crank-Nicolson 2(3) — Crank-Nicolson is a numerical solver based on the Runge-Kutta scheme providing an efficient and stable implicit method to solve Ordinary Differential Equations (ODEs) Initial Value Problems. Called by xcos.
  • IDA — Implicit Differential Algebraic equations system solver, providing an efficient and stable method to solve Differential Algebraic Equations system (DAEs) Initial Value Problems.
  • DDaskr — Double-precision Differential Algebraic equations system Solver with Krylov method and Rootfinding: numerical solver providing an efficient and stable method to solve Differential Algebraic Equations systems (DAEs) Initial Value Problems
  • Comparisons — This page compares solvers to determine which one is best fitted for the studied problem.
  • ode — программа решения обыкновенных дифференциальных уравнений
  • ode_discrete — программа решения обыкновенных дифференциальных уравнений, моделирование дискретного времени
  • ode_root — программа решения обыкновенных дифференциальных уравнений с поиском корней
  • odedc — программа решения дискретно-непрерывных ОДУ
  • impl — дифференциальное алгебраическое уравнение
Report an issue
<< LSodar Solvers Runge-Kutta 4(5) >>

Copyright (c) 2022-2024 (Dassault Systèmes)
Copyright (c) 2017-2022 (ESI Group)
Copyright (c) 2011-2017 (Scilab Enterprises)
Copyright (c) 1989-2012 (INRIA)
Copyright (c) 1989-2007 (ENPC)
with contributors
Last updated:
Tue Oct 24 14:37:12 CEST 2023