Overview

Version 2.

The package DELAUNAYSPARSE (ACM TOMS Algorithm 1012) contains serial and parallel codes, written in FORTRAN 2003 with OpenMP 4.5, for performing interpolation in medium to high dimensions via a sparse subset of the Delaunay triangulation. In addition to the original FORTRAN source code, this site can be used to download a Python 3.6+ wrapper and C/C++ bindings for DELAUNAYSPARSE. Note that each of the three downloads is self-contained, with the Python wrapper and C bindings each containing a subset of the FORTRAN as is needed to build and run their respective functionalities. Command line drivers, which accept formatted data files, are also available by downloading the original FORTRAN source code.

The serial driver subroutine is DELAUNAYSPARSES and the parallel driver subroutine is DELAUNAYSPARSEP. The subroutines DELAUNAYSPARSE{S|P} use the module REAL_PRECISION from HOMPACK90 (ACM TOMS Algorithm 777) for specifying the real data type, and a minimal copy of R. Fletcher's quadratic program solver BQPD with its critical dependencies. The master module DELSPARSE_MOD includes the REAL_PRECISION module and interface blocks for both DELAUNAYSPARSES and DELAUNAYSPARSEP, as well as an interface block for the new PROJECT subroutine, which wraps BQPD for convex hull projections and may be of separate interest. Comments at the beginning of the driver subroutines DELAUNAYSPARSE{S|P} document the arguments and usage, and examples demonstrating their usage are provided in the sample programs samples.f90 and samplep.f90.

Note that in version 1, the SLATEC subroutine DWNNLS (ACM TOMS Algorithm 587) for solving inequality constrained least squares problems was an additional requirement, but it has been replaced by BQPD starting in version 2.

Please note the BQPD carries a stricter LICENSE than DELAUNAYSPARSE

If needed, DELAUNAYSPARSE can be used without BQPD by setting the optional argument: EXTRAP=0.0_R8

However, this will not allow for extrapolation outside the convex hull. To enable extrapolation or to use the new PROJECT subroutine, please note the secondary LICENSE file in bqpd_min/LICENSE.

The physical organization of the main TOMS source download into files is as follows. Further details on using the Python wrapper and C bindings downloads are given in their respective README files, which are included in the downloads.

To check that the installation of DELAUNAYSPARSES and DELAUNAYSPARSEP is correct, assuming that your Fortran compiler allows mixing fixed format .f and free format .f90 files in the same compile command, use the command

$FORT $OPTS delsparse.f90 slatec.f lapack.f blas.f test_install.f90 \
  -o test_install $LIBS

where $FORT is a FORTRAN 2003 compliant compiler supporting OpenMP 4.5, $OPTS is a list of compiler options, and $LIBS is a list of flags to link the BLAS and LAPACK libraries, if those exist on your system (in which case the files blas.f and lapack.f can be omitted from the compile command). To run the parallel code, $OPTS must include the compiler option for OpenMP.

Then run the tests using

./test_install

To compile and link the sample main programs sample{s|p}.f90, use

$FORT $OPTS delsparse.f90 slatec.f lapack.f blas.f sample{s|p}.f90 \
  -o sample{s|p} $LIBS

similar to above. To run a sample main program, use

./sample{s|p} sample_input{2|4}d.dat

where sample_input{2|4}d.dat could be replaced by any other similarly formatted data file.

Usage Information for DELAUNAYSPARSE

DELAUNAYSPARSE solves the multivariate interpolation problem:

Given a set of N points PTS and a set of M interpolation points Q in R^D, for each interpolation point Q_i in Q, identify the set of D+1 data points in PTS that are the vertices of a Delaunay simplex containing Q_i.

These vertices can be used to calculate the Delaunay interpolant using a piecewise linear model.

For more information on the underlying algorithm, see

    @inproceedings{algorithm,
        author={Chang, Tyler H. and Watson, Layne T. and Lux, Thomas C. H. and
                Li, Bo and Xu, Li and Butt, Ali R. and Cameron, Kirk W. and
                Hong, Yili},
        title={A polynomial time algorithm for multivariate interpolation in
               arbitrary dimension via the {D}elaunay triangulation},
        year={2018},
        month={mar},
        booktitle={Proc. ACMSE 2018 Conference (ACMSE 18)},
        publisher={ACM},
        location={Richmond, KY},
        doi={10.1145/3190645.3190680}
    }

