ODEBND
ODEBND (version 2.0): Enclosing the Reachable Set of Parametric Nonlinear ODEs
Author:
Benoit C. Chachuat (b.chachuat@imperial.ac.uk), Boris Houska (bhouska@sjtu.edu.cn), Radoslav Paulen (Radoslav.Paulen@bci.tu-dortmund.de), Nikola Peric (nikola.peric09@imperial.ac.uk), Mario Villanueva (mario.villanueva10@imperial.ac.uk)
Version:
2.0
Date:
2012-2014

ODEBND is a collection of C++ class for computing enclosures of the reachable set of parametric nonlinear ordinary differential equations (ODEs) of the form

\[\forall t\in[t_0,t_{\rm f}]\,, \quad \dot{x}(t)=f(t, x(t),p) \quad\text{with}\quad x(t_0)=h(p), \]

where \(p\in P\subset\mathbb{R}^{n_p}\) denotes the parameter vector; \(x(t)\in\mathbb{R}^{n_x}\), the state vector at a given time \(t\); \(f: \mathbb{R}^{n_x}\times\mathbb{R}^{n_p}\to\mathbb{R}^{n_x}\), the right-hand side function; and \(h: \mathbb{R}^{n_p}\to\mathbb{R}^{n_x}\), the initial value function. It is assumed that \(f\) and \(h\) are factorable and sufficiently many times continuously differentiable.

Existing methods for set-valued ODE integration can be broadly classified into either continuous or discretized methods. Implementations for both methods are available in ODEBND:

  • Continuous-time set-valued integration construct a set of auxiliary ODEs whose solutions provide the desired enclosures. Two different approaches are implemented in the class mc::ODEBND_GSL (header file odebnd_gsl.hpp), which enable the propagation of convex enclosures and nonconvex enclosures based on interval arithmetic and polynomial model arithmetic, respectively. The types of polynomial models currently available are Taylor models and Chebyshev models, both through the library MC++ (ver 2.0). Moreover, the following contractors can be used to enhance the propagation:
    • differential inequalities contractor [Lakshmikantham & Leela, 1969; Walter, 1970; Harrison, 1977]
    • ellipsoidal contractor [Houska et al., 2012]
    Combination of these contractors with Taylor models are described, for instance, in [Chachuat & Villanueva, 2012; Villanueva et al., 2013a], and this combination is analogous for Chebyshev models. Moreover, the construction and convergence analysis of these methods are described in [Villanueva et al., in press].
  • Discretized set-valued integration proceed by first discretizing the integration horizon into finite steps and then propagating the reachable set enclosures between time steps by accounting for the discretization (truncation) errors. The class mc::ODEBND_VAL (header file odebnd_val.hpp) implements this approach based on polynomial model arithmetic. The types of polynomial models currently available are Taylor models and Chebyshev models, both through the library MC++ (ver 2.0). Besides propagation of polynomial models with interval remainders, mc::ODEBND_VAL enables ellipsoidal contractors [Houska et al., 2013]. A description and stability analysis of the discretized set-valued integrator can be found in [Villanueva et al, 2014; Houska et al., submitted].

The class mc::ODEBND_GSL and mc::ODEBND_VAL in ODEBND have three template parameters: T, PMT and PVT

  • T is the class that implements the underlying interval arithmetic. Both verified and non-verified types are currently supported. For verified interval computations [optional], the libraries PROFIL (header file mcprofil.hpp) and FILIB++ (header file mcfilib.hpp) are supported; the non-verified interval type mc::Interval in MC++ (header file interval.hpp) can also be used. More generally, any interval type can be used provided that the templated structure mc::Op<> is instantiated for that type - See MC++ documentation for details.
  • PMT and PVT are classes implementing polynomial model arithmetic. Currently available types through MC++ are: Taylor models (classes mc::TModel and mc::TVar, default template arguments); and Chebyshev models (classes mc::CModel and mc::CVar).

How Do I Obtain, Install and Run ODEBND?

ODEBND is released as open source code under the Eclipse Public License (EPL).

