changeset 22960:0d1422cb7e93

use F77_INT instead of octave_idx_type for liboctave ODE and DAE classes * DASPK.cc, DASRT.cc, DASSL.cc, LSODE.cc: Use F77_INT instead of octave_idx_type for integer data passed to Fortran subroutines. * DASPK.h, DASRT.h, DASSL.h, LSODE.h: Use octave_f77_int_type as needed in public interface files.
author John W. Eaton <jwe@octave.org>
date Tue, 27 Dec 2016 11:58:30 -0500
parents 12e4534d68b1
children 0e9606b04ae0
files liboctave/numeric/DASPK.cc liboctave/numeric/DASPK.h liboctave/numeric/DASRT.cc liboctave/numeric/DASRT.h liboctave/numeric/DASSL.cc liboctave/numeric/DASSL.h liboctave/numeric/LSODE.cc liboctave/numeric/LSODE.h
diffstat 8 files changed, 257 insertions(+), 247 deletions(-) [+]
line wrap: on
line diff
--- a/liboctave/numeric/DASPK.cc	Tue Dec 27 10:24:03 2016 -0500
+++ b/liboctave/numeric/DASPK.cc	Tue Dec 27 11:58:30 2016 -0500
@@ -34,46 +34,38 @@
 #include "lo-math.h"
 #include "quit.h"
 
-typedef octave_idx_type (*daspk_fcn_ptr) (const double&, const double*,
-                                          const double*, const double&,
-                                          double*, octave_idx_type&,
-                                          double*, octave_idx_type*);
+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 octave_idx_type (*daspk_jac_ptr) (const double&, const double*,
-                                          const double*, double*,
-                                          const double&, double*,
-                                          octave_idx_type*);
-
-typedef octave_idx_type (*daspk_psol_ptr) (const octave_idx_type&,
-                                           const double&, const double*,
-                                           const double*, const double*,
-                                           const double&, const double*,
-                                           double*, octave_idx_type*,
-                                           double*, const double&,
-                                           octave_idx_type&, double*,
-                                           octave_idx_type*);
+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*,
+  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_fun;
 static DAEFunc::DAEJacFunc user_jac;
-static octave_idx_type nn;
+static F77_INT nn;
 
