r2sundials-package {r2sundials}R Documentation

Wrapper for 'SUNDIALS' Solving ODE and Sensitivity Problem

Description

Wrapper for widely used 'SUNDIALS' software (SUite of Nonlinear and DIfferential/ALgebraic Equation Solvers) and more precisely to its 'CVODES' solver. It is aiming to solve ordinary differential equations (ODE) and optionally pending forward sensitivity problem. The wrapper is made 'R' friendly by allowing to pass custom parameters to user's callback functions. Such functions can be both written in 'R' and in 'C++' ('RcppArmadillo' flavor). In case of 'C++', performance is greatly improved so this option is highly advisable when performance matters. If provided, Jacobian matrix can be calculated either in dense or sparse format. In the latter case 'rmumps' package is used to solve corresponding linear systems. Root finding and pending event management are optional and can be specified as 'R' or 'C++' functions too. This makes them a very flexible tool for controlling the ODE system during the time course simulation. 'SUNDIALS' library was published in Hindmarsh et al. (2005) <doi:10.1145/1089014.1089020>. Interface for commonly used SUNDIALS library for solving ODE system. It is made R friendly by allowing to pass custom parameters to user's callback functions. They can be both written in R and in C++ (Rcpp flavor). In case of C++, performance is greatly improved so this option is highly advisable when performance matters. Jacobian matrix can be specified as either dense or sparse matrix. In the latter case rmumps package is used to solve corresponding linear systems. Root finding and corresponding event management are optional and can be too specified as R or C++ functions which makes them a very flexible tool for controlling the ODE system during the time course simulation.

Details

The DESCRIPTION file:

Package: r2sundials
Type: Package
Title: Wrapper for 'SUNDIALS' Solving ODE and Sensitivity Problem
Version: 6.5.0-4
Date: 2023-12-08
Authors@R: c( person("Serguei", "Sokol", role=c("cre", "aut"), email="sokol@insa-toulouse.fr"), person("Carol S.", "Woodward", role="ctb"), person("Daniel R.", "Reynolds", role="ctb"), person("Alan C.", "Hindmarsh", role="ctb"), person("David J.", "Gardner", role="ctb"), person("Cody J.", "Balos", role="ctb"), person("Radu", "Serban", role="ctb"), person("Scott D.", "Cohen", role="ctb"), person("Peter N.", "Brown", role="ctb"), person("George", "Byrne", role="ctb"), person("Allan G.", "Taylor", role="ctb"), person("Steven L.", "Lee", role="ctb"), person("Keith E.", "Grant", role="ctb"), person("Aaron", "Collier", role="ctb"), person("Lawrence E.", "Banks", role="ctb"), person("Steve", "Smith", role="ctb"), person("Cosmin", "Petra", role="ctb"), person("Slaven", "Peles", role="ctb"), person("John", "Loffeld", role="ctb"), person("Dan", "Shumaker", role="ctb"), person("Ulrike", "Yang", role="ctb"), person("James", "Almgren-Bell", role="ctb"), person("Shelby L.", "Lockhart", role="ctb"), person("Hilari C.", "Tiedeman", role="ctb"), person("Ting", "Yan", role="ctb"), person("Jean M.", "Sexton", role="ctb"), person("Chris", "White", role="ctb"), person("Lawrence Livermore National Security", role="cph"), person("Southern Methodist University", role="cph"), person("INSA", role="cph"), person("INRAE", role="cph"), person("CNRS", role="cph") )
Maintainer: Serguei Sokol <sokol@insa-toulouse.fr>
Description: Wrapper for widely used 'SUNDIALS' software (SUite of Nonlinear and DIfferential/ALgebraic Equation Solvers) and more precisely to its 'CVODES' solver. It is aiming to solve ordinary differential equations (ODE) and optionally pending forward sensitivity problem. The wrapper is made 'R' friendly by allowing to pass custom parameters to user's callback functions. Such functions can be both written in 'R' and in 'C++' ('RcppArmadillo' flavor). In case of 'C++', performance is greatly improved so this option is highly advisable when performance matters. If provided, Jacobian matrix can be calculated either in dense or sparse format. In the latter case 'rmumps' package is used to solve corresponding linear systems. Root finding and pending event management are optional and can be specified as 'R' or 'C++' functions too. This makes them a very flexible tool for controlling the ODE system during the time course simulation. 'SUNDIALS' library was published in Hindmarsh et al. (2005) <doi:10.1145/1089014.1089020>.
BugReports: https://github.com/sgsokol/r2sundials/issues
Depends: R (>= 3.0.2), rmumps (>= 5.2.1-6)
License: GPL (>=2)
Imports: Rcpp (>= 1.0.0)
LinkingTo: Rcpp, RcppArmadillo, rmumps (>= 5.2.1-6)
Suggests: RcppXPtrUtils, slam, RUnit, deSolve, RcppArmadillo
RoxygenNote: 7.2.3
Encoding: UTF-8
NeedsCompilation: yes
Biarch: FALSE
Author: Serguei Sokol [cre, aut], Carol S. Woodward [ctb], Daniel R. Reynolds [ctb], Alan C. Hindmarsh [ctb], David J. Gardner [ctb], Cody J. Balos [ctb], Radu Serban [ctb], Scott D. Cohen [ctb], Peter N. Brown [ctb], George Byrne [ctb], Allan G. Taylor [ctb], Steven L. Lee [ctb], Keith E. Grant [ctb], Aaron Collier [ctb], Lawrence E. Banks [ctb], Steve Smith [ctb], Cosmin Petra [ctb], Slaven Peles [ctb], John Loffeld [ctb], Dan Shumaker [ctb], Ulrike Yang [ctb], James Almgren-Bell [ctb], Shelby L. Lockhart [ctb], Hilari C. Tiedeman [ctb], Ting Yan [ctb], Jean M. Sexton [ctb], Chris White [ctb], Lawrence Livermore National Security [cph], Southern Methodist University [cph], INSA [cph], INRAE [cph], CNRS [cph]

