view liboctave/numeric/DASPK.cc @ 31605:e88a07dec498 stable

maint: Use macros to begin/end C++ namespaces. * oct-conf-post-public.in.h: Define two macros (OCTAVE_BEGIN_NAMESPACE, OCTAVE_END_NAMESPACE) that can be used to start/end a namespace. * mk-opts.pl, build-env.h, build-env.in.cc, __betainc__.cc, __contourc__.cc, __dsearchn__.cc, __eigs__.cc, __expint__.cc, __ftp__.cc, __gammainc__.cc, __ichol__.cc, __ilu__.cc, __isprimelarge__.cc, __lin_interpn__.cc, __magick_read__.cc, __pchip_deriv__.cc, __qp__.cc, amd.cc, auto-shlib.cc, auto-shlib.h, balance.cc, base-text-renderer.cc, base-text-renderer.h, besselj.cc, bitfcns.cc, bsxfun.cc, c-file-ptr-stream.cc, c-file-ptr-stream.h, call-stack.cc, call-stack.h, ccolamd.cc, cellfun.cc, chol.cc, colamd.cc, colloc.cc, conv2.cc, daspk.cc, dasrt.cc, dassl.cc, data.cc, data.h, debug.cc, defaults.cc, defaults.h, defun-int.h, defun.cc, det.cc, dirfns.cc, display.cc, display.h, dlmread.cc, dmperm.cc, dot.cc, dynamic-ld.cc, dynamic-ld.h, eig.cc, ellipj.cc, environment.cc, environment.h, error.cc, error.h, errwarn.h, event-manager.cc, event-manager.h, event-queue.cc, event-queue.h, fcn-info.cc, fcn-info.h, fft.cc, fft2.cc, fftn.cc, file-io.cc, filter.cc, find.cc, ft-text-renderer.cc, ft-text-renderer.h, gcd.cc, getgrent.cc, getpwent.cc, getrusage.cc, givens.cc, gl-render.cc, gl-render.h, gl2ps-print.cc, gl2ps-print.h, graphics-toolkit.cc, graphics-toolkit.h, graphics.cc, graphics.in.h, gsvd.cc, gtk-manager.cc, gtk-manager.h, hash.cc, help.cc, help.h, hess.cc, hex2num.cc, hook-fcn.cc, hook-fcn.h, input.cc, input.h, interpreter-private.cc, interpreter-private.h, interpreter.cc, interpreter.h, inv.cc, jsondecode.cc, jsonencode.cc, kron.cc, latex-text-renderer.cc, latex-text-renderer.h, load-path.cc, load-path.h, load-save.cc, load-save.h, lookup.cc, ls-ascii-helper.cc, ls-ascii-helper.h, ls-oct-text.cc, ls-utils.cc, ls-utils.h, lsode.cc, lu.cc, mappers.cc, matrix_type.cc, max.cc, mex-private.h, mex.cc, mgorth.cc, nproc.cc, oct-fstrm.cc, oct-fstrm.h, oct-hdf5-types.cc, oct-hdf5-types.h, oct-hist.cc, oct-hist.h, oct-iostrm.cc, oct-iostrm.h, oct-opengl.h, oct-prcstrm.cc, oct-prcstrm.h, oct-procbuf.cc, oct-procbuf.h, oct-process.cc, oct-process.h, oct-stdstrm.h, oct-stream.cc, oct-stream.h, oct-strstrm.cc, oct-strstrm.h, oct-tex-lexer.in.ll, oct-tex-parser.yy, ordqz.cc, ordschur.cc, pager.cc, pager.h, pinv.cc, pow2.cc, pr-flt-fmt.cc, pr-output.cc, procstream.cc, procstream.h, psi.cc, qr.cc, quad.cc, quadcc.cc, qz.cc, rand.cc, rcond.cc, regexp.cc, schur.cc, settings.cc, settings.h, sighandlers.cc, sighandlers.h, sparse-xdiv.cc, sparse-xdiv.h, sparse-xpow.cc, sparse-xpow.h, sparse.cc, spparms.cc, sqrtm.cc, stack-frame.cc, stack-frame.h, stream-euler.cc, strfind.cc, strfns.cc, sub2ind.cc, svd.cc, sylvester.cc, symbfact.cc, syminfo.cc, syminfo.h, symrcm.cc, symrec.cc, symrec.h, symscope.cc, symscope.h, symtab.cc, symtab.h, syscalls.cc, sysdep.cc, sysdep.h, text-engine.cc, text-engine.h, text-renderer.cc, text-renderer.h, time.cc, toplev.cc, tril.cc, tsearch.cc, typecast.cc, url-handle-manager.cc, url-handle-manager.h, urlwrite.cc, utils.cc, utils.h, variables.cc, variables.h, xdiv.cc, xdiv.h, xnorm.cc, xnorm.h, xpow.cc, xpow.h, __delaunayn__.cc, __fltk_uigetfile__.cc, __glpk__.cc, __init_fltk__.cc, __init_gnuplot__.cc, __ode15__.cc, __voronoi__.cc, audiodevinfo.cc, audioread.cc, convhulln.cc, fftw.cc, gzip.cc, mk-build-env-features.sh, mk-builtins.pl, cdef-class.cc, cdef-class.h, cdef-fwd.h, cdef-manager.cc, cdef-manager.h, cdef-method.cc, cdef-method.h, cdef-object.cc, cdef-object.h, cdef-package.cc, cdef-package.h, cdef-property.cc, cdef-property.h, cdef-utils.cc, cdef-utils.h, ov-base.cc, ov-base.h, ov-bool-mat.cc, ov-builtin.h, ov-cell.cc, ov-class.cc, ov-class.h, ov-classdef.cc, ov-classdef.h, ov-complex.cc, ov-fcn-handle.cc, ov-fcn-handle.h, ov-fcn.h, ov-java.cc, ov-java.h, ov-mex-fcn.h, ov-null-mat.cc, ov-oncleanup.cc, ov-struct.cc, ov-typeinfo.cc, ov-typeinfo.h, ov-usr-fcn.cc, ov-usr-fcn.h, ov.cc, ov.h, octave.cc, octave.h, mk-ops.sh, op-b-b.cc, op-b-bm.cc, op-b-sbm.cc, op-bm-b.cc, op-bm-bm.cc, op-bm-sbm.cc, op-cdm-cdm.cc, op-cell.cc, op-chm.cc, op-class.cc, op-cm-cm.cc, op-cm-cs.cc, op-cm-m.cc, op-cm-s.cc, op-cm-scm.cc, op-cm-sm.cc, op-cs-cm.cc, op-cs-cs.cc, op-cs-m.cc, op-cs-s.cc, op-cs-scm.cc, op-cs-sm.cc, op-dm-dm.cc, op-dm-scm.cc, op-dm-sm.cc, op-dm-template.cc, op-dms-template.cc, op-fcdm-fcdm.cc, op-fcm-fcm.cc, op-fcm-fcs.cc, op-fcm-fm.cc, op-fcm-fs.cc, op-fcn.cc, op-fcs-fcm.cc, op-fcs-fcs.cc, op-fcs-fm.cc, op-fcs-fs.cc, op-fdm-fdm.cc, op-fm-fcm.cc, op-fm-fcs.cc, op-fm-fm.cc, op-fm-fs.cc, op-fs-fcm.cc, op-fs-fcs.cc, op-fs-fm.cc, op-fs-fs.cc, op-i16-i16.cc, op-i32-i32.cc, op-i64-i64.cc, op-i8-i8.cc, op-int-concat.cc, op-m-cm.cc, op-m-cs.cc, op-m-m.cc, op-m-s.cc, op-m-scm.cc, op-m-sm.cc, op-mi.cc, op-pm-pm.cc, op-pm-scm.cc, op-pm-sm.cc, op-pm-template.cc, op-range.cc, op-s-cm.cc, op-s-cs.cc, op-s-m.cc, op-s-s.cc, op-s-scm.cc, op-s-sm.cc, op-sbm-b.cc, op-sbm-bm.cc, op-sbm-sbm.cc, op-scm-cm.cc, op-scm-cs.cc, op-scm-m.cc, op-scm-s.cc, op-scm-scm.cc, op-scm-sm.cc, op-sm-cm.cc, op-sm-cs.cc, op-sm-m.cc, op-sm-s.cc, op-sm-scm.cc, op-sm-sm.cc, op-str-m.cc, op-str-s.cc, op-str-str.cc, op-struct.cc, op-ui16-ui16.cc, op-ui32-ui32.cc, op-ui64-ui64.cc, op-ui8-ui8.cc, ops.h, anon-fcn-validator.cc, anon-fcn-validator.h, bp-table.cc, bp-table.h, comment-list.cc, comment-list.h, filepos.h, lex.h, lex.ll, oct-lvalue.cc, oct-lvalue.h, oct-parse.yy, parse.h, profiler.cc, profiler.h, pt-anon-scopes.cc, pt-anon-scopes.h, pt-arg-list.cc, pt-arg-list.h, pt-args-block.cc, pt-args-block.h, pt-array-list.cc, pt-array-list.h, pt-assign.cc, pt-assign.h, pt-binop.cc, pt-binop.h, pt-bp.cc, pt-bp.h, pt-cbinop.cc, pt-cbinop.h, pt-cell.cc, pt-cell.h, pt-check.cc, pt-check.h, pt-classdef.cc, pt-classdef.h, pt-cmd.h, pt-colon.cc, pt-colon.h, pt-const.cc, pt-const.h, pt-decl.cc, pt-decl.h, pt-eval.cc, pt-eval.h, pt-except.cc, pt-except.h, pt-exp.cc, pt-exp.h, pt-fcn-handle.cc, pt-fcn-handle.h, pt-id.cc, pt-id.h, pt-idx.cc, pt-idx.h, pt-jump.h, pt-loop.cc, pt-loop.h, pt-mat.cc, pt-mat.h, pt-misc.cc, pt-misc.h, pt-pr-code.cc, pt-pr-code.h, pt-select.cc, pt-select.h, pt-spmd.cc, pt-spmd.h, pt-stmt.cc, pt-stmt.h, pt-tm-const.cc, pt-tm-const.h, pt-unop.cc, pt-unop.h, pt-vm-eval.cc, pt-walk.cc, pt-walk.h, pt.cc, pt.h, token.cc, token.h, Range.cc, Range.h, idx-vector.cc, idx-vector.h, range-fwd.h, CollocWt.cc, CollocWt.h, aepbalance.cc, aepbalance.h, chol.cc, chol.h, gepbalance.cc, gepbalance.h, gsvd.cc, gsvd.h, hess.cc, hess.h, lo-mappers.cc, lo-mappers.h, lo-specfun.cc, lo-specfun.h, lu.cc, lu.h, oct-convn.cc, oct-convn.h, oct-fftw.cc, oct-fftw.h, oct-norm.cc, oct-norm.h, oct-rand.cc, oct-rand.h, oct-spparms.cc, oct-spparms.h, qr.cc, qr.h, qrp.cc, qrp.h, randgamma.cc, randgamma.h, randmtzig.cc, randmtzig.h, randpoisson.cc, randpoisson.h, schur.cc, schur.h, sparse-chol.cc, sparse-chol.h, sparse-lu.cc, sparse-lu.h, sparse-qr.cc, sparse-qr.h, svd.cc, svd.h, child-list.cc, child-list.h, dir-ops.cc, dir-ops.h, file-ops.cc, file-ops.h, file-stat.cc, file-stat.h, lo-sysdep.cc, lo-sysdep.h, lo-sysinfo.cc, lo-sysinfo.h, mach-info.cc, mach-info.h, oct-env.cc, oct-env.h, oct-group.cc, oct-group.h, oct-password.cc, oct-password.h, oct-syscalls.cc, oct-syscalls.h, oct-time.cc, oct-time.h, oct-uname.cc, oct-uname.h, action-container.cc, action-container.h, base-list.h, cmd-edit.cc, cmd-edit.h, cmd-hist.cc, cmd-hist.h, f77-fcn.h, file-info.cc, file-info.h, lo-array-errwarn.cc, lo-array-errwarn.h, lo-hash.cc, lo-hash.h, lo-ieee.h, lo-regexp.cc, lo-regexp.h, lo-utils.cc, lo-utils.h, oct-base64.cc, oct-base64.h, oct-glob.cc, oct-glob.h, oct-inttypes.h, oct-mutex.cc, oct-mutex.h, oct-refcount.h, oct-shlib.cc, oct-shlib.h, oct-sparse.cc, oct-sparse.h, oct-string.h, octave-preserve-stream-state.h, pathsearch.cc, pathsearch.h, quit.cc, quit.h, unwind-prot.cc, unwind-prot.h, url-transfer.cc, url-transfer.h : Use new macros to begin/end C++ namespaces.
author Rik <rik@octave.org>
date Thu, 01 Dec 2022 14:23:45 -0800
parents 51a3d3a69193
children 597f3ee61a48
line wrap: on
line source

