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 9db091b7aac5828c7ba2d185c0bed7895cf66b5b
parent 29f84c17e3fa7362a5e51c2b1a7fb5bf4082f655
Author: Christophe Coustet <christophe.coustet@meso-star.com>
Date:   Tue,  5 Nov 2019 10:38:39 +0100

Merge branch 'feature_segment_app_data' into develop

Diffstat:
Mcmake/CMakeLists.txt | 1+
Msrc/senc2d.h | 60++++++++++++++++++++++++++++++++++++------------------------
Msrc/senc2d_descriptor.c | 15---------------
Msrc/senc2d_enclosure.c | 47+++++++++++++----------------------------------
Msrc/senc2d_enclosure_data.h | 31++++++++++++++++++++-----------
Msrc/senc2d_internal_types.h | 25++++++++++---------------
Msrc/senc2d_scene.c | 35+++++++++++++++++++----------------
Msrc/senc2d_scene_analyze.c | 133++++++++++++++++++++++++-------------------------------------------------------
Msrc/senc2d_scene_c.h | 3---
Msrc/test_senc2d_descriptor.c | 23++++++-----------------
Msrc/test_senc2d_enclosure.c | 80++++++++++++++++++++++++++++++++++++++++++-------------------------------------
Msrc/test_senc2d_inconsistant_square.c | 20+++++++++-----------
Msrc/test_senc2d_many_enclosures.c | 2+-
Msrc/test_senc2d_many_segments.c | 2+-
Msrc/test_senc2d_sample_enclosure.c | 2+-
Msrc/test_senc2d_scene.c | 61+++++++++++++++++++++++--------------------------------------
Msrc/test_senc2d_square_behind_square.c | 12++++++------
Msrc/test_senc2d_square_in_square.c | 12++++++------
Msrc/test_senc2d_square_on_square.c | 6+++---
Msrc/test_senc2d_undefined_medium.c | 33++++++++++++++++-----------------
Asrc/test_senc2d_undefined_medium_attr.c | 366+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Msrc/test_senc2d_utils.h | 18++++--------------
22 files changed, 624 insertions(+), 363 deletions(-)

