changeset 29876:89bdb44db76f

move liboctave convn functions inside octave namespace * oct-convn.h (convn, convn_type): Move inside octave namespace. Also preserve global convn_type enum. Provide deprecated wrappers for old names. Change uses where needed. (convert_enum): New function. * oct-convn.cc: Move templates and all functions inside octave namespace.
author John W. Eaton <jwe@octave.org>
date Wed, 14 Jul 2021 00:18:49 -0400
parents 9f7132a05682
children 86393208b363
files libinterp/corefcn/conv2.cc liboctave/numeric/oct-convn.cc liboctave/numeric/oct-convn.h
diffstat 3 files changed, 366 insertions(+), 171 deletions(-) [+]
line wrap: on
line diff
--- a/libinterp/corefcn/conv2.cc	Tue Jul 13 23:14:33 2021 -0400
+++ b/libinterp/corefcn/conv2.cc	Wed Jul 14 00:18:49 2021 -0400
@@ -73,7 +73,7 @@
 
   std::string shape = "full";   // default
   bool separable = false;
-  convn_type ct = convn_full;
+  octave::convn_type ct = octave::convn_full;
 
   if (nargin == 3)
     {
@@ -92,11 +92,11 @@
     error ("conv2: A and B must be 1-D vectors or 2-D matrices");
 
   if (shape == "full")
-    ct = convn_full;
+    ct = octave::convn_full;
   else if (shape == "same")
-    ct = convn_same;
+    ct = octave::convn_same;
   else if (shape == "valid")
-    ct = convn_valid;
+    ct = octave::convn_valid;
   else
     error ("conv2: SHAPE type not valid");
 
@@ -120,13 +120,13 @@
                 {
                   FloatColumnVector v1 (args(0).float_vector_value ());
                   FloatRowVector v2 (args(1).float_vector_value ());
-                  retval = convn (a, v1, v2, ct);
+                  retval = octave::convn (a, v1, v2, ct);
                 }
               else
                 {
                   FloatComplexColumnVector v1 (args(0).float_complex_vector_value ());
                   FloatComplexRowVector v2 (args(1).float_complex_vector_value ());
-                  retval = convn (a, v1, v2, ct);
+                  retval = octave::convn (a, v1, v2, ct);
                 }
             }
           else
@@ -134,7 +134,7 @@
               FloatColumnVector v1 (args(0).float_vector_value ());
               FloatRowVector v2 (args(1).float_vector_value ());
               FloatMatrix a (args(2).float_matrix_value ());
-              retval = convn (a, v1, v2, ct);
+              retval = octave::convn (a, v1, v2, ct);
             }
         }
       else
@@ -147,13 +147,13 @@
                 {
                   ColumnVector v1 (args(0).vector_value ());
                   RowVector v2 (args(1).vector_value ());
-                  retval = convn (a, v1, v2, ct);
+                  retval = octave::convn (a, v1, v2, ct);
                 }
               else
                 {
                   ComplexColumnVector v1 (args(0).complex_vector_value ());
                   ComplexRowVector v2 (args(1).complex_vector_value ());
-                  retval = convn (a, v1, v2, ct);
+                  retval = octave::convn (a, v1, v2, ct);
                 }
             }
           else
@@ -161,7 +161,7 @@
               ColumnVector v1 (args(0).vector_value ());
               RowVector v2 (args(1).vector_value ());
               Matrix a (args(2).matrix_value ());
-              retval = convn (a, v1, v2, ct);
+              retval = octave::convn (a, v1, v2, ct);
             }
         }
     } // if (separable)
@@ -175,19 +175,19 @@
               if (args(1).isreal ())
                 {
                   FloatMatrix b (args(1).float_matrix_value ());
-                  retval = convn (a, b, ct);
+                  retval = octave::convn (a, b, ct);
                 }
               else
                 {
                   FloatComplexMatrix b (args(1).float_complex_matrix_value ());
-                  retval = convn (a, b, ct);
+                  retval = octave::convn (a, b, ct);
                 }
             }
           else
             {
               FloatMatrix a (args(0).float_matrix_value ());
               FloatMatrix b (args(1).float_matrix_value ());
-              retval = convn (a, b, ct);
+              retval = octave::convn (a, b, ct);
             }
         }
       else
