star-enclosures-2d

Extract enclosures from 2D geometry
git clone git://git.meso-star.fr/star-enclosures-2d.git
Log | Files | Refs | README | LICENSE

commit b09c9811dd3e557e0735a7a8b9f80860c9bea379
parent 8aba9d351ab42cc5db3da4c68346e627b45ac245
Author: Christophe Coustet <christophe.coustet@meso-star.com>
Date:   Wed, 19 Feb 2020 10:28:03 +0100

Add stuff to ease a possible change in API types

Diffstat:
Mcmake/CMakeLists.txt | 2++
Msrc/senc2d.h | 2+-
Msrc/senc2d_descriptor.c | 57+++++++++++++++++++++++++--------------------------------
Msrc/senc2d_enclosure.c | 23++++++++++-------------
Msrc/senc2d_enclosure_c.h | 3++-
Msrc/senc2d_enclosure_data.h | 17++++++++++-------
Msrc/senc2d_internal_types.h | 81++++++++++++++++++++++++++++++++++++++++++++-----------------------------------
Msrc/senc2d_scene.c | 113+++++++++++++++++++++++++++++++++++++------------------------------------------
Msrc/senc2d_scene_analyze.c | 60+++++++++++++++++++++++++++++++-----------------------------
Msrc/senc2d_scene_analyze_c.h | 6+++---
Msrc/senc2d_scene_c.h | 2+-
Msrc/test_senc2d_utils.h | 6++----
Msrc/test_senc2d_utils2.h | 18+++++++++---------
13 files changed, 194 insertions(+), 196 deletions(-)

