changeset 29629:93c8df989ea0

qhull: Use reentrant libqhull_r (bug #60016). * configure.ac, m4/acinclude.m4: Check for libqhull_r. Remove obsolete headers from test. * libinterp/dldfcn/oct-qhull.h: Include headers for reentrant QHull. * libinterp/dldfcn/module-files: Use QHULL_R_* flags. * libinterp/dldfcn/__delaunayn__.cc, libinterp/dldfcn/__voronoi__.cc, libinterp/dldfcn/convhulln.cc: Use functions from reentrant QHull. * contributors.in: Add Stefan Brüns to list of Octave contributors.
author Stefan Brüns <stefan.bruens@rwth-aachen.de>
date Fri, 02 Apr 2021 19:53:26 +0200
parents a3df57fa5c5d
children 60315909c884
files configure.ac doc/interpreter/contributors.in libinterp/dldfcn/__delaunayn__.cc libinterp/dldfcn/__voronoi__.cc libinterp/dldfcn/convhulln.cc libinterp/dldfcn/module-files libinterp/dldfcn/oct-qhull.h m4/acinclude.m4
diffstat 8 files changed, 121 insertions(+), 165 deletions(-) [+]
line wrap: on
line diff
--- a/configure.ac	Fri May 07 07:48:00 2021 -0700
+++ b/configure.ac	Fri Apr 02 19:53:26 2021 +0200
@@ -1344,15 +1344,15 @@
 
 ### Check for the Qhull library.
 
-OCTAVE_CHECK_LIB(qhull, QHull,
+OCTAVE_CHECK_LIB(qhull_r, QHull,
   [Qhull library not found.  This will result in loss of functionality for some geometry functions.],
-  [libqhull.h libqhull/libqhull.h qhull/libqhull.h qhull.h qhull/qhull.h],
+  [libqhull_r/libqhull_r.h libqhull_r.h],
   [qh_qhull], [], [],
-  [warn_qhull=
+  [warn_qhull_r=
   OCTAVE_CHECK_QHULL_VERSION
   OCTAVE_CHECK_LIB_QHULL_OK(
     [AC_DEFINE(HAVE_QHULL, 1, [Define to 1 if Qhull is available.])],
-    [warn_qhull="Qhull library found, but does not seem to work properly.  This will result in loss of functionality for some geometry functions.  Please try recompiling the library with -fno-strict-aliasing."])])
+    [warn_qhull_r="Qhull library found, but does not seem to work properly.  This will result in loss of functionality for some geometry functions.  Please try recompiling the library with -fno-strict-aliasing."])])
 
 ### Check for PCRE regex library.
 
@@ -3198,9 +3198,9 @@
   PortAudio libraries:           $PORTAUDIO_LIBS
   PTHREAD flags:                 $PTHREAD_CFLAGS
   PTHREAD libraries:             $PTHREAD_LIBS
-  QHULL CPPFLAGS:                $QHULL_CPPFLAGS
-  QHULL LDFLAGS:                 $QHULL_LDFLAGS
-  QHULL libraries:               $QHULL_LIBS
+  QHULL CPPFLAGS:                $QHULL_R_CPPFLAGS
+  QHULL LDFLAGS:                 $QHULL_R_LDFLAGS
+  QHULL libraries:               $QHULL_R_LIBS
   QRUPDATE CPPFLAGS:             $QRUPDATE_CPPFLAGS
   QRUPDATE LDFLAGS:              $QRUPDATE_LDFLAGS
   QRUPDATE libraries:            $QRUPDATE_LIBS
--- a/doc/interpreter/contributors.in	Fri May 07 07:48:00 2021 -0700
+++ b/doc/interpreter/contributors.in	Fri Apr 02 19:53:26 2021 +0200
@@ -41,6 +41,7 @@
 Marcus Brinkmann
 Max Brister
 Remy Bruno
+Stefan Brüns
 Clemens Buchacher
 Ansgar Burchard
 Antonius Burgers
--- a/libinterp/dldfcn/__delaunayn__.cc	Fri May 07 07:48:00 2021 -0700
+++ b/libinterp/dldfcn/__delaunayn__.cc	Fri Apr 02 19:53:26 2021 +0200
@@ -65,17 +65,17 @@
 
 #  include "oct-qhull.h"
 
-#  if defined (NEED_QHULL_VERSION)
+#  if defined (NEED_QHULL_R_VERSION)
 char qh_version[] = "__delaunayn__.oct 2007-08-21";
 #  endif
 
 static void