@@ -198,19 +198,19 @@
               if (args(1).isreal ())
                 {
                   Matrix b (args(1).matrix_value ());
-                  retval = convn (a, b, ct);
+                  retval = octave::convn (a, b, ct);
                 }
               else
                 {
                   ComplexMatrix b (args(1).complex_matrix_value ());
-                  retval = convn (a, b, ct);
+                  retval = octave::convn (a, b, ct);
                 }
             }
           else
             {
               Matrix a (args(0).matrix_value ());
               Matrix b (args(1).matrix_value ());
-              retval = convn (a, b, ct);
+              retval = octave::convn (a, b, ct);
             }
         }
 
@@ -332,17 +332,17 @@
     print_usage ();
 
   std::string shape = "full";   // default
-  convn_type ct = convn_full;
+  octave::convn_type ct = octave::convn_full;
 
   if (nargin == 3)
     shape = args(2).xstring_value ("convn: SHAPE must be a string");
 
   if (shape == "full")
-    ct = convn_full;
+    ct = octave::convn_full;
   else if (shape == "same")
-    ct = convn_same;
+    ct = octave::convn_same;
   else if (shape == "valid")
-    ct = convn_valid;
+    ct = octave::convn_valid;
   else
     error ("convn: SHAPE type not valid");
 
@@ -356,19 +356,19 @@
           if (args(1).isreal ())
             {
               FloatNDArray b (args(1).float_array_value ());
-              retval = convn (a, b, ct);
+              retval = octave::convn (a, b, ct);
             }
           else
             {
               FloatComplexNDArray b (args(1).float_complex_array_value ());
-              retval = convn (a, b, ct);
+              retval = octave::convn (a, b, ct);
             }
         }
       else
         {
           FloatNDArray a (args(0).float_array_value ());
           FloatNDArray b (args(1).float_array_value ());
-          retval = convn (a, b, ct);
+          retval = octave::convn (a, b, ct);
         }
     }
   else
@@ -379,19 +379,19 @@
           if (args(1).isreal ())
             {
               NDArray b (args(1).array_value ());
-              retval = convn (a, b, ct);
+              retval = octave::convn (a, b, ct);
             }
           else
             {
               ComplexNDArray b (args(1).complex_array_value ());
-              retval = convn (a, b, ct);
+              retval = octave::convn (a, b, ct);
             }
         }
       else
         {
           NDArray a (args(0).array_value ());
           NDArray b (args(1).array_value ());
-          retval = convn (a, b, ct);
+          retval = octave::convn (a, b, ct);
         }
     }
 
--- a/liboctave/numeric/oct-convn.cc	Tue Jul 13 23:14:33 2021 -0400
+++ b/liboctave/numeric/oct-convn.cc	Wed Jul 14 00:18:49 2021 -0400
@@ -50,14 +50,16 @@
 #include "fRowVector.h"
 #include "oct-convn.h"
 
-// 2d convolution with a matrix kernel.
-template <typename T, typename R>
-static void
-convolve_2d (const T *a, F77_INT ma, F77_INT na,
-             const R *b, F77_INT mb, F77_INT nb,
-             T *c, bool inner);
+namespace octave
+{
+  // 2d convolution with a matrix kernel.
+  template <typename T, typename R>
+  static void
+  convolve_2d (const T *a, F77_INT ma, F77_INT na,
+               const R *b, F77_INT mb, F77_INT nb,
+               T *c, bool inner);
 
-// Forward instances to our Fortran implementations.
+  // Forward instances to our Fortran implementations.
 #define FORWARD_IMPL(T_CXX, R_CXX, T, R, T_CAST, T_CONST_CAST,          \
                      R_CONST_CAST, f, F)                                \
   extern "C"                                                            \
@@ -87,112 +89,112 @@
                                        T_CAST (c)));                    \
   }
 
-FORWARD_IMPL (double, double, F77_DBLE, F77_DBLE, , , , d, D)
-FORWARD_IMPL (float, float, F77_REAL, F77_REAL, , , , s, S)
+  FORWARD_IMPL (double, double, F77_DBLE, F77_DBLE, , , , d, D)
+  FORWARD_IMPL (float, float, F77_REAL, F77_REAL, , , , s, S)
 