A number of third-party libraries are needed in order to compile and run mc::ODEBND_GSL and mc::ODEBND_VAL:

  • Library MC++ (version 2.0), for computing interval bounds and Taylor/Chebyshev models for factorable functions. Note that MC++ uses FADBAD++ (version 2.1) for automatic differentiation (forward, backward and Taylor methods of automatic differentiation) and BLAS/LAPACK for linear algebra computations. Moreover, MC++ is compatible with PROFIL and FILIB++ for verified interval arithmetic - See MC++ documentation for details.
  • Library GNU Scientific Library (GSL, version 1.15), for numerical integration of ODEs using the integrator gsl_odeiv2.

On distribution, the main directory should contain the files AUTHORS,CHANGELOG, INSTALL, LICENSE and README, as well as 3 subdirectories: src, test, and doc. The library comes with a Doxygen documentation, which can be accessed by opening the file doc/html/index.html with any html interpreter. Also, a test example is provided in the test subdirectory. Note that ODEBND consists of header files only, which can be found in the src directory, and we have only tested it under Linux/Ubuntu using the GCC compiler.

Running ODEBND first requires that you alter the makefile in order to point to the third-party libraries. A test example is provided in the subdirectory test. Change to this directory and type the following from the command line:

  $ make test1

This creates the executable file test1, which can be executed by doing:

  $ ./test1

The following output should be produced:


NON_VALIDATED INTEGRATION - APPROXIMATE ENCLOSURE OF REACHABLE SET:

 @t = 0.000000e+00 :
 x[0] = [  1.20000e+00 :  1.20000e+00 ]
 x[1] = [  1.10000e+00 :  1.10000e+00 ]
 @t = 5.000000e+00 :
 x[0] = [  8.03136e-01 :  8.18644e-01 ]
 x[1] = [  9.87776e-01 :  1.08750e+00 ]


DISCRETIZED SET-VALUED INTEGRATION - POLYNOMIAL MODEL ENCLOSURE OF REACHABLE SET:

 @t = 0.0000e+00 :
 x[0] = 
   a0    =  1.20000e+00
   R     = [  0.00000e+00 :  0.00000e+00 ]
   B     =  [  1.20000e+00 :  1.20000e+00 ]

 x[1] = 
   a0    =  1.10000e+00
   R     = [  0.00000e+00 :  0.00000e+00 ]
   B     =  [  1.10000e+00 :  1.10000e+00 ]

 Rx =
center:
 0.00000e+00
 0.00000e+00
shape:
 0.00000e+00 {0.00000e+00}
 0.00000e+00  0.00000e+00 

 @t = 5.0000e+00 :
 x[0] = 
   a0    =  8.08494e-01     0
   a1    = -7.58680e-03     1
   a2    =  2.54910e-03     2
   a3    = -1.57304e-05     3
   R     =  [ -2.73030e-05 :  2.73030e-05 ]
   B     =  [  8.02559e-01 :  8.19192e-01 ]

 x[1] = 
   a0    =  1.03753e+00     0
   a1    = -4.99778e-02     1
   a2    =  1.08644e-04     2
   a3    =  1.17340e-04     3
   R     =  [ -2.91695e-05 :  2.91695e-05 ]
   B     =  [  9.87098e-01 :  1.08818e+00 ]

 Rx =
center:
 0.00000e+00
 0.00000e+00
shape:
 7.45451e-10 {4.20329e-11}
 4.20329e-11  8.50861e-10 

 No STEPS  854
 No EVALATIONS   TRHS:  854 (IA)  854 (PM)   FTRHS: 0 (IA)  854 (PM)
 CPU TIME (SEC)     0.14775


CONTINUOUS SET-VALUED INTEGRATION - POLYNOMIAL MODEL ENCLOSURE OF REACHABLE SET:

 @t = 0.0000e+00 :
 x[0] = 
   a0    =  1.20000e+00
   R     = [  0.00000e+00 :  0.00000e+00 ]
   B     =  [  1.20000e+00 :  1.20000e+00 ]

 x[1] = 
   a0    =  1.10000e+00
   R     = [  0.00000e+00 :  0.00000e+00 ]
   B     =  [  1.10000e+00 :  1.10000e+00 ]

 Rx =
center:
 0.00000e+00
 0.00000e+00