////////////////////////////////////////////////////////////////////////
//
// Copyright (C) 1996-2022 The Octave Project Developers
//
// See the file COPYRIGHT.md in the top-level directory of this
// distribution or <https://octave.org/copyright/>.
//
// This file is part of Octave.
//
// Octave is free software: you can redistribute it and/or modify it
// under the terms of the GNU General Public License as published by
// the Free Software Foundation, either version 3 of the License, or
// (at your option) any later version.
//
// Octave is distributed in the hope that it will be useful, but
// WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
// GNU General Public License for more details.
//
// You should have received a copy of the GNU General Public License
// along with Octave; see the file COPYING.  If not, see
// <https://www.gnu.org/licenses/>.
//
////////////////////////////////////////////////////////////////////////

#if defined (HAVE_CONFIG_H)
#  include "config.h"
#endif

#include <cinttypes>
#include <sstream>

#include "DASPK.h"
#include "dMatrix.h"
#include "f77-fcn.h"
#include "lo-error.h"
#include "quit.h"

typedef F77_INT (*daspk_fcn_ptr) (const double&, const double *, const double *,
                                  const double&, double *, F77_INT&, double *,
                                  F77_INT *);

typedef F77_INT (*daspk_jac_ptr) (const double&, const double *, const double *,
                                  double *, const double&, double *, F77_INT *);