-FORWARD_IMPL (std::complex<double>, std::complex<double>,
-              F77_DBLE_CMPLX, F77_DBLE_CMPLX, F77_DBLE_CMPLX_ARG,
-              F77_CONST_DBLE_CMPLX_ARG, F77_CONST_DBLE_CMPLX_ARG, z, Z)
-FORWARD_IMPL (std::complex<float>, std::complex<float>,
-              F77_CMPLX, F77_CMPLX, F77_CMPLX_ARG,
-              F77_CONST_CMPLX_ARG, F77_CONST_CMPLX_ARG, c, C)
+  FORWARD_IMPL (std::complex<double>, std::complex<double>,
+                F77_DBLE_CMPLX, F77_DBLE_CMPLX, F77_DBLE_CMPLX_ARG,
+                F77_CONST_DBLE_CMPLX_ARG, F77_CONST_DBLE_CMPLX_ARG, z, Z)
+  FORWARD_IMPL (std::complex<float>, std::complex<float>,
+                F77_CMPLX, F77_CMPLX, F77_CMPLX_ARG,
+                F77_CONST_CMPLX_ARG, F77_CONST_CMPLX_ARG, c, C)
 
-FORWARD_IMPL (std::complex<double>, double,
-              F77_DBLE_CMPLX, F77_DBLE, F77_DBLE_CMPLX_ARG,
-              F77_CONST_DBLE_CMPLX_ARG, , zd, ZD)
-FORWARD_IMPL (std::complex<float>, float, F77_CMPLX, F77_REAL, F77_CMPLX_ARG,
-              F77_CONST_CMPLX_ARG, , cs, CS)
+  FORWARD_IMPL (std::complex<double>, double,
+                F77_DBLE_CMPLX, F77_DBLE, F77_DBLE_CMPLX_ARG,
+                F77_CONST_DBLE_CMPLX_ARG, , zd, ZD)
+  FORWARD_IMPL (std::complex<float>, float, F77_CMPLX, F77_REAL, F77_CMPLX_ARG,
+                F77_CONST_CMPLX_ARG, , cs, CS)
 
-template <typename T, typename R>
-void convolve_nd (const T *a, const dim_vector& ad, const dim_vector& acd,
-                  const R *b, const dim_vector& bd, const dim_vector& bcd,
-                  T *c, const dim_vector& ccd, int nd, bool inner)
-{
-  if (nd == 2)
-    {
-      F77_INT ad0 = octave::to_f77_int (ad(0));
-      F77_INT ad1 = octave::to_f77_int (ad(1));
+  template <typename T, typename R>
+  void convolve_nd (const T *a, const dim_vector& ad, const dim_vector& acd,
+                    const R *b, const dim_vector& bd, const dim_vector& bcd,
+                    T *c, const dim_vector& ccd, int nd, bool inner)
+  {
+    if (nd == 2)
+      {
+        F77_INT ad0 = octave::to_f77_int (ad(0));
+        F77_INT ad1 = octave::to_f77_int (ad(1));
 
-      F77_INT bd0 = octave::to_f77_int (bd(0));
-      F77_INT bd1 = octave::to_f77_int (bd(1));
+        F77_INT bd0 = octave::to_f77_int (bd(0));
+        F77_INT bd1 = octave::to_f77_int (bd(1));
 
-      convolve_2d<T, R> (a, ad0, ad1, b, bd0, bd1, c, inner);
-    }
-  else
-    {
-      octave_idx_type ma = acd(nd-2);
-      octave_idx_type na = ad(nd-1);
-      octave_idx_type mb = bcd(nd-2);
-      octave_idx_type nb = bd(nd-1);
-      octave_idx_type ldc = ccd(nd-2);
+        convolve_2d<T, R> (a, ad0, ad1, b, bd0, bd1, c, inner);
+      }
+    else
+      {
+        octave_idx_type ma = acd(nd-2);
+        octave_idx_type na = ad(nd-1);
+        octave_idx_type mb = bcd(nd-2);
+        octave_idx_type nb = bd(nd-1);
+        octave_idx_type ldc = ccd(nd-2);
 
-      if (inner)
-        {
-          for (octave_idx_type ja = 0; ja < na - nb + 1; ja++)
-            for (octave_idx_type jb = 0; jb < nb; jb++)
-              convolve_nd<T, R> (a + ma*(ja+jb), ad, acd,
-                                 b + mb*(nb-jb-1), bd, bcd,
-                                 c + ldc*ja, ccd, nd-1, inner);
-        }
-      else
-        {
-          for (octave_idx_type ja = 0; ja < na; ja++)
-            for (octave_idx_type jb = 0; jb < nb; jb++)
-              convolve_nd<T, R> (a + ma*ja, ad, acd, b + mb*jb, bd, bcd,
-                                 c + ldc*(ja+jb), ccd, nd-1, inner);
-        }
-    }
-}
+        if (inner)
+          {
+            for (octave_idx_type ja = 0; ja < na - nb + 1; ja++)
+              for (octave_idx_type jb = 0; jb < nb; jb++)
+                convolve_nd<T, R> (a + ma*(ja+jb), ad, acd,
+                                   b + mb*(nb-jb-1), bd, bcd,
+                                   c + ldc*ja, ccd, nd-1, inner);
+          }
+        else
+          {
+            for (octave_idx_type ja = 0; ja < na; ja++)
+              for (octave_idx_type jb = 0; jb < nb; jb++)
+                convolve_nd<T, R> (a + ma*ja, ad, acd, b + mb*jb, bd, bcd,
+                                   c + ldc*(ja+jb), ccd, nd-1, inner);
+          }
+      }
+  }
 