-free_qhull_memory ()
+free_qhull_memory (qhT *qh)
 {
-  qh_freeqhull (! qh_ALL);
+  qh_freeqhull (qh, ! qh_ALL);
 
   int curlong, totlong;
-  qh_memfreeshort (&curlong, &totlong);
+  qh_memfreeshort (qh, &curlong, &totlong);
 
   if (curlong || totlong)
     warning ("__delaunayn__: did not free %d bytes of long memory (%d pieces)",
@@ -155,10 +155,7 @@
       double *pt_array = p.fortran_vec ();
       boolT ismalloc = false;
 
-      // Qhull flags argument is not const char*
-      OCTAVE_LOCAL_BUFFER (char, flags, 9 + options.length ());
-
-      sprintf (flags, "qhull d %s", options.c_str ());
+      std::string cmd = "qhull d " + options;
 
       // Set the outfile pointer to stdout for status information.
       FILE *outfile = nullptr;
@@ -176,16 +173,19 @@
 
       octave::unwind_action close_errfile ([=] () { std::fclose (errfile); });
 
-      int exitcode = qh_new_qhull (dim, n, pt_array,
-                                   ismalloc, flags, outfile, errfile);
+      qhT context = { };
+      qhT *qh = &context;
 
-      octave::unwind_action free_memory ([] () { free_qhull_memory (); });
+      int exitcode = qh_new_qhull (qh, dim, n, pt_array, ismalloc, &cmd[0],
+                                   outfile, errfile);
+
+      octave::unwind_action free_memory ([qh] () { free_qhull_memory (qh); });
 
       if (exitcode)
         error ("__delaunayn__: qhull failed");
 
       // triangulate non-simplicial facets
-      qh_triangulate ();
+      qh_triangulate (qh);
 
       facetT *facet;
       vertexT *vertex, **vertexp;
@@ -212,7 +212,7 @@
 
               FOREACHvertex_ (facet->vertices)
                 {
-                  simpl(i, j++) = 1 + qh_pointid(vertex->point);
+                  simpl(i, j++) = 1 + qh_pointid(qh, vertex->point);
                 }
               i++;
             }
--- a/libinterp/dldfcn/__voronoi__.cc	Fri May 07 07:48:00 2021 -0700
+++ b/libinterp/dldfcn/__voronoi__.cc	Fri Apr 02 19:53:26 2021 +0200
@@ -60,17 +60,17 @@
 
 #  include "oct-qhull.h"
 
-#  if defined (NEED_QHULL_VERSION)
+#  if defined (NEED_QHULL_R_VERSION)
 char qh_version[] = "__voronoi__.oct 2007-07-24";
 #  endif
 
 static void
-free_qhull_memory ()
+free_qhull_memory (qhT *qh)
 {
-  qh_freeqhull (! qh_ALL);
+  qh_freeqhull (qh, ! qh_ALL);
 
   int curlong, totlong;
-  qh_memfreeshort (&curlong, &totlong);
+  qh_memfreeshort (qh, &curlong, &totlong);
 
   if (curlong || totlong)
     warning ("__voronoi__: did not free %d bytes of long memory (%d pieces)",
@@ -156,34 +156,32 @@
   FILE *outfile = nullptr;
   FILE *errfile = stderr;
 
-  // qh_new_qhull command and points arguments are not const...
+  qhT context = { };
+  qhT *qh = &context;
+
   std::string cmd = "qhull v" + options;
 
-  OCTAVE_LOCAL_BUFFER (char, cmd_str, cmd.length () + 1);
-
-  strcpy (cmd_str, cmd.c_str ());
+  int exitcode = qh_new_qhull (qh, dim, num_points, points.fortran_vec (),
+                               ismalloc, &cmd[0], outfile, errfile);
 
-  int exitcode = qh_new_qhull (dim, num_points, points.fortran_vec (),
-                               ismalloc, cmd_str, outfile, errfile);
-
-  octave::unwind_action free_memory ([] () { free_qhull_memory (); });
+  octave::unwind_action free_memory ([qh] () { free_qhull_memory (qh); });
 
   if (exitcode)
     error ("%s: qhull failed", caller.c_str ());
 
   // Calling findgood_all provides the number of Voronoi vertices
-  // (sets qh num_good).
+  // (sets qh->num_good).
 
-  qh_findgood_all (qh facet_list);
+  qh_findgood_all (qh, qh->facet_list);
 
   octave_idx_type num_voronoi_regions
-    = qh num_vertices - qh_setsize (qh del_vertices);
+    = qh->num_vertices - qh_setsize (qh, qh->del_vertices);
 
-  octave_idx_type num_voronoi_vertices = qh num_good;
+  octave_idx_type num_voronoi_vertices = qh->num_good;
 
   // Find the voronoi centers for all facets.
 
-  qh_setvoronoi_all ();
+  qh_setvoronoi_all (qh);
 
   facetT *facet;
   vertexT *vertex;
@@ -206,8 +204,8 @@
 
   FORALLvertices
     {
-      if (qh hull_dim == 3)
-        qh_order_vertexneighbors (vertex);
+      if (qh->hull_dim == 3)
+        qh_order_vertexneighbors (qh, vertex);
 
       bool infinity_seen = false;
 
@@ -271,12 +269,12 @@
 
   FORALLvertices
     {
-      if (qh hull_dim == 3)
-        qh_order_vertexneighbors (vertex);
+      if (qh->hull_dim == 3)
+        qh_order_vertexneighbors (qh, vertex);
 
       bool infinity_seen = false;
 
-      octave_idx_type idx = qh_pointid (vertex->point);
+      octave_idx_type idx = qh_pointid (qh, vertex->point);
 
       octave_idx_type num_vertices = ni[k++];
 
--- a/libinterp/dldfcn/convhulln.cc	Fri May 07 07:48:00 2021 -0700
+++ b/libinterp/dldfcn/convhulln.cc	Fri Apr 02 19:53:26 2021 +0200
@@ -53,17 +53,17 @@
 
 #  include "oct-qhull.h"
 
-#  if defined (NEED_QHULL_VERSION)
+#  if defined (NEED_QHULL_R_VERSION)
 char qh_version[] = "convhulln.oct 2007-07-24";
 #  endif
 
 static void
-free_qhull_memory ()
+free_qhull_memory (qhT *qh)
 {
-  qh_freeqhull (! qh_ALL);
+  qh_freeqhull (qh, ! qh_ALL);
 
   int curlong, totlong;
-  qh_memfreeshort (&curlong, &totlong);
+  qh_memfreeshort (qh, &curlong, &totlong);
 
   if (curlong || totlong)
     warning ("convhulln: did not free %d bytes of long memory (%d pieces)",
@@ -171,25 +171,22 @@
   FILE *outfile = nullptr;
   FILE *errfile = stderr;
 
-  // qh_new_qhull command and points arguments are not const...
+  qhT context = { };
+  qhT *qh = &context;
 
   std::string cmd = "qhull" + options;
 
-  OCTAVE_LOCAL_BUFFER (char, cmd_str, cmd.length () + 1);
-
-  strcpy (cmd_str, cmd.c_str ());
-
-  int exitcode = qh_new_qhull (dim, num_points, points.fortran_vec (),
-                               ismalloc, cmd_str, outfile, errfile);
-
-  octave::unwind_action free_memory ([] () { free_qhull_memory (); });
+  int exitcode = qh_new_qhull (qh, dim, num_points, points.fortran_vec (),
+                               ismalloc, &cmd[0], outfile, errfile);
+ 
+  octave::unwind_action free_memory ([qh] () { free_qhull_memory (qh); });
 
   if (exitcode)
     error ("convhulln: qhull failed");
 
   bool nonsimp_seen = false;
 
-  octave_idx_type nf = qh num_facets;
+  octave_idx_type nf = qh->num_facets;
 
   Matrix idx (nf, dim + 1);
 
@@ -201,7 +198,7 @@
     {
       octave_idx_type j = 0;
 
-      if (! (nonsimp_seen || facet->simplicial || qh hull_dim == 2))
+      if (! (nonsimp_seen || facet->simplicial || qh->hull_dim == 2))
         {
           nonsimp_seen = true;
 
@@ -212,14 +209,14 @@
 
       if (dim == 3)
         {
-          setT *vertices = qh_facet3vertex (facet);
+          setT *vertices = qh_facet3vertex (qh, facet);
 
           vertexT *vertex, **vertexp;
 
           FOREACHvertex_ (vertices)
-            idx(i, j++) = 1 + qh_pointid(vertex->point);
+            idx(i, j++) = 1 + qh_pointid(qh, vertex->point);
 
-          qh_settempfree (&vertices);
+          qh_settempfree (qh, &vertices);
         }
       else
         {
@@ -228,14 +225,14 @@
               vertexT *vertex, **vertexp;
 
               FOREACHvertex_ (facet->vertices)
-                idx(i, j++) = 1 + qh_pointid(vertex->point);
+                idx(i, j++) = 1 + qh_pointid(qh, vertex->point);
             }
           else
             {
               vertexT *vertex, **vertexp;
 
               FOREACHvertexreverse12_ (facet->vertices)
-                idx(i, j++) = 1 + qh_pointid(vertex->point);
+                idx(i, j++) = 1 + qh_pointid(qh, vertex->point);
             }
         }
       if (j < dim)
@@ -263,26 +260,26 @@
           if (! facet->normal)
             continue;
 
-          if (facet->upperdelaunay && qh ATinfinity)
+          if (facet->upperdelaunay && qh->ATinfinity)
             continue;
 
-          facet->f.area = area = qh_facetarea (facet);
+          facet->f.area = area = qh_facetarea (qh, facet);
           facet->isarea = True;
 
-          if (qh DELAUNAY)
+          if (qh->DELAUNAY)
             {
-              if (facet->upperdelaunay == qh UPPERdelaunay)
-                qh totarea += area;
+              if (facet->upperdelaunay == qh->UPPERdelaunay)
+                qh->totarea += area;
             }
           else
             {
-              qh totarea += area;
-              qh_distplane (qh interior_point, facet, &dist);
-              qh totvol += -dist * area/ qh hull_dim;
+              qh->totarea += area;
+              qh_distplane (qh, qh->interior_point, facet, &dist);
+              qh->totvol += -dist * area / qh->hull_dim;
             }
         }
 
-      retval(1) = octave_value (qh totvol);
+      retval(1) = octave_value (qh->totvol);
     }
 
   retval(0) = idx;
--- a/libinterp/dldfcn/module-files	Fri May 07 07:48:00 2021 -0700
+++ b/libinterp/dldfcn/module-files	Fri Apr 02 19:53:26 2021 +0200
@@ -1,13 +1,13 @@
 # FILE|CPPFLAGS|LDFLAGS|LIBRARIES
-__delaunayn__.cc|$(QHULL_CPPFLAGS)|$(QHULL_LDFLAGS)|$(QHULL_LIBS)
+__delaunayn__.cc|$(QHULL_R_CPPFLAGS)|$(QHULL_R_LDFLAGS)|$(QHULL_R_LIBS)
 __fltk_uigetfile__.cc|$(FLTK_CPPFLAGS) $(FT2_CPPFLAGS)|$(FLTK_LDFLAGS) $(FT2_LDFLAGS)|$(FLTK_LIBS) $(FT2_LIBS)
 __glpk__.cc|$(GLPK_CPPFLAGS)|$(GLPK_LDFLAGS)|$(GLPK_LIBS)
 __init_fltk__.cc|$(FLTK_CPPFLAGS) $(FT2_CPPFLAGS) $(FONTCONFIG_CPPFLAGS)|$(FLTK_LDFLAGS) $(FT2_LDFLAGS)|$(FLTK_LIBS) $(FT2_LIBS) $(OPENGL_LIBS)
 __init_gnuplot__.cc|$(FT2_CPPFLAGS) $(FONTCONFIG_CPPFLAGS)||
 __ode15__.cc|$(SUNDIALS_XCPPFLAGS)|$(SUNDIALS_XLDFLAGS)|$(SUNDIALS_XLIBS)
-__voronoi__.cc|$(QHULL_CPPFLAGS)|$(QHULL_LDFLAGS)|$(QHULL_LIBS)
+__voronoi__.cc|$(QHULL_R_CPPFLAGS)|$(QHULL_R_LDFLAGS)|$(QHULL_R_LIBS)
 audiodevinfo.cc|$(PORTAUDIO_CPPFLAGS)|$(PORTAUDIO_LDFLAGS)|$(PORTAUDIO_LIBS)
 audioread.cc|$(SNDFILE_CPPFLAGS)|$(SNDFILE_LDFLAGS)|$(SNDFILE_LIBS)
-convhulln.cc|$(QHULL_CPPFLAGS)|$(QHULL_LDFLAGS)|$(QHULL_LIBS)
+convhulln.cc|$(QHULL_R_CPPFLAGS)|$(QHULL_R_LDFLAGS)|$(QHULL_R_LIBS)
 fftw.cc|$(FFTW_XCPPFLAGS)|$(FFTW_XLDFLAGS)|$(FFTW_XLIBS)
 gzip.cc|$(Z_CPPFLAGS) $(BZ2_CPPFLAGS)|$(Z_LDFLAGS) $(BZ2_LDFLAGS)|$(Z_LIBS) $(BZ2_LIBS)
--- a/libinterp/dldfcn/oct-qhull.h	Fri May 07 07:48:00 2021 -0700
+++ b/libinterp/dldfcn/oct-qhull.h	Fri Apr 02 19:53:26 2021 +0200
@@ -32,32 +32,18 @@
 
 extern "C" {
 
-#if defined (HAVE_LIBQHULL_LIBQHULL_H)
-#  include <libqhull/libqhull.h>
-#  include <libqhull/qset.h>
-#  include <libqhull/geom.h>
-#  include <libqhull/poly.h>
-#  include <libqhull/io.h>
-#elif defined (HAVE_QHULL_LIBQHULL_H) || defined (HAVE_QHULL_QHULL_H)
-#  if defined (HAVE_QHULL_LIBQHULL_H)
-#    include <qhull/libqhull.h>
-#  else
-#    include <qhull/qhull.h>
-#  endif
-#  include <qhull/qset.h>
-#  include <qhull/geom.h>
-#  include <qhull/poly.h>
-#  include <qhull/io.h>
-#  elif defined (HAVE_LIBQHULL_H) || defined (HAVE_QHULL_H)
-#  if defined (HAVE_LIBQHULL_H)
-#    include <libqhull.h>
-#  else
-#    include <qhull.h>
-#  endif
-#  include <qset.h>
-#  include <geom.h>
-#  include <poly.h>
-#  include <io.h>
+#if defined (HAVE_LIBQHULL_R_LIBQHULL_R_H)
+#  include <libqhull_r/libqhull_r.h>
+#  include <libqhull_r/qset_r.h>
+#  include <libqhull_r/geom_r.h>
+#  include <libqhull_r/poly_r.h>
+#  include <libqhull_r/io_r.h>
+#elif defined (HAVE_LIBQHULL_R_H)
+#  include <libqhull_r.h>
+#  include <qset_r.h>
+#  include <geom_r.h>
+#  include <poly_r.h>
+#  include <io_r.h>
 #endif
 
 }
--- a/m4/acinclude.m4	Fri May 07 07:48:00 2021 -0700
+++ b/m4/acinclude.m4	Fri Apr 02 19:53:26 2021 +0200
@@ -1436,38 +1436,24 @@
 dnl Check whether Qhull works (does not crash).
 dnl
 AC_DEFUN([OCTAVE_CHECK_LIB_QHULL_OK], [
-  AC_CACHE_CHECK([whether the qhull library works],
-    [octave_cv_lib_qhull_ok],
+  AC_CACHE_CHECK([whether the qhull_r library works],
+    [octave_cv_lib_qhull_r_ok],
     [AC_RUN_IFELSE([AC_LANG_PROGRAM([[
         #include <stdio.h>
-        #if defined (HAVE_LIBQHULL_LIBQHULL_H)
-        # include <libqhull/libqhull.h>
-        # include <libqhull/qset.h>
-        # include <libqhull/geom.h>
-        # include <libqhull/poly.h>
-        # include <libqhull/io.h>
-        #elif defined (HAVE_QHULL_LIBQHULL_H) || defined (HAVE_QHULL_QHULL_H)
-        # if defined (HAVE_QHULL_LIBQHULL_H)
-        #  include <qhull/libqhull.h>
-        # else
-        #  include <qhull/qhull.h>
-        # endif
-        # include <qhull/qset.h>
-        # include <qhull/geom.h>
-        # include <qhull/poly.h>
-        # include <qhull/io.h>
-        #elif defined (HAVE_LIBQHULL_H) || defined (HAVE_QHULL_H)
-        # if defined (HAVE_LIBQHULL_H)
-        #  include <libqhull.h>
-        # else
-        #  include <qhull.h>
-        # endif
-        # include <qset.h>
-        # include <geom.h>
-        # include <poly.h>
-        # include <io.h>
+        #if defined (HAVE_LIBQHULL_R_LIBQHULL_R_H)
+        # include <libqhull_r/libqhull_r.h>
+        # include <libqhull_r/qset_r.h>
+        # include <libqhull_r/geom_r.h>
+        # include <libqhull_r/poly_r.h>
+        # include <libqhull_r/io_r.h>
+        #elif defined (HAVE_LIBQHULL_R_H)
+        # include <libqhull_r.h>
+        # include <qset_r.h>
+        # include <geom_r.h>
+        # include <poly_r.h>
+        # include <io_r.h>
         #endif
-        #if defined (NEED_QHULL_VERSION)
+        #if defined (NEED_QHULL_R_VERSION)
           char *qh_version = "version";
         #endif
         ]], [[
@@ -1475,13 +1461,15 @@
         int n = 4;
         coordT points[8] = { -0.5, -0.5, -0.5, 0.5, 0.5, -0.5, 0.5, 0.5 };
         boolT ismalloc = 0;
-        return qh_new_qhull (dim, n, points, ismalloc, "qhull ", 0, stderr);
+        qhT context;
+        qhT* qh = &context;
+        return qh_new_qhull (qh, dim, n, points, ismalloc, "qhull ", 0, stderr);
       ]])],
-      octave_cv_lib_qhull_ok=yes,
-      octave_cv_lib_qhull_ok=no,
-      octave_cv_lib_qhull_ok=yes)
+      octave_cv_lib_qhull_r_ok=yes,
+      octave_cv_lib_qhull_r_ok=no,
+      octave_cv_lib_qhull_r_ok=yes)
   ])
-  if test $octave_cv_lib_qhull_ok = yes; then
+  if test $octave_cv_lib_qhull_r_ok = yes; then
     $1
     :
   else
@@ -1578,44 +1566,30 @@
 dnl Check for the Qhull version.
 dnl
 AC_DEFUN([OCTAVE_CHECK_QHULL_VERSION], [
-  AC_CACHE_CHECK([for qh_version in $QHULL_LIBS],
-    [octave_cv_lib_qhull_version],
+  AC_CACHE_CHECK([for qh_version in $QHULL_R_LIBS],
+    [octave_cv_lib_qhull_r_version],
     [AC_LINK_IFELSE([AC_LANG_PROGRAM([[
         #include <stdio.h>
-        #if defined (HAVE_LIBQHULL_LIBQHULL_H)
-        # include <libqhull/libqhull.h>
-        # include <libqhull/qset.h>
-        # include <libqhull/geom.h>
-        # include <libqhull/poly.h>
-        # include <libqhull/io.h>
-        #elif defined (HAVE_QHULL_LIBQHULL_H) || defined (HAVE_QHULL_QHULL_H)
-        # if defined (HAVE_QHULL_LIBQHULL_H)
-        #  include <qhull/libqhull.h>
-        # else
-        #  include <qhull/qhull.h>
-        # endif
-        # include <qhull/qset.h>
-        # include <qhull/geom.h>
-        # include <qhull/poly.h>
-        # include <qhull/io.h>
-        #elif defined (HAVE_LIBQHULL_H) || defined (HAVE_QHULL_H)
-        # if defined (HAVE_LIBQHULL_H)
-        #  include <libqhull.h>
-        # else
-        #  include <qhull.h>
-        # endif
-        # include <qset.h>
-        # include <geom.h>
-        # include <poly.h>
-        # include <io.h>
+        #if defined (HAVE_LIBQHULL_R_LIBQHULL_R_H)
+        # include <libqhull_r/libqhull_r.h>
+        # include <libqhull_r/qset_r.h>
+        # include <libqhull_r/geom_r.h>
+        # include <libqhull_r/poly_r.h>
+        # include <libqhull_r/io_r.h>
+        #elif defined (HAVE_LIBQHULL_R_H)
+        # include <libqhull_r.h>
+        # include <qset_r.h>
+        # include <geom_r.h>
+        # include <poly_r.h>
+        # include <io_r.h>
         #endif
         ]], [[
         const char *tmp = qh_version;
       ]])],
-      octave_cv_lib_qhull_version=yes, octave_cv_lib_qhull_version=no)
+      octave_cv_lib_qhull_r_version=yes, octave_cv_lib_qhull_r_version=no)
   ])
-  if test $octave_cv_lib_qhull_version = no; then
-    AC_DEFINE(NEED_QHULL_VERSION, 1,
+  if test $octave_cv_lib_qhull_r_version = no; then
+    AC_DEFINE(NEED_QHULL_R_VERSION, 1,
       [Define to 1 if the Qhull library needs a qh_version variable defined.])
   fi
 ])