/*
  ---------------------------------

  user_eg.c
  sample code for calling qhull() from an application.  Uses non-reentrant libqhull

  call with:

     user_eg "cube/diamond options" "delaunay options" "halfspace options"

  for example:

     user_eg                             # return summaries

     user_eg "n" "o" "Fp"                # return normals, OFF, points

     user_eg "n Qt" "o" "Fp"             # triangulated cube

     user_eg "QR0 p" "QR0 v p" "QR0 Fp"  # rotate input and return points
                                         # 'v' returns Voronoi
                                         # transform is rotated for halfspaces

   main() makes three runs of qhull.

     1) compute the convex hull of a cube

     2a) compute the Delaunay triangulation of random points

     2b) find the Delaunay triangle closest to a point.

     3) compute the halfspace intersection of a diamond

 notes:

   For another example, see main() in unix.c and user_eg2.c.
   These examples, call qh_qhull() directly.  They allow
   tighter control on the code loaded with Qhull.

   For a C++ example, see user_eg3/user_eg3_r.cpp

   Summaries are sent to stderr if other output formats are used

   compiled by 'make bin/user_eg'

   see libqhull.h for data structures, macros, and user-callable functions.
*/

#include "libqhull/qhull_a.h"

/*-------------------------------------------------
-internal function prototypes
*/
void print_summary(void);
void makecube(coordT *points, int numpoints, int dim);
void makeDelaunay(coordT *points, int numpoints, int dim, int seed);
void findDelaunay(int dim);
void makehalf(coordT *points, int numpoints, int dim);

/*-------------------------------------------------
-print_summary( )
*/
void print_summary(void) {
  facetT *facet;
  int k;

  printf("\n%d vertices and %d facets with normals:\n",
                 qh num_vertices, qh num_facets);
  FORALLfacets {
    for (k=0; k < qh hull_dim; k++)
      printf("%6.2g ", facet->normal[k]);
    printf("\n");
  }
}

