Ordinary Differential Equations (ODEs) describe the rate of change of dependent variables with respect to a single independent variable and are used in many fields to model behavior of the system. There are many good C libraries available to solve (i.e., integrate systems of ODEs) and SUNDIALS available from the Lawrence Livermore National Laboratory is a one of the most popular and well-respected C library for solving non-stiff and stiff systems of ODEs.
Currently, this package provides an interface to the CVODE function (serial version) in the library which is used to solve ODEs (or Initial Value Problems). In future, we plan to provide interface for the other solvers in the library also. Right now, this package serves as a test case for providing an interface to the SUNDIALS library for R users.
One of the advantage of using this package is that all the source code of the SUNDIALS
library is bundled with the package itself, so it does not require the SUNDIALS
library to be installed on the machine separately (which is sometimes not trivial for Windows).
As described in the link above, the problem is from chemical kinetics, and consists of the following three rate equations:
\[ \begin{aligned} \frac{dy_1}{dt} &= -.04 \times y_1 + 10^4 \times y_2 \times y_3 \\ \frac{dy_2}{dt} &= .04 \times y_1 - 10^4 \times y_2 \times y_3 - 3 \times 10^7 \times y_2^2 \\ \frac{dy_3}{dt} &= 3 \times 10^7 \times y_2^2 \end{aligned} \]
on the interval from time
\[ t = 0.0 \] to
\[ t = 4 \times 10^{10} \], with initial conditions:
\[ y_1 = 1.0 , ~y_2 = y_3 = 0 \]
The problem is stiff.
The original example , While integrating the system, also uses the rootfinding feature to find the points at which
\[ y_1 = 1 \times 10^{-4} \]
or at which
\[ y_3 = 0.01 \],
but currently root-finding is not supported in this version. As, in the original example this package also solves the problem with the BDF method, Newton iteration with the SUNDENSE dense linear solver, however, without a user-supplied Jacobian routine (unlike the original example). It uses a scalar relative tolerance and a vector absolute tolerance (which can be provided as an input). Output is printed in decades from
\[ t = 0.4 \]
to
\[ t = 4 \times 10^{10} \].
Differential equations can be written as an R
function or as an Rcpp
function. Differential equations function must be written as
function(t, y){
# code to write differential equations
}
The key aspect to keep in mind is that the signature of the function must be function(t,y)
. As an example, we try to solve the cv_Roberts_dns.c
problem described above, the original code can be found here. An example of an R
function is as follows:
ODE_R <- function(t, y){
## initialize the derivative vector
ydot <- vector(mode = "numeric", length = length(y))
ydot[1] = -0.04 * y[1] + 1e04 * y[2] * y[3]
ydot[3] = 3e07 * y[2] * y[2]
ydot[2] = -ydot[1] - ydot[3]
ydot
}
Also, since this package using Rcpp
to bundle the C code, we can use the notation used in Rcpp
to describe the system of ODEs. The cv_Roberts_dns
problem describe above can be described in an Rcpp
function as follows (indices in C++
start from 0
, functions need to declare their return type, here NumericVector
and every expression ends in a semicolon, ;
) :
#include <Rcpp.h>
using namespace Rcpp;
// [[Rcpp::export]]
NumericVector ODE_Rcpp (double t, NumericVector y){
// Initialize ydot filled with zeros
NumericVector ydot(y.length());
ydot[0] = -0.04 * y[0] + 1e04 * y[1] * y[2];
ydot[2] = 3e07 * y[1] * y[1];
ydot[1] = -ydot[0] - ydot[2];
return ydot;
}
The above is a re-write of the cvRoberts_dns
example in the documentation of CVODE
. The original example can be found the document here.
The entire R
file to create right hand side of ODE function (which calculates rates of change) is as follows (also found in /inst/examples/cv_Roberts_dns.r
):
# ODEs described by an R function
ODE_R <- function(t, y){
# vector containing the right hand side gradients
ydot = vector(mode = "numeric", length = length(y))
# R indices start from 1
ydot[1] = -0.04 * y[1] + 1e04 * y[2] * y[3]
ydot[3] = 3e07 * y[2] * y[2]
ydot[2] = -ydot[1] - ydot[3]
ydot
}
# ODEs can also be described using Rcpp
Rcpp::sourceCpp(code = '
#include <Rcpp.h>
using namespace Rcpp;
// [[Rcpp::export]]
NumericVector ODE_Rcpp (double t, NumericVector y){
// Initialize ydot filled with zeros
NumericVector ydot(y.length());
ydot[0] = -0.04 * y[0] + 1e04 * y[1] * y[2];
ydot[2] = 3e07 * y[1] * y[1];
ydot[1] = -ydot[0] - ydot[2];
return ydot;
}')
# Generate time vector, IC and call cvode to solve the equations
# R code to genrate time vector, IC and solve the equations
time_vec <- c(0.0, 0.4, 4.0, 40.0, 4E2, 4E3, 4E4, 4E5, 4E6, 4E7, 4E8, 4E9, 4E10)
IC <- c(1,0,0)
reltol <- 1e-04
abstol <- c(1e-8,1e-14,1e-6)
## Solving the ODEs using cvode function
df1 <- cvode(time_vec, IC, ODE_R , reltol, abstol) ## using R
df2 <- cvode(time_vec, IC, ODE_Rcpp , reltol, abstol) ## using Rcpp
## Check that both solutions are identical
identical(df1, df2)
The final output is the df1
matrix in which first column is time, second, third and fourth column are the values of y1
, y2
and y3
respectively.
> df1
[,1] [,2] [,3] [,4]
[1,] 0e+00 1.000000e+00 0.000000e+00 0.00000000
[2,] 4e-01 9.851641e-01 3.386242e-05 0.01480205
[3,] 4e+00 9.055097e-01 2.240338e-05 0.09446793
[4,] 4e+01 7.158016e-01 9.185043e-06 0.28418924
[5,] 4e+02 4.505209e-01 3.222826e-06 0.54947590
[6,] 4e+03 1.832217e-01 8.943516e-07 0.81677741
[7,] 4e+04 3.898091e-02 1.621669e-07 0.96101893
[8,] 4e+05 4.936971e-03 1.984450e-08 0.99506301
[9,] 4e+06 5.170103e-04 2.069098e-09 0.99948299
[10,] 4e+07 5.204927e-05 2.082078e-10 0.99994795
[11,] 4e+08 5.184946e-06 2.073989e-11 0.99999482
[12,] 4e+09 5.246212e-07 2.098486e-12 0.99999948
[13,] 4e+10 6.043000e-08 2.417200e-13 0.99999994
The package sundialr
provides a way to interface with the famous ‘SUNDIALS’ C library to solver initial value problems. The package allows differential equation to be written in R
or using Rcpp
. Currently only a single initialization is supported, but we plan to add repeated initializations (e.g., repeated dosing) in near future. In addition, other solvers from the C library such as cvodes
and IDA
will be added soon.