view libinterp/corefcn/dasrt.cc @ 21200:fcac5dbbf9ed

maint: Indent #ifdef blocks in libinterp. * builtins.h, Cell.cc, __contourc__.cc, __dispatch__.cc, __dsearchn__.cc, __ichol__.cc, __ilu__.cc, __lin_interpn__.cc, __pchip_deriv__.cc, __qp__.cc, balance.cc, besselj.cc, betainc.cc, bitfcns.cc, bsxfun.cc, c-file-ptr-stream.cc, c-file-ptr-stream.h, cellfun.cc, colloc.cc, comment-list.cc, conv2.cc, daspk.cc, dasrt.cc, dassl.cc, data.cc, debug.cc, defaults.cc, defaults.in.h, defun-dld.h, defun.cc, defun.h, det.cc, dirfns.cc, display.cc, dlmread.cc, dot.cc, dynamic-ld.cc, eig.cc, ellipj.cc, error.cc, errwarn.cc, event-queue.cc, fft.cc, fft2.cc, fftn.cc, file-io.cc, filter.cc, find.cc, gammainc.cc, gcd.cc, getgrent.cc, getpwent.cc, getrusage.cc, givens.cc, gl-render.cc, gl2ps-print.cc, graphics.cc, graphics.in.h, gripes.cc, hash.cc, help.cc, hess.cc, hex2num.cc, input.cc, inv.cc, jit-ir.cc, jit-typeinfo.cc, jit-util.cc, jit-util.h, kron.cc, load-path.cc, load-save.cc, lookup.cc, ls-ascii-helper.cc, ls-hdf5.cc, ls-mat-ascii.cc, ls-mat4.cc, ls-mat5.cc, ls-oct-binary.cc, ls-oct-text.cc, ls-oct-text.h, ls-utils.cc, ls-utils.h, lsode.cc, lu.cc, luinc.cc, mappers.cc, matrix_type.cc, max.cc, mex.h, mexproto.h, mgorth.cc, nproc.cc, oct-errno.in.cc, oct-fstrm.cc, oct-hdf5-types.cc, oct-hdf5.h, oct-hist.cc, oct-iostrm.cc, oct-lvalue.cc, oct-map.cc, oct-prcstrm.cc, oct-procbuf.cc, oct-stream.cc, oct-strstrm.cc, octave-link.cc, ordschur.cc, pager.cc, pinv.cc, pr-output.cc, procstream.cc, profiler.cc, psi.cc, pt-jit.cc, quad.cc, quadcc.cc, qz.cc, rand.cc, rcond.cc, regexp.cc, schur.cc, sighandlers.cc, sparse-xdiv.cc, sparse-xpow.cc, sparse.cc, spparms.cc, sqrtm.cc, str2double.cc, strfind.cc, strfns.cc, sub2ind.cc, svd.cc, sylvester.cc, symtab.cc, syscalls.cc, sysdep.cc, sysdep.h, time.cc, toplev.cc, tril.cc, tsearch.cc, txt-eng-ft.cc, txt-eng.cc, typecast.cc, urlwrite.cc, utils.cc, variables.cc, xdiv.cc, xnorm.cc, xpow.cc, zfstream.cc, __delaunayn__.cc, __eigs__.cc, __fltk_uigetfile__.cc, __glpk__.cc, __init_fltk__.cc, __init_gnuplot__.cc, __magick_read__.cc, __osmesa_print__.cc, __voronoi__.cc, amd.cc, audiodevinfo.cc, audioread.cc, ccolamd.cc, chol.cc, colamd.cc, convhulln.cc, dmperm.cc, fftw.cc, oct-qhull.h, qr.cc, symbfact.cc, symrcm.cc, oct-conf.in.cc, ov-base-diag.cc, ov-base-int.cc, ov-base-mat.cc, ov-base-scalar.cc, ov-base-sparse.cc, ov-base.cc, ov-bool-mat.cc, ov-bool-sparse.cc, ov-bool.cc, ov-builtin.cc, ov-cell.cc, ov-ch-mat.cc, ov-class.cc, ov-classdef.cc, ov-colon.cc, ov-complex.cc, ov-cs-list.cc, ov-cx-diag.cc, ov-cx-mat.cc, ov-cx-sparse.cc, ov-dld-fcn.cc, ov-fcn-handle.cc, ov-fcn-inline.cc, ov-fcn.cc, ov-float.cc, ov-flt-complex.cc, ov-flt-cx-diag.cc, ov-flt-cx-mat.cc, ov-flt-re-diag.cc, ov-flt-re-mat.cc, ov-int16.cc, ov-int32.cc, ov-int64.cc, ov-int8.cc, ov-java.cc, ov-lazy-idx.cc, ov-mex-fcn.cc, ov-null-mat.cc, ov-oncleanup.cc, ov-perm.cc, ov-range.cc, ov-re-diag.cc, ov-re-mat.cc, ov-re-sparse.cc, ov-scalar.cc, ov-str-mat.cc, ov-struct.cc, ov-typeinfo.cc, ov-uint16.cc, ov-uint32.cc, ov-uint64.cc, ov-uint8.cc, ov-usr-fcn.cc, ov.cc, ovl.cc, octave.cc, 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-double-conv.cc, op-fcdm-fcdm.cc, op-fcdm-fdm.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-float-conv.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-int-conv.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-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, pt-arg-list.cc, pt-array-list.cc, pt-assign.cc, pt-binop.cc, pt-bp.cc, pt-cbinop.cc, pt-cell.cc, pt-check.cc, pt-classdef.cc, pt-cmd.cc, pt-colon.cc, pt-colon.h, pt-const.cc, pt-decl.cc, pt-eval.cc, pt-except.cc, pt-exp.cc, pt-fcn-handle.cc, pt-funcall.cc, pt-id.cc, pt-idx.cc, pt-jump.cc, pt-loop.cc, pt-mat.cc, pt-misc.cc, pt-pr-code.cc, pt-select.cc, pt-stmt.cc, pt-unop.cc, pt.cc, token.cc, Array-jit.cc, Array-os.cc, Array-sym.cc, Array-tc.cc, version.cc: Indent #ifdef blocks in libinterp.
author Rik <rik@octave.org>
date Fri, 05 Feb 2016 16:29:08 -0800
parents 90cd0f9442d5
children 40de9f8f23a6
line wrap: on
line source

