Torsten

Version 0.89.0

A pharmacokinetics/pharmacodynamics library for Stan

Coupled ODE Model Function

Table of Contents

1 Description

When the ODE system consists of two subsystems in form of

\begin{align*} y_1^\prime &= f_1(t, y_1), \\\
y_2^\prime &= f_2(t, y_1, y_2), \end{align*}

with \(y_1\), \(y_2\), \(f_1\), and \(f_2\) being vector-valued functions, and \(y_1^\prime\) independent of \(y_2\), the solution can be accelerated if \(y_1\) admits an analytical solution which can be introduced into the ODE for \(y_2\) for numerical integration. This structure arises in PK/PD models, where \(y_1\) describes a forcing PK function and \(y_2\) the PD effects. In the example of a Friberg-Karlsson semi-mechanistic model(see below), we observe an average speedup of \(\sim 47 \pm 18 %\) when using the mix solver in lieu of the numerical integrator. In the context, currently the couple solver supports one- & two-compartment for PK model, and rk45 & bdf integrator for nonlinear PD model.

2 Usage

matrix pmx_solve_onecpt_[ rk45 || bdf ](reduced_ODE_system, int nOde, time, amt, rate, ii, evid, cmt, addl, ss, theta, biovar, tlag [, real rel_tol, real abs_tol, int max_step, real as_rel_tol, real as_abs_tol, int as_max_step ] );
matrix pmx_solve_twocpt_[ rk45 || bdf ](reduced_ODE_system, int nOde, time, amt, rate, ii, evid, cmt, addl, ss, theta, biovar, tlag [, real rel_tol, real abs_tol, int max_step, real as_rel_tol, real as_abs_tol, int as_max_step ] );

3 Arguments

  • reduced_ODE_rhs The system numerically solve (\(y_2\) in the above discussion, also called the reduced system and nOde the number of equations in the reduced system. The function that defines a reduced system has an almost identical signature to that used for a full system, but takes one additional argument: \(y_1\), the PK states, i.e. solution to the PK ODEs.

    vector reduced_ODE_rhs(real t, vector y2, vector y1, real[] theta, real[] x_r, int[] x_i)
    
  • nCmt The number of compartments. Equivalently, the dimension of the ODE system.

  • rel_tol The relative tolerance for numerical integration, default to 1.0E-6.

  • abs_tol The absolute tolerance for numerical integration, default to 1.0E-6.

  • max_step The maximum number of steps in numerical integration, default to \(10^6\).

  • See Tables in Section Events specification for the rest of arguments.

4 Return value

An (nPk + nOde) × nt matrix, where nt is the size of time, and nPk equals to 2 in pmx_solve_onecpt_ functions and 3 in pmx_solve_twocpt_ functions.

Last updated on 30 Jun 2021
Published on 25 Jun 2021
 Edit on GitHub