/*--------------------------------------------------
-makecube- set points to vertices of cube
  points is numpoints X dim
*/
void makecube(coordT *points, int numpoints, int dim) {
  int j,k;
  coordT *point;

  for (j=0; jlocate a facet with qh_findbestfacet()
  calls qh_setdelaunay() to project the point to a parabaloid
warning:
  Errors if it finds a tricoplanar facet ('Qt').  The corresponding Delaunay triangle
  is in the set of tricoplanar facets or one of their neighbors.  This search
  is not implemented here.
*/
void findDelaunay(int dim) {
  int k;
  coordT point[ 100];
  boolT isoutside;
  realT bestdist;
  facetT *facet;
  vertexT *vertex, **vertexp;

  for (k=0; k < dim; k++)
    point[k]= 0.5;
  qh_setdelaunay(dim+1, 1, point);
  facet= qh_findbestfacet(point, qh_ALL, &bestdist, &isoutside);
  if (facet->tricoplanar) {
    fprintf(stderr, "findDelaunay: search not implemented for triangulated, non-simplicial Delaunay regions (tricoplanar facet, f%d).\n",
       facet->id);
    qh_errexit(qh_ERRqhull, facet, NULL);
  }
  FOREACHvertex_(facet->vertices) {
    for (k=0; k < dim; k++)
      printf("%5.2f ", vertex->point[k]);
    printf("\n");
  }
} /*.findDelaunay.*/

/*--------------------------------------------------
-makehalf- set points to halfspaces for a (dim)-dimensional diamond
  points is numpoints X dim+1

  each halfspace consists of dim coefficients followed by an offset
*/
void makehalf(coordT *points, int numpoints, int dim) {
  int j,k;
  coordT *point;

  for (j=0; j= 2 ? argv[1] : "");
  numpoints= SIZEcube;
  makecube(points, numpoints, DIM);
  for (i=numpoints; i--; )
    rows[i]= points+dim*i;
  qh_printmatrix(outfile, "input", rows, numpoints, dim);
  fflush(NULL);
  exitcode= qh_new_qhull(dim, numpoints, points, ismalloc,
                      flags, outfile, errfile);
  fflush(NULL);
  if (!exitcode) {                  /* if no error */
    /* 'qh facet_list' contains the convex hull */
    print_summary();
    FORALLfacets {
       /* ... your code ... */
    }
  }
#ifdef qh_NOmem
  qh_freeqhull(qh_ALL);
#else
  qh_freeqhull(!qh_ALL);                   /* free long memory  */
  qh_memfreeshort(&curlong, &totlong);    /* free short memory and memory allocator */
  if (curlong || totlong)
    fprintf(errfile, "qhull internal warning (user_eg, #1): did not free %d bytes of long memory (%d pieces)\n", totlong, curlong);
#endif

  /*
    Run 2: Delaunay triangulation, reusing the previous qh/qh_qh
  */

  printf( "\n========\ncompute %d-d Delaunay triangulation\n", dim);
  sprintf(flags, "qhull s d Tcv %s", argc >= 3 ? argv[2] : "");
  numpoints= SIZEcube;
  makeDelaunay(points, numpoints, dim, (int)time(NULL));
  for (i=numpoints; i--; )
    rows[i]= points+dim*i;
  qh_printmatrix(outfile, "input", rows, numpoints, dim);
  fflush(NULL);
  exitcode= qh_new_qhull(dim, numpoints, points, ismalloc,
                      flags, outfile, errfile);
  fflush(NULL);
  if (!exitcode) {                  /* if no error */
    /* 'qh facet_list' contains the convex hull */
    /* If you want a Voronoi diagram ('v') and do not request output (i.e., outfile=NULL),
       call qh_setvoronoi_all() after qh_new_qhull(). */
    print_summary();
    FORALLfacets {
       /* ... your code ... */
    }
    printf( "\nfind %d-d Delaunay triangle or adjacent triangle closest to [0.5, 0.5, ...]\n", dim);
    exitcode= setjmp(qh errexit);
    if (!exitcode) {
      /* Trap Qhull errors from findDelaunay().  Without the setjmp(), Qhull
         will exit() after reporting an error */
      qh NOerrexit= False;
      findDelaunay(DIM);
    }
    qh NOerrexit= True;
  }
#if qh_QHpointer  /* see user.h */
  {
    qhT *oldqhA, *oldqhB;
    coordT pointsB[DIM*TOTpoints]; /* array of coordinates for each point */

    printf( "\n========\nCompute a new triangulation as a separate instance of Qhull\n");
    oldqhA= qh_save_qhull();
    sprintf(flags, "qhull s d Tcv %s", argc >= 3 ? argv[2] : "");
    numpoints= SIZEcube;
    makeDelaunay(pointsB, numpoints, dim, (int)time(NULL)+1);
    for (i=numpoints; i--; )
      rows[i]= pointsB+dim*i;
    qh_printmatrix(outfile, "input", rows, numpoints, dim);
    fflush(NULL);
    exitcode= qh_new_qhull(dim, numpoints, pointsB, ismalloc,
                      flags, outfile, errfile);
    fflush(NULL);
    if (!exitcode)
      print_summary();
    printf( "\n========\nFree memory allocated by the new instance of Qhull, and redisplay the old results.\n");
    oldqhB= qh_save_qhull();
    qh_restore_qhull(&oldqhA);
    print_summary();
    printf( "\nfree first triangulation and restore second one.\n");
    qh_freeqhull(qh_ALL);               /* free short and long memory used by first call */
                                         /* do not use qh_memfreeshort */
    qh_restore_qhull(&oldqhB);
    printf( "\n\n");
    print_summary();
  }
#endif

#ifdef qh_NOmem
  qh_freeqhull(qh_ALL);
#else
  qh_freeqhull(!qh_ALL);                 /* free long memory */
  qh_memfreeshort(&curlong, &totlong);  /* free short memory and memory allocator */
  if (curlong || totlong)
    fprintf(errfile, "qhull internal warning (user_eg, #2): did not free %d bytes of long memory (%d pieces)\n", totlong, curlong);
#endif

  /*
    Run 3: halfspace intersection about the origin
  */
  printf( "\n========\ncompute halfspace intersection about the origin for a diamond\n");
  sprintf(flags, "qhull H0 s Tcv %s", argc >= 4 ? argv[3] : "Fp");
  numpoints= SIZEcube;
  makehalf(points, numpoints, dim);
  for (i=numpoints; i--; )
    rows[i]= points+(dim+1)*i;
  qh_printmatrix(outfile, "input as halfspace coefficients + offsets", rows, numpoints, dim+1);
  fflush(NULL);
  /* use qh_sethalfspace_all to transform the halfspaces yourself.
     If so, set 'qh feasible_point and do not use option 'Hn,...' [it would retransform the halfspaces]
  */
  exitcode= qh_new_qhull(dim+1, numpoints, points, ismalloc,
                      flags, outfile, errfile);
  fflush(NULL);
  if (!exitcode)
    print_summary();
#ifdef qh_NOmem
  qh_freeqhull(qh_ALL);
#else
  qh_freeqhull(!qh_ALL);
  qh_memfreeshort(&curlong, &totlong);
  if (curlong || totlong)  /* could also check previous runs */
    fprintf(stderr, "qhull internal warning (user_eg, #3): did not free %d bytes of long memory (%d pieces)\n",
       totlong, curlong);
#endif
  return exitcode;
} /* main */