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 bd25dbc2c166e41fbd7ce85dd1d2f7c6507b9fd7
parent 703a37a89bbef21f45329f06217c83e654c335b5
Author: Christophe Coustet <christophe.coustet@meso-star.com>
Date:   Wed,  1 Aug 2018 15:48:31 +0200

Seems to work, tests are OK; need to be checked again

Diffstat:
Msrc/senc2d_scene_analyze.c | 159+++++++++++++++++++++++++++++++++++++++----------------------------------------
Msrc/test_senc2d_utils.h | 2+-
2 files changed, 79 insertions(+), 82 deletions(-)

diff --git a/src/senc2d_scene_analyze.c b/src/senc2d_scene_analyze.c @@ -67,9 +67,9 @@ get_side_not_in_connex_component { side_id_t i = *first_side_not_in_component; while(i <= last_side - && (segsides[i].medium != medium || processed[i])) { + && (segsides[i].medium != medium + || (processed[SEGSIDE_2_SEG(i)] & SEGSIDE_2_SIDEFLAG(i)))) ++i; - } *first_side_not_in_component = i + 1; if(i > last_side) return SIDE_NULL__; return i; @@ -134,7 +134,7 @@ extract_connex_components ATOMIC* component_count, /* Shared error status. * We accept to overwritte an error with a different error */ - res_T* res) + res_T* p_res) { /* This function is called from an omp parallel block and executed * concurrently. */ @@ -143,6 +143,7 @@ extract_connex_components int64_t mm; struct darray_side_id stack; struct darray_side_id ids_of_sides_around_max_y_vertex; + struct htable_media tmp_media; const union double2* positions; const struct segment_tmp* segments_tmp; struct segment_comp* segments_comp; @@ -151,10 +152,9 @@ extract_connex_components /* An array to store the component being processed */ struct darray_side_id current_component; size_t ii, sz; - res_T res2 = RES_OK; ASSERT(segsides && desc && connex_components && segments_tmp_array - && segments_comp_array && s2d_view && component_count && res); + && segments_comp_array && s2d_view && component_count && p_res); alloc = descriptor_get_allocator(desc); scn = desc->scene; positions = darray_position_cdata_get(&scn->vertices); @@ -163,9 +163,10 @@ extract_connex_components darray_side_id_init(alloc, &stack); darray_side_id_init(alloc, &ids_of_sides_around_max_y_vertex); darray_side_id_init(alloc, &current_component); + htable_media_init(alloc, &tmp_media); processed = calloc(scn->nusegs, sizeof(unsigned char)); - if (!processed) { - *res = RES_MEM_ERR; + if(!processed) { + *p_res = RES_MEM_ERR; return; } @@ -191,8 +192,8 @@ extract_connex_components /* We loop on sides to build connex components. */ #pragma omp for schedule(dynamic) nowait - for (mm = 0; mm < (int64_t) scn->nmeds; mm++) { /* Process all media */ - const medium_id_t m = (medium_id_t) mm; + for(mm = 0; mm < (int64_t) scn->nmeds; mm++) { /* Process all media */ + const medium_id_t m = (medium_id_t)mm; /* Any not-already-used side is used as a starting point */ const struct side_range* media_use = darray_side_range_cdata_get(&scn->media_use) + m; @@ -200,17 +201,19 @@ extract_connex_components double max_y_ny; const seg_id_t last_side = media_use->last; char one = 1; + int component_canceled = 0; + res_T tmp_res = RES_OK; ATOMIC id; - if (*res != RES_OK) continue; - if (first_side_not_in_component == SIDE_NULL__) + if(*p_res != RES_OK) continue; + if(first_side_not_in_component == SIDE_NULL__) continue; /* Unused medium */ ASSERT(first_side_not_in_component < 2 * scn->nusegs); ASSERT(darray_side_id_size_get(&stack) == 0); ASSERT(darray_side_id_size_get(&current_component) == 0); - for (;;) { /* Process all components for this medium */ + for(;;) { /* Process all components for this medium */ const side_id_t start_side_id = get_side_not_in_connex_component - (last_side, segsides, processed, &first_side_not_in_component, m); + (last_side, segsides, processed, &first_side_not_in_component, m); side_id_t crt_side_id = start_side_id; side_id_t last_side_id = start_side_id; vrtx_id_t max_y_vrtx_id = VRTX_NULL__; @@ -218,22 +221,25 @@ extract_connex_components double max_y = -DBL_MAX; ASSERT(start_side_id == SIDE_NULL__ || start_side_id < 2 * scn->nusegs); + darray_side_id_clear(&current_component); + htable_media_clear(&tmp_media); - if (*res != RES_OK) break; - if (start_side_id == SIDE_NULL__) + if(*p_res != RES_OK) break; + if(start_side_id == SIDE_NULL__) break; /* start_side_id=SIDE_NULL__ => component done! */ #ifndef NDEBUG { - seg_id_t tid = SEGSIDE_2_SEG(start_side_id); + seg_id_t sid = SEGSIDE_2_SEG(start_side_id); enum side_id s = SEGSIDE_2_SIDE(start_side_id); medium_id_t side_med - = darray_segment_in_data_get(&desc->scene->segments_in)[tid].medium[s]; + = darray_segment_in_data_get(&desc->scene->segments_in)[sid].medium[s]; ASSERT(side_med == m); } #endif - for (;;) { /* Process all sides of this component */ + OK2(htable_media_set(&tmp_media, &m, &one)); + for(;;) { /* Process all sides of this component */ int i; enum side_flag crt_side_flag = SEGSIDE_2_SIDEFLAG(crt_side_id); struct segside* crt_side = segsides + crt_side_id; @@ -242,55 +248,47 @@ extract_connex_components darray_segment_in_cdata_get(&scn->segments_in) + crt_seg_id; unsigned char* seg_used = processed + crt_seg_id; const struct segment_tmp* const seg_tmp = segments_tmp + crt_seg_id; - res_T tmp_res; ASSERT(crt_seg_id < scn->nusegs); - if (*res != RES_OK) break; - + if(*p_res != RES_OK) break; /* Record Ymax information * Keep track of the appropriate vertex of the component in order * to cast a ray at the component grouping step of the algorithm. * The most appropriate vertex is (the) one with the greater Y * coordinate. */ - if (max_y < seg_tmp->max_y) { + if(max_y < seg_tmp->max_y) { /* New best vertex */ max_y = seg_tmp->max_y; /* New vertex: reset list of sides */ darray_side_id_clear(&ids_of_sides_around_max_y_vertex); /* Select a vertex with y == max_y */ - if (max_y == positions[seg_in->vertice_id[0]].pos.y) { + if(max_y == positions[seg_in->vertice_id[0]].pos.y) { max_y_vrtx_id = seg_in->vertice_id[0]; - } - else { + } else { ASSERT(max_y == positions[seg_in->vertice_id[1]].pos.y); max_y_vrtx_id = seg_in->vertice_id[1]; } /* List of sides using the vertex */ - tmp_res = darray_side_id_push_back(&ids_of_sides_around_max_y_vertex, - &crt_side_id); - if (tmp_res != RES_OK) *res = tmp_res; - } - else if (max_y == seg_tmp->max_y) { + OK2(darray_side_id_push_back(&ids_of_sides_around_max_y_vertex, + &crt_side_id)); + } else if(max_y == seg_tmp->max_y) { /* Does this segment use the currently selected max_y vertex? */ FOR_EACH(i, 0, 2) { - if (max_y_vrtx_id == seg_in->vertice_id[i]) { + if(max_y_vrtx_id == seg_in->vertice_id[i]) { /* List of sides using the vertex */ - tmp_res = darray_side_id_push_back(&ids_of_sides_around_max_y_vertex, - &crt_side_id); - if (tmp_res != RES_OK) *res = tmp_res; + OK2(darray_side_id_push_back(&ids_of_sides_around_max_y_vertex, + &crt_side_id)); break; } } } - if (*res != RES_OK) break; /* Record crt_side both as component and segment level */ - ASSERT((*seg_used & crt_side_flag) == 0); - *seg_used |= crt_side_flag; - tmp_res = darray_side_id_push_back(&current_component, &crt_side_id); - if (tmp_res != RES_OK) *res = tmp_res; - if (*res != RES_OK) break; + if((*seg_used & crt_side_flag) == 0) { + OK2(darray_side_id_push_back(&current_component, &crt_side_id)); + *seg_used |= crt_side_flag; + } /* Store neighbours' sides in a stack */ FOR_EACH(i, 0, 2) { @@ -299,36 +297,35 @@ extract_connex_components enum side_flag nbour_side_id = SEGSIDE_2_SIDEFLAG(neighbour_id); unsigned char* nbour_used = processed + nbour_seg_id; const struct segside* neighbour = segsides + neighbour_id; - if (neighbour->medium < m) { + if(*nbour_used & nbour_side_id) continue; /* Already processed */ + if(neighbour->medium < m) { /* Not the same medium. * Neighbour's medium id is less than current medium: the whole * component is to be processed by another thread (possibly the one * associated with neighbour's medium). */ - darray_side_id_clear(&current_component); - - - continue; + component_canceled = 1; + goto canceled; } - if (*nbour_used & nbour_side_id) continue; /* Already processed */ /* Mark neighbour as processed and stack it */ *nbour_used |= nbour_side_id; - tmp_res = darray_side_id_push_back(&stack, &neighbour_id); - if (tmp_res != RES_OK) *res = tmp_res; - if (*res != RES_OK) break; + OK2(darray_side_id_push_back(&stack, &neighbour_id)); + OK2(darray_side_id_push_back(&current_component, &neighbour_id)); + if(neighbour->medium != m) /* Can be a new medium */ + htable_media_set(&tmp_media, &neighbour->medium, &one); } - if (*res != RES_OK) break; sz = darray_side_id_size_get(&stack); - if (sz == 0) break; /* Empty stack => component is done! */ + if(sz == 0) break; /* Empty stack => component is done! */ crt_side_id = darray_side_id_cdata_get(&stack)[sz - 1]; darray_side_id_pop_back(&stack); last_side_id = MMAX(last_side_id, crt_side_id); } - if (*res != RES_OK) break; +canceled: + if(component_canceled) continue; /* Register the new component and get it initialized */ cc = MEM_ALLOC(alloc, sizeof(struct cc_descriptor)); - if (!cc) *res = RES_MEM_ERR; - if (*res != RES_OK) break; + if(!cc) *p_res = RES_MEM_ERR; + if(*p_res != RES_OK) break; ASSERT(m == segsides[start_side_id].medium); ASSERT(max_y_vrtx_id != VRTX_NULL__); @@ -339,21 +336,20 @@ extract_connex_components sz = darray_side_id_size_get(&current_component); ASSERT(sz <= SIDE_MAX__); cc->side_count = (side_id_t)sz; - cc->side_range.last = last_side_id; cc->side_range.first = start_side_id; + cc->side_range.last = last_side_id; cc->max_y_vrtx_id = max_y_vrtx_id; /* First medium of the component */ - htable_media_set(&cc->media, &m, &one); /* FIX FIX FIX */ + htable_media_copy_and_clear(&cc->media, &tmp_media); /* Write component membership in the global structure * No need for sync here as an unique thread write a given side */ - ASSERT(sizeof(cc->cc_id) >= 4); /* Atomic write */ - sz = darray_side_id_size_get(&current_component); + static_assert(sizeof(cc->cc_id) >= 4, "Cannot write IDs sync-free"); + ASSERT(IS_ALIGNED(segments_comp->component, sizeof(cc->cc_id))); FOR_EACH(ii, 0, sz) { const side_id_t s = darray_side_id_cdata_get(&current_component)[ii]; segments_comp[SEGSIDE_2_SEG(s)].component[SEGSIDE_2_SIDE(s)] = cc->cc_id; } - darray_side_id_clear(&current_component); /* Compute the normal at the max_y vertex. */ max_y_ny = 0; @@ -374,24 +370,22 @@ extract_connex_components * contribute (remember than x + y - x - y can be non-zero). */ ASSERT(seg_comp->component[SIDE_FRONT] == cc->cc_id || seg_comp->component[SIDE_BACK] == cc->cc_id); - if (seg_comp->component[SIDE_FRONT] == seg_comp->component[SIDE_BACK]) + if(seg_comp->component[SIDE_FRONT] == seg_comp->component[SIDE_BACK]) continue; d2_sub(edge, vertices[seg_in->vertice_id[1]].vec, vertices[seg_in->vertice_id[0]].vec); d2(normal, edge[1], -edge[0]); norm = d2_normalize(normal, normal); - ASSERT(norm); (void) norm; + ASSERT(norm); (void)norm; /* Geometrical normal points toward the back side */ - if (SEGSIDE_IS_FRONT(side_id)) { + if(SEGSIDE_IS_FRONT(side_id)) { max_y_ny -= normal[1]; - } - else { + } else { max_y_ny += normal[1]; } } - if (*res != RES_OK) break; cc->is_outer_border = (max_y_ny < 0); /* Need to synchronize connex_components growth as this global structure @@ -400,12 +394,12 @@ extract_connex_components { struct cc_descriptor** components; sz = darray_ptr_component_descriptor_size_get(connex_components); - if (sz <= cc->cc_id) { - res_T tmp_res = darray_ptr_component_descriptor_resize(connex_components, + if(sz <= cc->cc_id) { + tmp_res = darray_ptr_component_descriptor_resize(connex_components, 1 + cc->cc_id); - if (tmp_res != RES_OK) *res = tmp_res; + if(tmp_res != RES_OK) *p_res = tmp_res; } - if (*res == RES_OK) { + if(*p_res == RES_OK) { /* Don't set the pointer before resize as this can lead to move data */ components = darray_ptr_component_descriptor_data_get(connex_components); @@ -414,17 +408,20 @@ extract_connex_components } } } + tmp_error: + if(tmp_res != RES_OK) *p_res = tmp_res; } /* No barrier here */ free(processed); darray_side_id_release(&current_component); darray_side_id_release(&stack); darray_side_id_release(&ids_of_sides_around_max_y_vertex); + htable_media_release(&tmp_media); /* The first thread here creates the s2d view */ #pragma omp single nowait - if(*res == RES_OK) { - res_T tmp_res = RES_OK; + if(*p_res == RES_OK) { + res_T res = RES_OK; struct s2d_device* s2d = NULL; struct s2d_scene* s2d_scn = NULL; struct s2d_shape* s2d_shp = NULL; @@ -435,22 +432,22 @@ extract_connex_components attribs.get = get_scn_position; /* Put geometry in a 2D view */ - OK2(s2d_device_create(desc->scene->dev->logger, alloc, 0, &s2d)); - OK2(s2d_scene_create(s2d, &s2d_scn)); - OK2(s2d_shape_create_line_segments(s2d, &s2d_shp)); + OK(s2d_device_create(desc->scene->dev->logger, alloc, 0, &s2d)); + OK(s2d_scene_create(s2d, &s2d_scn)); + OK(s2d_shape_create_line_segments(s2d, &s2d_shp)); /* Back to API type for ntris and nverts */ ASSERT(desc->scene->nusegs < UINT_MAX); ASSERT(desc->scene->nuverts < UINT_MAX); - OK2(s2d_line_segments_setup_indexed_vertices(s2d_shp, + OK(s2d_line_segments_setup_indexed_vertices(s2d_shp, (unsigned)desc->scene->nusegs, get_scn_indices, (unsigned)desc->scene->nuverts, &attribs, 1, desc->scene)); s2d_line_segments_set_hit_filter_function(s2d_shp, self_hit_filter, - darray_segment_comp_data_get(segments_comp_array)); - OK2(s2d_scene_attach_shape(s2d_scn, s2d_shp)); - OK2(s2d_scene_view_create(s2d_scn, S2D_TRACE, s2d_view)); - tmp_error: - if(tmp_res != RES_OK) *res = tmp_res; + segments_comp_array); + OK(s2d_scene_attach_shape(s2d_scn, s2d_shp)); + OK(s2d_scene_view_create(s2d_scn, S2D_TRACE, s2d_view)); + error: + if(res != RES_OK) *p_res = res; if(s2d) S2D(device_ref_put(s2d)); if(s2d_scn) S2D(scene_ref_put(s2d_scn)); if(s2d_shp) S2D(shape_ref_put(s2d_shp)); @@ -459,7 +456,7 @@ extract_connex_components #ifndef NDEBUG /* Need to wait for all threads done to be able to check stuff */ #pragma omp barrier - if(*res != RES_OK) return; + if(*p_res != RES_OK) return; #pragma omp single { seg_id_t s_; diff --git a/src/test_senc2d_utils.h b/src/test_senc2d_utils.h @@ -246,7 +246,7 @@ static INLINE void check_desc(struct senc2d_descriptor* desc) size_t e_cpt = 0; CHK(senc2d_descriptor_get_max_medium(desc, &maxm) == RES_OK); CHK(senc2d_descriptor_get_enclosure_count(desc, &ecount) == RES_OK); - for (i = 0; i <= maxm; i++) { + for(i = 0; i <= maxm; i++) { unsigned j, ecount_bym; unsigned found = 0; CHK(senc2d_descriptor_get_enclosure_count_by_medium(desc, i, &ecount_bym) == RES_OK);