-static octave_idx_type
+static F77_INT
 ddaspk_f (const double& time, const double *state, const double *deriv,
-          const double&, double *delta, octave_idx_type& ires, double *,
-          octave_idx_type *)
+          const double&, double *delta, F77_INT& ires, double *, F77_INT *)
 {
   BEGIN_INTERRUPT_WITH_EXCEPTIONS;
 
@@ -81,13 +73,17 @@
   ColumnVector tmp_state (nn);
   ColumnVector tmp_delta (nn);
 
-  for (octave_idx_type i = 0; i < nn; i++)
+  for (F77_INT i = 0; i < nn; i++)
     {
       tmp_deriv.elem (i) = deriv[i];
       tmp_state.elem (i) = state[i];
     }
 
-  tmp_delta = user_fun (tmp_state, tmp_deriv, time, ires);
+  octave_idx_type tmp_ires = ires;
+
+  tmp_delta = user_fun (tmp_state, tmp_deriv, time, tmp_ires);
+
+  ires = to_f77_int (tmp_ires);
 
   if (ires >= 0)
     {
@@ -95,7 +91,7 @@
         ires = -2;
       else
         {
-          for (octave_idx_type i = 0; i < nn; i++)
+          for (F77_INT i = 0; i < nn; i++)
             delta[i] = tmp_delta.elem (i);
         }
     }
@@ -108,11 +104,11 @@
 //NEQ, T, Y, YPRIME, SAVR, WK, CJ, WGHT,
 //C                          WP, IWP, B, EPLIN, IER, RPAR, IPAR)
 
-static octave_idx_type
-ddaspk_psol (const octave_idx_type&, const double&, const double *,
+static F77_INT
+ddaspk_psol (const F77_INT&, const double&, const double *,
              const double *, const double *, const double&,
-             const double *, double *, octave_idx_type *, double *,
-             const double&, octave_idx_type&, double *, octave_idx_type*)
+             const double *, double *, F77_INT *, double *,
+             const double&, F77_INT&, double *, F77_INT*)
 {
   BEGIN_INTERRUPT_WITH_EXCEPTIONS;
 
@@ -123,9 +119,9 @@
   return 0;
 }
 
-static octave_idx_type
+static F77_INT
 ddaspk_j (const double& time, const double *state, const double *deriv,
-          double *pd, const double& cj, double *, octave_idx_type *)
+          double *pd, const double& cj, double *, F77_INT *)
 {
   BEGIN_INTERRUPT_WITH_EXCEPTIONS;
 
@@ -134,7 +130,7 @@
   ColumnVector tmp_state (nn);
   ColumnVector tmp_deriv (nn);
 
-  for (octave_idx_type i = 0; i < nn; i++)
+  for (F77_INT i = 0; i < nn; i++)
     {
       tmp_deriv.elem (i) = deriv[i];
       tmp_state.elem (i) = state[i];
@@ -142,8 +138,8 @@
 
   Matrix tmp_pd = user_jac (tmp_state, tmp_deriv, time, cj);
 
-  for (octave_idx_type j = 0; j < nn; j++)
-    for (octave_idx_type i = 0; i < nn; i++)
+  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);
 
   END_INTERRUPT_WITH_EXCEPTIONS;
@@ -167,10 +163,10 @@
 
       info.resize (dim_vector (20, 1));
 
-      for (octave_idx_type i = 0; i < 20; i++)
+      for (F77_INT i = 0; i < 20; i++)
         info(i) = 0;
 
-      octave_idx_type n = size ();
+      F77_INT n = to_f77_int (size ());
 
       nn = n;
 
@@ -241,8 +237,8 @@
       abs_tol = absolute_tolerance ();
       rel_tol = relative_tolerance ();
 
-      octave_idx_type abs_tol_len = abs_tol.numel ();
-      octave_idx_type rel_tol_len = rel_tol.numel ();
+      F77_INT abs_tol_len = to_f77_int (abs_tol.numel ());
+      F77_INT rel_tol_len = to_f77_int (rel_tol.numel ());
 
       if (abs_tol_len == 1 && rel_tol_len == 1)
         {
@@ -287,7 +283,7 @@
           if (maxord > 0 && maxord < 6)
             {
               info(8) = 1;
-              iwork(2) = maxord;
+              iwork(2) = to_f77_int (maxord);
             }
           else
             {
@@ -306,11 +302,13 @@
           {
             Array<octave_idx_type> ict = inequality_constraint_types ();
 
-            if (ict.numel () == n)
+            F77_INT ict_nel = to_f77_int (ict.numel ());
+
+            if (ict_nel == n)
               {
-                for (octave_idx_type i = 0; i < n; i++)
+                for (F77_INT i = 0; i < n; i++)
                   {
-                    octave_idx_type val = ict(i);
+                    F77_INT val = to_f77_int (ict(i));
                     if (val < -2 || val > 2)
                       {
                         // FIXME: Should this be a warning?
@@ -335,7 +333,7 @@
 
         case 0:
         case 2:
-          info(9) = eiq;
+          info(9) = to_f77_int (eiq);
           break;
 
         default:
@@ -354,9 +352,11 @@
 
               Array<octave_idx_type> av = algebraic_variables ();
 
-              if (av.numel () == n)
+              F77_INT av_nel = to_f77_int (av.numel ());
+
+              if (av_nel == n)
                 {
-                  octave_idx_type lid;
+                  F77_INT lid;
                   if (eiq == 0 || eiq == 2)
                     lid = 40;
                   else if (eiq == 1 || eiq == 3)
@@ -364,7 +364,7 @@
                   else
                     abort ();
 
-                  for (octave_idx_type i = 0; i < n; i++)
+                  for (F77_INT i = 0; i < n; i++)
                     iwork(lid+i) = av(i) ? -1 : 1;
                 }
               else
@@ -385,7 +385,7 @@
               return retval;
             }
 
-          info(10) = ccic;
+          info(10) = to_f77_int (ccic);
         }
 
       if (eavfet)
@@ -396,9 +396,11 @@
 
           Array<octave_idx_type> av = algebraic_variables ();
 
-          if (av.numel () == n)
+          F77_INT av_nel = to_f77_int (av.numel ());
+
+          if (av_nel == n)
             {
-              octave_idx_type lid;
+              F77_INT lid;
               if (eiq == 0 || eiq == 2)
                 lid = 40;
               else if (eiq == 1 || eiq == 3)
@@ -406,7 +408,7 @@
               else
                 abort ();
 
-              for (octave_idx_type i = 0; i < n; i++)
+              for (F77_INT i = 0; i < n; i++)
                 iwork(lid+i) = av(i) ? -1 : 1;
             }
         }
@@ -417,10 +419,10 @@
 
           if (ich.numel () == 6)
             {
-              iwork(31) = octave::math::nint_big (ich(0));
-              iwork(32) = octave::math::nint_big (ich(1));
-              iwork(33) = octave::math::nint_big (ich(2));
-              iwork(34) = octave::math::nint_big (ich(3));
+              iwork(31) = to_f77_int (octave::math::nint_big (ich(0)));
+              iwork(32) = to_f77_int (octave::math::nint_big (ich(1)));
+              iwork(33) = to_f77_int (octave::math::nint_big (ich(2)));
+              iwork(34) = to_f77_int (octave::math::nint_big (ich(3)));
 
               rwork(13) = ich(4);
               rwork(14) = ich(5);
@@ -443,7 +445,7 @@
         case 0:
         case 1:
         case 2:
-          info(17) = pici;
+          info(17) = to_f77_int (pici);
           break;
 
         default:
@@ -463,22 +465,26 @@
   double *px = x.fortran_vec ();
   double *pxdot = xdot.fortran_vec ();
 
-  octave_idx_type *pinfo = info.fortran_vec ();
+  F77_INT *pinfo = info.fortran_vec ();
 
   double *prel_tol = rel_tol.fortran_vec ();
   double *pabs_tol = abs_tol.fortran_vec ();
 
   double *prwork = rwork.fortran_vec ();
-  octave_idx_type *piwork = iwork.fortran_vec ();
+  F77_INT *piwork = iwork.fortran_vec ();
 
   double *dummy = 0;
-  octave_idx_type *idummy = 0;
+  F77_INT *idummy = 0;
+
+  F77_INT tmp_istate = to_f77_int (istate);
 
   F77_XFCN (ddaspk, DDASPK, (ddaspk_f, nn, t, px, pxdot, tout, pinfo,
-                             prel_tol, pabs_tol, istate, prwork, lrw,
+                             prel_tol, pabs_tol, tmp_istate, prwork, lrw,
                              piwork, liw, dummy, idummy, ddaspk_j,
                              ddaspk_psol));
 
+  istate = tmp_istate;
+
   switch (istate)
     {
     case 1: // A step was successfully taken in intermediate-output
@@ -549,14 +555,14 @@
   Matrix retval;
 
   octave_idx_type n_out = tout.numel ();
-  octave_idx_type n = size ();
+  F77_INT n = to_f77_int (size ());
 
   if (n_out > 0 && n > 0)
     {
       retval.resize (n_out, n);
       xdot_out.resize (n_out, n);
 
-      for (octave_idx_type i = 0; i < n; i++)
+      for (F77_INT i = 0; i < n; i++)
         {
           retval.elem (0, i) = x.elem (i);
           xdot_out.elem (0, i) = xdot.elem (i);
@@ -569,7 +575,7 @@
           if (integration_error)
             return retval;
 
-          for (octave_idx_type i = 0; i < n; i++)
+          for (F77_INT i = 0; i < n; i++)
             {
               retval.elem (j, i) = x_next.elem (i);
               xdot_out.elem (j, i) = xdot.elem (i);
@@ -594,14 +600,14 @@
   Matrix retval;
 
   octave_idx_type n_out = tout.numel ();
-  octave_idx_type n = size ();
+  F77_INT n = to_f77_int (size ());
 
   if (n_out > 0 && n > 0)
     {
       retval.resize (n_out, n);
       xdot_out.resize (n_out, n);
 
-      for (octave_idx_type i = 0; i < n; i++)
+      for (F77_INT i = 0; i < n; i++)
         {
           retval.elem (0, i) = x.elem (i);
           xdot_out.elem (0, i) = xdot.elem (i);
@@ -668,7 +674,7 @@
 
               if (save_output)
                 {
-                  for (octave_idx_type i = 0; i < n; i++)
+                  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) = xdot.elem (i);
--- a/liboctave/numeric/DASPK.h	Tue Dec 27 10:24:03 2016 -0500
+++ b/liboctave/numeric/DASPK.h	Tue Dec 27 11:58:30 2016 -0500
@@ -68,11 +68,11 @@
 
   bool initialized;
 
-  octave_idx_type liw;
-  octave_idx_type lrw;
+  octave_f77_int_type liw;
+  octave_f77_int_type lrw;
 
-  Array<octave_idx_type> info;
-  Array<octave_idx_type> iwork;
+  Array<octave_f77_int_type> info;
+  Array<octave_f77_int_type> iwork;
 
   Array<double> rwork;
 
--- a/liboctave/numeric/DASRT.cc	Tue Dec 27 10:24:03 2016 -0500
+++ b/liboctave/numeric/DASRT.cc	Tue Dec 27 11:58:30 2016 -0500
@@ -34,33 +34,25 @@
 #include "lo-math.h"
 #include "quit.h"
 
-typedef octave_idx_type (*dasrt_fcn_ptr) (const double&, const double*,
-                                          const double*, double*,
-                                          octave_idx_type&, double*,
-                                          octave_idx_type*);
+typedef F77_INT (*dasrt_fcn_ptr) (const double&, const double*, const double*,
+                                  double*, F77_INT&, double*, F77_INT*);
 
-typedef octave_idx_type (*dasrt_jac_ptr) (const double&, const double*,
-                                          const double*, double*,
-                                          const double&, double*,
-                                          octave_idx_type*);
+typedef F77_INT (*dasrt_jac_ptr) (const double&, const double*, const double*,
+                                  double*, const double&, double*, F77_INT*);
 
-typedef octave_idx_type (*dasrt_constr_ptr) (const octave_idx_type&,
-                                             const double&, const double*,
-                                             const octave_idx_type&,
-                                             double*, double*,
-                                             octave_idx_type*);
+typedef F77_INT (*dasrt_constr_ptr) (const F77_INT&, const double&,
+                                     const double*, const F77_INT&,
+                                     double*, double*, F77_INT*);
 
 extern "C"
 {
   F77_RET_T
-  F77_FUNC (ddasrt, DDASRT) (dasrt_fcn_ptr, const F77_INT&,
-                             F77_DBLE&, F77_DBLE*, F77_DBLE*, const F77_DBLE&,
-                             F77_INT*, const F77_DBLE*,
-                             const F77_DBLE*, F77_INT&, F77_DBLE*,
-                             const F77_INT&, F77_INT*,
-                             const F77_INT&, F77_DBLE*,
-                             F77_INT*, dasrt_jac_ptr,
-                             dasrt_constr_ptr, const F77_INT&,
+  F77_FUNC (ddasrt, DDASRT) (dasrt_fcn_ptr, const F77_INT&, F77_DBLE&,
+                             F77_DBLE*, F77_DBLE*, const F77_DBLE&, F77_INT*,
+                             const F77_DBLE*, const F77_DBLE*, F77_INT&,
+                             F77_DBLE*, const F77_INT&, F77_INT*,
+                             const F77_INT&, F77_DBLE*, F77_INT*,
+                             dasrt_jac_ptr, dasrt_constr_ptr, const F77_INT&,
                              F77_INT*);
 }
 
@@ -68,30 +60,34 @@
 static DAEFunc::DAEJacFunc user_jsub;
 static DAERTFunc::DAERTConstrFunc user_csub;
 
-static octave_idx_type nn;
+static F77_INT nn;
 
-static octave_idx_type
+static F77_INT
 ddasrt_f (const double& t, const double *state, const double *deriv,
-          double *delta, octave_idx_type& ires, double *, octave_idx_type *)
+          double *delta, F77_INT& ires, double *, F77_INT *)
 {
   BEGIN_INTERRUPT_WITH_EXCEPTIONS;
 
   ColumnVector tmp_state (nn);
   ColumnVector tmp_deriv (nn);
 
-  for (octave_idx_type i = 0; i < nn; i++)
+  for (F77_INT i = 0; i < nn; i++)
     {
       tmp_state(i) = state[i];
       tmp_deriv(i) = deriv[i];
     }
 
-  ColumnVector tmp_fval = (*user_fsub) (tmp_state, tmp_deriv, t, ires);
+  octave_idx_type tmp_ires = ires;
+
+  ColumnVector tmp_fval = (*user_fsub) (tmp_state, tmp_deriv, t, tmp_ires);
+
+  ires = to_f77_int (tmp_ires);
 
   if (tmp_fval.is_empty ())
     ires = -2;
   else
     {
-      for (octave_idx_type i = 0; i < nn; i++)
+      for (F77_INT i = 0; i < nn; i++)
         delta[i] = tmp_fval(i);
     }
 
@@ -100,9 +96,9 @@
   return 0;
 }
 
-octave_idx_type
+F77_INT
 ddasrt_j (const double& time, const double *state, const double *deriv,
-          double *pd, const double& cj, double *, octave_idx_type *)
+          double *pd, const double& cj, double *, F77_INT *)
 {
   BEGIN_INTERRUPT_WITH_EXCEPTIONS;
 
@@ -111,7 +107,7 @@
   ColumnVector tmp_state (nn);
   ColumnVector tmp_deriv (nn);
 
-  for (octave_idx_type i = 0; i < nn; i++)
+  for (F77_INT i = 0; i < nn; i++)
     {
       tmp_deriv.elem (i) = deriv[i];
       tmp_state.elem (i) = state[i];
@@ -119,8 +115,8 @@
 
   Matrix tmp_pd = (*user_jsub) (tmp_state, tmp_deriv, time, cj);
 
-  for (octave_idx_type j = 0; j < nn; j++)
-    for (octave_idx_type i = 0; i < nn; i++)
+  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);
 
   END_INTERRUPT_WITH_EXCEPTIONS;
@@ -128,21 +124,21 @@
   return 0;
 }
 
-static octave_idx_type
-ddasrt_g (const octave_idx_type& neq, const double& t, const double *state,
-          const octave_idx_type& ng, double *gout, double *, octave_idx_type *)
+static F77_INT
+ddasrt_g (const F77_INT& neq, const double& t, const double *state,
+          const F77_INT& ng, double *gout, double *, F77_INT *)
 {
   BEGIN_INTERRUPT_WITH_EXCEPTIONS;
 
-  octave_idx_type n = neq;
+  F77_INT n = neq;
 
   ColumnVector tmp_state (n);
-  for (octave_idx_type i = 0; i < n; i++)
+  for (F77_INT i = 0; i < n; i++)
     tmp_state(i) = state[i];
 
   ColumnVector tmp_fval = (*user_csub) (tmp_state, t);
 
-  for (octave_idx_type i = 0; i < ng; i++)
+  for (F77_INT i = 0; i < ng; i++)
     gout[i] = tmp_fval(i);
 
   END_INTERRUPT_WITH_EXCEPTIONS;
@@ -166,10 +162,10 @@
 
       info.resize (dim_vector (15, 1));
 
-      for (octave_idx_type i = 0; i < 15; i++)
+      for (F77_INT i = 0; i < 15; i++)
         info(i) = 0;
 
-      octave_idx_type n = size ();
+      F77_INT n = to_f77_int (size ());
 
       nn = n;
 
@@ -180,12 +176,12 @@
       if (user_csub)
         {
           ColumnVector tmp = (*user_csub) (x, t);
-          ng = tmp.numel ();
+          ng = to_f77_int (tmp.numel ());
         }
       else
         ng = 0;
 
-      octave_idx_type maxord = maximum_order ();
+      F77_INT maxord = to_f77_int (maximum_order ());
       if (maxord >= 0)
         {
           if (maxord > 0 && maxord < 6)
@@ -277,10 +273,11 @@
       else
         info(7) = 0;
 
-      if (step_limit () >= 0)
+      F77_INT sl = to_f77_int (step_limit ());
+      if (sl >= 0)
         {
           info(11) = 1;
-          iwork(20) = step_limit ();
+          iwork(20) = sl;
         }
       else
         info(11) = 0;
@@ -288,8 +285,8 @@
       abs_tol = absolute_tolerance ();
       rel_tol = relative_tolerance ();
 
-      octave_idx_type abs_tol_len = abs_tol.numel ();
-      octave_idx_type rel_tol_len = rel_tol.numel ();
+      F77_INT abs_tol_len = to_f77_int (abs_tol.numel ());
+      F77_INT rel_tol_len = to_f77_int (rel_tol.numel ());
 
       if (abs_tol_len == 1 && rel_tol_len == 1)
         {
@@ -314,24 +311,28 @@
   double *px = x.fortran_vec ();
   double *pxdot = xdot.fortran_vec ();
 
-  octave_idx_type *pinfo = info.fortran_vec ();
+  F77_INT *pinfo = info.fortran_vec ();
 
   double *prel_tol = rel_tol.fortran_vec ();
   double *pabs_tol = abs_tol.fortran_vec ();
 
   double *prwork = rwork.fortran_vec ();
-  octave_idx_type *piwork = iwork.fortran_vec ();
+  F77_INT *piwork = iwork.fortran_vec ();
 
-  octave_idx_type *pjroot = jroot.fortran_vec ();
+  F77_INT *pjroot = jroot.fortran_vec ();
 
   double *dummy = 0;
-  octave_idx_type *idummy = 0;
+  F77_INT *idummy = 0;
+
+  F77_INT tmp_istate = to_f77_int (istate);
 
   F77_XFCN (ddasrt, DDASRT, (ddasrt_f, nn, t, px, pxdot, tout, pinfo,
-                             prel_tol, pabs_tol, istate, prwork, lrw,
+                             prel_tol, pabs_tol, tmp_istate, prwork, lrw,
                              piwork, liw, dummy, idummy, ddasrt_j,
                              ddasrt_g, ng, pjroot));
 
+  istate = tmp_istate;
+
   switch (istate)
     {
     case 1: // A step was successfully taken in intermediate-output
@@ -392,14 +393,14 @@
   ColumnVector t_out = tout;
 
   octave_idx_type n_out = tout.numel ();
-  octave_idx_type n = size ();
+  F77_INT n = to_f77_int (size ());
 
   if (n_out > 0 && n > 0)
     {
       x_out.resize (n_out, n);
       xdot_out.resize (n_out, n);
 
-      for (octave_idx_type i = 0; i < n; i++)
+      for (F77_INT i = 0; i < n; i++)
         {
           x_out(0,i) = x(i);
           xdot_out(0,i) = xdot(i);
@@ -420,7 +421,7 @@
           else
             t_out(j) = tout(j);
 
-          for (octave_idx_type i = 0; i < n; i++)
+          for (F77_INT i = 0; i < n; i++)
             {
               x_out(j,i) = x(i);
               xdot_out(j,i) = xdot(i);
@@ -451,7 +452,7 @@
   ColumnVector t_outs = tout;
 
   octave_idx_type n_out = tout.numel ();
-  octave_idx_type n = size ();
+  F77_INT n = to_f77_int (size ());
 
   if (n_out > 0 && n > 0)
     {
@@ -474,14 +475,14 @@
               if (i_crit < n_crit)
                 next_crit = tcrit(i_crit);
 
-              octave_idx_type save_output;
+              bool save_output = false;
               double t_out;
 
               if (next_crit == next_out)
                 {
                   set_stop_time (next_crit);
                   t_out = next_out;
-                  save_output = 1;
+                  save_output = true;
                   i_out++;
                   i_crit++;
                   do_restart = true;
@@ -492,7 +493,7 @@
                     {
                       set_stop_time (next_crit);
                       t_out = next_crit;
-                      save_output = 0;
+                      save_output = false;
                       i_crit++;
                       do_restart = true;
                     }
@@ -500,7 +501,7 @@
                     {
                       clear_stop_time ();
                       t_out = next_out;
-                      save_output = 1;
+                      save_output = true;
                       i_out++;
                     }
                 }
@@ -508,7 +509,7 @@
                 {
                   set_stop_time (next_crit);
                   t_out = next_out;
-                  save_output = 1;
+                  save_output = true;
                   i_out++;
                 }
 
@@ -525,7 +526,7 @@
 
               if (save_output)
                 {
-                  for (octave_idx_type i = 0; i < n; i++)
+                  for (F77_INT i = 0; i < n; i++)
                     {
                       x_out(i_out-1,i) = x(i);
                       xdot_out(i_out-1,i) = xdot(i);
--- a/liboctave/numeric/DASRT.h	Tue Dec 27 10:24:03 2016 -0500
+++ b/liboctave/numeric/DASRT.h	Tue Dec 27 11:58:30 2016 -0500
@@ -105,14 +105,14 @@
 
   bool initialized;
 
-  octave_idx_type liw;
-  octave_idx_type lrw;
+  octave_f77_int_type liw;
+  octave_f77_int_type lrw;
 
-  octave_idx_type ng;
+  octave_f77_int_type ng;
 
-  Array<octave_idx_type> info;
-  Array<octave_idx_type> iwork;
-  Array<octave_idx_type> jroot;
+  Array<octave_f77_int_type> info;
+  Array<octave_f77_int_type> iwork;
+  Array<octave_f77_int_type> jroot;
 
   Array<double> rwork;
 
--- a/liboctave/numeric/DASSL.cc	Tue Dec 27 10:24:03 2016 -0500
+++ b/liboctave/numeric/DASSL.cc	Tue Dec 27 11:58:30 2016 -0500
@@ -34,37 +34,33 @@
 #include "lo-math.h"
 #include "quit.h"
 
-typedef octave_idx_type (*dassl_fcn_ptr) (const double&, const double*,
-                                          const double*, double*,
-                                          octave_idx_type&, double*,
-                                          octave_idx_type*);
+typedef F77_INT (*dassl_fcn_ptr) (const double&, const double*,
+                                  const double*, double*, F77_INT&,
+                                  double*, F77_INT*);
 
-typedef octave_idx_type (*dassl_jac_ptr) (const double&, const double*,
-                                          const double*, double*,
-                                          const double&, double*,
-                                          octave_idx_type*);
+typedef F77_INT (*dassl_jac_ptr) (const double&, const double*,
+                                  const double*, double*, const double&,
+                                  double*, F77_INT*);
 
 extern "C"
 {
   F77_RET_T
-  F77_FUNC (ddassl, DDASSL) (dassl_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*,
+  F77_FUNC (ddassl, DDASSL) (dassl_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*,
                              dassl_jac_ptr);
 }
 
 static DAEFunc::DAERHSFunc user_fun;
 static DAEFunc::DAEJacFunc user_jac;
 
-static octave_idx_type nn;
+static F77_INT nn;
 
-static octave_idx_type
+static F77_INT
 ddassl_f (const double& time, const double *state, const double *deriv,
-          double *delta, octave_idx_type& ires, double *, octave_idx_type *)
+          double *delta, F77_INT& ires, double *, F77_INT *)
 {
   BEGIN_INTERRUPT_WITH_EXCEPTIONS;
 
@@ -74,13 +70,17 @@
   ColumnVector tmp_state (nn);
   ColumnVector tmp_delta (nn);
 
-  for (octave_idx_type i = 0; i < nn; i++)
+  for (F77_INT i = 0; i < nn; i++)
     {
       tmp_deriv.elem (i) = deriv[i];
       tmp_state.elem (i) = state[i];
     }
 
-  tmp_delta = user_fun (tmp_state, tmp_deriv, time, ires);
+  octave_idx_type tmp_ires = ires;
+
+  tmp_delta = user_fun (tmp_state, tmp_deriv, time, tmp_ires);
+
+  ires = to_f77_int (tmp_ires);
 
   if (ires >= 0)
     {
@@ -88,7 +88,7 @@
         ires = -2;
       else
         {
-          for (octave_idx_type i = 0; i < nn; i++)
+          for (F77_INT i = 0; i < nn; i++)
             delta[i] = tmp_delta.elem (i);
         }
     }
@@ -98,9 +98,9 @@
   return 0;
 }
 
-static octave_idx_type
+static F77_INT
 ddassl_j (const double& time, const double *state, const double *deriv,
-          double *pd, const double& cj, double *, octave_idx_type *)
+          double *pd, const double& cj, double *, F77_INT *)
 {
   BEGIN_INTERRUPT_WITH_EXCEPTIONS;
 
@@ -109,7 +109,7 @@
   ColumnVector tmp_state (nn);
   ColumnVector tmp_deriv (nn);
 
-  for (octave_idx_type i = 0; i < nn; i++)
+  for (F77_INT i = 0; i < nn; i++)
     {
       tmp_deriv.elem (i) = deriv[i];
       tmp_state.elem (i) = state[i];
@@ -117,8 +117,8 @@
 
   Matrix tmp_pd = user_jac (tmp_state, tmp_deriv, time, cj);
 
-  for (octave_idx_type j = 0; j < nn; j++)
-    for (octave_idx_type i = 0; i < nn; i++)
+  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);
 
   END_INTERRUPT_WITH_EXCEPTIONS;
@@ -139,10 +139,10 @@
 
       info.resize (dim_vector (15, 1));
 
-      for (octave_idx_type i = 0; i < 15; i++)
+      for (F77_INT i = 0; i < 15; i++)
         info(i) = 0;
 
-      octave_idx_type n = size ();
+      F77_INT n = to_f77_int (size ());
 
       liw = 21 + n;
       lrw = 40 + 9*n + n*n;
@@ -217,15 +217,17 @@
       else
         info(7) = 0;
 
-      if (step_limit () >= 0)
+      F77_INT sl = to_f77_int (step_limit ());
+
+      if (sl >= 0)
         {
           info(11) = 1;
-          iwork(20) = step_limit ();
+          iwork(20) = sl;
         }
       else
         info(11) = 0;
 
-      octave_idx_type maxord = maximum_order ();
+      F77_INT maxord = to_f77_int (maximum_order ());
       if (maxord >= 0)
         {
           if (maxord > 0 && maxord < 6)
@@ -242,17 +244,17 @@
             }
         }
 
-      octave_idx_type enc = enforce_nonnegativity_constraints ();
+      F77_INT enc = to_f77_int (enforce_nonnegativity_constraints ());
       info(9) = enc ? 1 : 0;
 
-      octave_idx_type ccic = compute_consistent_initial_condition ();
+      F77_INT ccic = to_f77_int (compute_consistent_initial_condition ());
       info(10) = ccic ? 1 : 0;
 
       abs_tol = absolute_tolerance ();
       rel_tol = relative_tolerance ();
 
-      octave_idx_type abs_tol_len = abs_tol.numel ();
-      octave_idx_type rel_tol_len = rel_tol.numel ();
+      F77_INT abs_tol_len = to_f77_int (abs_tol.numel ());
+      F77_INT rel_tol_len = to_f77_int (rel_tol.numel ());
 
       if (abs_tol_len == 1 && rel_tol_len == 1)
         {
@@ -277,21 +279,25 @@
   double *px = x.fortran_vec ();
   double *pxdot = xdot.fortran_vec ();
 
-  octave_idx_type *pinfo = info.fortran_vec ();
+  F77_INT *pinfo = info.fortran_vec ();
 
   double *prel_tol = rel_tol.fortran_vec ();
   double *pabs_tol = abs_tol.fortran_vec ();
 
   double *prwork = rwork.fortran_vec ();
-  octave_idx_type *piwork = iwork.fortran_vec ();
+  F77_INT *piwork = iwork.fortran_vec ();
 
   double *dummy = 0;
-  octave_idx_type *idummy = 0;
+  F77_INT *idummy = 0;
+
+  F77_INT tmp_istate = to_f77_int (istate);
 
   F77_XFCN (ddassl, DDASSL, (ddassl_f, nn, t, px, pxdot, tout, pinfo,
-                             prel_tol, pabs_tol, istate, prwork, lrw,
+                             prel_tol, pabs_tol, tmp_istate, prwork, lrw,
                              piwork, liw, dummy, idummy, ddassl_j));
 
+  istate = tmp_istate;
+
   switch (istate)
     {
     case 1: // A step was successfully taken in intermediate-output
@@ -354,14 +360,14 @@
   Matrix retval;
 
   octave_idx_type n_out = tout.numel ();
-  octave_idx_type n = size ();
+  F77_INT n = to_f77_int (size ());
 
   if (n_out > 0 && n > 0)
     {
       retval.resize (n_out, n);
       xdot_out.resize (n_out, n);
 
-      for (octave_idx_type i = 0; i < n; i++)
+      for (F77_INT i = 0; i < n; i++)
         {
           retval.elem (0, i) = x.elem (i);
           xdot_out.elem (0, i) = xdot.elem (i);
@@ -374,7 +380,7 @@
           if (integration_error)
             return retval;
 
-          for (octave_idx_type i = 0; i < n; i++)
+          for (F77_INT i = 0; i < n; i++)
             {
               retval.elem (j, i) = x_next.elem (i);
               xdot_out.elem (j, i) = xdot.elem (i);
@@ -399,14 +405,14 @@
   Matrix retval;
 
   octave_idx_type n_out = tout.numel ();
-  octave_idx_type n = size ();
+  F77_INT n = to_f77_int (size ());
 
   if (n_out > 0 && n > 0)
     {
       retval.resize (n_out, n);
       xdot_out.resize (n_out, n);
 
-      for (octave_idx_type i = 0; i < n; i++)
+      for (F77_INT i = 0; i < n; i++)
         {
           retval.elem (0, i) = x.elem (i);
           xdot_out.elem (0, i) = xdot.elem (i);
@@ -473,7 +479,7 @@
 
               if (save_output)
                 {
-                  for (octave_idx_type i = 0; i < n; i++)
+                  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) = xdot.elem (i);
--- a/liboctave/numeric/DASSL.h	Tue Dec 27 10:24:03 2016 -0500
+++ b/liboctave/numeric/DASSL.h	Tue Dec 27 11:58:30 2016 -0500
@@ -68,11 +68,11 @@
 
   bool initialized;
 
-  octave_idx_type liw;
-  octave_idx_type lrw;
+  octave_f77_int_type liw;
+  octave_f77_int_type lrw;
 
-  Array<octave_idx_type> info;
-  Array<octave_idx_type> iwork;
+  Array<octave_f77_int_type> info;
+  Array<octave_f77_int_type> iwork;
 
   Array<double> rwork;
 
--- a/liboctave/numeric/LSODE.cc	Tue Dec 27 10:24:03 2016 -0500
+++ b/liboctave/numeric/LSODE.cc	Tue Dec 27 11:58:30 2016 -0500
@@ -34,25 +34,20 @@
 #include "lo-math.h"
 #include "quit.h"
 
-typedef octave_idx_type (*lsode_fcn_ptr) (const octave_idx_type&,
-                                          const double&, double*,
-                                          double*, octave_idx_type&);
+typedef F77_INT (*lsode_fcn_ptr) (const F77_INT&, const double&, double*,
+                                  double*, F77_INT&);
 
-typedef octave_idx_type (*lsode_jac_ptr) (const octave_idx_type&,
-                                          const double&, double*,
-                                          const octave_idx_type&,
-                                          const octave_idx_type&,
-                                          double*, const octave_idx_type&);
+typedef F77_INT (*lsode_jac_ptr) (const F77_INT&, const double&, double*,
+                                  const F77_INT&, const F77_INT&, double*,
+                                  const F77_INT&);
 
 extern "C"
 {
   F77_RET_T
-  F77_FUNC (dlsode, DLSODE) (lsode_fcn_ptr, F77_INT&, F77_DBLE*,
-                             F77_DBLE&, F77_DBLE&, F77_INT&, F77_DBLE&,
-                             const F77_DBLE*, F77_INT&,
-                             F77_INT&, F77_INT&,
-                             F77_DBLE*, F77_INT&, F77_INT*,
-                             F77_INT&, lsode_jac_ptr,
+  F77_FUNC (dlsode, DLSODE) (lsode_fcn_ptr, F77_INT&, F77_DBLE*, F77_DBLE&,
+                             F77_DBLE&, F77_INT&, F77_DBLE&, const F77_DBLE*,
+                             F77_INT&, F77_INT&, F77_INT&, F77_DBLE*,
+                             F77_INT&, F77_INT*, F77_INT&, lsode_jac_ptr,
                              F77_INT&);
 }
 
@@ -60,9 +55,9 @@
 static ODEFunc::ODEJacFunc user_jac;
 static ColumnVector *tmp_x;
 
-static octave_idx_type
-lsode_f (const octave_idx_type& neq, const double& time, double *,
-         double *deriv, octave_idx_type& ierr)
+static F77_INT
+lsode_f (const F77_INT& neq, const double& time, double *, double *deriv,
+         F77_INT& ierr)
 {
   BEGIN_INTERRUPT_WITH_EXCEPTIONS;
 
@@ -78,7 +73,7 @@
     ierr = -1;
   else
     {
-      for (octave_idx_type i = 0; i < neq; i++)
+      for (F77_INT i = 0; i < neq; i++)
         deriv[i] = tmp_deriv.elem (i);
     }
 
@@ -87,10 +82,9 @@
   return 0;
 }
 
-static octave_idx_type
-lsode_j (const octave_idx_type& neq, const double& time, double *,
-         const octave_idx_type&, const octave_idx_type&, double *pd,
-         const octave_idx_type& nrowpd)
+static F77_INT
+lsode_j (const F77_INT& neq, const double& time, double *, const F77_INT&,
+         const F77_INT&, double *pd, const F77_INT& nrowpd)
 {
   BEGIN_INTERRUPT_WITH_EXCEPTIONS;
 
@@ -102,8 +96,8 @@
 
   tmp_jac = (*user_jac) (*tmp_x, time);
 
-  for (octave_idx_type j = 0; j < neq; j++)
-    for (octave_idx_type i = 0; i < neq; i++)
+  for (F77_INT j = 0; j < neq; j++)
+    for (F77_INT i = 0; i < neq; i++)
       pd[nrowpd * j + i] = tmp_jac (i, j);
 
   END_INTERRUPT_WITH_EXCEPTIONS;
@@ -116,7 +110,7 @@
 {
   ColumnVector retval;
 
-  static octave_idx_type nn = 0;
+  static F77_INT nn = 0;
 
   if (! initialized || restart || ODEFunc::reset || LSODE_options::reset)
     {
@@ -126,7 +120,7 @@
 
       istate = 1;
 
-      octave_idx_type n = size ();
+      F77_INT n = to_f77_int (size ());
 
       nn = n;
 
@@ -154,23 +148,23 @@
           lrw = 22 + 16 * n;
         }
 
-      maxord = maximum_order ();
-
       iwork.resize (dim_vector (liw, 1));
 
-      for (octave_idx_type i = 4; i < 9; i++)
+      for (F77_INT i = 4; i < 9; i++)
         iwork(i) = 0;
 
       rwork.resize (dim_vector (lrw, 1));
 
-      for (octave_idx_type i = 4; i < 9; i++)
+      for (F77_INT i = 4; i < 9; i++)
         rwork(i) = 0;
 
+      octave_idx_type maxord = maximum_order ();
+
       if (maxord >= 0)
         {
           if (maxord > 0 && maxord <= max_maxord)
             {
-              iwork(4) = maxord;
+              iwork(4) = to_f77_int (maxord);
               iopt = 1;
             }
           else
@@ -226,7 +220,7 @@
       rel_tol = relative_tolerance ();
       abs_tol = absolute_tolerance ();
 
-      octave_idx_type abs_tol_len = abs_tol.numel ();
+      F77_INT abs_tol_len = to_f77_int (abs_tol.numel ());
 
       if (abs_tol_len == 1)
         itol = 1;
@@ -263,7 +257,7 @@
           iopt = 1;
         }
 
-      octave_idx_type sl = step_limit ();
+      F77_INT sl = to_f77_int (step_limit ());
       if (sl > 0)
         {
           iwork(5) = sl;
@@ -277,13 +271,17 @@
 
   double *pabs_tol = abs_tol.fortran_vec ();
 
-  octave_idx_type *piwork = iwork.fortran_vec ();
+  F77_INT *piwork = iwork.fortran_vec ();
   double *prwork = rwork.fortran_vec ();
 
+  F77_INT tmp_istate = to_f77_int (istate);
+
   F77_XFCN (dlsode, DLSODE, (lsode_f, nn, px, t, tout, itol, rel_tol,
-                             pabs_tol, itask, istate, iopt, prwork, lrw,
+                             pabs_tol, itask, tmp_istate, iopt, prwork, lrw,
                              piwork, liw, lsode_j, method_flag));
 
+  istate = tmp_istate;
+
   switch (istate)
     {
     case 1:  // prior to initial integration step.
@@ -386,13 +384,13 @@
   Matrix retval;
 
   octave_idx_type n_out = tout.numel ();
-  octave_idx_type n = size ();
+  F77_INT n = to_f77_int (size ());
 
   if (n_out > 0 && n > 0)
     {
       retval.resize (n_out, n);
 
-      for (octave_idx_type i = 0; i < n; i++)
+      for (F77_INT i = 0; i < n; i++)
         retval.elem (0, i) = x.elem (i);
 
       for (octave_idx_type j = 1; j < n_out; j++)
@@ -402,7 +400,7 @@
           if (integration_error)
             return retval;
 
-          for (octave_idx_type i = 0; i < n; i++)
+          for (F77_INT i = 0; i < n; i++)
             retval.elem (j, i) = x_next.elem (i);
         }
     }
@@ -416,13 +414,13 @@
   Matrix retval;
 
   octave_idx_type n_out = tout.numel ();
-  octave_idx_type n = size ();
+  F77_INT n = to_f77_int (size ());
 
   if (n_out > 0 && n > 0)
     {
       retval.resize (n_out, n);
 
-      for (octave_idx_type i = 0; i < n; i++)
+      for (F77_INT i = 0; i < n; i++)
         retval.elem (0, i) = x.elem (i);
 
       octave_idx_type n_crit = tcrit.numel ();
@@ -441,14 +439,14 @@
               if (i_crit < n_crit)
                 next_crit = tcrit.elem (i_crit);
 
-              octave_idx_type save_output;
+              bool save_output = false;
               double t_out;
 
               if (next_crit == next_out)
                 {
                   set_stop_time (next_crit);
                   t_out = next_out;
-                  save_output = 1;
+                  save_output = true;
                   i_out++;
                   i_crit++;
                   do_restart = true;
@@ -459,7 +457,7 @@
                     {
                       set_stop_time (next_crit);
                       t_out = next_crit;
-                      save_output = 0;
+                      save_output = false;
                       i_crit++;
                       do_restart = true;
                     }
@@ -467,7 +465,7 @@
                     {
                       clear_stop_time ();
                       t_out = next_out;
-                      save_output = 1;
+                      save_output = true;
                       i_out++;
                     }
                 }
@@ -475,7 +473,7 @@
                 {
                   set_stop_time (next_crit);
                   t_out = next_out;
-                  save_output = 1;
+                  save_output = true;
                   i_out++;
                 }
 
@@ -486,7 +484,7 @@
 
               if (save_output)
                 {
-                  for (octave_idx_type i = 0; i < n; i++)
+                  for (F77_INT i = 0; i < n; i++)
                     retval.elem (i_out-1, i) = x_next.elem (i);
                 }
 
--- a/liboctave/numeric/LSODE.h	Tue Dec 27 10:24:03 2016 -0500
+++ b/liboctave/numeric/LSODE.h	Tue Dec 27 11:58:30 2016 -0500
@@ -37,12 +37,12 @@
 
   LSODE (void)
     : ODE (), LSODE_options (), initialized (false), method_flag (0),
-      maxord (0), itask (0), iopt (0), itol (0), liw (0), lrw (0),
+      itask (0), iopt (0), itol (0), liw (0), lrw (0),
       iwork (), rwork (), rel_tol (0.0), abs_tol () { }
 
   LSODE (const ColumnVector& s, double tm, const ODEFunc& f)
     : ODE (s, tm, f), LSODE_options (), initialized (false), method_flag (0),
-      maxord (0), itask (0), iopt (0), itol (0), liw (0), lrw (0),
+      itask (0), iopt (0), itol (0), liw (0), lrw (0),
       iwork (), rwork (), rel_tol (0.0), abs_tol () { }
 
   ~LSODE (void) = default;
@@ -59,16 +59,15 @@
 
   bool initialized;
 
-  octave_idx_type method_flag;
-  octave_idx_type maxord;
-  octave_idx_type itask;
-  octave_idx_type iopt;
-  octave_idx_type itol;
+  octave_f77_int_type method_flag;
+  octave_f77_int_type itask;
+  octave_f77_int_type iopt;
+  octave_f77_int_type itol;
 
-  octave_idx_type liw;
-  octave_idx_type lrw;
+  octave_f77_int_type liw;
+  octave_f77_int_type lrw;
 
-  Array<octave_idx_type> iwork;
+  Array<octave_f77_int_type> iwork;
   Array<double> rwork;
 
   double rel_tol;