-// Arbitrary convolutor.
-// The 2nd array is assumed to be the smaller one.
-template <typename T, typename R>
-static MArray<T>
-convolve (const MArray<T>& a, const MArray<R>& b,
-          convn_type ct)
-{
-  if (a.isempty () || b.isempty ())
-    return MArray<T> ();
+  // Arbitrary convolutor.
+  // The 2nd array is assumed to be the smaller one.
+  template <typename T, typename R>
+  static MArray<T>
+  convolve (const MArray<T>& a, const MArray<R>& b,
+            convn_type ct)
+  {
+    if (a.isempty () || b.isempty ())
+      return MArray<T> ();
 
-  int nd = std::max (a.ndims (), b.ndims ());
-  const dim_vector adims = a.dims ().redim (nd);
-  const dim_vector bdims = b.dims ().redim (nd);
-  dim_vector cdims = dim_vector::alloc (nd);
+    int nd = std::max (a.ndims (), b.ndims ());
+    const dim_vector adims = a.dims ().redim (nd);
+    const dim_vector bdims = b.dims ().redim (nd);
+    dim_vector cdims = dim_vector::alloc (nd);
 
-  for (int i = 0; i < nd; i++)
-    {
-      if (ct == convn_valid)
-        cdims(i) = std::max (adims(i) - bdims(i) + 1,
-                             static_cast<octave_idx_type> (0));
-      else
-        cdims(i) = std::max (adims(i) + bdims(i) - 1,
-                             static_cast<octave_idx_type> (0));
-    }
+    for (int i = 0; i < nd; i++)
+      {
+        if (ct == convn_valid)
+          cdims(i) = std::max (adims(i) - bdims(i) + 1,
+                               static_cast<octave_idx_type> (0));
+        else
+          cdims(i) = std::max (adims(i) + bdims(i) - 1,
+                               static_cast<octave_idx_type> (0));
+      }
 
-  MArray<T> c (cdims, T ());
+    MArray<T> c (cdims, T ());
 
-  // "valid" shape can sometimes result in empty matrices which must avoid
-  // calling Fortran code which does not expect this (bug #52067)
-  if (c.isempty ())
-    return c;
+    // "valid" shape can sometimes result in empty matrices which must avoid
+    // calling Fortran code which does not expect this (bug #52067)
+    if (c.isempty ())
+      return c;
 
-  convolve_nd<T, R> (a.fortran_vec (), adims, adims.cumulative (),
-                     b.fortran_vec (), bdims, bdims.cumulative (),
-                     c.fortran_vec (), cdims.cumulative (),
-                     nd, ct == convn_valid);
+    convolve_nd<T, R> (a.fortran_vec (), adims, adims.cumulative (),
+                       b.fortran_vec (), bdims, bdims.cumulative (),
+                       c.fortran_vec (), cdims.cumulative (),
+                       nd, ct == convn_valid);
 
-  if (ct == convn_same)
-    {
-      // Pick the relevant part.
-      Array<octave::idx_vector> sidx (dim_vector (nd, 1));
+    if (ct == convn_same)
+      {
+        // Pick the relevant part.
+        Array<octave::idx_vector> sidx (dim_vector (nd, 1));
 
-      for (int i = 0; i < nd; i++)
-        sidx(i) = octave::idx_vector::make_range (bdims(i)/2, 1, adims(i));
-      c = c.index (sidx);
-    }
+        for (int i = 0; i < nd; i++)
+          sidx(i) = octave::idx_vector::make_range (bdims(i)/2, 1, adims(i));
+        c = c.index (sidx);
+      }
 
-  return c;
-}
+    return c;
+  }
 
 #define CONV_DEFS(TPREF, RPREF)                                         \
   TPREF ## NDArray                                                      \
