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 4b3b9adbe93b97236f280c968e7a7a258fd858ac
parent 9d5aec88857ad39dedc0e70e10b1d4c37eab49c6
Author: Christophe Coustet <christophe.coustet@meso-star.com>
Date:   Tue, 27 Feb 2018 13:56:05 +0100

Move openmp block from outside to inside a fn.

Diffstat:
Msrc/senc_scene_analyze.c | 435+++++++++++++++++++++++++++++++++++++++----------------------------------------
1 file changed, 216 insertions(+), 219 deletions(-)

diff --git a/src/senc_scene_analyze.c b/src/senc_scene_analyze.c @@ -172,11 +172,9 @@ find_component_Zmax if(cc->max_vrtx[2] < trg_tmp->max_z) { change = 1; /* Try first to improve z */ - } - else if(cc->max_z_nz <= 0 && fabs(cc->max_z_nz) < fabs(side_nz)) { + } else if(cc->max_z_nz <= 0 && fabs(cc->max_z_nz) < fabs(side_nz)) { change = 1; /* If nz <= 0, the more negative the better */ - } - else if(cc->max_z_nz > 0 && cc->max_z_nz < side_nz) { + } else if(cc->max_z_nz > 0 && cc->max_z_nz < side_nz) { change = 1; /* If nz > 0, the more positive the better */ } if(change) { @@ -666,206 +664,210 @@ collect_and_link_neighbours struct trgside* trgsides, struct darray_triangle_tmp* triangles_tmp_array) { - trg_id_t t; + res_T res = RES_OK; const struct triangle_in *triangles_in; struct triangle_tmp *triangles_tmp; - edge_id_t approx_nbedges, e, edge_count; const union double3* vertices; - size_t tmp; - const int num_thread = omp_get_num_threads(); - const int rank = omp_get_thread_num(); - struct darray_neighbourhood neighbourhood_by_edge; - /* Htable used to give an id to edges */ - struct htable_edge_id edge_ids; - res_T res = RES_OK; + ASSERT(scn && trgsides && triangles_tmp_array); ASSERT((size_t) scn->nuverts + (size_t) scn->nutris + 2 <= EDGE_MAX__); - darray_neighbourhood_init(scn->dev->allocator, &neighbourhood_by_edge); - - /* Make some room for edges. */ - approx_nbedges = 4 + (num_thread == 1 ? (edge_id_t)(scn->nuverts + scn->nutris) - : (edge_id_t)((scn->nuverts + scn->nutris + num_thread) / (0.75 * num_thread))); - htable_edge_id_init(scn->dev->allocator, &edge_ids); - OK(htable_edge_id_reserve(&edge_ids, approx_nbedges)); - OK(darray_neighbourhood_reserve(&neighbourhood_by_edge, approx_nbedges)); - - /* Loop on triangles to register edges. */ triangles_in = darray_triangle_in_cdata_get(&scn->triangles_in); triangles_tmp = darray_triangle_tmp_data_get(triangles_tmp_array); + vertices = darray_position_cdata_get(&scn->vertices); + ASSERT(scn->nutris == darray_triangle_tmp_size_get(triangles_tmp_array)); - FOR_EACH(t, 0, scn->nutris) { - struct trg_edge edge; - char ee; - ASSERT(t < scn->nutris); - 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 = v0 + v1 + MMIN(v0,v1); /* v0,v1 and v1,v0 must give the same hash!!! */ - if (h % num_thread != rank) continue; - /* Create edge. */ - set_edge(v0, v1, &edge, &triangles_tmp[t].reversed_edge[ee]); - /* Find edge id; create it if not already done. */ - p_id = htable_edge_id_find(&edge_ids, &edge); - if (p_id) { - struct edge_neighbourhood* n = - darray_neighbourhood_data_get(&neighbourhood_by_edge) + *p_id; - struct neighbour_info* info; - size_t sz; - /* Add neighbour info to existing edge's neighbour list */ - ASSERT(n->edge.vrtx0 == edge.vrtx0 && n->edge.vrtx1 == edge.vrtx1); - sz = darray_neighbour_size_get(&n->neighbours); - OK(darray_neighbour_resize(&n->neighbours, 1 + sz)); - info = darray_neighbour_data_get(&n->neighbours) + sz; - info->trg_id = t; - info->edge_rank = ee; - } - else { - /* Create id */ - edge_id_t id; - struct edge_neighbourhood* n; - struct neighbour_info* info; - size_t sz = htable_edge_id_size_get(&edge_ids); - ASSERT(sz <= EDGE_MAX__); - id = (edge_id_t)sz; - ASSERT(htable_edge_id_size_get(&edge_ids) - == darray_neighbourhood_size_get(&neighbourhood_by_edge)); - OK(htable_edge_id_set(&edge_ids, &edge, &id)); - OK(darray_neighbourhood_resize(&neighbourhood_by_edge, 1 + sz)); - n = darray_neighbourhood_data_get(&neighbourhood_by_edge) + sz; - /* Add neighbour info to a newly created edge's neighbour list */ - n->edge = edge; - ASSERT(darray_neighbour_size_get(&n->neighbours) == 0); - OK(darray_neighbour_reserve(&n->neighbours, 2)); - OK(darray_neighbour_resize(&n->neighbours, 1)); - info = darray_neighbour_data_get(&n->neighbours); - info->trg_id = t; - info->edge_rank = ee; +#pragma omp parallel + { + const int thread_count = omp_get_num_threads(); + const int rank = omp_get_thread_num(); + /* Htable used to give an id to edges */ + struct htable_edge_id edge_ids; + /* Array to keep neighbourhood of edges */ + struct darray_neighbourhood neighbourhood_by_edge; + edge_id_t e, edge_count; + edge_id_t nbedges_guess; + trg_id_t t; + size_t tmp; + + htable_edge_id_init(scn->dev->allocator, &edge_ids); + darray_neighbourhood_init(scn->dev->allocator, &neighbourhood_by_edge); + + /* Make some room for edges. */ + nbedges_guess = 4 + (thread_count == 1 + ? (edge_id_t) (scn->nuverts + scn->nutris) + : (edge_id_t) ((scn->nuverts + scn->nutris) / (0.75 * thread_count))); + OK(darray_neighbourhood_reserve(&neighbourhood_by_edge, nbedges_guess)); + OK(htable_edge_id_reserve(&edge_ids, nbedges_guess)); + + /* Loop on triangles to register edges. */ + FOR_EACH(t, 0, scn->nutris) { + struct trg_edge edge; + 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 = + v0 + v1 + MMIN(v0, v1); /* v0,v1 and v1,v0 must give the same hash!!! */ + if(h % thread_count != rank) continue; + /* Create edge. */ + set_edge(v0, v1, &edge, &triangles_tmp[t].reversed_edge[ee]); + /* Find edge id; create it if not already done. */ + p_id = htable_edge_id_find(&edge_ids, &edge); + if(p_id) { + struct edge_neighbourhood* n = + darray_neighbourhood_data_get(&neighbourhood_by_edge) + *p_id; + struct neighbour_info* info; + size_t sz; + /* Add neighbour info to existing edge's neighbour list */ + ASSERT(n->edge.vrtx0 == edge.vrtx0 && n->edge.vrtx1 == edge.vrtx1); + sz = darray_neighbour_size_get(&n->neighbours); + OK(darray_neighbour_resize(&n->neighbours, 1 + sz)); + info = darray_neighbour_data_get(&n->neighbours) + sz; + info->trg_id = t; + info->edge_rank = ee; + } else { + /* Create id */ + edge_id_t id; + struct edge_neighbourhood* n; + struct neighbour_info* info; + size_t sz = htable_edge_id_size_get(&edge_ids); + ASSERT(sz <= EDGE_MAX__); + id = (edge_id_t) sz; + ASSERT(htable_edge_id_size_get(&edge_ids) + == darray_neighbourhood_size_get(&neighbourhood_by_edge)); + OK(htable_edge_id_set(&edge_ids, &edge, &id)); + OK(darray_neighbourhood_resize(&neighbourhood_by_edge, 1 + sz)); + n = darray_neighbourhood_data_get(&neighbourhood_by_edge) + sz; + /* Add neighbour info to a newly created edge's neighbour list */ + n->edge = edge; + ASSERT(darray_neighbour_size_get(&n->neighbours) == 0); + OK(darray_neighbour_reserve(&n->neighbours, 2)); + OK(darray_neighbour_resize(&n->neighbours, 1)); + info = darray_neighbour_data_get(&n->neighbours); + info->trg_id = t; + info->edge_rank = ee; + } } } - } - /* Loop on collected edges. - * For each edge sort triangle sides by rotation angle - * and connect neighbours. */ - vertices = darray_position_cdata_get(&scn->vertices); - tmp = darray_neighbourhood_size_get(&neighbourhood_by_edge); - ASSERT(tmp <= EDGE_MAX__); - 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; - struct edge_neighbourhood* neighbourhood - = 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; - tmp = darray_neighbour_size_get(neighbour_list); - ASSERT(tmp <= SIDE_MAX__); - neighbour_count = (side_id_t)tmp; - ASSERT(neighbour_count); - v0 = neighbourhood->edge.vrtx0; - v1 = neighbourhood->edge.vrtx1; - d3_sub(common_edge, vertices[v1].vec, vertices[v0].vec); - if (vertices[v0].pos.z > vertices[v1].pos.z) { - mz_edge = vertices[v0].pos.z; - mz_vid_edge = 0; - } - else { - mz_edge = vertices[v1].pos.z; - mz_vid_edge = 1; - } - norm = d3_normalize(common_edge, common_edge); - ASSERT(norm); - d33_basis(basis, common_edge); - d33_inverse(basis, basis); - FOR_EACH(i, 0, neighbour_count) { - struct neighbour_info* neighbour_info - = 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; - v2 = trg_in->vertice_id[(neighbour_info->edge_rank + 2) % 3]; - if (vertices[v2].pos.z > mz_edge) { - mz = vertices[v2].pos.z; - mz_vid = 2; - } - else { - mz = mz_edge; - mz_vid = mz_vid_edge; + /* Loop on collected edges. + * For each edge sort triangle sides by rotation angle + * and connect neighbours. */ + tmp = darray_neighbourhood_size_get(&neighbourhood_by_edge); + ASSERT(tmp <= EDGE_MAX__); + 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; + struct edge_neighbourhood* neighbourhood + = 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; + tmp = darray_neighbour_size_get(neighbour_list); + ASSERT(tmp <= SIDE_MAX__); + neighbour_count = (side_id_t) tmp; + ASSERT(neighbour_count); + v0 = neighbourhood->edge.vrtx0; + v1 = neighbourhood->edge.vrtx1; + d3_sub(common_edge, vertices[v1].vec, vertices[v0].vec); + if(vertices[v0].pos.z > vertices[v1].pos.z) { + mz_edge = vertices[v0].pos.z; + mz_vid_edge = 0; + } else { + mz_edge = vertices[v1].pos.z; + mz_vid_edge = 1; } - /* Compute the actual vertex id - * 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; + norm = d3_normalize(common_edge, common_edge); + ASSERT(norm); + d33_basis(basis, common_edge); + d33_inverse(basis, basis); + FOR_EACH(i, 0, neighbour_count) { + struct neighbour_info* neighbour_info + = 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; + v2 = trg_in->vertice_id[(neighbour_info->edge_rank + 2) % 3]; + if(vertices[v2].pos.z > mz_edge) { + mz = vertices[v2].pos.z; + mz_vid = 2; + } else { + mz = mz_edge; + mz_vid = mz_vid_edge; + } + /* Compute the actual vertex id + * 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; + } 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; + } + ASSERT(neighbour->max_z <= mz); + neighbour->max_z = mz; + ASSERT(0 <= actual_vid && 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); + d33_muld3(edge, basis, edge); + ASSERT(d3_len(edge) && (edge[0] || edge[1])); + neighbour_info->angle = atan2(edge[1], edge[0]); + if(neighbour_info->angle < 0) neighbour_info->angle += 2 * PI; + /* Due to catastrophic cancelation, -eps+2pi translates to 2pi */ + ASSERT(0 <= neighbour_info->angle && neighbour_info->angle <= 2 * PI); } - 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; + /* Sort triangles by rotation angle */ + qsort(darray_neighbour_data_get(neighbour_list), neighbour_count, + sizeof(struct neighbour_info), neighbour_cmp); + /* Link sides. + * Create cycles of sides by neighbourhood around common edge. */ + FOR_EACH(i, 0, neighbour_count) { + /* Neighbourhood info for current pair of triangles */ + const struct neighbour_info* current + = darray_neighbour_cdata_get(neighbour_list) + i; + 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; + /* User id of current triangles */ + const trg_id_t crt_id = current->trg_id; + const trg_id_t ccw_id = ccw_neighbour->trg_id; + /* Facing sides of triangles */ + const enum side_id crt_side + = triangles_tmp[crt_id].reversed_edge[crt_edge] ? SIDE_BACK : SIDE_FRONT; + const enum side_id ccw_side + = triangles_tmp[ccw_id].reversed_edge[ccw_edge] ? SIDE_FRONT : SIDE_BACK; + /* Index of sides in trgsides */ + const side_id_t crt_side_idx = TRGIDxSIDE_2_TRGSIDE(crt_id, crt_side); + const side_id_t ccw_side_idx = TRGIDxSIDE_2_TRGSIDE(ccw_id, ccw_side); + /* Side ptrs */ + struct trgside* const p_crt_side = trgsides + crt_side_idx; + struct trgside* const p_ccw_side = trgsides + ccw_side_idx; + /* Link sides */ + p_crt_side->facing_side_id[crt_edge] = ccw_side_idx; + p_ccw_side->facing_side_id[ccw_edge] = crt_side_idx; + /* Record media */ + p_crt_side->medium = triangles_in[crt_id].medium[crt_side]; + p_ccw_side->medium = triangles_in[ccw_id].medium[ccw_side]; + ASSERT(p_crt_side->medium < scn->nmeds); + ASSERT(p_ccw_side->medium < scn->nmeds); + p_crt_side->list_id = FLAG_LIST_SIDE_LIST; + p_ccw_side->list_id = FLAG_LIST_SIDE_LIST; } - ASSERT(neighbour->max_z <= mz); - neighbour->max_z = mz; - ASSERT(0 <= actual_vid && 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); - d33_muld3(edge, basis, edge); - ASSERT(d3_len(edge) && (edge[0] || edge[1])); - neighbour_info->angle = atan2(edge[1], edge[0]); - if (neighbour_info->angle < 0) neighbour_info->angle += 2 * PI; - /* Due to catastrophic cancelation, -eps+2pi translates to 2pi */ - ASSERT(0 <= neighbour_info->angle && neighbour_info->angle <= 2 * PI); - } - /* Sort triangles by rotation angle */ - qsort(darray_neighbour_data_get(neighbour_list), neighbour_count, - sizeof(struct neighbour_info), neighbour_cmp); - /* Link sides. - * Create cycles of sides by neighbourhood around common edge. */ - FOR_EACH(i, 0, neighbour_count) { - /* Neighbourhood info for current pair of triangles */ - const struct neighbour_info* current - = darray_neighbour_cdata_get(neighbour_list) + i; - 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; - /* User id of current triangles */ - const trg_id_t crt_id = current->trg_id; - const trg_id_t ccw_id = ccw_neighbour->trg_id; - /* Facing sides of triangles */ - const enum side_id crt_side - = triangles_tmp[crt_id].reversed_edge[crt_edge] ? SIDE_BACK : SIDE_FRONT; - const enum side_id ccw_side - = triangles_tmp[ccw_id].reversed_edge[ccw_edge] ? SIDE_FRONT : SIDE_BACK; - /* Index of sides in trgsides */ - const side_id_t crt_side_idx = TRGIDxSIDE_2_TRGSIDE(crt_id, crt_side); - const side_id_t ccw_side_idx = TRGIDxSIDE_2_TRGSIDE(ccw_id, ccw_side); - /* Side ptrs */ - struct trgside* const p_crt_side = trgsides + crt_side_idx; - struct trgside* const p_ccw_side = trgsides + ccw_side_idx; - /* Link sides */ - p_crt_side->facing_side_id[crt_edge] = ccw_side_idx; - p_ccw_side->facing_side_id[ccw_edge] = crt_side_idx; - /* Record media */ - p_crt_side->medium = triangles_in[crt_id].medium[crt_side]; - p_ccw_side->medium = triangles_in[ccw_id].medium[ccw_side]; - ASSERT(p_crt_side->medium < scn->nmeds); - ASSERT(p_ccw_side->medium < scn->nmeds); - p_crt_side->list_id = FLAG_LIST_SIDE_LIST; - p_ccw_side->list_id = FLAG_LIST_SIDE_LIST; } +error: + darray_neighbourhood_release(&neighbourhood_by_edge); + htable_edge_id_release(&edge_ids); } -exit: - darray_neighbourhood_release(&neighbourhood_by_edge); - htable_edge_id_release(&edge_ids); + return res; -error: - goto exit; } static res_T @@ -879,39 +881,37 @@ build_result struct enclosure_data* enclosures; const struct triangle_in* triangles_in; struct triangle_enc* triangles_enc; - struct htable_vrtx_id vtable; - size_t tmp; - component_id_t cc_count; - enclosure_id_t e; - trg_id_t fst_ixd = 0; - trg_id_t sgd_ixd; + int ee; ASSERT(desc && connex_components); alloc = descriptor_get_allocator(desc); - tmp = darray_cc_descriptor_size_get(connex_components); - ASSERT(tmp < COMPONENT_MAX__); - cc_count = (component_id_t)tmp; + ASSERT(darray_cc_descriptor_size_get(connex_components) < COMPONENT_MAX__); cc_descriptors = darray_cc_descriptor_data_get(connex_components); enclosures = darray_enclosure_data_get(&desc->enclosures); triangles_in = darray_triangle_in_cdata_get(&desc->scene->triangles_in); - OK(darray_triangle_enc_resize(&desc->triangles_enc, desc->scene->nutris)); + OK2(darray_triangle_enc_resize(&desc->triangles_enc, desc->scene->nutris), + error_); triangles_enc = darray_triangle_enc_data_get(&desc->triangles_enc); - htable_vrtx_id_init(alloc, &vtable); - /* Merge enclosures' membership information */ - FOR_EACH(e, 0, desc->enclosures_count) { +#pragma omp parallel for schedule(dynamic) + FOR_EACH(ee, 0, desc->enclosures_count) { + 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; uint64_t* enc_memb64; + trg_id_t fst_ixd = 0; + trg_id_t sgd_ixd; trg_id_t t; trg_id_t first_member_id = desc->scene->nutris; trg_id_t last_member_id = 0; side_id_t side_count = 0; - ASSERT(enc->cc_count != 0 && enc->cc_count <= cc_count); - ASSERT(enc->first_component <= cc_count); + ASSERT(ee <= ENCLOSURE_MAX__); ASSERT(current != NULL); + + if(res != RES_OK) continue; enc->header.enclosure_id = (unsigned)e; /* Back to API type */ enc->header.is_infinite = (e == 0); enc->header.enclosed_medium @@ -927,16 +927,16 @@ build_result % sizeof(*enc_memb64) == 0); ASSERT(darray_char_size_get(&enc->side_membership) == desc->scene->nutris); + /* Merge enclosures' membership information */ for(;;) { uint64_t* cc_memb64; - trg_id_t t64; + trg_id_t t64, tmin64, tmax64; size_t tmin, tmax; ASSERT(enc->header.enclosed_medium == (unsigned)current->medium); side_count += current->side_count; first_member_id = MMIN(first_member_id, current->first_member_id); last_member_id = MMAX(last_member_id, current->last_member_id); if(current->next_component == COMPONENT_NULL__) break; - ASSERT(current->next_component <= cc_count); current = cc_descriptors + current->next_component; ASSERT(darray_char_size_get(&current->side_membership) == desc->scene->nutris); @@ -948,9 +948,9 @@ build_result / sizeof(*enc_memb64); ASSERT(tmin <= tmax); ASSERT(tmax < TRG_MAX__); - FOR_EACH(t64, (trg_id_t)tmin, (trg_id_t)tmax) { - enc_memb64[t64] |= cc_memb64[t64]; - } + 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); } @@ -958,7 +958,7 @@ build_result OK(darray_triangle_in_resize(&enc->sides, side_count)); OK(darray_vrtx_id_reserve(&enc->vertices, side_count / 2)); /* New vertex numbering scheme local to the enclosure */ - htable_vrtx_id_clear(&vtable); + htable_vrtx_id_init(alloc, &vtable); ASSERT(desc->scene->nutris == darray_triangle_in_size_get(&desc->scene->triangles_in)); fst_ixd = 0; @@ -979,7 +979,7 @@ build_result vertice_id[i] = *id; /* Known vertex */ } else { /* Create new association */ - tmp = htable_vrtx_id_size_get(&vtable); + size_t tmp = htable_vrtx_id_size_get(&vtable); ASSERT(tmp == darray_vrtx_id_size_get(&enc->vertices)); ASSERT(tmp < VRTX_MAX__); vertice_id[i] = (vrtx_id_t)tmp; @@ -1014,13 +1014,14 @@ build_result if(fst_ixd == sgd_ixd) break; } darray_char_purge(&enc->side_membership); + htable_vrtx_id_release(&vtable); + error: + /* Cannot exit openmp block */ + continue; } -exit: - htable_vrtx_id_release(&vtable); +error_: return res; -error: - goto exit; } #include <rsys/clock_time.h> @@ -1080,19 +1081,15 @@ senc_scene_analyze(struct senc_scene* scn, struct senc_descriptor** out_desc) OK(darray_triangle_tmp_resize(&triangles_tmp, scn->nutris)); trgsides = MEM_CALLOC(scn->dev->allocator, 2 * scn->nutris, sizeof(struct trgside)); - if (!trgsides) { + if(!trgsides) { res = RES_MEM_ERR; goto error; } /* Step 1: build neighbourhoods */ -#pragma omp parallel - { - res_T tmp_res = - collect_and_link_neighbours(scn, trgsides, &triangles_tmp); - if(tmp_res != RES_OK) res = tmp_res; - } - if (res != RES_OK) { + res = collect_and_link_neighbours(scn, trgsides, &triangles_tmp); + + if(res != RES_OK) { log_err(scn->dev, "%s: could not build neighbourhoods from scene.\n", FUNC_NAME); goto error;