For more information on the DELAUNAYSPARSE software (version 1), see

    @article{toms1012,
        author={Chang, Tyler H. and Watson, Layne T. and Lux, Thomas C. H.
                and Butt, Ali R. and Cameron, Kirk W. and Hong, Yili},
        title={Algorithm 1012: {DELAUNAYSPARSE}: Interpolation via a sparse
               subset of the {D}elaunay triangulation in medium to high
               dimensions},
        year={2020},
        journal={ACM Trans. Math. Softw.},
        volume={46},
        number={4},
        articleno={38},
        nopages={20},
        doi={10.1145/3422818}
    }

For more information on the changes in version 2, see

    @article{toms1012remark,
        author={Chang, Tyler H. and Watson, Layne T. and Leyffer, Sven
                and Lux, Thomas C. H. and Almohri, Hussain M. J.},
        title={Remark on {Algorithm 1012}: Computing projections with large
        data sets},
        year={2024},
        journal={ACM Trans. Math. Softw.},
        nopages={8},
        note={to appear}
    }

DELAUNAYSPARSE contains a Fortran module

as well as C bindings two command-line drivers and a Python 3 wrapper

These interfaces are described in the following sections.

Fortran interface

DELAUNAYSPARSE is written in Fortran 2003, and this is its native interface. The Fortran interface contains two drivers:

DELAUNAYSPARSES

The interface for DELAUNAYSPARSES is

SUBROUTINE DELAUNAYSPARSES( D, N, PTS, M, Q, SIMPS, WEIGHTS, IERR,     &
                            INTERP_IN, INTERP_OUT, EPS, EXTRAP, RNORM, &
                            IBUDGET, CHAIN, EXACT                      )

Each of the above parameters is described below.

On input:

On output:

Optional arguments:

Subroutines and functions directly referenced from BLAS are

and from LAPACK are

The quadratic programming solver BQPD is also used. For more information, see Annals of Operations Research, 46 : 307--334 (1993). The module REAL_PRECISION from HOMPACK90 (ACM TOMS Algorithm 777) is used for the real data type. The REAL_PRECISION module, DELAUNAYSPARSES, and BQPD and its dependencies comply with the Fortran 2008 standard.

DELAUNAYSPARSEP

The interface for DELAUNAYSPARSEP is

    SUBROUTINE DELAUNAYSPARSEP( D, N, PTS, M, Q, SIMPS, WEIGHTS, IERR,     &
                                INTERP_IN, INTERP_OUT, EPS, EXTRAP, RNORM, &
                                IBUDGET, CHAIN, EXACT, PMODE               )

Each of the above parameters is described below.

On input:

On output:

Optional arguments:

Subroutines and functions directly referenced from BLAS are

and from LAPACK are

The quadratic programming solver BQPD is also used. For more information, see Annals of Operations Research, 46 : 307--334 (1993). The module REAL_PRECISION from HOMPACK90 (ACM TOMS Algorithm 777) is used for the real data type. The REAL_PRECISION module, DELAUNAYSPARSEP, and BQPD and its dependencies comply with the Fortran 2008 standard.

PROJECT

The interface for PROJECT is

SUBROUTINE PROJECT(D, N, PTS, Z, RNORM, IERR, EPS, WEIGHTS)

Project a point outside the convex hull of the point set onto the convex hull by solving a non negatively constrained least squares problem with 1 equality constraint (an instance of the WNNLS problem):

   min_X   || MATMUL(PTS, X) - Z ||   s.t.   X >= 0, SUM(X) == 1
The solution to the WNNLS problem stated above gives the projection Z_hat as a convex combination of the data points:
  Z_hat = MATMUL(PTS, X).
The above WNNLS problem is solved via R. Fletcher's QP solver BQPD. Compared to other existing (D)WNNLS solvers, BQPD's flexible nature allows us to exploit the sparsity in the solution X, which should contain at most D positive entries (inactive constraints).

On input:

On output:

Optional arguments:

Subroutines and functions directly referenced from BLAS are

BQPD, its utility functions, and its sparse linear algebra library are also referenced.


Notes: DELAUNAYSPARSE is available free of charge via a permissive MIT LICENSE. However, BQPD uses a slightly more restrictive license located in the bqpd_min subdirectory.