@@ -214,9 +216,10 @@
     return convolve (a, c * r, ct);                                     \
   }
 
-CONV_DEFS ( , )
-CONV_DEFS (Complex, )
-CONV_DEFS (Complex, Complex)
-CONV_DEFS (Float, Float)
-CONV_DEFS (FloatComplex, Float)
-CONV_DEFS (FloatComplex, FloatComplex)
+  CONV_DEFS ( , )
+  CONV_DEFS (Complex, )
+  CONV_DEFS (Complex, Complex)
+  CONV_DEFS (Float, Float)
+  CONV_DEFS (FloatComplex, Float)
+  CONV_DEFS (FloatComplex, FloatComplex)
+}
--- a/liboctave/numeric/oct-convn.h	Tue Jul 13 23:14:33 2021 -0400
+++ b/liboctave/numeric/oct-convn.h	Wed Jul 14 00:18:49 2021 -0400
@@ -48,85 +48,277 @@
 class FloatComplexMatrix;
 class FloatComplexNDArray;
 
+// The remaining includes can be removed when the global enum
+// declaration, the convert_enum function, and the deprecated functions
+// at the end of this file are removed.
+
+#include <cstdlib>
+#include "dMatrix.h"
+#include "dNDArray.h"
+#include "CMatrix.h"
+#include "CNDArray.h"
+#include "fMatrix.h"
+#include "fNDArray.h"
+#include "fCMatrix.h"
+#include "fCNDArray.h"
+
+// FIXME: Is there any sane way to move a global enum to a namespace and
+// tag the global one as deprecated when it is used as a parameter in
+// public functions that also need to be tagged as deprecated?
+
 enum convn_type
+  {
+    convn_full,
+    convn_same,
+    convn_valid
+  };
+
+namespace octave
 {
-  convn_full,
-  convn_same,
-  convn_valid
-};
+  enum convn_type
+    {
+      convn_full,
+      convn_same,
+      convn_valid
+    };
+
+  // double real X double real
+
+  extern OCTAVE_API NDArray
+  convn (const NDArray& a, const NDArray& b, convn_type ct);
+
+  extern OCTAVE_API Matrix
+  convn (const Matrix& a, const Matrix& b, convn_type ct);
+
+  extern OCTAVE_API Matrix
+  convn (const Matrix& a, const ColumnVector& c, const RowVector& r,
+         convn_type ct);
+
+  // double complex X double real
 
-// double real X double real
+  extern OCTAVE_API ComplexNDArray
+  convn (const ComplexNDArray& a, const NDArray& b, convn_type ct);
+
+  extern OCTAVE_API ComplexMatrix
+  convn (const ComplexMatrix& a, const Matrix& b, convn_type ct);
+
+  extern OCTAVE_API ComplexMatrix
+  convn (const ComplexMatrix& a, const ColumnVector& c, const RowVector& r,
+         convn_type ct);
+
+  // double complex X double complex
+
+  extern OCTAVE_API ComplexNDArray
+  convn (const ComplexNDArray& a, const ComplexNDArray& b, convn_type ct);
+
+  extern OCTAVE_API ComplexMatrix
+  convn (const ComplexMatrix& a, const ComplexMatrix& b, convn_type ct);
+
+  extern OCTAVE_API ComplexMatrix
+  convn (const ComplexMatrix& a, const ComplexColumnVector& c,
+         const ComplexRowVector& r, convn_type ct);
+
+  // float real X float real
 
