commit d5d860d7610ab1447e6cc0c769a186a4fd7a5949
parent 40830282abcb264403f2725f64250b479a53fc03
Author: Vincent Forest <vincent.forest@meso-star.com>
Date: Wed, 26 Sep 2018 12:38:22 +0200
Provide "voxels" to the "challenge_merge" functor
Diffstat:
8 files changed, 450 insertions(+), 63 deletions(-)
diff --git a/src/svx.h b/src/svx.h
@@ -88,10 +88,11 @@ struct svx_voxel_desc {
const size_t nvoxels, /* #submitted data */
void* ctx); /* Pointer toward user data */
- /* Check if the voxel's data can be merged */
+ /* Check if the voxel's data can be merged. Note that the `id' field of the
+ * submitted voxels is undefined since these voxels are temporaries. */
int
(*challenge_merge)
- (const void* voxels[], /* Data candidates to the merge */
+ (const struct svx_voxel voxels[], /* Voxels candidate to the merge */
const size_t nvoxels, /* #candidates */
void* ctx); /* Pointer toward user data */
diff --git a/src/svx_bintree.c b/src/svx_bintree.c
@@ -71,8 +71,21 @@ svx_bintree_create
/* The binary tree definition that must be a power of two */
bintree->definition = round_up_pow2(nvoxels);
+ bintree->frame[0] = axis;
+
+ /* Setup the binary tree AABB in world space */
+ bintree->lower[axis] = lower;
+ bintree->upper[axis] = upper;
+
+ /* Compute the world space range of the binary tree */
+ vox_sz = (upper - lower) / (double)nvoxels;
+ bintree->tree_low[axis] = lower;
+ bintree->tree_upp[axis] = lower + (double)bintree->definition * vox_sz;
+ bintree->tree_size[axis] = bintree->tree_upp[axis] - bintree->tree_low[axis];
+
/* Initialize the bintree builder */
- res = bintree_builder_init(&bldr, bintree->definition, desc, &bintree->buffer);
+ res = bintree_builder_init(&bldr, bintree->definition, bintree->tree_low,
+ bintree->tree_upp, bintree->frame, desc, &bintree->buffer);
if(res != RES_OK) goto error;
/* Allocate the memory space of a temporary voxel */
@@ -105,20 +118,14 @@ svx_bintree_create
ASSERT(bldr.tree_depth > bldr.non_empty_lvl);
#ifndef NDEBUG
- bintree_check(&bintree->buffer, bintree->root);
+ {
+ size_t nleaves = 0;
+ bintree_check(&bintree->buffer, bintree->root, &nleaves);
+ CHK(nleaves == bintree->nleaves);
+ }
#endif
- bintree->lower[axis] = lower;
- bintree->upper[axis] = upper;
- bintree->frame[0] = axis;
-
- /* Compute the world space range of the binary tree */
- vox_sz = (upper - lower) / (double)nvoxels;
- bintree->tree_low[axis] = lower;
- bintree->tree_upp[axis] = lower + (double)bintree->definition * vox_sz;
- bintree->tree_size[axis] = bintree->tree_upp[axis] - bintree->tree_low[axis];
-
- /* Adjuste the binary tree definition with respect to the binary tree depth
+ /* Adjust the binary tree definition with respect to the binary tree depth
* since finest levels might be removed during the build due to the merging
* process */
bintree->definition = bintree->definition / (1ull<<bldr.non_empty_lvl);
diff --git a/src/svx_octree.c b/src/svx_octree.c
@@ -77,14 +77,39 @@ svx_octree_create
res = tree_create(dev, SVX_OCTREE, desc->size, &oct);
if(res != RES_OK) goto error;
+ oct->frame[0] = SVX_AXIS_X;
+ oct->frame[1] = SVX_AXIS_Y;
+ oct->frame[2] = SVX_AXIS_Z;
+
/* Compute the octree definition */
oct->definition = MMAX(nvoxels[0], 2);
oct->definition = MMAX(nvoxels[1], oct->definition);
oct->definition = MMAX(nvoxels[2], oct->definition);
oct->definition = round_up_pow2(oct->definition);
+ /* Setup the octree AABB in world space */
+ d3_set(oct->lower, lower);
+ d3_set(oct->upper, upper);
+
+ /* Compute the voxel size in world space */
+ vox_sz[0] = (upper[0] - lower[0]) / (double)nvoxels[0];
+ vox_sz[1] = (upper[1] - lower[1]) / (double)nvoxels[1];
+ vox_sz[2] = (upper[2] - lower[2]) / (double)nvoxels[2];
+
+ /* Compute the octree AABB in world space */
+ oct->tree_low[0] = lower[0];
+ oct->tree_low[1] = lower[1];
+ oct->tree_low[2] = lower[2];
+ oct->tree_upp[0] = oct->tree_low[0] + (double)oct->definition * vox_sz[0];
+ oct->tree_upp[1] = oct->tree_low[1] + (double)oct->definition * vox_sz[1];
+ oct->tree_upp[2] = oct->tree_low[2] + (double)oct->definition * vox_sz[2];
+ oct->tree_size[0] = oct->tree_upp[0] - oct->tree_low[0];
+ oct->tree_size[1] = oct->tree_upp[1] - oct->tree_low[1];
+ oct->tree_size[2] = oct->tree_upp[2] - oct->tree_low[2];
+
/* Intialize the octree builder */
- res = octree_builder_init(&bldr, oct->definition, desc, &oct->buffer);
+ res = octree_builder_init(&bldr, oct->definition, oct->tree_low,
+ oct->tree_upp, oct->frame, desc, &oct->buffer);
if(res != RES_OK) goto error;
vox.data = MEM_CALLOC(dev->allocator, 1, desc->size);
@@ -127,31 +152,12 @@ svx_octree_create
ASSERT(bldr.tree_depth > bldr.non_empty_lvl);
#ifndef NDEBUG
- octree_check(&oct->buffer, oct->root);
+ {
+ size_t nleaves = 0;
+ octree_check(&oct->buffer, oct->root, &nleaves);
+ CHK(nleaves == oct->nleaves);
+ }
#endif
- oct->frame[0] = SVX_AXIS_X;
- oct->frame[1] = SVX_AXIS_Y;
- oct->frame[2] = SVX_AXIS_Z;
-
- /* Setup the octree AABB in world space */
- d3_set(oct->lower, lower);
- d3_set(oct->upper, upper);
-
- /* Compute the voxel size in world space */
- vox_sz[0] = (upper[0] - lower[0]) / (double)nvoxels[0];
- vox_sz[1] = (upper[1] - lower[1]) / (double)nvoxels[1];
- vox_sz[2] = (upper[2] - lower[2]) / (double)nvoxels[2];
-
- /* Compute the octree AABB in world space */
- oct->tree_low[0] = lower[0];
- oct->tree_low[1] = lower[1];
- oct->tree_low[2] = lower[2];
- oct->tree_upp[0] = oct->tree_low[0] + (double)oct->definition * vox_sz[0];
- oct->tree_upp[1] = oct->tree_low[1] + (double)oct->definition * vox_sz[1];
- oct->tree_upp[2] = oct->tree_low[2] + (double)oct->definition * vox_sz[2];
- oct->tree_size[0] = oct->tree_upp[0] - oct->tree_low[0];
- oct->tree_size[1] = oct->tree_upp[1] - oct->tree_low[1];
- oct->tree_size[2] = oct->tree_upp[2] - oct->tree_low[2];
/* Adjust the octree definition with respect to the octree depth since finest
* levels might be removed during the build due to the merging process */
diff --git a/src/svx_tree_builder.h b/src/svx_tree_builder.h
@@ -19,6 +19,7 @@
#define SVX_TREE_BUILDER_H
#include "svx_buffer.h"
+#include <rsys/double3.h>
#define TREE_DEPTH_MAX 16 /* Maximum depth of a tree */
@@ -55,7 +56,7 @@ struct XD(node) {
ALIGN(16) char data[NCHILDREN][SVX_MAX_SIZEOF_VOXEL]; /* Data of the leaves */
};
-/* Stacked children of an tree node */
+/* Stacked children of a tree node */
struct XD(stack) {
struct XD(node) nodes[NCHILDREN]; /* List of registered children */
uint8_t mask; /* Mask of valid children nodes (0 = empty) */
@@ -65,6 +66,11 @@ struct XD(builder) {
struct XD(stack) stacks[TREE_DEPTH_MAX];
struct buffer* buffer;
+ double lower[3]; /* Tree lower bound in world space */
+ double upper[3]; /* Tree upper bound in world space */
+ double voxsz[3]; /* Size of the finest voxel in world space */
+ enum svx_axis frame[3];
+
const struct svx_voxel_desc* desc;
size_t nleaves; /* Number of emitted leaves */
@@ -79,7 +85,7 @@ struct XD(builder) {
******************************************************************************/
#ifndef NDEBUG
static void
-XD(check)(struct buffer* buf, const struct buffer_index root)
+XD(check)(struct buffer* buf, const struct buffer_index root, size_t* nleaves)
{
const struct buffer_xnode* node;
int ichild;
@@ -94,10 +100,11 @@ XD(check)(struct buffer* buf, const struct buffer_index root)
struct buffer_index iattr;
iattr = buffer_get_child_attr_index(buf, root, ichild);
ASSERT(buffer_get_attr(buf, iattr) != NULL);
+ *nleaves += 1;
} else {
struct buffer_index child;
child = buffer_get_child_node_index(buf, root, ichild);
- XD(check)(buf, child);
+ XD(check)(buf, child, nleaves);
}
}
}
@@ -125,10 +132,20 @@ static INLINE void
XD(stack_setup_node)
(struct XD(stack)* stack,
struct XD(node)* node,
+ const double node_low[3], /* World space node lower bound */
+ const double node_sz[3], /* World space node size */
+ const size_t node_depth, /* Depth of the node */
+ const enum svx_axis frame[3],
const struct svx_voxel_desc* desc)
{
+ double child_sz[3];
+ double grandchild_sz[3];
int ichild;
+ int i;
ASSERT(stack && node && check_svx_voxel_desc(desc));
+ ASSERT(node_low && node_sz);
+ ASSERT(node_sz[0] > 0 && node_sz[1] > 0 && node_sz[2] > 0);
+ ASSERT(frame);
node->ichild_node = BUFFER_INDEX_NULL;
node->ichild_attr = BUFFER_INDEX_NULL;
@@ -137,22 +154,55 @@ XD(stack_setup_node)
if(stack->mask == 0) return; /* Empty stack */
+ d3_splat(child_sz, INF);
+ d3_splat(grandchild_sz, INF);
+ FOR_EACH(i, 0, TREE_DIMENSION) {
+ child_sz[frame[i]] = node_sz[frame[i]] * 0.5f;
+ grandchild_sz[frame[i]] = node_sz[frame[i]] * 0.25f;
+ }
+
/* Try to merge the child's leaves */
FOR_EACH(ichild, 0, NCHILDREN) {
const void* data[NCHILDREN];
+ struct svx_voxel voxels[NCHILDREN];
+ double child_low[3];
const uint8_t ichild_flag = (uint8_t)BIT(ichild);
struct XD(node)* child = stack->nodes + ichild;
int igrandchild;
+ d3_splat(child_low, -INF);
+ FOR_EACH(i, 0, TREE_DIMENSION) {
+ const int iaxis = frame[i];
+ child_low[iaxis] = node_low[iaxis];
+ if(ichild & (NCHILDREN >> (1+i))) child_low[iaxis] += child_sz[iaxis];
+ }
+
/* Fetch the grandchildren data */
FOR_EACH(igrandchild, 0, NCHILDREN) {
+ struct svx_voxel* vox = &voxels[igrandchild];
+
+ voxels[igrandchild].data = child->data[igrandchild];
+ voxels[igrandchild].depth = node_depth + 2;
+ voxels[igrandchild].id = SIZE_MAX;
+ voxels[igrandchild].is_leaf = (child->is_leaf & BIT(igrandchild))!=0;
+ voxels[igrandchild].data = child->data[igrandchild];
data[igrandchild] = child->data[igrandchild];
+
+ d3_splat(vox->lower,-INF);
+ d3_splat(vox->upper, INF);
+ FOR_EACH(i, 0, TREE_DIMENSION) {
+ const int iaxis = frame[i];
+ if(igrandchild & (NCHILDREN >> (1+i))) {
+ vox->lower[iaxis] += grandchild_sz[iaxis];
+ }
+ vox->upper[iaxis] = vox->lower[iaxis] + grandchild_sz[iaxis];
+ }
}
desc->merge(node->data[ichild], data, NCHILDREN, desc->context);
if(child->is_leaf == (BIT(NCHILDREN)-1)/*all active bitmask*/
- && desc->challenge_merge(data, NCHILDREN, desc->context)) {
+ && desc->challenge_merge(voxels, NCHILDREN, desc->context)) {
/* The node becomes a leaf : the children does not exist anymore */
node->is_leaf |= ichild_flag;
stack->mask ^= ichild_flag;
@@ -278,12 +328,16 @@ static res_T
XD(builder_init)
(struct XD(builder)* bldr,
const size_t definition,
+ const double lower[3], /* Lower bound of the tree */
+ const double upper[3], /* Upper bound of the tree */
+ const enum svx_axis frame[3],
const struct svx_voxel_desc* desc,
struct buffer* buffer)
{
- int ilvl;
+ int ilvl, i;
res_T res = RES_OK;
ASSERT(bldr && IS_POW2(definition) && check_svx_voxel_desc(desc));
+ ASSERT(lower && upper && frame);
memset(bldr, 0, sizeof(struct XD(builder)));
/* Compute the maximum depth of the tree */
@@ -303,6 +357,17 @@ XD(builder_init)
bldr->desc = desc;
bldr->buffer = buffer;
bldr->non_empty_lvl = bldr->tree_depth - 1;
+ bldr->frame[0] = frame[0];
+ bldr->frame[1] = frame[1];
+ bldr->frame[2] = frame[2];
+ d3_splat(bldr->lower,-INF);
+ d3_splat(bldr->upper, INF);
+ d3_splat(bldr->voxsz, INF);
+
+ FOR_EACH(i, 0, TREE_DIMENSION) {
+ const int iaxis = frame[i];
+ bldr->voxsz[iaxis] = (upper[iaxis] - lower[iaxis]) / (double)definition;
+ }
exit:
return res;
@@ -339,9 +404,15 @@ XD(builder_add_voxel)
* while the penultimate describes the root node itself. These 2 levels
* contain all the voxels and can be skipped */
FOR_EACH(ilvl, 0, bldr->tree_depth-2/*The 2 last voxels contain all voxels*/) {
- struct XD(node)* stack_node;
+ uint32_t parent_coords[3] = {UINT32_MAX, UINT32_MAX, UINT32_MAX};
+ double parent_sz[3];
+ double parent_low[3];
+ size_t parent_depth;
+ uint64_t parent_mcode;
+ struct XD(node)* parent_node;
uint64_t mcode_max_lvl;
size_t nleaves;
+ int i;
/* Compute the maximum morton code value for the current tree level */
mcode_max_lvl =
@@ -350,18 +421,35 @@ XD(builder_add_voxel)
if(mcode_xor < mcode_max_lvl) break;
+ /* Compute the parent node attribute */
+ parent_mcode = bldr->mcode >> (TREE_DIMENSION * (ilvl+2));
+ d3_splat(parent_sz, INF);
+ d3_splat(parent_low, -INF);
+ FOR_EACH(i, 0, TREE_DIMENSION) {
+ const int iaxis = bldr->frame[i];
+ parent_sz[iaxis] = bldr->voxsz[iaxis] * (double)(1 << (ilvl+2));
+ parent_coords[iaxis] = morton3D_decode_u21
+ (parent_mcode >> (TREE_DIMENSION-1-i));
+ parent_low[iaxis] =
+ parent_coords[iaxis] * parent_sz[iaxis] + bldr->lower[iaxis];
+ }
+
+ parent_depth = (size_t)bldr->tree_depth - 2 - ilvl;
+ ASSERT(parent_depth <= (size_t)bldr->tree_depth-2);
+
/* Retrieve the node index of the next level */
- inode = (int)(bldr->mcode >> (TREE_DIMENSION*(ilvl+2))) & (NCHILDREN-1);
+ inode = (int)parent_mcode & (NCHILDREN-1);
/* The next voxel is not in the ilvl^th stack. Setup the parent node of the
* nodes registered into the stack */
- stack_node = &bldr->stacks[ilvl+1].nodes[inode];
- XD(stack_setup_node)(&bldr->stacks[ilvl], stack_node, bldr->desc);
+ parent_node = &bldr->stacks[ilvl+1].nodes[inode];
+ XD(stack_setup_node)(&bldr->stacks[ilvl], parent_node, parent_low,
+ parent_sz, parent_depth, bldr->frame, bldr->desc);
bldr->stacks[ilvl+1].mask |= (uint8_t)BIT(inode);
/* Write the nodes of the stack of the current tree level into the buf */
res = XD(stack_write)
- (&bldr->stacks[ilvl], bldr->buffer, &stack_node->ichild_node, &nleaves);
+ (&bldr->stacks[ilvl], bldr->buffer, &parent_node->ichild_node, &nleaves);
if(res != RES_OK) goto error;
bldr->nleaves += nleaves;
@@ -410,16 +498,38 @@ XD(builder_finalize)
/* Flush the stacked nodes */
FOR_EACH(ilvl, 0, bldr->tree_depth-1) {
+ uint32_t parent_coords[3];
+ double parent_sz[3];
+ double parent_low[3];
+ size_t parent_depth;
+ uint64_t parent_mcode;
struct XD(node)* parent_node;
+ int i;
if(bldr->stacks[ilvl].mask == 0) continue;
+ /* Compute the parent node attribute */
+ parent_mcode = bldr->mcode >> (TREE_DIMENSION * (ilvl+2));
+ d3_splat(parent_sz, INF);
+ d3_splat(parent_low,-INF);
+ FOR_EACH(i, 0, TREE_DIMENSION) {
+ const int iaxis = bldr->frame[i];
+ parent_coords[iaxis] = morton3D_decode_u21
+ (parent_mcode >> (TREE_DIMENSION-1-i));
+ parent_sz[iaxis] = bldr->voxsz[iaxis] * (double)(1 << (ilvl+2));
+ parent_low[iaxis] = parent_coords[iaxis] * parent_sz[iaxis];
+ }
+
+ parent_depth = (size_t)(bldr->tree_depth - 2 - ilvl);
+ ASSERT(parent_depth <= (size_t)bldr->tree_depth-2);
+
/* Retrieve the node index of the next level */
- inode = (bldr->mcode >> (TREE_DIMENSION*(ilvl+2))) & (NCHILDREN-1);
+ inode = (int)parent_mcode & (NCHILDREN-1);
/* Setup the parent node of the nodes registered into the current stack */
parent_node = &bldr->stacks[ilvl+1].nodes[inode]; /* Fetch the parent node */
- XD(stack_setup_node)(&bldr->stacks[ilvl], parent_node, bldr->desc);
+ XD(stack_setup_node)(&bldr->stacks[ilvl], parent_node, parent_low,
+ parent_sz, parent_depth, bldr->frame, bldr->desc);
bldr->stacks[ilvl+1].mask |= (uint8_t)BIT(inode);
/* Write the stacked nodes of the current level */
diff --git a/src/test_svx_bintree.c b/src/test_svx_bintree.c
@@ -36,6 +36,20 @@ struct at_context {
size_t depth;
};
+struct aabb {
+ double lower[3];
+ double upper[3];
+ size_t depth;
+};
+
+struct build_context {
+ double voxsz[3];
+ double lower[3];
+ double upper[3];
+ size_t max_depth;
+};
+
+
static void
get(const size_t xyz[3], void* dst, void* ctx)
{
@@ -47,7 +61,7 @@ get(const size_t xyz[3], void* dst, void* ctx)
}
static int
-no_merge(const void* voxels[], const size_t nvoxels, void* ctx)
+no_merge(const struct svx_voxel voxels[], const size_t nvoxels, void* ctx)
{
CHK(voxels != NULL);
CHK(nvoxels != 0);
@@ -122,6 +136,95 @@ check_leaves
}
static void
+get_aabb(const size_t xyz[3], void* dst, void* ctx)
+{
+ const struct build_context* build_ctx = ctx;
+ struct aabb* aabb = dst;
+
+ d3_splat(aabb->lower,-INF);
+ d3_splat(aabb->upper, INF);
+ aabb->lower[AXIS] =
+ (double)xyz[AXIS] * build_ctx->voxsz[AXIS] + build_ctx->lower[AXIS];
+ aabb->upper[AXIS] = aabb->lower[AXIS] + build_ctx->voxsz[AXIS];
+}
+
+
+static void
+merge_aabb(void* dst, const void* voxels[], const size_t nvoxels, void* ctx)
+{
+ const struct build_context* build_ctx = ctx;
+ double voxsz[3];
+ struct aabb* aabb = dst;
+ size_t depth = SIZE_MAX;
+ size_t i;
+ CHK(dst && voxels && nvoxels && ctx);
+
+ d3_splat(aabb->lower,-INF);
+ d3_splat(aabb->upper, INF);
+ aabb->depth = 0;
+ aabb->lower[AXIS] = DBL_MAX;
+ aabb->upper[AXIS] =-DBL_MAX;
+
+ FOR_EACH(i, 0, nvoxels) {
+ const struct aabb* vox_aabb = voxels[i];
+ if(depth == SIZE_MAX) {
+ depth = vox_aabb->depth;
+ } else if(vox_aabb->depth != 0) { /* Not invalid */
+ CHK(depth == vox_aabb->depth);
+ aabb->lower[AXIS] = MMIN(vox_aabb->lower[AXIS], aabb->lower[AXIS]);
+ }
+ }
+
+ if(depth == 0) { /* Invalid voxel */
+ aabb->lower[AXIS] = 0;
+ aabb->upper[AXIS] = 0;
+ aabb->depth = 0;
+ } else {
+ double upper[3];
+ CHK(build_ctx->max_depth >= depth);
+ CHK(depth > 0);
+ aabb->depth = depth - 1;
+
+ i = build_ctx->max_depth - aabb->depth;
+ voxsz[AXIS] = build_ctx->voxsz[AXIS] * (double)(1<<i);
+
+ /* Clamp voxel to grid size */
+ upper[AXIS] = MMIN(aabb->lower[AXIS] + voxsz[AXIS], build_ctx->upper[AXIS]);
+
+ /* Adjust the voxel size from the clampd voxel */
+ voxsz[AXIS] = upper[AXIS] - aabb->lower[AXIS];
+
+ CHK(eq_eps(voxsz[AXIS], aabb->upper[AXIS] - aabb->lower[AXIS], 1.e-6));
+ }
+}
+
+static int
+challenge_aabb(const struct svx_voxel voxels[], const size_t nvoxels, void* ctx)
+{
+ const struct build_context* build_ctx = ctx;
+ size_t depth = SIZE_MAX;
+ size_t i;
+ CHK(voxels && nvoxels && ctx);
+
+ FOR_EACH(i, 0, nvoxels) {
+ double voxsz[3];
+ const struct aabb* aabb = voxels[i].data;
+ const size_t n = build_ctx->max_depth - aabb->depth;
+ CHK(build_ctx->max_depth >= aabb->depth);
+
+ if(depth == SIZE_MAX) {
+ depth = aabb->depth;
+ } else if(aabb->depth != 0) { /* Not invalid */
+ CHK(depth == aabb->depth);
+ CHK(depth == voxels[i].depth);
+ voxsz[AXIS] = build_ctx->voxsz[AXIS] * (double)(1<<n);
+ CHK(eq_eps(voxsz[AXIS], aabb->upper[AXIS] - aabb->lower[AXIS], 1.e-6));
+ }
+ }
+ return 1;
+}
+
+static void
test_at_accessor(struct svx_tree* tree, const size_t nvoxels)
{
struct svx_tree_desc desc;
@@ -185,6 +288,7 @@ main(int argc, char** argv)
struct mem_allocator allocator;
struct svx_voxel_desc vox_desc = SVX_VOXEL_DESC_NULL;
struct svx_tree_desc tree_desc = SVX_TREE_DESC_NULL;
+ struct build_context build_ctx;
struct leaves_context ctx;
enum svx_axis axis;
double low, upp;
@@ -267,8 +371,17 @@ main(int argc, char** argv)
CHK(tree_desc.upper[axis] == upp);
test_at_accessor(tree, nvxls);
+ CHK(svx_tree_ref_put(tree) == RES_OK);
+
+ vox_desc.get = get_aabb;
+ vox_desc.merge = merge_aabb;
+ vox_desc.challenge_merge = challenge_aabb;
+ vox_desc.context = &build_ctx;
+ vox_desc.size = sizeof(struct aabb);
+
+ CHK(NEW_TREE(dev, low, upp, nvxls, axis, &vox_desc, &tree) == RES_OK);
- #undef NEW_SCN
+ #undef NEW_TREE
CHK(svx_tree_ref_put(tree) == RES_OK);
CHK(svx_device_ref_put(dev) == RES_OK);
diff --git a/src/test_svx_bintree_trace_ray.c b/src/test_svx_bintree_trace_ray.c
@@ -82,7 +82,8 @@ voxels_merge(void* dst, const void* voxels[], const size_t nvoxels, void* ctx)
}
static int
-voxels_challenge_merge(const void* voxels[], const size_t nvoxels, void* ctx)
+voxels_challenge_merge
+ (const struct svx_voxel voxels[], const size_t nvoxels, void* ctx)
{
size_t ivoxel;
int merge = 1;
@@ -90,10 +91,10 @@ voxels_challenge_merge(const void* voxels[], const size_t nvoxels, void* ctx)
CHK(voxels && nvoxels && ctx);
- ref = *(char*)(voxels[0]);
+ ref = *(char*)(voxels[0].data);
for(ivoxel=1; merge && ivoxel<nvoxels; ++ivoxel) {
- const char* voxel_data = voxels[ivoxel];
+ const char* voxel_data = voxels[ivoxel].data;
merge = (ref == *voxel_data);
}
return merge;
diff --git a/src/test_svx_octree.c b/src/test_svx_octree.c
@@ -19,6 +19,8 @@
#include <rsys/double3.h>
+#include <string.h>
+
struct leaves_context {
double* lower;
double* upper;
@@ -31,8 +33,21 @@ struct at_context {
size_t depth;
};
+struct aabb {
+ double lower[3];
+ double upper[3];
+ size_t depth;
+};
+
+struct build_context {
+ double voxsz[3];
+ double lower[3];
+ double upper[3];
+ size_t max_depth;
+};
+
static int
-no_merge(const void* voxels[], const size_t nvoxels, void* ctx)
+no_merge(const struct svx_voxel voxels[], const size_t nvoxels, void* ctx)
{
CHK(voxels != NULL);
CHK(nvoxels != 0);
@@ -41,7 +56,7 @@ no_merge(const void* voxels[], const size_t nvoxels, void* ctx)
}
static int
-merge_level0(const void** voxels, const size_t nvoxels, void* ctx)
+merge_level0(const struct svx_voxel voxels[], const size_t nvoxels, void* ctx)
{
double min_val = DBL_MAX;
double max_val =-DBL_MAX;
@@ -51,11 +66,17 @@ merge_level0(const void** voxels, const size_t nvoxels, void* ctx)
CHK((intptr_t)ctx == 0xDECAFBAD);
FOR_EACH(i, 0, nvoxels) {
- const double* val = voxels[i];
+ const double* val = voxels[i].data;
min_val = MMIN(min_val, *val);
max_val = MMAX(max_val, *val);
}
+ FOR_EACH(i, 0, nvoxels) {
+ if(max_val - min_val < 8) {
+ CHK(voxels[i].depth == 5);
+ }
+ }
+
return (max_val - min_val) < 8;
}
@@ -171,6 +192,109 @@ max_lod
}
static void
+get_aabb(const size_t xyz[3], void* dst, void* ctx)
+{
+ const struct build_context* build_ctx = ctx;
+ struct aabb* aabb = dst;
+
+ aabb->lower[0] = (double)xyz[0] * build_ctx->voxsz[0] + build_ctx->lower[0];
+ aabb->lower[1] = (double)xyz[1] * build_ctx->voxsz[1] + build_ctx->lower[1];
+ aabb->lower[2] = (double)xyz[2] * build_ctx->voxsz[2] + build_ctx->lower[2];
+ d3_add(aabb->upper, aabb->lower, build_ctx->voxsz);
+ aabb->depth = build_ctx->max_depth;
+}
+
+static void
+merge_aabb(void* dst, const void* voxels[], const size_t nvoxels, void* ctx)
+{
+ const struct build_context* build_ctx = ctx;
+ double voxsz[3];
+ struct aabb* aabb = dst;
+ size_t depth = SIZE_MAX;
+ size_t i;
+ CHK(dst && voxels && nvoxels && ctx);
+
+ d3_splat(aabb->lower, DBL_MAX);
+ d3_splat(aabb->upper,-DBL_MAX);
+ aabb->depth = 0;
+
+ FOR_EACH(i, 0, nvoxels) {
+ const struct aabb* vox_aabb = voxels[i];
+ if(depth == SIZE_MAX) {
+ depth = vox_aabb->depth;
+ } else if(vox_aabb->depth != 0) { /* Not invalid */
+ CHK(depth == vox_aabb->depth);
+ aabb->lower[0] = MMIN(vox_aabb->lower[0], aabb->lower[0]);
+ aabb->lower[1] = MMIN(vox_aabb->lower[1], aabb->lower[1]);
+ aabb->lower[2] = MMIN(vox_aabb->lower[2], aabb->lower[2]);
+ aabb->upper[0] = MMAX(vox_aabb->upper[0], aabb->upper[0]);
+ aabb->upper[1] = MMAX(vox_aabb->upper[1], aabb->upper[1]);
+ aabb->upper[2] = MMAX(vox_aabb->upper[2], aabb->upper[2]);
+ }
+ }
+
+ if(depth == 0) { /* Invalid voxel */
+ memset(aabb, 0, sizeof(*aabb));
+ } else {
+ double upper[3];
+ CHK(build_ctx->max_depth >= depth);
+ CHK(depth > 0);
+ aabb->depth = depth - 1;
+
+ i = build_ctx->max_depth - aabb->depth;
+ voxsz[0] = build_ctx->voxsz[0] * (double)(1<<i);
+ voxsz[1] = build_ctx->voxsz[1] * (double)(1<<i);
+ voxsz[2] = build_ctx->voxsz[2] * (double)(1<<i);
+
+ /* Clamp voxel to grid size */
+ upper[0] = MMIN(aabb->lower[0] + voxsz[0], build_ctx->upper[0]);
+ upper[1] = MMIN(aabb->lower[1] + voxsz[1], build_ctx->upper[1]);
+ upper[2] = MMIN(aabb->lower[2] + voxsz[2], build_ctx->upper[2]);
+
+ /* Adjust the voxel size from the clampd voxel */
+ voxsz[0] = upper[0] - aabb->lower[0];
+ voxsz[1] = upper[1] - aabb->lower[1];
+ voxsz[2] = upper[2] - aabb->lower[2];
+
+ CHK(eq_eps(voxsz[0], aabb->upper[0] - aabb->lower[0], 1.e-6));
+ CHK(eq_eps(voxsz[1], aabb->upper[1] - aabb->lower[1], 1.e-6));
+ CHK(eq_eps(voxsz[2], aabb->upper[2] - aabb->lower[2], 1.e-6));
+ }
+}
+
+static int
+challenge_aabb(const struct svx_voxel voxels[], const size_t nvoxels, void* ctx)
+{
+ const struct build_context* build_ctx = ctx;
+ size_t depth = SIZE_MAX;
+ size_t i;
+ CHK(voxels && nvoxels && ctx);
+
+ FOR_EACH(i, 0, nvoxels) {
+ double voxsz[3];
+ const struct aabb* aabb = voxels[i].data;
+ const size_t n = build_ctx->max_depth - aabb->depth;
+ CHK(build_ctx->max_depth >= aabb->depth);
+
+ if(depth == SIZE_MAX) {
+ depth = aabb->depth;
+ } else if(aabb->depth != 0) { /* Not invalid */
+ CHK(depth == aabb->depth);
+ CHK(depth == voxels[i].depth);
+
+ voxsz[0] = build_ctx->voxsz[0] * (double)(1<<n);
+ voxsz[1] = build_ctx->voxsz[1] * (double)(1<<n);
+ voxsz[2] = build_ctx->voxsz[2] * (double)(1<<n);
+ CHK(eq_eps(voxsz[0], aabb->upper[0] - aabb->lower[0], 1.e-6));
+ CHK(eq_eps(voxsz[1], aabb->upper[1] - aabb->lower[1], 1.e-6));
+ CHK(eq_eps(voxsz[2], aabb->upper[2] - aabb->lower[2], 1.e-6));
+ }
+ }
+
+ return 1;
+}
+
+static void
test_at_accessor(struct svx_tree* oct, const size_t nvoxels[3])
{
struct svx_tree_desc tree_desc;
@@ -268,6 +392,7 @@ main(int argc, char** argv)
struct mem_allocator allocator;
struct svx_voxel_desc voxdesc = SVX_VOXEL_DESC_NULL;
struct svx_tree_desc tree_desc = SVX_TREE_DESC_NULL;
+ struct build_context build_ctx;
double low[3];
double upp[3];
size_t nvxls[3];
@@ -369,6 +494,29 @@ main(int argc, char** argv)
dump_data(stdout, oct, TYPE_FLOAT, 1, write_scalars);
+ CHK(svx_tree_ref_put(oct) == RES_OK);
+
+ nvxls[0] = 32;
+ nvxls[1] = 16;
+ nvxls[2] = 33;
+
+ build_ctx.max_depth = (size_t)log2i
+ ((int)round_up_pow2(MMAX(MMAX(nvxls[0], nvxls[1]), nvxls[2])));
+
+ d3_set(build_ctx.lower, low);
+ d3_set(build_ctx.upper, upp);
+ build_ctx.voxsz[0] = (upp[0]-low[0])/(double)nvxls[0];
+ build_ctx.voxsz[1] = (upp[1]-low[1])/(double)nvxls[1];
+ build_ctx.voxsz[2] = (upp[2]-low[2])/(double)nvxls[2];
+
+ voxdesc.get = get_aabb;
+ voxdesc.merge = merge_aabb;
+ voxdesc.challenge_merge = challenge_aabb;
+ voxdesc.context = &build_ctx;
+ voxdesc.size = sizeof(struct aabb);
+
+ CHK(NEW_SCN(dev, low, upp, nvxls, &voxdesc, &oct) == RES_OK);
+
#undef NEW_SCN
CHK(svx_device_ref_put(dev) == RES_OK);
diff --git a/src/test_svx_octree_trace_ray.c b/src/test_svx_octree_trace_ray.c
@@ -97,7 +97,8 @@ voxels_merge(void* dst, const void* voxels[], const size_t nvoxels, void* ctx)
}
static int
-voxels_challenge_merge(const void* voxels[], const size_t nvoxels, void* ctx)
+voxels_challenge_merge
+ (const struct svx_voxel voxels[], const size_t nvoxels, void* ctx)
{
size_t ivoxel;
int merge = 1;
@@ -105,10 +106,10 @@ voxels_challenge_merge(const void* voxels[], const size_t nvoxels, void* ctx)
CHK(voxels && nvoxels && ctx);
- ref = *(char*)(voxels[0]);
+ ref = *(char*)(voxels[0].data);
for(ivoxel=1; merge && ivoxel<nvoxels; ++ivoxel) {
- const char* voxel_data = voxels[ivoxel];
+ const char* voxel_data = voxels[ivoxel].data;
merge = (ref == *voxel_data);
}
return merge;