star-line

Structure for accelerating line importance sampling
git clone git://git.meso-star.fr/star-line.git
Log | Files | Refs | README | LICENSE

commit 0ea298a9db29a1fe241ab3df4a341d8c9e40c938
parent 444d1e73ceedb1ba93b60861213c8618e3bdb6e1
Author: Vincent Forest <vincent.forest@meso-star.com>
Date:   Fri, 10 Apr 2026 11:06:05 +0200

Remove the legacy binary tree construction code

This code was temporarily retained to test for regressions and
unexpected behavior during the deployment of the new construction
procedure, in which the tree's arity is now variable. For binary trees,
the two procedures should produce the same results.

Diffstat:
Msrc/sln_tree_build.c | 303+------------------------------------------------------------------------------
1 file changed, 2 insertions(+), 301 deletions(-)

diff --git a/src/sln_tree_build.c b/src/sln_tree_build.c @@ -87,98 +87,6 @@ eval_ka } static res_T -merge_node_polylines - (struct sln_tree* tree, - const struct sln_node* node0, - const struct sln_node* node1, - struct sln_node* merged_node) -{ - struct sln_vertex* vertices = NULL; - size_t vertices_range[2]; - size_t ivtx; - size_t ivtx0, ivtx0_max; - size_t ivtx1, ivtx1_max; - res_T res = RES_OK; - ASSERT(tree && node0 && node1 && merged_node); - - /* Define the vertices range of the merged polyline */ - vertices_range[0] = darray_vertex_size_get(&tree->vertices); - vertices_range[1] = vertices_range[0] + node0->nvertices + node1->nvertices - 1; - - /* Allocate the memory space to store the new polyline */ - res = darray_vertex_resize(&tree->vertices, vertices_range[1] + 1); - if(res != RES_OK) { - ERROR(tree->sln, "Error in merging polylines -- %s.\n", res_to_cstr(res)); - goto error; - } - vertices = darray_vertex_data_get(&tree->vertices); - - ivtx0 = node0->ivertex; - ivtx1 = node1->ivertex; - ivtx0_max = ivtx0 + node0->nvertices - 1; - ivtx1_max = ivtx1 + node1->nvertices - 1; - FOR_EACH(ivtx, vertices_range[0], vertices_range[1]+1) { - const double nu0 = ivtx0 <= ivtx0_max ? vertices[ivtx0].wavenumber : INF; - const double nu1 = ivtx1 <= ivtx1_max ? vertices[ivtx1].wavenumber : INF; - double ka; - float nu; - - if(nu0 < nu1) { - /* The vertex comes from the node0 */ - nu = vertices[ivtx0].wavenumber; - ka = vertices[ivtx0].ka + eval_ka(tree, node1, nu); - ++ivtx0; - } else if(nu0 > nu1) { - /* The vertex comes from the node1 */ - nu = vertices[ivtx1].wavenumber; - ka = vertices[ivtx1].ka + eval_ka(tree, node0, nu); - ++ivtx1; - } else { - /* The vertex is shared by node0 and node1 */ - nu = vertices[ivtx0].wavenumber; - ka = vertices[ivtx0].ka + vertices[ivtx1].ka; - --vertices_range[1]; /* Remove duplicate */ - ++ivtx0; - ++ivtx1; - } - vertices[ivtx].wavenumber = nu; - vertices[ivtx].ka = (float)ka; - } - - /* Decimate the resulting polyline */ - res = polyline_decimate(tree->sln, darray_vertex_data_get(&tree->vertices), - vertices_range, (float)tree->args.mesh_decimation_err, tree->args.mesh_type); - if(res != RES_OK) goto error; - - /* Shrink the size of the vertices */ - darray_vertex_resize(&tree->vertices, vertices_range[1] + 1); - - /* Setup the node polyline */ - merged_node->ivertex = vertices_range[0]; - merged_node->nvertices = ui64_to_ui32(vertices_range[1] - vertices_range[0] + 1); - - /* It is necessary to ensure that the recorded vertices define a polyline - * along which any value (calculated by linear interpolation) is well above - * the sum of the corresponding values of the polylines to be merged. However, - * although this is guaranteed by definition for the vertices of the polyline, - * numerical uncertainty may nevertheless introduce errors that violate this - * criterion. Hence the following adjustment, which slightly increases the ka - * of the mesh so as to guarantee this constraint between the mesh of a node - * and that of its children */ - if(tree->args.mesh_type == SLN_MESH_UPPER) { - FOR_EACH(ivtx, vertices_range[0], vertices_range[1]+1/*inclusive bound*/) { - const double ka = vertices[ivtx].ka; - vertices[ivtx].ka = (float)(ka + MMAX(ka*1e-2, 1e-6)); - } - } - -exit: - return res; -error: - goto exit; -} - -static res_T merge_children_polylines (struct sln_tree* tree, const size_t inode) @@ -333,102 +241,6 @@ static res_T build_polylines(struct sln_tree* tree) { /* Stack */ - #define STACK_SIZE (SLN_TREE_DEPTH_MAX*2/*tree arity*/) - size_t stack[STACK_SIZE]; - size_t istack = 0; - - /* Progress */ - size_t nnodes_total = 0; - size_t nnodes_processed = 0; - int progress = 0; - - /* Miscellaneous */ - size_t inode = 0; - res_T res = RES_OK; - - ASSERT(tree); /* Pre-condition */ - - #define LOG_MSG "Meshing: %3d%%\n" - #define NODE(Id) (darray_node_data_get(&tree->nodes) + (Id)) - #define IS_LEAF(Id) (NODE(Id)->offset == 0) - - nnodes_total = darray_node_size_get(&tree->nodes); - - /* Print that nothing has been done yet */ - INFO(tree->sln, LOG_MSG, progress); - - /* Push back SIZE_MAX which, once pop up, will mark the end of recursion */ - stack[istack++] = SIZE_MAX; - - inode = 0; /* Root node */ - while(inode != SIZE_MAX) { - const size_t istack_saved = istack; - - if(IS_LEAF(inode)) { - res = build_leaf_polyline(tree, NODE(inode)); - if(res != RES_OK) goto error; - - inode = stack[--istack]; /* Pop the next node */ - - } else { - const size_t ichild0 = inode + NODE(inode)->offset + 0; - const size_t ichild1 = inode + NODE(inode)->offset + 1; - - /* Child nodes have their polyline created */ - if(NODE(ichild0)->nvertices) { - ASSERT(NODE(ichild1)->nvertices != 0); - - /* Build the polyline of the current node by merging the polylines of - * its children */ - res = merge_node_polylines - (tree, NODE(ichild0), NODE(ichild1), NODE(inode)); - if(res != RES_OK) goto error; - inode = stack[--istack]; - - } else { - ASSERT(NODE(ichild1)->nvertices == 0); - - ASSERT(istack < STACK_SIZE - 2); - stack[istack++] = inode; /* Push the current node */ - stack[istack++] = ichild1; /* Push child1 */ - - /* Recursively build the polyline of child0 */ - inode = ichild0; - } - } - - /* Handle progression bar */ - if(istack < istack_saved) { - int pcent = 0; - - nnodes_processed += istack_saved - istack; - ASSERT(nnodes_processed <= nnodes_total); - - /* Print progress update */ - pcent = (int)((double)nnodes_processed*100.0/(double)nnodes_total+0.5); - if(pcent/10 > progress/10) { - progress = pcent; - INFO(tree->sln, LOG_MSG, progress); - } - } - } - ASSERT(nnodes_processed == nnodes_total); - - #undef NODE - #undef IS_LEAF - #undef LOG_MSG - #undef STACK_SIZE - -exit: - return res; -error: - goto exit; -} - -static res_T -build_polylines2(struct sln_tree* tree) -{ - /* Stack */ #define STACK_SIZE (SLN_TREE_DEPTH_MAX*SLN_TREE_ARITY_MAX) size_t stack[STACK_SIZE]; size_t istack = 0; @@ -533,117 +345,6 @@ static res_T partition_lines(struct sln_tree* tree) { /* Stack */ - #define STACK_SIZE SLN_TREE_DEPTH_MAX - size_t stack[STACK_SIZE]; - size_t istack = 0; - - /* Progress */ - size_t nlines = 0; - size_t nlines_total = 0; - size_t nlines_processed = 0; - int progress = 0; - - /* Miscellaneous */ - size_t inode = 0; - res_T res = RES_OK; - - ASSERT(tree); /* Pre-condition */ - - SHTR(line_list_get_size(tree->args.lines, &nlines)); - - nlines_total = nlines; - nlines_processed = 0; - - #define LOG_MSG "Partitioning: %3d%%\n" - #define NODE(Id) (darray_node_data_get(&tree->nodes) + (Id)) - #define CREATE_NODE { \ - res = darray_node_push_back(&tree->nodes, &SLN_NODE_NULL); \ - if(res != RES_OK) goto error; \ - } (void)0 - - CREATE_NODE; /* Alloc the root node */ - - /* Setup the root node */ - NODE(0)->range[0] = 0; - NODE(0)->range[1] = nlines - 1; - - /* Push back SIZE_MAX which, once pop up, will mark the end of recursion */ - stack[istack++] = SIZE_MAX; - - /* Print that although nothing has been done yet, - * the calculation is nevertheless in progress */ - INFO(tree->sln, LOG_MSG, progress); - - inode = 0; /* Root node */ - while(inode != SIZE_MAX) { - /* #lines into the node */ - nlines = NODE(inode)->range[1] - NODE(inode)->range[0] + 1; - - /* Make a leaf */ - if(nlines <= MAX_NLINES_PER_LEAF) { - int pcent = 0; - - NODE(inode)->offset = 0; - inode = stack[--istack]; /* Pop the next node */ - - ASSERT(nlines_processed + nlines <= nlines_total); - nlines_processed += nlines; - - /* Print progress update */ - pcent = (int)((double)nlines_processed*100.0/(double)nlines_total+0.5); - if(pcent/10 > progress/10) { - progress = pcent; - INFO(tree->sln, LOG_MSG, progress); - } - - /* Split the node */ - } else { - /* Median split */ - const size_t split = NODE(inode)->range[0] + (nlines + 1/*ceil*/)/2; - - /* Compute the index toward the 2 children (first is the left child, - * followed by the right child) */ - size_t ichildren = darray_node_size_get(&tree->nodes); - - ASSERT(NODE(inode)->range[0] <= split); - ASSERT(NODE(inode)->range[1] >= split); - ASSERT(ichildren > inode); - - /* Define the offset from the current node to its children */ - NODE(inode)->offset = ui64_to_ui32((uint64_t)(ichildren - inode)); - - CREATE_NODE; /* Alloc left child */ - CREATE_NODE; /* Alloc right child */ - - /* Setup the left child */ - NODE(ichildren+0)->range[0] = NODE(inode)->range[0]; - NODE(ichildren+0)->range[1] = split-1; - /* Setup the right child */ - NODE(ichildren+1)->range[0] = split; - NODE(ichildren+1)->range[1] = NODE(inode)->range[1]; - - inode = ichildren + 0; /* Make the left child current node */ - - ASSERT(istack < STACK_SIZE - 1); - stack[istack++] = ichildren + 1; /* Push the right child */ - } - } - ASSERT(nlines_processed == nlines_total); - - #undef NODE - #undef CREATE_NODE - #undef STACK_SIZE -exit: - return res; -error: - darray_node_clear(&tree->nodes); - goto exit; -} - -static res_T -partition_lines2(struct sln_tree* tree) -{ - /* Stack */ #define STACK_SIZE (SLN_TREE_DEPTH_MAX*(SLN_TREE_ARITY_MAX-1)) size_t stack[STACK_SIZE]; size_t istack = 0; @@ -782,8 +483,8 @@ tree_build(struct sln_tree* tree) { res_T res = RES_OK; - if((res = partition_lines2(tree)) != RES_OK) goto error; - if((res = build_polylines2(tree)) != RES_OK) goto error; + if((res = partition_lines(tree)) != RES_OK) goto error; + if((res = build_polylines(tree)) != RES_OK) goto error; exit: return res;