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 c5b914b2ed26fff04309fbe124f0c5a3dc873189
parent cbb5805453ccc42394c74676e6f4d39dad1db97d
Author: Christophe Coustet <christophe.coustet@meso-star.com>
Date:   Thu, 12 Apr 2018 17:39:53 +0200

Some performance mprovement by easing omp thread syncs.

Diffstat:
Msrc/senc_internal_types.h | 23++++++++++++++++-------
Msrc/senc_scene_analyze.c | 1511++++++++++++++++++++++++++++++++++++++++++-------------------------------------
Msrc/senc_scene_analyze_c.h | 23+++++++++++++++--------
3 files changed, 826 insertions(+), 731 deletions(-)

diff --git a/src/senc_internal_types.h b/src/senc_internal_types.h @@ -22,18 +22,27 @@ /* Utility macros */ #ifdef NDEBUG -#define OK2(Expr, Label)\ - if((res = (Expr)) != RES_OK) goto Label; +#define OK2(Expr)\ + if((tmp_res = (Expr)) != RES_OK) goto tmp_error; + +#define OK(Expr)\ + if((res = (Expr)) != RES_OK) goto error; #else -#define OK2(Expr, Label)\ +#define OK2(Expr)\ + if((tmp_res = (Expr)) != RES_OK) {\ + fprintf(stderr, "%s: error code set to %d at line %d\n", FUNC_NAME,\ + tmp_res, __LINE__);\ + goto tmp_error;\ + } + +#define OK(Expr)\ if((res = (Expr)) != RES_OK) {\ - fprintf(stderr, "%s: error code set to %d at line %d\n", FUNC_NAME, res, __LINE__);\ - goto Label;\ + fprintf(stderr, "%s: error code set to %d at line %d\n", FUNC_NAME,\ + res, __LINE__);\ + goto error;\ } #endif -#define OK(Expr) OK2((Expr), error) - /* Side IDs are uint32_t */ typedef uint32_t side_id_t; #define SIDE_MAX__ (UINT32_MAX-1) diff --git a/src/senc_scene_analyze.c b/src/senc_scene_analyze.c @@ -155,7 +155,7 @@ get_scn_indices(const unsigned itri, unsigned ids[3], void* ctx) { darray_triangle_in_cdata_get(&scene->triangles_in) + itri; FOR_EACH(i, 0, 3) { ASSERT(trg->vertice_id[i] < scene->nverts); - ids[i] = (unsigned) trg->vertice_id[i]; /* Back to API type */ + ids[i] = (unsigned)trg->vertice_id[i]; /* Back to API type */ } } @@ -194,350 +194,352 @@ self_hit_filter return (hit_component == *origin_component); } -static res_T +static void extract_connex_components (struct senc_descriptor* desc, struct trgside* trgsides, struct darray_ptr_component_descriptor* connex_components, const struct darray_triangle_tmp* triangles_tmp_array, struct darray_triangle_comp* triangles_comp, - struct s3d_scene_view** s3d_view) + struct s3d_scene_view** s3d_view, + ATOMIC* component_count, + /* Shared error status. + * We accept to overwritte an error with a different error */ + res_T* res) { - res_T res = RES_OK; + /* This function is called from an omp parallel block and executed + * concurrently. */ const struct senc_scene* scn; struct mem_allocator* alloc; - ATOMIC component_count = 0; - volatile int exit_for = 0; int64_t mm; + struct darray_side_id stack; #ifndef NDEBUG trg_id_t t; component_id_t c; #endif - ASSERT(trgsides && desc && connex_components && triangles_tmp_array); - ASSERT(darray_ptr_component_descriptor_size_get(connex_components) == 0); + ASSERT(trgsides && desc && connex_components && triangles_tmp_array + && component_count && res); alloc = descriptor_get_allocator(desc); scn = desc->scene; - /* Just a hint; to avoid contention on first loop */ - OK2(darray_ptr_component_descriptor_reserve(connex_components, 2 * scn->nmeds), - error_); /* Cannot goto into openmp block */ - #ifndef NDEBUG - FOR_EACH(t, 0, scn->nutris) { - const struct triangle_in* trg_in = - darray_triangle_in_cdata_get(&scn->triangles_in) + t; - const struct side_range* media_use - = darray_side_range_cdata_get(&scn->media_use); - FOR_EACH(mm, 0, 2) { - const side_id_t side = TRGIDxSIDE_2_TRGSIDE(t, mm); - const medium_id_t medium = trg_in->medium[mm]; - ASSERT(media_use[medium].first <= side && side <= media_use[medium].last); + #pragma omp single + { + ASSERT(darray_ptr_component_descriptor_size_get(connex_components) == 0); + FOR_EACH(t, 0, scn->nutris) { + const struct triangle_in* trg_in = + darray_triangle_in_cdata_get(&scn->triangles_in) + t; + const struct side_range* media_use + = darray_side_range_cdata_get(&scn->media_use); + FOR_EACH(mm, 0, 2) { + const side_id_t side = TRGIDxSIDE_2_TRGSIDE(t, mm); + const medium_id_t medium = trg_in->medium[mm]; + ASSERT(media_use[medium].first <= side && side + <= media_use[medium].last); + } } } #endif - #pragma omp parallel num_threads(scn->dev->nthreads) - { - struct darray_side_id stack; - 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; - - if(exit_for) continue; - first_side_not_in_component - = darray_side_range_cdata_get(&scn->media_use)[m].first; - if(first_side_not_in_component == SIDE_NULL__) - continue; /* Unused medium */ - ASSERT(first_side_not_in_component < 2 * scn->nutris); - ASSERT(darray_side_id_size_get(&stack) == 0); - for(;;) { /* Process all components for this medium */ - const side_id_t start_side_id = get_side_not_in_connex_component(scn, - trgsides, &first_side_not_in_component, m); - side_id_t crt_side_id = start_side_id; - side_id_t last_side_id = start_side_id; - ASSERT(start_side_id == SIDE_NULL__ || start_side_id < 2 * scn->nutris); - if(start_side_id == SIDE_NULL__) - break; /* start_side_id=SIDE_NULL__ => done! */ - ASSERT(trgsides[start_side_id].list_id == FLAG_LIST_SIDE_LIST); + 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; + + if(*res != RES_OK) continue; + first_side_not_in_component + = darray_side_range_cdata_get(&scn->media_use)[m].first; + if(first_side_not_in_component == SIDE_NULL__) + continue; /* Unused medium */ + ASSERT(first_side_not_in_component < 2 * scn->nutris); + ASSERT(darray_side_id_size_get(&stack) == 0); + for(;;) { /* Process all components for this medium */ + const side_id_t start_side_id = get_side_not_in_connex_component(scn, + trgsides, &first_side_not_in_component, m); + side_id_t crt_side_id = start_side_id; + side_id_t last_side_id = start_side_id; + ASSERT(start_side_id == SIDE_NULL__ || start_side_id < 2 * scn->nutris); + if(start_side_id == SIDE_NULL__) + break; /* start_side_id=SIDE_NULL__ => done! */ + ASSERT(trgsides[start_side_id].list_id == FLAG_LIST_SIDE_LIST); #ifndef NDEBUG - { - trg_id_t tid = TRGSIDE_2_TRG(start_side_id); - enum side_flag s = TRGSIDE_2_SIDE(start_side_id); - medium_id_t side_med - = darray_triangle_in_data_get(&desc->scene->triangles_in)[tid].medium[s]; - ASSERT(side_med == m); - } + { + trg_id_t tid = TRGSIDE_2_TRG(start_side_id); + enum side_flag s = TRGSIDE_2_SIDE(start_side_id); + medium_id_t side_med + = darray_triangle_in_data_get(&desc->scene->triangles_in)[tid].medium[s]; + ASSERT(side_med == m); + } #endif - /* Create and init a new component */ - cc = MEM_ALLOC(alloc, sizeof(struct cc_descriptor)); - if(!cc) { - res = RES_MEM_ERR; - goto error1; - } - cc_descriptor_init(alloc, cc); - ASSERT(m == trgsides[start_side_id].medium); - cc->cc_id = (component_id_t)(ATOMIC_INCR(&component_count) - 1); - cc->medium = m; - cc->side_range.first = start_side_id; - - for(;;) { /* Process all sides of this component */ - int i; - enum side_flag crt_side_flag = TRGSIDE_2_SIDE(crt_side_id); - struct trgside* crt_side = trgsides + crt_side_id; - const trg_id_t crt_trg_id = TRGSIDE_2_TRG(crt_side_id); - const struct triangle_in* trg_in = - darray_triangle_in_cdata_get(&scn->triangles_in) + crt_trg_id; - struct triangle_comp* trg_comp = - darray_triangle_comp_data_get(triangles_comp) + crt_trg_id; - const struct triangle_tmp* const trg_tmp = - darray_triangle_tmp_cdata_get(triangles_tmp_array) + crt_trg_id; - const union double3* vertices = - darray_position_cdata_get(&scn->vertices); - ASSERT(crt_trg_id < scn->nutris); - ASSERT(trgsides[crt_side_id].medium == m); - - /* Record Zmax 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. - * The most appropriate vertex is (the) one with the greater Z - * coordinate. - * If more than one vertex/side has the same Z, we want the side that - * most faces Z (that is the one with the greater nz). - * This is mandatory to select the correct side when both sides of a - * triangle are candidate. */ - if(cc->max_vrtx[2] <= trg_tmp->max_z) { - /* Can either improve z or nz */ - - if(cc->max_z_side_id == TRGSIDE_OPPOSITE(crt_side_id)) { - /* Both sides are in cc and the opposite side is currently selected - * Just keep the side with nz>0 */ - if(cc->max_z_nz < 0) { - /* Change side! */ - cc->max_z_nz = -cc->max_z_nz; - cc->max_z_side_id = crt_side_id; - } + /* Create and init a new component */ + cc = MEM_ALLOC(alloc, sizeof(struct cc_descriptor)); + if(!cc) { + *res = RES_MEM_ERR; + continue; + } + cc_descriptor_init(alloc, cc); + ASSERT(m == trgsides[start_side_id].medium); + cc->cc_id = (component_id_t)(ATOMIC_INCR(component_count) - 1); + cc->medium = m; + cc->side_range.first = start_side_id; + + for(;;) { /* Process all sides of this component */ + int i; + enum side_flag crt_side_flag = TRGSIDE_2_SIDE(crt_side_id); + struct trgside* crt_side = trgsides + crt_side_id; + const trg_id_t crt_trg_id = TRGSIDE_2_TRG(crt_side_id); + const struct triangle_in* trg_in = + darray_triangle_in_cdata_get(&scn->triangles_in) + crt_trg_id; + struct triangle_comp* trg_comp = + darray_triangle_comp_data_get(triangles_comp) + crt_trg_id; + const struct triangle_tmp* const trg_tmp = + darray_triangle_tmp_cdata_get(triangles_tmp_array) + crt_trg_id; + const union double3* vertices = + darray_position_cdata_get(&scn->vertices); + ASSERT(crt_trg_id < scn->nutris); + ASSERT(trgsides[crt_side_id].medium == m); + + /* Record Zmax 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. + * The most appropriate vertex is (the) one with the greater Z + * coordinate. + * If more than one vertex/side has the same Z, we want the side that + * most faces Z (that is the one with the greater nz). + * This is mandatory to select the correct side when both sides of a + * triangle are candidate. */ + if(cc->max_vrtx[2] <= trg_tmp->max_z) { + /* Can either improve z or nz */ + + if(cc->max_z_side_id == TRGSIDE_OPPOSITE(crt_side_id)) { + /* Both sides are in cc and the opposite side is currently selected + * Just keep the side with nz>0 */ + if(cc->max_z_nz < 0) { + /* Change side! */ + cc->max_z_nz = -cc->max_z_nz; + cc->max_z_side_id = crt_side_id; + } + } else { + double edge0[3], edge1[3], normal[3], norm, side_nz; + int process = 0; + + d3_sub(edge0, vertices[trg_in->vertice_id[1]].vec, + vertices[trg_in->vertice_id[0]].vec); + d3_sub(edge1, vertices[trg_in->vertice_id[2]].vec, + vertices[trg_in->vertice_id[0]].vec); + d3_cross(normal, edge0, edge1); + norm = d3_normalize(normal, normal); + ASSERT(norm); (void) norm; + + if(TRGSIDE_IS_FRONT(crt_side_id)) { + side_nz = normal[2]; + process = 1; } else { - double edge0[3], edge1[3], normal[3], norm, side_nz; - int process = 0; - - d3_sub(edge0, vertices[trg_in->vertice_id[1]].vec, - vertices[trg_in->vertice_id[0]].vec); - d3_sub(edge1, vertices[trg_in->vertice_id[2]].vec, - vertices[trg_in->vertice_id[0]].vec); - d3_cross(normal, edge0, edge1); - norm = d3_normalize(normal, normal); - ASSERT(norm); (void) norm; - - if(TRGSIDE_IS_FRONT(crt_side_id)) { - side_nz = normal[2]; - process = 1; - } else { - side_nz = -normal[2]; - process = 1; + side_nz = -normal[2]; + process = 1; + } + + if(process) { + int change = 0; + 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)) { + change = 1; /* If nz <= 0, the more negative the better */ + } else if(cc->max_z_nz > 0 && cc->max_z_nz < side_nz) { + change = 1; /* If nz > 0, the more positive the better */ } - if(process) { - int change = 0; - 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)) { - change = 1; /* If nz <= 0, the more negative the better */ - } 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) { - cc->max_z_nz = side_nz; - cc->max_z_side_id = crt_side_id; - ASSERT(trg_tmp->max_z_vrtx_rank < 3); - ASSERT(trg_in->vertice_id[trg_tmp->max_z_vrtx_rank] < scn->nverts); - cc->max_z_vrtx_id = trg_in->vertice_id[trg_tmp->max_z_vrtx_rank]; - ASSERT(trg_tmp->max_z == vertices[cc->max_z_vrtx_id].pos.z); - d3_set(cc->max_vrtx, vertices[cc->max_z_vrtx_id].vec); - } + if(change) { + cc->max_z_nz = side_nz; + cc->max_z_side_id = crt_side_id; + ASSERT(trg_tmp->max_z_vrtx_rank < 3); + ASSERT(trg_in->vertice_id[trg_tmp->max_z_vrtx_rank] < scn->nverts); + cc->max_z_vrtx_id = trg_in->vertice_id[trg_tmp->max_z_vrtx_rank]; + ASSERT(trg_tmp->max_z == vertices[cc->max_z_vrtx_id].pos.z); + d3_set(cc->max_vrtx, vertices[cc->max_z_vrtx_id].vec); } } } + } - /* Record crt_side both as component and triangle level */ - cc->side_count++; - trgsides[crt_side_id].list_id = FLAG_LIST_COMPONENT; - ASSERT(trg_comp->component[crt_side_flag] == COMPONENT_NULL__); - trg_comp->component[crt_side_flag] = cc->cc_id; + /* Record crt_side both as component and triangle level */ + cc->side_count++; + trgsides[crt_side_id].list_id = FLAG_LIST_COMPONENT; + ASSERT(trg_comp->component[crt_side_flag] == COMPONENT_NULL__); + trg_comp->component[crt_side_flag] = cc->cc_id; #ifndef NDEBUG - crt_side->member_of_cc = cc->cc_id; + crt_side->member_of_cc = cc->cc_id; #endif - /* Store neighbour sides in a waiting stack */ - FOR_EACH(i, 0, 3) { - side_id_t neighbour_id = crt_side->facing_side_id[i]; - trg_id_t nbour_trg_id = TRGSIDE_2_TRG(neighbour_id); - const struct trgside* neighbour = trgsides + neighbour_id; - ASSERT(m == crt_side->medium); - if(neighbour->medium != crt_side->medium) { - /* Found medium discontinuity! Model topology is broken. */ - const struct triangle_in* triangles_in - = darray_triangle_in_cdata_get(&scn->triangles_in); - const union double3* positions - = darray_position_cdata_get(&scn->vertices); - log_err(scn->dev, - "Medium mismatch found between neighbour triangles %lu %s" - " side and %u %s side.\n", - (unsigned long)triangles_in[crt_trg_id].global_id, - TRGSIDE_IS_FRONT(crt_side_id) ? "front" : "back", - triangles_in[nbour_trg_id].global_id, - TRGSIDE_IS_FRONT(neighbour_id) ? "front" : "back"); - log_err(scn->dev, - "Triangle %lu:\n (%g %g %g)\n (%g %g %g)\n (%g %g %g)\n", - (unsigned long)triangles_in[crt_trg_id].global_id, - SPLIT3(positions[triangles_in[crt_trg_id].vertice_id[0]].vec), - SPLIT3(positions[triangles_in[crt_trg_id].vertice_id[1]].vec), - SPLIT3(positions[triangles_in[crt_trg_id].vertice_id[2]].vec)); - log_err(scn->dev, - "Triangle %lu:\n (%g %g %g)\n (%g %g %g)\n (%g %g %g)\n", - (unsigned long)triangles_in[nbour_trg_id].global_id, - SPLIT3(positions[triangles_in[nbour_trg_id].vertice_id[0]].vec), - 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", - (unsigned long)neighbour->medium, (unsigned long)crt_side->medium); - res = RES_BAD_ARG; - goto error1; - } - if(neighbour->list_id == FLAG_LIST_COMPONENT) { - /* Already processed */ + /* Store neighbour sides in a waiting stack */ + FOR_EACH(i, 0, 3) { + side_id_t neighbour_id = crt_side->facing_side_id[i]; + trg_id_t nbour_trg_id = TRGSIDE_2_TRG(neighbour_id); + const struct trgside* neighbour = trgsides + neighbour_id; + ASSERT(m == crt_side->medium); + if(neighbour->medium != crt_side->medium) { + /* Found medium discontinuity! Model topology is broken. */ + const struct triangle_in* triangles_in + = darray_triangle_in_cdata_get(&scn->triangles_in); + const union double3* positions + = darray_position_cdata_get(&scn->vertices); + log_err(scn->dev, + "Medium mismatch found between neighbour triangles %lu %s" + " side and %u %s side.\n", + (unsigned long)triangles_in[crt_trg_id].global_id, + TRGSIDE_IS_FRONT(crt_side_id) ? "front" : "back", + triangles_in[nbour_trg_id].global_id, + TRGSIDE_IS_FRONT(neighbour_id) ? "front" : "back"); + log_err(scn->dev, + "Triangle %lu:\n (%g %g %g)\n (%g %g %g)\n (%g %g %g)\n", + (unsigned long)triangles_in[crt_trg_id].global_id, + SPLIT3(positions[triangles_in[crt_trg_id].vertice_id[0]].vec), + SPLIT3(positions[triangles_in[crt_trg_id].vertice_id[1]].vec), + SPLIT3(positions[triangles_in[crt_trg_id].vertice_id[2]].vec)); + log_err(scn->dev, + "Triangle %lu:\n (%g %g %g)\n (%g %g %g)\n (%g %g %g)\n", + (unsigned long)triangles_in[nbour_trg_id].global_id, + SPLIT3(positions[triangles_in[nbour_trg_id].vertice_id[0]].vec), + 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", + (unsigned long)neighbour->medium, (unsigned long)crt_side->medium); + *res = RES_BAD_ARG; + continue; + } + if(neighbour->list_id == FLAG_LIST_COMPONENT) { + /* Already processed */ #ifndef NDEBUG - ASSERT(neighbour->member_of_cc == cc->cc_id); + ASSERT(neighbour->member_of_cc == cc->cc_id); #endif - continue; - } - if(neighbour->list_id == FLAG_WAITING_STACK) { - continue; /* Already processed */ - } - add_side_to_stack(scn, &stack, trgsides, neighbour_id); + continue; } - if(darray_side_id_size_get(&stack) == 0) - break; /* Empty stack => connex component is done! */ - crt_side_id = get_side_from_stack(trgsides, &stack); - last_side_id = MMAX(last_side_id, crt_side_id); - } - /* Keep track of this new connex component */ - cc->side_range.last = last_side_id; - /* 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); - if(sz <= cc->cc_id) { - res_T tmp_res = darray_ptr_component_descriptor_resize(connex_components, - 1 + cc->cc_id); - if(tmp_res != RES_OK) res = tmp_res; - } - if(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); - ASSERT(components[cc->cc_id] == NULL); - components[cc->cc_id] = cc; + if(neighbour->list_id == FLAG_WAITING_STACK) { + continue; /* Already processed */ } + add_side_to_stack(scn, &stack, trgsides, neighbour_id); + } + if(darray_side_id_size_get(&stack) == 0) + break; /* Empty stack => connex component is done! */ + crt_side_id = get_side_from_stack(trgsides, &stack); + last_side_id = MMAX(last_side_id, crt_side_id); + } + /* Keep track of this new connex component */ + cc->side_range.last = last_side_id; + /* 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); + if(sz <= cc->cc_id) { + res_T tmp_res = darray_ptr_component_descriptor_resize(connex_components, + 1 + cc->cc_id); + if(tmp_res != RES_OK) *res = tmp_res; + } + if(*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); + ASSERT(components[cc->cc_id] == NULL); + components[cc->cc_id] = cc; } - OK2(res, error1); } - continue; - error1: - /* Cannot goto out of openmp block */ - exit_for = 1; - continue; - } - /* No barrier here (nowait clause). - * The first thread executes the single block */ - darray_side_id_release(&stack); - #pragma omp single nowait - if(res == RES_OK) { - struct s3d_device* s3d = NULL; - struct s3d_scene* s3d_scn = NULL; - struct s3d_shape* s3d_shp = NULL; - struct s3d_vertex_data attribs; - - attribs.type = S3D_FLOAT3; - attribs.usage = S3D_POSITION; - attribs.get = get_scn_position; - - /* Put geometry in a 3D view */ - OK(s3d_device_create(desc->scene->dev->logger, alloc, 0, &s3d)); - OK(s3d_scene_create(s3d, &s3d_scn)); - OK(s3d_shape_create_mesh(s3d, &s3d_shp)); - - /* Back to API type for ntris and nverts */ - ASSERT(desc->scene->nutris < UINT_MAX); - ASSERT(desc->scene->nuverts < UINT_MAX); - OK(s3d_mesh_setup_indexed_vertices(s3d_shp, (unsigned) desc->scene->nutris, - get_scn_indices, (unsigned) desc->scene->nuverts, &attribs, 1, desc->scene)); - s3d_mesh_set_hit_filter_function(s3d_shp, self_hit_filter, triangles_comp); - OK(s3d_scene_attach_shape(s3d_scn, s3d_shp)); - OK(s3d_scene_view_create(s3d_scn, S3D_TRACE, s3d_view)); - error: - if(s3d) S3D(device_ref_put(s3d)); - if(s3d_scn) S3D(scene_ref_put(s3d_scn)); - if(s3d_shp) S3D(shape_ref_put(s3d_shp)); } + } /* No barrier here */ + + darray_side_id_release(&stack); + + /* The first thread here creates the s3d view */ + #pragma omp single nowait + if(*res == RES_OK) { + res_T tmp_res = RES_OK; + struct s3d_device* s3d = NULL; + struct s3d_scene* s3d_scn = NULL; + struct s3d_shape* s3d_shp = NULL; + struct s3d_vertex_data attribs; + + attribs.type = S3D_FLOAT3; + attribs.usage = S3D_POSITION; + attribs.get = get_scn_position; + + /* Put geometry in a 3D view */ + OK2(s3d_device_create(desc->scene->dev->logger, alloc, 0, &s3d)); + OK2(s3d_scene_create(s3d, &s3d_scn)); + OK2(s3d_shape_create_mesh(s3d, &s3d_shp)); + + /* Back to API type for ntris and nverts */ + ASSERT(desc->scene->nutris < UINT_MAX); + ASSERT(desc->scene->nuverts < UINT_MAX); + OK2(s3d_mesh_setup_indexed_vertices(s3d_shp, + (unsigned)desc->scene->nutris, get_scn_indices, + (unsigned)desc->scene->nuverts, &attribs, 1, desc->scene)); + s3d_mesh_set_hit_filter_function(s3d_shp, self_hit_filter, triangles_comp); + OK2(s3d_scene_attach_shape(s3d_scn, s3d_shp)); + OK2(s3d_scene_view_create(s3d_scn, S3D_TRACE, s3d_view)); + tmp_error: + if(tmp_res != RES_OK) *res = tmp_res; + if(s3d) S3D(device_ref_put(s3d)); + if(s3d_scn) S3D(scene_ref_put(s3d_scn)); + if(s3d_shp) S3D(shape_ref_put(s3d_shp)); } - OK2(res, error_); + if(*res != RES_OK) return; - ASSERT(component_count == - (int)darray_ptr_component_descriptor_size_get(connex_components)); #ifndef NDEBUG - FOR_EACH(t, 0, scn->nutris) { - struct triangle_comp* trg_comp = - darray_triangle_comp_data_get(triangles_comp) + t; - ASSERT(trg_comp->component[SIDE_FRONT] != COMPONENT_NULL__); - ASSERT(trg_comp->component[SIDE_BACK] != COMPONENT_NULL__); - } - FOR_EACH(c, 0, component_count) { - struct cc_descriptor** components = - darray_ptr_component_descriptor_data_get(connex_components); - ASSERT(components[c] != NULL && - components[c]->cc_id == c); + /* Need to wait for all threads done to be able to check stuff */ + #pragma omp barrier + #pragma omp single + { + ASSERT(*component_count == + (int)darray_ptr_component_descriptor_size_get(connex_components)); + FOR_EACH(t, 0, scn->nutris) { + struct triangle_comp* trg_comp = + darray_triangle_comp_data_get(triangles_comp) + t; + ASSERT(trg_comp->component[SIDE_FRONT] != COMPONENT_NULL__); + ASSERT(trg_comp->component[SIDE_BACK] != COMPONENT_NULL__); + } + FOR_EACH(c, 0, *component_count) { + struct cc_descriptor** components = + darray_ptr_component_descriptor_data_get(connex_components); + ASSERT(components[c] != NULL && + components[c]->cc_id == c); + } } #endif - -exit: - ASSERT(desc->triangle_count - == darray_triangle_comp_size_get(triangles_comp)); - /* triangles_enc is still unused: no size to assert */ - return res; -error_: - goto exit; } -static res_T +static void group_connex_components (struct senc_descriptor* desc, struct trgside* trgsides, struct darray_triangle_comp* triangles_comp, struct darray_ptr_component_descriptor* connex_components, - struct s3d_scene_view* s3d_view) + struct s3d_scene_view* s3d_view, + ATOMIC* next_enclosure_id, + struct cc_descriptor** infinity_first_cc, + /* Shared error status. + * We accept to overwritte an error with a different error */ + res_T* res) { - res_T res = RES_OK; + /* This function is called from an omp parallel block and executed + * concurrently. */ struct cc_descriptor** descriptors; size_t tmp; component_id_t cc_count; int64_t ccc; - volatile int exit_for = 0; - ATOMIC next_enclosure_id = desc->enclosures_count; - struct cc_descriptor* infinity_first_cc = NULL; + (void)trgsides; - ASSERT(desc && trgsides && triangles_comp && connex_components); + ASSERT(desc && trgsides && triangles_comp && connex_components + && s3d_view && next_enclosure_id &&infinity_first_cc && res); + ASSERT(desc->enclosures_count == 1); descriptors = darray_ptr_component_descriptor_data_get(connex_components); tmp = darray_ptr_component_descriptor_size_get(connex_components); @@ -545,8 +547,9 @@ group_connex_components cc_count = (component_id_t)tmp; /* Cast rays to find links between connex components */ - #pragma omp parallel for num_threads(desc->scene->dev->nthreads) + #pragma omp for for(ccc = 0; ccc < (int64_t)cc_count; ccc++) { + res_T tmp_res = RES_OK; component_id_t c = (component_id_t)ccc; struct s3d_hit hit = S3D_HIT_NULL; float origin[3]; @@ -558,7 +561,7 @@ group_connex_components component_id_t self_hit_component = origin_trg->component[1 - TRGSIDE_2_SIDE(cc->max_z_side_id)]; - if(exit_for) continue; + if(*res != RES_OK) continue; ASSERT(cc->cc_id == c); ASSERT(cc->cc_group_root == CC_GROUP_ID_NONE); @@ -566,7 +569,7 @@ group_connex_components int64_t id; /* Don't need to cast a ray */ cc->cc_group_root = cc->cc_id; /* New group with self as root */ - id = ATOMIC_INCR(&next_enclosure_id) - 1; + id = ATOMIC_INCR(next_enclosure_id) - 1; ASSERT(id < ENCLOSURE_MAX__); cc->enclosure_id = (enclosure_id_t)id; continue; @@ -584,17 +587,21 @@ group_connex_components == TRGSIDE_OPPOSITE(cc->max_z_side_id)))); f3_set_d3(origin, cc->max_vrtx); /* Self-hit data: self hit if hit this component "on the other side" */ - OK2(s3d_scene_view_trace_ray(s3d_view, origin, dir, range, - &self_hit_component, &hit), error_); + tmp_res = s3d_scene_view_trace_ray(s3d_view, origin, dir, range, + &self_hit_component, &hit); + if(tmp_res != RES_OK) { + *res = tmp_res; + continue; + } /* If no hit, the component is facing an infinite medium */ if(S3D_HIT_NONE(&hit)) { cc->cc_group_root = CC_GROUP_ROOT_INFINITE; cc->enclosure_id = 0; /* 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_z_side_id; - const medium_id_t infinity_medium = infinity_first_cc->medium; + 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_z_side_id; + const medium_id_t infinity_medium = (*infinity_first_cc)->medium; /* Medium mismatch! Model topology is broken. */ const trg_id_t t1 = TRGSIDE_2_TRG(infinity_first_side); const trg_id_t t2 = TRGSIDE_2_TRG(cc->max_z_side_id); @@ -623,8 +630,7 @@ group_connex_components SPLIT3(positions[triangles_in[t2].vertice_id[2]].vec)); log_err(desc->scene->dev, "Media: %lu VS %lu\n", (unsigned long)infinity_medium, (unsigned long)cc->medium); - res = RES_BAD_ARG; - goto error_; + *res = RES_BAD_ARG; } /* Same medium as previous members of the group: OK */ continue; @@ -682,428 +688,426 @@ group_connex_components log_err(desc->scene->dev, "Media: %lu VS %lu\n", (unsigned long)hit_trg_in->medium[hit_side], (unsigned long)cc->medium); - - res = RES_BAD_ARG; - goto error_; + *res = RES_BAD_ARG; } } - continue; - error_: - /* Cannot goto out of openmp block */ - exit_for = 1; } - ASSERT(next_enclosure_id < ENCLOSURE_MAX__); - desc->enclosures_count = (enclosure_id_t)next_enclosure_id; - OK(res); + /* Implicit barrier here */ + ASSERT(*next_enclosure_id < ENCLOSURE_MAX__); + if(*res != RES_OK) return; - /* Post-process links to group connex components */ - OK(darray_enclosure_resize(&desc->enclosures, desc->enclosures_count)); - FOR_EACH(ccc, 0, cc_count) { - component_id_t c = (component_id_t)ccc; - struct cc_descriptor* const cc = descriptors[c]; - 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) - + other_desc->cc_group_root); + /* One thread post-processes links to group connex components */ + #pragma omp single + { + res_T tmp_res = RES_OK; + desc->enclosures_count = (enclosure_id_t)*next_enclosure_id; + tmp_res = darray_enclosure_resize(&desc->enclosures, desc->enclosures_count); + if(tmp_res != RES_OK) { + *res = tmp_res; + } else { + FOR_EACH(ccc, 0, cc_count) { + component_id_t c = (component_id_t)ccc; + struct cc_descriptor* const cc = descriptors[c]; + 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) + + other_desc->cc_group_root); + } + ASSERT(other_desc->cc_group_root != CC_GROUP_ROOT_NONE); + ASSERT(other_desc->enclosure_id != CC_GROUP_ID_NONE); + ASSERT(cc->medium == other_desc->medium); + cc->cc_group_root = other_desc->cc_group_root; + 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); + enclosures[cc->enclosure_id].side_range.last + = MMAX(enclosures[cc->enclosure_id].side_range.last, cc->side_range.last); + enclosures[cc->enclosure_id].side_count += cc->side_count; + } } - ASSERT(other_desc->cc_group_root != CC_GROUP_ROOT_NONE); - ASSERT(other_desc->enclosure_id != CC_GROUP_ID_NONE); - ASSERT(cc->medium == other_desc->medium); - cc->cc_group_root = other_desc->cc_group_root; - 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); - enclosures[cc->enclosure_id].side_range.last - = MMAX(enclosures[cc->enclosure_id].side_range.last, cc->side_range.last); - enclosures[cc->enclosure_id].side_count += cc->side_count; } - -exit: - return res; -error: - goto exit; + /* Implicit barrier here */ } -static res_T +static void collect_and_link_neighbours (struct senc_scene* scn, struct trgside* trgsides, - struct darray_triangle_tmp* triangles_tmp_array) + struct darray_triangle_tmp* triangles_tmp_array, + /* Shared error status. + * We accept to overwritte an error with a different error */ + res_T* res) { - res_T res = RES_OK; + /* This function is called from an omp parallel block and executed + * concurrently. + * Resize / Push operations on neighbourhood_by_vertex are valid + * because each neighbourhood is processes by an unique thread */ const struct triangle_in *triangles_in; struct triangle_tmp *triangles_tmp; const union double3* vertices; + 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 + * Resize/Push operations on neighbourhood_by_edge are valid in the + * openmp block because it is thread local data */ + struct darray_neighbourhood neighbourhood_by_edge; + edge_id_t edge_count; + edge_id_t nbedges_guess; + edge_id_t e; + trg_id_t t; + size_t sz; + res_T tmp_res; - ASSERT(scn && trgsides && triangles_tmp_array); + ASSERT(scn && trgsides && triangles_tmp_array && res); ASSERT((size_t)scn->nuverts + (size_t)scn->nutris + 2 <= EDGE_MAX__); + htable_edge_id_init(scn->dev->allocator, &edge_ids); + darray_neighbourhood_init(scn->dev->allocator, &neighbourhood_by_edge); + 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)); - - #pragma omp parallel num_threads(scn->dev->nthreads) - { - 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 - * Resize/Push operations on neighbourhood_by_edge are valid in the - * openmp block because it is thread local data */ - 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. - * All threads considering all the edges and processing some */ - FOR_EACH(t, 0, scn->nutris) { - struct trg_edge edge; - 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 int64_t h = - /* v0,v1 and v1,v0 must give the same hash!!! */ - v0 + v1 + (int64_t)MMIN(v0, v1); - 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->common_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->common_edge_rank = ee; - } + + /* 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))); + OK2(darray_neighbourhood_reserve(&neighbourhood_by_edge, nbedges_guess)); + OK2(htable_edge_id_reserve(&edge_ids, nbedges_guess)); + + /* Loop on triangles to register edges. + * All threads considering all the edges and processing some */ + FOR_EACH(t, 0, scn->nutris) { + struct trg_edge edge; + unsigned char ee; + FOR_EACH(ee, 0, 3) { + edge_id_t* p_id; + size_t n_sz; + struct edge_neighbourhood* neighbourhood; + struct neighbour_info* info; + 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 int64_t h = + /* v0,v1 and v1,v0 must give the same hash!!! */ + v0 + v1 + (int64_t)MMIN(v0, v1); + 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) { + neighbourhood = + darray_neighbourhood_data_get(&neighbourhood_by_edge) + *p_id; + ASSERT(n->edge.vrtx0 == edge.vrtx0 && n->edge.vrtx1 == edge.vrtx1); + } else { + /* Create id */ + edge_id_t id; + 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)); + OK2(htable_edge_id_set(&edge_ids, &edge, &id)); + OK2(darray_neighbourhood_resize(&neighbourhood_by_edge, 1 + sz)); + neighbourhood = darray_neighbourhood_data_get(&neighbourhood_by_edge) + sz; + /* Add neighbour info to a newly created edge's neighbour list */ + neighbourhood->edge = edge; + ASSERT(darray_neighbour_size_get(&neighbourhood->neighbours) == 0); + /* Just a guess: few edges will have less than 2 neighbours */ + OK2(darray_neighbour_reserve(&neighbourhood->neighbours, 2)); } + /* Add neighbour info to neighbourhood */ + n_sz = darray_neighbour_size_get(&neighbourhood->neighbours); + OK2(darray_neighbour_resize(&neighbourhood->neighbours, 1 + n_sz)); + info = darray_neighbour_data_get(&neighbourhood->neighbours) + n_sz; + info->trg_id = t; + info->common_edge_rank = ee; + } + } /* No barrier here. */ + + /* Loop on collected edges. + * For each edge sort triangle sides by rotation angle + * and connect neighbours. */ + sz = darray_neighbourhood_size_get(&neighbourhood_by_edge); + ASSERT(sz <= EDGE_MAX__); + edge_count = (edge_id_t)sz; + FOR_EACH(e, 0, edge_count) { + double edge[3], common_edge[3], basis[9], norm, max_z, maxz_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; + unsigned char maxz_vrank, maxz_vrank_edge; + sz = darray_neighbour_size_get(neighbour_list); + ASSERT(sz <= SIDE_MAX__); + neighbour_count = (side_id_t)sz; + 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) { + maxz_edge = vertices[v0].pos.z; + maxz_vrank_edge = 0; + } else { + maxz_edge = vertices[v1].pos.z; + maxz_vrank_edge = 1; } - /* No implicit barrier here. */ - - /* 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, max_z, maxz_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; - unsigned char maxz_vrank, maxz_vrank_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) { - maxz_edge = vertices[v0].pos.z; - maxz_vrank_edge = 0; + norm = d3_normalize(common_edge, common_edge); + ASSERT(norm); (void)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; + unsigned char actual_vrank; + v2 = trg_in->vertice_id[(neighbour_info->common_edge_rank + 2) % 3]; + if(vertices[v2].pos.z > maxz_edge) { + max_z = vertices[v2].pos.z; + maxz_vrank = 2; } else { - maxz_edge = vertices[v1].pos.z; - maxz_vrank_edge = 1; + max_z = maxz_edge; + maxz_vrank = maxz_vrank_edge; } - norm = d3_normalize(common_edge, common_edge); - ASSERT(norm); (void)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; - unsigned char actual_vrank; - v2 = trg_in->vertice_id[(neighbour_info->common_edge_rank + 2) % 3]; - if(vertices[v2].pos.z > maxz_edge) { - max_z = vertices[v2].pos.z; - maxz_vrank = 2; - } else { - max_z = maxz_edge; - maxz_vrank = maxz_vrank_edge; - } - /* Compute the actual vertex id - * as vertices are not in the v0 v1 v2 order in the actual triangle */ - if(maxz_vrank == 2) { - actual_vrank = - (unsigned char)((neighbour_info->common_edge_rank + 2) % 3); - } else { - unsigned char is_r = - neighbour->reversed_edge[neighbour_info->common_edge_rank]; - ASSERT(maxz_vrank == 0 || maxz_vrank == 1); - actual_vrank = (unsigned char)((is_r ? 1 - maxz_vrank : maxz_vrank) - + neighbour_info->common_edge_rank) % 3; - } - ASSERT(neighbour->max_z <= max_z); - neighbour->max_z = max_z; - ASSERT(actual_vrank <= 2); - neighbour->max_z_vrtx_rank = actual_vrank; - /* 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 unsigned char crt_edge = current->common_edge_rank; - const unsigned char ccw_edge = ccw_neighbour->common_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; + /* Compute the actual vertex id + * as vertices are not in the v0 v1 v2 order in the actual triangle */ + if(maxz_vrank == 2) { + actual_vrank = + (unsigned char)((neighbour_info->common_edge_rank + 2) % 3); + } else { + unsigned char is_r = + neighbour->reversed_edge[neighbour_info->common_edge_rank]; + ASSERT(maxz_vrank == 0 || maxz_vrank == 1); + actual_vrank = (unsigned char)((is_r ? 1 - maxz_vrank : maxz_vrank) + + neighbour_info->common_edge_rank) % 3; } + ASSERT(neighbour->max_z <= max_z); + neighbour->max_z = max_z; + ASSERT(actual_vrank <= 2); + neighbour->max_z_vrtx_rank = actual_vrank; + /* 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 unsigned char crt_edge = current->common_edge_rank; + const unsigned char ccw_edge = ccw_neighbour->common_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); } - - return res; +tmp_error: + if(tmp_res != RES_OK) *res = tmp_res; + /* Threads are allowed to return whitout sync. */ + htable_edge_id_release(&edge_ids); + darray_neighbourhood_release(&neighbourhood_by_edge); } -static res_T +static void build_result (struct senc_descriptor* desc, const struct darray_ptr_component_descriptor* connex_components, - const struct darray_triangle_comp* triangles_comp_array) + const struct darray_triangle_comp* triangles_comp_array, + /* Shared error status. + * We accept to overwritte an error with a different error */ + res_T* res) { - res_T res = RES_OK; + /* This function is called from an omp parallel block and executed + * concurrently. */ struct mem_allocator* alloc; struct cc_descriptor* const* cc_descriptors; struct enclosure_data* enclosures; const struct triangle_in* triangles_in; struct triangle_enc* triangles_enc; const struct triangle_comp* triangles_comp; - volatile int exit_for = 0; + struct htable_vrtx_id vtable; + int64_t tt; + int64_t ee; ASSERT(desc && connex_components && triangles_comp_array); alloc = descriptor_get_allocator(desc); - ASSERT(darray_ptr_component_descriptor_size_get(connex_components) < COMPONENT_MAX__); + ASSERT(darray_ptr_component_descriptor_size_get(connex_components) + < COMPONENT_MAX__); cc_descriptors = darray_ptr_component_descriptor_cdata_get(connex_components); enclosures = darray_enclosure_data_get(&desc->enclosures); triangles_in = darray_triangle_in_cdata_get(&desc->scene->triangles_in); triangles_comp = darray_triangle_comp_cdata_get(triangles_comp_array); - OK2(darray_triangle_enc_resize(&desc->triangles_enc, desc->scene->nutris), - error_); - triangles_enc = darray_triangle_enc_data_get(&desc->triangles_enc); - #pragma omp parallel num_threads(desc->scene->dev->nthreads) + #pragma omp single { - struct htable_vrtx_id vtable; - int64_t tt; - int64_t ee; - - /* Build global enclosure information */ - #pragma omp for - for(tt = 0; tt < (int64_t) desc->scene->nutris; tt++) { - trg_id_t t = (trg_id_t)tt; - const component_id_t cf_id = triangles_comp[t].component[SIDE_FRONT]; - const component_id_t cb_id = triangles_comp[t].component[SIDE_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(triangles_enc[t].enclosure[SIDE_FRONT] == ENCLOSURE_NULL__); - triangles_enc[t].enclosure[SIDE_FRONT] = ef_id; - ASSERT(triangles_enc[t].enclosure[SIDE_BACK] == ENCLOSURE_NULL__); - triangles_enc[t].enclosure[SIDE_BACK] = eb_id; - } - /* Implicit barrier here */ - - /* Resize/push operations on enclosure's fields are valid in the - * openmp block as a given enclosure is processed by a single thread */ - htable_vrtx_id_init(alloc, &vtable); - - ASSERT(desc->enclosures_count <= ENCLOSURE_MAX__); - #pragma omp for schedule(dynamic) nowait - for(ee = 0; ee < (int64_t)desc->enclosures_count; ee++) { - const enclosure_id_t e = (enclosure_id_t)ee; - struct enclosure_data* enc = enclosures + e; - const struct cc_descriptor* current = cc_descriptors[enc->first_component]; - trg_id_t fst_idx = 0; - trg_id_t sgd_idx = enc->side_count; - trg_id_t t; - ASSERT(enc->first_component - < darray_ptr_component_descriptor_size_get(connex_components)); - ASSERT(current->cc_id == enc->first_component); - - if(exit_for) continue; - ASSERT(e <= UINT_MAX); - enc->header.enclosure_id = (unsigned)e; /* Back to API type */ - ASSERT(current->enclosure_id == enc->header.enclosure_id); - enc->header.is_infinite = (e == 0); - enc->header.enclosed_medium - = (unsigned)current->medium; /* Back to API type */ - ASSERT(enc->header.enclosed_medium < desc->scene->nmeds); - - /* Build side and vertex lists. */ - OK(darray_triangle_in_resize(&enc->sides, enc->side_count)); - /* Size is just a int */ - OK(darray_vrtx_id_reserve(&enc->vertices, - (size_t)(enc->side_count * 0.6))); - /* New vertex numbering scheme local to the enclosure */ - htable_vrtx_id_clear(&vtable); - ASSERT(desc->scene->nutris - == darray_triangle_in_size_get(&desc->scene->triangles_in)); - /* Put at the end the back-faces of triangles that also have their - * front-face in the list. */ - for(t = TRGSIDE_2_TRG(enc->side_range.first); - t <= TRGSIDE_2_TRG(enc->side_range.last); - t++) - { - const struct triangle_in* trg_in = triangles_in + t; - struct triangle_in* trg; - unsigned vertice_id[3]; - int i; - if(triangles_enc[t].enclosure[SIDE_FRONT] != e - && triangles_enc[t].enclosure[SIDE_BACK] != e) - continue; - ++enc->header.unique_triangle_count; + res_T tmp_res = + darray_triangle_enc_resize(&desc->triangles_enc, desc->scene->nutris); + if(tmp_res != RES_OK) *res = tmp_res; + }/* Implicit barrier here. */ + if(*res != RES_OK) return; + triangles_enc = darray_triangle_enc_data_get(&desc->triangles_enc); - FOR_EACH(i, 0, 3) { - vrtx_id_t* id = htable_vrtx_id_find(&vtable, trg_in->vertice_id + i); - if(id) { - vertice_id[i] = *id; /* Known vertex */ - } else { - /* Create new association */ - 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; - OK(htable_vrtx_id_set(&vtable, trg_in->vertice_id + i, - vertice_id + i)); - OK(darray_vrtx_id_push_back(&enc->vertices, trg_in->vertice_id + i)); - ++enc->header.vertices_count; - } - } - ASSERT(triangles_enc[t].enclosure[SIDE_FRONT] == e - || triangles_enc[t].enclosure[SIDE_BACK] == e); - if(triangles_enc[t].enclosure[SIDE_FRONT] == e) { - ++enc->header.triangle_count; - trg = darray_triangle_in_data_get(&enc->sides) + fst_idx++; - - FOR_EACH(i, 0, 2) trg->medium[i] = trg_in->medium[i]; - trg->global_id = trg_in->global_id; - FOR_EACH(i, 0, 3) trg->vertice_id[i] = vertice_id[i]; - } - if(triangles_enc[t].enclosure[SIDE_BACK] == e) { - ++enc->header.triangle_count; - trg = darray_triangle_in_data_get(&enc->sides) + - ((triangles_enc[t].enclosure[SIDE_FRONT] == e) ? --sgd_idx : fst_idx++); - - FOR_EACH(i, 0, 2) trg->medium[i] = trg_in->medium[1 - i]; - trg->global_id = trg_in->global_id; - FOR_EACH(i, 0, 3) trg->vertice_id[i] = vertice_id[2 - i]; + /* Build global enclosure information */ + #pragma omp for + for(tt = 0; tt < (int64_t)desc->scene->nutris; tt++) { + trg_id_t t = (trg_id_t)tt; + const component_id_t cf_id = triangles_comp[t].component[SIDE_FRONT]; + const component_id_t cb_id = triangles_comp[t].component[SIDE_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(triangles_enc[t].enclosure[SIDE_FRONT] == ENCLOSURE_NULL__); + triangles_enc[t].enclosure[SIDE_FRONT] = ef_id; + ASSERT(triangles_enc[t].enclosure[SIDE_BACK] == ENCLOSURE_NULL__); + triangles_enc[t].enclosure[SIDE_BACK] = eb_id; + } + /* Implicit barrier here */ + + /* Resize/push operations on enclosure's fields are valid in the + * openmp block as a given enclosure is processed by a single thread */ + htable_vrtx_id_init(alloc, &vtable); + + ASSERT(desc->enclosures_count <= ENCLOSURE_MAX__); + #pragma omp for schedule(dynamic) nowait + for(ee = 0; ee < (int64_t)desc->enclosures_count; ee++) { + const enclosure_id_t e = (enclosure_id_t)ee; + struct enclosure_data* enc = enclosures + e; + const struct cc_descriptor* current = cc_descriptors[enc->first_component]; + trg_id_t fst_idx = 0; + trg_id_t sgd_idx = enc->side_count; + trg_id_t t; + res_T tmp_res = RES_OK; + ASSERT(enc->first_component + < darray_ptr_component_descriptor_size_get(connex_components)); + ASSERT(current->cc_id == enc->first_component); + + if(*res != RES_OK) continue; + ASSERT(e <= UINT_MAX); + enc->header.enclosure_id = (unsigned)e; /* Back to API type */ + ASSERT(current->enclosure_id == enc->header.enclosure_id); + enc->header.is_infinite = (e == 0); + enc->header.enclosed_medium + = (unsigned)current->medium; /* Back to API type */ + ASSERT(enc->header.enclosed_medium < desc->scene->nmeds); + + /* Build side and vertex lists. */ + OK2(darray_triangle_in_resize(&enc->sides, enc->side_count)); + /* Size is just a int */ + OK2(darray_vrtx_id_reserve(&enc->vertices, + (size_t)(enc->side_count * 0.6))); + /* New vertex numbering scheme local to the enclosure */ + htable_vrtx_id_clear(&vtable); + ASSERT(desc->scene->nutris + == darray_triangle_in_size_get(&desc->scene->triangles_in)); + /* Put at the end the back-faces of triangles that also have their + * front-face in the list. */ + for(t = TRGSIDE_2_TRG(enc->side_range.first); + t <= TRGSIDE_2_TRG(enc->side_range.last); + t++) + { + const struct triangle_in* trg_in = triangles_in + t; + struct triangle_in* trg; + unsigned vertice_id[3]; + int i; + if(triangles_enc[t].enclosure[SIDE_FRONT] != e + && triangles_enc[t].enclosure[SIDE_BACK] != e) + continue; + ++enc->header.unique_triangle_count; + + FOR_EACH(i, 0, 3) { + vrtx_id_t* p_id = htable_vrtx_id_find(&vtable, trg_in->vertice_id + i); + if(p_id) { + vertice_id[i] = *p_id; /* Known vertex */ + } else { + /* Create new association */ + 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; + OK2(htable_vrtx_id_set(&vtable, trg_in->vertice_id + i, + vertice_id + i)); + OK2(darray_vrtx_id_push_back(&enc->vertices, trg_in->vertice_id + i)); + ++enc->header.vertices_count; } - if(fst_idx == sgd_idx) break; } - continue; - error: - /* Cannot goto out of openmp block */ - exit_for = 1; - continue; + ASSERT(triangles_enc[t].enclosure[SIDE_FRONT] == e + || triangles_enc[t].enclosure[SIDE_BACK] == e); + if(triangles_enc[t].enclosure[SIDE_FRONT] == e) { + ++enc->header.triangle_count; + trg = darray_triangle_in_data_get(&enc->sides) + fst_idx++; + + FOR_EACH(i, 0, 2) trg->medium[i] = trg_in->medium[i]; + trg->global_id = trg_in->global_id; + FOR_EACH(i, 0, 3) trg->vertice_id[i] = vertice_id[i]; + } + if(triangles_enc[t].enclosure[SIDE_BACK] == e) { + ++enc->header.triangle_count; + trg = darray_triangle_in_data_get(&enc->sides) + + ((triangles_enc[t].enclosure[SIDE_FRONT] == e) ? --sgd_idx : fst_idx++); + + FOR_EACH(i, 0, 2) trg->medium[i] = trg_in->medium[1 - i]; + trg->global_id = trg_in->global_id; + FOR_EACH(i, 0, 3) trg->vertice_id[i] = vertice_id[2 - i]; + } + if(fst_idx == sgd_idx) break; } - htable_vrtx_id_release(&vtable); - } - OK2(res, error_); - -exit: - return res; -error_: - goto exit; + continue; + tmp_error: + ASSERT(tmp_res != RES_OK); + *res = tmp_res; + } /* No barrier here */ + htable_vrtx_id_release(&vtable); } /******************************************************************************* @@ -1112,7 +1116,6 @@ error_: res_T senc_scene_analyze(struct senc_scene* scn, struct senc_descriptor** out_desc) { - res_T res = RES_OK; struct senc_descriptor* desc = NULL; /* By triangle tmp data */ struct darray_triangle_tmp triangles_tmp; @@ -1121,15 +1124,18 @@ senc_scene_analyze(struct senc_scene* scn, struct senc_descriptor** out_desc) * They are refered to by arrays of ids. */ struct darray_ptr_component_descriptor connex_components; char connex_components_initialized = 0; - /* Triangle neighbourhood by edge. */ - struct darray_neighbourhood neighbourhood_by_edge; - char neighbourhood_by_edge_initialized = 0; /* Store by-triangle components */ struct darray_triangle_comp triangles_comp; char triangles_comp_initialized = 0; /* Array of triangle sides. */ struct trgside* trgsides = NULL; struct s3d_scene_view* s3d_view = NULL; + /* Atomic counters to share beetwen threads */ + ATOMIC component_count = 0; + ATOMIC next_enclosure_id = 1; + /* Used as args to have shared vars between threads in functions */ + struct cc_descriptor* infinity_first_cc = NULL; + res_T res = RES_OK; if(!scn || !out_desc) return RES_BAD_ARG; @@ -1152,66 +1158,139 @@ senc_scene_analyze(struct senc_scene* scn, struct senc_descriptor** out_desc) goto error; } - /* Step 1: build neighbourhoods */ - res = collect_and_link_neighbours(scn, trgsides, &triangles_tmp); + /* The end of the analyze is multithreaded */ + ASSERT(scn->dev->nthreads > 0); + #pragma omp parallel num_threads(scn->dev->nthreads) + { + /* Step 1: build neighbourhoods */ + collect_and_link_neighbours(scn, trgsides, &triangles_tmp, &res); + /* No barrier at the end of step 1: data used in step 1 cannot be + * released / data produced by step 1 cannot be used + * until next sync point */ - if(res != RES_OK) { - log_err(scn->dev, - "%s: could not build neighbourhoods from scene.\n", FUNC_NAME); - goto error; - } + if(res != RES_OK) { + { + #pragma omp single nowait + { + log_err(scn->dev, + "%s: could not build neighbourhoods from scene.\n", FUNC_NAME); + } /* No barrier here */ + } + } + if(res != RES_OK) goto error_; + + /* The first thread here allocates some data. + * Barrier needed at the end to ensure data created before any use. */ + #pragma omp single + { + res_T tmp_res = RES_OK; + darray_ptr_component_descriptor_init(scn->dev->allocator, + &connex_components); + connex_components_initialized = 1; + /* Just a hint; to limit contention */ + OK2(darray_ptr_component_descriptor_reserve(&connex_components, + 2 * scn->nmeds)); + darray_triangle_comp_init(scn->dev->allocator, &triangles_comp); + triangles_comp_initialized = 1; + OK2(darray_triangle_comp_resize(&triangles_comp, scn->nutris)); + tmp_error: + if(tmp_res != RES_OK) res = tmp_res; + } + /* Implicit barrier here: constraints on step 1 data are now met */ + if(res != RES_OK) goto error_; + + /* Step 2: extract triangle connex components */ + extract_connex_components(desc, trgsides, &connex_components, + &triangles_tmp, &triangles_comp, &s3d_view, &component_count, &res); + /* No barrier at the end of step 2: data used in step 2 cannot be + * released / data produced by step 2 cannot be used + * until next sync point */ + + if(res != RES_OK) { + #pragma omp single nowait + { + log_err(scn->dev, + "%s: could not extract connex components from scene.\n", FUNC_NAME); + } /* No barrier here */ + } + if(res != RES_OK) goto error_; - darray_ptr_component_descriptor_init(scn->dev->allocator, &connex_components); - connex_components_initialized = 1; - darray_triangle_comp_init(scn->dev->allocator, &triangles_comp); - triangles_comp_initialized = 1; - OK(darray_triangle_comp_resize(&triangles_comp, scn->nutris)); - - /* Step 2: extract triangle connex components */ - res = extract_connex_components(desc, trgsides, &connex_components, - &triangles_tmp, &triangles_comp, &s3d_view); - if(res != RES_OK) { - log_err(scn->dev, - "%s: could not extract connex components from scene.\n", FUNC_NAME); - goto error; - } + #pragma omp barrier + /* Constraints on step 2 data are now met */ - darray_triangle_tmp_release(&triangles_tmp); - triangles_tmp_initialized = 0; + /* One thread releases some data before going to step 3, + * the others go to step 3 without sync */ + #pragma omp single nowait + { + darray_triangle_tmp_release(&triangles_tmp); + triangles_tmp_initialized = 0; + } /* No barrier here */ + + /* Step 3: group components */ + group_connex_components(desc, trgsides, &triangles_comp, + &connex_components, s3d_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 */ + + /* One thread releases some data before going to step 4, + * the others go to step 4 without sync */ + #pragma omp single nowait + { + if(s3d_view) S3D(scene_view_ref_put(s3d_view)); + s3d_view = NULL; + } /* No barrier here */ - /* Step 3: group components */ - res = group_connex_components(desc, trgsides, &triangles_comp, - &connex_components, s3d_view); - if (s3d_view) S3D(scene_view_ref_put(s3d_view)); - if(res != RES_OK) { - log_err(scn->dev, - "%s: could not group connex components from scene.\n", FUNC_NAME); - goto error; - } + if(res != RES_OK) { + #pragma omp single nowait + { + log_err(scn->dev, + "%s: could not group connex components from scene.\n", FUNC_NAME); + } /* No barrier here */ + } + if(res != RES_OK) goto error_; - /* Build result. */ - res = build_result(desc, &connex_components, &triangles_comp); - if(res != RES_OK) { - log_err(scn->dev, "%s: could not build result.\n", FUNC_NAME); - goto error; - } + /* Step 4: Build result */ + build_result(desc, &connex_components, &triangles_comp, &res); + /* No barrier at the end of step 4: data used in step 4 cannot be + * released / data produced by step 4 cannot be used + * until next sync point */ + + if(res != RES_OK) { + #pragma omp single nowait + { + log_err(scn->dev, "%s: could not build result.\n", FUNC_NAME); + } /* No barrier here */ + } + if(res != RES_OK) goto error_; + + #pragma omp barrier + /* Constraints on step 4 data are now met */ - darray_triangle_comp_release(&triangles_comp); - triangles_comp_initialized = 0; + /* Some threads release data */ + #pragma omp sections nowait + { + #pragma omp section + { + ASSERT(connex_components_initialized); + custom_darray_ptr_component_descriptor_release(&connex_components); + connex_components_initialized = 0; + } + #pragma omp section + { + darray_triangle_comp_release(&triangles_comp); + triangles_comp_initialized = 0; + } + } /* No barrier here */ +error_: + ; + } /* Implicit barrier here */ + if(res != RES_OK) goto error; exit: - if(connex_components_initialized) { - size_t c, cc_count = - darray_ptr_component_descriptor_size_get(&connex_components); - struct cc_descriptor** components = - darray_ptr_component_descriptor_data_get(&connex_components); - FOR_EACH(c, 0, cc_count) { - ptr_component_descriptor_release(scn->dev->allocator, components + c); - } - darray_ptr_component_descriptor_release(&connex_components); - } - if(neighbourhood_by_edge_initialized) - darray_neighbourhood_release(&neighbourhood_by_edge); + if(connex_components_initialized) + custom_darray_ptr_component_descriptor_release(&connex_components); + if(s3d_view) S3D(scene_view_ref_put(s3d_view)); if(triangles_tmp_initialized) darray_triangle_tmp_release(&triangles_tmp); if(triangles_comp_initialized) darray_triangle_comp_release(&triangles_comp); if(trgsides) MEM_RM(scn->dev->allocator, trgsides); diff --git a/src/senc_scene_analyze_c.h b/src/senc_scene_analyze_c.h @@ -146,20 +146,27 @@ ptr_component_descriptor_init ASSERT(data); *data = NULL; } -static FINLINE void -ptr_component_descriptor_release - (struct mem_allocator* allocator, - struct cc_descriptor** data) -{ - ASSERT(allocator && data); - MEM_RM(allocator, *data); -} #define DARRAY_NAME ptr_component_descriptor #define DARRAY_DATA struct cc_descriptor* #define DARRAY_FUNCTOR_INIT ptr_component_descriptor_init #include <rsys/dynamic_array.h> +/* Need allocator to free array elts: cannot rely on standard + * darray release stuff */ +static FINLINE void +custom_darray_ptr_component_descriptor_release + (struct darray_ptr_component_descriptor* array) +{ + size_t c, cc_count; + struct cc_descriptor** components; + if(!array) return; + cc_count = darray_ptr_component_descriptor_size_get(array); + components = darray_ptr_component_descriptor_data_get(array); + FOR_EACH(c, 0, cc_count) MEM_RM(array->allocator, components[c]); + darray_ptr_component_descriptor_release(array); +} + /* Triangle information. * Depending on lifespan, information is kept in different places: * - triangle_in for user provided information (kept in scene)