typedef F77_INT (*daspk_psol_ptr) (const F77_INT&, const double&,
                                   const double *, const double *,
                                   const double *, const double&,
                                   const double *, double *, F77_INT *,
                                   double *, const double&, F77_INT&,
                                   double *, F77_INT *);

extern "C"
{
  F77_RET_T
  F77_FUNC (ddaspk, DDASPK) (daspk_fcn_ptr, const F77_INT&, F77_DBLE&,
                             F77_DBLE *, F77_DBLE *, F77_DBLE&, const F77_INT *,
                             const F77_DBLE *, const F77_DBLE *, F77_INT&,
                             F77_DBLE *, const F77_INT&, F77_INT *,
                             const F77_INT&, const F77_DBLE *, const F77_INT *,
                             daspk_jac_ptr, daspk_psol_ptr);
}

static DAEFunc::DAERHSFunc user_fcn;
static DAEFunc::DAEJacFunc user_jac;
static F77_INT nn;

static F77_INT
ddaspk_f (const double& time, const double *state, const double *deriv,
          const double&, double *delta, F77_INT& ires, double *, F77_INT *)
{
  ColumnVector tmp_deriv (nn);
  ColumnVector tmp_state (nn);
  ColumnVector tmp_delta (nn);

  for (F77_INT i = 0; i < nn; i++)
    {
      tmp_deriv.elem (i) = deriv[i];
      tmp_state.elem (i) = state[i];
    }

  octave_idx_type tmp_ires = ires;

  tmp_delta = user_fcn (tmp_state, tmp_deriv, time, tmp_ires);

  ires = octave::to_f77_int (tmp_ires);

  if (ires >= 0)
    {
      if (tmp_delta.isempty ())
        ires = -2;
      else
        {
          for (F77_INT i = 0; i < nn; i++)
            delta[i] = tmp_delta.elem (i);
        }
    }

  return 0;
}