-extern OCTAVE_API NDArray
-convn (const NDArray& a, const NDArray& b, convn_type ct);
+  extern OCTAVE_API FloatNDArray
+  convn (const FloatNDArray& a, const FloatNDArray& b, convn_type ct);
+
+  extern OCTAVE_API FloatMatrix
+  convn (const FloatMatrix& a, const FloatMatrix& b, convn_type ct);
+
+  extern OCTAVE_API FloatMatrix
+  convn (const FloatMatrix& a, const FloatColumnVector& c,
+         const FloatRowVector& r, convn_type ct);
+
+  // float complex X float real
+
+  extern OCTAVE_API FloatComplexNDArray
+  convn (const FloatComplexNDArray& a, const FloatNDArray& b, convn_type ct);
+
+  extern OCTAVE_API FloatComplexMatrix
+  convn (const FloatComplexMatrix& a, const FloatMatrix& b, convn_type ct);
+
+  extern OCTAVE_API FloatComplexMatrix
+  convn (const FloatComplexMatrix& a, const FloatColumnVector& c,
+         const FloatRowVector& r, convn_type ct);
+
+  // float complex X float complex
+
+  extern OCTAVE_API FloatComplexNDArray
+  convn (const FloatComplexNDArray& a, const FloatComplexNDArray& b,
+         convn_type ct);
+
+  extern OCTAVE_API FloatComplexMatrix
+  convn (const FloatComplexMatrix& a, const FloatComplexMatrix& b,
+         convn_type ct);
+
+  extern OCTAVE_API FloatComplexMatrix
+  convn (const FloatComplexMatrix& a, const FloatComplexColumnVector& c,
+         const FloatComplexRowVector& r, convn_type ct);
 
-extern OCTAVE_API Matrix
-convn (const Matrix& a, const Matrix& b, convn_type ct);
+  convn_type convert_enum (::convn_type ct)
+  {
+    switch (ct)
+      {
+      case ::convn_full:
+        return convn_full;
+
+      case ::convn_same:
+        return convn_same;
+
+      case ::convn_valid:
+        return convn_valid;
+
+      default:
+        abort ();
+      }
+  }
+}
 
-extern OCTAVE_API Matrix
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline NDArray
+convn (const NDArray& a, const NDArray& b, convn_type ct)
+{
+  return octave::convn (a, b, static_cast<octave::convn_type> (ct));
+}
+
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline Matrix
+convn (const Matrix& a, const Matrix& b, convn_type ct)
+{
+  return octave::convn (a, b, octave::convert_enum (ct));
+}
+
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline Matrix
 convn (const Matrix& a, const ColumnVector& c, const RowVector& r,
-       convn_type ct);
+       convn_type ct)
+{
+  return octave::convn (a, c, r, octave::convert_enum (ct));
+}
 
 // double complex X double real
 
-extern OCTAVE_API ComplexNDArray
-convn (const ComplexNDArray& a, const NDArray& b, convn_type ct);
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline ComplexNDArray
+convn (const ComplexNDArray& a, const NDArray& b, convn_type ct)
+{
+  return octave::convn (a, b, octave::convert_enum (ct));
+}
 
-extern OCTAVE_API ComplexMatrix
-convn (const ComplexMatrix& a, const Matrix& b, convn_type ct);
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline ComplexMatrix
+convn (const ComplexMatrix& a, const Matrix& b, convn_type ct)
+{
+  return octave::convn (a, b, octave::convert_enum (ct));
+}
 
-extern OCTAVE_API ComplexMatrix
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline ComplexMatrix
 convn (const ComplexMatrix& a, const ColumnVector& c, const RowVector& r,
-       convn_type ct);
+       convn_type ct)
+{
+  return octave::convn (a, c, r, octave::convert_enum (ct));
+}
 
 // double complex X double complex
 
