MALOC 0.1
|
00001 00030 #ifndef _VPRED_H_ 00031 #define _VPRED_H_ 00032 00033 #include <maloc/maloc_base.h> 00034 00035 /* random() prototype seems to be missing in <stdlib.h> */ 00036 /* 00037 * if !defined(VOSF1) 00038 * extern long int random(void); 00039 * endif 00040 */ 00041 00042 /* On some machines, the exact arithmetic routines might be defeated by the */ 00043 /* use of internal extended precision floating-point registers. Sometimes */ 00044 /* this problem can be fixed by defining certain values to be volatile, */ 00045 /* thus forcing them to be stored to memory and rounded off. This isn't */ 00046 /* a great solution, though, as it slows the arithmetic down. */ 00047 /* */ 00048 /* To try this out, write "#define INEXACT volatile" below. Normally, */ 00049 /* however, INEXACT should be defined to be nothing. ("#define INEXACT".) */ 00050 00052 #define INEXACT /* Nothing */ 00053 /* #define INEXACT volatile */ 00054 00056 #define REAL double 00057 00059 #define REALPRINT doubleprint 00060 00063 #define REALRAND doublerand 00064 00067 #define NARROWRAND narrowdoublerand 00068 00070 #define UNIFORMRAND uniformdoublerand 00071 00084 void Vpred_exactinit(void); 00085 00097 REAL Vpred_orient2d(REAL *pa, REAL *pb, REAL *pc); 00098 00110 REAL Vpred_orient2dfast(REAL *pa, REAL *pb, REAL *pc); 00111 00123 REAL Vpred_orient2dexact(REAL *pa, REAL *pb, REAL *pc); 00124 00139 REAL Vpred_orient3d(REAL *pa, REAL *pb, REAL *pc, REAL *pd); 00140 00155 REAL Vpred_orient3dfast(REAL *pa, REAL *pb, REAL *pc, REAL *pd); 00156 00171 REAL Vpred_orient3dexact(REAL *pa, REAL *pb, REAL *pc, REAL *pd); 00172 00185 REAL Vpred_incircle(REAL *pa, REAL *pb, REAL *pc, REAL *pd); 00186 00199 REAL Vpred_incirclefast(REAL *pa, REAL *pb, REAL *pc, REAL *pd); 00200 00213 REAL Vpred_incircleexact(REAL *pa, REAL *pb, REAL *pc, REAL *pd); 00214 00228 REAL Vpred_insphere(REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL *pe); 00229 00243 REAL Vpred_inspherefast(REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL *pe); 00244 00258 REAL Vpred_insphereexact(REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL *pe); 00259 00260 #endif /* _VPRED_H_ */ 00261 00262