//NEQ, T, Y, YPRIME, SAVR, WK, CJ, WGHT,
//C                          WP, IWP, B, EPLIN, IER, RPAR, IPAR)

static F77_INT
ddaspk_psol (const F77_INT&, const double&, const double *,
             const double *, const double *, const double&,
             const double *, double *, F77_INT *, double *,
             const double&, F77_INT&, double *, F77_INT *)
{
  (*current_liboctave_error_handler) ("daspk: PSOL is not implemented");

  return 0;
}

static F77_INT
ddaspk_j (const double& time, const double *state, const double *deriv,
          double *pd, const double& cj, double *, F77_INT *)
{
  // FIXME: would be nice to avoid copying the data.

  ColumnVector tmp_state (nn);
  ColumnVector tmp_deriv (nn);

  for (F77_INT i = 0; i < nn; i++)
    {
      tmp_deriv.elem (i) = deriv[i];
      tmp_state.elem (i) = state[i];
    }

  Matrix tmp_pd = user_jac (tmp_state, tmp_deriv, time, cj);

  for (F77_INT j = 0; j < nn; j++)
    for (F77_INT i = 0; i < nn; i++)
      pd[nn * j + i] = tmp_pd.elem (i, j);

  return 0;
}

ColumnVector
DASPK::do_integrate (double tout)
{
  // FIXME: should handle all this option stuff just once for each new problem.

  ColumnVector retval;

  if (! m_initialized || m_restart || DAEFunc::m_reset
      || DASPK_options::m_reset)
    {
      m_integration_error = false;

      m_initialized = true;

      m_info.resize (dim_vector (20, 1));

      for (F77_INT i = 0; i < 20; i++)
        m_info(i) = 0;

      F77_INT n = octave::to_f77_int (size ());

      nn = n;

      m_info(0) = 0;

      if (m_stop_time_set)
        {
          m_rwork(0) = m_stop_time;
          m_info(3) = 1;
        }
      else
        m_info(3) = 0;

      // DAEFunc

      user_fcn = DAEFunc::function ();
      user_jac = DAEFunc::jacobian_function ();

      if (user_fcn)
        {
          octave_idx_type ires = 0;

          ColumnVector res = (*user_fcn) (m_x, m_xdot, m_t, ires);

          if (res.numel () != m_x.numel ())
            {
              // FIXME: Should this be a warning?
              (*current_liboctave_error_handler)
                ("daspk: inconsistent sizes for state and residual vectors");

              m_integration_error = true;
              return retval;
            }
        }
      else
        {
          // FIXME: Should this be a warning?
          (*current_liboctave_error_handler)
            ("daspk: no user supplied RHS subroutine!");

          m_integration_error = true;
          return retval;
        }

      m_info(4) = (user_jac ? 1 : 0);

      DAEFunc::m_reset = false;

      octave_idx_type eiq = enforce_inequality_constraints ();
      octave_idx_type ccic = compute_consistent_initial_condition ();
      octave_idx_type eavfet = exclude_algebraic_variables_from_error_test ();

      m_liw = 40 + n;
      if (eiq == 1 || eiq == 3)
        m_liw += n;
      if (ccic == 1 || eavfet == 1)
        m_liw += n;

      m_lrw = 50 + 9*n + n*n;
      if (eavfet == 1)
        m_lrw += n;

      m_iwork.resize (dim_vector (m_liw, 1));
      m_rwork.resize (dim_vector (m_lrw, 1));

      // DASPK_options

      m_abs_tol = absolute_tolerance ();
      m_rel_tol = relative_tolerance ();

      F77_INT abs_tol_len = octave::to_f77_int (m_abs_tol.numel ());
      F77_INT rel_tol_len = octave::to_f77_int (m_rel_tol.numel ());

      if (abs_tol_len == 1 && rel_tol_len == 1)
        {
          m_info(1) = 0;
        }
      else if (abs_tol_len == n && rel_tol_len == n)
        {
          m_info(1) = 1;
        }
      else
        {
          // FIXME: Should this be a warning?
          (*current_liboctave_error_handler)
            ("daspk: inconsistent sizes for tolerance arrays");

          m_integration_error = true;
          return retval;
        }

      double hmax = maximum_step_size ();
      if (hmax >= 0.0)
        {
          m_rwork(1) = hmax;
          m_info(6) = 1;
        }
      else
        m_info(6) = 0;

      double h0 = initial_step_size ();
      if (h0 >= 0.0)
        {
          m_rwork(2) = h0;
          m_info(7) = 1;
        }
      else
        m_info(7) = 0;

      octave_idx_type maxord = maximum_order ();
      if (maxord >= 0)
        {
          if (maxord > 0 && maxord < 6)
            {
              m_info(8) = 1;
              m_iwork(2) = octave::to_f77_int (maxord);
            }
          else
            {
              // FIXME: Should this be a warning?
              (*current_liboctave_error_handler)
                ("daspk: invalid value for maximum order");
              m_integration_error = true;
              return retval;
            }
        }

      switch (eiq)
        {
        case 1:
        case 3:
          {
            Array<octave_idx_type> ict = inequality_constraint_types ();

            F77_INT ict_nel = octave::to_f77_int (ict.numel ());

            if (ict_nel == n)
              {
                for (F77_INT i = 0; i < n; i++)
                  {
                    F77_INT val = octave::to_f77_int (ict(i));
                    if (val < -2 || val > 2)
                      {
                        // FIXME: Should this be a warning?
                        (*current_liboctave_error_handler)
                          ("daspk: invalid value for inequality constraint type");
                        m_integration_error = true;
                        return retval;
                      }
                    m_iwork(40+i) = val;
                  }
              }
            else
              {
                // FIXME: Should this be a warning?
                (*current_liboctave_error_handler)
                  ("daspk: inequality constraint types size mismatch");
                m_integration_error = true;
                return retval;
              }
          }
          OCTAVE_FALLTHROUGH;

        case 0:
        case 2:
          m_info(9) = octave::to_f77_int (eiq);
          break;

        default:
          // FIXME: Should this be a warning?
          (*current_liboctave_error_handler)
            ("daspk: invalid value for enforce inequality constraints option");
          m_integration_error = true;
          return retval;
        }

      if (ccic)
        {
          if (ccic == 1)
            {
              // FIXME: this code is duplicated below.

              Array<octave_idx_type> av = algebraic_variables ();

              F77_INT av_nel = octave::to_f77_int (av.numel ());

              if (av_nel == n)
                {
                  F77_INT lid;
                  if (eiq == 0 || eiq == 2)
                    lid = 40;
                  else if (eiq == 1 || eiq == 3)
                    lid = 40 + n;
                  else
                    (*current_liboctave_error_handler)
                      ("daspk: invalid value for eiq: "
                       "%" OCTAVE_IDX_TYPE_FORMAT, eiq);

                  for (F77_INT i = 0; i < n; i++)
                    m_iwork(lid+i) = (av(i) ? -1 : 1);
                }
              else
                {
                  // FIXME: Should this be a warning?
                  (*current_liboctave_error_handler)
                    ("daspk: algebraic variables size mismatch");
                  m_integration_error = true;
                  return retval;
                }
            }
          else if (ccic != 2)
            {
              // FIXME: Should this be a warning?
              (*current_liboctave_error_handler)
                ("daspk: invalid value for compute consistent initial condition option");
              m_integration_error = true;
              return retval;
            }

          m_info(10) = octave::to_f77_int (ccic);
        }

      if (eavfet)
        {
          m_info(15) = 1;

          // FIXME: this code is duplicated above.

          Array<octave_idx_type> av = algebraic_variables ();

          F77_INT av_nel = octave::to_f77_int (av.numel ());

          if (av_nel == n)
            {
              F77_INT lid;
              if (eiq == 0 || eiq == 2)
                lid = 40;
              else if (eiq == 1 || eiq == 3)
                lid = 40 + n;
              else
                (*current_liboctave_error_handler)
                  ("daspk: invalid value for eiq: %" OCTAVE_IDX_TYPE_FORMAT,
                   eiq);

              for (F77_INT i = 0; i < n; i++)
                m_iwork(lid+i) = (av(i) ? -1 : 1);
            }
        }

      if (use_initial_condition_heuristics ())
        {
          Array<double> ich = initial_condition_heuristics ();

          if (ich.numel () == 6)
            {
              m_iwork(31) = octave::to_f77_int (octave::math::nint_big (ich(0)));
              m_iwork(32) = octave::to_f77_int (octave::math::nint_big (ich(1)));
              m_iwork(33) = octave::to_f77_int (octave::math::nint_big (ich(2)));
              m_iwork(34) = octave::to_f77_int (octave::math::nint_big (ich(3)));

              m_rwork(13) = ich(4);
              m_rwork(14) = ich(5);
            }
          else
            {
              // FIXME: Should this be a warning?
              (*current_liboctave_error_handler)
                ("daspk: invalid initial condition heuristics option");
              m_integration_error = true;
              return retval;
            }

          m_info(16) = 1;
        }

      octave_idx_type pici = print_initial_condition_info ();
      switch (pici)
        {
        case 0:
        case 1:
        case 2:
          m_info(17) = octave::to_f77_int (pici);
          break;

        default:
          // FIXME: Should this be a warning?
          (*current_liboctave_error_handler)
            ("daspk: invalid value for print initial condition info option");
          m_integration_error = true;
          return retval;
          break;
        }

      DASPK_options::m_reset = false;

      m_restart = false;
    }

  double *px = m_x.fortran_vec ();
  double *pxdot = m_xdot.fortran_vec ();

  F77_INT *pinfo = m_info.fortran_vec ();

  double *prel_tol = m_rel_tol.fortran_vec ();
  double *pabs_tol = m_abs_tol.fortran_vec ();

  double *prwork = m_rwork.fortran_vec ();
  F77_INT *piwork = m_iwork.fortran_vec ();

  double *dummy = nullptr;
  F77_INT *idummy = nullptr;

  F77_INT tmp_istate = octave::to_f77_int (m_istate);

  F77_XFCN (ddaspk, DDASPK, (ddaspk_f, nn, m_t, px, pxdot, tout, pinfo,
                             prel_tol, pabs_tol, tmp_istate, prwork, m_lrw,
                             piwork, m_liw, dummy, idummy, ddaspk_j,
                             ddaspk_psol));

  m_istate = tmp_istate;

  switch (m_istate)
    {
    case 1: // A step was successfully taken in intermediate-output
            // mode.  The code has not yet reached TOUT.
    case 2: // The integration to TSTOP was successfully completed
            // (T=TSTOP) by stepping exactly to TSTOP.
    case 3: // The integration to TOUT was successfully completed
            // (T=TOUT) by stepping past TOUT.  Y(*) is obtained by
            // interpolation.  YPRIME(*) is obtained by interpolation.
    case 4: // The initial condition calculation, with
            // INFO(11) > 0, was successful, and INFO(14) = 1.
            // No integration steps were taken, and the solution
            // is not considered to have been started.
      retval = m_x;
      m_t = tout;
      break;

    case -1: // A large amount of work has been expended.  (~500 steps).
    case -2: // The error tolerances are too stringent.
    case -3: // The local error test cannot be satisfied because you
             // specified a zero component in ATOL and the
             // corresponding computed solution component is zero.
             // Thus, a pure relative error test is impossible for
             // this component.
    case -6: // DDASPK had repeated error test failures on the last
             // attempted step.
    case -7: // The corrector could not converge.
    case -8: // The matrix of partial derivatives is singular.
    case -9: // The corrector could not converge.  There were repeated
             // error test failures in this step.
    case -10: // The corrector could not converge because IRES was
              // equal to minus one.
    case -11: // IRES equal to -2 was encountered and control is being
              // returned to the calling program.
    case -12: // DDASPK failed to compute the initial YPRIME.
    case -13: // Unrecoverable error encountered inside user's
              // PSOL routine, and control is being returned to
              // the calling program.
    case -14: // The Krylov linear system solver could not
              // achieve convergence.
    case -33: // The code has encountered trouble from which it cannot
              // recover.  A message is printed explaining the trouble
              // and control is returned to the calling program.  For
              // example, this occurs when invalid input is detected.
      m_integration_error = true;
      break;

    default:
      m_integration_error = true;
      (*current_liboctave_error_handler)
        ("unrecognized value of m_istate (= %" OCTAVE_IDX_TYPE_FORMAT ") "
         "returned from ddaspk", m_istate);
      break;
    }

  return retval;
}