diff --git a/cmake/CMakeLists.txt b/cmake/CMakeLists.txt @@ -131,6 +131,7 @@ if(NOT NO_TEST) new_test(test_senc2d_sample_enclosure) new_test(test_senc2d_scene) new_test(test_senc2d_undefined_medium) + new_test(test_senc2d_undefined_medium_attr) target_link_libraries(test_senc2d_sample_enclosure StarSP) rcmake_copy_runtime_libraries(test_senc2d_sample_enclosure) diff --git a/src/senc2d.h b/src/senc2d.h @@ -18,6 +18,8 @@ #include <rsys/rsys.h> +#include <limits.h> + /* Library symbol management */ #if defined(SENC2D_SHARED_BUILD) #define SENC2D_API extern EXPORT_SYM @@ -59,6 +61,12 @@ struct senc2d_device; struct senc2d_scene; struct senc2d_enclosure; +/* A type to discriminate segment sides */ +enum senc2d_side { + SENC2D_FRONT, + SENC2D_BACK +}; + /* Enclosure2D header type */ struct senc2d_enclosure_header { /* The ID of the enclosure; 0, 1, ... */ @@ -169,20 +177,37 @@ senc2d_scene_reserve * of SENC_UNDEFINED_MEDIUM that causes no error) and are deduplicated. * The special value SENC2D_UNDEFINED_MEDIUM denotes an undefined medium. * It can be used to define the 2 sides of a segment at different times. - * When deduplicating segments, the first occurence is kept (with its original - * global_id). Users can provide their own global ids for segments; these ids - * are not used by the library but are returned as-is by some API calls. */ + * When deduplicating segments, the first occurence remains. + * The add_segment and merge_segment callbacks can be used for attributes + * management including segment IDs; they allow the client app to store + * its own data. They can also fail and stop the add_geometry call. */ SENC2D_API res_T senc2d_scene_add_geometry (struct senc2d_scene* scene, + /* Number of added segments */ const unsigned segments_count, + /* User function that provides vertices ids for added segments */ void(*indices)(const unsigned iseg, unsigned ids[2], void* context), + /* User function that provides media ids for added segments */ void(*media) /* Can be NULL <=> SENC2D_UNDEFINED_MEDIUM medium used */ (const unsigned iseg, unsigned med[2], void* context), - void(*global_id) /* Can be NULL <=> use segment rank */ - (const unsigned iseg, unsigned* gid, void* context), + /* Number of added vertices */ const unsigned vertices_count, + /* User function that provides coordinates for added vertices */ void(*position)(const unsigned ivert, double pos[2], void* context), + /* Called for each new segment so that the client app can manage its own + * segment data/properties/attributes. + * If return is not RES_OK, add_geometry stops and fails. */ + res_T(*add_segment) /* Can be NULL */ + (const unsigned global_id, const unsigned iseg, void* context), + /* Called if the iseg_th segment of the current add_geometry is equal to + * the global_id_th global segment so that the client app can try to merge + * its own segment data. The reversed_segment arg indicates if the segment + * vertices' order is the same it was when the segment was first added. + * If return is not RES_OK, add_geometry stops and fails. */ + res_T(*merge_segment) /* Can be NULL */ + (const unsigned global_id, const unsigned iseg, const int reversed_segment, + void* context), void* context); /* Returns a descriptor of the scene that holds the analysis' result. */ @@ -347,14 +372,6 @@ senc2d_descriptor_get_global_segment_enclosures const unsigned iseg, unsigned enclosures[2]); -/* Returns the global id of the iseg_th global unique segment, either the user - * provided one or the default one. */ -SENC2D_API res_T -senc2d_descriptor_get_global_segment_global_id - (const struct senc2d_descriptor* descriptor, - const unsigned iseg, - unsigned* gid); - /* Returns the number of vertices that are frontier vertices: * - that have arity 1 (single segment using the vertex) * - that connect 2 different media */ @@ -399,7 +416,8 @@ senc2d_enclosure_get_header (const struct senc2d_enclosure* enclosure, struct senc2d_enclosure_header* header); -/* Returns the iseg_th segment of an enclosure. */ +/* Returns the iseg_th segment of an enclosure. + * Indices are local to the enclosure. */ SENC2D_API res_T senc2d_enclosure_get_segment (const struct senc2d_enclosure* enclosure, @@ -413,20 +431,14 @@ senc2d_enclosure_get_vertex const unsigned ivert, double coord[2]); -/* Returns the front and back side media ids of the iseg_th segment of an - * enclosure. */ -SENC2D_API res_T -senc2d_enclosure_get_segment_media - (const struct senc2d_enclosure* enclosure, - const unsigned iseg, - unsigned medium[2]); - -/* Returns the global id of the iseg_th segment of an enclosure. */ +/* Returns the global id of the iseg_th segment of an enclosure + * and the involved side. */ SENC2D_API res_T senc2d_enclosure_get_segment_global_id (const struct senc2d_enclosure* enclosure, const unsigned iseg, - unsigned* gid); + unsigned* gid, + enum senc2d_side* side); /* Returns the id of the imed_th medium of an enclosure. */ SENC2D_API res_T diff --git a/src/senc2d_descriptor.c b/src/senc2d_descriptor.c @@ -272,21 +272,6 @@ senc2d_descriptor_get_global_segment_enclosures } res_T -senc2d_descriptor_get_global_segment_global_id - (const struct senc2d_descriptor* desc, - const unsigned iseg, - unsigned* gid) -{ - const struct segment_in* seg; - if(!gid || !desc - || iseg >= darray_segment_in_size_get(&desc->scene->segments_in)) - return RES_BAD_ARG; - seg = darray_segment_in_cdata_get(&desc->scene->segments_in) + iseg; - *gid = seg->global_id; - return RES_OK; -} - -res_T senc2d_descriptor_get_frontier_vertices_count (const struct senc2d_descriptor* desc, unsigned* count) diff --git a/src/senc2d_enclosure.c b/src/senc2d_enclosure.c @@ -81,17 +81,17 @@ senc2d_enclosure_get_segment const unsigned iseg, unsigned indices[2]) { - const struct segment_in* segment; + const struct side_enc* side; int i; if(!enclosure || !indices || iseg >= enclosure->data->header.segment_count) return RES_BAD_ARG; - ASSERT(darray_segment_in_size_get(&enclosure->data->sides) + ASSERT(darray_sides_enc_size_get(&enclosure->data->sides) == enclosure->data->header.segment_count); - segment = darray_segment_in_cdata_get(&enclosure->data->sides) + iseg; + side = darray_sides_enc_cdata_get(&enclosure->data->sides) + iseg; FOR_EACH(i, 0, 2) { - ASSERT(segment->vertice_id[i] < UINT_MAX); - indices[i] = (unsigned)segment->vertice_id[i]; /* Back to API type */ + ASSERT(side->vertice_id[i] < UINT_MAX); + indices[i] = (unsigned)side->vertice_id[i]; /* Back to API type */ } return RES_OK; } @@ -118,42 +118,21 @@ senc2d_enclosure_get_vertex } res_T -senc2d_enclosure_get_segment_media - (const struct senc2d_enclosure* enclosure, - const unsigned iseg, - unsigned medium[2]) -{ - const struct segment_in* segment; - int i; - if(!enclosure || !medium - || iseg >= enclosure->data->header.segment_count) - return RES_BAD_ARG; - ASSERT(darray_segment_in_size_get(&enclosure->data->sides) - == enclosure->data->header.segment_count); - segment = darray_segment_in_cdata_get(&enclosure->data->sides) + iseg; - FOR_EACH(i, 0, 2) { -#if (UINT_MAX < MEDIUM_MAX__) - ASSERT(segment->medium[i] < UINT_MAX); -#endif - medium[i] = (unsigned)segment->medium[i]; /* Back to API type */ - } - return RES_OK; -} - -res_T senc2d_enclosure_get_segment_global_id (const struct senc2d_enclosure* enclosure, const unsigned iseg, - unsigned* gid) + unsigned* gid, + enum senc2d_side* sde) { - const struct segment_in* segment; - if(!enclosure || !gid + const struct side_enc* side; + if(!enclosure || !gid || !sde || iseg >= enclosure->data->header.segment_count) return RES_BAD_ARG; - ASSERT(darray_segment_in_size_get(&enclosure->data->sides) + ASSERT(darray_sides_enc_size_get(&enclosure->data->sides) == enclosure->data->header.segment_count); - segment = darray_segment_in_cdata_get(&enclosure->data->sides) + iseg; - *gid = segment->global_id; + side = darray_sides_enc_cdata_get(&enclosure->data->sides) + iseg; + *gid = (unsigned)SEGSIDE_2_SEG(side->side_id); + *sde = SEGSIDE_2_SIDE(side->side_id); return RES_OK; } diff --git a/src/senc2d_enclosure_data.h b/src/senc2d_enclosure_data.h @@ -21,6 +21,21 @@ #include <rsys/hash_table.h> #include <rsys/dynamic_array.h> +#include "senc2d.h" +#include "senc2d_scene_c.h" +#include "senc2d_internal_types.h" + +#include <limits.h> + +struct side_enc { + vrtx_id_t vertice_id[2]; + side_id_t side_id; +}; + +#define DARRAY_NAME sides_enc +#define DARRAY_DATA struct side_enc +#include <rsys/dynamic_array.h> + /* unsigned char array with init to zero */ static FINLINE void zero_init_uchar @@ -32,12 +47,6 @@ zero_init_uchar #define DARRAY_FUNCTOR_INIT zero_init_uchar #include <rsys/dynamic_array_uchar.h> -#include "senc2d.h" -#include "senc2d_scene_c.h" -#include "senc2d_internal_types.h" - -#include <limits.h> - static void init_header(struct senc2d_enclosure_header* header) { @@ -112,7 +121,7 @@ error: struct enclosure_data { struct senc2d_enclosure_header header; /* Same segment can appear twice if both sides */ - struct darray_segment_in sides; + struct darray_sides_enc sides; /* Index of vertices in scene's unique vertices */ struct darray_vrtx_id vertices; /* List of the enclosed media */ @@ -137,7 +146,7 @@ enclosure_data_init(struct mem_allocator* alloc, struct enclosure_data* enc) { enc->side_range.first = SIDE_NULL__; enc->side_range.last = 0; enc->side_count = 0; - darray_segment_in_init(alloc, &enc->sides); + darray_sides_enc_init(alloc, &enc->sides); darray_vrtx_id_init(alloc, &enc->vertices); darray_uchar_init(alloc, &enc->tmp_enclosed_media); darray_media_init(alloc, &enc->enclosed_media); @@ -155,7 +164,7 @@ enclosure_data_copy dst->first_component = src->first_component; dst->side_range = src->side_range; dst->side_count = src->side_count; - OK(darray_segment_in_copy(&dst->sides, &src->sides)); + OK(darray_sides_enc_copy(&dst->sides, &src->sides)); OK(darray_vrtx_id_copy(&dst->vertices, &src->vertices)); OK(darray_uchar_copy(&dst->tmp_enclosed_media, &src->tmp_enclosed_media)); OK(darray_media_copy(&dst->enclosed_media, &src->enclosed_media)); @@ -166,7 +175,7 @@ error: static FINLINE void enclosure_data_release(struct enclosure_data* n) { ASSERT(n); - darray_segment_in_release(&n->sides); + darray_sides_enc_release(&n->sides); darray_vrtx_id_release(&n->vertices); darray_uchar_release(&n->tmp_enclosed_media); darray_media_release(&n->enclosed_media); @@ -184,7 +193,7 @@ enclosure_data_copy_and_release dst->first_component = src->first_component; dst->side_range = src->side_range; dst->side_count = src->side_count; - OK(darray_segment_in_copy_and_release(&dst->sides, &src->sides)); + OK(darray_sides_enc_copy_and_release(&dst->sides, &src->sides)); OK(darray_vrtx_id_copy_and_release(&dst->vertices, &src->vertices)); OK(darray_uchar_copy_and_release(&dst->tmp_enclosed_media, &src->tmp_enclosed_media)); diff --git a/src/senc2d_internal_types.h b/src/senc2d_internal_types.h @@ -95,12 +95,6 @@ enum side_flag { FLAG_BACK = BIT(1) }; -/* This one is used as an index to arrays */ -enum side_id { - SIDE_FRONT = 0, - SIDE_BACK = 1 -}; - /* Utility macros */ static FINLINE seg_id_t SEGSIDE_2_SEG(side_id_t s) { @@ -113,9 +107,9 @@ SEGSIDE_IS_FRONT(side_id_t s) { return (s & 1) == 0; } -static FINLINE enum side_id +static FINLINE enum senc2d_side SEGSIDE_2_SIDE(side_id_t s) { - return (s & 1) ? SIDE_BACK : SIDE_FRONT; + return (s & 1) ? SENC2D_BACK : SENC2D_FRONT; } static FINLINE enum side_flag @@ -123,22 +117,23 @@ SEGSIDE_2_SIDEFLAG(side_id_t s) { return (s & 1) ? FLAG_BACK : FLAG_FRONT; } -static FINLINE enum side_flag +static FINLINE unsigned char SIDE_CANCELED_FLAG(enum side_flag f) { - return ((unsigned char)f) << 4; + ASSERT((((unsigned)f) << 4) <= UCHAR_MAX); + return (unsigned char)(f << 4); } static FINLINE side_id_t -SEGIDxSIDE_2_SEGSIDE(seg_id_t s, enum side_id i) { - ASSERT((((size_t)s << 1) | (i == SIDE_BACK)) < SIDE_MAX__); - ASSERT(i == SIDE_FRONT || i == SIDE_BACK); - return (side_id_t)((s << 1) | (i == SIDE_BACK)); +SEGIDxSIDE_2_SEGSIDE(seg_id_t s, enum senc2d_side i) { + ASSERT((((size_t)s << 1) | (i == SENC2D_BACK)) < SIDE_MAX__); + ASSERT(i == SENC2D_FRONT || i == SENC2D_BACK); + return (side_id_t)((s << 1) | (i == SENC2D_BACK)); } static FINLINE side_id_t SEGSIDE_OPPOSITE(side_id_t s) { return SEGIDxSIDE_2_SEGSIDE(SEGSIDE_2_SEG(s), - SEGSIDE_IS_FRONT(s) ? SIDE_BACK : SIDE_FRONT); + SEGSIDE_IS_FRONT(s) ? SENC2D_BACK : SENC2D_FRONT); } #endif /* SENC2D_INTERNAL_TYPES_H */ diff --git a/src/senc2d_scene.c b/src/senc2d_scene.c @@ -130,11 +130,12 @@ res_T senc2d_scene_add_geometry (struct senc2d_scene* scn, const unsigned nsegs, - void(*indices)(const unsigned iseg, unsigned ids[2], void* ctx), - void(*media)(const unsigned iseg, unsigned med[2], void* ctx), - void(*global_id)(const unsigned iseg, unsigned* gid, void* ctx), + void(*indices)(const unsigned, unsigned* ids, void*), + void(*media)(const unsigned, unsigned* med, void*), const unsigned nverts, - void(*position)(const unsigned ivert, double pos[2], void* ctx), + void(*position)(const unsigned, double* pos, void*), + res_T(*add_segment)(const unsigned, const unsigned, void*), + res_T(*merge_segment)(const unsigned, const unsigned, const int, void*), void* ctx) { struct darray_vrtx_id unique_vertice_ids; @@ -199,11 +200,6 @@ senc2d_scene_add_geometry struct segment_in tmp, *range_adjust_ptr = NULL; seg_id_t* p_seg; char reversed; - if(global_id) { - global_id(i, &tmp.global_id, ctx); - } else { - tmp.global_id = (unsigned)(scn->nsegs + i); - } indices(i, ind, ctx); /* API: indices needs an unsigned */ FOR_EACH(j, 0, 2) { if(ind[j] >= nverts) { @@ -218,7 +214,7 @@ senc2d_scene_add_geometry const union double2* positions = darray_position_cdata_get(&scn->vertices); log_err(scn->dev, "%s: segment %lu is degenerate.\n", - FUNC_NAME, (unsigned long)tmp.global_id); + FUNC_NAME, (unsigned long)scn->nsegs + i); log_err(scn->dev, " (%g %g) (%g %g)\n", SPLIT2(positions[seg[i].vertice_id[0]].vec), SPLIT2(positions[seg[i].vertice_id[1]].vec)); @@ -255,11 +251,10 @@ senc2d_scene_add_geometry = darray_position_cdata_get(&scn->vertices); log_err(scn->dev, "%s: segment %lu is a duplicate" " of segment %lu with incoherent media.\n", - FUNC_NAME, (unsigned long)tmp.global_id, - (unsigned long)seg[*p_seg].global_id); + FUNC_NAME, (unsigned long)(scn->nsegs + i), (unsigned long)*p_seg); log_err(scn->dev, "Segment %lu:\n (%g %g ) (%g %g)\n", - (unsigned long)seg[*p_seg].global_id, + (unsigned long)*p_seg, SPLIT2(positions[seg[*p_seg].vertice_id[0]].vec), SPLIT2(positions[seg[*p_seg].vertice_id[1]].vec)); log_err(scn->dev, "Media: (%lu, %lu) VS (%lu, %lu)\n", @@ -271,9 +266,6 @@ senc2d_scene_add_geometry goto error; } else { /* Legit duplicate */ - log_warn(scn->dev, "%s: segment %lu is a duplicate of segment %lu.\n", - FUNC_NAME, (unsigned long)tmp.global_id, - (unsigned long)seg[*p_seg].global_id); range_adjust_ptr = darray_segment_in_data_get(&scn->segments_in) + *p_seg; /* Replace possible undefined media */ FOR_EACH(j, 0, 2) { @@ -283,6 +275,14 @@ senc2d_scene_add_geometry scn->sides_with_defined_medium_count++; } } + if (merge_segment) { + OK(merge_segment((unsigned)*p_seg, i, same, ctx)); + } + else { + log_warn(scn->dev, + "%s: segment %lu is a duplicate of segment %lu.\n", + FUNC_NAME, (unsigned long)(scn->nsegs + i), (unsigned long)*p_seg); + } } } else { /* New segment */ @@ -295,6 +295,9 @@ senc2d_scene_add_geometry if(tmp.medium[j] != SENC2D_UNDEFINED_MEDIUM) scn->sides_with_defined_medium_count++; } + if (add_segment) { + OK(add_segment((unsigned)u, i, ctx)); + } ++actual_nusegs; } if(range_adjust_ptr) { diff --git a/src/senc2d_scene_analyze.c b/src/senc2d_scene_analyze.c @@ -116,7 +116,7 @@ self_hit_filter const struct darray_segment_comp* segments_comp = filter_data; const component_id_t* origin_component = ray_data; const struct segment_comp* hit_seg_comp; - enum side_id hit_side; + enum senc2d_side hit_side; component_id_t hit_component; (void)ray_org; (void)ray_dir; @@ -124,7 +124,7 @@ self_hit_filter ASSERT(hit->prim.prim_id < darray_segment_comp_size_get(segments_comp)); hit_seg_comp = darray_segment_comp_cdata_get(segments_comp) + hit->prim.prim_id; - hit_side = (hit->normal[1] > 0) ? SIDE_FRONT : SIDE_BACK; + hit_side = (hit->normal[1] > 0) ? SENC2D_FRONT : SENC2D_BACK; hit_component = hit_seg_comp->component[hit_side]; /* If self hit, distance should be small */ @@ -257,7 +257,7 @@ extract_connex_components #ifndef NDEBUG { seg_id_t sid = SEGSIDE_2_SEG(start_side_id); - enum side_id s = SEGSIDE_2_SIDE(start_side_id); + enum senc2d_side s = SEGSIDE_2_SIDE(start_side_id); medium_id_t side_med = darray_segment_in_data_get(&desc->scene->segments_in)[sid].medium[s]; ASSERT(side_med == m); @@ -323,7 +323,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 |= (unsigned char)crt_side_flag; + *seg_used = *seg_used |(unsigned char)crt_side_flag; } /* Store neighbours' sides in a stack */ @@ -336,7 +336,7 @@ extract_connex_components medium_id_t nbour_med_idx = (neighbour->medium == SENC2D_UNDEFINED_MEDIUM) ? scn->next_medium_idx : neighbour->medium; if(neighbour->medium < m - || (*nbour_used & (unsigned char)SIDE_CANCELED_FLAG(nbour_side_id))) + || (*nbour_used & SIDE_CANCELED_FLAG(nbour_side_id))) { /* 1) Not the same medium. * Neighbour's medium id is less than current medium: the whole @@ -361,7 +361,7 @@ extract_connex_components ASSERT(*used & (unsigned char)used_side_flag); /* Set the used flag for sides in cancelled component as leading * to further cancellations */ - *used |= (unsigned char)SIDE_CANCELED_FLAG(used_side_flag); + *used |= SIDE_CANCELED_FLAG(used_side_flag); } goto canceled; @@ -411,7 +411,7 @@ canceled: FOR_EACH(ii, 0, sz) { const side_id_t s = darray_side_id_cdata_get(&current_component)[ii]; seg_id_t segid = SEGSIDE_2_SEG(s); - enum side_id sid = SEGSIDE_2_SIDE(s); + enum senc2d_side sid = SEGSIDE_2_SIDE(s); component_id_t* cmp = segments_comp[segid].component; ASSERT(cmp[sid] == COMPONENT_NULL__); ASSERT(segsides[s].medium >= m); @@ -426,7 +426,7 @@ canceled: const side_id_t side_id = darray_side_id_cdata_get(&ids_of_sides_around_max_y_vertex)[ii]; const seg_id_t seg_id = SEGSIDE_2_SEG(side_id); - enum side_id s = SEGSIDE_2_SIDE(side_id); + enum senc2d_side s = SEGSIDE_2_SIDE(side_id); const struct segment_in* seg_in = darray_segment_in_cdata_get(&scn->segments_in) + seg_id; struct segment_comp* seg_comp = segments_comp + seg_id; @@ -438,7 +438,7 @@ canceled: * regardless of numeric accuracy, we need to prevent them to * contribute (remember than x + y - y == x can be false). */ ASSERT(seg_comp->component[s] == cc->cc_id); (void)s; - if(seg_comp->component[SIDE_FRONT] == seg_comp->component[SIDE_BACK]) + if(seg_comp->component[SENC2D_FRONT] == seg_comp->component[SENC2D_BACK]) continue; d2_sub(edge, vertices[seg_in->vertice_id[1]].vec, @@ -542,8 +542,8 @@ canceled: FOR_EACH(s_, 0, scn->nusegs) { struct segment_comp* seg_comp = darray_segment_comp_data_get(segments_comp_array) + s_; - ASSERT(seg_comp->component[SIDE_FRONT] != COMPONENT_NULL__); - ASSERT(seg_comp->component[SIDE_BACK] != COMPONENT_NULL__); + ASSERT(seg_comp->component[SENC2D_FRONT] != COMPONENT_NULL__); + ASSERT(seg_comp->component[SENC2D_BACK] != COMPONENT_NULL__); } FOR_EACH(c, 0, ATOMIC_GET(component_count)) { struct cc_descriptor** components = @@ -632,13 +632,13 @@ group_connex_components const seg_id_t hit_seg_id = (seg_id_t)hit.prim.prim_id; const struct segment_comp* hit_seg_comp = darray_segment_comp_cdata_get(segments_comp) + hit_seg_id; - enum side_id hit_side = + enum senc2d_side hit_side = ((hit.normal[1] < 0) /* Facing geometrical normal of hit */ == ((desc->scene->convention & SENC2D_CONVENTION_NORMAL_FRONT) != 0)) /* Warning: following Embree 2 convention for geometrical normals, * the Star2D hit normal is left-handed while star-enclosure2D uses * right-handed convention */ - ? SIDE_BACK : SIDE_FRONT; + ? SENC2D_BACK : SENC2D_FRONT; ASSERT(hit.normal[1] != 0); ASSERT(hit_seg_id < desc->scene->nusegs); @@ -859,11 +859,11 @@ collect_and_link_neighbours const seg_id_t ccw_id = ccw_neighbour->seg_id; /* Facing sides of segments */ const int front = ((scn->convention & SENC2D_CONVENTION_NORMAL_FRONT) != 0); - const enum side_id crt_side - = current->normal_toward_next_neighbour == front ? SIDE_FRONT : SIDE_BACK; - const enum side_id ccw_side + const enum senc2d_side crt_side + = current->normal_toward_next_neighbour == front ? SENC2D_FRONT : SENC2D_BACK; + const enum senc2d_side ccw_side = ccw_neighbour->normal_toward_next_neighbour == front ? - SIDE_BACK : SIDE_FRONT; + SENC2D_BACK : SENC2D_FRONT; /* Index of sides in segsides */ const side_id_t crt_side_idx = SEGIDxSIDE_2_SEGSIDE(crt_id, crt_side); const side_id_t ccw_side_idx = SEGIDxSIDE_2_SEGSIDE(ccw_id, ccw_side); @@ -881,10 +881,8 @@ collect_and_link_neighbours prev_id = previous->seg_id; log_err(scn->dev, "%s: found 2 overlying segments (%lu & %lu).\n", FUNC_NAME, - (unsigned long)segments_in[crt_id].global_id, - (unsigned long)segments_in[prev_id].global_id); + (unsigned long)crt_id, (unsigned long)prev_id); - /* TODO */ *res = RES_BAD_OP; return; } @@ -914,7 +912,7 @@ collect_and_link_neighbours { log_warn(scn->dev, "%s: found possible frontier involving segment %lu.\n", - FUNC_NAME, (unsigned long)segments_in[crt_id].global_id); + FUNC_NAME, (unsigned long)crt_id); darray_vrtx_id_push_back(frontiers, &v); } } @@ -975,16 +973,16 @@ build_result #pragma omp for for(sg = 0; sg < (int64_t)scn->nusegs; sg++) { seg_id_t s = (seg_id_t)sg; - const component_id_t cf_id = segments_comp[s].component[SIDE_FRONT]; - const component_id_t cb_id = segments_comp[s].component[SIDE_BACK]; + const component_id_t cf_id = segments_comp[s].component[SENC2D_FRONT]; + const component_id_t cb_id = segments_comp[s].component[SENC2D_BACK]; const struct cc_descriptor* cf = cc_descriptors[cf_id]; const struct cc_descriptor* cb = cc_descriptors[cb_id]; const enclosure_id_t ef_id = cf->enclosure_id; const enclosure_id_t eb_id = cb->enclosure_id; - ASSERT(segments_enc[s].enclosure[SIDE_FRONT] == ENCLOSURE_NULL__); - segments_enc[s].enclosure[SIDE_FRONT] = ef_id; - ASSERT(segments_enc[s].enclosure[SIDE_BACK] == ENCLOSURE_NULL__); - segments_enc[s].enclosure[SIDE_BACK] = eb_id; + ASSERT(segments_enc[s].enclosure[SENC2D_FRONT] == ENCLOSURE_NULL__); + segments_enc[s].enclosure[SENC2D_FRONT] = ef_id; + ASSERT(segments_enc[s].enclosure[SENC2D_BACK] == ENCLOSURE_NULL__); + segments_enc[s].enclosure[SENC2D_BACK] = eb_id; } /* Implicit barrier here */ @@ -1045,7 +1043,7 @@ build_result if(*res != RES_OK) continue; /* Build side and vertex lists. */ - OK2(darray_segment_in_resize(&enc->sides, enc->side_count)); + OK2(darray_sides_enc_resize(&enc->sides, enc->side_count)); /* Size is just a int */ OK2(darray_vrtx_id_reserve(&enc->vertices, enc->side_count + 1)); /* New vertex numbering scheme local to the enclosure */ @@ -1058,11 +1056,11 @@ build_result s++) { const struct segment_in* seg_in = segments_in + s; - struct segment_in* seg; + struct side_enc* side_enc; unsigned vertice_id[2]; int i; - if(segments_enc[s].enclosure[SIDE_FRONT] != e - && segments_enc[s].enclosure[SIDE_BACK] != e) + if(segments_enc[s].enclosure[SENC2D_FRONT] != e + && segments_enc[s].enclosure[SENC2D_BACK] != e) continue; ++enc->header.unique_segment_count; @@ -1082,41 +1080,33 @@ build_result ++enc->header.vertices_count; } } - ASSERT(segments_enc[s].enclosure[SIDE_FRONT] == e - || segments_enc[s].enclosure[SIDE_BACK] == e); - if(segments_enc[s].enclosure[SIDE_FRONT] == e) { + ASSERT(segments_enc[s].enclosure[SENC2D_FRONT] == e + || segments_enc[s].enclosure[SENC2D_BACK] == e); + if(segments_enc[s].enclosure[SENC2D_FRONT] == e) { /* Front side of the original segment is member of the enclosure */ int input_normal_in = normals_front; int revert_segment = (input_normal_in != output_normal_in); ++enc->header.segment_count; - seg = darray_segment_in_data_get(&enc->sides) + fst_idx++; + side_enc = darray_sides_enc_data_get(&enc->sides) + fst_idx++; FOR_EACH(i, 0, 2) { int ii = revert_segment ? 1 - i : i; - seg->medium[i] = seg_in->medium[ii]; - } - seg->global_id = seg_in->global_id; - FOR_EACH(i, 0, 2) { - int ii = revert_segment ? 1 - i : i; - seg->vertice_id[i] = vertice_id[ii]; + side_enc->vertice_id[i] = vertice_id[ii]; } + side_enc->side_id = SEGIDxSIDE_2_SEGSIDE(s, SENC2D_FRONT); } - if(segments_enc[s].enclosure[SIDE_BACK] == e) { + if(segments_enc[s].enclosure[SENC2D_BACK] == e) { /* Back side of the original segment is member of the enclosure */ int input_normal_in = normals_back; int revert_segment = (input_normal_in != output_normal_in); ++enc->header.segment_count; /* If both sides are in the enclosure, put the second side at the end */ - seg = darray_segment_in_data_get(&enc->sides) + - ((segments_enc[s].enclosure[SIDE_FRONT] == e) ? --sgd_idx : fst_idx++); - FOR_EACH(i, 0, 2) { - int ii = revert_segment ? 1 - i : i; - seg->medium[i] = seg_in->medium[ii]; - } - seg->global_id = seg_in->global_id; + side_enc = darray_sides_enc_data_get(&enc->sides) + + ((segments_enc[s].enclosure[SENC2D_FRONT] == e) ? --sgd_idx : fst_idx++); FOR_EACH(i, 0, 2) { int ii = revert_segment ? 1 - i : i; - seg->vertice_id[i] = vertice_id[ii]; + side_enc->vertice_id[i] = vertice_id[ii]; } + side_enc->side_id = SEGIDxSIDE_2_SEGSIDE(s, SENC2D_BACK); } if(fst_idx == sgd_idx) break; } @@ -1163,8 +1153,6 @@ senc2d_scene_analyze /* Atomic counters to share beetwen threads */ ATOMIC component_count = 0; ATOMIC next_enclosure_id = 1; - char dump[64]; - struct time t0, t1; res_T res = RES_OK; res_T res2 = RES_OK; @@ -1179,11 +1167,6 @@ senc2d_scene_analyze if(!scn->nusegs) goto exit; - #pragma omp single nowait - if(scn->dev->verbose) { - time_current(&t0); - } - darray_segment_tmp_init(scn->dev->allocator, &segments_tmp); segments_tmp_initialized = 1; darray_vrtx_id_init(scn->dev->allocator, &frontiers); @@ -1265,15 +1248,6 @@ senc2d_scene_analyze neighbourhood_by_vertex_initialized = 0; } /* No barrier here */ - #pragma omp single nowait - if(scn->dev->verbose) { - time_sub(&t1, time_current(&t1), &t0); - time_current(&t0); - time_dump(&t1, TIME_MSEC | TIME_SEC | TIME_MIN, NULL, dump, sizeof(dump)); - log_info(scn->dev, - "senc2d_scene_analyze: collect_and_link_neighbours step in %s\n", dump); - } - /* Step 2: extract segment connex components */ extract_connex_components(desc, segsides, &connex_components, &segments_tmp, &segments_comp, &s2d_view, &component_count, &res); @@ -1301,15 +1275,6 @@ senc2d_scene_analyze segments_tmp_initialized = 0; } /* No barrier here */ - #pragma omp single nowait - if(scn->dev->verbose) { - time_sub(&t1, time_current(&t1), &t0); - time_current(&t0); - time_dump(&t1, TIME_MSEC | TIME_SEC | TIME_MIN, NULL, dump, sizeof(dump)); - log_info(scn->dev, - "senc2d_scene_analyze: extract_connex_components step in %s\n", dump); - } - /* Step 3: group components */ group_connex_components(desc, segsides, &segments_comp, &connex_components, s2d_view, &next_enclosure_id, &res); @@ -1333,15 +1298,6 @@ senc2d_scene_analyze s2d_view = NULL; } /* No barrier here */ - #pragma omp single nowait - if(scn->dev->verbose) { - time_sub(&t1, time_current(&t1), &t0); - time_current(&t0); - time_dump(&t1, TIME_MSEC | TIME_SEC | TIME_MIN, NULL, dump, sizeof(dump)); - log_info(scn->dev, - "senc2d_scene_analyze: group_connex_components step in %s\n", dump); - } - /* Step 4: Build result */ build_result(desc, &connex_components, &segments_comp, &frontiers, &res); @@ -1375,15 +1331,6 @@ senc2d_scene_analyze } } /* No barrier here */ - #pragma omp single nowait - if(scn->dev->verbose) { - time_sub(&t1, time_current(&t1), &t0); - time_current(&t0); - time_dump(&t1, TIME_MSEC | TIME_SEC | TIME_MIN, NULL, dump, sizeof(dump)); - log_info(scn->dev, - "senc2d_scene_analyze: build_result step in %s\n", dump); - } - error_: ; } /* Implicit barrier here */ diff --git a/src/senc2d_scene_c.h b/src/senc2d_scene_c.h @@ -48,8 +48,6 @@ struct segment_in { vrtx_id_t vertice_id[2]; /* Ids of this segment's media */ medium_id_t medium[2]; - /* Segment index in user world (that is regardless of deduplication). */ - unsigned global_id; }; static FINLINE void @@ -59,7 +57,6 @@ segment_in_init(struct mem_allocator* alloc, struct segment_in* seg) { ASSERT(seg); FOR_EACH(i, 0, 2) seg->vertice_id[i] = VRTX_NULL__; FOR_EACH(i, 0, 2) seg->medium[i] = SENC2D_UNDEFINED_MEDIUM; - seg->global_id = 0; } #define DARRAY_NAME segment_in diff --git a/src/test_senc2d_descriptor.c b/src/test_senc2d_descriptor.c @@ -53,8 +53,8 @@ main(int argc, char** argv) ctx.front_media = medium0; ctx.back_media = medium1; - CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL, - nvertices, get_position, &ctx) == RES_OK); + CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); CHK(senc2d_scene_analyze(scn, &desc) == RES_OK); @@ -167,22 +167,11 @@ main(int argc, char** argv) == RES_OK); CHK(enclosures[0] == 0 && enclosures[1] == 1); - CHK(senc2d_descriptor_get_global_segment_global_id(NULL, 0, indices) - == RES_BAD_ARG); - CHK(senc2d_descriptor_get_global_segment_global_id(NULL, nvertices, indices) - == RES_BAD_ARG); - CHK(senc2d_descriptor_get_global_segment_global_id(desc, 0, NULL) - == RES_BAD_ARG); - CHK(senc2d_descriptor_get_global_segment_global_id(desc, 0, indices) - == RES_OK); - /* No duplicates and no custom id: user id is unique vertex id */ - CHK(indices[0] == 0); - /* Add valid duplicate geometry */ CHK(senc2d_descriptor_ref_put(desc) == RES_OK); desc = NULL; CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, - NULL, nvertices, get_position, &ctx) == RES_OK); + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); CHK(senc2d_scene_analyze(scn, &desc) == RES_OK); CHK(senc2d_descriptor_get_global_vertices_count(desc, &count) == RES_OK); @@ -198,8 +187,8 @@ main(int argc, char** argv) desc = NULL; ctx.front_media = medium1; ctx.back_media = medium0; - CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL, - nvertices, get_position, &ctx) == RES_BAD_ARG); + CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, NULL, NULL, &ctx) == RES_BAD_ARG); CHK(senc2d_scene_ref_put(scn) == RES_OK); if(desc) CHK(senc2d_descriptor_ref_put(desc) == RES_OK); @@ -212,7 +201,7 @@ main(int argc, char** argv) == RES_OK); CHK(senc2d_scene_add_geometry(scn, nsegments - 1, get_indices, get_media, - NULL, nvertices, get_position, &ctx) == RES_OK); + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); CHK(senc2d_scene_analyze(scn, &desc) == RES_OK); diff --git a/src/test_senc2d_enclosure.c b/src/test_senc2d_enclosure.c @@ -22,7 +22,7 @@ #include <star/s2d.h> static void -test(int convention) +test(const int convention) { struct mem_allocator allocator; struct senc2d_descriptor* desc = NULL; @@ -36,13 +36,15 @@ test(int convention) struct s2d_shape* s2d_shp = NULL; struct s2d_vertex_data s2d_attribs; unsigned indices[2][2]; - unsigned medium, media[2]; + unsigned medium; unsigned gid; + enum senc2d_side side; double vrtx[2]; struct context ctx; - unsigned i, n, t, ecount; + unsigned i, n, s, ecount; int conv; - int is_front, is_in; + const int conv_front = (convention & SENC2D_CONVENTION_NORMAL_FRONT) != 0; + const int conv_in = (convention & SENC2D_CONVENTION_NORMAL_INSIDE) != 0; CHK(mem_init_proxy_allocator(&allocator, &mem_default_allocator) == RES_OK); CHK(senc2d_device_create(NULL, &allocator, SENC2D_NTHREADS_DEFAULT, 1, &dev) @@ -52,8 +54,6 @@ test(int convention) CHK(senc2d_scene_get_convention(scn, &conv) == RES_OK); CHK(conv == convention); - is_front = (conv & SENC2D_CONVENTION_NORMAL_FRONT) != 0; - is_in = (conv & SENC2D_CONVENTION_NORMAL_INSIDE) != 0; s2d_attribs.type = S2D_FLOAT2; s2d_attribs.usage = S2D_POSITION; @@ -76,7 +76,7 @@ test(int convention) ctx.back_media = medium1; CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, - NULL, nvertices, get_position, &ctx) == RES_OK); + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); CHK(senc2d_scene_analyze(scn, &desc) == RES_OK); @@ -112,32 +112,38 @@ test(int convention) CHK(senc2d_enclosure_get_vertex(NULL, nvertices, NULL) == RES_BAD_ARG); CHK(senc2d_enclosure_get_vertex(enclosure, 0, vrtx) == RES_OK); - CHK(senc2d_enclosure_get_segment_media(NULL, 0, media) == RES_BAD_ARG); - CHK(senc2d_enclosure_get_segment_media(enclosure, nsegments, media) + CHK(senc2d_enclosure_get_segment_global_id(NULL, 0, &gid, NULL) == RES_BAD_ARG); - CHK(senc2d_enclosure_get_segment_media(enclosure, 0, NULL) == RES_BAD_ARG); - CHK(senc2d_enclosure_get_segment_media(NULL, nsegments, media) + CHK(senc2d_enclosure_get_segment_global_id(enclosure, nsegments, &gid, NULL) == RES_BAD_ARG); - CHK(senc2d_enclosure_get_segment_media(NULL, 0, NULL) == RES_BAD_ARG); - CHK(senc2d_enclosure_get_segment_media(enclosure, nsegments, NULL) + CHK(senc2d_enclosure_get_segment_global_id(enclosure, 0, NULL, NULL) == RES_BAD_ARG); - CHK(senc2d_enclosure_get_segment_media(NULL, nsegments, NULL) + CHK(senc2d_enclosure_get_segment_global_id(NULL, nsegments, &gid, NULL) == RES_BAD_ARG); - CHK(senc2d_enclosure_get_segment_media(enclosure, 0, media) == RES_OK); - - CHK(senc2d_enclosure_get_segment_global_id(NULL, 0, &gid) == RES_BAD_ARG); - CHK(senc2d_enclosure_get_segment_global_id(enclosure, nsegments, &gid) + CHK(senc2d_enclosure_get_segment_global_id(NULL, 0, NULL, NULL) + == RES_BAD_ARG); + CHK(senc2d_enclosure_get_segment_global_id(enclosure, nsegments, NULL, NULL) + == RES_BAD_ARG); + CHK(senc2d_enclosure_get_segment_global_id(NULL, nsegments, NULL, NULL) + == RES_BAD_ARG); + CHK(senc2d_enclosure_get_segment_global_id(enclosure, 0, &gid, NULL) + == RES_BAD_ARG); + CHK(senc2d_enclosure_get_segment_global_id(NULL, 0, &gid, &side) == RES_BAD_ARG); - CHK(senc2d_enclosure_get_segment_global_id(enclosure, 0, NULL) + CHK(senc2d_enclosure_get_segment_global_id(enclosure, nsegments, &gid, &side) == RES_BAD_ARG); - CHK(senc2d_enclosure_get_segment_global_id(NULL, nsegments, &gid) + CHK(senc2d_enclosure_get_segment_global_id(enclosure, 0, NULL, &side) == RES_BAD_ARG); - CHK(senc2d_enclosure_get_segment_global_id(NULL, 0, NULL) == RES_BAD_ARG); - CHK(senc2d_enclosure_get_segment_global_id(enclosure, nsegments, NULL) + CHK(senc2d_enclosure_get_segment_global_id(NULL, nsegments, &gid, &side) == RES_BAD_ARG); - CHK(senc2d_enclosure_get_segment_global_id(NULL, nsegments, NULL) + CHK(senc2d_enclosure_get_segment_global_id(NULL, 0, NULL, &side) == RES_BAD_ARG); - CHK(senc2d_enclosure_get_segment_global_id(enclosure, 0, &gid) == RES_OK); + CHK(senc2d_enclosure_get_segment_global_id(enclosure, nsegments, NULL, &side) + == RES_BAD_ARG); + CHK(senc2d_enclosure_get_segment_global_id(NULL, nsegments, NULL, &side) + == RES_BAD_ARG); + CHK(senc2d_enclosure_get_segment_global_id(enclosure, 0, &gid, &side) + == RES_OK); CHK(senc2d_enclosure_get_medium(NULL, 0, &medium) == RES_BAD_ARG); CHK(senc2d_enclosure_get_medium(enclosure, 2, &medium) == RES_BAD_ARG); @@ -164,14 +170,15 @@ test(int convention) /* Geometrical normals point outside the square in input segments: * if convention is front, front medium (0) is outside, * that is medium 0's enclosure is infinite */ - CHK(is_front == ((medium == 0) == header.is_infinite)); + CHK(conv_front == ((medium == 0) == header.is_infinite)); CHK(header.segment_count == nsegments); CHK(header.unique_segment_count == nsegments); CHK(header.vertices_count == nvertices); - FOR_EACH(t, 0, header.segment_count) { - CHK(senc2d_enclosure_get_segment_global_id(enclosure, t, &gid) == RES_OK); - CHK(gid == t); + FOR_EACH(s, 0, header.segment_count) { + CHK(senc2d_enclosure_get_segment_global_id(enclosure, s, &gid, &side) == RES_OK); + CHK(gid == s); + CHK(side == (medium == 0) ? SENC2D_FRONT : SENC2D_BACK); } CHK(senc2d_enclosure_ref_put(enclosure) == RES_OK); @@ -191,9 +198,9 @@ test(int convention) * of input segments for enclosure 0 iff convention is inside. * The opposite holds for enclosure 1. */ cmp_seg(n, enclosures[0], box_indices + 2 * n, box_vertices, &same, &reversed); - CHK((same && !reversed) == is_in); + CHK((same && !reversed) == conv_in); cmp_seg(n, enclosures[1], box_indices + 2 * n, box_vertices, &same, &reversed); - CHK(same && reversed == is_in); + CHK(same && reversed == conv_in); } FOR_EACH(i, 0, 2) CHK(senc2d_enclosure_ref_put(enclosures[i]) == RES_OK); @@ -204,9 +211,7 @@ test(int convention) /* Same 2D square, but with a hole (incomplete). * 1 single enclosure including both sides of segments */ - CHK(senc2d_scene_create(dev, - SENC2D_CONVENTION_NORMAL_BACK | SENC2D_CONVENTION_NORMAL_INSIDE, &scn) - == RES_OK); + CHK(senc2d_scene_create(dev, convention, &scn) == RES_OK); ctx.positions = box_vertices; ctx.indices = box_indices; @@ -216,7 +221,7 @@ test(int convention) ctx.back_media = medium1; CHK(senc2d_scene_add_geometry(scn, nsegments - 1, get_indices, get_media, - NULL, nvertices, get_position, &ctx) == RES_OK); + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); CHK(senc2d_scene_analyze(scn, &desc) == RES_OK); @@ -239,11 +244,12 @@ test(int convention) CHK(header.vertices_count == nvertices); CHK(header.is_infinite == 1); - FOR_EACH(t, 0, header.unique_segment_count) { + FOR_EACH(s, 0, header.segment_count) { + CHK(senc2d_enclosure_get_segment_global_id(enclosure, s, &gid, &side) == RES_OK); /* The first unique_segment_count segments of an enclosure * are unique segments */ - CHK(senc2d_enclosure_get_segment_global_id(enclosure, t, &gid) == RES_OK); - CHK(gid == t); + if (s < header.unique_segment_count) CHK(gid == s); + CHK(side == (s < header.unique_segment_count) ? SENC2D_FRONT : SENC2D_BACK); } FOR_EACH(n, 0, header.unique_segment_count) { diff --git a/src/test_senc2d_inconsistant_square.c b/src/test_senc2d_inconsistant_square.c @@ -44,7 +44,7 @@ static const unsigned inconsistant_medium_back[4] = { 0, 1, 1, 1 }; static void -test(int convention) +test(const int convention) { struct mem_allocator allocator; struct senc2d_descriptor* desc = NULL; @@ -78,7 +78,7 @@ test(int convention) ctx.back_media = inconsistant_medium_back; CHK(senc2d_scene_add_geometry(scn, inconsistant_box_nsegments, get_indices, - get_media, NULL, nvertices, get_position, &ctx) == RES_OK); + get_media, nvertices, get_position, NULL, NULL, &ctx) == RES_OK); CHK(senc2d_scene_analyze(scn, &desc) == RES_OK); @@ -93,8 +93,8 @@ test(int convention) FOR_EACH(e, 0, ecount) { unsigned medium, expected_external_medium, expected_medium; char name[128]; - int front_inside; - int expected_side; + enum senc2d_side side, expected_side; + unsigned gid; CHK(senc2d_descriptor_get_enclosure(desc, e, &enclosure) == RES_OK); CHK(senc2d_enclosure_get_header(enclosure, &header) == RES_OK); CHK(header.enclosed_media_count == 1); @@ -105,9 +105,6 @@ test(int convention) expected_medium = (header.is_infinite ? expected_external_medium : !expected_external_medium); CHK(medium == expected_medium); - /* Common media, for non input segments */ - front_inside = (conv_front == conv_in); - expected_side = front_inside ? 0 : 1; sprintf(name, "test_inconsistant_square_%s_%s_%u.obj", conv_front ? "front" : "back", conv_in ? "in" : "out", e); @@ -115,16 +112,17 @@ test(int convention) FOR_EACH(i, 0, header.segment_count) { int same, reversed, fst_reversed; - unsigned med[2]; - fst_reversed = (e == 0) == conv_in; - CHK(senc2d_enclosure_get_segment_media(enclosure, i, med) == RES_OK); - CHK(med[expected_side] == medium); + fst_reversed = ((e == 0) == conv_in); + expected_side = (inconsistant_medium_front[i] == expected_medium) + ? SENC2D_FRONT : SENC2D_BACK; cmp_seg(i, enclosure, inconsistant_box_indices + (2 * i), box_vertices, &same, &reversed); /* Should be made of the same segments */ CHK(same); CHK(i ? reversed != fst_reversed : reversed == fst_reversed); + CHK(senc2d_enclosure_get_segment_global_id(enclosure, i, &gid, &side) == RES_OK); + CHK(side == expected_side); } SENC2D(enclosure_ref_put(enclosure)); } diff --git a/src/test_senc2d_many_enclosures.c b/src/test_senc2d_many_enclosures.c @@ -106,7 +106,7 @@ main(int argc, char** argv) ctx.scale = k + 1; d2(ctx.offset, center_x, center_y); CHK(senc2d_scene_add_geometry(scn, circ_seg_count, get_ctx_indices, - get_ctx_media, NULL, circ_vrtx_count, get_ctx_position, &ctx) + get_ctx_media, circ_vrtx_count, get_ctx_position, NULL, NULL, &ctx) == RES_OK); } } diff --git a/src/test_senc2d_many_segments.c b/src/test_senc2d_many_segments.c @@ -98,7 +98,7 @@ main(int argc, char** argv) m1 = i; d2(ctx.offset, 0, i * 10); CHK(senc2d_scene_add_geometry(scn, circ_seg_count, get_ctx_indices, - get_ctx_media, NULL, circ_vrtx_count, get_ctx_position, &ctx) + get_ctx_media, circ_vrtx_count, get_ctx_position, NULL, NULL, &ctx) == RES_OK); } circle_release(&ctx); diff --git a/src/test_senc2d_sample_enclosure.c b/src/test_senc2d_sample_enclosure.c @@ -73,7 +73,7 @@ main(int argc, char** argv) ctx.back_media = medium0; CHK(senc2d_scene_add_geometry(scn, nsegments - 1, get_indices, - get_media, NULL, nvertices, get_position, &ctx) == RES_OK); + get_media, nvertices, get_position, NULL, NULL, &ctx) == RES_OK); CHK(senc2d_scene_analyze(scn, &desc) == RES_OK); diff --git a/src/test_senc2d_scene.c b/src/test_senc2d_scene.c @@ -111,20 +111,19 @@ main(int argc, char** argv) d2(ctx.offset, 0, 0); ctx.front_media = medium0; ctx.back_media = medium1; - ctx.global_ids = gid_face; CHK(senc2d_scene_add_geometry(NULL, nsegments, get_indices, get_media, - get_global_id, nvertices, get_position, &ctx) == RES_BAD_ARG); - CHK(senc2d_scene_add_geometry(scn, 0, get_indices, get_media, get_global_id, - nvertices, get_position, &ctx) == RES_BAD_ARG); + nvertices, get_position, NULL, NULL, &ctx) == RES_BAD_ARG); + CHK(senc2d_scene_add_geometry(scn, 0, get_indices, get_media, + nvertices, get_position, NULL, NULL, &ctx) == RES_BAD_ARG); CHK(senc2d_scene_add_geometry(scn, nsegments, NULL, get_media, - get_global_id, nvertices, get_position, &ctx) == RES_BAD_ARG); + nvertices, get_position, NULL, NULL, &ctx) == RES_BAD_ARG); CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, - get_global_id, 0, get_position, &ctx) == RES_BAD_ARG); + 0, get_position, NULL, NULL, &ctx) == RES_BAD_ARG); CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, - get_global_id, nvertices, NULL, &ctx) == RES_BAD_ARG); + nvertices, NULL, NULL, NULL, &ctx) == RES_BAD_ARG); CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, - get_global_id, nvertices, get_position, &ctx) == RES_OK); + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); CHK(senc2d_scene_get_segments_count(scn, &count) == RES_OK); CHK(count == nsegments); @@ -182,20 +181,13 @@ main(int argc, char** argv) CHK(convention == (SENC2D_CONVENTION_NORMAL_BACK | SENC2D_CONVENTION_NORMAL_INSIDE)); CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, - get_global_id, nvertices, get_position, &ctx) == RES_OK); + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); CHK(senc2d_scene_analyze(scn, &desc) == RES_OK); /* Check that medium 0 is inside */ CHK(senc2d_descriptor_get_enclosure_by_medium(desc, 0, 0, &enc) == RES_OK); CHK(senc2d_enclosure_get_header(enc, &header) == RES_OK); CHK(!header.is_infinite); CHK(senc2d_enclosure_ref_put(enc) == RES_OK); - - FOR_EACH(i, 0, nsegments) { - unsigned gid; - CHK(senc2d_descriptor_get_global_segment_global_id(desc, i, &gid) == RES_OK); - /* gid has been set to gid_face. */ - CHK(gid == gid_face[i]); - } CHK(senc2d_descriptor_get_global_segment_media(desc, 0, medback) == RES_OK); ctx.front_media = medium1_3; @@ -205,7 +197,7 @@ main(int argc, char** argv) SENC2D_CONVENTION_NORMAL_FRONT | SENC2D_CONVENTION_NORMAL_INSIDE, &scn) == RES_OK); CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, - get_global_id, nvertices, get_position, &ctx) == RES_OK); + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); /* Medium mismatch between neighbour segments, but OK */ CHK(senc2d_scene_analyze(scn, &desc) == RES_OK); @@ -227,55 +219,48 @@ main(int argc, char** argv) CHK(senc2d_scene_create(dev, SENC2D_CONVENTION_NORMAL_FRONT | SENC2D_CONVENTION_NORMAL_INSIDE, &scn) == RES_OK); - CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL, - nvertices, get_position, &ctx) == RES_OK); + CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); CHK(senc2d_scene_analyze(scn, &desc) == RES_OK); /* Check that medium 0 is outside */ CHK(senc2d_descriptor_get_enclosure_by_medium(desc, 0, 0, &enc) == RES_OK); CHK(senc2d_enclosure_get_header(enc, &header) == RES_OK); CHK(header.is_infinite); CHK(senc2d_enclosure_ref_put(enc) == RES_OK); - - FOR_EACH(i, 0, nsegments) { - unsigned gid; - CHK(senc2d_descriptor_get_global_segment_global_id(desc, i, &gid) == RES_OK); - /* Default gid: segments rank. */ - CHK(gid == i); - } CHK(senc2d_descriptor_get_global_segment_media(desc, 0, medfront) == RES_OK); FOR_EACH(i, 0, 2) CHK(medback[i] == medfront[i]); /* Invalid vertex ID */ - CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL, - nvertices - 1, get_position, &ctx) == RES_BAD_ARG); + CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices - 1, get_position, NULL, NULL, &ctx) == RES_BAD_ARG); /* Incoherent medium on a duplicate segment */ ctx.back_media = medium1_3; - CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL, - nvertices, get_position, &ctx) == RES_BAD_ARG); + CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, NULL, NULL, &ctx) == RES_BAD_ARG); /* It is OK add geometry after a failed add */ ctx.back_media = medium1; - CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL, - nvertices, get_position, &ctx) == RES_OK); + CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); /* Coherent medium on duplicate segment */ ctx.back_media = medium1; - CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL, - nvertices, get_position, &ctx) == RES_OK); + CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); /* Coherent medium on duplicate segment V2 */ ctx.reverse_med = 1; ctx.front_media = medium1; ctx.back_media = medium0; - CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL, - nvertices, get_position, &ctx) == RES_OK); + CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); /* Coherent medium on duplicate segment V3 */ ctx.reverse_med = 0; ctx.reverse_vrtx = 1; - CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL, - nvertices, get_position, &ctx) == RES_OK); + CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); CHK(senc2d_scene_ref_put(scn) == RES_OK); CHK(senc2d_device_ref_put(dev) == RES_OK); diff --git a/src/test_senc2d_square_behind_square.c b/src/test_senc2d_square_behind_square.c @@ -48,8 +48,8 @@ main(int argc, char** argv) ctx.back_media = medium1; /* First square */ - CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL, - nvertices, get_position, &ctx) == RES_OK); + CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); /* +Y from the first square, * big enough to prevent rays from the first square to miss this one */ @@ -57,8 +57,8 @@ main(int argc, char** argv) ctx.scale = 5; /* Second square */ - CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL, - nvertices, get_position, &ctx) == RES_OK); + CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); CHK(senc2d_scene_analyze(scn, &desc) == RES_OK); @@ -86,8 +86,8 @@ main(int argc, char** argv) ctx.back_media = medium0; /* Third square */ - CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL, - nvertices, get_position, &ctx) == RES_OK); + CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); if(desc) CHK(senc2d_descriptor_ref_put(desc) == RES_OK); desc = NULL; diff --git a/src/test_senc2d_square_in_square.c b/src/test_senc2d_square_in_square.c @@ -50,8 +50,8 @@ main(int argc, char** argv) ctx.back_media = medium1; /* First square */ - CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL, - nvertices, get_position, &ctx) == RES_OK); + CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); d2(ctx.offset, -1, -1); ctx.scale = 3; @@ -60,8 +60,8 @@ main(int argc, char** argv) ctx.reverse_vrtx = 1; /* Second square */ - CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL, - nvertices, get_position, &ctx) == RES_OK); + CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); CHK(senc2d_scene_analyze(scn, &desc) == RES_OK); @@ -92,8 +92,8 @@ main(int argc, char** argv) ctx.back_media = medium0; /* Third square */ - CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL, - nvertices, get_position, &ctx) == RES_OK); + CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); if(desc) CHK(senc2d_descriptor_ref_put(desc) == RES_OK); desc = NULL; diff --git a/src/test_senc2d_square_on_square.c b/src/test_senc2d_square_on_square.c @@ -75,7 +75,7 @@ main(int argc, char** argv) /* Add square #1 */ CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, - NULL, nvertices, get_position, &ctx) == RES_OK); + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); d2(ctx.offset, 1, 2); ctx.scale = 1; @@ -87,7 +87,7 @@ main(int argc, char** argv) /* Add square #2 (has a duplicate face with square #1) */ CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, - NULL, nvertices, get_position, &ctx) == RES_OK); + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); ctx.positions = box_vertices; /* Can use distorded square for square #3 */ d2(ctx.offset, 0, 0); @@ -101,7 +101,7 @@ main(int argc, char** argv) /* Add square #3 */ CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, - NULL, nvertices, get_position, &ctx) == RES_OK); + nvertices, get_position, NULL, NULL, &ctx) == RES_OK); CHK(senc2d_scene_analyze(scn, &desc) == RES_OK); diff --git a/src/test_senc2d_undefined_medium.c b/src/test_senc2d_undefined_medium.c @@ -20,7 +20,7 @@ #include <rsys/double2.h> static void -test(int convention) +test(const int convention) { struct mem_allocator allocator; struct senc2d_descriptor* desc = NULL; @@ -30,17 +30,16 @@ test(int convention) struct senc2d_enclosure_header header; unsigned medium, expected_external_medium, expected_internal_medium; unsigned gid; + enum senc2d_side side; struct context ctx; unsigned i, s, ecount, vcount, scount; - int is_front, is_in; unsigned media[4]; unsigned rev_box_indices[sizeof(box_indices) / sizeof(*box_indices)]; - int conv_front, conv_in; + int conv_front; conv_front = (convention & SENC2D_CONVENTION_NORMAL_FRONT) != 0; - conv_in = (convention & SENC2D_CONVENTION_NORMAL_INSIDE) != 0; - /* Create the box with reversed segments */ + /* Create a box with reversed segments */ FOR_EACH(i, 0, sizeof(rev_box_indices) / sizeof(*rev_box_indices)) { switch (i % 2) { case 0: rev_box_indices[i] = box_indices[i + 1]; break; @@ -51,8 +50,6 @@ test(int convention) OK(senc2d_device_create(NULL, &allocator, SENC2D_NTHREADS_DEFAULT, 1, &dev)); OK(senc2d_scene_create(dev, convention, &scn)); - is_front = (convention & SENC2D_CONVENTION_NORMAL_FRONT) != 0; - is_in = (convention & SENC2D_CONVENTION_NORMAL_INSIDE) != 0; /* A 2D square. * 2 enclosures (inside, outside) sharing the same segments, @@ -72,13 +69,13 @@ test(int convention) for(i = 0; i < sizeof(media) / sizeof(*media); i++) media[i] = SENC2D_UNDEFINED_MEDIUM; OK(senc2d_scene_add_geometry(scn, nsegments, get_indices, NULL, - NULL, nvertices, get_position, &ctx)); + nvertices, get_position, NULL, NULL, &ctx)); OK(senc2d_scene_get_unique_sides_without_medium_count(scn, &scount)); CHK(scount == 2 * nsegments); /* Add geometry with no media information on the front sides */ OK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, - NULL, nvertices, get_position, &ctx)); + nvertices, get_position, NULL, NULL, &ctx)); OK(senc2d_scene_get_unique_sides_without_medium_count(scn, &scount)); CHK(scount == nsegments); @@ -120,7 +117,8 @@ test(int convention) FOR_EACH(s, 0, header.segment_count) { unsigned ind[2]; - OK(senc2d_enclosure_get_segment_global_id(enclosure, s, &gid)); + OK(senc2d_enclosure_get_segment_global_id(enclosure, s, &gid, &side)); + CHK(side == (medium == 1) ? SENC2D_BACK : SENC2D_FRONT); CHK(gid == s); OK(senc2d_enclosure_get_segment(enclosure, s, ind)); } @@ -133,7 +131,7 @@ test(int convention) for(i = 0; i < sizeof(media) / sizeof(*media); i++) media[i] = (i % 2) ? 0 : SENC2D_UNDEFINED_MEDIUM; OK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, - NULL, nvertices, get_position, &ctx)); + nvertices, get_position, NULL, NULL, &ctx)); OK(senc2d_scene_get_unique_sides_without_medium_count(scn, &scount)); CHK(scount == nsegments / 2); @@ -163,7 +161,7 @@ test(int convention) ctx.indices = rev_box_indices; SWAP(const unsigned*, ctx.front_media, ctx.back_media); OK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, - NULL, nvertices, get_position, &ctx)); + nvertices, get_position, NULL, NULL, &ctx)); OK(senc2d_scene_get_unique_sides_without_medium_count(scn, &scount)); CHK(scount == nsegments / 2); @@ -175,7 +173,7 @@ test(int convention) for(i = 0; i < sizeof(media) / sizeof(*media); i++) media[i] = (i % 2) ? SENC2D_UNDEFINED_MEDIUM : 0; OK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, - NULL, nvertices, get_position, &ctx)); + nvertices, get_position, NULL, NULL, &ctx)); OK(senc2d_scene_get_unique_sides_without_medium_count(scn, &scount)); CHK(scount == 0); @@ -205,7 +203,7 @@ test(int convention) for(i = 0; i < sizeof(media) / sizeof(*media); i++) media[i] = 0; OK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, - NULL, nvertices, get_position, &ctx)); + nvertices, get_position, NULL, NULL, &ctx)); OK(senc2d_scene_get_unique_sides_without_medium_count(scn, &scount)); CHK(scount == 0); @@ -213,7 +211,7 @@ test(int convention) for(i = 0; i < sizeof(media) / sizeof(*media); i++) media[i] = (i % 2); BA(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, - NULL, nvertices, get_position, &ctx)); + nvertices, get_position, NULL, NULL, &ctx)); OK(senc2d_scene_get_unique_sides_without_medium_count(scn, &scount)); CHK(scount == 0); @@ -239,7 +237,7 @@ test(int convention) /* Geometrical normals point outside the cube in input segments: * if convention is front, front medium (0) is outside, * that is medium 0's enclosure is infinite */ - CHK(is_front == ((medium == 0) == header.is_infinite)); + CHK(conv_front == ((medium == 0) == header.is_infinite)); CHK(header.segment_count == nsegments); CHK(header.unique_segment_count == nsegments); CHK(header.vertices_count == nvertices); @@ -253,7 +251,8 @@ test(int convention) OK(senc2d_enclosure_ref_put(ee)); FOR_EACH(s, 0, header.segment_count) { - OK(senc2d_enclosure_get_segment_global_id(enclosure, s, &gid)); + OK(senc2d_enclosure_get_segment_global_id(enclosure, s, &gid, &side)); + CHK(side == (medium == 1) ? SENC2D_BACK : SENC2D_FRONT); CHK(gid == s); } OK(senc2d_enclosure_ref_put(enclosure)); diff --git a/src/test_senc2d_undefined_medium_attr.c b/src/test_senc2d_undefined_medium_attr.c @@ -0,0 +1,366 @@ +/* Copyright (C) |Meso|Star> 2016-2018 (contact@meso-star.com) +* +* This program is free software: you can redistribute it and/or modify +* it under the terms of the GNU General Public License as published by +* the Free Software Foundation, either version 3 of the License, or +* (at your option) any later version. +* +* This program is distributed in the hope that it will be useful, +* but WITHOUT ANY WARRANTY; without even the implied warranty of +* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the +* GNU General Public License for more details. +* +* You should have received a copy of the GNU General Public License +* along with this program. If not, see <http://www.gnu.org/licenses/>. */ + +#include "senc2d.h" +#include "senc2d_s2d_wrapper.h" +#include "test_senc2d_utils.h" + +#include <rsys/double2.h> + +#include <limits.h> + +#define INVALID_INTFACE_ID UINT_MAX + +static FINLINE void +set_null_id(struct mem_allocator* alloc, unsigned* data) +{ + ASSERT(data); (void)alloc; + *data = INVALID_INTFACE_ID; +} + +#include <rsys/dynamic_array.h> +#define DARRAY_NAME intface_id +#define DARRAY_FUNCTOR_INIT set_null_id +#define DARRAY_DATA unsigned +#include <rsys/dynamic_array.h> + +/* Manage interface properties */ +struct merge_ctx { + struct darray_intface_id global_interf_data; + const unsigned* current_add_interf_data; +}; + +static res_T +add_seg + (const unsigned global_id, + const unsigned iseg, + void* context) +{ + res_T res = RES_OK; + struct context* ctx = context; + struct merge_ctx* merge_ctx; + unsigned interf; + ASSERT(ctx); + merge_ctx = ctx->custom; + /* Get interface information from ctx */ + interf = merge_ctx->current_add_interf_data[iseg]; + /* Keep data */ + res = darray_intface_id_resize(&merge_ctx->global_interf_data, global_id + 1); + if (res != RES_OK) return res; + darray_intface_id_data_get(&merge_ctx->global_interf_data)[global_id] = interf; + return res; +} + +static res_T +merge_seg + (const unsigned global_id, + const unsigned iseg, + const int reversed_segment, + void* context) +{ + res_T res = RES_OK; + struct context* ctx = context; + struct merge_ctx* merge_ctx; + int need_merge; + unsigned interf; + unsigned* interf_data; + ASSERT(ctx); (void)reversed_segment; + merge_ctx = ctx->custom; + res = darray_intface_id_resize(&merge_ctx->global_interf_data, global_id + 1); + if (res != RES_OK) return res; + /* Get interface information from ctx */ + interf = merge_ctx->current_add_interf_data[iseg]; + interf_data = darray_intface_id_data_get(&merge_ctx->global_interf_data); + + need_merge = (interf_data[global_id] != INVALID_INTFACE_ID); + if (need_merge) { + if (interf_data[global_id] != interf) + /* Previous interface id is different: no possible merge */ + return RES_BAD_ARG; + } else { + /* Triangle is known, but without interface information: create */ + interf_data[global_id] = interf; + } + /* Same data: no merge needed */ + return RES_OK; +} + +static void +test(const int convention) +{ + struct mem_allocator allocator; + struct senc2d_descriptor* desc = NULL; + struct senc2d_device* dev = NULL; + struct senc2d_scene* scn = NULL; + struct senc2d_enclosure* enclosure; + struct senc2d_enclosure_header header; + unsigned medium, expected_external_medium, expected_internal_medium; + unsigned gid; + enum senc2d_side side; + struct context ctx; + unsigned i, s, ecount, vcount, scount; + unsigned media[4], interface_ids[4] = { 0, 0, 1, 1, }; + unsigned rev_box_indices[sizeof(box_indices) / sizeof(*box_indices)]; + struct merge_ctx merge_ctx; + int conv_front = (convention & SENC2D_CONVENTION_NORMAL_FRONT) != 0; + + darray_intface_id_init(&allocator, &merge_ctx.global_interf_data); + merge_ctx.current_add_interf_data = interface_ids; + + /* Create a box with reversed segments */ + FOR_EACH(i, 0, sizeof(rev_box_indices) / sizeof(*rev_box_indices)) { + switch (i % 2) { + case 0: rev_box_indices[i] = box_indices[i + 1]; break; + case 1: rev_box_indices[i] = box_indices[i - 1]; break; + } + } + OK(mem_init_proxy_allocator(&allocator, &mem_default_allocator)); + OK(senc2d_device_create(NULL, &allocator, SENC2D_NTHREADS_DEFAULT, 1, &dev)); + + OK(senc2d_scene_create(dev, convention, &scn)); + + /* A 2D square. + * 2 enclosures (inside, outside) sharing the same segments, + * but opposite sides */ + ctx.positions = box_vertices; + ctx.indices = box_indices; + ctx.scale = 1; + ctx.reverse_vrtx = 0; + ctx.reverse_med = 0; + d2(ctx.offset, 0, 0); + ctx.front_media = media; + ctx.back_media = medium1; + ctx.custom = &merge_ctx; + + /* Can add the same segments again defined/undefined media in any order */ + + /* Add geometry with no media information on both sides */ + for (i = 0; i < sizeof(media) / sizeof(*media); i++) + media[i] = SENC2D_UNDEFINED_MEDIUM; + OK(senc2d_scene_add_geometry(scn, nsegments, get_indices, NULL, + nvertices, get_position, add_seg, merge_seg, &ctx)); + OK(senc2d_scene_get_unique_sides_without_medium_count(scn, &scount)); + CHK(scount == 2 * nsegments); + + /* If merge fails, add geometry fails */ + interface_ids[0] = 6; + BA(senc2d_scene_add_geometry(scn, nsegments, get_indices, NULL, + nvertices, get_position, add_seg, merge_seg, &ctx)); + interface_ids[0] = 0; + + /* Add geometry with no media information on the front sides */ + OK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, add_seg, merge_seg, &ctx)); + OK(senc2d_scene_get_unique_sides_without_medium_count(scn, &scount)); + CHK(scount == nsegments); + + /* Analyze with undefined media on the front sides */ + OK(senc2d_scene_analyze(scn, &desc)); + OK(senc2d_descriptor_get_enclosure_count(desc, &ecount)); + CHK(ecount == 2); + + FOR_EACH(i, 0, ecount) { + struct senc2d_enclosure* ee; + struct senc2d_enclosure_header hh; + unsigned cc; + OK(senc2d_descriptor_get_enclosure(desc, i, &enclosure)); + OK(senc2d_enclosure_get_header(enclosure, &header)); + + CHK(header.enclosure_id == i); + CHK(header.enclosed_media_count == 1); + + OK(senc2d_enclosure_get_medium(enclosure, 0, &medium)); + /* Geometrical normals point outside the cube in input segments: + * if convention is front, front medium (undef) is outside, + * that is medium 0's enclosure is infinite */ + expected_external_medium = conv_front ? SENC2D_UNDEFINED_MEDIUM : 1; + expected_internal_medium = conv_front ? 1 : SENC2D_UNDEFINED_MEDIUM; + + CHK(medium == (header.is_infinite + ? expected_external_medium : expected_internal_medium)); + CHK(header.segment_count == nsegments); + CHK(header.unique_segment_count == nsegments); + CHK(header.vertices_count == nvertices); + CHK(header.is_infinite == (i == 0)); + + OK(senc2d_descriptor_get_enclosure_count_by_medium(desc, medium, &cc)); + CHK(cc == 1); + OK(senc2d_descriptor_get_enclosure_by_medium(desc, medium, 0, &ee)); + OK(senc2d_enclosure_get_header(ee, &hh)); + CHK(header.enclosure_id == hh.enclosure_id); + OK(senc2d_enclosure_ref_put(ee)); + + FOR_EACH(s, 0, header.segment_count) { + unsigned ind[2]; + OK(senc2d_enclosure_get_segment_global_id(enclosure, s, &gid, &side)); + CHK(gid == s); + CHK(side == (medium == 1) ? SENC2D_BACK : SENC2D_FRONT); + OK(senc2d_enclosure_get_segment(enclosure, s, ind)); + } + OK(senc2d_enclosure_ref_put(enclosure)); + } + OK(senc2d_descriptor_ref_put(desc)); + + /* Same geometry, front media are defined for odd segments */ + for (i = 0; i < sizeof(media) / sizeof(*media); i++) + media[i] = (i % 2) ? 0 : SENC2D_UNDEFINED_MEDIUM; + OK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, add_seg, merge_seg, &ctx)); + OK(senc2d_scene_get_unique_sides_without_medium_count(scn, &scount)); + CHK(scount == nsegments / 2); + + /* Analyze with undefined media */ + OK(senc2d_scene_analyze(scn, &desc)); + OK(senc2d_descriptor_ref_put(desc)); + + /* Get the deduplicated geometry without (successful) analyze */ + OK(senc2d_scene_get_unique_vertices_count(scn, &vcount)); + CHK(vcount == nvertices); + OK(senc2d_scene_get_unique_segments_count(scn, &scount)); + CHK(scount == nsegments); + FOR_EACH(i, 0, scount) { + int j; + unsigned med[2], ids[2]; + OK(senc2d_scene_get_unique_segment(scn, i, ids)); + OK(senc2d_scene_get_unique_segment_media(scn, i, med)); + CHK(med[0] == ((i % 2) ? 0 : SENC2D_UNDEFINED_MEDIUM) && med[1] == 1); + FOR_EACH(j, 0, 2) { + double pos[2]; + CHK(ids[j] < vcount); + OK(senc2d_scene_get_unique_vertex(scn, ids[j], pos)); + } + } + + /* Same information again, using a reversed box */ + ctx.indices = rev_box_indices; + SWAP(const unsigned*, ctx.front_media, ctx.back_media); + OK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, add_seg, merge_seg, &ctx)); + OK(senc2d_scene_get_unique_sides_without_medium_count(scn, &scount)); + CHK(scount == nsegments / 2); + + /* Analyze with undefined media */ + OK(senc2d_scene_analyze(scn, &desc)); + OK(senc2d_descriptor_ref_put(desc)); + + /* Define media for remaining segments, using reversed box */ + for (i = 0; i < sizeof(media) / sizeof(*media); i++) + media[i] = (i % 2) ? SENC2D_UNDEFINED_MEDIUM : 0; + OK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, add_seg, merge_seg, &ctx)); + OK(senc2d_scene_get_unique_sides_without_medium_count(scn, &scount)); + CHK(scount == 0); + + /* Get the deduplicated geometry without (successful) analyze */ + OK(senc2d_scene_get_unique_vertices_count(scn, &vcount)); + CHK(vcount == nvertices); + OK(senc2d_scene_get_unique_segments_count(scn, &scount)); + CHK(scount == nsegments); + FOR_EACH(i, 0, scount) { + int j; + unsigned med[2], ids[2]; + OK(senc2d_scene_get_unique_segment(scn, i, ids)); + OK(senc2d_scene_get_unique_segment_media(scn, i, med)); + CHK(med[0] == 0 && med[1] == 1); + FOR_EACH(j, 0, 2) { + double pos[2]; + CHK(ids[j] < vcount); + OK(senc2d_scene_get_unique_vertex(scn, ids[j], pos)); + } + } + + /* Analyze with all media defined */ + OK(senc2d_scene_analyze(scn, &desc)); + OK(senc2d_descriptor_ref_put(desc)); + + /* Define media for all segments, nothing new here */ + for (i = 0; i < sizeof(media) / sizeof(*media); i++) + media[i] = 0; + OK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, add_seg, merge_seg, &ctx)); + OK(senc2d_scene_get_unique_sides_without_medium_count(scn, &scount)); + CHK(scount == 0); + + /* Define incoherent media for some segments */ + for (i = 0; i < sizeof(media) / sizeof(*media); i++) + media[i] = (i % 2); + BA(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, + nvertices, get_position, add_seg, merge_seg, &ctx)); + OK(senc2d_scene_get_unique_sides_without_medium_count(scn, &scount)); + CHK(scount == 0); + + /* Scene is still OK and can be analyzed */ + OK(senc2d_scene_analyze(scn, &desc)); + + OK(senc2d_descriptor_get_global_segments_count(desc, &scount)); + CHK(scount == sizeof(media) / sizeof(*media)); + + OK(senc2d_descriptor_get_enclosure_count(desc, &ecount)); + CHK(ecount == 2); + + FOR_EACH(i, 0, ecount) { + struct senc2d_enclosure* ee; + struct senc2d_enclosure_header hh; + unsigned cc; + OK(senc2d_descriptor_get_enclosure(desc, i, &enclosure)); + OK(senc2d_enclosure_get_header(enclosure, &header)); + + CHK(header.enclosure_id == i); + CHK(header.enclosed_media_count == 1); + OK(senc2d_enclosure_get_medium(enclosure, 0, &medium)); + /* Geometrical normals point outside the cube in input segments: + * if convention is front, front medium (0) is outside, + * that is medium 0's enclosure is infinite */ + CHK(conv_front == ((medium == 0) == header.is_infinite)); + CHK(header.segment_count == nsegments); + CHK(header.unique_segment_count == nsegments); + CHK(header.vertices_count == nvertices); + CHK(header.is_infinite == (i == 0)); + + OK(senc2d_descriptor_get_enclosure_count_by_medium(desc, medium, &cc)); + CHK(cc == 1); + OK(senc2d_descriptor_get_enclosure_by_medium(desc, medium, 0, &ee)); + OK(senc2d_enclosure_get_header(ee, &hh)); + CHK(header.enclosure_id == hh.enclosure_id); + OK(senc2d_enclosure_ref_put(ee)); + + FOR_EACH(s, 0, header.segment_count) { + OK(senc2d_enclosure_get_segment_global_id(enclosure, s, &gid, &side)); + CHK(gid == s); + CHK(side == (medium == 1) ? SENC2D_BACK : SENC2D_FRONT); + } + OK(senc2d_enclosure_ref_put(enclosure)); + } + + SENC2D(scene_ref_put(scn)); + SENC2D(device_ref_put(dev)); + SENC2D(descriptor_ref_put(desc)); + darray_intface_id_release(&merge_ctx.global_interf_data); + + check_memory_allocator(&allocator); + mem_shutdown_proxy_allocator(&allocator); + CHK(mem_allocated_size() == 0); +} + +int +main(int argc, char** argv) +{ + (void)argc, (void)argv; + test(SENC2D_CONVENTION_NORMAL_FRONT | SENC2D_CONVENTION_NORMAL_INSIDE); + test(SENC2D_CONVENTION_NORMAL_BACK | SENC2D_CONVENTION_NORMAL_INSIDE); + test(SENC2D_CONVENTION_NORMAL_FRONT | SENC2D_CONVENTION_NORMAL_OUTSIDE); + test(SENC2D_CONVENTION_NORMAL_BACK | SENC2D_CONVENTION_NORMAL_OUTSIDE); + return 0; +} diff --git a/src/test_senc2d_utils.h b/src/test_senc2d_utils.h @@ -71,7 +71,7 @@ struct context { unsigned* indices; const unsigned* front_media; const unsigned* back_media; - const unsigned* global_ids; + void* custom; double offset[2]; double scale; char reverse_vrtx, reverse_med; @@ -84,10 +84,8 @@ static const unsigned medium1_3[4] = { 1, 1, 3, 1 }; static const unsigned medium1_back0[4] = { 1, 1, 1, 0 }; static const unsigned medium1_front0[4] = { 1, 0, 1, 1 }; -static const unsigned gid_face[4] = { 0, 1, 2, 3 }; - static INLINE void -get_indices(const unsigned iseg, unsigned ids[2], const void* context) +get_indices(const unsigned iseg, unsigned ids[2], void* context) { const struct context* ctx = context; ASSERT(ids && ctx); @@ -96,7 +94,7 @@ get_indices(const unsigned iseg, unsigned ids[2], const void* context) } static INLINE void -get_position(const unsigned ivert, double pos[2], const void* context) +get_position(const unsigned ivert, double pos[2], void* context) { const struct context* ctx = context; ASSERT(pos && ctx); @@ -105,7 +103,7 @@ get_position(const unsigned ivert, double pos[2], const void* context) } static INLINE void -get_media(const unsigned iseg, unsigned medium[2], const void* context) +get_media(const unsigned iseg, unsigned medium[2], void* context) { const struct context* ctx = context; ASSERT(medium && ctx); @@ -113,14 +111,6 @@ get_media(const unsigned iseg, unsigned medium[2], const void* context) medium[ctx->reverse_med ? 0 : 1] = ctx->back_media[iseg]; } -static INLINE void -get_global_id(const unsigned iseg, unsigned* gid, const void* context) -{ - const struct context* ctx = context; - ASSERT(gid && context); - *gid = ctx->global_ids[iseg]; -} - /******************************************************************************* * Miscellaneous ******************************************************************************/