star-enclosures-3d

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

commit 2169ed6ff92ad461dcff2eb14dd0f88400c72f57
parent 6c6304774b8b0711ab3af0289df320dac6f6422a
Author: Christophe Coustet <christophe.coustet@meso-star.com>
Date:   Wed, 28 Feb 2018 14:26:08 +0100

Build Linux.

Diffstat:
Mcmake/CMakeLists.txt | 3++-
Msrc/senc_descriptor.c | 4++++
Msrc/senc_enclosure.c | 2++
Msrc/senc_enclosure_data.h | 13++++++++-----
Msrc/senc_s3d_wrapper.h | 16++++++++--------
Msrc/senc_scene.c | 8+++++---
Msrc/senc_scene_analyze.c | 117++++++++++++++++++++++++++++++++-----------------------------------------------
Msrc/senc_scene_analyze_c.h | 34+++++++++++++++++-----------------
Msrc/senc_scene_c.h | 6+++---
Msrc/test_senc_cube_on_cube.c | 4++++
Msrc/test_senc_many_triangles.c | 6+++---
Msrc/test_senc_utils.h | 14+++++++-------
12 files changed, 110 insertions(+), 117 deletions(-)

diff --git a/cmake/CMakeLists.txt b/cmake/CMakeLists.txt @@ -87,7 +87,8 @@ set_target_properties(senc PROPERTIES rcmake_copy_runtime_libraries(senc) if(CMAKE_COMPILER_IS_GNUCC) - set_target_properties(senc PROPERTIES LINK_FLAGS "${OpenMP_C_FLAGS} -lm") + set_target_properties(senc PROPERTIES LINK_FLAGS "${OpenMP_C_FLAGS}") + target_link_libraries(senc m) endif() rcmake_setup_devel(senc StarEnc ${VERSION} senc_version.h) diff --git a/src/senc_descriptor.c b/src/senc_descriptor.c @@ -180,7 +180,9 @@ senc_descriptor_get_global_triangle_media return RES_BAD_ARG; trg = darray_triangle_in_cdata_get(&desc->scene->triangles_in) + itri; FOR_EACH(i, 0, 2) { +#if (UINT_MAX < MEDIUM_MAX__) ASSERT(trg->medium[i] < UINT_MAX); +#endif media[i] = (unsigned)trg->medium[i]; /* Back to API type */ } return RES_OK; @@ -199,7 +201,9 @@ senc_descriptor_get_global_triangle_enclosures return RES_BAD_ARG; trg = darray_triangle_enc_cdata_get(&desc->triangles_enc) + itri; FOR_EACH(i, 0, 2) { +#if (UINT_MAX < ENCLOSURE_MAX__) ASSERT(trg->enclosure[i] < UINT_MAX); +#endif enclosures[i] = (unsigned)trg->enclosure[i]; /* Back to API type */ } return RES_OK; diff --git a/src/senc_enclosure.c b/src/senc_enclosure.c @@ -130,7 +130,9 @@ senc_enclosure_get_triangle_media == enclosure->data->header.triangle_count); triangle = darray_triangle_in_cdata_get(&enclosure->data->sides) + itri; FOR_EACH(i, 0, 2) { +#if (UINT_MAX < MEDIUM_MAX__) ASSERT(triangle->medium[i] < UINT_MAX); +#endif medium[i] = (unsigned)triangle->medium[i]; /* Back to API type */ } return RES_OK; diff --git a/src/senc_enclosure_data.h b/src/senc_enclosure_data.h @@ -18,11 +18,14 @@ #include <rsys/rsys.h> #include <rsys/ref_count.h> +#include <rsys/dynamic_array_uchar.h> #include "senc.h" #include "senc_scene_c.h" #include "senc_internal_types.h" +#include <limits.h> + static void init_header(struct enclosure_header* header) { @@ -42,7 +45,7 @@ struct enclosure_data { /* Index of vertices in scene's unique vertices */ struct darray_vrtx_id vertices; /* TODO: use only 1 bit per side (now 1 char per triangle) */ - struct darray_char side_membership; + struct darray_uchar side_membership; /* Number of components involved in this enclosure */ component_id_t cc_count; /* Linked list of the components */ @@ -57,7 +60,7 @@ enclosure_data_init(struct mem_allocator* alloc, struct enclosure_data* enc) { enc->first_component = COMPONENT_NULL__; darray_triangle_in_init(alloc, &enc->sides); darray_vrtx_id_init(alloc, &enc->vertices); - darray_char_init(alloc, &enc->side_membership); + darray_uchar_init(alloc, &enc->side_membership); } static FINLINE res_T @@ -72,7 +75,7 @@ enclosure_data_copy dst->first_component = src->first_component; OK(darray_triangle_in_copy(&dst->sides, &src->sides)); OK(darray_vrtx_id_copy(&dst->vertices, &src->vertices)); - OK(darray_char_copy(&dst->side_membership, &src->side_membership)); + OK(darray_uchar_copy(&dst->side_membership, &src->side_membership)); error: return res; } @@ -82,7 +85,7 @@ enclosure_data_release(struct enclosure_data* n) { ASSERT(n); darray_triangle_in_release(&n->sides); darray_vrtx_id_release(&n->vertices); - darray_char_release(&n->side_membership); + darray_uchar_release(&n->side_membership); } static FINLINE res_T @@ -97,7 +100,7 @@ enclosure_data_copy_and_release dst->first_component = src->first_component; OK(darray_triangle_in_copy_and_release(&dst->sides, &src->sides)); OK(darray_vrtx_id_copy_and_release(&dst->vertices, &src->vertices)); - OK(darray_char_copy(&dst->side_membership, &src->side_membership)); + OK(darray_uchar_copy(&dst->side_membership, &src->side_membership)); error: return res; } diff --git a/src/senc_s3d_wrapper.h b/src/senc_s3d_wrapper.h @@ -21,7 +21,7 @@ #include <rsys/rsys.h> #include <rsys/float3.h> -void FINLINE +FINLINE void senc_descriptor_get_global_indices__ (const unsigned itri, unsigned indices[3], @@ -31,10 +31,10 @@ senc_descriptor_get_global_indices__ res_T r; ASSERT(indices && ctx); r = senc_descriptor_get_global_triangle(descriptor, itri, indices); - ASSERT(r == RES_OK); + ASSERT(r == RES_OK); (void)r; } -void FINLINE +FINLINE void senc_descriptor_get_global_vertices__ (const unsigned ivert, float coord[3], @@ -45,11 +45,11 @@ senc_descriptor_get_global_vertices__ res_T r; ASSERT(coord && ctx); r = senc_descriptor_get_global_vertex(descriptor, ivert, tmp); - ASSERT(r == RES_OK); + ASSERT(r == RES_OK); (void)r; f3_set_d3(coord, tmp); } -void FINLINE +FINLINE void senc_enclosure_get_triangle__ (const unsigned itri, unsigned indices[3], @@ -59,10 +59,10 @@ senc_enclosure_get_triangle__ res_T r; ASSERT(indices && ctx); r = senc_enclosure_get_triangle(enclosure, itri, indices); - ASSERT(r == RES_OK); + ASSERT(r == RES_OK); (void)r; } -void FINLINE +FINLINE void senc_enclosure_get_vertex__ (const unsigned ivert, float coord[3], @@ -73,7 +73,7 @@ senc_enclosure_get_vertex__ res_T r; ASSERT(coord && ctx); r = senc_enclosure_get_vertex(enclosure, ivert, tmp); - ASSERT(r == RES_OK); + ASSERT(r == RES_OK); (void)r; f3_set_d3(coord, tmp); } diff --git a/src/senc_scene.c b/src/senc_scene.c @@ -180,7 +180,7 @@ senc_scene_add_geometry || tmp.vertice_id[1] == tmp.vertice_id[2]) { const union double3* positions = darray_position_cdata_get(&scn->vertices); - log_err(scn->dev, "%s: triangle %u is degenerate.\n", + log_err(scn->dev, "%s: triangle %lu is degenerate.\n", FUNC_NAME, (unsigned long)tmp.global_id); log_err(scn->dev, " (%g %g %g)\n (%g %g %g)\n (%g %g %g)\n", SPLIT3(positions[trg[i].vertice_id[0]].vec), @@ -232,8 +232,10 @@ senc_scene_add_geometry SPLIT3(positions[trg[*p_trg].vertice_id[1]].vec), SPLIT3(positions[trg[*p_trg].vertice_id[2]].vec)); log_err(scn->dev, "Media: (%lu, %lu) VS (%lu, %lu)\n", - umed[ureversed? 1 : 0], umed[ureversed ? 0 : 1], - med[reversed ? 1 : 0], med[reversed ? 0 : 1]); + (unsigned long)umed[ureversed? 1 : 0], + (unsigned long)umed[ureversed ? 0 : 1], + (unsigned long)med[reversed ? 1 : 0], + (unsigned long)med[reversed ? 0 : 1]); res = RES_BAD_ARG; goto error; } else { diff --git a/src/senc_scene_analyze.c b/src/senc_scene_analyze.c @@ -36,7 +36,7 @@ #define CC_DESCRIPTOR_NULL__ {\ {0,0,-DBL_MAX}, -1, SIDE_NULL__, VRTX_NULL__, 0, MEDIUM_NULL__,\ CC_ID_NONE, CC_GROUP_ROOT_NONE, CC_GROUP_ID_NONE, CC_ID_NONE, TRG_NULL__,\ - TRG_NULL__\ + TRG_NULL__, { NULL, 0, 0, NULL }\ } const struct cc_descriptor CC_DESCRIPTOR_NULL = CC_DESCRIPTOR_NULL__; @@ -51,32 +51,7 @@ const struct cc_descriptor CC_DESCRIPTOR_NULL = CC_DESCRIPTOR_NULL__; /******************************************************************************* * Helper function ******************************************************************************/ -static void -dumplist - (const struct trgside* trgsides, - const struct darray_side_id* side_ids, - const enum list_id list_id) -{ - side_id_t i; - size_t tmp; - (void)list_id; - ASSERT(trgsides && side_ids); - printf("\n"); - tmp = darray_side_id_size_get(side_ids); - ASSERT(tmp <= SIDE_MAX__); - FOR_EACH(i, 0, (side_id_t)tmp) { - const side_id_t id = darray_side_id_cdata_get(side_ids)[i]; - const struct trgside* const side = trgsides + id; - printf("Side %lu (%lu %lu %lu)\n", - (unsigned long)id, - (unsigned long)side->facing_side_id[0], - (unsigned long)side->facing_side_id[1], - (unsigned long)side->facing_side_id[2]); - ASSERT(side->list_id & list_id); - } -} - -static int +static INLINE int find_side_in_list (const struct trgside* trgsides, const struct darray_side_id* side_ids, @@ -117,13 +92,13 @@ find_component_Zmax struct triangle_in* triangles_in; struct triangle_tmp* triangles_tmp; const union double3* vertices; - const char* side_membership; + const unsigned char* side_membership; ASSERT(scn && triangles_tmp_array && cc); vertices = darray_position_cdata_get(&scn->vertices); triangles_in = darray_triangle_in_data_get(&scn->triangles_in); triangles_tmp = darray_triangle_tmp_data_get(triangles_tmp_array); - side_membership = darray_char_cdata_get(&cc->side_membership); + side_membership = darray_uchar_cdata_get(&cc->side_membership); /* Build the sorted list of side ids */ ASSERT(darray_triangle_tmp_size_get(triangles_tmp_array) @@ -131,7 +106,7 @@ find_component_Zmax ASSERT(scn->nutris == darray_triangle_in_size_get(&scn->triangles_in)); FOR_EACH(trid, 0, scn->nutris) { - const char member = side_membership[trid]; + const unsigned char member = side_membership[trid]; struct triangle_in* const trg_in = triangles_in + trid; struct triangle_tmp* const trg_tmp = triangles_tmp + trid; enum side_id side; @@ -155,7 +130,7 @@ find_component_Zmax vertices[trg_in->vertice_id[0]].vec); d3_cross(normal, edge0, edge1); norm = d3_normalize(normal, normal); - ASSERT(norm); + ASSERT(norm); (void)norm; if((member & FLAG_FRONT) && (member & FLAG_BACK)) { /* Select the side with nz>0 */ @@ -268,7 +243,7 @@ extract_connex_components trgsides, &first_side_not_in_component); side_id_t crt_side_id = start_side_id; side_id_t last_side_id = start_side_id; - char* side_membership; + unsigned char* side_membership; size_t sz; ASSERT(start_side_id == SIDE_NULL__ || start_side_id < 2 * scn->nutris); @@ -280,12 +255,12 @@ extract_connex_components sz = darray_cc_descriptor_size_get(connex_components); OK(darray_cc_descriptor_resize(connex_components, 1+ sz)); current_cc = darray_cc_descriptor_data_get(connex_components) + sz; - /* 1 char per triangle (2 sides); allow uint64_t* access */ - OK(darray_char_reserve(&current_cc->side_membership, + /* 1 uchar per triangle (2 sides); allow uint64_t* access */ + OK(darray_uchar_reserve(&current_cc->side_membership, ((scn->nutris + 7) / 8) * 8)); - OK(darray_char_resize(&current_cc->side_membership, scn->nutris)); - side_membership = darray_char_data_get(&current_cc->side_membership); - memset(side_membership, 0, scn->nutris * sizeof(char)); + OK(darray_uchar_resize(&current_cc->side_membership, scn->nutris)); + side_membership = darray_uchar_data_get(&current_cc->side_membership); + memset(side_membership, 0, scn->nutris * sizeof(unsigned char)); m = trgsides[start_side_id].medium; ASSERT(sz <= COMPONENT_MAX__); current_cc->cc_id = (component_id_t)sz; @@ -350,7 +325,7 @@ extract_connex_components SPLIT3(positions[triangles_in[nbour_trg_id].vertice_id[1]].vec), SPLIT3(positions[triangles_in[nbour_trg_id].vertice_id[2]].vec)); log_err(desc->scene->dev, "Media: %lu VS %lu\n", - neighbour->medium, crt_side->medium); + (unsigned long)neighbour->medium, (unsigned long)crt_side->medium); res = RES_BAD_ARG; goto error; } @@ -386,7 +361,7 @@ error: goto exit; } -void +static void get_scn_indices(const unsigned itri, unsigned ids[3], void* ctx) { int i; const struct senc_scene* scene = ctx; @@ -398,7 +373,7 @@ get_scn_indices(const unsigned itri, unsigned ids[3], void* ctx) { } } -void +static void get_scn_position(const unsigned ivert, float pos[3], void* ctx) { const struct senc_scene* scene = ctx; const union double3* pt = @@ -452,7 +427,7 @@ group_connex_components component_id_t cc_count, c; medium_id_t infinite_medium = MEDIUM_NULL__; side_id_t infinite_medium_first_side = SIDE_MAX__; - char infinite_medium_is_known = 0; + unsigned char infinite_medium_is_known = 0; (void)trgsides; ASSERT(desc && trgsides && triangles_comp && connex_components); @@ -556,7 +531,7 @@ group_connex_components SPLIT3(positions[triangles_in[t2].vertice_id[1]].vec), SPLIT3(positions[triangles_in[t2].vertice_id[2]].vec)); log_err(desc->scene->dev, "Media: %lu VS %lu\n", - infinite_medium, cc->medium); + (unsigned long)infinite_medium, (unsigned long)cc->medium); res = RES_BAD_ARG; goto error; } @@ -584,7 +559,7 @@ group_connex_components hit_desc = darray_cc_descriptor_cdata_get(connex_components) + cc->cc_group_root; ASSERT(hit_desc->medium == hit_trg_in->medium[hit_side]); - ASSERT(darray_char_cdata_get(&hit_desc->side_membership)[hit_trg_id] + ASSERT(darray_uchar_cdata_get(&hit_desc->side_membership)[hit_trg_id] & TRGSIDE_IS_FRONT(hit_side) ? FLAG_FRONT : FLAG_BACK); } #endif @@ -598,7 +573,7 @@ group_connex_components = darray_position_cdata_get(&desc->scene->vertices); log_err(desc->scene->dev, "Medium mismatch found between triangle %lu %s side and triangle" - " %u %s side facing each other.\n", + " %lu %s side facing each other.\n", (unsigned long)triangles_in[t1].global_id, TRGSIDE_IS_FRONT(hit_side) ? "front" : "back", (unsigned long)triangles_in[t2].global_id, @@ -616,7 +591,8 @@ group_connex_components SPLIT3(positions[triangles_in[t2].vertice_id[1]].vec), SPLIT3(positions[triangles_in[t2].vertice_id[2]].vec)); log_err(desc->scene->dev, "Media: %lu VS %lu\n", - hit_trg_in->medium[hit_side], cc->medium); + (unsigned long)hit_trg_in->medium[hit_side], + (unsigned long)cc->medium); res = RES_BAD_ARG; goto error; } @@ -704,13 +680,13 @@ collect_and_link_neighbours /* Loop on triangles to register edges. */ FOR_EACH(t, 0, scn->nutris) { struct trg_edge edge; - char ee; + unsigned char ee; FOR_EACH(ee, 0, 3) { edge_id_t* p_id; const vrtx_id_t v0 = triangles_in[t].vertice_id[ee]; const vrtx_id_t v1 = triangles_in[t].vertice_id[(ee + 1) % 3]; /* Process only "my" edges! */ - const int h = + const int64_t h = v0 + v1 + MMIN(v0, v1); /* v0,v1 and v1,v0 must give the same hash!!! */ if(h % thread_count != rank) continue; /* Create edge. */ @@ -759,7 +735,7 @@ collect_and_link_neighbours * and connect neighbours. */ tmp = darray_neighbourhood_size_get(&neighbourhood_by_edge); ASSERT(tmp <= EDGE_MAX__); - edge_count = (edge_id_t) tmp; + edge_count = (edge_id_t)tmp; FOR_EACH(e, 0, edge_count) { double edge[3], common_edge[3], basis[9], norm, mz, mz_edge; vrtx_id_t v0, v1, v2; @@ -767,10 +743,10 @@ collect_and_link_neighbours = darray_neighbourhood_data_get(&neighbourhood_by_edge) + e; struct darray_neighbour* neighbour_list = &neighbourhood->neighbours; side_id_t i, neighbour_count; - char mz_vid, mz_vid_edge; + unsigned char mz_vid, mz_vid_edge; tmp = darray_neighbour_size_get(neighbour_list); ASSERT(tmp <= SIDE_MAX__); - neighbour_count = (side_id_t) tmp; + neighbour_count = (side_id_t)tmp; ASSERT(neighbour_count); v0 = neighbourhood->edge.vrtx0; v1 = neighbourhood->edge.vrtx1; @@ -783,7 +759,7 @@ collect_and_link_neighbours mz_vid_edge = 1; } norm = d3_normalize(common_edge, common_edge); - ASSERT(norm); + ASSERT(norm); (void)norm; d33_basis(basis, common_edge); d33_inverse(basis, basis); FOR_EACH(i, 0, neighbour_count) { @@ -791,7 +767,7 @@ collect_and_link_neighbours = darray_neighbour_data_get(neighbour_list) + i; const struct triangle_in *trg_in = triangles_in + neighbour_info->trg_id; struct triangle_tmp *neighbour = triangles_tmp + neighbour_info->trg_id; - char actual_vid; + unsigned char actual_vid; v2 = trg_in->vertice_id[(neighbour_info->edge_rank + 2) % 3]; if(vertices[v2].pos.z > mz_edge) { mz = vertices[v2].pos.z; @@ -801,17 +777,18 @@ collect_and_link_neighbours mz_vid = mz_vid_edge; } /* Compute the actual vertex id - * as vertices are not in the v0 v1 v2 order in the actual triangle */ + * as vertices are not in the v0 v1 v2 order in the actual triangle */ if(mz_vid == 2) { - actual_vid = (2 + neighbour_info->edge_rank) % 3; + actual_vid = (unsigned char)((2 + neighbour_info->edge_rank) % 3); } else { int is_r = neighbour->reversed_edge[neighbour_info->edge_rank]; ASSERT(mz_vid == 0 || mz_vid == 1); - actual_vid = ((is_r ? 1 - mz_vid : mz_vid) + neighbour_info->edge_rank) % 3; + actual_vid = (unsigned char)(((is_r ? 1 - mz_vid : mz_vid) + + neighbour_info->edge_rank) % 3); } ASSERT(neighbour->max_z <= mz); neighbour->max_z = mz; - ASSERT(0 <= actual_vid && actual_vid <= 2); + ASSERT(actual_vid <= 2); neighbour->max_z_vrtx_id = actual_vid; /* Compute rotation angle around common edge */ d3_sub(edge, vertices[v2].vec, vertices[v0].vec); @@ -834,8 +811,8 @@ collect_and_link_neighbours const struct neighbour_info* ccw_neighbour = darray_neighbour_cdata_get(neighbour_list) + (i + 1) % neighbour_count; /* Rank of the edge of interest in triangles */ - const char crt_edge = current->edge_rank; - const char ccw_edge = ccw_neighbour->edge_rank; + const unsigned char crt_edge = current->edge_rank; + const unsigned char ccw_edge = ccw_neighbour->edge_rank; /* User id of current triangles */ const trg_id_t crt_id = current->trg_id; const trg_id_t ccw_id = ccw_neighbour->trg_id; @@ -895,12 +872,12 @@ build_result triangles_enc = darray_triangle_enc_data_get(&desc->triangles_enc); #pragma omp parallel for schedule(dynamic) - FOR_EACH(ee, 0, desc->enclosures_count) { + for(ee = 0; ee < (int)desc->enclosures_count; ee++) { const enclosure_id_t e = (enclosure_id_t)ee; struct enclosure_data* enc = enclosures + e; struct cc_descriptor* current = cc_descriptors + enc->first_component; struct htable_vrtx_id vtable; - const char* enc_memb; + const unsigned char* enc_memb; uint64_t* enc_memb64; trg_id_t fst_ixd = 0; trg_id_t sgd_ixd; @@ -918,14 +895,14 @@ build_result = (unsigned)current->medium; /* Back to API type */ ASSERT(enc->header.enclosed_medium < desc->scene->nmeds); /* Get side_membership information from component */ - darray_char_copy_and_clear(&enc->side_membership, + darray_uchar_copy_and_clear(&enc->side_membership, &current->side_membership); - enc_memb = darray_char_cdata_get(&enc->side_membership); + enc_memb = darray_uchar_cdata_get(&enc->side_membership); enc_memb64 = (uint64_t*)enc_memb; /* Check we can use uint64_t to merge bit vectors */ - ASSERT(darray_char_capacity(&enc->side_membership) + ASSERT(darray_uchar_capacity(&enc->side_membership) % sizeof(*enc_memb64) == 0); - ASSERT(darray_char_size_get(&enc->side_membership) + ASSERT(darray_uchar_size_get(&enc->side_membership) == desc->scene->nutris); /* Merge enclosures' membership information */ for(;;) { @@ -938,20 +915,20 @@ build_result last_member_id = MMAX(last_member_id, current->last_member_id); if(current->next_component == COMPONENT_NULL__) break; current = cc_descriptors + current->next_component; - ASSERT(darray_char_size_get(&current->side_membership) + ASSERT(darray_uchar_size_get(&current->side_membership) == desc->scene->nutris); - cc_memb64 = (uint64_t*)darray_char_cdata_get(&current->side_membership); - ASSERT(darray_char_capacity(&current->side_membership) + cc_memb64 = (uint64_t*)darray_uchar_cdata_get(&current->side_membership); + ASSERT(darray_uchar_capacity(&current->side_membership) % sizeof(*enc_memb64) == 0); tmin = current->first_member_id / sizeof(*enc_memb64); - tmax = (darray_char_size_get(&enc->side_membership) + sizeof(*enc_memb64) - 1) + tmax = (darray_uchar_size_get(&enc->side_membership) + sizeof(*enc_memb64) - 1) / sizeof(*enc_memb64); ASSERT(tmin <= tmax); ASSERT(tmax < TRG_MAX__); tmin64 = (trg_id_t)tmin; tmax64 = (trg_id_t)tmax; FOR_EACH(t64, tmin64, tmax64) enc_memb64[t64] |= cc_memb64[t64]; - darray_char_purge(&current->side_membership); + darray_uchar_purge(&current->side_membership); } /* Translate membership into side and vertex lists. */ @@ -1013,7 +990,7 @@ build_result } if(fst_ixd == sgd_ixd) break; } - darray_char_purge(&enc->side_membership); + darray_uchar_purge(&enc->side_membership); htable_vrtx_id_release(&vtable); error: /* Cannot exit openmp block */ diff --git a/src/senc_scene_analyze_c.h b/src/senc_scene_analyze_c.h @@ -20,7 +20,7 @@ #include "senc_internal_types.h" #include <rsys/mem_allocator.h> -#include <rsys/dynamic_array_char.h> +#include <rsys/dynamic_array_uchar.h> #include <rsys/hash_table.h> #include <rsys/double3.h> @@ -63,10 +63,10 @@ set_edge (const vrtx_id_t vrtx0, const vrtx_id_t vrtx1, struct trg_edge* edge, - char* reversed) + unsigned char* reversed) { ASSERT(edge && reversed && vrtx0 != vrtx1); - ASSERT(*reversed == CHAR_MAX); /* Should not be already set. */ + ASSERT(*reversed == UCHAR_MAX); /* Should not be already set. */ if(vrtx0 < vrtx1) { edge->vrtx0 = vrtx0; edge->vrtx1 = vrtx1; @@ -93,7 +93,7 @@ struct trgside { /* Id of this trgside's medium */ medium_id_t medium; /* The list containing the trgside; made of enum list_id flags */ - char list_id; + unsigned char list_id; /* Implicit information that we don't need to store: * - triangle_id @@ -155,10 +155,10 @@ struct cc_descriptor { enclosure_id_t enclosure_id; /* To create linked lists of componnents */ component_id_t next_component; - /* TODO: use only 1 bit per side (now 1 char per triangle) */ + /* TODO: use only 1 bit per side (now 1 uchar per triangle) */ trg_id_t first_member_id; trg_id_t last_member_id; - struct darray_char side_membership; + struct darray_uchar side_membership; }; extern const struct cc_descriptor CC_DESCRIPTOR_NULL; @@ -168,14 +168,14 @@ cc_descriptor_init(struct mem_allocator* alloc, struct cc_descriptor* data) { ASSERT(data); *data = CC_DESCRIPTOR_NULL; - darray_char_init(alloc, &data->side_membership); + darray_uchar_init(alloc, &data->side_membership); } static FINLINE void cc_descriptor_release(struct cc_descriptor* data) { ASSERT(data); - darray_char_release(&data->side_membership); + darray_uchar_release(&data->side_membership); } static FINLINE void @@ -193,7 +193,7 @@ cc_descriptor_clear(struct cc_descriptor* data) data->next_component = CC_DESCRIPTOR_NULL.next_component; data->first_member_id = CC_DESCRIPTOR_NULL.first_member_id; data->last_member_id = CC_DESCRIPTOR_NULL.last_member_id; - darray_char_clear(&data->side_membership); + darray_uchar_clear(&data->side_membership); } static FINLINE void @@ -211,7 +211,7 @@ cc_descriptor_purge(struct cc_descriptor* data) data->next_component = CC_DESCRIPTOR_NULL.next_component; data->first_member_id = CC_DESCRIPTOR_NULL.first_member_id; data->last_member_id = CC_DESCRIPTOR_NULL.last_member_id; - darray_char_purge(&data->side_membership); + darray_uchar_purge(&data->side_membership); } static FINLINE res_T @@ -230,7 +230,7 @@ cc_descriptor_copy(struct cc_descriptor* dst, const struct cc_descriptor* src) dst->next_component = src->next_component; dst->first_member_id = src->first_member_id; dst->last_member_id = src->last_member_id; - return darray_char_copy(&dst->side_membership, &src->side_membership); + return darray_uchar_copy(&dst->side_membership, &src->side_membership); } static FINLINE res_T @@ -251,7 +251,7 @@ cc_descriptor_copy_and_release dst->next_component = src->next_component; dst->first_member_id = src->first_member_id; dst->last_member_id = src->last_member_id; - return darray_char_copy_and_release(&dst->side_membership, + return darray_uchar_copy_and_release(&dst->side_membership, &src->side_membership); } @@ -274,9 +274,9 @@ cc_descriptor_copy_and_release struct triangle_tmp { /* Are the edges of the triangle defined in the same order than * the edges they are linked to? */ - char reversed_edge[3]; + unsigned char reversed_edge[3]; /* tmp data used to find the +Z-most vertex of components */ - char max_z_vrtx_id; + unsigned char max_z_vrtx_id; double max_z; }; @@ -286,8 +286,8 @@ triangle_tmp_init(struct mem_allocator* alloc, struct triangle_tmp* trg) { int i; (void) alloc; ASSERT(trg); - FOR_EACH(i, 0, 3) trg->reversed_edge[i] = CHAR_MAX; - trg->max_z_vrtx_id = CHAR_MAX; + FOR_EACH(i, 0, 3) trg->reversed_edge[i] = UCHAR_MAX; + trg->max_z_vrtx_id = UCHAR_MAX; trg->max_z = -DBL_MAX; } #define DARRAY_FUNCTOR_INIT triangle_tmp_init @@ -306,7 +306,7 @@ triangle_tmp_init(struct mem_allocator* alloc, struct triangle_tmp* trg) { struct neighbour_info { trg_id_t trg_id; /* Rank of the edge in the triangle (in [0 2]) */ - char edge_rank; + unsigned char edge_rank; double angle; }; #define DARRAY_NAME neighbour diff --git a/src/senc_scene_c.h b/src/senc_scene_c.h @@ -82,7 +82,7 @@ triangle_in_flip(struct triangle_in* trg) { trg->medium[1] = m; } -FINLINE int +static FINLINE int vrtx_eq(const union double3* v1, const union double3* v2) { ASSERT(v1 && v2); @@ -108,7 +108,7 @@ union vrtx_id3 { vrtx_id_t vec[3]; }; -FINLINE char /* Return 1 if reversed */ +static FINLINE char /* Return 1 if reversed */ trg_make_key(union vrtx_id3* k, const vrtx_id_t t[3]) { ASSERT(t); @@ -162,7 +162,7 @@ trg_make_key(union vrtx_id3* k, const vrtx_id_t t[3]) } } -FINLINE int +static FINLINE int trg_key_eq(const union vrtx_id3* k1, const union vrtx_id3* k2) { ASSERT(k1 && k2); diff --git a/src/test_senc_cube_on_cube.c b/src/test_senc_cube_on_cube.c @@ -13,11 +13,15 @@ * You should have received a copy of the GNU General Public License * along with this program. If not, see <http://www.gnu.org/licenses/>. */ +#define _POSIX_C_SOURCE 200112L /* snprintf */ + #include "senc.h" #include "test_senc_utils.h" #include <rsys/double3.h> +#include <stdio.h> + /* Z ^ diff --git a/src/test_senc_many_triangles.c b/src/test_senc_many_triangles.c @@ -17,11 +17,11 @@ #include "test_senc_utils.h" #include <star/s3dut.h> - #include <rsys/double3.h> - #include <rsys/clock_time.h> +#include <limits.h> + struct s3dut_context { struct s3dut_mesh_data data; struct context ctx; @@ -127,7 +127,7 @@ main(int argc, char** argv) FOR_EACH(i, 0, count) { struct senc_enclosure* enclosure; - struct enclosure_header* header; + const struct enclosure_header* header; CHK(senc_descriptor_get_enclosure(desc, i, &enclosure) == RES_OK); CHK(senc_enclosure_get_header(enclosure, &header) == RES_OK); CHK(header->triangle_count == diff --git a/src/test_senc_utils.h b/src/test_senc_utils.h @@ -77,7 +77,7 @@ static const unsigned medium1_front0[12] = { 0, 0, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1 static const unsigned gid_face[12] = { 0, 0, 1, 1, 2, 2, 3, 3, 4, 4, 5, 5 }; -static void +static INLINE void get_indices(const unsigned itri, unsigned ids[3], void* context) { struct context* ctx = context; @@ -87,7 +87,7 @@ get_indices(const unsigned itri, unsigned ids[3], void* context) ids[ctx->reverse_vrtx ? 1 : 2] = ctx->indices[itri * 3 + 2]; } -static void +static INLINE void get_position(const unsigned ivert, double pos[3], void* context) { struct context* ctx = context; @@ -97,7 +97,7 @@ get_position(const unsigned ivert, double pos[3], void* context) pos[2] = ctx->positions[ivert * 3 + 2] * ctx->scale + ctx->offset[2]; } -static void +static INLINE void get_media(const unsigned itri, unsigned medium[2], void* context) { struct context* ctx = context; @@ -106,7 +106,7 @@ get_media(const unsigned itri, unsigned medium[2], void* context) medium[ctx->reverse_med ? 0 : 1] = ctx->back_media[itri]; } -static void +static INLINE void get_global_id(const unsigned itri, unsigned* gid, void* context) { struct context* ctx = context; @@ -139,7 +139,7 @@ dump_mesh } } -static void +static INLINE void dump_enclosure (struct senc_descriptor* desc, const unsigned enc, @@ -167,8 +167,8 @@ dump_enclosure FOR_EACH(i, 0, header->triangle_count) { unsigned indices[3]; CHK(senc_enclosure_get_triangle(enclosure, i, indices) == RES_OK); - fprintf(stream, "f %lu %lu %lu\n", - 1+indices[0], 1+indices[1], 1+indices[2]); + fprintf(stream, "f %lu %lu %lu\n", (unsigned long)(1+indices[0]), + (unsigned long)(1+indices[1]), (unsigned long)(1+indices[2])); } CHK(senc_enclosure_ref_put(enclosure) == RES_OK); fclose(stream);