commit 55a940ef222b54638d0655213ad556a139f380b7
parent b3ead38f781671a746209ece4088070d41709a84
Author: Christophe Coustet <christophe.coustet@meso-star.com>
Date: Tue, 7 Feb 2023 16:41:36 +0100
Change internal representation to int64
We experienced difficulties due to internal representation in double and
conversions to int64 for offsetting and clipping using different
scaling/unscaling recipes. The simplest fix was to scale once at
creation time and to work on int64 data. Unscaling is done when data
access is requested.
Diffstat:
7 files changed, 213 insertions(+), 152 deletions(-)
diff --git a/src/scpr.h b/src/scpr.h
@@ -84,10 +84,14 @@ scpr_polygon_ref_put
/* To ensure constant precision, vertice coordinates are truncated, and have a
* limited range. Range checking along with the associated truncation process
* occur at vertices setup.
- * The range of the coordinates is [-INT64_MAX*0.25E-6 +INT64_MAX*0.25E-6], that
- * is approximately [-2.3E18 +2.3E18], and vertex coordinates are truncated past
- * the 6th decimal place. It is an error to use out-of-range values.
- * Truncated coordinates can be retrieved using the appropriate getters. */
+ * The allowed range is [-1e12 +1e12], and vertex coordinates are truncated
+ * approximately past the 6th decimal place. It is an error to use out-of-range
+ * values.
+ * Also, the number of components and vertices can be changed due to a
+ * simplification process and one should not take for granted that the number of
+ * components and vertices stay as provided at construction time.
+ * The actual counts and coordinates should be retrieved using the appropriate
+ * getters. */
SCPR_API res_T
scpr_polygon_setup_indexed_vertices
(struct scpr_polygon* polygon,
@@ -143,10 +147,7 @@ scpr_mesh_ref_put
(struct scpr_mesh* mesh);
/* In a way similar to polygons, meshes have their vertices coordinates range
- * limited and truncated.
- * However, there is currently a difference between the 2 setup operations that
- * prevent to share the same coordinate space.
- * Fixing this will require a change in the Clipper2 library. */
+ * limited to [-1e12 1e12] and truncated past the 6th decimal place. */
SCPR_API res_T
scpr_mesh_setup_indexed_vertices
(struct scpr_mesh* mesh,
diff --git a/src/scpr_c.h b/src/scpr_c.h
@@ -21,7 +21,7 @@
#include <rsys/rsys.h>
#include <rsys/ref_count.h>
#include <rsys/double2.h>
-#include <rsys/dynamic_array_double.h>
+#include <rsys/dynamic_array.h>
#include <rsys/dynamic_array_size_t.h>
#include <rsys/hash_table.h>
@@ -30,6 +30,10 @@
#undef PI
#include <clipper2/clipper.h>
+#define DARRAY_NAME int64
+#define DARRAY_DATA int64_t
+#include <rsys/dynamic_array.h>
+
#define ERR(Expr) { if((res = (Expr)) != RES_OK) goto error; } (void)0
#define TRY(Expr) \
@@ -56,12 +60,16 @@
* precision and divides the coordinate range by 10. */
#define PRECISION 6
+/* Maximum range is somewhat arbitrarily set at [-1e12 + 1e12].
+ * The rationale is that beyond this range accuracy decreases sharply. */
+#define RANGE_MAX 1e12
+
struct mem_allocator;
-struct vertex { double pos[2]; };
+struct vertex { int64_t pos[2]; };
static FINLINE int
vertex_eq(const struct vertex* a, const struct vertex* b)
-{ ASSERT(a && b); return d2_eq(a->pos, b->pos); }
+{ ASSERT(a && b); return (a->pos[0] == b->pos[0]) && (a->pos[1] == b->pos[1]); }
/* Define the vertex to index hash table */
#define HTABLE_NAME vertex
@@ -71,15 +79,15 @@ vertex_eq(const struct vertex* a, const struct vertex* b)
#include <rsys/hash_table.h>
struct scpr_polygon {
- Clipper2Lib::PathsD paths;
- double lower[2], upper[2]; /* Polygon AABB */
+ Clipper2Lib::Paths64 paths;
+ int64_t lower[2], upper[2]; /* Polygon AABB */
ref_T ref;
struct mem_allocator* allocator;
};
struct scpr_mesh {
- struct darray_double coords;
+ struct darray_int64 coords;
struct darray_size_t indices;
ref_T ref;
@@ -87,31 +95,36 @@ struct scpr_mesh {
};
res_T INLINE
-check_and_truncate_vertex
- (double pt[2])
+round_and_set_vertex
+ (const double in[2],
+ int64_t out[2])
{
- double scale = pow(10, PRECISION);
- const int64_t MAX_COORD = INT64_MAX >> 2;
- const int64_t MIN_COORD = -MAX_COORD;
- const double max_coord_d = (double)MAX_COORD;
- const double min_coord_d = (double)MIN_COORD;
- double tmp[2];
- int64_t tmp2[2];
- ASSERT(pt);
- /* Check value is in-range like Clipper2Lib::ScalePath */
- if((pt[0] * scale < min_coord_d) || (pt[0] * scale > max_coord_d)
- || (pt[1] * scale < min_coord_d) || (pt[1] * scale > max_coord_d))
+ /* to optimize scaling / descaling precision
+ * set the scale to a power of double's radix (2) (#25) */
+ double scale = pow(2, (int)log2(pow(10, PRECISION)) + 1);
+ ASSERT(in && out);
+ /* Check value is in-range */
+ if(in[0] < -RANGE_MAX || in[0] > RANGE_MAX
+ || in[1] < -RANGE_MAX || in[1] > RANGE_MAX)
{
return RES_BAD_ARG;
}
- /* Truncate precision to ensure further consistency */
- tmp2[0] = std::llround(pt[0] * scale);
- tmp2[1] = std::llround(pt[1] * scale);
- /* Store truncated vertex */
- tmp[0] = (double)tmp2[0] / scale;
- tmp[1] = (double)tmp2[1] / scale;
- d2_set(pt, tmp);
+ out[0] = std::llround(in[0] * scale);
+ out[1] = std::llround(in[1] * scale);
return RES_OK;
}
+void INLINE
+get_vertex
+ (const int64_t in[2],
+ double out[2])
+{
+ /* to optimize scaling / descaling precision
+ * set the scale to a power of double's radix (2) (#25) */
+ double inv_scale = pow(2, -(int)log2(pow(10, PRECISION)) - 1);
+ ASSERT(in && out);
+ out[0] = (double)in[0] * inv_scale;
+ out[1] = (double)in[1] * inv_scale;
+}
+
#endif
diff --git a/src/scpr_mesh.c b/src/scpr_mesh.c
@@ -39,10 +39,10 @@ scpr_operation_to_clip_type(const enum scpr_operation op)
static INLINE int
aabb_intersect
- (const double lower0[2],
- const double upper0[2],
- const double lower1[2],
- const double upper1[2])
+ (const int64_t lower0[2],
+ const int64_t upper0[2],
+ const int64_t lower1[2],
+ const int64_t upper1[2])
{
ASSERT(lower0 && upper0 && lower1 && upper1);
return lower0[0] < upper1[0]
@@ -52,7 +52,7 @@ aabb_intersect
}
static INLINE int
-aabb_is_degenerated(const double lower[2], const double upper[2])
+aabb_is_degenerated(const int64_t lower[2], const int64_t upper[2])
{
ASSERT(lower && upper);
return lower[0] >= upper[0] || lower[1] >= upper[1];
@@ -60,50 +60,53 @@ aabb_is_degenerated(const double lower[2], const double upper[2])
static void
triangle_compute_aabb
- (double tri[3][2],
- double lower[2],
- double upper[2])
+ (int64_t tri[3][2],
+ int64_t lower[2],
+ int64_t upper[2])
{
size_t ivert;
ASSERT(tri && lower && upper);
- d2_splat(lower, DBL_MAX);
- d2_splat(upper,-DBL_MAX);
+ lower[0] = lower[1] = INT64_MAX;
+ upper[0] = upper[1] = -INT64_MAX;
FOR_EACH(ivert, 0, 3) {
- d2_min(lower, lower, tri[ivert]);
- d2_max(upper, upper, tri[ivert]);
+ int i;
+ for(i = 0; i < 2; i++) {
+ if(lower[i] > tri[ivert][i]) lower[i] = tri[ivert][i];
+ if(upper[i] < tri[ivert][i]) upper[i] = tri[ivert][i];
+ }
}
}
static res_T
register_vertex
- (const double pos[2],
- struct darray_double* coords,
+ (const int64_t pos[2],
+ struct darray_int64* coords,
struct darray_size_t* indices,
struct htable_vertex* vertices)
{
struct vertex v;
size_t* pid, id;
+ unsigned i;
res_T res = RES_OK;
ASSERT(pos && coords && indices && vertices);
- d2_set(v.pos, pos);
+ for(i = 0; i < 2; i++) v.pos[i] = pos[i];
- /* FIXME dirty hack */
- v.pos[0] = (float)v.pos[0];
- v.pos[1] = (float)v.pos[1];
+ v.pos[0] = v.pos[0];
+ v.pos[1] = v.pos[1];
pid = htable_vertex_find(vertices, &v);
if(pid) {
id = *pid;
} else {
- const size_t ivert = darray_double_size_get(coords);
+ const size_t count = darray_int64_size_get(coords);
- ERR(darray_double_resize(coords, ivert+2/*#coords*/));
- d2_set(darray_double_data_get(coords)+ivert, pos);
+ ERR(darray_int64_resize(coords, count+2/*#coords*/));
+ for(i = 0; i < 2; i++) darray_int64_data_get(coords)[count+i] = pos[i];
- id = ivert / 2;
+ id = count / 2;
ERR(htable_vertex_set(vertices, &v, &id));
}
@@ -117,8 +120,8 @@ error:
static FINLINE res_T
register_triangle
- (double tri[3][2],
- struct darray_double* coords, /* Vertex buffer */
+ (int64_t tri[3][2],
+ struct darray_int64* coords, /* Vertex buffer */
struct darray_size_t* indices, /* Index buffer */
struct htable_vertex* vertices) /* Map a vertex to its index */
{
@@ -137,8 +140,8 @@ error:
static res_T
register_paths
- (const Clipper2Lib::PathsD& paths,
- struct darray_double* coords, /* Vertex buffer */
+ (const Clipper2Lib::Paths64& paths,
+ struct darray_int64* coords, /* Vertex buffer */
struct darray_size_t* indices, /* Index buffer */
struct htable_vertex* vertices, /* Map a vertex to its index */
struct polygon* polygon) /* Use to triangulate the clipped polygons */
@@ -151,7 +154,7 @@ register_paths
FOR_EACH(ipath, 0, paths.size()) {
if(paths[ipath].size() == 3) {
FOR_EACH(ivert, 0, 3) {
- double pos[2];
+ int64_t pos[2];
pos[0] = paths[ipath][ivert].x;
pos[1] = paths[ipath][ivert].y;
ERR(register_vertex(pos, coords, indices, vertices));
@@ -165,23 +168,27 @@ register_paths
POLYGON(clear(polygon));
FOR_EACH(ivert, 0, paths[ipath].size()) {
float fpos[3] = {0.f, 0.f, 0.f};
- double pos[2];
+ double posd[2];
+ int64_t pos64[2];
- pos[0] = paths[ipath][ivert].x;
- pos[1] = paths[ipath][ivert].y;
+ pos64[0] = paths[ipath][ivert].x;
+ pos64[1] = paths[ipath][ivert].y;
+ get_vertex(pos64, posd);
- fpos[0] = (float)pos[0], fpos[1] = (float)pos[1];
+ fpos[0] = (float)posd[0], fpos[1] = (float)posd[1];
ERR(polygon_vertex_add(polygon, fpos));
}
ERR(polygon_triangulate(polygon, &ids, &nids));
FOR_EACH(ivert, 0, nids) {
float fpos[3];
+ int64_t pos64[2];
double pos[2];
POLYGON(vertex_get(polygon, ids[ivert], fpos));
- pos[0] = (float)fpos[0];
- pos[1] = (float)fpos[1];
- ERR(register_vertex(pos, coords, indices, vertices));
+ pos[0] = (double)fpos[0];
+ pos[1] = (double)fpos[1];
+ ERR(round_and_set_vertex(pos, pos64));
+ ERR(register_vertex(pos64, coords, indices, vertices));
}
}
}
@@ -192,22 +199,26 @@ error:
}
static void
-mesh_compute_aabb(const struct scpr_mesh* mesh, double lower[2], double upper[2])
+mesh_compute_aabb(const struct scpr_mesh* mesh, int64_t lower[2], int64_t upper[2])
{
- size_t itri, ntris;
+ size_t itri, ntris, i;
SCPR(mesh_get_triangles_count(mesh, &ntris));
- d2_splat(lower, DBL_MAX);
- d2_splat(upper,-DBL_MAX);
+ lower[0] = lower[1] = INT64_MAX;
+ upper[0] = upper[1] = -INT64_MAX;
FOR_EACH(itri, 0, ntris) {
size_t ids[3], ivert;
SCPR(mesh_get_indices(mesh, itri, ids));
FOR_EACH(ivert, 0, 3) {
double pos[2];
+ int64_t pos64[2];
SCPR(mesh_get_position(mesh, ids[ivert], pos));
- d2_min(lower, lower, pos);
- d2_max(upper, upper, pos);
+ CHK(RES_OK == round_and_set_vertex(pos, pos64));
+ for(i = 0; i < 2; i++) {
+ if(lower[i] > pos64[i]) lower[i] = pos64[i];
+ if(upper[i] < pos64[i]) upper[i] = pos64[i];
+ }
}
}
}
@@ -218,7 +229,7 @@ mesh_release(ref_T* ref)
struct scpr_mesh* mesh;
ASSERT(ref);
mesh = CONTAINER_OF(ref, struct scpr_mesh, ref);
- darray_double_release(&mesh->coords);
+ darray_int64_release(&mesh->coords);
darray_size_t_release(&mesh->indices);
MEM_RM(mesh->allocator, mesh);
}
@@ -245,7 +256,7 @@ scpr_mesh_create(struct mem_allocator* mem_allocator, struct scpr_mesh** out_mes
res = RES_MEM_ERR;
goto error;
}
- darray_double_init(mesh->allocator, &mesh->coords);
+ darray_int64_init(mesh->allocator, &mesh->coords);
darray_size_t_init(mesh->allocator, &mesh->indices);
ref_init(&mesh->ref);
mesh->allocator = allocator;
@@ -287,7 +298,7 @@ scpr_mesh_setup_indexed_vertices
void (*get_position)(const size_t ivert, double pos[2], void* ctx),
void* data)
{
- double* pos = NULL;
+ int64_t* pos = NULL;
size_t* ids = NULL;
size_t i;
res_T res = RES_OK;
@@ -297,21 +308,21 @@ scpr_mesh_setup_indexed_vertices
goto error;
}
- ERR(darray_double_resize(&mesh->coords, nverts*2/*#coords per vertex*/));
+ ERR(darray_int64_resize(&mesh->coords, nverts*2/*#coords per vertex*/));
ERR(darray_size_t_resize(&mesh->indices, ntris*3/*#vertices per triangle*/));
/* Fetch mesh positions */
- pos = darray_double_data_get(&mesh->coords);
+ pos = darray_int64_data_get(&mesh->coords);
FOR_EACH(i, 0, nverts) {
- get_position(i, pos+i*2, data);
- /* Check range and truncate precision to ensure further consistency */
- ERR(check_and_truncate_vertex(pos+i*2));
+ double posd[2];
+ get_position(i, posd, data);
+ ERR(round_and_set_vertex(posd, pos + i*2));
}
/* Fetch mesh indices */
ids = darray_size_t_data_get(&mesh->indices);
FOR_EACH(i, 0, ntris) {
- get_indices(i, ids +i*3, data);
+ get_indices(i, ids + i*3, data);
if(ids[i*3+0] >= nverts
|| ids[i*3+1] >= nverts
|| ids[i*3+2] >= nverts) {
@@ -324,7 +335,7 @@ exit:
return res;
error:
if(mesh) {
- darray_double_clear(&mesh->coords);
+ darray_int64_clear(&mesh->coords);
darray_size_t_clear(&mesh->indices);
}
goto exit;
@@ -334,8 +345,8 @@ res_T
scpr_mesh_get_vertices_count(const struct scpr_mesh* mesh, size_t* nverts)
{
if(!mesh || !nverts) return RES_BAD_ARG;
- ASSERT((darray_double_size_get(&mesh->coords) % 2/*#coords per vertex*/)==0);
- *nverts = darray_double_size_get(&mesh->coords) / 2/*#coords per vertex*/;
+ ASSERT((darray_int64_size_get(&mesh->coords) % 2/*#coords per vertex*/)==0);
+ *nverts = darray_int64_size_get(&mesh->coords) / 2/*#coords per vertex*/;
return RES_OK;
}
@@ -368,12 +379,14 @@ scpr_mesh_get_position
(const struct scpr_mesh* mesh, const size_t ivert, double pos[2])
{
size_t nverts;
+ int64_t pos64[2];
const size_t i = ivert * 2/*#coords per vertex*/;
if(!mesh || !pos) return RES_BAD_ARG;
SCPR(mesh_get_vertices_count(mesh, &nverts));
if(ivert >= nverts) return RES_BAD_ARG;
- pos[0] = darray_double_cdata_get(&mesh->coords)[i+0];
- pos[1] = darray_double_cdata_get(&mesh->coords)[i+1];
+ pos64[0] = darray_int64_cdata_get(&mesh->coords)[i+0];
+ pos64[1] = darray_int64_cdata_get(&mesh->coords)[i+1];
+ get_vertex(pos64, pos);
return RES_OK;
}
@@ -383,17 +396,18 @@ scpr_mesh_clip
const enum scpr_operation op,
struct scpr_polygon* poly_desc)
{
- double lower[2], upper[2];
+ int64_t lower[2], upper[2];
struct polygon* polygon = NULL; /* Use to triangulate clipped polygons */
- struct darray_double coords; /* Coordinates of the clipped mesh */
+ struct darray_int64 coords; /* Coordinates of the clipped mesh */
struct darray_size_t indices; /* Indices of the clipped mesh */
struct htable_vertex vertices; /* Map a coordinate to its index */
- Clipper2Lib::ClipperD clipper(PRECISION);
- Clipper2Lib::PathsD output; /* Contour of the clipped polgyon */
- Clipper2Lib::PathsD cand_path; /* Contour of the candidate polygon */
- Clipper2Lib::PathD tmp;
+ Clipper2Lib::Clipper64 clipper;
+ Clipper2Lib::Paths64 output; /* Contour of the clipped polgyon */
+ Clipper2Lib::Paths64 cand_path; /* Contour of the candidate polygon */
+ Clipper2Lib::Path64 tmp;
Clipper2Lib::ClipType clip_type; /* Type of clipping to perform */
size_t itri, ntris, ivert;
+ int i;
res_T res = RES_OK;
@@ -402,7 +416,7 @@ scpr_mesh_clip
clip_type = scpr_operation_to_clip_type(op);
- darray_double_init(mesh->allocator, &coords);
+ darray_int64_init(mesh->allocator, &coords);
darray_size_t_init(mesh->allocator, &indices);
htable_vertex_init(mesh->allocator, &vertices);
@@ -412,8 +426,10 @@ scpr_mesh_clip
if(aabb_is_degenerated(lower, upper)) goto exit;
/* Compute the overall aabb of the candidate and the clip polygon */
- d2_min(lower, lower, poly_desc->lower);
- d2_max(upper, upper, poly_desc->upper);
+ for(i = 0; i < 2; i++) {
+ if(lower[i] > poly_desc->lower[i]) lower[i] = poly_desc->lower[i];
+ if(upper[i] < poly_desc->upper[i]) upper[i] = poly_desc->upper[i];
+ }
/* Create the polygon structure used to triangulate the clipped polygons */
ERR(polygon_create(mesh->allocator, &polygon));
@@ -423,18 +439,22 @@ scpr_mesh_clip
FOR_EACH(itri, 0, ntris) {
size_t ids[3];
double tri[3][2];
+ int64_t tri64[3][2];
/* Fetch the triangle vertices and compute its AABB */
SCPR(mesh_get_indices(mesh, itri, ids));
SCPR(mesh_get_position(mesh, ids[0], tri[0]));
SCPR(mesh_get_position(mesh, ids[1], tri[1]));
SCPR(mesh_get_position(mesh, ids[2], tri[2]));
- triangle_compute_aabb(tri, lower, upper);
+ ERR(round_and_set_vertex(tri[0], tri64[0]));
+ ERR(round_and_set_vertex(tri[1], tri64[1]));
+ ERR(round_and_set_vertex(tri[2], tri64[2]));
+ triangle_compute_aabb(tri64, lower, upper);
/* Do not clip triangles that do not intersect the clip AABB */
if(!aabb_intersect(lower, upper, poly_desc->lower, poly_desc->upper)) {
if(op != SCPR_AND) {
- ERR(register_triangle(tri, &coords, &indices, &vertices));
+ ERR(register_triangle(tri64, &coords, &indices, &vertices));
}
continue;
}
@@ -443,9 +463,9 @@ scpr_mesh_clip
cand_path.clear();
tmp.clear();
FOR_EACH(ivert, 0, 3) {
- Clipper2Lib::PointD pt;
- pt.x = tri[ivert][0];
- pt.y = tri[ivert][1];
+ Clipper2Lib::Point64 pt;
+ pt.x = tri64[ivert][0];
+ pt.y = tri64[ivert][1];
tmp.push_back(pt);
}
cand_path.push_back(tmp);
@@ -460,12 +480,12 @@ scpr_mesh_clip
ERR(register_paths (output, &coords, &indices, &vertices, polygon));
}
- ERR(darray_double_copy_and_clear(&mesh->coords, &coords));
+ ERR(darray_int64_copy_and_clear(&mesh->coords, &coords));
ERR(darray_size_t_copy_and_clear(&mesh->indices, &indices));
exit:
if(polygon) POLYGON(ref_put(polygon));
- darray_double_release(&coords);
+ darray_int64_release(&coords);
darray_size_t_release(&indices);
htable_vertex_release(&vertices);
return res;
diff --git a/src/scpr_polygon.c b/src/scpr_polygon.c
@@ -50,14 +50,14 @@ polygon_release(ref_T* ref)
ASSERT(ref);
polygon = CONTAINER_OF(ref, struct scpr_polygon, ref);
/* Call destructor for paths */
- polygon->paths.Clipper2Lib::PathsD::~PathsD();
+ polygon->paths.Clipper2Lib::Paths64::~Paths64();
MEM_RM(polygon->allocator, polygon);
}
static int
path_is_eq
- (const Clipper2Lib::PathD* p1,
- const Clipper2Lib::PathD* p2)
+ (const Clipper2Lib::Path64* p1,
+ const Clipper2Lib::Path64* p2)
{
size_t i, first_vtx, sz;
int opposite_cw;
@@ -88,8 +88,8 @@ path_is_eq
static int
one_path_is_eq
- (const Clipper2Lib::PathsD* pp,
- const Clipper2Lib::PathD* p,
+ (const Clipper2Lib::Paths64* pp,
+ const Clipper2Lib::Path64* p,
char* matched)
{
size_t i, sz;
@@ -133,8 +133,8 @@ scpr_polygon_create
polygon->allocator = allocator;
/* Allocate paths the C++ way (placement new) */
new (&polygon->paths) Clipper2Lib::PathsD;
- d2_splat(polygon->lower, DBL_MAX);
- d2_splat(polygon->upper,-DBL_MAX);
+ polygon->lower[0] = polygon->lower[1] = INT64_MAX;
+ polygon->upper[0] = polygon->upper[1] = -INT64_MAX;
exit:
if(out_polygon) *out_polygon = polygon;
@@ -176,8 +176,10 @@ scpr_polygon_setup_indexed_vertices
void* data)
{
size_t c;
- double scale = pow(10, PRECISION);
- Clipper2Lib::PathsD paths;
+ /* to optimize scaling / descaling precision
+ * set the scale to a power of double's radix (2) (#25) */
+ double inv_scale = pow(2, -(int)log2(pow(10, PRECISION)) - 1);
+ Clipper2Lib::Paths64 paths;
res_T res = RES_OK;
if(!polygon || !get_nverts || !get_position || !data) {
@@ -197,32 +199,36 @@ scpr_polygon_setup_indexed_vertices
/* Fetch polygon positions for connex component c */
FOR_EACH(i, 0, nverts) {
double tmp[2];
- Clipper2Lib::PointD pt;
+ int64_t tmp64[2];
+ Clipper2Lib::Point64 pt;
get_position(c, i, tmp, data);
/* Check range and truncate precision to ensure further consistency */
- ERR(check_and_truncate_vertex(tmp));
- pt.x = tmp[0];
- pt.y = tmp[1];
+ ERR(round_and_set_vertex(tmp, tmp64));
+ pt.x = tmp64[0];
+ pt.y = tmp64[1];
paths[c][i] = pt;
}
}
/* Merge vertices, ... */
- TRY(paths = Clipper2Lib::SimplifyPaths(paths, 1/scale));
+ TRY(paths = Clipper2Lib::SimplifyPaths(paths, inv_scale));
polygon->paths = std::move(paths);
/* Build bounding box */
- d2_splat(polygon->lower, DBL_MAX);
- d2_splat(polygon->upper,-DBL_MAX);
+ polygon->lower[0] = polygon->lower[1] = INT64_MAX;
+ polygon->upper[0] = polygon->upper[1] = -INT64_MAX;
FOR_EACH(c, 0, ncomponents) {
size_t i, nverts;
get_nverts(c, &nverts, data);
FOR_EACH(i, 0, nverts) {
- double tmp[2];
+ int64_t tmp[2];
+ int j;
tmp[0] = polygon->paths[c][i].x;
tmp[1] = polygon->paths[c][i].y;
- d2_min(polygon->lower, polygon->lower, tmp);
- d2_max(polygon->upper, polygon->upper, tmp);
+ for(j = 0; j < 2; j++) {
+ if(polygon->lower[j] > tmp[j]) polygon->lower[j] = tmp[j];
+ if(polygon->upper[j] < tmp[j]) polygon->upper[j] = tmp[j];
+ }
}
}
@@ -231,8 +237,8 @@ exit:
error:
if(polygon) {
polygon->paths.clear();
- d2_splat(polygon->lower, DBL_MAX);
- d2_splat(polygon->upper,-DBL_MAX);
+ polygon->lower[0] = polygon->lower[1] = INT64_MAX;
+ polygon->upper[0] = polygon->upper[1] = -INT64_MAX;
}
goto exit;
}
@@ -255,8 +261,13 @@ scpr_polygon_create_copy
copy_created = 1;
copy->paths = src_polygon->paths;
- d2_set(copy->lower, src_polygon->lower);
- d2_set(copy->upper, src_polygon->upper);
+ int i;
+ for(i = 0; i < 2; i++) {
+ if(copy->lower[i] > src_polygon->lower[i])
+ copy->lower[i] = src_polygon->lower[i];
+ if(copy->upper[i] < src_polygon->upper[i])
+ copy->upper[i] = src_polygon->upper[i];
+ }
exit:
if(out_polygon) *out_polygon = copy;
@@ -301,15 +312,17 @@ scpr_polygon_get_position
{
size_t nverts;
const size_t i = ivert;
- const Clipper2Lib::PointD* pt;
+ const Clipper2Lib::Point64* pt;
+ int64_t pos64[2];
if(!polygon || !pos || icomponent >= polygon->paths.size()) {
return RES_BAD_ARG;
}
SCPR(polygon_get_vertices_count(polygon, icomponent, &nverts));
if(ivert >= nverts) return RES_BAD_ARG;
pt = &polygon->paths[icomponent][i];
- pos[0] = pt->x;
- pos[1] = pt->y;
+ pos64[0] = pt->x;
+ pos64[1] = pt->y;
+ get_vertex(pos64, pos);
return RES_OK;
}
@@ -320,9 +333,12 @@ scpr_offset_polygon
const enum scpr_join_type join_type)
{
size_t c;
- Clipper2Lib::PathD tmp;
+ Clipper2Lib::Path64 tmp;
Clipper2Lib::JoinType cjt;
- Clipper2Lib::PathsD polygon;
+ Clipper2Lib::Paths64 polygon;
+ /* to optimize scaling / descaling precision
+ * set the scale to a power of double's radix (2) (#25) */
+ double scale = pow(2, (int)log2(pow(10, PRECISION)) + 1);
res_T res = RES_OK;
if(!poly_desc) {
@@ -332,11 +348,15 @@ scpr_offset_polygon
/* Check offset polygon will be in-range */
if(offset > 0) {
- double pt[2];
- d2_sub(pt, poly_desc->upper, d2(pt, -offset, -offset));
- ERR(check_and_truncate_vertex(pt));
- d2_add(pt, poly_desc->upper, d2(pt, offset, offset));
- ERR(check_and_truncate_vertex(pt));
+ int i;
+ for(i = 0; i < 2; i++) {
+ if((double)poly_desc->lower[i] - offset < -RANGE_MAX
+ || (double)poly_desc->upper[i] + offset > RANGE_MAX)
+ {
+ res = RES_BAD_ARG;
+ goto error;
+ }
+ }
}
/* Check join type */
@@ -350,23 +370,26 @@ scpr_offset_polygon
goto error;
}
- /* Some known problems when offset=0, not leaving polygon unchanged */
cjt = scpr_join_type_to_clipper_join_type(join_type);
- TRY(poly_desc->paths = Clipper2Lib::InflatePaths(poly_desc->paths, offset, cjt,
- Clipper2Lib::EndType::Polygon, 2, PRECISION));
+ TRY(poly_desc->paths = Clipper2Lib::InflatePaths(poly_desc->paths,
+ offset * scale, cjt, Clipper2Lib::EndType::Polygon));
/* Rebuild AABB */
- d2_splat(poly_desc->lower, DBL_MAX);
- d2_splat(poly_desc->upper,-DBL_MAX);
+ poly_desc->lower[0] = poly_desc->lower[1] = INT64_MAX;
+ poly_desc->upper[0] = poly_desc->upper[1] = -INT64_MAX;
FOR_EACH(c, 0, poly_desc->paths.size()) {
size_t i, nverts;
nverts = poly_desc->paths[c].size();
FOR_EACH(i, 0, nverts) {
- double pos[2];
- ERR(scpr_polygon_get_position(poly_desc, c, i, pos));
- d2_min(poly_desc->lower, poly_desc->lower, pos);
- d2_max(poly_desc->upper, poly_desc->upper, pos);
+ int64_t pos64[2];
+ int j;
+ pos64[0] = poly_desc->paths[c][i].x;
+ pos64[1] = poly_desc->paths[c][i].y;
+ for(j = 0; j < 2; j++) {
+ if(poly_desc->lower[j] > pos64[j]) poly_desc->lower[j] = pos64[j];
+ if(poly_desc->upper[j] < pos64[j]) poly_desc->upper[j] = pos64[j];
+ }
}
}
@@ -396,8 +419,10 @@ scpr_polygon_eq
*is_eq = 0;
goto exit;
}
- if(!d2_eq(polygon1->lower, polygon2->lower)
- || !d2_eq(polygon1->upper, polygon2->upper))
+ if(polygon1->lower[0] != polygon2->lower[0]
+ || polygon1->lower[1] != polygon2->lower[1]
+ || polygon1->upper[0] != polygon2->upper[0]
+ || polygon1->upper[1] != polygon2->upper[1])
{
*is_eq = 0;
goto exit;
diff --git a/src/test_scpr_offset.c b/src/test_scpr_offset.c
@@ -38,13 +38,13 @@ test_single(void)
const double coords2[] = {
0.12345678901234, 0.0,
0.0, 1.0,
- 1.0, 2305843009213,
+ 1.0, 1000000000000,
1.0, 0.0
};
const double coords2_reverse[] = {
0.12345678901234, 0.0,
1.0, 0.0,
- 1.0, 2305843009213,
+ 1.0, 1000000000000,
0.0, 1.0
};
double** coords;
diff --git a/src/test_scpr_polygon.c b/src/test_scpr_polygon.c
@@ -43,7 +43,7 @@ main(int argc, char** argv)
0.5, 0.0,
0.0, 1.0,
1.0, 1.0,
- 1.0, 2305843009214,
+ 1.0, 1000000000001,
1.0, 0.0
};
double** coords;
diff --git a/src/test_scpr_utils.h b/src/test_scpr_utils.h
@@ -104,7 +104,9 @@ check_stability
(const struct scpr_polygon* polygon)
{
res_T res = RES_OK;
- double scale = 1e6;
+ /* to optimize scaling / descaling precision
+ * set the scale to a power of double's radix (2) (#25) */
+ double scale = pow(2, (int)log2(pow(10, 6)) + 1);
size_t i, j, ccount, vcount;
ASSERT(polygon);