1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273
|
/*****************************************************************************/
/* */
/* (tricall.c) */
/* */
/* Example program that demonstrates how to call Triangle. */
/* */
/* Accompanies Triangle Version 1.6 */
/* July 19, 1996 */
/* */
/* This file is placed in the public domain (but the file that it calls */
/* is still copyrighted!) by */
/* Jonathan Richard Shewchuk */
/* 2360 Woolsey #H */
/* Berkeley, California 94705-1927 */
/* jrs@cs.berkeley.edu */
/* */
/*****************************************************************************/
/* If SINGLE is defined when triangle.o is compiled, it should also be */
/* defined here. If not, it should not be defined here. */
/* #define SINGLE */
#ifdef SINGLE
#define REAL float
#else /* not SINGLE */
#define REAL double
#endif /* not SINGLE */
#include <stdio.h>
#include <stdlib.h>
#include "triangle.h"
/*****************************************************************************/
/* */
/* report() Print the input or output. */
/* */
/*****************************************************************************/
void report(io, markers, reporttriangles, reportneighbors, reportsegments,
reportedges, reportnorms)
struct triangulateio *io;
int markers;
int reporttriangles;
int reportneighbors;
int reportsegments;
int reportedges;
int reportnorms;
{
int i, j;
for (i = 0; i < io->numberofpoints; i++) {
printf("Point %4d:", i);
for (j = 0; j < 2; j++) {
printf(" %.6g", io->pointlist[i * 2 + j]);
}
if (io->numberofpointattributes > 0) {
printf(" attributes");
}
for (j = 0; j < io->numberofpointattributes; j++) {
printf(" %.6g",
io->pointattributelist[i * io->numberofpointattributes + j]);
}
if (markers) {
printf(" marker %d\n", io->pointmarkerlist[i]);
} else {
printf("\n");
}
}
printf("\n");
if (reporttriangles || reportneighbors) {
for (i = 0; i < io->numberoftriangles; i++) {
if (reporttriangles) {
printf("Triangle %4d points:", i);
for (j = 0; j < io->numberofcorners; j++) {
printf(" %4d", io->trianglelist[i * io->numberofcorners + j]);
}
if (io->numberoftriangleattributes > 0) {
printf(" attributes");
}
for (j = 0; j < io->numberoftriangleattributes; j++) {
printf(" %.6g", io->triangleattributelist[i *
io->numberoftriangleattributes + j]);
}
printf("\n");
}
if (reportneighbors) {
printf("Triangle %4d neighbors:", i);
for (j = 0; j < 3; j++) {
printf(" %4d", io->neighborlist[i * 3 + j]);
}
printf("\n");
}
}
printf("\n");
}
if (reportsegments) {
for (i = 0; i < io->numberofsegments; i++) {
printf("Segment %4d points:", i);
for (j = 0; j < 2; j++) {
printf(" %4d", io->segmentlist[i * 2 + j]);
}
if (markers) {
printf(" marker %d\n", io->segmentmarkerlist[i]);
} else {
printf("\n");
}
}
printf("\n");
}
if (reportedges) {
for (i = 0; i < io->numberofedges; i++) {
printf("Edge %4d points:", i);
for (j = 0; j < 2; j++) {
printf(" %4d", io->edgelist[i * 2 + j]);
}
if (reportnorms && (io->edgelist[i * 2 + 1] == -1)) {
for (j = 0; j < 2; j++) {
printf(" %.6g", io->normlist[i * 2 + j]);
}
}
if (markers) {
printf(" marker %d\n", io->edgemarkerlist[i]);
} else {
printf("\n");
}
}
printf("\n");
}
}
/*****************************************************************************/
/* */
/* main() Create and refine a mesh. */
/* */
/*****************************************************************************/
int main()
{
struct triangulateio in, mid, out, vorout;
/* Define input points. */
in.numberofpoints = 4;
in.numberofpointattributes = 1;
in.pointlist = (REAL *) malloc(in.numberofpoints * 2 * sizeof(REAL));
in.pointlist[0] = 0.0;
in.pointlist[1] = 0.0;
in.pointlist[2] = 1.0;
in.pointlist[3] = 0.0;
in.pointlist[4] = 1.0;
in.pointlist[5] = 10.0;
in.pointlist[6] = 0.0;
in.pointlist[7] = 10.0;
in.pointattributelist = (REAL *) malloc(in.numberofpoints *
in.numberofpointattributes *
sizeof(REAL));
in.pointattributelist[0] = 0.0;
in.pointattributelist[1] = 1.0;
in.pointattributelist[2] = 11.0;
in.pointattributelist[3] = 10.0;
in.pointmarkerlist = (int *) malloc(in.numberofpoints * sizeof(int));
in.pointmarkerlist[0] = 0;
in.pointmarkerlist[1] = 2;
in.pointmarkerlist[2] = 0;
in.pointmarkerlist[3] = 0;
in.numberofsegments = 0;
in.numberofholes = 0;
in.numberofregions = 1;
in.regionlist = (REAL *) malloc(in.numberofregions * 4 * sizeof(REAL));
in.regionlist[0] = 0.5;
in.regionlist[1] = 5.0;
in.regionlist[2] = 7.0; /* Regional attribute (for whole mesh). */
in.regionlist[3] = 0.1; /* Area constraint that will not be used. */
printf("Input point set:\n\n");
report(&in, 1, 0, 0, 0, 0, 0);
/* Make necessary initializations so that Triangle can return a */
/* triangulation in `mid' and a voronoi diagram in `vorout'. */
mid.pointlist = (REAL *) NULL; /* Not needed if -N switch used. */
/* Not needed if -N switch used or number of point attributes is zero: */
mid.pointattributelist = (REAL *) NULL;
mid.pointmarkerlist = (int *) NULL; /* Not needed if -N or -B switch used. */
mid.trianglelist = (int *) NULL; /* Not needed if -E switch used. */
/* Not needed if -E switch used or number of triangle attributes is zero: */
mid.triangleattributelist = (REAL *) NULL;
mid.neighborlist = (int *) NULL; /* Needed only if -n switch used. */
/* Needed only if segments are output (-p or -c) and -P not used: */
mid.segmentlist = (int *) NULL;
/* Needed only if segments are output (-p or -c) and -P and -B not used: */
mid.segmentmarkerlist = (int *) NULL;
mid.edgelist = (int *) NULL; /* Needed only if -e switch used. */
mid.edgemarkerlist = (int *) NULL; /* Needed if -e used and -B not used. */
vorout.pointlist = (REAL *) NULL; /* Needed only if -v switch used. */
/* Needed only if -v switch used and number of attributes is not zero: */
vorout.pointattributelist = (REAL *) NULL;
vorout.edgelist = (int *) NULL; /* Needed only if -v switch used. */
vorout.normlist = (REAL *) NULL; /* Needed only if -v switch used. */
/* Triangulate the points. Switches are chosen to read and write a */
/* PSLG (p), preserve the convex hull (c), number everything from */
/* zero (z), assign a regional attribute to each element (A), and */
/* produce an edge list (e), a Voronoi diagram (v), and a triangle */
/* neighbor list (n). */
triangulate("pczAevn", &in, &mid, &vorout);
printf("Initial triangulation:\n\n");
report(&mid, 1, 1, 1, 1, 1, 0);
printf("Initial Voronoi diagram:\n\n");
report(&vorout, 0, 0, 0, 0, 1, 1);
/* Attach area constraints to the triangles in preparation for */
/* refining the triangulation. */
/* Needed only if -r and -a switches used: */
mid.trianglearealist = (REAL *) malloc(mid.numberoftriangles * sizeof(REAL));
mid.trianglearealist[0] = 3.0;
mid.trianglearealist[1] = 1.0;
/* Make necessary initializations so that Triangle can return a */
/* triangulation in `out'. */
out.pointlist = (REAL *) NULL; /* Not needed if -N switch used. */
/* Not needed if -N switch used or number of attributes is zero: */
out.pointattributelist = (REAL *) NULL;
out.trianglelist = (int *) NULL; /* Not needed if -E switch used. */
/* Not needed if -E switch used or number of triangle attributes is zero: */
out.triangleattributelist = (REAL *) NULL;
/* Refine the triangulation according to the attached */
/* triangle area constraints. */
triangulate("prazBP", &mid, &out, (struct triangulateio *) NULL);
printf("Refined triangulation:\n\n");
report(&out, 0, 1, 0, 0, 0, 0);
/* Free all allocated arrays, including those allocated by Triangle. */
free(in.pointlist);
free(in.pointattributelist);
free(in.pointmarkerlist);
free(in.regionlist);
free(mid.pointlist);
free(mid.pointattributelist);
free(mid.pointmarkerlist);
free(mid.trianglelist);
free(mid.triangleattributelist);
free(mid.trianglearealist);
free(mid.neighborlist);
free(mid.segmentlist);
free(mid.segmentmarkerlist);
free(mid.edgelist);
free(mid.edgemarkerlist);
free(vorout.pointlist);
free(vorout.pointattributelist);
free(vorout.edgelist);
free(vorout.normlist);
free(out.pointlist);
free(out.pointattributelist);
free(out.trianglelist);
free(out.triangleattributelist);
return 0;
}
|