/*

Copyright (C) 2002-2015 John W. Eaton

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
<http://www.gnu.org/licenses/>.

*/

#ifdef HAVE_CONFIG_H
#  include <config.h>
#endif

#include <iostream>
#include <string>

#include "DASRT.h"
#include "lo-mappers.h"

#include "defun.h"
#include "error.h"
#include "errwarn.h"
#include "ovl.h"
#include "ov-fcn.h"
#include "ov-cell.h"
#include "pager.h"
#include "parse.h"
#include "unwind-prot.h"
#include "utils.h"
#include "variables.h"

#include "DASRT-opts.cc"

// Global pointers for user defined function required by dasrt.
static octave_function *dasrt_f;
static octave_function *dasrt_j;
static octave_function *dasrt_cf;

// Have we warned about imaginary values returned from user function?
static bool warned_fcn_imaginary = false;
static bool warned_jac_imaginary = false;
static bool warned_cf_imaginary = false;

// Is this a recursive call?
static int call_depth = 0;

static ColumnVector
dasrt_user_f (const ColumnVector& x, const ColumnVector& xdot,
              double t, octave_idx_type&)
{
  ColumnVector retval;

  assert (x.numel () == xdot.numel ());

  octave_value_list args;

  args(2) = t;
  args(1) = xdot;
  args(0) = x;

  if (dasrt_f)
    {
      octave_value_list tmp;

      try
        {
          tmp = dasrt_f->do_multi_index_op (1, args);
        }
      catch (octave_execution_exception& e)
        {
          err_user_supplied_eval (e, "dasrt");
        }

      if (tmp.length () == 0 || ! tmp(0).is_defined ())
        err_user_supplied_eval ("dasrt");

      if (! warned_fcn_imaginary && tmp(0).is_complex_type ())
        {
          warning ("dasrt: ignoring imaginary part returned from user-supplied function");
          warned_fcn_imaginary = true;
        }

      retval = tmp(0).vector_value ();

      if (retval.is_empty ())
        err_user_supplied_eval ("dasrt");
    }

  return retval;
}

