commit a5c9569f03edf67d648f848cd1a6170f29116bc8
parent fa6104395702ed57a162490fd6152d2e0b43077f
Author: Vincent Forest <vincent.forest@meso-star.com>
Date: Fri, 11 Dec 2020 18:51:32 +0100
Add the suvm_polyhedron_intersect_aabb function
Refactoring of the suvm_primitive related code
Diffstat:
8 files changed, 493 insertions(+), 349 deletions(-)
diff --git a/cmake/CMakeLists.txt b/cmake/CMakeLists.txt
@@ -47,6 +47,7 @@ set(VERSION ${VERSION_MAJOR}.${VERSION_MINOR}.${VERSION_PATCH})
set(SUVM_FILES_SRC
suvm_device.c
+ suvm_primitive.c
suvm_volume.c
suvm_volume_at.c
suvm_volume_intersect_aabb.c)
diff --git a/src/suvm.h b/src/suvm.h
@@ -17,6 +17,7 @@
#define SUVM_H
#include <rsys/rsys.h>
+#include <float.h>
/* Library symbol management */
#if defined(SUVM_SHARED_BUILD) /* Build shared library */
@@ -44,19 +45,48 @@ struct suvm_data {
size_t size; /* Size of the data in bytes */
size_t alignment; /* Alignment of the data */
};
-static const struct suvm_data SUVM_DATA_NULL;
+#define SUVM_DATA_NULL__ {NULL, 0, 0}
+static const struct suvm_data SUVM_DATA_NULL = SUVM_DATA_NULL__;
struct suvm_primitive {
const void* data; /* Data of the primitive */
const void* vertex_data[SUVM_PRIMITIVE_MAX_VERTICES_COUNT]; /* Vertex data */
size_t iprim; /* Identifier of the primitive */
size_t nvertices; /* #vertices of the primitive */
+
+ /* Internal data */
+ const struct suvm_volume* volume__;
};
-#define SUVM_PRIMITIVE_NULL__ {NULL, {NULL}, SIZE_MAX, SIZE_MAX}
+#define SUVM_PRIMITIVE_NULL__ { \
+ NULL, /* Primitive data */ \
+ {NULL}, /* Vertex data */ \
+ SIZE_MAX, /* Primitive id */ \
+ SIZE_MAX, /* #vertices */ \
+ NULL /* Pointer toward its associated volume */ \
+}
static const struct suvm_primitive SUVM_PRIMITIVE_NULL = SUVM_PRIMITIVE_NULL__;
#define SUVM_PRIMITIVE_NONE(Prim) ((Prim)->iprim == SIZE_MAX)
+/* Precomputed data from a suvm_primitive used to speed up AABB/primitive
+ * intersection tests */
+struct suvm_polyhedron {
+ float v[4/*#vertices*/][3/*#coords*/]; /* Vertices */
+ float N[4/*#vertices*/][3/*#coords*/]; /* Normals */
+ float D[4/*#facets*/]; /* Slope parameter of the plane (D = -(Ax+By+Cz)) */
+
+ /* 2D equation of the silhouette edges projected along the X, Y and Z axis of
+ * the form Ax + By + C = 0, with (A, B) the normals of the edge */
+ float Ep[3/*#axes*/][4/*#edges max*/][3/*#equation parameters*/];
+
+ /* Number of silhouette edges in the ZY, XZ and XY planes (3 or 4) */
+ int nEp[3/*#coords*/];
+
+ /* Tetrahedron axis aligned bounding box */
+ float lower[3/*#coords*/];
+ float upper[3/*#coords*/];
+};
+
struct suvm_tetrahedral_mesh_args {
size_t ntetrahedra; /* #tetrahedra */
size_t nvertices; /* #vertices */
@@ -78,7 +108,11 @@ struct suvm_tetrahedral_mesh_args {
void* context; /* Client data set as the last param of the callbacks */
};
-static const struct suvm_tetrahedral_mesh_args SUVM_TETRAHEDRAL_MESH_ARGS_NULL;
+#define SUVM_TETRAHEDRAL_MESH_ARGS_NULL__ { \
+ 0, 0, NULL, NULL, SUVM_DATA_NULL__, SUVM_DATA_NULL__, 0, NULL \
+}
+static const struct suvm_tetrahedral_mesh_args SUVM_TETRAHEDRAL_MESH_ARGS_NULL =
+ SUVM_TETRAHEDRAL_MESH_ARGS_NULL__;
/* Callback invoked on suvm_volume_intersect_aabb invocation on each primitives
* intersected by the submitted Axis Aligned Bounding Box */
@@ -162,6 +196,31 @@ suvm_volume_intersect_aabb
suvm_primitive_intersect_aabb_T clbk,
void* context); /* User data sent as the last argument of clbk */
+SUVM_API res_T
+suvm_volume_get_primitives_count
+ (const struct suvm_volume* volume,
+ size_t* nprims);
+
+SUVM_API res_T
+suvm_volume_get_primitive
+ (const struct suvm_volume* volume,
+ const size_t iprim, /* In [0, suvm_volume_get_primitives_count[ */
+ struct suvm_primitive* prim);
+
+/*******************************************************************************
+ * Primitive API
+ ******************************************************************************/
+SUVM_API res_T
+suvm_primitive_setup_polyhedron
+ (const struct suvm_primitive* prim,
+ struct suvm_polyhedron* poly);
+
+SUVM_API int
+suvm_polyhedron_intersect_aabb
+ (const struct suvm_polyhedron* poly,
+ const float low[3],
+ const float upp[3]);
+
END_DECLS
#endif /* SUVM_H */
diff --git a/src/suvm_primitive.c b/src/suvm_primitive.c
@@ -0,0 +1,339 @@
+/* Copyright (C) 2020 |Meso|Star> (contact@meso-star.com)
+ *
+ * This program is free software: you can redistribute it and/or modify
+ * it under the terms of the GNU General Public License as published by
+ * the Free Software Foundation, either version 3 of the License, or
+ * (at your option) any later version.
+ *
+ * This program is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
+ * GNU General Public License for more details.
+ *
+ * You should have received a copy of the GNU General Public License
+ * along with this program. If not, see <http://www.gnu.org/licenses/>. */
+
+#include "suvm.h"
+#include "suvm_volume.h"
+
+#include <rsys/float2.h>
+#include <rsys/float3.h>
+
+/*******************************************************************************
+ * Exported functions
+ ******************************************************************************/
+res_T
+suvm_primitive_setup_polyhedron
+ (const struct suvm_primitive* prim,
+ struct suvm_polyhedron* poly)
+{
+ if(!prim || !poly) return RES_BAD_ARG;
+ tetrahedron_setup(prim->volume__, NULL, NULL, prim->iprim, poly);
+ return RES_OK;
+}
+
+/* Define if an axis aligned bounding box and a tetrahedron are intersecting.
+ * Its implementation follows the algorithms proposed by N. Greene in
+ * "Detecting intersection of a rectangular solid and a convex polyhedron"
+ * (Graphics Gems IV, 1994, p74--82) */
+int
+suvm_polyhedron_intersect_aabb
+ (const struct suvm_polyhedron* tetra,
+ const float low[3],
+ const float upp[3])
+{
+ enum { X, Y, Z, AXES_COUNT }; /* Syntactic sugar */
+ float intersect_aabb_low[AXES_COUNT];
+ float intersect_aabb_upp[AXES_COUNT];
+ int i;
+ int nplanes_including_aabb;
+
+ ASSERT(tetra && low && !upp);
+ ASSERT(low[0] < upp[0]);
+ ASSERT(low[1] < upp[1]);
+ ASSERT(low[2] < upp[2]);
+
+ /* Check the intersection between the AABB and the tetra AABB */
+ FOR_EACH(i, 0, AXES_COUNT) {
+ intersect_aabb_low[i] = MMAX(low[i], tetra->lower[i]);
+ intersect_aabb_upp[i] = MMIN(upp[i], tetra->upper[i]);
+ if(intersect_aabb_low[i] > intersect_aabb_upp[i]) /* Do not intersect */
+ return 0;
+ }
+
+ /* Check if the tetrahedron is included into the aabb */
+ if(intersect_aabb_low[X] == low[X]
+ && intersect_aabb_low[Y] == low[Y]
+ && intersect_aabb_low[Z] == low[Z]) {
+ return 1;
+ }
+
+ /* #planes that totally includes the aabb on its positive side */
+ nplanes_including_aabb = 0;
+
+ /* Check the tetrahedron planes against the aabb */
+ FOR_EACH(i, 0, 4) {
+ float p[3], n[3];
+
+ /* Define the aabb vertices that is the farthest in the positive/negative
+ * direction of the face normal. */
+ f3_set(p, upp);
+ f3_set(n, low);
+ if(tetra->N[i][X] < 0) SWAP(float, p[X], n[X]);
+ if(tetra->N[i][Y] < 0) SWAP(float, p[Y], n[Y]);
+ if(tetra->N[i][Z] < 0) SWAP(float, p[Z], n[Z]);
+
+ /* Check that the box is totally outside the plane, To do this, check that
+ * 'p', aka the farthest aabb vertex in the positive direction of the plane
+ * normal, is not in the positive half space. In this case, the aabb is
+ * entirely outside the tetrahedron */
+ if((f3_dot(tetra->N[i], p) + tetra->D[i]) < 0) return 0;
+
+ /* Check if the box is totally inside the given plane. To do this, check
+ * that 'n', aka the farthest aabb vertex in the negative direction of the
+ * plane normal, is inside the positive half space */
+ if((f3_dot(tetra->N[i], n) + tetra->D[i]) > 0) {
+ /* Register this plane as a plane including the aabb */
+ nplanes_including_aabb += 1;
+ }
+ }
+
+ /* Check if the aabb is entirely included into the tetrahedron */
+ if(nplanes_including_aabb == 4) return 1;
+
+ /* For each silhouette edge in each projection plane, check if it totally
+ * exclude the projected aabb */
+ FOR_EACH(i, 0, AXES_COUNT) {
+ /* Compute the remaining axes.
+ * On 1st iteration 'j' and 'k' define the YZ plane
+ * On 2nd iteration 'j' and 'k' define the ZX plane
+ * On 3rd ietration 'j' and 'k' define the XY plane */
+ const int j = (i + 1) % AXES_COUNT;
+ const int k = (j + 1) % AXES_COUNT;
+
+ int iedge;
+ FOR_EACH(iedge, 0, tetra->nEp[i]) {
+ float p[2];
+
+ /* Define the aabb vertices that is the farthest in the positive/negative
+ * direction of the normal of the silhouette edge. */
+ p[0] = tetra->Ep[i][iedge][0] > 0 ? upp[j] : low[j];
+ p[1] = tetra->Ep[i][iedge][1] > 0 ? upp[k] : low[k];
+
+ /* Check that the projected aabb along the 'i'th axis is totally outside
+ * the current silhouette edge. That means that the aabb is entirely
+ * outside the tetrahedron and thus that they do not intersect */
+ if((f2_dot(tetra->Ep[i][iedge], p) + tetra->Ep[i][iedge][2]) < 0) return 0;
+ }
+ }
+
+ /* The tetrahedron and the aabb are partially intersecting */
+ return 1;
+}
+
+
+/*******************************************************************************
+ * Local functions
+ ******************************************************************************/
+void
+tetrahedron_setup
+ (const struct suvm_volume* vol,
+ const float lower[3],
+ const float upper[3],
+ const size_t itetra,
+ struct suvm_polyhedron* tetra)
+{
+ enum { X, Y, Z, AXES_COUNT }; /* Syntactic sugar */
+ const float* low = lower;
+ const float* upp = upper;
+ float low__[3];
+ float upp__[3];
+#ifndef NDEBUG
+ float center[AXES_COUNT]; /* Center of the tetrahedron */
+#endif
+ float e[6][AXES_COUNT];
+ float (*v)[AXES_COUNT];
+ float (*N)[AXES_COUNT];
+ float (*Ep)[4][3];
+ int i;
+ ASSERT(vol && tetra);
+
+ /* Fetch tetrahedron vertices */
+ volume_primitive_get_vertex_position(vol, itetra, 0, tetra->v[0]);
+ volume_primitive_get_vertex_position(vol, itetra, 1, tetra->v[1]);
+ volume_primitive_get_vertex_position(vol, itetra, 2, tetra->v[2]);
+ volume_primitive_get_vertex_position(vol, itetra, 4, tetra->v[3]);
+ v = tetra->v;
+
+#ifndef NDEBUG
+ /* Compute the center of the tetrahedron */
+ center[X] = (v[0][X] + v[1][X] + v[2][X] + v[3][X]) * 0.25f;
+ center[Y] = (v[0][Y] + v[1][Y] + v[2][Y] + v[3][Y]) * 0.25f;
+ center[Z] = (v[0][Z] + v[1][Z] + v[2][Z] + v[3][Z]) * 0.25f;
+#endif
+
+ /* Define the primitive AABB if necessary */
+ if(!low) {
+ low__[0] = MMIN(MMIN(v[0][X], v[1][X]), MMIN(v[2][X], v[3][X]));
+ low__[1] = MMIN(MMIN(v[0][Y], v[1][Y]), MMIN(v[2][Y], v[3][Y]));
+ low__[2] = MMIN(MMIN(v[0][Z], v[1][Z]), MMIN(v[2][Z], v[3][Z]));
+ low = low__;
+ }
+ if(!upp) {
+ upp__[0] = MMAX(MMAX(v[0][X], v[1][X]), MMAX(v[2][X], v[3][X]));
+ upp__[1] = MMAX(MMAX(v[0][Y], v[1][Y]), MMAX(v[2][Y], v[3][Y]));
+ upp__[2] = MMAX(MMAX(v[0][Z], v[1][Z]), MMAX(v[2][Z], v[3][Z]));
+ upp = upp__;
+ }
+
+#ifndef NDEBUG
+ /* Check argument consistency */
+ if(lower) {
+ ASSERT(MMIN(MMIN(v[0][X], v[1][X]), MMIN(v[2][X], v[3][X])) == lower[X]);
+ ASSERT(MMIN(MMIN(v[0][Y], v[1][Y]), MMIN(v[2][Y], v[3][Y])) == lower[Y]);
+ ASSERT(MMIN(MMIN(v[0][Z], v[1][Z]), MMIN(v[2][Z], v[3][Z])) == lower[Z]);
+ }
+ if(upper) {
+ ASSERT(MMAX(MMAX(v[0][X], v[1][X]), MMAX(v[2][X], v[3][X])) == upper[X]);
+ ASSERT(MMAX(MMAX(v[0][Y], v[1][Y]), MMAX(v[2][Y], v[3][Y])) == upper[Y]);
+ ASSERT(MMAX(MMAX(v[0][Z], v[1][Z]), MMAX(v[2][Z], v[3][Z])) == upper[Z]);
+ }
+#endif
+
+ /* Setup tetrahedron AABB */
+ f3_set(tetra->lower, low);
+ f3_set(tetra->upper, upp);
+
+ /* Fetch tetrahedron normals */
+ volume_primitive_get_facet_normal(vol, itetra, 0, tetra->N[0]);
+ volume_primitive_get_facet_normal(vol, itetra, 1, tetra->N[1]);
+ volume_primitive_get_facet_normal(vol, itetra, 2, tetra->N[2]);
+ volume_primitive_get_facet_normal(vol, itetra, 3, tetra->N[3]);
+ N = tetra->N;
+
+ /* Compute the slope of the planes */
+ tetra->D[0] = -f3_dot(tetra->N[0], v[0]);
+ tetra->D[1] = -f3_dot(tetra->N[1], v[1]);
+ tetra->D[2] = -f3_dot(tetra->N[2], v[2]);
+ tetra->D[3] = -f3_dot(tetra->N[3], v[3]);
+
+ /* Compute tetrahedron edges */
+ f3_sub(e[0], v[1], v[0]);
+ f3_sub(e[1], v[2], v[1]);
+ f3_sub(e[2], v[0], v[2]);
+ f3_sub(e[3], v[0], v[3]);
+ f3_sub(e[4], v[1], v[3]);
+ f3_sub(e[5], v[2], v[3]);
+
+ /* Detect the silhouette edges of the tetrahedron once projected in the YZ,
+ * ZX and XY planes, and compute their 2D equation in the projected plane. The
+ * variable 'i' defines the projection axis which is successively X (YZ
+ * plane), Y (ZX plane) and Z (XY plane). The axes of the corresponding 2D
+ * repairs are defined by the 'j' and 'k' variables. */
+ Ep = tetra->Ep;
+ FOR_EACH(i, 0, AXES_COUNT) {
+ /* On 1st iteration 'j' and 'k' define the YZ plane
+ * On 2nd iteration 'j' and 'k' define the ZX plane
+ * On 3rd ietration 'j' and 'k' define the XY plane */
+ const int j = (i + 1) % AXES_COUNT;
+ const int k = (j + 1) % AXES_COUNT;
+
+ /* Register the number of detected silhouette edges */
+ int n = 0;
+
+ /* To detect the silhouette edges, check the sign of the normals of two
+ * adjacent facets for the coordinate of the projection axis. If the signs
+ * are the same, the facets look at the same direction regarding the
+ * projection axis. The other case means that one facet points toward the
+ * projection axis whole the other one points to the opposite direction:
+ * this is the definition of a silhouette edge.
+ *
+ * Once detected, we compute the equation of the silhouette edge:
+ *
+ * Ax + By + C = 0
+ *
+ * with A and B the coordinates of the edge normals and C the slope of the
+ * edge. Computing the normal is actually as simple as swizzling the
+ * coordinates of the edges projected in 2D and sign tuning. Let the
+ * projection axis X (i==X) and thus the YZ projection plane (j==Y &&
+ * k==Z). Assuming that the first edge 'e[0]' is a silhouette edge, ie the
+ * edge between the facet0 and the facet1 (refers to the suvm_volume.h for
+ * the geometric layout of a tetrahedron) once projected in the YZ plane
+ * the edge is {0, e[0][Y], e[0][Z]}. Its normal N can is defined as the
+ * cross product between this projected edge ands the projection axis:
+ *
+ * | 0 | | 1 | | 0 |
+ * N = | e[0][Y] | X | 0 | = | e[0][Z] |
+ * | e[0][Z] | | 0 | |-e[0][Y] |
+ *
+ * In the projection plane, the _unormalized_ normal of the {e[0][Y],
+ * e[0][Z]} edge is thus {e[0][Z], -e[0][Y]}. More generally, the normal of
+ * an edge 'I' {e[I][j], e[I][k]} projected along the 'i' axis in the 'jk'
+ * plane is {e[I][k],-e[I]j]}. Anyway, one has to revert the normal
+ * orientation to ensure a specific convention wrt the tetrahedron
+ * footprint. In our case we ensure that the normal points toward the
+ * footprint. Finalle the edge slope 'C' can be computed as:
+ * | A |
+ * C =- | B | . | Pj, Pk |
+ *
+ * with {Pj, Pk} a point belonging to the edge as for instance one of its
+ * vertices. */
+ if(signf(N[0][i]) != signf(N[1][i])) { /* The edge 0 is silhouette */
+ f2_normalize(Ep[i][n], f2(Ep[i][n], e[0][k], -e[0][j]));
+ if(N[0][i] > 0) f2_minus(Ep[i][n], Ep[i][n]);
+ Ep[i][n][2] = -f2_dot(Ep[i][n], v[0]);
+ ++n;
+ }
+ if(signf(N[0][i]) != signf(N[2][i])) { /* The edge 1 is silhouette */
+ f2_normalize(Ep[i][n], f2(Ep[i][n], e[1][k], -e[1][j]));
+ if(N[0][i] > 0) f2_minus(Ep[i][n], Ep[i][n]);
+ Ep[i][n][2] = -f2_dot(Ep[i][n], v[1]);
+ ++n;
+ }
+ if(signf(N[0][i]) != signf(N[3][i])) { /* The edge 2 is silhouette */
+ f2_normalize(Ep[i][n], f2(Ep[i][n], e[2][k], -e[2][j]));
+ if(N[0][i] > 0) f2_minus(Ep[i][n], Ep[i][n]);
+ Ep[i][n][2] = -f2_dot(Ep[i][n], v[2]);
+ ++n;
+ }
+ if(signf(N[1][i]) != signf(N[3][i])) { /* The edge 3 is silhouette */
+ f2_normalize(Ep[i][n], f2(Ep[i][n], e[3][k], -e[3][j]));
+ if(N[1][i] > 0) f2_minus(Ep[i][n], Ep[i][n]);
+ Ep[i][n][2] = -f2_dot(Ep[i][n], v[0]);
+ ++n;
+ }
+ if(signf(N[1][i]) != signf(N[2][i])) { /* The edge 4 is silhouette */
+ f2_normalize(Ep[X][n], f2(Ep[i][n], e[4][k], -e[4][j]));
+ if(N[2][i] > 0) f2_minus(Ep[i][n], Ep[i][n]);
+ Ep[i][n][2] = -f2_dot(Ep[i][n], v[1]);
+ ++n;
+ }
+ if(signf(N[2][i]) != signf(N[3][i])) { /* The edge 5 is silhouette */
+ f2_normalize(Ep[i][n], f2(Ep[i][n], e[5][k], -e[5][j]));
+ if(N[3][i] > 0) f2_minus(Ep[i][n], Ep[i][n]);
+ Ep[i][n][2] = -f2_dot(Ep[i][n], v[2]);
+ ++n;
+ }
+ /* A tetrahedron can have 3 or 4 silhouette edges once projected in 2D */
+ ASSERT(n == 3 || n == 4);
+ tetra->nEp[i] = n; /* Register the #silouhette edges for this project axis */
+
+#ifndef NDEBUG
+ /* Assert that the normals points toward the expected half space */
+ {
+ int iedge;
+ FOR_EACH(iedge, 0, n) {
+ /* Evaluate the edge equation regarding the projected polyhedron center
+ * and check that its sign is positive */
+ const float dst =
+ Ep[i][iedge][0] * center[j]
+ + Ep[i][iedge][1] * center[k]
+ + Ep[i][iedge][2];
+ ASSERT(dst > 0);
+ }
+ }
+#endif
+ }
+}
+
+
diff --git a/src/suvm_volume.c b/src/suvm_volume.c
@@ -179,10 +179,10 @@ compute_tetrahedral_mesh_normals(struct suvm_volume* vol)
float* n1 = darray_float_data_get(&vol->normals) + (itetra*4 + 1)*3;
float* n2 = darray_float_data_get(&vol->normals) + (itetra*4 + 2)*3;
float* n3 = darray_float_data_get(&vol->normals) + (itetra*4 + 3)*3;
- volume_compute_tetrahedron_normal(vol, itetra, 0, n0);
- volume_compute_tetrahedron_normal(vol, itetra, 1, n1);
- volume_compute_tetrahedron_normal(vol, itetra, 2, n2);
- volume_compute_tetrahedron_normal(vol, itetra, 3, n3);
+ volume_primitive_compute_facet_normal(vol, itetra, 0, n0);
+ volume_primitive_compute_facet_normal(vol, itetra, 1, n1);
+ volume_primitive_compute_facet_normal(vol, itetra, 2, n2);
+ volume_primitive_compute_facet_normal(vol, itetra, 3, n3);
}
exit:
@@ -456,10 +456,10 @@ check_tetrahedra_normals(const struct suvm_volume* vol)
verts[3] = darray_float_cdata_get(&vol->positions) + ids[3]*3/*#coords*/;
/* Fetch the tetrahedron normals */
- volume_get_tetrahedron_normal(vol, itetra, 0, normals[0]);
- volume_get_tetrahedron_normal(vol, itetra, 1, normals[1]);
- volume_get_tetrahedron_normal(vol, itetra, 2, normals[2]);
- volume_get_tetrahedron_normal(vol, itetra, 3, normals[3]);
+ volume_primitive_get_facet_normal(vol, itetra, 0, normals[0]);
+ volume_primitive_get_facet_normal(vol, itetra, 1, normals[1]);
+ volume_primitive_get_facet_normal(vol, itetra, 2, normals[2]);
+ volume_primitive_get_facet_normal(vol, itetra, 3, normals[3]);
/* Compute the position at the center of the tetrahedron */
pt[0] = (verts[0][0] + verts[1][0] + verts[2][0] + verts[3][0]) * 0.25f;
@@ -607,3 +607,23 @@ suvm_volume_get_aabb
upper[2] = volume->upp[2];
return RES_OK;
}
+
+res_T
+suvm_volume_get_primitives_count(const struct suvm_volume* vol, size_t* nprims)
+{
+ if(!vol || !nprims) return RES_BAD_ARG;
+ *nprims = volume_get_primitives_count(vol);
+ return RES_OK;
+}
+
+res_T
+suvm_volume_get_primitive
+ (const struct suvm_volume* vol,
+ const size_t iprim,
+ struct suvm_primitive* prim)
+{
+ if(!vol|| !prim || iprim >= volume_get_primitives_count(vol))
+ return RES_BAD_ARG;
+ volume_primitive_setup(vol, iprim, prim);
+ return RES_OK;
+}
diff --git a/src/suvm_volume.h b/src/suvm_volume.h
@@ -122,7 +122,7 @@ volume_get_vertices_count(const struct suvm_volume* vol)
}
static FINLINE float*
-volume_get_primitive_vertex_position
+volume_primitive_get_vertex_position
(const struct suvm_volume* vol,
const size_t iprim,
const size_t ivert, /* In [0, 3] */
@@ -138,9 +138,9 @@ volume_get_primitive_vertex_position
}
static INLINE float*
-volume_compute_tetrahedron_normal
+volume_primitive_compute_facet_normal
(const struct suvm_volume* vol,
- const size_t itetra,
+ const size_t iprim,
const int ifacet/*In [0,3]*/,
float normal[3])
{
@@ -154,10 +154,11 @@ volume_compute_tetrahedron_normal
/* Facet indices, Ensure CCW ordering */
const int facet[4][3] = {{0,1,2}, {0,3,1}, {1,3,2}, {2,3,0}};
- ASSERT(vol && itetra < volume_get_primitives_count(vol) && normal);
+ ASSERT(vol && iprim < volume_get_primitives_count(vol) && normal);
+ ASSERT(ifacet < 4);
/* Fetch tetrahedron indices */
- ids = darray_u32_cdata_get(&vol->indices) + itetra*4;
+ ids = darray_u32_cdata_get(&vol->indices) + iprim*4;
/* Fetch facet vertices */
v0 = darray_float_cdata_get(&vol->positions) + ids[facet[ifacet][0]]*3;
@@ -174,19 +175,19 @@ volume_compute_tetrahedron_normal
}
static FINLINE float*
-volume_get_tetrahedron_normal
+volume_primitive_get_facet_normal
(const struct suvm_volume* vol,
- const size_t itetra,
+ const size_t iprim,
const int ifacet, /* In [0, 3] */
float normal[3])
{
- ASSERT(vol && itetra < volume_get_primitives_count(vol) && normal);
+ ASSERT(vol && iprim < volume_get_primitives_count(vol) && normal);
if(!darray_float_size_get(&vol->normals)) {
- volume_compute_tetrahedron_normal(vol, itetra, ifacet, normal);
+ volume_primitive_compute_facet_normal(vol, iprim, ifacet, normal);
} else {
const size_t id =
- (itetra*4/*#facets per tetra*/ + (size_t)ifacet)*3/*#coords per normal*/;
+ (iprim*4/*#facets per tetra*/ + (size_t)ifacet)*3/*#coords per normal*/;
ASSERT(id + 3 <= darray_float_size_get(&vol->normals));
normal[0] = darray_float_cdata_get(&vol->normals)[id+0];
normal[1] = darray_float_cdata_get(&vol->normals)[id+1];
@@ -206,7 +207,7 @@ buffer_at(const struct buffer* buf, const size_t i)
}
static INLINE const void*
-volume_get_primitive_data(const struct suvm_volume* vol, const size_t iprim)
+volume_primitive_get_data(const struct suvm_volume* vol, const size_t iprim)
{
ASSERT(vol && vol->has_prim_data && iprim < volume_get_primitives_count(vol));
return buffer_at(&vol->prim_data, iprim);
@@ -225,5 +226,34 @@ volume_primitive_get_vertex_data
return buffer_at(&vol->vert_data, ids[ivert]);
}
+static INLINE void
+volume_primitive_setup
+ (const struct suvm_volume* vol,
+ const size_t iprim,
+ struct suvm_primitive* prim)
+{
+ ASSERT(vol && prim && iprim < volume_get_primitives_count(vol));
+ prim->nvertices = 4;
+ prim->iprim = iprim;
+ prim->volume__ = vol;
+ if(vol->has_prim_data) {
+ prim->data = volume_primitive_get_data(vol, prim->iprim);
+ }
+ if(vol->has_vert_data) {
+ prim->vertex_data[0] = volume_primitive_get_vertex_data(vol, iprim, 0);
+ prim->vertex_data[1] = volume_primitive_get_vertex_data(vol, iprim, 1);
+ prim->vertex_data[2] = volume_primitive_get_vertex_data(vol, iprim, 2);
+ prim->vertex_data[3] = volume_primitive_get_vertex_data(vol, iprim, 3);
+ }
+}
+
+extern LOCAL_SYM void
+tetrahedron_setup
+ (const struct suvm_volume* vol,
+ const float low[3], /* NULL <=> compute on the fly */
+ const float upp[3], /* NULL <=> compute on the fly */
+ const size_t iprim,
+ struct suvm_polyhedron* polyhedron);
+
#endif /* SUVM_VOLUME_H */
diff --git a/src/suvm_volume_at.c b/src/suvm_volume_at.c
@@ -56,7 +56,7 @@ intersect_leaf
v0 = darray_float_cdata_get(&vol->positions) + tetra[0]*3;
/* Compute the distance from pos to the bottom facet */
- volume_get_tetrahedron_normal(vol, itetra, 0, n0);
+ volume_primitive_get_facet_normal(vol, itetra, 0, n0);
d0 = f3_dot(n0, f3_sub(p, pos, v0));
if(d0 < 0) return 0;
@@ -65,13 +65,13 @@ intersect_leaf
f3_sub(p, pos, v3);
/* Compute the distance from pos to the 'side' facets */
- volume_get_tetrahedron_normal(vol, itetra, 1, n1);
+ volume_primitive_get_facet_normal(vol, itetra, 1, n1);
d1 = f3_dot(n1, p);
if(d1 < 0) return 0;
- volume_get_tetrahedron_normal(vol, itetra, 2, n2);
+ volume_primitive_get_facet_normal(vol, itetra, 2, n2);
d2 = f3_dot(n2, p);
if(d2 < 0) return 0;
- volume_get_tetrahedron_normal(vol, itetra, 3, n3);
+ volume_primitive_get_facet_normal(vol, itetra, 3, n3);
d3 = f3_dot(n3, p);
if(d3 < 0) return 0;
@@ -184,17 +184,7 @@ suvm_volume_at
/* pos lies into a tetrahedron */
if(iprim != SIZE_MAX) {
- prim->iprim = iprim;
- prim->nvertices = 4;
- if(vol->has_prim_data) {
- prim->data = volume_get_primitive_data(vol, iprim);
- }
- if(vol->has_vert_data) {
- prim->vertex_data[0] = volume_primitive_get_vertex_data(vol, iprim, 0);
- prim->vertex_data[1] = volume_primitive_get_vertex_data(vol, iprim, 1);
- prim->vertex_data[2] = volume_primitive_get_vertex_data(vol, iprim, 2);
- prim->vertex_data[3] = volume_primitive_get_vertex_data(vol, iprim, 3);
- }
+ volume_primitive_setup(vol, iprim, prim);
}
exit:
diff --git a/src/suvm_volume_intersect_aabb.c b/src/suvm_volume_intersect_aabb.c
@@ -20,301 +20,9 @@
#include <rsys/float2.h>
#include <rsys/float3.h>
-enum { X, Y, Z, AXES_COUNT }; /* Syntactic sugar */
-
-struct tetrahedron {
- float v[4][AXES_COUNT]; /* Tetrahedron vertices */
- float N[4][AXES_COUNT]; /* Tetrahedron normals */
- float D[4]; /* Slope parameter of the plane equation Ax + By + Cz + D = 0 */
-
- /* 2D equation of the silhouette edges projected along the X, Y and Z axis of
- * the form Ax + By + C = 0, with (A, B) the normals of the edge */
- float Ep[AXES_COUNT][4/*#edges max*/][3];
-
- /* Number of silhouette edges in the ZY, XZ and XY planes (3 or 4) */
- int nEp[AXES_COUNT];
-
- /* Tetrahedron axis aligned bounding box */
- float lower[AXES_COUNT];
- float upper[AXES_COUNT];
-};
-
/*******************************************************************************
* Helper functions
******************************************************************************/
-static INLINE void
-tetrahedron_setup
- (const struct suvm_volume* vol,
- const float low[3],
- const float upp[3],
- const size_t itetra,
- struct tetrahedron* tetra)
-{
-#ifndef NDEBUG
- float center[AXES_COUNT]; /* Center of the tetrahedron */
-#endif
- float e[6][AXES_COUNT];
- float (*v)[AXES_COUNT];
- float (*N)[AXES_COUNT];
- float (*Ep)[4][3];
- int i;
- ASSERT(vol && low && upp && tetra);
-
- /* Fetch tetrahedron vertices */
- volume_get_primitive_vertex_position(vol, itetra, 0, tetra->v[0]);
- volume_get_primitive_vertex_position(vol, itetra, 1, tetra->v[1]);
- volume_get_primitive_vertex_position(vol, itetra, 2, tetra->v[2]);
- volume_get_primitive_vertex_position(vol, itetra, 4, tetra->v[3]);
- v = tetra->v;
-
-#ifndef NDEBUG
- /* Compute the center of the tetrahedron */
- center[X] = (v[0][X] + v[1][X] + v[2][X] + v[3][X]) * 0.25f;
- center[Y] = (v[0][Y] + v[1][Y] + v[2][Y] + v[3][Y]) * 0.25f;
- center[Z] = (v[0][Z] + v[1][Z] + v[2][Z] + v[3][Z]) * 0.25f;
-#endif
-
- /* Check argument consistency */
- ASSERT(MMIN(MMIN(v[0][X], v[1][X]), MMIN(v[2][X], v[3][X])) == low[X]);
- ASSERT(MMIN(MMIN(v[0][Y], v[1][Y]), MMIN(v[2][Y], v[3][Y])) == low[Y]);
- ASSERT(MMIN(MMIN(v[0][Z], v[1][Z]), MMIN(v[2][Z], v[3][Z])) == low[Z]);
- ASSERT(MMAX(MMAX(v[0][X], v[1][X]), MMAX(v[2][X], v[3][X])) == upp[X]);
- ASSERT(MMAX(MMAX(v[0][Y], v[1][Y]), MMAX(v[2][Y], v[3][Y])) == upp[Y]);
- ASSERT(MMAX(MMAX(v[0][Z], v[1][Z]), MMAX(v[2][Z], v[3][Z])) == upp[Z]);
-
- /* Setup tetrahedron AABB */
- f3_set(tetra->lower, low);
- f3_set(tetra->upper, upp);
-
- /* Fetch tetrahedron normals */
- volume_get_tetrahedron_normal(vol, itetra, 0, tetra->N[0]);
- volume_get_tetrahedron_normal(vol, itetra, 1, tetra->N[1]);
- volume_get_tetrahedron_normal(vol, itetra, 2, tetra->N[2]);
- volume_get_tetrahedron_normal(vol, itetra, 3, tetra->N[3]);
- N = tetra->N;
-
- /* Compute the slope of the planes */
- tetra->D[0] = -f3_dot(tetra->N[0], v[0]);
- tetra->D[1] = -f3_dot(tetra->N[1], v[1]);
- tetra->D[2] = -f3_dot(tetra->N[2], v[2]);
- tetra->D[3] = -f3_dot(tetra->N[3], v[3]);
-
- /* Compute tetrahedron edges */
- f3_sub(e[0], v[1], v[0]);
- f3_sub(e[1], v[2], v[1]);
- f3_sub(e[2], v[0], v[2]);
- f3_sub(e[3], v[0], v[3]);
- f3_sub(e[4], v[1], v[3]);
- f3_sub(e[5], v[2], v[3]);
-
- /* Detect the silhouette edges of the tetrahedron once projected in the YZ,
- * ZX and XY planes, and compute their 2D equation in the projected plane. The
- * variable 'i' defines the projection axis which is successively X (YZ
- * plane), Y (ZX plane) and Z (XY plane). The axes of the corresponding 2D
- * repairs are defined by the 'j' and 'k' variables. */
- Ep = tetra->Ep;
- FOR_EACH(i, 0, AXES_COUNT) {
- /* On 1st iteration 'j' and 'k' define the YZ plane
- * On 2nd iteration 'j' and 'k' define the ZX plane
- * On 3rd ietration 'j' and 'k' define the XY plane */
- const int j = (i + 1) % AXES_COUNT;
- const int k = (j + 1) % AXES_COUNT;
-
- /* Register the number of detected silhouette edges */
- int n = 0;
-
- /* To detect the silhouette edges, check the sign of the normals of two
- * adjacent facets for the coordinate of the projection axis. If the signs
- * are the same, the facets look at the same direction regarding the
- * projection axis. The other case means that one facet points toward the
- * projection axis whole the other one points to the opposite direction:
- * this is the definition of a silhouette edge.
- *
- * Once detected, we compute the equation of the silhouette edge:
- *
- * Ax + By + C = 0
- *
- * with A and B the coordinates of the edge normals and C the slope of the
- * edge. Computing the normal is actually as simple as swizzling the
- * coordinates of the edges projected in 2D and sign tuning. Let the
- * projection axis X (i==X) and thus the YZ projection plane (j==Y &&
- * k==Z). Assuming that the first edge 'e[0]' is a silhouette edge, ie the
- * edge between the facet0 and the facet1 (refers to the suvm_volume.h for
- * the geometric layout of a tetrahedron) once projected in the YZ plane
- * the edge is {0, e[0][Y], e[0][Z]}. Its normal N can is defined as the
- * cross product between this projected edge ands the projection axis:
- *
- * | 0 | | 1 | | 0 |
- * N = | e[0][Y] | X | 0 | = | e[0][Z] |
- * | e[0][Z] | | 0 | |-e[0][Y] |
- *
- * In the projection plane, the _unormalized_ normal of the {e[0][Y],
- * e[0][Z]} edge is thus {e[0][Z], -e[0][Y]}. More generally, the normal of
- * an edge 'I' {e[I][j], e[I][k]} projected along the 'i' axis in the 'jk'
- * plane is {e[I][k],-e[I]j]}. Anyway, one has to revert the normal
- * orientation to ensure a specific convention wrt the tetrahedron
- * footprint. In our case we ensure that the normal points toward the
- * footprint. Finalle the edge slope 'C' can be computed as:
- * | A |
- * C =- | B | . | Pj, Pk |
- *
- * with {Pj, Pk} a point belonging to the edge as for instance one of its
- * vertices. */
- if(signf(N[0][i]) != signf(N[1][i])) { /* The edge 0 is silhouette */
- f2_normalize(Ep[i][n], f2(Ep[i][n], e[0][k], -e[0][j]));
- if(N[0][i] > 0) f2_minus(Ep[i][n], Ep[i][n]);
- Ep[i][n][2] = -f2_dot(Ep[i][n], v[0]);
- ++n;
- }
- if(signf(N[0][i]) != signf(N[2][i])) { /* The edge 1 is silhouette */
- f2_normalize(Ep[i][n], f2(Ep[i][n], e[1][k], -e[1][j]));
- if(N[0][i] > 0) f2_minus(Ep[i][n], Ep[i][n]);
- Ep[i][n][2] = -f2_dot(Ep[i][n], v[1]);
- ++n;
- }
- if(signf(N[0][i]) != signf(N[3][i])) { /* The edge 2 is silhouette */
- f2_normalize(Ep[i][n], f2(Ep[i][n], e[2][k], -e[2][j]));
- if(N[0][i] > 0) f2_minus(Ep[i][n], Ep[i][n]);
- Ep[i][n][2] = -f2_dot(Ep[i][n], v[2]);
- ++n;
- }
- if(signf(N[1][i]) != signf(N[3][i])) { /* The edge 3 is silhouette */
- f2_normalize(Ep[i][n], f2(Ep[i][n], e[3][k], -e[3][j]));
- if(N[1][i] > 0) f2_minus(Ep[i][n], Ep[i][n]);
- Ep[i][n][2] = -f2_dot(Ep[i][n], v[0]);
- ++n;
- }
- if(signf(N[1][i]) != signf(N[2][i])) { /* The edge 4 is silhouette */
- f2_normalize(Ep[X][n], f2(Ep[i][n], e[4][k], -e[4][j]));
- if(N[2][i] > 0) f2_minus(Ep[i][n], Ep[i][n]);
- Ep[i][n][2] = -f2_dot(Ep[i][n], v[1]);
- ++n;
- }
- if(signf(N[2][i]) != signf(N[3][i])) { /* The edge 5 is silhouette */
- f2_normalize(Ep[i][n], f2(Ep[i][n], e[5][k], -e[5][j]));
- if(N[3][i] > 0) f2_minus(Ep[i][n], Ep[i][n]);
- Ep[i][n][2] = -f2_dot(Ep[i][n], v[2]);
- ++n;
- }
- /* A tetrahedron can have 3 or 4 silhouette edges once projected in 2D */
- ASSERT(n == 3 || n == 4);
- tetra->nEp[i] = n; /* Register the #silouhette edges for this project axis */
-
-#ifndef NDEBUG
- /* Assert that the normals points toward the expected half space */
- {
- int iedge;
- FOR_EACH(iedge, 0, n) {
- /* Evaluate the edge equation regarding the projected polyhedron center
- * and check that its sign is positive */
- const float dst =
- Ep[i][iedge][0] * center[j]
- + Ep[i][iedge][1] * center[k]
- + Ep[i][iedge][2];
- ASSERT(dst > 0);
- }
- }
-#endif
- }
-}
-
-/* Define if an axis aligned bounding box and a tetrahedron are intersecting.
- * Its implementation follows the algorithms proposed by N. Greene in
- * "Detecting intersection of a rectangular solid and a convex polyhedron"
- * (Graphics Gems IV, 1994, p74--82) */
-static INLINE int
-aabb_intersect_tetrahedron
- (const float low[3],
- const float upp[3],
- const struct tetrahedron* tetra)
-{
- float intersect_aabb_low[3];
- float intersect_aabb_upp[3];
- int i;
- int nplanes_including_aabb;
- ASSERT(low && upp && tetra);
- ASSERT(low[X] < upp[X]);
- ASSERT(low[Y] < upp[Y]);
- ASSERT(low[Z] < upp[Z]);
-
- /* Check the intersection between the AABB and the tetra AABB */
- FOR_EACH(i, 0, AXES_COUNT) {
- intersect_aabb_low[i] = MMAX(low[i], tetra->lower[i]);
- intersect_aabb_upp[i] = MMIN(upp[i], tetra->upper[i]);
- if(intersect_aabb_low[i] > intersect_aabb_upp[i]) /* Do not intersect */
- return 0;
- }
-
- /* Check if the tetrahedron is included into the aabb */
- if(intersect_aabb_low[X] == low[X]
- && intersect_aabb_low[Y] == low[Y]
- && intersect_aabb_low[Z] == low[Z]) {
- return 1;
- }
-
- /* #planes that totally includes the aabb on its positive side */
- nplanes_including_aabb = 0;
-
- /* Check the tetrahedron planes against the aabb */
- FOR_EACH(i, 0, 4) {
- float p[3], n[3];
-
- /* Define the aabb vertices that is the farthest in the positive/negative
- * direction of the face normal. */
- f3_set(p, upp);
- f3_set(n, low);
- if(tetra->N[i][X] < 0) SWAP(float, p[X], n[X]);
- if(tetra->N[i][Y] < 0) SWAP(float, p[Y], n[Y]);
- if(tetra->N[i][Z] < 0) SWAP(float, p[Z], n[Z]);
-
- /* Check that the box is totally outside the plane, To do this, check that
- * 'p', aka the farthest aabb vertex in the positive direction of the plane
- * normal, is not in the positive half space. In this case, the aabb is
- * entirely outside the tetrahedron */
- if((f3_dot(tetra->N[i], p) + tetra->D[i]) < 0) return 0;
-
- /* Check if the box is totally inside the given plane. To do this, check
- * that 'n', aka the farthest aabb vertex in the negative direction of the
- * plane normal, is inside the positive half space */
- if((f3_dot(tetra->N[i], n) + tetra->D[i]) > 0) {
- /* Register this plane as a plane including the aabb */
- nplanes_including_aabb += 1;
- }
- }
-
- /* Check if the aabb is entirely included into the tetrahedron */
- if(nplanes_including_aabb == 4) return 1;
-
- /* For each silhouette edge in each projection plane, check if it totally
- * exclude the projected aabb */
- FOR_EACH(i, 0, AXES_COUNT) {
- /* Compute the remaining axes.
- * On 1st iteration 'j' and 'k' define the YZ plane
- * On 2nd iteration 'j' and 'k' define the ZX plane
- * On 3rd ietration 'j' and 'k' define the XY plane */
- const int j = (i + 1) % AXES_COUNT;
- const int k = (j + 1) % AXES_COUNT;
-
- int iedge;
- FOR_EACH(iedge, 0, tetra->nEp[i]) {
- float p[2];
-
- /* Define the aabb vertices that is the farthest in the positive/negative
- * direction of the normal of the silhouette edge. */
- p[0] = tetra->Ep[i][iedge][0] > 0 ? upp[j] : low[j];
- p[1] = tetra->Ep[i][iedge][1] > 0 ? upp[k] : low[k];
-
- /* Check that the projected aabb along the 'i'th axis is totally outside
- * the current silhouette edge. That means that the aabb is entirely
- * outside the tetrahedron and thus that they do not intersect */
- if((f2_dot(tetra->Ep[i][iedge], p) + tetra->Ep[i][iedge][2]) < 0) return 0;
- }
- }
-
- /* The tetrahedron and the aabb are partially intersecting */
- return 1;
-}
-
static INLINE int
aabb_intersect_aabb
(const float low0[3],
@@ -392,17 +100,7 @@ suvm_volume_intersect_aabb
if(aabb_intersect_aabb(lowf, uppf, leaf->low, leaf->upp)) {
struct suvm_primitive prim = SUVM_PRIMITIVE_NULL;
/* Setup the intersected primitive */
- prim.nvertices = 4;
- prim.iprim = leaf->prim_id;
- if(vol->has_prim_data) {
- prim.data = volume_get_primitive_data(vol, prim.iprim);
- }
- if(vol->has_vert_data) {
- prim.vertex_data[0] = volume_primitive_get_vertex_data(vol, prim.iprim, 0);
- prim.vertex_data[1] = volume_primitive_get_vertex_data(vol, prim.iprim, 1);
- prim.vertex_data[2] = volume_primitive_get_vertex_data(vol, prim.iprim, 2);
- prim.vertex_data[3] = volume_primitive_get_vertex_data(vol, prim.iprim, 3);
- }
+ volume_primitive_setup(vol, leaf->prim_id, &prim);
/* Invoke the user callback on the intersected primitive */
clbk(&prim, low, upp, context);
}
diff --git a/src/test_suvm_volume.c b/src/test_suvm_volume.c
@@ -225,7 +225,10 @@ test_tetrahedral_mesh_creation(struct suvm_device* dev)
}
static void
-check_prim(const struct suvm_primitive* prim, struct mesh* msh)
+check_prim
+ (const struct suvm_primitive* prim,
+ const struct suvm_volume* vol,
+ struct mesh* msh)
{
size_t ids[4];
double verts[4][3];
@@ -234,6 +237,8 @@ check_prim(const struct suvm_primitive* prim, struct mesh* msh)
/* Check primitive data */
CHK(prim != NULL);
+ CHK(vol != NULL);
+ CHK(prim->volume__ == vol);
CHK(!SUVM_PRIMITIVE_NONE(prim));
CHK(prim->data != NULL);
CHK(prim->data != NULL);
@@ -271,6 +276,7 @@ check_prim(const struct suvm_primitive* prim, struct mesh* msh)
static void
check_prim_at
(const struct suvm_primitive* prim,
+ const struct suvm_volume* vol,
const double pos[3],
const double bcoords[4],
struct mesh* msh)
@@ -283,7 +289,7 @@ check_prim_at
double p[3];
const double* vert_data[4];
- check_prim(prim, msh);
+ check_prim(prim, vol, msh);
/* Fetch tetrahedron vertices */
get_indices(prim->iprim, ids, msh);
@@ -345,6 +351,7 @@ check_prim_at
static void
check_prims_intersect_aabb
(struct darray_prim* primitives,
+ const struct suvm_volume* vol,
const double low[3],
const double upp[3],
struct mesh* msh)
@@ -363,7 +370,7 @@ check_prims_intersect_aabb
nprims = darray_prim_size_get(primitives);
FOR_EACH(iprim, 0, nprims) {
- check_prim(prims + iprim, msh);
+ check_prim(prims + iprim, vol, msh);
/* Fetch tetrahedron vertices */
get_indices(prims[iprim].iprim, ids, msh);
@@ -485,7 +492,7 @@ test_volume_at_cube(struct suvm_device* dev)
CHK(suvm_volume_at(vol, pos, &prim, bcoords) == RES_OK);
CHK(!SUVM_PRIMITIVE_NONE(&prim));
CHK(prim.iprim == 2);
- check_prim_at(&prim, pos, bcoords, &msh);
+ check_prim_at(&prim, vol, pos, bcoords, &msh);
CHK(suvm_volume_get_aabb(NULL, low, upp) == RES_BAD_ARG);
CHK(suvm_volume_get_aabb(vol, NULL, upp) == RES_BAD_ARG);
@@ -499,7 +506,7 @@ test_volume_at_cube(struct suvm_device* dev)
samp[1] = low[1] + rand_canonic()*(upp[1] - low[1]);
samp[2] = low[2] + rand_canonic()*(upp[2] - low[2]);
CHK(suvm_volume_at(vol, samp, &prim, bcoords) == RES_OK);
- check_prim_at(&prim, samp, bcoords, &msh);
+ check_prim_at(&prim, vol, samp, bcoords, &msh);
}
pos[0] = upp[0] + 0.1;
@@ -523,7 +530,7 @@ test_volume_at_cube(struct suvm_device* dev)
(vol, low, upp, prim_intersect_aabb, &prims) == RES_OK);
CHK(darray_prim_size_get(&prims) == box_ntetras);
- check_prims_intersect_aabb(&prims, low, upp, &msh);
+ check_prims_intersect_aabb(&prims, vol, low, upp, &msh);
CHK(suvm_volume_intersect_aabb
(vol, upp, low, prim_intersect_aabb, &prims) == RES_BAD_ARG);
@@ -548,7 +555,7 @@ test_volume_at_cube(struct suvm_device* dev)
CHK(suvm_volume_intersect_aabb
(vol, low, upp, prim_intersect_aabb, &prims) == RES_OK);
CHK(darray_prim_size_get(&prims) == 10);
- check_prims_intersect_aabb(&prims, low, upp, &msh);
+ check_prims_intersect_aabb(&prims, vol, low, upp, &msh);
CHK(suvm_volume_ref_put(vol) == RES_OK);
@@ -626,7 +633,7 @@ test_volume_at_ball(struct suvm_device* dev)
if(!SUVM_PRIMITIVE_NONE(&prim)) {
check_at += 1;
- check_prim_at(&prim, samp, bcoords, &msh);
+ check_prim_at(&prim, vol, samp, bcoords, &msh);
}
}
CHK(check_at != 0);
@@ -657,7 +664,7 @@ test_volume_at_ball(struct suvm_device* dev)
(vol, aabb_low, aabb_upp, prim_intersect_aabb, &prims) == RES_OK);
check_intersect_aabb = darray_prim_size_get(&prims) != 0;
- check_prims_intersect_aabb(&prims, aabb_low, aabb_upp, &msh);
+ check_prims_intersect_aabb(&prims, vol, aabb_low, aabb_upp, &msh);
darray_prim_clear(&prims);
}
CHK(check_intersect_aabb != 0);