star-line

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

sln_tree_build.c (22419B)


      1 /* Copyright (C) 2022, 2026 |Méso|Star> (contact@meso-star.com)
      2  * Copyright (C) 2026 Université de Lorraine
      3  * Copyright (C) 2022 Centre National de la Recherche Scientifique
      4  * Copyright (C) 2022 Université Paul Sabatier
      5  *
      6  * This program is free software: you can redistribute it and/or modify
      7  * it under the terms of the GNU General Public License as published by
      8  * the Free Software Foundation, either version 3 of the License, or
      9  * (at your option) any later version.
     10  *
     11  * This program is distributed in the hope that it will be useful,
     12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
     13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
     14  * GNU General Public License for more details.
     15  *
     16  * You should have received a copy of the GNU General Public License
     17  * along with this program. If not, see <http://www.gnu.org/licenses/>. */
     18 
     19 #include "sln.h"
     20 #include "sln_device_c.h"
     21 #include "sln_polyline.h"
     22 #include "sln_tree_c.h"
     23 
     24 #include <star/shtr.h>
     25 
     26 #include <rsys/cstr.h>
     27 
     28 /* Maximum number of lines per leaf */
     29 #define MAX_NLINES_PER_LEAF 1
     30 
     31 /* STACK_SIZE is set to the maximum depth of the partitioning tree generated by
     32  * the tree_build function. This is actually the maximum stack size where tree
     33  * nodes are pushed by the recursive build process
     34  *
     35  * FIXME the previous comment is not true with a tree whose arity is not 2 */
     36 #define STACK_SIZE 64
     37 
     38 /*******************************************************************************
     39  * Helper functions
     40  ******************************************************************************/
     41 static FINLINE uint32_t
     42 ui64_to_ui32(const uint64_t ui64)
     43 {
     44   if(ui64 > UINT32_MAX)
     45     VFATAL("%s: overflow %lu.\n", ARG2(FUNC_NAME, ((unsigned long)ui64)));
     46   return (uint32_t)ui64;
     47 }
     48 
     49 static INLINE res_T
     50 build_leaf_polyline(struct sln_tree* tree, struct sln_node* leaf)
     51 {
     52   size_t vertices_range[2];
     53   res_T res = RES_OK;
     54 
     55   /* Currently assume that we have only one line per leaf */
     56   ASSERT(leaf->range[0] == leaf->range[1]);
     57 
     58   /* Line meshing */
     59   res = line_mesh(tree, leaf->range[0], tree->args.nvertices_hint, vertices_range);
     60   if(res != RES_OK) goto error;
     61 
     62   /* Decimate the line mesh */
     63   res = polyline_decimate(tree->sln, darray_vertex_data_get(&tree->vertices),
     64      vertices_range, (float)tree->args.mesh_decimation_err, tree->args.mesh_type);
     65   if(res != RES_OK) goto error;
     66 
     67   /* Shrink the size of the vertices */
     68   darray_vertex_resize(&tree->vertices, vertices_range[1] + 1);
     69 
     70   /* Setup the leaf polyline  */
     71   leaf->ivertex = vertices_range[0];
     72   leaf->nvertices = ui64_to_ui32(vertices_range[1] - vertices_range[0] + 1);
     73 
     74 exit:
     75   return res;
     76 error:
     77   goto exit;
     78 }
     79 
     80 static INLINE double
     81 eval_ka
     82   (const struct sln_tree* tree,
     83    const struct sln_node* node,
     84    const double wavenumber)
     85 {
     86   struct sln_mesh mesh = SLN_MESH_NULL;
     87   double ka = 0;
     88   ASSERT(tree && node);
     89 
     90   /* Whether the mesh to be constructed corresponds to the spectrum or its upper
     91    * limit, use the node mesh to calculate the value of ka at a given wave
     92    * number. Calculating the value from the node lines would take far too long*/
     93   SLN(node_get_mesh(tree, node, &mesh));
     94   ka = sln_mesh_eval(&mesh, wavenumber);
     95 
     96   return ka;
     97 }
     98 
     99 static res_T
    100 merge_node_polylines
    101   (struct sln_tree* tree,
    102    const struct sln_node* node0,
    103    const struct sln_node* node1,
    104    struct sln_node* merged_node)
    105 {
    106   struct sln_vertex* vertices = NULL;
    107   size_t vertices_range[2];
    108   size_t ivtx;
    109   size_t ivtx0, ivtx0_max;
    110   size_t ivtx1, ivtx1_max;
    111   res_T res = RES_OK;
    112   ASSERT(tree && node0 && node1 && merged_node);
    113 
    114   /* Define the vertices range of the merged polyline */
    115   vertices_range[0] = darray_vertex_size_get(&tree->vertices);
    116   vertices_range[1] = vertices_range[0] + node0->nvertices + node1->nvertices - 1;
    117 
    118   /* Allocate the memory space to store the new polyline */
    119   res = darray_vertex_resize(&tree->vertices, vertices_range[1] + 1);
    120   if(res != RES_OK) {
    121     ERROR(tree->sln, "Error in merging polylines -- %s.\n", res_to_cstr(res));
    122     goto error;
    123   }
    124   vertices = darray_vertex_data_get(&tree->vertices);
    125 
    126   ivtx0 = node0->ivertex;
    127   ivtx1 = node1->ivertex;
    128   ivtx0_max = ivtx0 + node0->nvertices - 1;
    129   ivtx1_max = ivtx1 + node1->nvertices - 1;
    130   FOR_EACH(ivtx, vertices_range[0], vertices_range[1]+1) {
    131     const double nu0 = ivtx0 <= ivtx0_max ? vertices[ivtx0].wavenumber : INF;
    132     const double nu1 = ivtx1 <= ivtx1_max ? vertices[ivtx1].wavenumber : INF;
    133     double ka;
    134     float nu;
    135 
    136     if(nu0 < nu1) {
    137       /* The vertex comes from the node0 */
    138       nu = vertices[ivtx0].wavenumber;
    139       ka = vertices[ivtx0].ka + eval_ka(tree, node1, nu);
    140       ++ivtx0;
    141     } else if(nu0 > nu1) {
    142       /* The vertex comes from the node1 */
    143       nu = vertices[ivtx1].wavenumber;
    144       ka = vertices[ivtx1].ka + eval_ka(tree, node0, nu);
    145       ++ivtx1;
    146     } else {
    147       /* The vertex is shared by node0 and node1 */
    148       nu = vertices[ivtx0].wavenumber;
    149       ka = vertices[ivtx0].ka + vertices[ivtx1].ka;
    150       --vertices_range[1]; /* Remove duplicate */
    151       ++ivtx0;
    152       ++ivtx1;
    153     }
    154     vertices[ivtx].wavenumber = nu;
    155     vertices[ivtx].ka = (float)ka;
    156   }
    157 
    158   /* Decimate the resulting polyline */
    159   res = polyline_decimate(tree->sln, darray_vertex_data_get(&tree->vertices),
    160      vertices_range, (float)tree->args.mesh_decimation_err, tree->args.mesh_type);
    161   if(res != RES_OK) goto error;
    162 
    163   /* Shrink the size of the vertices */
    164   darray_vertex_resize(&tree->vertices, vertices_range[1] + 1);
    165 
    166   /* Setup the node polyline */
    167   merged_node->ivertex = vertices_range[0];
    168   merged_node->nvertices = ui64_to_ui32(vertices_range[1] - vertices_range[0] + 1);
    169 
    170   /* It is necessary to ensure that the recorded vertices define a polyline
    171    * along which any value (calculated by linear interpolation) is well above
    172    * the sum of the corresponding values of the polylines to be merged. However,
    173    * although this is guaranteed by definition for the vertices of the polyline,
    174    * numerical uncertainty may nevertheless introduce errors that violate this
    175    * criterion. Hence the following adjustment, which slightly increases the ka
    176    * of the mesh so as to guarantee this constraint between the mesh of a node
    177    * and that of its children */
    178   if(tree->args.mesh_type == SLN_MESH_UPPER) {
    179     FOR_EACH(ivtx, vertices_range[0], vertices_range[1]+1/*inclusive bound*/) {
    180       const double ka = vertices[ivtx].ka;
    181       vertices[ivtx].ka = (float)(ka + MMAX(ka*1e-2, 1e-6));
    182     }
    183   }
    184 
    185 exit:
    186   return res;
    187 error:
    188   goto exit;
    189 }
    190 
    191 static res_T
    192 merge_children_polylines
    193   (struct sln_tree* tree,
    194    const size_t inode)
    195 {
    196   /* Helper constant */
    197   static const size_t NO_MORE_VERTEX = SIZE_MAX;
    198 
    199   /* Polyline vertices */
    200   size_t children_ivtx[SLN_TREE_ARITY_MAX] = {0};
    201   size_t vertices_range[2] = {0,0};
    202   struct sln_vertex* vertices = NULL;
    203   size_t ivtx = 0;
    204   size_t nvertices = 0;
    205 
    206   /* Miscellaneous */
    207   unsigned i = 0;
    208   res_T res = RES_OK;
    209 
    210   /* Pre-conditions */
    211   ASSERT(tree && inode < darray_node_size_get(&tree->nodes));
    212   ASSERT(tree->args.arity >= 2);
    213   ASSERT(tree->args.arity <= SLN_TREE_ARITY_MAX);
    214 
    215   #define NODE(Id) (darray_node_data_get(&tree->nodes) + (Id))
    216 
    217   /* Compute the number of vertices to be merged,
    218    * i.e., the sum of vertices of the children */
    219   nvertices = 0;
    220   FOR_EACH(i, 0, tree->args.arity) {
    221     const size_t ichild = inode + NODE(inode)->offset + i;
    222     nvertices += NODE(ichild)->nvertices;
    223   }
    224 
    225   /* Define the vertices range of the merged polyline */
    226   vertices_range[0] = darray_vertex_size_get(&tree->vertices);
    227   vertices_range[1] = vertices_range[0] + nvertices - 1/*inclusive bound*/;
    228 
    229   /* Allocate the memory space to store the new polyline */
    230   res = darray_vertex_resize(&tree->vertices, vertices_range[1]+1);
    231   if(res != RES_OK) {
    232     ERROR(tree->sln, "Error in merging polylines -- %s.\n", res_to_cstr(res));
    233     goto error;
    234   }
    235   vertices = darray_vertex_data_get(&tree->vertices);
    236 
    237   /* Initialize the vertex index list. For each child, the initial value
    238    * corresponds to the index of its first vertex. This index will be
    239    * incremented as vertices are merged into the parent polyline. */
    240   FOR_EACH(i, 0, tree->args.arity) {
    241     const size_t ichild = inode + NODE(inode)->offset + i;
    242     children_ivtx[i] = NODE(ichild)->ivertex;
    243   }
    244 
    245   FOR_EACH(ivtx, vertices_range[0], vertices_range[1]+1/*inclusive bound*/) {
    246     double ka = 0;
    247     double nu = INF;
    248 
    249     /* The number of vertices corresponding to the current wave number for which
    250      * the parent ka is calculated. It is at least equal to one, since this nu
    251      * is defined by the child vertices, but may be greater if multiple children
    252      * share the same vertex, i.e., a ka value calculated for the same nu */
    253     unsigned nvertices_merged = 0;
    254 
    255     /* Find the minimum wave number among the vertices of the child vertices
    256      * that are candidates for merging */
    257     FOR_EACH(i, 0, tree->args.arity) {
    258       const size_t child_ivtx = children_ivtx[i];
    259       if(child_ivtx != NO_MORE_VERTEX) {
    260         nu = MMIN(nu, vertices[child_ivtx].wavenumber);
    261       }
    262     }
    263     ASSERT(nu != INF); /* At least one vertex must have been found */
    264 
    265     /* Compute the value of ka at the wave number determined above */
    266     FOR_EACH(i, 0, tree->args.arity) {
    267       const size_t child_ivtx = children_ivtx[i];
    268       const struct sln_node* child = NODE(inode + NODE(inode)->offset + i);
    269 
    270       if(child_ivtx == NO_MORE_VERTEX
    271       || nu != vertices[child_ivtx].wavenumber) {
    272         /* The wave number does not correspond to a vertex in the current
    273          * child's mesh. Therefore, its contribution to the parent node's ka is
    274          * computed  */
    275         ka += eval_ka(tree, child, nu);
    276 
    277       } else {
    278         /* The wave number is the one for which the child node stores a ka
    279          * value. Add it to the parent node's ka value and designate the child's
    280          * next vertex as a candidate for merging into the parent. The exception
    281          * is when all vertices of the child have already been merged. In this
    282          * case, report that the child no longer has any candidate vertices */
    283         ka += vertices[child_ivtx].ka;
    284         ++children_ivtx[i];
    285         if(children_ivtx[i] >= child->ivertex + child->nvertices) {
    286           children_ivtx[i] = NO_MORE_VERTEX;
    287         }
    288         ++nvertices_merged; /* Record that a vertex has been merged */
    289       }
    290     }
    291 
    292     /* Setup the parent vertex */
    293     vertices[ivtx].wavenumber = (float)nu;
    294     vertices[ivtx].ka = (float)ka;
    295 
    296     /* If multiple child vertices have been merged, then a single wave number
    297      * corresponds to a vertex with multiple children. The number of parent
    298      * vertices is therefore no longer the sum of the number of its children's
    299      * vertices, since some vertices are duplicated. Hence the following
    300      * adjustment, which removes the duplicate vertices. */
    301     vertices_range[1] -= (nvertices_merged-1);
    302   }
    303 
    304   /* Decimate the resulting polyline */
    305   res = polyline_decimate(tree->sln, darray_vertex_data_get(&tree->vertices),
    306      vertices_range, (float)tree->args.mesh_decimation_err, tree->args.mesh_type);
    307   if(res != RES_OK) goto error;
    308 
    309   /* Shrink the size of the vertices */
    310   darray_vertex_resize(&tree->vertices, vertices_range[1] + 1);
    311 
    312   /* Setup the node polyline */
    313   NODE(inode)->ivertex = vertices_range[0];
    314   NODE(inode)->nvertices = ui64_to_ui32(vertices_range[1] - vertices_range[0] + 1);
    315 
    316   /* It is necessary to ensure that the recorded vertices define a polyline
    317    * along which any value (calculated by linear interpolation) is well above
    318    * the sum of the corresponding values of the polylines to be merged. However,
    319    * although this is guaranteed by definition for the vertices of the polyline,
    320    * numerical uncertainty may nevertheless introduce errors that violate this
    321    * criterion. Hence the following adjustment, which slightly increases the ka
    322    * of the mesh so as to guarantee this constraint between the mesh of a node
    323    * and that of its children */
    324   if(tree->args.mesh_type == SLN_MESH_UPPER) {
    325     FOR_EACH(ivtx, vertices_range[0], vertices_range[1]+1/*inclusive bound*/) {
    326       const double ka = vertices[ivtx].ka;
    327       vertices[ivtx].ka = (float)(ka + MMAX(ka*1e-2, 1e-6));
    328     }
    329   }
    330 
    331   #undef NODE
    332 
    333 exit:
    334   return res;
    335 error:
    336   goto exit;
    337 }
    338 
    339 static res_T
    340 build_polylines(struct sln_tree* tree)
    341 {
    342   size_t stack[STACK_SIZE*2];
    343   size_t istack = 0;
    344   size_t inode = 0;
    345   size_t nnodes_total = 0;
    346   size_t nnodes_processed = 0;
    347   int progress = 0;
    348   res_T res = RES_OK;
    349   ASSERT(tree);
    350 
    351   #define LOG_MSG "Meshing: %3d%%\n"
    352   #define NODE(Id) (darray_node_data_get(&tree->nodes) + (Id))
    353   #define IS_LEAF(Id) (NODE(Id)->offset == 0)
    354 
    355   nnodes_total = darray_node_size_get(&tree->nodes);
    356 
    357   /* Print that nothing has been done yet */
    358   INFO(tree->sln, LOG_MSG, progress);
    359 
    360   /* Push back SIZE_MAX which, once pop up, will mark the end of recursion */
    361   stack[istack++] = SIZE_MAX;
    362 
    363   inode = 0; /* Root node */
    364   while(inode != SIZE_MAX) {
    365     const size_t istack_saved = istack;
    366 
    367     if(IS_LEAF(inode)) {
    368       res = build_leaf_polyline(tree, NODE(inode));
    369       if(res != RES_OK) goto error;
    370 
    371       inode = stack[--istack]; /* Pop the next node */
    372 
    373     } else {
    374       const size_t ichild0 = inode + NODE(inode)->offset + 0;
    375       const size_t ichild1 = inode + NODE(inode)->offset + 1;
    376 
    377       /* Child nodes have their polyline created */
    378       if(NODE(ichild0)->nvertices) {
    379         ASSERT(NODE(ichild1)->nvertices != 0);
    380 
    381         /* Build the polyline of the current node by merging the polylines of
    382          * its children */
    383 #if 0
    384         res = merge_node_polylines
    385           (tree, NODE(ichild0), NODE(ichild1), NODE(inode));
    386 #else
    387         res = merge_children_polylines(tree, inode);
    388 #endif
    389         if(res != RES_OK) goto error;
    390         inode = stack[--istack];
    391 
    392       } else {
    393         ASSERT(NODE(ichild1)->nvertices == 0);
    394 
    395         ASSERT(istack < STACK_SIZE - 2);
    396         stack[istack++] = inode; /* Push the current node */
    397         stack[istack++] = ichild1; /* Push child1 */
    398 
    399         /* Recursively build the polyline of child0 */
    400         inode = ichild0;
    401       }
    402     }
    403 
    404     /* Handle progression bar */
    405     if(istack < istack_saved) {
    406       int pcent = 0;
    407 
    408       nnodes_processed += istack_saved - istack;
    409       ASSERT(nnodes_processed <= nnodes_total);
    410 
    411       /* Print progress update */
    412       pcent = (int)((double)nnodes_processed*100.0/(double)nnodes_total+0.5);
    413       if(pcent/10 > progress/10) {
    414         progress = pcent;
    415         INFO(tree->sln, LOG_MSG, progress);
    416       }
    417     }
    418   }
    419   ASSERT(nnodes_processed == nnodes_total);
    420 
    421   #undef NODE
    422   #undef IS_LEAF
    423   #undef LOG_MSG
    424 
    425 exit:
    426   return res;
    427 error:
    428   goto exit;
    429 }
    430 
    431 static res_T
    432 partition_lines(struct sln_tree* tree)
    433 {
    434   size_t stack[STACK_SIZE];
    435   size_t istack = 0;
    436   size_t inode = 0;
    437   size_t nlines = 0;
    438   size_t nlines_total = 0;
    439   size_t nlines_processed = 0;
    440   int progress = 0;
    441   res_T res = RES_OK;
    442 
    443   ASSERT(tree);
    444   SHTR(line_list_get_size(tree->args.lines, &nlines));
    445 
    446   nlines_total = nlines;
    447   nlines_processed = 0;
    448 
    449   #define LOG_MSG "Partitioning: %3d%%\n"
    450   #define NODE(Id) (darray_node_data_get(&tree->nodes) + (Id))
    451   #define CREATE_NODE {                                                        \
    452     res = darray_node_push_back(&tree->nodes, &SLN_NODE_NULL);                 \
    453     if(res != RES_OK) goto error;                                              \
    454   } (void)0
    455 
    456   CREATE_NODE; /* Alloc the root node */
    457 
    458   /* Setup the root node */
    459   NODE(0)->range[0] = 0;
    460   NODE(0)->range[1] = nlines - 1;
    461 
    462   /* Push back SIZE_MAX which, once pop up, will mark the end of recursion */
    463   stack[istack++] = SIZE_MAX;
    464 
    465   /* Print that although nothing has been done yet,
    466    * the calculation is nevertheless in progress */
    467   INFO(tree->sln, LOG_MSG, progress);
    468 
    469   inode = 0; /* Root node */
    470   while(inode != SIZE_MAX) {
    471     /* #lines into the node */
    472     nlines = NODE(inode)->range[1] - NODE(inode)->range[0] + 1;
    473 
    474     /* Make a leaf */
    475     if(nlines <= MAX_NLINES_PER_LEAF) {
    476       int pcent = 0;
    477 
    478       NODE(inode)->offset = 0;
    479       inode = stack[--istack]; /* Pop the next node */
    480 
    481       ASSERT(nlines_processed + nlines <= nlines_total);
    482       nlines_processed += nlines;
    483 
    484       /* Print progress update */
    485       pcent = (int)((double)nlines_processed*100.0/(double)nlines_total+0.5);
    486       if(pcent/10 > progress/10) {
    487         progress = pcent;
    488         INFO(tree->sln, LOG_MSG, progress);
    489       }
    490 
    491     /* Split the node  */
    492     } else {
    493       /* Median split */
    494       const size_t split = NODE(inode)->range[0] + (nlines + 1/*ceil*/)/2;
    495 
    496       /* Compute the index toward the 2 children (first is the left child,
    497        * followed by the right child) */
    498       size_t ichildren = darray_node_size_get(&tree->nodes);
    499 
    500       ASSERT(NODE(inode)->range[0] <= split);
    501       ASSERT(NODE(inode)->range[1] >= split);
    502       ASSERT(ichildren > inode);
    503 
    504       /* Define the offset from the current node to its children */
    505       NODE(inode)->offset = ui64_to_ui32((uint64_t)(ichildren - inode));
    506 
    507       CREATE_NODE; /* Alloc left child */
    508       CREATE_NODE; /* Alloc right child */
    509 
    510       /* Setup the left child  */
    511       NODE(ichildren+0)->range[0] = NODE(inode)->range[0];
    512       NODE(ichildren+0)->range[1] = split-1;
    513       /* Setup the right child */
    514       NODE(ichildren+1)->range[0] = split;
    515       NODE(ichildren+1)->range[1] = NODE(inode)->range[1];
    516 
    517       inode = ichildren + 0; /* Make the left child current node */
    518 
    519       ASSERT(istack < STACK_SIZE - 1);
    520       stack[istack++] = ichildren + 1; /* Push the right child */
    521     }
    522   }
    523   ASSERT(nlines_processed == nlines_total);
    524 
    525   #undef NODE
    526   #undef CREATE_NODE
    527 
    528 exit:
    529   return res;
    530 error:
    531   darray_node_clear(&tree->nodes);
    532   goto exit;
    533 }
    534 
    535 static res_T
    536 partition_lines2(struct sln_tree* tree)
    537 {
    538   size_t stack[STACK_SIZE];
    539   size_t istack = 0;
    540   size_t inode = 0;
    541   size_t nlines = 0;
    542   size_t nlines_total = 0;
    543   size_t nlines_processed = 0;
    544   int progress = 0;
    545   unsigned arity = 0;
    546   res_T res = RES_OK;
    547 
    548   ASSERT(tree);
    549   SHTR(line_list_get_size(tree->args.lines, &nlines));
    550 
    551   arity = tree->args.arity;
    552   nlines_total = nlines;
    553   nlines_processed = 0;
    554 
    555   #define LOG_MSG "Partitioning: %3d%%\n"
    556   #define NODE(Id) (darray_node_data_get(&tree->nodes) + (Id))
    557   #define CREATE_NODE {                                                        \
    558     res = darray_node_push_back(&tree->nodes, &SLN_NODE_NULL);                 \
    559     if(res != RES_OK) goto error;                                              \
    560   } (void)0
    561 
    562   CREATE_NODE; /* Alloc the root node */
    563 
    564   /* Setup the root node */
    565   NODE(0)->range[0] = 0;
    566   NODE(0)->range[1] = nlines - 1;
    567 
    568   /* Push back SIZE_MAX which, once pop up, will mark the end of recursion */
    569   stack[istack++] = SIZE_MAX;
    570 
    571   /* Print that although nothing has been done yet,
    572    * the calculation is nevertheless in progress */
    573   INFO(tree->sln, LOG_MSG, progress);
    574 
    575   inode = 0; /* Root node */
    576   while(inode != SIZE_MAX) {
    577     /* #lines into the node */
    578     nlines = NODE(inode)->range[1] - NODE(inode)->range[0] + 1;
    579 
    580     /* Make a leaf */
    581     if(nlines <= MAX_NLINES_PER_LEAF) {
    582       int pcent = 0;
    583 
    584       NODE(inode)->offset = 0;
    585       inode = stack[--istack]; /* Pop the next node */
    586 
    587       ASSERT(nlines_processed + nlines <= nlines_total);
    588       nlines_processed += nlines;
    589 
    590       /* Print progress update */
    591       pcent = (int)((double)nlines_processed*100.0/(double)nlines_total+0.5);
    592       if(pcent/10 > progress/10) {
    593         progress = pcent;
    594         INFO(tree->sln, LOG_MSG, progress);
    595       }
    596 
    597     /* Split the node  */
    598     } else {
    599       size_t node_range[2] = {0,0};
    600       size_t i = 0;
    601 
    602       /* Compute how the lines of a node are distributed among its children,
    603        * as well as the number of children that node has */
    604       const size_t nlines_child = (nlines + arity/2/*ceil*/)/arity;
    605       const size_t nchildren = (nlines + nlines_child/2/*ceil*/)/nlines_child;
    606 
    607       /* Calculate the index of the first child */
    608       size_t ichildren = darray_node_size_get(&tree->nodes);
    609 
    610       ASSERT(nchildren <= arity);
    611       ASSERT(ichildren > inode);
    612 
    613       node_range[0] = NODE(inode)->range[0];
    614       node_range[1] = NODE(inode)->range[1];
    615 
    616       /* Define the offset from the current node to its children */
    617       NODE(inode)->offset = ui64_to_ui32((uint64_t)(ichildren - inode));
    618 
    619       FOR_EACH(i, 0, nchildren) {
    620         CREATE_NODE;
    621 
    622         /* Set the range of lines line for the newly created child. Note that
    623          * the boundaries of the range are inclusive, which is why 1 is
    624          * subtracted to the upper bound */
    625         NODE(ichildren+i)->range[0] = node_range[0] + nlines_child*(i+0);
    626         NODE(ichildren+i)->range[1] = node_range[0] + nlines_child*(i+1)-1;
    627 
    628         /* Check that the child's lines are a subset of the parent's lines. Note
    629          * that the upper bound of the child's line interval is not checked, as
    630          * it may temporarily exceed the parent's interval since, up to this
    631          * point, it is assumed that all children contain the same number of
    632          * lines. This is no longer true if the number of lines in the parent is
    633          * not a multiple of the computed number of lines per child. The upper
    634          * bound of the last child's range is therefore readjusted at the end of
    635          * the loop (see below) */
    636         ASSERT(NODE(ichildren+i)->range[0] <= node_range[1]);
    637       }
    638       /* Readjust the upper bound of the last child to ensure that its line
    639        * interval is fully contained within its parent node. See above for an
    640        * explanation of why */
    641       NODE(ichildren + nchildren - 1)->range[1] = node_range[1];
    642 
    643       inode = ichildren; /* Make the first child the current node */
    644       FOR_EACH(i, 1, nchildren) {
    645         stack[istack++] = ichildren + i; /* Push the other children */
    646       }
    647     }
    648   }
    649   ASSERT(nlines_processed == nlines_total);
    650 
    651   #undef NODE
    652   #undef CREATE_NODE
    653 exit:
    654   return res;
    655 error:
    656   darray_node_clear(&tree->nodes);
    657   goto exit;
    658 }
    659 
    660 /*******************************************************************************
    661  * Local functions
    662  ******************************************************************************/
    663 res_T
    664 tree_build(struct sln_tree* tree)
    665 {
    666   res_T res = RES_OK;
    667 
    668   if((res = partition_lines2(tree)) != RES_OK) goto error;
    669   if((res = build_polylines(tree)) != RES_OK) goto error;
    670 
    671 exit:
    672   return res;
    673 error:
    674   goto exit;
    675 }