Matrix
DASPK::do_integrate (const ColumnVector& tout)
{
  Matrix dummy;
  return integrate (tout, dummy);
}

Matrix
DASPK::integrate (const ColumnVector& tout, Matrix& xdot_out)
{
  Matrix retval;

  octave_idx_type n_out = tout.numel ();
  F77_INT n = octave::to_f77_int (size ());

  if (n_out > 0 && n > 0)
    {
      retval.resize (n_out, n);
      xdot_out.resize (n_out, n);

      for (F77_INT i = 0; i < n; i++)
        {
          retval.elem (0, i) = m_x.elem (i);
          xdot_out.elem (0, i) = m_xdot.elem (i);
        }

      for (octave_idx_type j = 1; j < n_out; j++)
        {
          ColumnVector x_next = do_integrate (tout.elem (j));

          if (m_integration_error)
            return retval;

          for (F77_INT i = 0; i < n; i++)
            {
              retval.elem (j, i) = x_next.elem (i);
              xdot_out.elem (j, i) = m_xdot.elem (i);
            }
        }
    }

  return retval;
}

Matrix
DASPK::do_integrate (const ColumnVector& tout, const ColumnVector& tcrit)
{
  Matrix dummy;
  return integrate (tout, dummy, tcrit);
}

Matrix
DASPK::integrate (const ColumnVector& tout, Matrix& xdot_out,
                  const ColumnVector& tcrit)
{
  Matrix retval;

  octave_idx_type n_out = tout.numel ();
  F77_INT n = octave::to_f77_int (size ());

  if (n_out > 0 && n > 0)
    {
      retval.resize (n_out, n);
      xdot_out.resize (n_out, n);

      for (F77_INT i = 0; i < n; i++)
        {
          retval.elem (0, i) = m_x.elem (i);
          xdot_out.elem (0, i) = m_xdot.elem (i);
        }

      octave_idx_type n_crit = tcrit.numel ();

      if (n_crit > 0)
        {
          octave_idx_type i_crit = 0;
          octave_idx_type i_out = 1;
          double next_crit = tcrit.elem (0);
          double next_out;
          while (i_out < n_out)
            {
              bool do_restart = false;

              next_out = tout.elem (i_out);
              if (i_crit < n_crit)
                next_crit = tcrit.elem (i_crit);

              bool save_output;
              double t_out;

              if (next_crit == next_out)
                {
                  set_stop_time (next_crit);
                  t_out = next_out;
                  save_output = true;
                  i_out++;
                  i_crit++;
                  do_restart = true;
                }
              else if (next_crit < next_out)
                {
                  if (i_crit < n_crit)
                    {
                      set_stop_time (next_crit);
                      t_out = next_crit;
                      save_output = false;
                      i_crit++;
                      do_restart = true;
                    }
                  else
                    {
                      clear_stop_time ();
                      t_out = next_out;
                      save_output = true;
                      i_out++;
                    }
                }
              else
                {
                  set_stop_time (next_crit);
                  t_out = next_out;
                  save_output = true;
                  i_out++;
                }

              ColumnVector x_next = do_integrate (t_out);

              if (m_integration_error)
                return retval;

              if (save_output)
                {
                  for (F77_INT i = 0; i < n; i++)
                    {
                      retval.elem (i_out-1, i) = x_next.elem (i);
                      xdot_out.elem (i_out-1, i) = m_xdot.elem (i);
                    }
                }

              if (do_restart)
                force_restart ();
            }
        }
      else
        {
          retval = integrate (tout, xdot_out);

          if (m_integration_error)
            return retval;
        }
    }

  return retval;
}