static ColumnVector
dasrt_user_cf (const ColumnVector& x, double t)
{
  ColumnVector retval;

  octave_value_list args;

  args(1) = t;
  args(0) = x;

  if (dasrt_cf)
    {
      octave_value_list tmp;

      try
        {
          tmp = dasrt_cf->do_multi_index_op (1, args);
        }
      catch (octave_execution_exception& e)
        {
          err_user_supplied_eval (e, "dasrt");
        }

      if (tmp.length () == 0 || ! tmp(0).is_defined ())
        err_user_supplied_eval ("dasrt");

      if (! warned_cf_imaginary && tmp(0).is_complex_type ())
        {
          warning ("dasrt: ignoring imaginary part returned from user-supplied constraint function");
          warned_cf_imaginary = true;
        }

      retval = tmp(0).vector_value ();

      if (retval.is_empty ())
        err_user_supplied_eval ("dasrt");
    }

  return retval;
}

static Matrix
dasrt_user_j (const ColumnVector& x, const ColumnVector& xdot,
              double t, double cj)
{
  Matrix retval;

  assert (x.numel () == xdot.numel ());

  octave_value_list args;

  args(3) = cj;
  args(2) = t;
  args(1) = xdot;
  args(0) = x;

  if (dasrt_j)
    {
      octave_value_list tmp;

      try
        {
          tmp = dasrt_j->do_multi_index_op (1, args);
        }
      catch (octave_execution_exception& e)
        {
          err_user_supplied_eval (e, "dasrt");
        }

      int tlen = tmp.length ();
      if (tlen == 0 || ! tmp(0).is_defined ())
        err_user_supplied_eval ("dasrt");

      if (! warned_jac_imaginary && tmp(0).is_complex_type ())
        {
          warning ("dasrt: ignoring imaginary part returned from user-supplied jacobian function");
          warned_jac_imaginary = true;
        }

      retval = tmp(0).matrix_value ();

      if (retval.is_empty ())
        err_user_supplied_eval ("dasrt");
    }

  return retval;
}