-extern OCTAVE_API ComplexNDArray
-convn (const ComplexNDArray& a, const ComplexNDArray& b, convn_type ct);
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline ComplexNDArray
+convn (const ComplexNDArray& a, const ComplexNDArray& b, convn_type ct)
+{
+  return octave::convn (a, b, octave::convert_enum (ct));
+}
 
-extern OCTAVE_API ComplexMatrix
-convn (const ComplexMatrix& a, const ComplexMatrix& b, convn_type ct);
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline ComplexMatrix
+convn (const ComplexMatrix& a, const ComplexMatrix& b, convn_type ct)
+{
+  return octave::convn (a, b, octave::convert_enum (ct));
+}
 
-extern OCTAVE_API ComplexMatrix
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline ComplexMatrix
 convn (const ComplexMatrix& a, const ComplexColumnVector& c,
-       const ComplexRowVector& r, convn_type ct);
+       const ComplexRowVector& r, convn_type ct)
+{
+  return octave::convn (a, c, r, octave::convert_enum (ct));
+}
 
 // float real X float real
 
-extern OCTAVE_API FloatNDArray
-convn (const FloatNDArray& a, const FloatNDArray& b, convn_type ct);
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline FloatNDArray
+convn (const FloatNDArray& a, const FloatNDArray& b, convn_type ct)
+{
+  return octave::convn (a, b, octave::convert_enum (ct));
+}
 
-extern OCTAVE_API FloatMatrix
-convn (const FloatMatrix& a, const FloatMatrix& b, convn_type ct);
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline FloatMatrix
+convn (const FloatMatrix& a, const FloatMatrix& b, convn_type ct)
+{
+  return octave::convn (a, b, octave::convert_enum (ct));
+}
 
-extern OCTAVE_API FloatMatrix
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline FloatMatrix
 convn (const FloatMatrix& a, const FloatColumnVector& c,
-       const FloatRowVector& r, convn_type ct);
+       const FloatRowVector& r, convn_type ct)
+{
+  return octave::convn (a, c, r, octave::convert_enum (ct));
+}
 
 // float complex X float real
 
-extern OCTAVE_API FloatComplexNDArray
-convn (const FloatComplexNDArray& a, const FloatNDArray& b, convn_type ct);
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline FloatComplexNDArray
+convn (const FloatComplexNDArray& a, const FloatNDArray& b,
+       convn_type ct)
+{
+  return octave::convn (a, b, octave::convert_enum (ct));
+}
 
-extern OCTAVE_API FloatComplexMatrix
-convn (const FloatComplexMatrix& a, const FloatMatrix& b, convn_type ct);
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline FloatComplexMatrix
+convn (const FloatComplexMatrix& a, const FloatMatrix& b,
+       convn_type ct)
+{
+  return octave::convn (a, b, octave::convert_enum (ct));
+}
 
-extern OCTAVE_API FloatComplexMatrix
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline FloatComplexMatrix
 convn (const FloatComplexMatrix& a, const FloatColumnVector& c,
-       const FloatRowVector& r, convn_type ct);
+       const FloatRowVector& r, convn_type ct)
+{
+  return octave::convn (a, c, r, octave::convert_enum (ct));
+}
 
 // float complex X float complex
 
-extern OCTAVE_API FloatComplexNDArray
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline FloatComplexNDArray
 convn (const FloatComplexNDArray& a, const FloatComplexNDArray& b,
-       convn_type ct);
+       convn_type ct)
+{
+  return octave::convn (a, b, octave::convert_enum (ct));
+}
 
-extern OCTAVE_API FloatComplexMatrix
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline FloatComplexMatrix
 convn (const FloatComplexMatrix& a, const FloatComplexMatrix& b,
-       convn_type ct);
+       convn_type ct)
+{
+  return octave::convn (a, b, octave::convert_enum (ct));
+}
 
-extern OCTAVE_API FloatComplexMatrix
+OCTAVE_DEPRECATED (7, "use 'octave::convn' instead")
+inline FloatComplexMatrix
 convn (const FloatComplexMatrix& a, const FloatComplexColumnVector& c,
-       const FloatComplexRowVector& r, convn_type ct);
+       const FloatComplexRowVector& r, convn_type ct)
+{
+  return octave::convn (a, c, r, octave::convert_enum (ct));
+}
 
 #endif