std::string
DASPK::error_message (void) const
{
  std::string retval;

  std::ostringstream buf;
  buf << m_t;
  std::string t_curr = buf.str ();

  switch (m_istate)
    {
    case 1:
      retval = "a step was successfully taken in intermediate-output mode.";
      break;

    case 2:
      retval = "integration completed by stepping exactly to TOUT";
      break;

    case 3:
      retval = "integration to tout completed by stepping past TOUT";
      break;

    case 4:
      retval = "initial condition calculation completed successfully";
      break;

    case -1:
      retval = "a large amount of work has been expended (t =" + t_curr + ')';
      break;

    case -2:
      retval = "the error tolerances are too stringent";
      break;

    case -3:
      retval = "error weight became zero during problem. (t = " + t_curr +
               "; solution component i vanished, and atol or atol(i) == 0)";
      break;

    case -6:
      retval = "repeated error test failures on the last attempted step (t = "
               + t_curr + ')';
      break;

    case -7:
      retval = "the corrector could not converge (t = " + t_curr + ')';
      break;

    case -8:
      retval = "the matrix of partial derivatives is singular (t = " + t_curr +
               ')';
      break;

    case -9:
      retval = "the corrector could not converge (t = " + t_curr +
               "; repeated test failures)";
      break;

    case -10:
      retval = "corrector could not converge because IRES was -1 (t = "
               + t_curr + ')';
      break;

    case -11:
      retval = "return requested in user-supplied function (t = " + t_curr +
               ')';
      break;

    case -12:
      retval = "failed to compute consistent initial conditions";
      break;

    case -13:
      retval = "unrecoverable error encountered inside user's PSOL function (t = "
               + t_curr + ')';
      break;

    case -14:
      retval = "the Krylov linear system solver failed to converge (t = "
               + t_curr + ')';
      break;

    case -33:
      retval = "unrecoverable error (see printed message)";
      break;

    default:
      retval = "unknown error state";
      break;
    }

  return retval;
}