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 b053263689c5644d5e475a5a0bda4ce550b72a1d
parent 702181ede1b14fdf97e3755a31af4c97cec5e94e
Author: Christophe Coustet <christophe.coustet@meso-star.com>
Date:   Fri, 20 Apr 2018 14:25:44 +0200

BugFix: fix interior/exterior detection for components.

Diffstat:
Msrc/senc2d_descriptor_c.h | 4+---
Msrc/senc2d_scene_analyze.c | 269++++++++++++++++++++++++++++++++++++++-----------------------------------------
Msrc/senc2d_scene_analyze_c.h | 11++---------
3 files changed, 133 insertions(+), 151 deletions(-)

diff --git a/src/senc2d_descriptor_c.h b/src/senc2d_descriptor_c.h @@ -32,7 +32,6 @@ struct segment_comp { component_id_t component[2]; }; -#ifndef NDEBUG static void segment_comp_init(struct mem_allocator* alloc, struct segment_comp* seg) { int i; @@ -40,11 +39,10 @@ segment_comp_init(struct mem_allocator* alloc, struct segment_comp* seg) { ASSERT(seg); FOR_EACH(i, 0, 2) seg->component[i] = COMPONENT_NULL__; } -#define DARRAY_FUNCTOR_INIT segment_comp_init -#endif #define DARRAY_NAME segment_comp #define DARRAY_DATA struct segment_comp +#define DARRAY_FUNCTOR_INIT segment_comp_init #include <rsys/dynamic_array.h> struct segment_enc { diff --git a/src/senc2d_scene_analyze.c b/src/senc2d_scene_analyze.c @@ -43,8 +43,8 @@ #include <stdlib.h> #define CC_DESCRIPTOR_NULL__ {\ - {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,\ + CHAR_MAX, VRTX_NULL__, 0, MEDIUM_NULL__,\ + CC_ID_NONE, CC_GROUP_ROOT_NONE, ENCLOSURE_NULL__,\ { SEG_NULL__, 0}\ } const struct cc_descriptor CC_DESCRIPTOR_NULL = CC_DESCRIPTOR_NULL__; @@ -90,7 +90,7 @@ neighbour_cmp(const void* w1, const void* w2) return (a1 > a2) - (a1 < a2); } -static FINLINE void +static FINLINE res_T add_side_to_stack (const struct senc2d_scene* scn, struct darray_side_id* stack, @@ -102,8 +102,8 @@ add_side_to_stack && side_id < SIDE_MAX__ && side_id < 2 * scn->nusegs); ASSERT((darray_side_id_size_get(stack) > 128) || !find_side_in_list(segsides, stack, side_id, FLAG_WAITING_STACK)); - darray_side_id_push_back(stack, &side_id); segsides[side_id].list_id = FLAG_WAITING_STACK; + return darray_side_id_push_back(stack, &side_id); } static side_id_t @@ -180,7 +180,7 @@ self_hit_filter enum side_id hit_side; component_id_t hit_component; - (void) ray_org; (void) ray_dir; + (void)ray_org; (void)ray_dir; ASSERT(hit && segments_comp && origin_component); ASSERT(hit->prim.prim_id < darray_segment_comp_size_get(segments_comp)); hit_seg_comp = darray_segment_comp_cdata_get(segments_comp) @@ -212,6 +212,9 @@ extract_connex_components struct mem_allocator* alloc; int64_t mm; struct darray_side_id stack; + struct darray_side_id ids_of_sides_around_max_y_vertex; + const union double2* positions; + size_t sz, ii; #ifndef NDEBUG seg_id_t s_; component_id_t c; @@ -221,6 +224,9 @@ extract_connex_components && s2d_view && component_count && res); alloc = descriptor_get_allocator(desc); scn = desc->scene; + positions = darray_position_cdata_get(&scn->vertices); + darray_side_id_init(alloc, &stack); + darray_side_id_init(alloc, &ids_of_sides_around_max_y_vertex); #ifndef NDEBUG #pragma omp single @@ -241,14 +247,13 @@ extract_connex_components } /* Implicit barrier here */ #endif - darray_side_id_init(alloc, &stack); - #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; struct cc_descriptor* cc; /* Any not-already-used side is used as a starting point */ side_id_t first_side_not_in_component; + double max_y_ny; if(*res != RES_OK) continue; first_side_not_in_component @@ -262,7 +267,11 @@ extract_connex_components segsides, &first_side_not_in_component, m); side_id_t crt_side_id = start_side_id; side_id_t last_side_id = start_side_id; + double max_y = -DBL_MAX; + ASSERT(start_side_id == SIDE_NULL__ || start_side_id < 2 * scn->nusegs); + + if(*res != RES_OK) break; if(start_side_id == SIDE_NULL__) break; /* start_side_id=SIDE_NULL__ => done! */ ASSERT(segsides[start_side_id].list_id == FLAG_LIST_SIDE_LIST); @@ -300,73 +309,44 @@ extract_connex_components darray_segment_comp_data_get(segments_comp) + crt_seg_id; const struct segment_tmp* const seg_tmp = darray_segment_tmp_cdata_get(segments_tmp_array) + crt_seg_id; - const union double2* vertices = - darray_position_cdata_get(&scn->vertices); ASSERT(crt_seg_id < scn->nusegs); ASSERT(segsides[crt_side_id].medium == m); + if(*res != RES_OK) break; + /* Record Ymax information - * Keep track of the appropriate vertex/side of the connex component - * in order to cast a ray at the component grouping step of the - * algorithm. + * Keep track of the appropriate vertex of the connex 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 more than one vertex/side has the same Y, we want the side that - * most faces Y (that is the one with the greater ny). - * This is mandatory to select the correct side when both sides of a - * segment are candidate. */ - if(cc->max_vrtx[1] <= seg_tmp->max_y) { - /* Can either improve y or ny */ - - if(cc->max_y_side_id == SEGSIDE_OPPOSITE(crt_side_id)) { - /* Both sides are in cc and the opposite side is currently selected - * Just keep the side with ny>0 */ - if(cc->max_y_ny < 0) { - /* Change side! */ - cc->max_y_ny = -cc->max_y_ny; - cc->max_y_side_id = crt_side_id; - } + * coordinate. */ + 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 one vertex with y == max_y */ + if(max_y == positions[seg_in->vertice_id[0]].pos.y) { + cc->max_y_vrtx_id = seg_in->vertice_id[0]; } else { - double edge[2], normal[2], norm, side_ny; - int process = 0; - - 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; - - /* Geometrical normal points toward the right */ - if(SEGSIDE_IS_FRONT(crt_side_id)) { - side_ny = -normal[1]; - process = 1; - } else { - side_ny = normal[1]; - process = 1; - } - - if(process) { - int change = 0; - if(cc->max_vrtx[1] < seg_tmp->max_y) { - change = 1; /* Try first to improve y */ - } else if(cc->max_y_ny <= 0 && fabs(cc->max_y_ny) < fabs(side_ny)) { - change = 1; /* If ny <= 0, the more negative the better */ - } else if(cc->max_y_ny > 0 && cc->max_y_ny < side_ny) { - change = 1; /* If ny > 0, the more positive the better */ - } - - if(change) { - cc->max_y_ny = side_ny; - cc->max_y_side_id = crt_side_id; - ASSERT(seg_tmp->max_y_vrtx_rank < 2); - ASSERT(seg_in->vertice_id[seg_tmp->max_y_vrtx_rank] < scn->nverts); - cc->max_y_vrtx_id = seg_in->vertice_id[seg_tmp->max_y_vrtx_rank]; - ASSERT(seg_tmp->max_y == vertices[cc->max_y_vrtx_id].pos.y); - d2_set(cc->max_vrtx, vertices[cc->max_y_vrtx_id].vec); - } + ASSERT(max_y == positions[seg_in->vertice_id[1]].pos.y); + cc->max_y_vrtx_id = seg_in->vertice_id[1]; + } + /* List of sides using the vertex */ + *res = 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(cc->max_y_vrtx_id == seg_in->vertice_id[i]) { + /* List of sides using the vertex */ + *res = 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 */ cc->side_count++; @@ -387,8 +367,6 @@ extract_connex_components /* Found medium discontinuity! Model topology is broken. */ const struct segment_in* segments_in = darray_segment_in_cdata_get(&scn->segments_in); - const union double2* positions - = darray_position_cdata_get(&scn->vertices); log_err(scn->dev, "Medium mismatch found between neighbour segments %lu %s" " side and %u %s side.\n", @@ -409,7 +387,7 @@ extract_connex_components log_err(desc->scene->dev, "Media: %lu VS %lu\n", (unsigned long)neighbour->medium, (unsigned long)crt_side->medium); *res = RES_BAD_ARG; - continue; + break; } if(neighbour->list_id == FLAG_LIST_COMPONENT) { /* Already processed */ @@ -421,21 +399,64 @@ extract_connex_components if(neighbour->list_id == FLAG_WAITING_STACK) { continue; /* Already processed */ } - add_side_to_stack(scn, &stack, segsides, neighbour_id); + *res = add_side_to_stack(scn, &stack, segsides, neighbour_id); + if(*res != RES_OK) break; } + if(*res != RES_OK) break; if(darray_side_id_size_get(&stack) == 0) break; /* Empty stack => connex component is done! */ crt_side_id = get_side_from_stack(segsides, &stack); last_side_id = MMAX(last_side_id, crt_side_id); } + if(*res != RES_OK) break; /* Keep track of this new connex component */ cc->side_range.last = last_side_id; + + /* Compute the normal at the max_y vertex. */ + max_y_ny = 0; + sz = darray_side_id_size_get(&ids_of_sides_around_max_y_vertex); + FOR_EACH(ii, 0, sz) { + 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); + const struct segment_in* seg_in = + darray_segment_in_cdata_get(&scn->segments_in) + seg_id; + struct segment_comp* seg_comp = + darray_segment_comp_data_get(segments_comp) + seg_id; + const union double2* vertices = + darray_position_cdata_get(&scn->vertices); + double edge[2], normal[2], norm; + + /* To garanty that segments with 2 sides in the component total to 0 + * regardless of numeric accuracy, we need to prevent them to + * 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]) + 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; + + /* Geometrical normal points toward the back side */ + if(SEGSIDE_IS_FRONT(side_id)) { + max_y_ny -= normal[1]; + } 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 * is accessed by multipe threads */ #pragma omp critical { struct cc_descriptor** components; - size_t sz = darray_ptr_component_descriptor_size_get(connex_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, 1 + cc->cc_id); @@ -453,6 +474,7 @@ extract_connex_components } /* No barrier here */ darray_side_id_release(&stack); + darray_side_id_release(&ids_of_sides_around_max_y_vertex); /* The first thread here creates the s2d view */ #pragma omp single nowait @@ -531,6 +553,7 @@ group_connex_components /* This function is called from an omp parallel block and executed * concurrently. */ struct cc_descriptor** descriptors; + const union double2* positions; size_t tmp; component_id_t cc_count; int64_t ccc; @@ -544,6 +567,7 @@ group_connex_components tmp = darray_ptr_component_descriptor_size_get(connex_components); ASSERT(tmp <= COMPONENT_MAX__); cc_count = (component_id_t)tmp; + positions = darray_position_cdata_get(&desc->scene->vertices); /* Cast rays to find links between connex components */ #pragma omp for @@ -555,16 +579,14 @@ group_connex_components const float dir[2] = { 0, 1 }; const float range[2] = { 0, FLT_MAX }; struct cc_descriptor* const cc = descriptors[c]; - const struct segment_comp* origin_seg = - darray_segment_comp_cdata_get(segments_comp) + cc->max_y_vrtx_id; - component_id_t self_hit_component - = origin_seg->component[1 - SEGSIDE_2_SIDE(cc->max_y_side_id)]; + component_id_t self_hit_component = cc->cc_id; + const double* max_vrtx = positions[cc->max_y_vrtx_id].vec; if(*res != RES_OK) continue; ASSERT(cc->cc_id == c); ASSERT(cc->cc_group_root == CC_GROUP_ID_NONE); - if(cc->max_y_ny < 0) { + if(cc->is_outer_border) { int64_t id; /* Don't need to cast a ray */ cc->cc_group_root = cc->cc_id; /* New group with self as root */ @@ -574,15 +596,7 @@ group_connex_components continue; } - ASSERT(cc->max_y_ny != 0 - /* The only situation with ny==0 we can think of is this one: */ - || (segsides[cc->max_y_side_id].medium - == segsides[SEGSIDE_OPPOSITE(cc->max_y_side_id)].medium - && (segsides[cc->max_y_side_id].facing_side_id[0] - == SEGSIDE_OPPOSITE(cc->max_y_side_id) - || segsides[cc->max_y_side_id].facing_side_id[1] - == SEGSIDE_OPPOSITE(cc->max_y_side_id)))); - f2_set_d2(origin, cc->max_vrtx); + f2_set_d2(origin, max_vrtx); /* Self-hit data: self hit if hit this component "on the other side" */ tmp_res = s2d_scene_view_trace_ray(s2d_view, origin, dir, range, &self_hit_component, &hit); @@ -597,34 +611,24 @@ group_connex_components /* Keep track of the first component facing infinity */ ATOMIC_CAS_PTR(infinity_first_cc, cc, NULL); if((*infinity_first_cc)->medium != cc->medium) { - const side_id_t infinity_first_side = (*infinity_first_cc)->max_y_side_id; - const medium_id_t infinity_medium = (*infinity_first_cc)->medium; /* Medium mismatch! Model topology is broken. */ - const seg_id_t t1 = SEGSIDE_2_SEG(infinity_first_side); - const seg_id_t t2 = SEGSIDE_2_SEG(cc->max_y_side_id); - const struct segment_in* segments_in - = darray_segment_in_cdata_get(&desc->scene->segments_in); - const union double2* positions - = darray_position_cdata_get(&desc->scene->vertices); + const double* infinity_max_vrtx = + positions[(*infinity_first_cc)->max_y_vrtx_id].vec; log_err(desc->scene->dev, - "Medium mismatch found between segment %lu %s side and segment" - " %lu %s side, both facing infinity.\n", - (unsigned long)segments_in[t1].global_id, - SEGSIDE_IS_FRONT(infinity_first_side) ? "front" : "back", - (unsigned long)segments_in[t2].global_id, - SEGSIDE_IS_FRONT(cc->max_y_side_id) ? "front" : "back"); + "Medium mismatch found between vertex %lu and vertex %lu," + " both facing infinity.\n", + (unsigned long)(*infinity_first_cc)->max_y_vrtx_id, + (unsigned long)cc->max_y_vrtx_id); log_err(desc->scene->dev, - "Segment %lu:\n (%g %g) (%g %g)\n", - (unsigned long)segments_in[t1].global_id, - SPLIT2(positions[segments_in[t1].vertice_id[0]].vec), - SPLIT2(positions[segments_in[t1].vertice_id[1]].vec)); + "Vertex %lu: (%g %g)\n", + (unsigned long)(*infinity_first_cc)->max_y_vrtx_id, + SPLIT2(infinity_max_vrtx)); log_err(desc->scene->dev, - "Segment %lu:\n (%g %g) (%g %g)\n", - (unsigned long)segments_in[t2].global_id, - SPLIT2(positions[segments_in[t2].vertice_id[0]].vec), - SPLIT2(positions[segments_in[t2].vertice_id[1]].vec)); + "Vertex %lu: (%g %g)\n", + (unsigned long)cc->max_y_vrtx_id, + SPLIT2(max_vrtx)); log_err(desc->scene->dev, "Media: %lu VS %lu\n", - (unsigned long)infinity_medium, (unsigned long)cc->medium); + (unsigned long)(*infinity_first_cc)->medium, (unsigned long)cc->medium); *res = RES_BAD_ARG; } /* Same medium as previous members of the group: OK */ @@ -654,32 +658,27 @@ group_connex_components #endif if(hit_seg_in->medium[hit_side] != cc->medium) { /* Medium mismatch! Model topology is broken. */ - const seg_id_t t1 = SEGSIDE_2_SEG(hit_side_id); - const seg_id_t t2 = SEGSIDE_2_SEG(cc->max_y_side_id); + const seg_id_t seg = SEGSIDE_2_SEG(hit_side_id); const struct segment_in* segments_in = darray_segment_in_cdata_get(&desc->scene->segments_in); - const union double2* positions - = darray_position_cdata_get(&desc->scene->vertices); log_err(desc->scene->dev, - "Medium mismatch found between segment %lu %s side and segment" + "Medium mismatch found between vertex %lu and segment" " %lu %s side facing each other.\n", - (unsigned long)segments_in[t1].global_id, - SEGSIDE_IS_FRONT(hit_side) ? "front" : "back", - (unsigned long)segments_in[t2].global_id, - SEGSIDE_IS_FRONT(cc->max_y_side_id) ? "front" : "back"); + (unsigned long)cc->max_y_vrtx_id, + (unsigned long)segments_in[seg].global_id, + SEGSIDE_IS_FRONT(hit_side_id) ? "front" : "back"); log_err(desc->scene->dev, - "Segment %lu:\n (%g %g) (%g %g)\n", - (unsigned long)segments_in[t1].global_id, - SPLIT2(positions[segments_in[t1].vertice_id[0]].vec), - SPLIT2(positions[segments_in[t1].vertice_id[1]].vec)); + "Vertex %lu: (%g %g)\n", + (unsigned long)cc->max_y_vrtx_id, + SPLIT2(max_vrtx)); log_err(desc->scene->dev, "Segment %lu:\n (%g %g) (%g %g)\n", - (unsigned long)segments_in[t2].global_id, - SPLIT2(positions[segments_in[t2].vertice_id[0]].vec), - SPLIT2(positions[segments_in[t2].vertice_id[1]].vec)); + (unsigned long)segments_in[seg].global_id, + SPLIT2(positions[segments_in[seg].vertice_id[0]].vec), + SPLIT2(positions[segments_in[seg].vertice_id[1]].vec)); log_err(desc->scene->dev, "Media: %lu VS %lu\n", - (unsigned long)hit_seg_in->medium[hit_side], - (unsigned long)cc->medium); + (unsigned long)cc->medium, + (unsigned long)hit_seg_in->medium[hit_side]); *res = RES_BAD_ARG; } } @@ -703,7 +702,6 @@ group_connex_components const struct cc_descriptor* other_desc = cc; struct enclosure_data* enclosures = darray_enclosure_data_get(&desc->enclosures); - component_id_t fst; while(other_desc->enclosure_id == CC_GROUP_ID_NONE) { other_desc = *(darray_ptr_component_descriptor_cdata_get(connex_components) @@ -716,8 +714,6 @@ group_connex_components cc->enclosure_id = other_desc->enclosure_id; ++enclosures[cc->enclosure_id].cc_count; /* Linked list of componnents */ - fst = enclosures[cc->enclosure_id].first_component; - cc->enclosure_next_component = fst; enclosures[cc->enclosure_id].first_component = cc->cc_id; enclosures[cc->enclosure_id].side_range.first = MMIN(enclosures[cc->enclosure_id].side_range.first, cc->side_range.first); @@ -828,23 +824,15 @@ collect_and_link_neighbours neighbour_count = (side_id_t)sz; FOR_EACH(i, 0, neighbour_count) { double max_y, disp[2]; - unsigned char max_y_vrank; struct neighbour_info* neighbour_info = darray_neighbour_data_get(neighbourhood) + i; const struct segment_in* seg_in = segments_in + neighbour_info->seg_id; struct segment_tmp* neighbour = segments_tmp + neighbour_info->seg_id; other_vrtx = seg_in->vertice_id[(neighbour_info->common_vertex_rank + 1) % 2]; - if(vertices[other_vrtx].pos.y > vertices[common_vrtx].pos.y) { - max_y = vertices[other_vrtx].pos.y; - max_y_vrank = (unsigned char)(1 - neighbour_info->common_vertex_rank); - } else { - max_y = vertices[common_vrtx].pos.y; - max_y_vrank = neighbour_info->common_vertex_rank; - } + max_y = MMAX(vertices[other_vrtx].pos.y, vertices[common_vrtx].pos.y); ASSERT(neighbour->max_y <= max_y); neighbour->max_y = max_y; - neighbour->max_y_vrtx_rank = max_y_vrank; /* Compute rotation angle around common vertex (in world system) */ d2_sub(disp, vertices[other_vrtx].vec, vertices[common_vrtx].vec); ASSERT(disp[0] || disp[1]); @@ -866,6 +854,9 @@ collect_and_link_neighbours = darray_neighbour_cdata_get(neighbourhood) + (i + 1) % neighbour_count; /* Rank of the end of interest in segments */ const unsigned char crt_end = current->common_vertex_rank; + /* Here ccw refers to the rotation around the common vertex + * and has nothing to do with vertices order in segment definition + * nor Front/Back side convention */ const unsigned char ccw_end = ccw_neighbour->common_vertex_rank; /* User id of current segments */ const seg_id_t crt_id = current->seg_id; @@ -984,7 +975,7 @@ build_result { tmp_res = darray_uint_push_back(ids_by_medium, &e); } - if (tmp_res != RES_OK) *res = tmp_res; + if(tmp_res != RES_OK) *res = tmp_res; if(*res != RES_OK) continue; /* Build side and vertex lists. */ @@ -1195,7 +1186,7 @@ senc2d_scene_analyze group_connex_components(desc, segsides, &segments_comp, &connex_components, s2d_view, &next_enclosure_id, &infinity_first_cc, &res); /* Barrier at the end of step 3: data used in step 3 can be released / - * data produced by step 2 can be used */ + * data produced by step 3 can be used */ /* One thread releases some data before going to step 4, * the others go to step 4 without sync */ diff --git a/src/senc2d_scene_analyze_c.h b/src/senc2d_scene_analyze_c.h @@ -68,9 +68,8 @@ struct segside { #define CC_GROUP_ROOT_INFINITE (COMPONENT_MAX__-1) #define CC_GROUP_ID_NONE COMPONENT_MAX__ struct cc_descriptor { - double max_vrtx[2]; - double max_y_ny; - side_id_t max_y_side_id; + /* Does this component is an outer border of an enclosure or an inner one? */ + char is_outer_border; vrtx_id_t max_y_vrtx_id; side_id_t side_count; medium_id_t medium; @@ -78,11 +77,8 @@ struct cc_descriptor { component_id_t cc_id; component_id_t cc_group_root; enclosure_id_t enclosure_id; - /* To create by-medium linked lists of componnents */ - component_id_t enclosure_next_component; /* Range of sides member of this component */ struct side_range side_range; - }; extern const struct cc_descriptor CC_DESCRIPTOR_NULL; @@ -135,8 +131,6 @@ custom_darray_ptr_component_descriptor_release * - segment_enc for information describing enclosures (kept in * senc2d_descriptor). */ struct segment_tmp { - /* tmp data used to find the +Y-most vertex of components */ - unsigned char max_y_vrtx_rank; double max_y; }; @@ -145,7 +139,6 @@ static FINLINE void segment_tmp_init(struct mem_allocator* alloc, struct segment_tmp* seg) { (void)alloc; ASSERT(seg); - seg->max_y_vrtx_rank = UCHAR_MAX; seg->max_y = -DBL_MAX; } #define DARRAY_FUNCTOR_INIT segment_tmp_init