GNU Octave 10.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 
Loading...
Searching...
No Matches
__delaunayn__.cc
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////
2//
3// Copyright (C) 2000-2025 The Octave Project Developers
4//
5// See the file COPYRIGHT.md in the top-level directory of this
6// distribution or <https://octave.org/copyright/>.
7//
8// This file is part of Octave.
9//
10// Octave is free software: you can redistribute it and/or modify it
11// under the terms of the GNU General Public License as published by
12// the Free Software Foundation, either version 3 of the License, or
13// (at your option) any later version.
14//
15// Octave is distributed in the hope that it will be useful, but
16// WITHOUT ANY WARRANTY; without even the implied warranty of
17// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18// GNU General Public License for more details.
19//
20// You should have received a copy of the GNU General Public License
21// along with Octave; see the file COPYING. If not, see
22// <https://www.gnu.org/licenses/>.
23//
24////////////////////////////////////////////////////////////////////////
25
26/*
27 16. July 2000 - Kai Habel: first release
28
29 25. September 2002 - Changes by Rafael Laboissiere <rafael@laboissiere.net>
30
31 * Added Qbb option to normalize the input and avoid crashes in Octave.
32 * delaunayn accepts now a second (optional) argument that must be a string
33 containing extra options to the qhull command.
34 * Fixed doc string. The dimension of the result matrix is [m, dim+1], and
35 not [n, dim-1].
36
37 6. June 2006: Changes by Alexander Barth <abarth@marine.usf.edu>
38
39 * triangulate non-simplicial facets
40 * allow options to be specified as cell array of strings
41 * change the default options (for compatibility with matlab)
42*/
43
44#if defined (HAVE_CONFIG_H)
45# include "config.h"
46#endif
47
48#include <cstdio>
49
50#include <limits>
51#include <string>
52
53#include "Array.h"
54#include "dMatrix.h"
55#include "dRowVector.h"
56#include "oct-locbuf.h"
57#include "unwind-prot.h"
58
59#include "defun-dld.h"
60#include "error.h"
61#include "errwarn.h"
62#include "ovl.h"
63
64#if defined (HAVE_QHULL)
65# include "oct-qhull.h"
66# if defined (NEED_QHULL_R_VERSION)
67char qh_version[] = "__delaunayn__.oct 2007-08-21";
68# endif
69#endif
70
72
73#if defined (HAVE_QHULL)
74
75static void
76free_qhull_memory (qhT *qh)
77{
78 qh_freeqhull (qh, ! qh_ALL);
79
80 int curlong, totlong;
81 qh_memfreeshort (qh, &curlong, &totlong);
82
83 if (curlong || totlong)
84 warning ("__delaunayn__: did not free %d bytes of long memory (%d pieces)",
85 totlong, curlong);
86}
87
88static bool
89octave_qhull_dims_ok (octave_idx_type dim, octave_idx_type n, const char *who)
90{
91 if (sizeof (octave_idx_type) > sizeof (int))
92 {
93 int maxval = std::numeric_limits<int>::max ();
94
95 if (dim > maxval || n > maxval)
96 error ("%s: dimension too large for Qhull", who);
97 }
98
99 return true;
100}
101
102#endif
103
104DEFUN_DLD (__delaunayn__, args, ,
105 doc: /* -*- texinfo -*-
106@deftypefn {} {@var{T} =} __delaunayn__ (@var{pts})
107@deftypefnx {} {@var{T} =} __delaunayn__ (@var{pts}, @var{options})
108Undocumented internal function.
109@end deftypefn */)
110
111{
112#if defined (HAVE_QHULL)
113
114 int nargin = args.length ();
115
116 if (nargin < 1 || nargin > 2)
117 print_usage ();
118
119 octave_value_list retval;
120
121 retval(0) = 0.0;
122
123 Matrix p (args(0).matrix_value ());
124 const octave_idx_type dim = p.columns ();
125 const octave_idx_type n = p.rows ();
126
127 if (! octave_qhull_dims_ok (dim, n, "__delaynayn__"))
128 return retval;
129
130 // Default options
131 std::string options;
132 if (dim <= 3)
133 options = "Qt Qbb Qc";
134 else
135 options = "Qt Qbb Qc Qx";
136
137 if (nargin == 2)
138 {
139 if (args(1).is_string ())
140 options = args(1).string_value ();
141 else if (args(1).isempty ())
142 ; // Use default options
143 else if (args(1).iscellstr ())
144 {
145 options = "";
146 Array<std::string> tmp = args(1).cellstr_value ();
147
148 for (octave_idx_type i = 0; i < tmp.numel (); i++)
149 options += tmp(i) + ' ';
150 }
151 else
152 error ("__delaunayn__: OPTIONS argument must be a string, cell array of strings, or empty");
153 }
154
155 if (n > dim + 1)
156 {
157 p = p.transpose ();
158 double *pt_array = p.rwdata ();
159 boolT ismalloc = false;
160
161 std::string cmd = "qhull d " + options;
162
163 // Set the outfile pointer to stdout for status information.
164 FILE *outfile = nullptr;
165
166 // Set the errfile pointer to stderr for errors and summary information.
167 // Note: pointer cannot be NULL to disable reporting, unlike outfile.
168#if defined (OCTAVE_HAVE_WINDOWS_FILESYSTEM) && ! defined (OCTAVE_HAVE_POSIX_FILESYSTEM)
169 FILE *errfile = std::fopen ("NUL", "w");
170#else
171 FILE *errfile = std::fopen ("/dev/null", "w");
172#endif
173
174 if (! errfile)
175 error ("__delaunayn__: unable to redirect Qhull errors to /dev/null");
176
177 unwind_action close_errfile ([errfile] () { std::fclose (errfile); });
178
179 qhT context = { };
180 qhT *qh = &context;
181
182 int exitcode = qh_new_qhull (qh, dim, n, pt_array, ismalloc, &cmd[0],
183 outfile, errfile);
184
185 unwind_action free_memory ([qh] () { free_qhull_memory (qh); });
186
187 if (exitcode)
188 error ("__delaunayn__: qhull failed");
189
190 // triangulate non-simplicial facets
191 qh_triangulate (qh);
192
193 facetT *facet;
194 vertexT *vertex, **vertexp;
195 octave_idx_type nf = 0;
196 octave_idx_type i = 0;
197
198 FORALLfacets
199 {
200 if (! facet->upperdelaunay)
201 nf++;
202
203 // Double check. Non-simplicial facets will cause segfault below
204 if (! facet->simplicial)
205 error ("__delaunayn__: Qhull returned non-simplicial facets -- try delaunayn with different options");
206 }
207
208 Matrix simpl (nf, dim+1);
209
210 FORALLfacets
211 {
212 if (! facet->upperdelaunay)
213 {
214 octave_idx_type j = 0;
215
216 FOREACHvertex_ (facet->vertices)
217 {
218 simpl(i, j++) = 1 + qh_pointid(qh, vertex->point);
219 }
220 i++;
221 }
222 }
223
224 retval(0) = simpl;
225 }
226 else if (n == dim + 1)
227 {
228 // FIXME: One should check if nx points span a simplex.
229 // I will look at this later.
230 RowVector vec (n);
231 for (octave_idx_type i = 0; i < n; i++)
232 vec(i) = i + 1.0;
233
234 retval(0) = vec;
235 }
236
237 return retval;
238
239#else
240
241 octave_unused_parameter (args);
242
243 err_disabled_feature ("__delaunayn__", "Qhull");
244
245#endif
246}
247
248/*
249## No test needed for internal helper function.
250%!assert (1)
251*/
252
253OCTAVE_END_NAMESPACE(octave)
N Dimensional Array with copy-on-write semantics.
Definition Array.h:130
octave_idx_type rows() const
Definition Array.h:463
octave_idx_type columns() const
Definition Array.h:475
T * rwdata()
Size of the specified dimension.
octave_idx_type numel() const
Number of elements in the array.
Definition Array.h:418
Matrix transpose() const
Definition dMatrix.h:138
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
#define DEFUN_DLD(name, args_name, nargout_name, doc)
Macro to define an at run time dynamically loadable builtin function.
Definition defun-dld.h:61
void print_usage()
Definition defun-int.h:72
void warning(const char *fmt,...)
Definition error.cc:1078
void error(const char *fmt,...)
Definition error.cc:1003
void err_disabled_feature(const std::string &fcn, const std::string &feature, const std::string &pkg)
Definition errwarn.cc:53