DEFUN (dasrt, args, nargout,
       "-*- texinfo -*-\n\
@deftypefn  {} {[@var{x}, @var{xdot}, @var{t_out}, @var{istat}, @var{msg}] =} dasrt (@var{fcn}, [], @var{x_0}, @var{xdot_0}, @var{t})\n\
@deftypefnx {} {@dots{} =} dasrt (@var{fcn}, @var{g}, @var{x_0}, @var{xdot_0}, @var{t})\n\
@deftypefnx {} {@dots{} =} dasrt (@var{fcn}, [], @var{x_0}, @var{xdot_0}, @var{t}, @var{t_crit})\n\
@deftypefnx {} {@dots{} =} dasrt (@var{fcn}, @var{g}, @var{x_0}, @var{xdot_0}, @var{t}, @var{t_crit})\n\
Solve the set of differential-algebraic equations\n\
@tex\n\
$$ 0 = f (x, \\dot{x}, t) $$\n\
with\n\
$$ x(t_0) = x_0, \\dot{x}(t_0) = \\dot{x}_0 $$\n\
@end tex\n\
@ifnottex\n\
\n\
@example\n\
0 = f (x, xdot, t)\n\
@end example\n\
\n\
@noindent\n\
with\n\
\n\
@example\n\
x(t_0) = x_0, xdot(t_0) = xdot_0\n\
@end example\n\
\n\
@end ifnottex\n\
with functional stopping criteria (root solving).\n\
\n\
The solution is returned in the matrices @var{x} and @var{xdot},\n\
with each row in the result matrices corresponding to one of the\n\
elements in the vector @var{t_out}.  The first element of @var{t}\n\
should be @math{t_0} and correspond to the initial state of the\n\
system @var{x_0} and its derivative @var{xdot_0}, so that the first\n\
row of the output @var{x} is @var{x_0} and the first row\n\
of the output @var{xdot} is @var{xdot_0}.\n\
\n\
The vector @var{t} provides an upper limit on the length of the\n\
integration.  If the stopping condition is met, the vector\n\
@var{t_out} will be shorter than @var{t}, and the final element of\n\
@var{t_out} will be the point at which the stopping condition was met,\n\
and may not correspond to any element of the vector @var{t}.\n\
\n\
The first argument, @var{fcn}, is a string, inline, or function handle\n\
that names the function @math{f} to call to compute the vector of\n\
residuals for the set of equations.  It must have the form\n\
\n\
@example\n\
@var{res} = f (@var{x}, @var{xdot}, @var{t})\n\
@end example\n\
\n\
@noindent\n\
in which @var{x}, @var{xdot}, and @var{res} are vectors, and @var{t} is a\n\
scalar.\n\
\n\
If @var{fcn} is a two-element string array or a two-element cell array\n\
of strings, inline functions, or function handles, the first element names\n\
the function @math{f} described above, and the second element names a\n\
function to compute the modified Jacobian\n\
\n\
@tex\n\
$$\n\
J = {\\partial f \\over \\partial x}\n\
  + c {\\partial f \\over \\partial \\dot{x}}\n\
$$\n\
@end tex\n\
@ifnottex\n\
\n\
@example\n\
@group\n\
      df       df\n\
jac = -- + c ------\n\
      dx     d xdot\n\
@end group\n\
@end example\n\
\n\
@end ifnottex\n\
\n\
The modified Jacobian function must have the form\n\
\n\
@example\n\
@group\n\
\n\
@var{jac} = j (@var{x}, @var{xdot}, @var{t}, @var{c})\n\
\n\
@end group\n\
@end example\n\
\n\
The optional second argument names a function that defines the\n\
constraint functions whose roots are desired during the integration.\n\
This function must have the form\n\
\n\
@example\n\
@var{g_out} = g (@var{x}, @var{t})\n\
@end example\n\
\n\
@noindent\n\
and return a vector of the constraint function values.\n\
If the value of any of the constraint functions changes sign, @sc{dasrt}\n\
will attempt to stop the integration at the point of the sign change.\n\
\n\
If the name of the constraint function is omitted, @code{dasrt} solves\n\
the same problem as @code{daspk} or @code{dassl}.\n\
\n\
Note that because of numerical errors in the constraint functions\n\
due to round-off and integration error, @sc{dasrt} may return false\n\
roots, or return the same root at two or more nearly equal values of\n\
@var{T}.  If such false roots are suspected, the user should consider\n\
smaller error tolerances or higher precision in the evaluation of the\n\
constraint functions.\n\
\n\
If a root of some constraint function defines the end of the problem,\n\
the input to @sc{dasrt} should nevertheless allow integration to a\n\
point slightly past that root, so that @sc{dasrt} can locate the root\n\
by interpolation.\n\
\n\
The third and fourth arguments to @code{dasrt} specify the initial\n\
condition of the states and their derivatives, and the fourth argument\n\
specifies a vector of output times at which the solution is desired,\n\
including the time corresponding to the initial condition.\n\
\n\
The set of initial states and derivatives are not strictly required to\n\
be consistent.  In practice, however, @sc{dassl} is not very good at\n\
determining a consistent set for you, so it is best if you ensure that\n\
the initial values result in the function evaluating to zero.\n\
\n\
The sixth argument is optional, and may be used to specify a set of\n\
times that the DAE solver should not integrate past.  It is useful for\n\
avoiding difficulties with singularities and points where there is a\n\
discontinuity in the derivative.\n\
\n\
After a successful computation, the value of @var{istate} will be\n\
greater than zero (consistent with the Fortran version of @sc{dassl}).\n\
\n\
If the computation is not successful, the value of @var{istate} will be\n\
less than zero and @var{msg} will contain additional information.\n\
\n\
You can use the function @code{dasrt_options} to set optional\n\
parameters for @code{dasrt}.\n\
@seealso{dasrt_options, daspk, dasrt, lsode}\n\
@end deftypefn")
{
  int nargin = args.length ();

  if (nargin < 4 || nargin > 6)
    print_usage ();

  warned_fcn_imaginary = false;
  warned_jac_imaginary = false;
  warned_cf_imaginary = false;

  octave_value_list retval (5);

  unwind_protect frame;

  frame.protect_var (call_depth);
  call_depth++;

  if (call_depth > 1)
    error ("dasrt: invalid recursive call");

  int argp = 0;
  std::string fcn_name, fname, jac_name, jname;
  dasrt_f = 0;
  dasrt_j = 0;
  dasrt_cf = 0;

  // Check all the arguments.  Are they the right animals?

  // Here's where I take care of f and j in one shot:

  octave_value f_arg = args(0);

  if (f_arg.is_cell ())
    {
      Cell c = f_arg.cell_value ();
      if (c.numel () == 1)
        f_arg = c(0);
      else if (c.numel () == 2)
        {
          if (c(0).is_function_handle () || c(0).is_inline_function ())
            dasrt_f = c(0).function_value ();
          else
            {
              fcn_name = unique_symbol_name ("__dasrt_fcn__");
              fname = "function y = ";
              fname.append (fcn_name);
              fname.append (" (x, xdot, t) y = ");
              dasrt_f = extract_function (c(0), "dasrt", fcn_name, fname,
                                          "; endfunction");
            }

          if (dasrt_f)
            {
              if (c(1).is_function_handle () || c(1).is_inline_function ())
                dasrt_j = c(1).function_value ();
              else
                {
                  jac_name = unique_symbol_name ("__dasrt_jac__");
                  jname = "function jac = ";
                  jname.append (jac_name);
                  jname.append (" (x, xdot, t, cj) jac = ");
                  dasrt_j = extract_function (c(1), "dasrt", jac_name, jname,
                                              "; endfunction");

                  if (! dasrt_j)
                    {
                      if (fcn_name.length ())
                        clear_function (fcn_name);
                      dasrt_f = 0;
                    }
                }
            }
        }
      else
        error ("dasrt: incorrect number of elements in cell array");
    }

  if (! dasrt_f && ! f_arg.is_cell ())
    {
      if (f_arg.is_function_handle () || f_arg.is_inline_function ())
        dasrt_f = f_arg.function_value ();
      else
        {
          switch (f_arg.rows ())
            {
            case 1:
              fcn_name = unique_symbol_name ("__dasrt_fcn__");
              fname = "function y = ";
              fname.append (fcn_name);
              fname.append (" (x, xdot, t) y = ");
              dasrt_f = extract_function (f_arg, "dasrt", fcn_name, fname,
                                          "; endfunction");
              break;

            case 2:
              {
                string_vector tmp = args(0).string_vector_value ();

                fcn_name = unique_symbol_name ("__dasrt_fcn__");
                fname = "function y = ";
                fname.append (fcn_name);
                fname.append (" (x, xdot, t) y = ");
                dasrt_f = extract_function (tmp(0), "dasrt", fcn_name,
                                            fname, "; endfunction");

                if (dasrt_f)
                  {
                    jac_name = unique_symbol_name ("__dasrt_jac__");
                    jname = "function jac = ";
                    jname.append (jac_name);
                    jname.append (" (x, xdot, t, cj) jac = ");
                    dasrt_j = extract_function (tmp(1), "dasrt", jac_name,
                                                jname, "; endfunction");

                    if (! dasrt_j)
                      dasrt_f = 0;
                  }
              }
              break;

            default:
              error ("dasrt: first arg should be a string or 2-element string array");
            }
        }
    }

  if (! dasrt_f)
    return retval;

  DAERTFunc func (dasrt_user_f);

  argp++;

  if (args(1).is_function_handle () || args(1).is_inline_function ())
    {
      dasrt_cf = args(1).function_value ();

      if (! dasrt_cf)
        error ("dasrt: invalid constraint function G");

      argp++;

      func.set_constraint_function (dasrt_user_cf);
    }
  else if (args(1).is_string ())
    {
      dasrt_cf = is_valid_function (args(1), "dasrt", true);
      if (! dasrt_cf)
        error ("dasrt: invalid constraint function G");

      argp++;

      func.set_constraint_function (dasrt_user_cf);
    }

  ColumnVector state = args(argp).xvector_value ("dasrt: initial state X_0 must be a vector");

  ColumnVector stateprime = args(argp).xvector_value ("dasrt: initial derivatives XDOT_0 must be a vector");
  argp++;

  ColumnVector out_times = args(argp).xvector_value ("dasrt: output time variable T must be a vector");
  argp++;

  double tzero = out_times (0);

  ColumnVector crit_times;

  bool crit_times_set = false;

  if (argp < nargin)
    {
      crit_times = args(argp).xvector_value ("dasrt: list of critical times T_CRIT must be a vector");
      argp++;

      crit_times_set = true;
    }

  if (dasrt_j)
    func.set_jacobian_function (dasrt_user_j);

  DASRT_result output;

  DASRT dae = DASRT (state, stateprime, tzero, func);

  dae.set_options (dasrt_opts);

  if (crit_times_set)
    output = dae.integrate (out_times, crit_times);
  else
    output = dae.integrate (out_times);

  if (fcn_name.length ())
    clear_function (fcn_name);
  if (jac_name.length ())
    clear_function (jac_name);

  std::string msg = dae.error_message ();

  if (dae.integration_ok ())
    {
      retval(0) = output.state ();
      retval(1) = output.deriv ();
      retval(2) = output.times ();
    }
  else
    {
      if (nargout < 4)
        error ("dasrt: %s", msg.c_str ());

      retval(0) = Matrix ();
      retval(1) = Matrix ();
      retval(2) = Matrix ();
    }

  retval(3) = static_cast<double> (dae.integration_state ());
  retval(4) = msg;

  return retval;
}