shape:
 0.00000e+00 {0.00000e+00}
 0.00000e+00  0.00000e+00 

 @t = 5.0000e+00 :
 x[0] = 
   a0    =  8.08494e-01     0
   a1    = -7.58680e-03     1
   a2    =  2.54910e-03     2
   a3    = -1.57304e-05     3
   R     =  [ -2.46530e-05 :  2.46530e-05 ]
   B     =  [  8.03082e-01 :  8.18670e-01 ]

 x[1] = 
   a0    =  1.03753e+00     0
   a1    = -4.99778e-02     1
   a2    =  1.08644e-04     2
   a3    =  1.17340e-04     3
   R     =  [ -2.65044e-05 :  2.65044e-05 ]
   B     =  [  9.87519e-01 :  1.08776e+00 ]

 Rx =
center:
 0.00000e+00
 0.00000e+00
shape:
 6.07771e-10 {4.14806e-11}
 4.14806e-11  7.02486e-10 

 No STEPS  266
 No EVALATIONS   RHS: 1874   JAC: 0
 CPU TIME (SEC)     0.01465

How Do I Specify my ODE Problem in ODEBND?

Suppose we want to bound the reachable set of the following Lotka-Volterra problem:

\[\left\{\begin{array}{rl} \dot{x}_1(t)=px_1(t)\,[1-x_2(t)] & \quad \text{with}\quad x_1(0)=1.2\\ \dot{x}_2(t)=px_2(t)\,[x_1(t)-1] & \quad \text{with}\quad x_2(0)=1.1, \end{array}\right.\]

that depends on the (scalar) parameter \(p\).

The implementation starts by defining the variables and functions in the form of a DAG using the classes mc::FFGraph and mc::FFVar in MC++:

      mc::FFGraph IVP;    // DAG describing the problem

      const unsigned int NP = 1;  // Parameter dimension
      const unsigned int NX = 2;  // State dimension

      mc::FFVar P[NP];    // Parameter array
      for( unsigned int i=0; i<NP; i++ ) P[i].set( &IVP );

      mc::FFVar X[NX];    // State array
      for( unsigned int i=0; i<NX; i++ ) X[i].set( &IVP );

      mc::FFVar RHS[NX];  // Right-hand-side function
      RHS[0] = P[0] * X[0] * ( 1. - X[1] );
      RHS[1] = P[0] * X[1] * ( X[0] - 1. );

      mc::FFVar IC[NX];   // Initial-value function
      IC[0] = 1.2;
      IC[1] = 1.1;

How Do I Compute an Interval Enclosure of the Reachable Set of my ODE Problem?

The focus here is on computing an enclosure of the reachable using continuous-time set-valued integration in mc::ODEBND_GSL, but discretized set-valued integration in mc::ODEBND_VAL proceeds in the exact same way.

We start by defining an instance of mc::ODEBND_GSL and set the variables \(x\), the parameters \(p\), the right-hand-side function \(f\) and the initial-value function \(h\) for the Lotka-Volterra problem defined above:

      typedef mc::Interval I;    // Default interval type in MC++
      mc::ODEBND_GSL<I> LV;      // Continuous-time set-valued integrator

      LV.set_dag( &IVP );
      LV.set_state( NX, X );
      LV.set_parameter( NP, P );
      LV.set_dynamic( RHS );
      LV.set_initial( IC );

For simplicity, we have defined the template type T of mc::ODEBND_GSL as the default interval type mc::Interval in MC++. The other two template parameters are not specified and therefore default to the Taylor model types mc::TModel and mc::TVar here.

Now, suppose we want to compute an interval enclosure of the reachable set of this problem, with \(p\in [2.95,3.05]\) and final time \(t_{\rm f}=2\). This is done by invoking the method mc::ODEBND_GSL::bounds, after defining the parameter set \(P\) and integration timespan \([t_0,t_{\rm f}]\):

      I Ip[NP];
      Ip[0] = I(2.95,3.05);

      I Ix0[NX], Ixf[NX];
      I* Ixk[2] = {Ix0, Ixf};

      double tk[2] = {0., 5.};

      LV.bounds( 1, tk, Ip, Ixk );

Note that enclosures at intermediate times can also be computed by modifying the arrays tk and Ixk accordingly.

The following result is displayed during the computations:

 @t = 0.0000e+00 :
 x[0] = [  1.20000e+00 :  1.20000e+00 ]
 x[1] = [  1.10000e+00 :  1.10000e+00 ]
 @t = 5.0000e+00 :
 x[0] = [  7.95292e-01 :  8.21697e-01 ]
 x[1] = [  9.81849e-01 :  1.09322e+00 ]
 No STEPS  96
 No EVALATIONS   RHS: 714   JAC: 0
 CPU TIME (SEC)     0.00484

With the default options, this enclosure was generated using an ellipsoidal contractor with wrapping mitigation order of 1 and the integration method is an explicit embedded Runge-Kutta-Fehlberg (4,5) method with both relative and absolute tolerances set to 1e-7.

The default options can be modified via the public member mc::ODEBND_GSL::options; see How Are the Options Set in ODEBND? for a complete list of options. For instance, display can be suppressed and intermediate results can be recorded by setting the DISPLAY and RESRECORD options as:

      LV.options.DISPLAY = 0;
      LV.options.RESRECORD = true;

When intermediate results are recorded during the enclosure propagation, the computed enclosures at stage times can be exported to a file by using the method mc::ODEBND_GSL::record after an integration was completed:

      #include <fstream>
      std::ofstream ofile( "LV.out", std::ios_base::out );
      LV.record( ofile );

This creates (or updates) the file LV.out in the current directory as follows -- This requires the option RESRECORD to be set to `true' before calling mc::ODEBND_GSL::bounds:

   0.00000e+00   1.20000e+00   1.20000e+00   1.10000e+00   1.10000e+00
   5.00000e+00   7.95292e-01   8.21697e-01   9.81849e-01   1.09322e+00

Estimate of the Hausdorff distance between the computed enclosure and the actual reachable set - after projection onto each variable - can be computed by using the method mc::ODEBND_GSL::hausdorff as follows:

      const unsigned int NSAMP = 100;  // Number of sampling points
      double Hx0[NX], Hxf[NX];
      double* Hxk[2] = {Hx0, Hxf};
      LV.hausdorff( 1, tk, Ip, Hxk, NSAMP );

This method simply evaluates the ODEs by sampling parameter set at NSAMP equally spaced points, and keeping whichever response is minimal or maximal. Nonetheless, such a brute force approach quickly becomes computationally intractable when the ODE problem contains more than just a few parameters, especially for a large number of sample points for each parameter.

The following results are displayed in this case:

 @t = 0.0000e+00 :
 dH[0] = 0.0000e+00
 dH[1] = 0.0000e+00
 @t = 5.0000e+00 :
 dH[0] = 7.8393e-03
 dH[1] = 5.9276e-03

In particular, such estimates of the Hausdorff metric are useful to assess the convergence rate of the bounding techniques.

How Do I Compute a Polynomial Model Enclosure of the Reachable Set of my ODE Problem?

We now consider computing a 4th-order Taylor model enclosure of the reachable set of the same Lotka-Volterra system. Note that computing a Chebyshev model can be done in the exact same way, by using the template arguments mc::CModel and mc::CVar for the instantiation of mc::ODEBND_GSL.

A Taylor model of the parameter set \(P\) is defined as follows:

      typedef mc::TModel<I> TM;
      typedef mc::TVar<I>   TV;

      const unsigned int NTM = 4; // Order of Taylor model
      TM TMenv( NP, NTM );        // Taylor model environment
      TV TMp[NP];
      TMp[0] = TV( &TMenv, 0, I(2.95,3.05) );

The desired Taylor model enclosure is computed exactly as previously - the method mc::ODEBND_GSL::bounds being overloaded for polynomial model types:

      TV TMx0[NX], TMxf[NX];
      TV* TMxk[2] = {TMx0, TMxf};

      LV.bounds( 1, tk, TMp, TMxk );

The following results are obtained:

 @t = 0.0000e+00 :
 x[0] = 
   a0    =  1.20000e+00
   R     = [  0.00000e+00 :  0.00000e+00 ]
   B     =  [  1.20000e+00 :  1.20000e+00 ]

 x[1] = 
   a0    =  1.10000e+00
   R     = [  0.00000e+00 :  0.00000e+00 ]
   B     =  [  1.10000e+00 :  1.10000e+00 ]

 Rx =
center:
 0.00000e+00
 0.00000e+00
shape:
 0.00000e+00 {0.00000e+00}
 0.00000e+00  0.00000e+00 

 @t = 5.0000e+00 :
 x[0] = 
   a0    =  8.05943e-01     0
   a1    = -1.50791e-01     1
   a2    =  2.04224e+00     2
   a3    = -5.10978e-01     3
   a4    = -1.16864e+00     4
   R     =  [ -2.91408e-05 :  2.91408e-05 ]
   B     =  [  8.03060e-01 :  8.18682e-01 ]

 x[1] = 
   a0    =  1.03742e+00     0
   a1    = -1.00659e+00     1
   a2    =  9.72546e-02     2
   a3    =  3.75263e+00     3
   a4    = -4.16086e+00     4
   R     =  [ -3.09151e-05 :  3.09151e-05 ]
   B     =  [  9.86807e-01 :  1.08849e+00 ]

 Rx =
center:
 0.00000e+00
 0.00000e+00
shape:
 8.49186e-10 {6.10540e-11}
 6.10540e-11  9.55741e-10 

 No STEPS  136
 No EVALATIONS   RHS: 1024   JAC: 0
 CPU TIME (SEC)     0.00940

Estimates of the Hausdorff distance between the Taylor model remainder bound and the actual range of the remainder function - after projection onto each variable - can be computed by using the method mc::ODEBND_GSL::hausdorff as follows:

      double Hx0[NX], Hxf[NX];
      double* Hxk[2] = {Hx0, Hxf};
      LV.hausdorff( 1, tk, TMp, Hxk, NSAMP );

The following result is displayed after the computations:

 @t = 0.0000e+00 :
 dH[0] = 0.0000e+00
 dH[1] = 0.0000e+00
 @t = 5.0000e+00 :
 dH[0] = 2.9298e-05
 dH[1] = 3.1074e-05

The overestimation is much reduced compared to the interval bounding case above.

How Are the Options Set in ODEBND?

The options are defined in the structures mc::ODEBND_GSL::Options and mc::ODEBND_VAL::Options, and they can be modified through their public (non-static) members mc::ODEBND_GSL::options and mc::ODEBND_VAL::options:

      mc::ODEBND_GSL<T> LV;
      LV.options.RTOL    = LV.options.ATOL = 1e-8;
      LV.options.INTMETH = mc::ODEBND_GSL<T>::Options::RK8PD;
      LV.options.WRAPMIT = mc::ODEBND_GSL<T>::Options::DINEQ;

Full lists of options and corresponding default values are reported in the following tables.

Options in mc::ODEBND_GSL: name, type and description
Name TypeDefault Description
INTMETH mc::ODEBND_GSL::Options::INTEGRATION_METHOD mc::ODEBND_GSL::Options::RKF45 Numerical integration method
WRAPMIT mc::ODEBND_GSL::Options::WRAPPING_STRATEGY mc::ODEBND_GSL::Options::ELLIPS Wrapping mitigation strategy
ORDMIT unsigned int 1 Order of nonlinearity bounder in ellipsoidal contractor
PMOPT PMT::Options PMT::Options() Options of polynomial model arithmetic for wrapping mitigation with ellipsoidal contractor
ODESLVOPT PMT::Options mc::ODESLV_GSL::Options() Options of non-validated ODE solver for reachable set inner-approximation
USEINV bool true Use the specified invariants for bounds contraction
H0 double 1e-2 Initial step-size
HMIN double 0e0 (unrestricted) Minimal step-size
HMAX double 0e0 (unrestricted) Maximal step-size
NMAX unsigned int 0 (unlimited) Maximal number of steps per time stage
ATOL double 1e-6 Absolute tolerance for gsl_odeiv2
RTOL double 1e-6 Relative tolerance for gsl_odeiv2
QTOL double mc::machprec() Tolerance when dividing by trace of shape matrix in ellipsoidal bounds
DISPLAY int 1 Display option
RESRECORD bool false Keep track of intermediate results


Options in mc::ODEBND_VAL: name, type and description
Name TypeDefault Description
TSORDER unsigned int 7 Taylor series expansion order
WRAPMIT mc::ODEBND_VAL::Options::WRAPPING_STRATEGY mc::ODEBND_VAL::Options::ELLIPS Wrapping mitigation strategy
ORDMIT unsigned int 1 Order of nonlinearity bounder in ellipsoidal contractor
DMAX double 1e20 Maximum enclosure diameter
HMIN double 1e-8 Minimal step-size
HMAX double 1e8 Maximal step-size
HREDUC double 0.8 Reduction factor for step-size validation (between 0-1)
HSTAB bool true Consider stepsize dependence in Taylor truncation term for stabilizing properties
PMVALID bool false Use polynomial models for stepsize validation
TSTOP bool true Stop and reinitialize integrator at time steps
ODESLVOPT PMT::Options mc::ODESLV_GSL::Options() Options of non-validated ODE solver for reachable set inner-approximation
USEINV bool true Use the specified invariants for bounds contraction
TOL double 1e-7 Tolerance on truncation error for step size selection
SCALING bool true Adjust scaling for state components
ATOL double 1e-10 Relative tolerance for gsl_odeiv2
QTOL double mc::machprec() Tolerance when dividing by trace of shape matrix in ellipsoidal bounds
DISPLAY int 1 Display option
RESRECORD bool false Keep track of intermediate results

What Errors Can I Encounter with ODEBND?

Errors are managed based on the exception handling mechanism of the C++ language. Each time an error is encountered during a set-valued integration, a class object mc::ODEBND_GSL::Exceptions or mc::ODEBND_VAL::Exceptions is thrown, which contains the type and description of the error. It is the user's responsibility to test whether an exception was thrown during the integration of parametric ODEs, and make the appropriate changes accordingly. Should an exception be thrown and not caught by the calling program, the execution will stop.

Possible errors encountered during set-valued integration of parametric ODEs are reported in the following tables.

Exceptions in mc::ODEBND_GSL: number and description
Number Description
1 State reinitialization/discontinuity at intermediate stage not allowed
-33 Calling a feature not yet implemented


Exceptions in mc::ODEBND_VAL: number and description
Number Description
1 State reinitialization/discontinuity at intermediate stage not allowed
2 Invalid stepsize reduction parameter
-33 Calling a feature not yet implemented

Naturally, exceptions may be thrown by the embedded classes T, PMT or PVT themsevelves as well as by third-party programs.

References

Chachuat, B. and M.E. Villanueva, "Bounding the solutions of parametric ODEs: When Taylor models meet differential inequalities", Proceedings of the 22nd European Conference on Computer Aided Process Engineering (ESCAPE 22), I.D.L. Bogle and M. Fairweather (Eds), 2012

Harrison, G. W., "Dynamic Models with Uncertain Parameters", Proceedings of the 1st International Conference on Mathematical Modeling, X. Avula (Ed), vol. 1, pp. 295-304, 1977

Houska, B., F. Logist, J. Van Impe, and M. Diehl, "Robust optimization of nonlinear dynamic systems with application to a jacketed tubular reactor", Journal of Process Control, 22(6):1152-1160, 2012

Houska, B., M.E. Villanueva, and B. Chachuat, "A validated integration algorithm for nonlinear ODEs using Taylor models and ellipsoidal calculus", 52nd IEEE Conference on Decision and Control (CDC), 10-13 December 2013, Florence, Italy

Houska, B., M.E. Villanueva, and B. Chachuat, "Stable set-valued integration of nonlinear dynamic systems using affine set parameterizations", SIAM Journal on Numerical Analysis, submitted

Lakshmikantham, V., S. Leela, "Differential and Integral Inequalities, Theory and Applications: Volume I, Ordinary Differential Equations", Academic Press, 1969

Villanueva, M.E., B. Houska, and B. Chachuat, "Unified framework for the propagation of continuous-time enclosures for parametric nonlinear ODEs", Journal of Global Optimization, in press

Villanueva, M.E., R. Paulen, B. Houska, and B. Chachuat, "Enclosing the reachable set of parametric ODEs using Taylor models and ellipsoidal calculus", Proceedings of the 23rd European Conference on Computer Aided Process Engineering (ESCAPE 23), A. Kraslawski and I. Turunen (Eds), 2013

Villanueva, M.E., B. Houska, and B. Chachuat, "On the stability of set-valued integration for nonlinear parametric ODEs", Proceedings of the 24th European Conference on Computer Aided Process Engineering (ESCAPE 24), J.J. Klemes, P.S. Varbanov and P.Y. Liew (Eds), 2014

Walter, W., "Differential and Integral Inequalities", Springer-Verlag, Berlin, 1970