Index of help topics:

get_cnst                Get Predefined Constant
r2cvodes                Solving ODE System and Sensitivity Equations
r2sundials-package      Wrapper for 'SUNDIALS' Solving ODE and
                        Sensitivity Problem

The main function of the package is r2cvodes() which wraps and converts its arguments in data structures and parameters convenient for calling cvodes() from SUNDIALS library.

When using r2sundials::r2cvodes(), some callback functions have to be written by the user either in R or in C++ (RcppArmadillo flavor). One of them is mandatory and defines the rhs of the ODE system to solve. Others are optional and can be used to

Note that if 'SUNDIALS' 'CVODES' is compiled with parameter SUNDIALS_INDEX_SIZE set to 32, some memory copying is spared if a C++ function calculating sparse Jacobian is provided by user.

The version number of this package if made of CVODES original version, e.g. 5.0.0 followed by one digit for R wrapper release. This can give 5.0.0-1.

Author(s)

Serguei Sokol [cre, aut], Carol S. Woodward [ctb], Daniel R. Reynolds [ctb], Alan C. Hindmarsh [ctb], David J. Gardner [ctb], Cody J. Balos [ctb], Radu Serban [ctb], Scott D. Cohen [ctb], Peter N. Brown [ctb], George Byrne [ctb], Allan G. Taylor [ctb], Steven L. Lee [ctb], Keith E. Grant [ctb], Aaron Collier [ctb], Lawrence E. Banks [ctb], Steve Smith [ctb], Cosmin Petra [ctb], Slaven Peles [ctb], John Loffeld [ctb], Dan Shumaker [ctb], Ulrike Yang [ctb], James Almgren-Bell [ctb], Shelby L. Lockhart [ctb], Hilari C. Tiedeman [ctb], Ting Yan [ctb], Jean M. Sexton [ctb], Chris White [ctb], Lawrence Livermore National Security [cph], Southern Methodist University [cph], INSA [cph], INRAE [cph], CNRS [cph] Maintainer: Serguei Sokol <sokol@insa-toulouse.fr>

References

Original SUNDIALS CVODES user documentation (search for cvs_guide.pdf)

See Also

deSolve

Examples

  # a very simple ODE for exponential transition from 0 to 1: y'=1-y, y(0)=0
  frhs=function(t, y, param, psens) 1-y
  ti=seq(0, 5, length.out=101)
  y0=0
  res_exp=r2sundials::r2cvodes(y0, ti, frhs)
  plot(ti, res_exp[1,], t="l", xlab="Time", ylab="Y")

[Package r2sundials version 6.5.0-4 Index]