Comparison of ODE solvers: Example NAODE

Alfonso R. Reyes

2017-11-08

In this vignette we will use the example from the book “Numerical Solutions of Ordinary Differential Equations” by Atkinson, Han and Stewart. Google Books

Build the ODE class (without class accumulator)

Comparison of solutions: RK45 vs analytical solution

For the differential equation:

\[\dfrac{dy}{dt} = - y\] the analytical solution is: \[y(t) = e^{-t}\]

Find the errors if: \[y(0) = 1\]

library(rODE)

# ODETest.R

setClass("ODETest", slots = c(
    stack = "environment"           # environment object inside the class
    ),
    contains = c("ODE")
    )

setMethod("initialize", "ODETest", function(.Object, ...) {
    .Object@stack$n <-  0               # "n" belongs to the class environment
    .Object@state   <- c(1.0, 0.0)
    return(.Object)
})

setMethod("getExactSolution", "ODETest", function(object, t, ...) {
    # analytical solution
    return(exp(-t))
})

setMethod("getState", "ODETest", function(object, ...) {
    object@state
})

setMethod("getRate", "ODETest", function(object, state, ...) {
    object@rate[1] <- -state[1]
    object@rate[2] <-  1                        # rate of change of time, dt/dt
    object@stack$n <-  object@stack$n + 1       # add 1 to the rate count
    object@rate
})

# constructor
ODETest <- function() {
    odetest <- new("ODETest")
    odetest
}
## [1] "initialize"
## [1] "getExactSolution"
## [1] "getState"
## [1] "getRate"

Euler

ComparisonEulerApp <- function(stepSize) {
    ode <- new("ODETest")
    ode_solver <- Euler(ode)
    ode_solver <- setStepSize(ode_solver, stepSize)
    time <-  0
    rowVector <- vector("list")
    i <-  1
    while (time < 5) {
        state <- getState(ode_solver@ode)
        time <- state[2]
        error <- getExactSolution(ode_solver@ode, time) - state[1]
        rowVector[[i]] <- list(step_size = stepSize, 
                               t = time, 
                               y_t = state[1], 
                               error = error, 
                               rel_err = error / getExactSolution(ode_solver@ode, time),
                               steps = ode_solver@ode@stack$n
                               )
        ode_solver <- step(ode_solver)
        i <- i + 1
    }
    data.table::rbindlist(rowVector)
}
dt <- ComparisonEulerApp(stepSize = 0.2)
dt[round(t,1) %in% c(1.0, 2.0, 3.0, 4.0, 5.0)]
##    step_size t         y_t       error   rel_err steps
## 1:       0.2 1 0.327680000 0.040199441 0.1092734     5
## 2:       0.2 2 0.107374182 0.027961101 0.2066061    10
## 3:       0.2 3 0.035184372 0.014602696 0.2933030    15
## 4:       0.2 4 0.011529215 0.006786424 0.3705262    20
## 5:       0.2 5 0.003777893 0.002960054 0.4393109    25
# get a summary table for different step sizes
get_table <- function(stepSize) {
    dt <- ComparisonEulerApp(stepSize)
    dt[round(t,2) %in% c(1.0, 2.0, 3.0, 4.0, 5.0)]
}

step_sizes <- c(0.2, 0.1, 0.05)
dt_li <- lapply(step_sizes, get_table)
data.table::rbindlist(dt_li)
##     step_size t         y_t        error    rel_err steps
##  1:      0.20 1 0.327680000 0.0401994412 0.10927341     5
##  2:      0.20 2 0.107374182 0.0279611008 0.20660614    10
##  3:      0.20 3 0.035184372 0.0146026963 0.29330300    15
##  4:      0.20 4 0.011529215 0.0067864238 0.37052619    20
##  5:      0.20 5 0.003777893 0.0029600538 0.43931094    25
##  6:      0.10 1 0.348678440 0.0192010011 0.05219373    10
##  7:      0.10 2 0.121576655 0.0137586286 0.10166328    20
##  8:      0.10 3 0.042391158 0.0073959101 0.14855083    30
##  9:      0.10 4 0.014780883 0.0035347559 0.19299114    40
## 10:      0.10 5 0.005153775 0.0015841718 0.23511194    50
## 11:      0.05 1 0.358485922 0.0093935188 0.02553423    20
## 12:      0.05 2 0.128512157 0.0068231267 0.05041647    40
## 13:      0.05 3 0.046069799 0.0037172694 0.07466335    60
## 14:      0.05 4 0.016515374 0.0018002645 0.09829111    80
## 15:      0.05 5 0.005920529 0.0008174178 0.12131555   100

Using RK4

ComparisonRK4App <- function(stepSize) {
    ode <- new("ODETest")
    ode_solver <- RK4(ode)
    ode_solver <- setStepSize(ode_solver, stepSize)
    time <-  0
    rowVector <- vector("list")
    i <-  1
    while (time < 5) {
        state <- getState(ode_solver@ode)
        time <- state[2]
        error <- getExactSolution(ode_solver@ode, time) - state[1]
        rowVector[[i]] <- list(step_size = stepSize, 
                               t = time, 
                               y_t = state[1], 
                               error = error, 
                               rel_err = error / getExactSolution(ode_solver@ode, time),
                               steps = ode_solver@ode@stack$n
                               )
        ode_solver <- step(ode_solver)
        i <- i + 1
    }
    data.table::rbindlist(rowVector)
}

# get a summary table for different step sizes
get_table <- function(stepSize) {
    dt <- ComparisonRK4App(stepSize)
    dt[round(t,2) %in% c(1.0, 2.0, 3.0, 4.0, 5.0)]
}

step_sizes <- c(0.2, 0.1, 0.05)
dt_li <- lapply(step_sizes, get_table)
data.table::rbindlist(dt_li)
##     step_size t         y_t         error       rel_err steps
##  1:      0.20 1 0.367885238 -5.796954e-06 -1.575775e-05    20
##  2:      0.20 2 0.135339548 -4.265194e-06 -3.151576e-05    40
##  3:      0.20 3 0.049789422 -2.353634e-06 -4.727401e-05    60
##  4:      0.20 4 0.018316793 -1.154481e-06 -6.303251e-05    80
##  5:      0.20 5 0.006738478 -5.308913e-07 -7.879125e-05   100
##  6:      0.10 1 0.367879774 -3.332411e-07 -9.058431e-07    40
##  7:      0.10 2 0.135335528 -2.451852e-07 -1.811687e-06    80
##  8:      0.10 3 0.049787204 -1.352979e-07 -2.717532e-06   120
##  9:      0.10 4 0.018315705 -6.636447e-08 -3.623377e-06   160
## 10:      0.10 5 0.006737978 -3.051767e-08 -4.529224e-06   200
## 11:      0.05 1 0.367879461 -1.997610e-08 -5.430066e-08    80
## 12:      0.05 2 0.135335298 -1.469759e-08 -1.086013e-07   160
## 13:      0.05 3 0.049787076 -8.110413e-09 -1.629020e-07   240
## 14:      0.05 4 0.018315643 -3.978206e-09 -2.172027e-07   320
## 15:      0.05 5 0.006737949 -1.829375e-09 -2.715033e-07   400