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 428f3d13e76bec797dbbf9624ea8ed7513d1cbff
parent 695c5ff90841d5493b64f0331bdd634d80ef0201
Author: Christophe Coustet <christophe.coustet@meso-star.com>
Date:   Fri, 14 Aug 2020 10:58:09 +0200

BugFix: grouping of components in contact.

The recent algorithm based on closest point queries introduced some
regression on this. But notice that for components in contact at a
single vertex, no previous version has ever been correct.

Diffstat:
Msrc/senc3d_scene_analyze.c | 452++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++-----------------
1 file changed, 359 insertions(+), 93 deletions(-)

diff --git a/src/senc3d_scene_analyze.c b/src/senc3d_scene_analyze.c @@ -53,8 +53,14 @@ const struct cc_descriptor CC_DESCRIPTOR_NULL = CC_DESCRIPTOR_NULL__; #define DARRAY_DATA component_id_t #include <rsys/dynamic_array.h> -struct filter_ctx { +#define HTABLE_NAME component_id +#define HTABLE_KEY component_id_t +#define HTABLE_DATA char +#include <rsys/hash_table.h> + +struct filter_ctx1 { struct senc3d_scene* scn; + struct s3d_scene_view* view; component_id_t origin_component; struct darray_triangle_comp* triangles_comp; struct darray_ptr_component_descriptor* components; @@ -62,6 +68,25 @@ struct filter_ctx { component_id_t hit_component; }; +struct filter_ctx2 { + struct darray_triangle_comp* triangles_comp; + int cpt; + component_id_t component; +}; + +enum fctx_type { + FCTX1, + FCTX2 +}; + +struct filter_ctx { + enum fctx_type type; + union { + struct filter_ctx1 ctx1; + struct filter_ctx2 ctx2; + } c; +}; + #define HTABLE_NAME overlap #define HTABLE_KEY trg_id_t #define HTABLE_DATA char @@ -131,79 +156,201 @@ self_hit_filter void* ray_data, void* filter_data) { - struct filter_ctx* filter_ctx = ray_data; - const struct triangle_comp* hit_trg_comp; - struct s3d_attrib pos; - float s, dir[3]; + struct filter_ctx* fctx_ = ray_data; + struct filter_ctx1* fctx; + const struct triangle_comp* components; + const component_id_t* hit_comp; + float s = 0; enum senc3d_side hit_side; + int i; + double mz = -INF; + const struct triangle_in* triangles; + const struct triangle_in* trg = NULL; + const union double3* vertices; (void)ray_org; (void)ray_dir; (void)filter_data; - ASSERT(hit && filter_ctx); + ASSERT(hit && fctx_); + + if(fctx_->type == FCTX2) { + /* The filter is used to count the hits on some component along + * an infinite ray */ + struct filter_ctx2* ctx2 = &fctx_->c.ctx2; + ASSERT(hit->prim.prim_id + < darray_triangle_comp_size_get(ctx2->triangles_comp)); + components = darray_triangle_comp_cdata_get(ctx2->triangles_comp); + hit_comp = components[hit->prim.prim_id].component; + if(hit_comp[SENC3D_FRONT] == ctx2->component + || hit_comp[SENC3D_BACK] == ctx2->component) + ctx2->cpt++; + return 1; /* Reject to continue counting */ + } + + /* The filter is called from a point query on successive hits found from + * ray_org, that belongs to origin_component. It can keep or reject the hit. + * Hits are only submitted inside a certain radius from ray_org, that is + * decreased to the hit distance for every hit that is kept. + * At the end, the last kept hit (= the closest), determines a component to + * which origin_component is linked. At a later stage the algorithm proceeds + * linked components to determine their relative inclusions. + * + * For each hit, the filter computes if the hit is on a component above + * origin_component (that is with >= Z). + * If the hit is distant (dist>0), we just keep the hit as a valid candidate, + * but things get more tricky when dist==0 (ray_org is a vertex where some + * other components are in contact with origin_component). + * In this case, one of the other components can include the origin_component + * (greater volume needed), or they can be disjoint, with (at least) ray_org + * as a common vertex (they can also partially intersect, but this is invalid + * and remains undetected by star enclosures). */ + ASSERT(fctx_->type == FCTX1); + fctx = &fctx_->c.ctx1; + components = darray_triangle_comp_cdata_get(fctx->triangles_comp); + hit_comp = components[hit->prim.prim_id].component; + triangles = darray_triangle_in_cdata_get(&fctx->scn->triangles_in); + vertices = darray_position_cdata_get(&fctx->scn->vertices); ASSERT(hit->prim.prim_id - < darray_triangle_comp_size_get(filter_ctx->triangles_comp)); + < darray_triangle_comp_size_get(fctx->triangles_comp)); ASSERT(hit->uv[0] == CLAMP(hit->uv[0], 0, 1)); ASSERT(hit->uv[1] == CLAMP(hit->uv[1], 0, 1)); - hit_trg_comp = darray_triangle_comp_cdata_get(filter_ctx->triangles_comp) - + hit->prim.prim_id; + /* No self hit */ - if(hit_trg_comp->component[SENC3D_FRONT] == filter_ctx->origin_component - || hit_trg_comp->component[SENC3D_BACK] == filter_ctx->origin_component) + if(hit_comp[SENC3D_FRONT] == fctx->origin_component + || hit_comp[SENC3D_BACK] == fctx->origin_component) return 1; /* Reject */ if(hit->distance == 0) { - /* dir = + z; s = dir . n; */ - s = hit->normal[2]; - } else { - CHK(s3d_primitive_get_attrib(&hit->prim, S3D_POSITION, hit->uv, &pos) - == RES_OK); - /* Cannot go downwards */ - if(pos.value[2] <= ray_org[2]) - return 1; /* Reject */ - /* In closest_point queries ray_dir is not informed */ - f3_sub(dir, pos.value, ray_org); - s = f3_dot(dir, hit->normal); + /* origin_component is in contact with some other components + * We will need further exploration to know if they should be considered */ + struct cc_descriptor* const* descriptors + = darray_ptr_component_descriptor_cdata_get(fctx->components); + int n; + + /* If same component, process only once */ + FOR_EACH(n, 0, (hit_comp[SENC3D_FRONT] == hit_comp[SENC3D_BACK] ? 1 : 2)) { + const enum senc3d_side sides[2] = { SENC3D_FRONT, SENC3D_BACK }; + component_id_t c = hit_comp[sides[n]]; + side_id_t side; + double pt[3] = { 0, 0, 0 }; + float org[3], dir[3] = { 0, 0, 1 }, rg[2] = { 0, FLT_MAX }; + struct s3d_hit hit2 = S3D_HIT_NULL; + struct filter_ctx fctx2; + ASSERT(c < darray_ptr_component_descriptor_size_get(fctx->components)); + if(descriptors[c]->is_outer_border) { + if(fabs(descriptors[c]->_6volume) + <= fabs(descriptors[fctx->origin_component]->_6volume)) + /* Component is not large enough to include origin_component */ + continue; + } else { + vrtx_id_t c_z_id = descriptors[c]->max_z_vrtx_id; + vrtx_id_t o_z_id = descriptors[fctx->origin_component]->max_z_vrtx_id; + ASSERT(c_z_id < darray_position_size_get(&fctx->scn->vertices)); + ASSERT(o_z_id < darray_position_size_get(&fctx->scn->vertices)); + if(vertices[c_z_id].pos.z <= vertices[o_z_id].pos.z) + /* Component is not above origin_component */ + continue; + } + /* Check if component c includes origin_component; as components cannot + * overlap, testing a single point is OK, as long as the point is not on + * c boundary (can be on origin_component boundary, though). + * As this case is supposed to be rare, we go for a basic algorithm */ + for(side = descriptors[fctx->origin_component]->side_range.first; + side <= descriptors[fctx->origin_component]->side_range.last; + side++) + { + /* Find a triangle on origin_component boundary that is not on c + * boundary (the 2 components cannot share all their triangles) */ + trg_id_t t = TRGSIDE_2_TRG(side); + const component_id_t* candidate_comp = components[t].component; + if(candidate_comp[SENC3D_FRONT] != fctx->origin_component + && candidate_comp[SENC3D_BACK] != fctx->origin_component) + continue; + if(candidate_comp[SENC3D_FRONT] == c + || candidate_comp[SENC3D_BACK] == c) + continue; + /* This triangle is OK */ + trg = triangles + t; + break; + } + ASSERT(trg != NULL); + /* Any point on trg not on an edge can do the trick: use the barycenter */ + FOR_EACH(i, 0, 3) { + vrtx_id_t v = trg->vertice_id[i]; + ASSERT(v < darray_position_size_get(&fctx->scn->vertices)); + d3_add(pt, pt, vertices[v].vec); + } + d3_divd(pt, pt, 3); + f3_set_d3(org, pt); + /* Trace a ray and count intersections with components of trg */ + fctx2.type = FCTX2; + fctx2.c.ctx2.triangles_comp = fctx->triangles_comp; + fctx2.c.ctx2.cpt = 0; + fctx2.c.ctx2.component = c; + S3D(scene_view_trace_ray(fctx->view, org, dir, rg, &fctx2, &hit2)); + ASSERT(S3D_HIT_NONE(&hit2)); /* The ray is supposed to go to infinity */ + /* origin_component is linked_to an outer component if cpt is odd, + * linked_to an inner component if cpt is even */ + if(descriptors[c]->is_outer_border == (fctx2.c.ctx2.cpt % 2)) { + fctx->hit_component = fctx2.c.ctx2.component; + return 0 /* origin_component is linked to c: keep */; + } + } + return 1; /* Reject */ + } + + /* Reject hits with < Z */ + ASSERT(hit->prim.prim_id < + darray_triangle_in_size_get(&fctx->scn->triangles_in)); + trg = triangles + hit->prim.prim_id; + /* Check if the hit is above ray_org (hit[2] > ray_org[2]) + * As we cannot rely on numerical accuracy when computing hit positions, + * we use the triangle vertices to check if some part of the hit triangle + * is above ray_org */ + FOR_EACH(i, 0, 3) { + vrtx_id_t v = trg->vertice_id[i]; + ASSERT(v < darray_position_size_get(&fctx->scn->vertices)); + const union double3* p = vertices + v; + if(i == 0 || mz < p->pos.z) mz = p->pos.z; } + if(mz <= ray_org[2]) + return 1; /* Hit triangle is below ray_org: reject */ - if(s == 0) return 1; /* Reject */ + if(hit_comp[SENC3D_FRONT] == hit_comp[SENC3D_BACK]) { + /* Easy case and hit component is known */ + fctx->hit_component = hit_comp[SENC3D_FRONT]; + return 0; /* Keep */ + } + + /* Compute hit side by using the vertex with best accuracy */ + FOR_EACH(i, 0, 3) { + float tmps, tmp[3], dir[3]; + vrtx_id_t v = trg->vertice_id[i]; + ASSERT(v < darray_position_size_get(&fctx->scn->vertices)); + const union double3* p = vertices + v; + f3_sub(dir, f3_set_d3(tmp, p->vec), ray_org); + tmps = f3_dot(dir, hit->normal); + if(i == 0 || fabsf(s) < fabsf(tmps)) s = tmps; + } + + /* We cannot know which side to consider if s==0. + * As hit distance > 0 and the 2 sides belong to 2 different components, + * another candidate must selected afterwards (can be at greater distance, + * disallowing to restrict the search distance here) */ + if(s == 0) { + fctx->hit_component = COMPONENT_NULL__; + return 1; /* Reject */ + } /* Determine which side was hit */ hit_side = ((s < 0) /* Facing geometrical normal of hit */ - == ((filter_ctx->scn->convention & SENC3D_CONVENTION_NORMAL_FRONT) != 0)) + == ((fctx->scn->convention & SENC3D_CONVENTION_NORMAL_FRONT) != 0)) /* Warning: following Embree 2 convention for geometrical normals, * the Star3D hit normal is left-handed while star-enclosures-3d uses * right-handed convention */ ? SENC3D_BACK : SENC3D_FRONT; - filter_ctx->hit_component = hit_trg_comp->component[hit_side]; + fctx->hit_component = hit_comp[hit_side]; - if(hit->distance == 0) { - /* To avoid create loops, allow dist=0 only if the hit component has - * a higher max_z than the origin component */ - struct cc_descriptor* const* descriptors - = darray_ptr_component_descriptor_cdata_get(filter_ctx->components); - const struct cc_descriptor* origin_desc; - const struct cc_descriptor* hit_desc; - const union double3* positions - = darray_position_cdata_get(&filter_ctx->scn->vertices); - ASSERT(filter_ctx->origin_component - < darray_ptr_component_descriptor_size_get(filter_ctx->components)); - ASSERT(filter_ctx->hit_component - < darray_ptr_component_descriptor_size_get(filter_ctx->components)); - origin_desc = descriptors[filter_ctx->origin_component]; - hit_desc = descriptors[filter_ctx->hit_component]; (void)hit_desc; - ASSERT(origin_desc->max_z_vrtx_id - < darray_position_size_get(&filter_ctx->scn->vertices)); - ASSERT(hit_desc->max_z_vrtx_id - < darray_position_size_get(&filter_ctx->scn->vertices)); - ASSERT(positions[origin_desc->max_z_vrtx_id].pos.z - <= positions[hit_desc->max_z_vrtx_id].pos.z); - if(positions[origin_desc->max_z_vrtx_id].pos.z - == positions[hit_desc->max_z_vrtx_id].pos.z) - { - return 1; - } - } return 0; /* Keep */ } @@ -413,7 +560,6 @@ extract_connex_components * to further cancellations */ *used |= SIDE_CANCELED_FLAG(used_side_flag); } - goto canceled; } if(*nbour_used & nbour_side_id) continue; /* Already processed */ @@ -640,6 +786,73 @@ extract_connex_components #endif } +#ifndef NDEBUG +static int +compare_components + (const void* ptr1, const void* ptr2) +{ + const struct cc_descriptor* const* cc1 = ptr1; + const struct cc_descriptor* const* cc2 = ptr2; + ASSERT((*cc1)->side_range.first != (*cc2)->side_range.first); + return (int)(*cc1)->side_range.first - (int)(*cc2)->side_range.first; +} + +static res_T +reorder_components + (struct senc3d_scene* scn, + const struct darray_ptr_component_descriptor* connex_components, + struct darray_triangle_comp* triangles_comp_array) +{ + struct cc_descriptor* const* cc_descriptors; + struct triangle_comp* triangles_comp; + component_id_t c; + component_id_t* new_ids = NULL; + size_t t, cc_count; + res_T res = RES_OK; + + ASSERT(connex_components && triangles_comp_array); + + cc_count = darray_ptr_component_descriptor_size_get(connex_components); + ASSERT(cc_count <= COMPONENT_MAX__); + cc_descriptors = darray_ptr_component_descriptor_cdata_get(connex_components); + + /* Single thread code: can allocate here */ + new_ids = MEM_ALLOC(scn->dev->allocator, sizeof(*new_ids) * cc_count); + if(!new_ids) { + res = RES_MEM_ERR; + goto error; + } + + FOR_EACH(c, 0, cc_count) ASSERT(cc_descriptors[c]->cc_id == c); + /* Sort components by first side */ + qsort(cc_descriptors, cc_count, sizeof(*cc_descriptors), compare_components); + /* Make conversion table */ + FOR_EACH(c, 0, cc_count) { + component_id_t rank = cc_descriptors[c]->cc_id; + new_ids[rank] = c; + } + triangles_comp = darray_triangle_comp_data_get(triangles_comp_array); + FOR_EACH(c, 0, cc_count) { + ASSERT(new_ids[cc_descriptors[c]->cc_id] == c); + /* Update the enclosure ID in cc_descriptor */ + cc_descriptors[c]->cc_id = c; + } + FOR_EACH(t, 0, scn->ntris) { + struct triangle_comp* comp = triangles_comp + t; + component_id_t fc = comp->component[SENC3D_FRONT]; + component_id_t bc = comp->component[SENC3D_BACK]; + ASSERT(new_ids[fc] < cc_count); + comp->component[SENC3D_FRONT] = new_ids[fc]; + ASSERT(new_ids[bc] < cc_count); + comp->component[SENC3D_BACK] = new_ids[bc]; + } + +error: + if(new_ids) MEM_RM(scn->dev->allocator, new_ids); + return res; +} +#endif + static void group_connex_components (struct senc3d_scene* scn, @@ -658,23 +871,34 @@ group_connex_components size_t tmp; component_id_t cc_count; int64_t ccc; - struct filter_ctx filter_ctx; + struct filter_ctx fctx; float lower[3], upper[3]; ASSERT(scn && triangles_comp && connex_components && s3d_view && next_enclosure_id && res); ASSERT(scn->analyze.enclosures_count == 1); +#ifndef NDEBUG +#pragma omp single + { + /* Ensure reproducible cc_ids for debugging purpose */ + *res = reorder_components(scn, connex_components, triangles_comp); + if(*res != RES_OK) return; + } +#endif + descriptors = darray_ptr_component_descriptor_data_get(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(&scn->vertices); - filter_ctx.scn = scn; - filter_ctx.triangles_comp = triangles_comp; - filter_ctx.components = connex_components; + fctx.type = FCTX1; + fctx.c.ctx1.scn = scn; + fctx.c.ctx1.view = s3d_view; + fctx.c.ctx1.triangles_comp = triangles_comp; + fctx.c.ctx1.components = connex_components; *res = s3d_scene_view_get_aabb(s3d_view, lower, upper); - if(*res != RES_OK) return; + if(*res != RES_OK) goto end; /* Cast rays to find links between connex components */ #pragma omp for schedule(dynamic) for(ccc = 0; ccc < (int64_t)cc_count; ccc++) { @@ -704,8 +928,9 @@ group_connex_components f3_set_d3(origin, max_vrtx); /* Self-hit data: self hit if hit this component "on the other side" */ - filter_ctx.origin_component = cc->cc_id; - /* Limit search radius. Only +Z moves are permitted */ + fctx.c.ctx1.origin_component = cc->cc_id; + /* Limit search radius. Only upwards moves (+Z) will be considered. + * Ue a radius allowing to reach the closest top vertex of scene's AABB */ FOR_EACH(i, 0, 2) { ASSERT(lower[i] <= origin[i] && origin[i] <= upper[i]); rrr[i] = MMIN(origin[i] - lower[i], upper[i] - origin[i]); @@ -713,24 +938,24 @@ group_connex_components ASSERT(lower[2] <= origin[2] && origin[2] <= upper[2]); rrr[2] = upper[2] - origin[2]; r = f3_len(rrr) + FLT_EPSILON; /* Ensure r > 0 */ - tmp_res = s3d_scene_view_closest_point(s3d_view, origin, r, &filter_ctx, - &hit); + tmp_res = s3d_scene_view_closest_point(s3d_view, origin, r, &fctx, &hit); if(tmp_res != RES_OK) { *res = tmp_res; continue; } + ASSERT(fctx.c.ctx1.hit_component != COMPONENT_NULL__); /* 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; } else { /* If hit, group this component */ - cc->cc_group_root = filter_ctx.hit_component; + cc->cc_group_root = fctx.c.ctx1.hit_component; ASSERT(cc->cc_group_root < cc_count); } } /* Implicit barrier here */ - if(*res != RES_OK) return; + if(*res != RES_OK) goto end; /* One thread post-processes links to group connex components */ #pragma omp single @@ -782,6 +1007,8 @@ group_connex_components } } /* Implicit barrier here */ +end: + return; } static void @@ -1055,13 +1282,62 @@ compare_enclosures return (int)e1->side_range.first - (int)e2->side_range.first; } +static res_T +reorder_enclosures + (struct senc3d_scene* scn, + const struct darray_ptr_component_descriptor* connex_components) +{ + struct enclosure_data* enclosures; + struct cc_descriptor* const* cc_descriptors; + enclosure_id_t* new_ids; + enclosure_id_t e; + size_t c, cc_count; + res_T res = RES_OK; + + /* Single thread code: can allocate here */ + new_ids = MEM_ALLOC(scn->dev->allocator, + sizeof(*new_ids) * scn->analyze.enclosures_count); + if(!new_ids) { + res = RES_MEM_ERR; + goto error; + } + cc_count = darray_ptr_component_descriptor_size_get(connex_components); + ASSERT(cc_count <= COMPONENT_MAX__); + cc_descriptors = darray_ptr_component_descriptor_cdata_get(connex_components); + enclosures = darray_enclosure_data_get(&scn->analyze.enclosures); + /* Store initial enclosure order */ + FOR_EACH(e, 0, scn->analyze.enclosures_count) + enclosures[e].header.enclosure_id = e; + /* Sort enclosures by first side while keeping enclosure 0 + * at rank 0 (its a convention) */ + qsort(enclosures + 1, scn->analyze.enclosures_count - 1, sizeof(*enclosures), + compare_enclosures); + /* Make conversion table */ + FOR_EACH(e, 0, scn->analyze.enclosures_count) { + enclosure_id_t rank = enclosures[e].header.enclosure_id; + new_ids[rank] = e; + } + FOR_EACH(e, 0, scn->analyze.enclosures_count) { + ASSERT(new_ids[enclosures[e].header.enclosure_id] == e); + enclosures[e].header.enclosure_id = e; + } + FOR_EACH(c, 0, cc_count) { + /* Update the enclosure ID in cc_descriptor */ + enclosure_id_t new_id = new_ids[cc_descriptors[c]->enclosure_id]; + ASSERT(new_id < scn->analyze.enclosures_count); + cc_descriptors[c]->enclosure_id = new_id; + } +error: + if(new_ids) MEM_RM(scn->dev->allocator, new_ids); + return res; +} + static void build_result (struct senc3d_scene* scn, const struct darray_ptr_component_descriptor* connex_components, const struct darray_triangle_comp* triangles_comp_array, struct darray_frontier_edge* frontiers, - enclosure_id_t* ordered_ids, /* Shared error status. * We accept to overwrite an error with a different error */ res_T* res) @@ -1098,33 +1374,28 @@ build_result #pragma omp single { + size_t c; enclosure_id_t e; - size_t d; res_T tmp_res = darray_triangle_enc_resize(&scn->analyze.triangles_enc, scn->ntris); - if(tmp_res != RES_OK) { + if(tmp_res != RES_OK) *res = tmp_res; - goto single_err; + /* Sort enclosures by first side to ensure reproducible results */ + else *res = reorder_enclosures(scn, connex_components); + if(*res != RES_OK) goto single_err; + /* Compute area and volume */ + FOR_EACH(c, 0, cc_count) { + struct enclosure_data* enc; + ASSERT(cc_descriptors[c]->enclosure_id < scn->analyze.enclosures_count); + enc = enclosures + cc_descriptors[c]->enclosure_id; + /* Sum up areas and volumes of components int enclosures */ + enc->header.area += cc_descriptors[c]->_2area; + enc->header.volume += cc_descriptors[c]->_6volume; } - /* Store initial enclosure order */ - FOR_EACH(e, 0, scn->analyze.enclosures_count) - enclosures[e].header.enclosure_id = e; - /* Move enclosures by first side while keeping enclosure 0 - * at rank 0 (its a convention) */ - qsort(enclosures + 1, scn->analyze.enclosures_count - 1, - sizeof(*enclosures), compare_enclosures); - /* Make conversion table */ FOR_EACH(e, 0, scn->analyze.enclosures_count) { - enclosure_id_t rank = enclosures[e].header.enclosure_id; - ordered_ids[rank] = e; - } - FOR_EACH(d, 0, cc_count) { - enclosure_id_t new_id = ordered_ids[cc_descriptors[d]->enclosure_id]; - /* Update the enclosure ID in cc_descriptor */ - cc_descriptors[d]->enclosure_id = new_id; - /* Sum up areas and volumes of components int oenclosures */ - enclosures[new_id].header.area += cc_descriptors[d]->_2area * 0.5; - enclosures[new_id].header.volume += cc_descriptors[d]->_6volume / 6; + struct enclosure_data* enc = enclosures + e; + enc->header.area /= 2; + enc->header.volume /= 6; } single_err: (void)0; }/* Implicit barrier here. */ @@ -1317,7 +1588,6 @@ scene_analyze /* Atomic counters to share beetwen threads */ ATOMIC component_count = 0; ATOMIC next_enclosure_id = 1; - enclosure_id_t* ordered_ids = NULL; res_T res = RES_OK; res_T res2 = RES_OK; @@ -1395,6 +1665,7 @@ scene_analyze if(tmp_res != RES_OK) goto tmp_error2; htable_overlap_iterator_next(&it); } + /* Sort overlapping triangle ids */ qsort(darray_trg_id_data_get(&scn->analyze.overlapping_ids), darray_trg_id_size_get(&scn->analyze.overlapping_ids), sizeof(*darray_trg_id_cdata_get(&scn->analyze.overlapping_ids)), @@ -1474,15 +1745,11 @@ scene_analyze { if(s3d_view) S3D(scene_view_ref_put(s3d_view)); s3d_view = NULL; - ordered_ids = MEM_ALLOC(scn->dev->allocator, - sizeof(*ordered_ids) * scn->analyze.enclosures_count); - if(!ordered_ids) res = RES_MEM_ERR; } /* Implicit barrier here */ if(res != RES_OK) goto error_; /* Step 4: Build result */ - build_result(scn, &connex_components, &triangles_comp, &frontiers, - ordered_ids, &res); + build_result(scn, &connex_components, &triangles_comp, &frontiers, &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 */ @@ -1504,7 +1771,6 @@ scene_analyze { #pragma omp section { - MEM_RM(scn->dev->allocator, ordered_ids); custom_darray_ptr_component_descriptor_release(&connex_components); connex_components_initialized = 0; }