diff --git a/cmake/CMakeLists.txt b/cmake/CMakeLists.txt @@ -108,6 +108,8 @@ add_library(senc2d SHARED target_link_libraries(senc2d RSys Star2D ${MATH_LIB}) set_target_properties(senc2d PROPERTIES +# C99 needed in case of printf %zu used +# C_STANDARD 99 DEFINE_SYMBOL SENC2D_SHARED_BUILD VERSION ${VERSION} COMPILE_FLAGS ${OpenMP_C_FLAGS} diff --git a/src/senc2d.h b/src/senc2d.h @@ -87,7 +87,7 @@ struct senc2d_enclosure_header { unsigned enclosed_media_count; /* Is the enclosure open/infinite? * Only the outermost enclosure is infinite. */ - char is_infinite; + int is_infinite; }; /* We consider the geometrical normal Ng to a segment V0 V1 that verifies diff --git a/src/senc2d_descriptor.c b/src/senc2d_descriptor.c @@ -26,34 +26,32 @@ *****************************************************************************/ res_T senc2d_scene_get_max_medium - (const struct senc2d_scene* scn, unsigned* max_medium_id) + (const struct senc2d_scene* scn, medium_id_t* max_medium_id) { if(!scn || !max_medium_id) return RES_BAD_ARG; - ASSERT(scn->next_medium_idx < UINT_MAX); /* API type */ - *max_medium_id = (unsigned)scn->next_medium_idx - 1; + *max_medium_id = scn->next_medium_idx - 1; return RES_OK; } res_T senc2d_scene_get_enclosure_count - (const struct senc2d_scene* scn, unsigned* count) + (const struct senc2d_scene* scn, enclosure_id_t* count) { - size_t tmp; if(!scn || !count) return RES_BAD_ARG; - tmp = darray_enclosure_size_get(&scn->analyze.enclosures); - ASSERT(tmp < UINT_MAX); /* API type */ - ASSERT(scn->analyze.enclosures_count == tmp); - *count = (unsigned)tmp; + ASSERT(scn->analyze.enclosures_count == + darray_enclosure_size_get(&scn->analyze.enclosures)); + *count = scn->analyze.enclosures_count; return RES_OK; } res_T senc2d_scene_get_enclosure_count_by_medium (const struct senc2d_scene* scn, - const unsigned imed, - unsigned* count) + const medium_id_t imed, + enclosure_id_t* count) { - size_t tmp, m_idx; + size_t tmp; + size_t m_idx; const struct darray_enc_id* enc_ids; if(!scn || !count || (imed != SENC2D_UNSPECIFIED_MEDIUM && imed >= scn->next_medium_idx)) @@ -64,15 +62,15 @@ senc2d_scene_get_enclosure_count_by_medium enc_ids = darray_enc_ids_array_cdata_get(&scn->analyze.enc_ids_array_by_medium) + m_idx; tmp = darray_enc_id_size_get(enc_ids); - ASSERT(tmp < UINT_MAX); /* API type */ - *count = (unsigned)tmp; + ASSERT(tmp <= ENCLOSURE_MAX__); /* API type */ + *count = (enclosure_id_t)tmp; return RES_OK; } FINLINE res_T senc2d_scene_get_enclosure (struct senc2d_scene* scn, - const unsigned idx, + const enclosure_id_t idx, struct senc2d_enclosure** out_enc) { struct senc2d_enclosure* enc; @@ -88,13 +86,13 @@ senc2d_scene_get_enclosure res_T senc2d_scene_get_enclosure_by_medium (struct senc2d_scene* scn, - const unsigned imed, - const unsigned idx, + const medium_id_t imed, + const enclosure_id_t idx, struct senc2d_enclosure** out_enc) { size_t m_idx; const struct darray_enc_id* enc_ids; - unsigned index; + enclosure_id_t index; if(!scn || !out_enc || (imed != SENC2D_UNSPECIFIED_MEDIUM && imed >= scn->next_medium_idx)) return RES_BAD_ARG; @@ -111,8 +109,8 @@ senc2d_scene_get_enclosure_by_medium res_T senc2d_scene_get_segment_enclosures (const struct senc2d_scene* scn, - const unsigned iseg, - unsigned enclosures[2]) + const seg_id_t iseg, + enclosure_id_t enclosures[2]) { const struct segment_enc* seg; int i; @@ -120,40 +118,35 @@ senc2d_scene_get_segment_enclosures || iseg >= darray_segment_enc_size_get(&scn->analyze.segments_enc)) return RES_BAD_ARG; seg = darray_segment_enc_cdata_get(&scn->analyze.segments_enc) + iseg; - FOR_EACH(i, 0, 2) { -#if (UINT_MAX < ENCLOSURE_MAX__) - ASSERT(seg->enclosure[i] < UINT_MAX); -#endif - enclosures[i] = (unsigned)seg->enclosure[i]; /* Back to API type */ - } + FOR_EACH(i, 0, 2) enclosures[i] = seg->enclosure[i]; return RES_OK; } res_T senc2d_scene_get_frontier_vertice_count (const struct senc2d_scene* scn, - unsigned* count) + vrtx_id_t* count) { size_t tmp; if(!scn || !count) return RES_BAD_ARG; tmp = darray_vrtx_id_size_get(&scn->analyze.frontiers); - ASSERT(tmp < UINT_MAX); - *count = (unsigned)tmp; + ASSERT(tmp <= VRTX_MAX__); + *count = (vrtx_id_t)tmp; /* Back to API type */ return RES_OK; } res_T senc2d_scene_get_frontier_vertex (const struct senc2d_scene* scn, - const unsigned iver, - unsigned* vrtx_id) + const vrtx_id_t iver, + vrtx_id_t* vrtx_id) { vrtx_id_t vrtx; if(!vrtx_id || !scn || iver >= darray_vrtx_id_size_get(&scn->analyze.frontiers)) return RES_BAD_ARG; vrtx = darray_vrtx_id_cdata_get(&scn->analyze.frontiers)[iver]; - *vrtx_id = (unsigned)vrtx; /* Back to API type */ + *vrtx_id = vrtx; return RES_OK; } diff --git a/src/senc2d_enclosure.c b/src/senc2d_enclosure.c @@ -45,7 +45,7 @@ enclosure_release(ref_T * ref) struct senc2d_enclosure* enclosure_create (struct senc2d_scene* scn, - const unsigned idx) + const enclosure_id_t idx) { struct senc2d_enclosure* enc; ASSERT(scn && idx < darray_enclosure_size_get(&scn->analyze.enclosures)); @@ -77,8 +77,8 @@ senc2d_enclosure_get_header res_T senc2d_enclosure_get_segment (const struct senc2d_enclosure* enclosure, - const unsigned iseg, - unsigned indices[2]) + const seg_id_t iseg, + vrtx_id_t indices[2]) { const struct side_enc* side; int i; @@ -88,17 +88,14 @@ senc2d_enclosure_get_segment ASSERT(darray_sides_enc_size_get(&enclosure->data->sides) == enclosure->data->header.primitives_count); side = darray_sides_enc_cdata_get(&enclosure->data->sides) + iseg; - FOR_EACH(i, 0, 2) { - ASSERT(side->vertice_id[i] < UINT_MAX); - indices[i] = (unsigned)side->vertice_id[i]; /* Back to API type */ - } + FOR_EACH(i, 0, 2) indices[i] = side->vertice_id[i]; return RES_OK; } res_T senc2d_enclosure_get_vertex (const struct senc2d_enclosure* enclosure, - const unsigned ivert, + const vrtx_id_t ivert, double coord[2]) { if(!enclosure || !coord @@ -119,8 +116,8 @@ senc2d_enclosure_get_vertex res_T senc2d_enclosure_get_segment_id (const struct senc2d_enclosure* enclosure, - const unsigned iseg, - unsigned* gid, + const seg_id_t iseg, + seg_id_t* gid, enum senc2d_side* sde) { const struct side_enc* side; @@ -130,7 +127,7 @@ senc2d_enclosure_get_segment_id ASSERT(darray_sides_enc_size_get(&enclosure->data->sides) == enclosure->data->header.primitives_count); side = darray_sides_enc_cdata_get(&enclosure->data->sides) + iseg; - *gid = (unsigned)SEGSIDE_2_SEG(side->side_id); + *gid = SEGSIDE_2_SEG(side->side_id); *sde = SEGSIDE_2_SIDE(side->side_id); return RES_OK; } @@ -138,8 +135,8 @@ senc2d_enclosure_get_segment_id res_T senc2d_enclosure_get_medium (const struct senc2d_enclosure* enclosure, - const unsigned imed, - unsigned* medium) + const medium_id_t imed, + medium_id_t* medium) { if(!enclosure || !medium || imed >= enclosure->data->header.enclosed_media_count) diff --git a/src/senc2d_enclosure_c.h b/src/senc2d_enclosure_c.h @@ -19,6 +19,7 @@ #include <rsys/ref_count.h> #include "senc2d.h" +#include "senc2d_internal_types.h" struct enclosure_data; struct senc2d_scene; @@ -32,6 +33,6 @@ struct senc2d_enclosure { struct senc2d_enclosure* enclosure_create (struct senc2d_scene* scene, - const unsigned idx); + const enclosure_id_t idx); #endif /* SENC2D_ENCLOSURE_C_H */ diff --git a/src/senc2d_enclosure_data.h b/src/senc2d_enclosure_data.h @@ -46,10 +46,10 @@ struct side_enc { #define DARRAY_DATA struct side_enc #include <rsys/dynamic_array.h> -/* unsigned char array with init to zero */ +/* uchar array with init to zero */ static FINLINE void zero_init_uchar - (struct mem_allocator* alloc, unsigned char* data) + (struct mem_allocator* alloc, uchar* data) { ASSERT(data); (void) alloc; @@ -77,12 +77,12 @@ init_header(struct senc2d_enclosure_header* header) static FINLINE res_T bool_array_of_media_merge (struct darray_uchar* dst, - const unsigned char* src, + const uchar* src, const medium_id_t sz) { res_T res = RES_OK; medium_id_t i; - unsigned char* data_dst; + uchar* data_dst; ASSERT(src && dst); @@ -107,8 +107,8 @@ bool_array_of_media_to_darray_media const medium_id_t next_medium_idx) { res_T res = RES_OK; - unsigned m_idx; - const unsigned char* data; + size_t m_idx; + const uchar* data; ASSERT(src && dst); @@ -118,8 +118,11 @@ bool_array_of_media_to_darray_media darray_media_clear(dst); if(res != RES_OK) goto error; FOR_EACH(m_idx, 0, next_medium_idx + 1) { - medium_id_t medium = m_idx ? (medium_id_t)(m_idx - 1) : SENC2D_UNSPECIFIED_MEDIUM; + medium_id_t medium; + size_t mm = m_idx ? (medium_id_t)(m_idx - 1) : SENC2D_UNSPECIFIED_MEDIUM; if(!data[m_idx]) continue; + ASSERT(mm <= MEDIUM_MAX__); + medium = (medium_id_t)mm; res = darray_media_push_back(dst, &medium); if(res != RES_OK) goto error; } diff --git a/src/senc2d_internal_types.h b/src/senc2d_internal_types.h @@ -46,51 +46,60 @@ } #endif -/* Side IDs are uint32_t */ -typedef uint32_t side_id_t; -#define SIDE_MAX__ (UINT32_MAX-1) -#define SIDE_NULL__ UINT32_MAX - -/* Seg IDs use internally side_id_t */ -/* Cannot be larger than unsigned, as the API uses it. */ -typedef side_id_t seg_id_t; -/* SEG_MAX__ is limited to allow to count sides */ -#define SEG_MAX__ (SIDE_MAX__ / 2) -#define SEG_NULL__ UINT32_MAX - -/* Vertex IDs are internally uint32_t */ -/* Cannot be larger than unsigned, as the API uses it. */ -typedef uint32_t vrtx_id_t; -#define VRTX_MAX__ (UINT32_MAX-1) -#define VRTX_NULL__ UINT32_MAX +/* Helper type */ +typedef unsigned char uchar; + +/* The following types must be defined accordingly with the types + * used in senc2d.h */ + +/* Segment IDs use the same type than Side IDs */ +typedef unsigned seg_id_t; +/* SEG_MAX__ is limited to half the max of the base type to allow to count + * sides */ +#define SEG_MAX__ (UINT_MAX/2) +#define SEG_NULL__ UINT_MAX +#define PRTF_SEG "%u" + +/* Side IDs type use the same base type than Segment IDs */ +typedef seg_id_t side_id_t; +#define SIDE_MAX__ (2*SEG_MAX__) +#define SIDE_NULL__ SEG_NULL__ + +/* Vertex IDs type */ +typedef unsigned vrtx_id_t; +#define VRTX_MAX__ (UINT_MAX-1) +#define VRTX_NULL__ UINT_MAX +#define PRTF_VRTX "%u" /* Edge IDs use the same type than vertex IDs */ -/* Cannot be larger than unsigned, as the API uses it. */ typedef vrtx_id_t edge_id_t; #define EDGE_MAX__ VRTX_MAX__ #define EDGE_NULL__ VRTX_NULL__ -/* Medium IDs are internally uint32_t */ -/* Should nnot be larger than unsigned, as the API uses it. */ -typedef uint32_t medium_id_t; -#define MEDIUM_MAX__ INT32_MAX -#define MEDIUM_NULL__ UINT32_MAX +/* Medium IDs type */ +typedef unsigned medium_id_t; +#define MEDIUM_MAX__ (UINT_MAX-1) /* MAX is for unspecified medium */ +#define MEDIUM_NULL__ UINT_MAX +#define PRTF_MDM "%u" -/* Enclosure IDs are internally uint32_t */ -/* Cannot be larger than unsigned, as the API uses it. */ -typedef uint32_t enclosure_id_t; -#define ENCLOSURE_MAX__ UINT32_MAX -#define ENCLOSURE_NULL__ UINT32_MAX +/* Enclosure IDs type */ +typedef unsigned enclosure_id_t; +#define ENCLOSURE_MAX__ (UINT_MAX-1) +#define ENCLOSURE_NULL__ UINT_MAX /* Component IDs use the same type than enclosure IDs */ typedef enclosure_id_t component_id_t; -#define COMPONENT_MAX__ (UINT32_MAX - 2) /* To allow special values */ -#define COMPONENT_NULL__ UINT32_MAX +#define COMPONENT_MAX__ (UINT_MAX-2) /* To allow special values */ +#define COMPONENT_NULL__ UINT_MAX /* Special values */ -#define CC_GROUP_ROOT_NONE UINT32_MAX -#define CC_GROUP_ROOT_INFINITE (UINT32_MAX - 1) -#define CC_GROUP_ID_NONE UINT32_MAX -#define CC_ID_NONE UINT32_MAX +#define CC_GROUP_ROOT_NONE UINT_MAX +#define CC_GROUP_ROOT_INFINITE (UINT_MAX-1) +#define CC_GROUP_ID_NONE UINT_MAX +#define CC_ID_NONE UINT_MAX + +#if (MEDIUM_MAX__+1 != SENC2D_UNSPECIFIED_MEDIUM) +#error "Inconsistant values" +#endif /* This one is used as flag */ enum side_flag { @@ -120,10 +129,10 @@ SEGSIDE_2_SIDEFLAG(side_id_t s) { return (s & 1) ? FLAG_BACK : FLAG_FRONT; } -static FINLINE unsigned char +static FINLINE uchar SIDE_CANCELED_FLAG(enum side_flag f) { ASSERT((f << 4) <= UCHAR_MAX); - return (unsigned char)(f << 4); + return (uchar)(f << 4); } static FINLINE side_id_t diff --git a/src/senc2d_scene.c b/src/senc2d_scene.c @@ -53,7 +53,8 @@ compatible_medium (const medium_id_t m1, const medium_id_t m2) { - if(m1 == SENC2D_UNSPECIFIED_MEDIUM || m2 == SENC2D_UNSPECIFIED_MEDIUM) return 1; + if(m1 == SENC2D_UNSPECIFIED_MEDIUM || m2 == SENC2D_UNSPECIFIED_MEDIUM) + return 1; return (m1 == m2); } @@ -64,11 +65,11 @@ res_T senc2d_scene_create (struct senc2d_device* dev, const int conv, - const unsigned nsegs, - void(*indices)(const unsigned, unsigned*, void*), - void(*media)(const unsigned, unsigned*, void*), - const unsigned nverts, - void(*position)(const unsigned, double*, void* ctx), + const seg_id_t nsegs, + void(*indices)(const seg_id_t, vrtx_id_t*, void*), + void(*media)(const seg_id_t, medium_id_t*, void*), + const vrtx_id_t nverts, + void(*position)(const vrtx_id_t, double*, void* ctx), void* ctx, struct senc2d_scene** out_scn) { @@ -76,10 +77,11 @@ senc2d_scene_create /* Tables to detect duplicates */ struct htable_vrtx unique_vertices; struct htable_seg unique_segments; - unsigned i; + vrtx_id_t nv; + seg_id_t ns; res_T res = RES_OK; - if(!dev || !out_scn || !indices || !position + if(!dev || !out_scn || !indices || !position || !nverts || !nsegs /* Convention must be set both regarding FRONT/BACK and INSIDE/OUTSIDE */ || !(conv & (SENC2D_CONVENTION_NORMAL_FRONT | SENC2D_CONVENTION_NORMAL_BACK)) || !(conv & (SENC2D_CONVENTION_NORMAL_INSIDE | SENC2D_CONVENTION_NORMAL_OUTSIDE))) @@ -113,90 +115,90 @@ senc2d_scene_create OK(darray_enclosure_resize(&scn->analyze.enclosures, 1)); scn->analyze.enclosures_count = 1; - if(!scn || !indices || !position || !nverts || !nsegs) { - res = RES_BAD_ARG; - goto error; - } - OK(darray_position_reserve(&scn->vertices, scn->nverts)); OK(darray_segment_in_reserve(&scn->segments_in, scn->nsegs)); OK(htable_vrtx_reserve(&unique_vertices, scn->nverts)); OK(htable_seg_reserve(&unique_segments, scn->nsegs)); /* Get vertices */ - FOR_EACH(i, 0, nverts) { + FOR_EACH(nv, 0, nverts) { vrtx_id_t* p_vrtx; union double2 tmp; - /* API: position needs an unsigned */ - position(i, tmp.vec, ctx); + position(nv, tmp.vec, ctx); p_vrtx = htable_vrtx_find(&unique_vertices, &tmp); if(p_vrtx) { /* Duplicate vertex */ - log_err(scn->dev, "%s: vertex %u is a duplicate.\n", - FUNC_NAME, i); + log_err(scn->dev, "%s: vertex "PRTF_VRTX" is a duplicate.\n", + FUNC_NAME, nv); res = RES_BAD_ARG; goto error; } /* New vertex */ - ASSERT(i == htable_vrtx_size_get(&unique_vertices)); + ASSERT(nv == htable_vrtx_size_get(&unique_vertices)); OK(darray_position_push_back(&scn->vertices, &tmp)); - OK(htable_vrtx_set(&unique_vertices, &tmp, &i)); + OK(htable_vrtx_set(&unique_vertices, &tmp, &nv)); } /* Get segments */ - FOR_EACH(i, 0, nsegs) { + FOR_EACH(ns, 0, nsegs) { int j; - unsigned med[2] = { SENC2D_UNSPECIFIED_MEDIUM, SENC2D_UNSPECIFIED_MEDIUM }; - unsigned ind[2]; + seg_id_t s; + medium_id_t med[2] + = { SENC2D_UNSPECIFIED_MEDIUM, SENC2D_UNSPECIFIED_MEDIUM }; + vrtx_id_t ind[2]; union vrtx_id2 seg_key; struct segment_in tmp; seg_id_t* p_seg; - indices(i, ind, ctx); /* API: indices need unsigneds */ + indices(ns, ind, ctx); FOR_EACH(j, 0, 2) { if(ind[j] >= nverts) { - log_err(scn->dev, "%s: segment %u uses invalid vertex id %u.\n", - FUNC_NAME, i, ind[j]); + log_err(scn->dev, + "%s: segment "PRTF_SEG" uses invalid vertex id "PRTF_VRTX".\n", + FUNC_NAME, ns, ind[j]); res = RES_BAD_ARG; goto error; } ASSERT(ind[j] <= VRTX_MAX__); - tmp.vertice_id[j] = (vrtx_id_t)ind[j]; + tmp.vertice_id[j] = ind[j]; } - if(tmp.vertice_id[0] == tmp.vertice_id[1]) - { - log_err(scn->dev, "%s: segment %u is degenerated.\n", - FUNC_NAME, i); + if(tmp.vertice_id[0] == tmp.vertice_id[1]) { + log_err(scn->dev, "%s: segment "PRTF_SEG" is degenerated.\n", + FUNC_NAME, ns); res = RES_BAD_ARG; goto error; } /* Get media */ - if(media) media(i, med, ctx); /* API: media needs an unsigned */ + if(media) media(ns, med, ctx); seg_make_key(&seg_key, tmp.vertice_id); p_seg = htable_seg_find(&unique_segments, &seg_key); if(p_seg) { /* Duplicate segment */ - log_err(scn->dev, "%s: segment %u is a duplicate.\n", - FUNC_NAME, i); + log_err(scn->dev, "%s: segment "PRTF_SEG" is a duplicate.\n", + FUNC_NAME, ns); res = RES_BAD_ARG; goto error; } /* New segment */ - ASSERT(i == htable_seg_size_get(&unique_segments)); - OK(htable_seg_set(&unique_segments, &seg_key, &i)); - FOR_EACH(j, 0, 2) { + ASSERT(ns == htable_seg_size_get(&unique_segments)); + OK(htable_seg_set(&unique_segments, &seg_key, &ns)); + for(s = SENC2D_FRONT; s <= SENC2D_BACK; s += (SENC2D_BACK - SENC2D_FRONT)) { struct side_range* media_use; - unsigned m_idx = (med[j] == SENC2D_UNSPECIFIED_MEDIUM) ? 0 : med[j] + 1; - tmp.medium[j] = (medium_id_t)med[j]; + size_t m_idx = (med[s] == SENC2D_UNSPECIFIED_MEDIUM) ? 0 : med[s] + 1; + ASSERT(m_idx <= MEDIUM_MAX__); + tmp.medium[s] = med[s]; if(m_idx >= scn->next_medium_idx) { - scn->next_medium_idx = m_idx; + medium_id_t medium; + ASSERT(m_idx <= MEDIUM_MAX__); + medium = (medium_id_t)m_idx; + scn->next_medium_idx = medium; darray_side_range_resize(&scn->media_use, 1 + m_idx); } /* media_use 0 is for SENC2D_UNSPECIFIED_MEDIUM */ media_use = darray_side_range_data_get(&scn->media_use) + m_idx; media_use->first = - MMIN(media_use->first, SEGIDxSIDE_2_SEGSIDE((seg_id_t)i, j)); + MMIN(media_use->first, SEGIDxSIDE_2_SEGSIDE(ns, s)); ASSERT(media_use->first < 2 * (scn->nsegs + 1)); media_use->last = - MMAX(media_use->last, SEGIDxSIDE_2_SEGSIDE((seg_id_t)i, j)); + MMAX(media_use->last, SEGIDxSIDE_2_SEGSIDE(ns, s)); ASSERT(media_use->last < 2 * (scn->nsegs + 1)); ASSERT(media_use->first <= media_use->last); } @@ -235,7 +237,7 @@ senc2d_scene_get_convention res_T senc2d_scene_get_segments_count (const struct senc2d_scene* scn, - unsigned* count) + seg_id_t* count) { if(!scn || !count) return RES_BAD_ARG; *count = scn->nsegs; @@ -245,8 +247,8 @@ senc2d_scene_get_segments_count res_T senc2d_scene_get_segment (const struct senc2d_scene* scn, - const unsigned iseg, - unsigned indices[2]) + const seg_id_t iseg, + vrtx_id_t indices[2]) { const struct segment_in* seg; int i; @@ -254,19 +256,15 @@ senc2d_scene_get_segment || iseg >= darray_segment_in_size_get(&scn->segments_in)) return RES_BAD_ARG; seg = darray_segment_in_cdata_get(&scn->segments_in) + iseg; - - FOR_EACH(i, 0, 2) { - ASSERT(seg->vertice_id[i] < UINT_MAX); - indices[i] = (unsigned)seg->vertice_id[i]; /* Back to API type */ - } + FOR_EACH(i, 0, 2) indices[i] = seg->vertice_id[i]; return RES_OK; } res_T senc2d_scene_get_segment_media (const struct senc2d_scene* scn, - const unsigned iseg, - unsigned media[2]) + const seg_id_t iseg, + medium_id_t media[2]) { const struct segment_in* seg; int i; @@ -274,18 +272,14 @@ senc2d_scene_get_segment_media || iseg >= darray_segment_in_size_get(&scn->segments_in)) return RES_BAD_ARG; seg = darray_segment_in_cdata_get(&scn->segments_in) + iseg; - - FOR_EACH(i, 0, 2) { - ASSERT(seg->vertice_id[i] < UINT_MAX); - media[i] = (unsigned)seg->medium[i]; /* Back to API type */ - } + FOR_EACH(i, 0, 2) media[i] = seg->medium[i]; return RES_OK; } res_T senc2d_scene_get_vertices_count (const struct senc2d_scene* scn, - unsigned* count) + vrtx_id_t* count) { if(!scn || !count) return RES_BAD_ARG; *count = scn->nverts; @@ -295,14 +289,13 @@ senc2d_scene_get_vertices_count res_T senc2d_scene_get_vertex (const struct senc2d_scene* scn, - const unsigned ivert, + const vrtx_id_t ivert, double coord[2]) { const union double2* v; if(!scn || !coord || ivert >= darray_position_size_get(&scn->vertices)) return RES_BAD_ARG; - v = darray_position_cdata_get(&scn->vertices) + ivert; d2_set(coord, v->vec); return RES_OK; diff --git a/src/senc2d_scene_analyze.c b/src/senc2d_scene_analyze.c @@ -69,7 +69,7 @@ static side_id_t get_side_not_in_connex_component (const side_id_t last_side, const struct segside* segsides, - const unsigned char* processed, + const uchar* processed, side_id_t* first_side_not_in_component, const medium_id_t medium) { @@ -87,6 +87,7 @@ get_side_not_in_connex_component } } +/* Here unsigned are required by s2d API */ static void get_scn_indices(const unsigned iseg, unsigned ids[2], void* ctx) { int i; @@ -95,10 +96,12 @@ get_scn_indices(const unsigned iseg, unsigned ids[2], void* ctx) { darray_segment_in_cdata_get(&scene->segments_in) + iseg; FOR_EACH(i, 0, 2) { ASSERT(seg->vertice_id[i] < scene->nverts); - ids[i] = (unsigned)seg->vertice_id[i]; /* Back to API type */ + ASSERT(seg->vertice_id[i] <= UINT_MAX); + ids[i] = (unsigned)seg->vertice_id[i]; /* Back to s2d API type */ } } +/* Here unsigned are required by s2d API */ static void get_scn_position(const unsigned ivert, float pos[2], void* ctx) { const struct senc2d_scene* scene = ctx; @@ -126,7 +129,6 @@ self_hit_filter + hit->prim.prim_id; return (hit_seg_comp->component[SENC2D_FRONT] == *origin_component || hit_seg_comp->component[SENC2D_BACK] == *origin_component); - } static void @@ -145,18 +147,18 @@ extract_connex_components /* This function is called from an omp parallel block and executed * concurrently. */ struct mem_allocator* alloc; - int64_t m_idx; + int64_t m_idx; /* OpenMP requires a signed type for the for loop variable */ struct darray_side_id stack; struct darray_side_id ids_of_sides_around_max_y_vertex; const union double2* positions; const struct segment_tmp* segments_tmp; struct segment_comp* segments_comp; /* An array to flag sides when processed */ - unsigned char* processed; + uchar* processed; /* An array to store the component being processed */ struct darray_side_id current_component; /* A bool array to store media of the component being processed */ - unsigned char* current_media = NULL; + uchar* current_media = NULL; size_t sz, ii; ASSERT(scn && segsides && connex_components && segments_tmp_array @@ -168,7 +170,7 @@ extract_connex_components darray_side_id_init(alloc, &stack); darray_side_id_init(alloc, &ids_of_sides_around_max_y_vertex); darray_side_id_init(alloc, &current_component); - processed = MEM_CALLOC(alloc, scn->nsegs, sizeof(unsigned char)); + processed = MEM_CALLOC(alloc, scn->nsegs, sizeof(uchar)); if(!processed) { *p_res = RES_MEM_ERR; return; @@ -265,7 +267,7 @@ extract_connex_components const seg_id_t crt_seg_id = SEGSIDE_2_SEG(crt_side_id); const struct segment_in* seg_in = darray_segment_in_cdata_get(&scn->segments_in) + crt_seg_id; - unsigned char* seg_used = processed + crt_seg_id; + uchar* seg_used = processed + crt_seg_id; const struct segment_tmp* const seg_tmp = segments_tmp + crt_seg_id; ASSERT(crt_seg_id < scn->nsegs); @@ -308,7 +310,7 @@ extract_connex_components /* Record crt_side both as component and segment level */ if((*seg_used & crt_side_flag) == 0) { OK2(darray_side_id_push_back(&current_component, &crt_side_id)); - *seg_used = *seg_used | (unsigned char)crt_side_flag; + *seg_used = *seg_used | (uchar)crt_side_flag; } /* Store neighbour's sides in a waiting stack */ @@ -316,7 +318,7 @@ extract_connex_components side_id_t neighbour_id = crt_side->facing_side_id[i]; seg_id_t nbour_seg_id = SEGSIDE_2_SEG(neighbour_id); enum side_flag nbour_side_id = SEGSIDE_2_SIDEFLAG(neighbour_id); - unsigned char* nbour_used = processed + nbour_seg_id; + uchar* nbour_used = processed + nbour_seg_id; const struct segside* neighbour = segsides + neighbour_id; medium_id_t nbour_med_idx = (neighbour->medium == SENC2D_UNSPECIFIED_MEDIUM) ? 0 : neighbour->medium + 1; @@ -342,8 +344,8 @@ extract_connex_components seg_id_t used_seg_id = SEGSIDE_2_SEG(used_side); enum side_flag used_side_flag = SEGSIDE_2_SIDEFLAG(used_side); - unsigned char* used = processed + used_seg_id; - ASSERT(*used & (unsigned char)used_side_flag); + uchar* used = processed + used_seg_id; + ASSERT(*used & (uchar)used_side_flag); /* Set the used flag for sides in cancelled component as leading * to further cancellations */ *used |= SIDE_CANCELED_FLAG(used_side_flag); @@ -353,7 +355,7 @@ extract_connex_components } if(*nbour_used & nbour_side_id) continue; /* Already processed */ /* Mark neighbour as processed and stack it */ - *nbour_used |= (unsigned char)nbour_side_id; + *nbour_used |= (uchar)nbour_side_id; OK2(darray_side_id_push_back(&stack, &neighbour_id)); OK2(darray_side_id_push_back(&current_component, &neighbour_id)); current_media[nbour_med_idx] = 1; @@ -497,7 +499,7 @@ extract_connex_components OK(s2d_scene_create(s2d, &s2d_scn)); OK(s2d_shape_create_line_segments(s2d, &s2d_shp)); - /* Back to API type for nsegs and nverts */ + /* Back to s2d API type for nsegs and nverts */ ASSERT(scn->nsegs < UINT_MAX); ASSERT(scn->nverts < UINT_MAX); OK(s2d_line_segments_setup_indexed_vertices(s2d_shp, @@ -530,7 +532,7 @@ extract_connex_components ASSERT(seg_comp->component[SENC2D_FRONT] != COMPONENT_NULL__); ASSERT(seg_comp->component[SENC2D_BACK] != COMPONENT_NULL__); } - FOR_EACH(c, 0, ATOMIC_GET(component_count)) { + FOR_EACH(c, 0, (component_id_t)ATOMIC_GET(component_count)) { struct cc_descriptor** components = darray_ptr_component_descriptor_data_get(connex_components); ASSERT(components[c] != NULL && components[c]->cc_id == c); @@ -647,7 +649,7 @@ group_connex_components } else { struct enclosure_data* enclosures = darray_enclosure_data_get(&scn->analyze.enclosures); - FOR_EACH(ccc, 0, cc_count) { + FOR_EACH(ccc, 0, (int64_t)cc_count) { component_id_t c = (component_id_t)ccc; struct cc_descriptor* const cc = descriptors[c]; const struct cc_descriptor* other_desc = cc; @@ -726,7 +728,7 @@ collect_and_link_neighbours /* Loop on segments to collect edges' neighbours. * All threads considering all the vertices and processing some */ FOR_EACH(s, 0, scn->nsegs) { - unsigned char vv; + uchar vv; FOR_EACH(vv, 0, 2) { struct darray_neighbour* neighbourhood; struct neighbour_info* info; @@ -840,11 +842,11 @@ collect_and_link_neighbours const struct neighbour_info* ccw_neighbour = darray_neighbour_cdata_get(neighbourhood) + (i + 1) % neighbour_count; /* Rank of the end of interest in segments */ - const unsigned char crt_end = current->common_vertex_rank; + const uchar crt_end = current->common_vertex_rank; /* Here ccw refers to the rotation around the common vertex * and has nothing to do with vertices order in segment definition * nor Front/Back side convention */ - const unsigned char ccw_end = ccw_neighbour->common_vertex_rank; + const uchar ccw_end = ccw_neighbour->common_vertex_rank; /* User id of current segments */ const seg_id_t crt_id = current->seg_id; const seg_id_t ccw_id = ccw_neighbour->seg_id; @@ -871,8 +873,8 @@ collect_and_link_neighbours previous = darray_neighbour_cdata_get(neighbourhood) + i - 1; prev_id = previous->seg_id; log_err(scn->dev, - "%s: found 2 overlying segments (%u & %u).\n", FUNC_NAME, - (unsigned)crt_id, (unsigned)prev_id); + "%s: found 2 overlying segments ("PRTF_SEG" & "PRTF_SEG").\n", + FUNC_NAME, crt_id, prev_id); tmp_res = RES_BAD_OP; goto tmp_error; } @@ -901,8 +903,8 @@ collect_and_link_neighbours #pragma omp critical { log_warn(scn->dev, - "%s: found frontier involving segment %lu.\n", - FUNC_NAME, (unsigned long)crt_id); + "%s: found frontier involving segment "PRTF_SEG".\n", + FUNC_NAME, crt_id); darray_vrtx_id_push_back(frontiers, &v); } } @@ -998,24 +1000,24 @@ build_result == enc->first_component); if(*res != RES_OK) continue; - ASSERT(e <= UINT_MAX); - enc->header.enclosure_id = (unsigned)e; /* Back to API type */ + ASSERT(e <= ENCLOSURE_MAX__); + enc->header.enclosure_id = (enclosure_id_t)e; /* Back to API type */ ASSERT(cc_descriptors[enc->first_component]->enclosure_id == enc->header.enclosure_id); enc->header.is_infinite = (e == 0); - ASSERT(darray_uchar_size_get(&enc->tmp_enclosed_media) <= UINT_MAX); ASSERT(enc->header.enclosed_media_count < 1 + scn->next_medium_idx); OK2(bool_array_of_media_to_darray_media (&enc->enclosed_media, &enc->tmp_enclosed_media, scn->next_medium_idx)); + ASSERT(darray_media_size_get(&enc->enclosed_media) <= MEDIUM_MAX__); enc->header.enclosed_media_count - = (unsigned)darray_media_size_get(&enc->enclosed_media); + = (medium_id_t)darray_media_size_get(&enc->enclosed_media); darray_uchar_purge(&enc->tmp_enclosed_media); /* Add this enclosure in relevant by-medium lists */ FOR_EACH(m, 0, enc->header.enclosed_media_count) { medium_id_t medium = darray_media_cdata_get(&enc->enclosed_media)[m]; - unsigned m_idx = (medium == SENC2D_UNSPECIFIED_MEDIUM) ? 0 : medium + 1; + size_t m_idx = (medium == SENC2D_UNSPECIFIED_MEDIUM) ? 0 : medium + 1; struct darray_enc_id* enc_ids_array_by_medium; ASSERT(medium == SENC2D_UNSPECIFIED_MEDIUM || medium < scn->next_medium_idx); ASSERT(darray_enc_ids_array_size_get(&scn->analyze.enc_ids_array_by_medium) @@ -1049,7 +1051,7 @@ build_result { const struct segment_in* seg_in = segments_in + s; struct side_enc* side_enc; - unsigned vertice_id[2]; + vrtx_id_t vertice_id[2]; int i; if(segments_enc[s].enclosure[SENC2D_FRONT] != e && segments_enc[s].enclosure[SENC2D_BACK] != e) diff --git a/src/senc2d_scene_analyze_c.h b/src/senc2d_scene_analyze_c.h @@ -54,7 +54,7 @@ struct cc_descriptor { /* Range of sides member of this component */ struct side_range side_range; /* Media used by this component */ - unsigned char* media; + uchar* media; }; extern const struct cc_descriptor CC_DESCRIPTOR_NULL; @@ -132,10 +132,10 @@ struct neighbour_info { double angle; seg_id_t seg_id; /* Rank of the vertex in the segment (in [0 1]) */ - unsigned char common_vertex_rank; + uchar common_vertex_rank; /* Does the geometrical normal point towards the next neighbour * (if not, it points towards the previous one)? */ - unsigned char normal_toward_next_neighbour; + uchar normal_toward_next_neighbour; }; #define DARRAY_NAME neighbour #define DARRAY_DATA struct neighbour_info diff --git a/src/senc2d_scene_c.h b/src/senc2d_scene_c.h @@ -84,7 +84,7 @@ set_edge (const vrtx_id_t vrtx0, const vrtx_id_t vrtx1, struct seg_edge* edge, - unsigned char* reversed) + uchar* reversed) { ASSERT(edge && reversed && vrtx0 != vrtx1); ASSERT(*reversed == UCHAR_MAX); /* Should not be already set. */ diff --git a/src/test_senc2d_utils.h b/src/test_senc2d_utils.h @@ -152,8 +152,7 @@ dump_global FOR_EACH(i, 0, segments_count) { unsigned indices[2]; OK(senc2d_scene_get_segment(scn, i, indices)); - fprintf(stream, "l %lu %lu\n", - (unsigned long)(1 + indices[0]), (unsigned long)(1 + indices[1])); + fprintf(stream, "l %u %u\n", 1 + indices[0], 1 + indices[1]); } fclose(stream); } @@ -186,8 +185,7 @@ dump_enclosure FOR_EACH(i, 0, header.primitives_count) { unsigned indices[2]; OK(senc2d_enclosure_get_segment(enclosure, i, indices)); - fprintf(stream, "l %lu %lu\n", - (unsigned long)(1+indices[0]), (unsigned long)(1+indices[1])); + fprintf(stream, "l %u %u\n", 1+indices[0], 1+indices[1]); } OK(senc2d_enclosure_ref_put(enclosure)); fclose(stream); diff --git a/src/test_senc2d_utils2.h b/src/test_senc2d_utils2.h @@ -29,25 +29,25 @@ static void get_ctx_indices(const unsigned iseg, unsigned ids[2], void* context) { struct context* ctx = context; - unsigned s3dut_iseg, circ_idx, v_offset; + unsigned template_iseg, circ_idx, v_offset; unsigned circ_seg_count, circ_vrtx_count; ASSERT(ids && ctx); - /* Get circ_idx along with the s3dut vertice index */ + /* Get circ_idx along with the template vertice index */ circ_seg_count = (unsigned)sa_size(ctx->indices) / 2; circ_vrtx_count = (unsigned)sa_size(ctx->positions) / 2; ASSERT(iseg < NB_CIRC * circ_seg_count); - s3dut_iseg = iseg % circ_seg_count; + template_iseg = iseg % circ_seg_count; circ_idx = iseg / circ_seg_count; - ASSERT(ctx->indices[s3dut_iseg * 2 + 0] <= UINT_MAX - && ctx->indices[s3dut_iseg * 2 + 1] <= UINT_MAX); + ASSERT(ctx->indices[template_iseg * 2 + 0] <= UINT_MAX + && ctx->indices[template_iseg * 2 + 1] <= UINT_MAX); /* Compute the vertex index in the user numbering - * from circ_idx and s3dut data; vertex related getters - * will have to get the s3dut index back */ + * from circ_idx and template data; vertex related getters + * will have to get the template index back */ v_offset = circ_idx * circ_vrtx_count; ids[ctx->reverse_vrtx ? 1 : 0] - = v_offset + (unsigned)ctx->indices[s3dut_iseg * 2 + 0]; + = v_offset + (unsigned)ctx->indices[template_iseg * 2 + 0]; ids[ctx->reverse_vrtx ? 0 : 1] - = v_offset + (unsigned)ctx->indices[s3dut_iseg * 2 + 1]; + = v_offset + (unsigned)ctx->indices[template_iseg * 2 + 1]; } static void