commit 8436345d9432589c2dbcf9b064448d0c24dcaaf3
parent 930e5955e962b945526176511d921e15aa20316e
Author: Vincent Forest <vincent.forest@meso-star.com>
Date: Wed, 9 Jan 2019 15:39:07 +0100
Replace Embree2 by Embree3 in Ray-Tracing backend
No test is performed yet. Only compilation is OK on GNU/Linux.
Diffstat:
10 files changed, 538 insertions(+), 361 deletions(-)
diff --git a/cmake/CMakeLists.txt b/cmake/CMakeLists.txt
@@ -30,7 +30,7 @@
# knowledge of the CeCILL license and that you accept its terms.
cmake_minimum_required(VERSION 2.8)
-project(star-3d C CXX)
+project(star-3d C)
enable_testing()
set(S3D_SOURCE_DIR ${PROJECT_SOURCE_DIR}/../src)
@@ -39,10 +39,7 @@ option(NO_TEST "Disable the test" OFF)
################################################################################
# Check dependencies
################################################################################
-get_filename_component(_current_source_dir ${CMAKE_CURRENT_LIST_FILE} PATH)
-set(Embree_DIR ${_current_source_dir}/)
-
-find_package(Embree REQUIRED)
+find_package(Embree 3 REQUIRED)
find_package(RCMake 0.2.2 REQUIRED)
find_package(RSys 0.6 REQUIRED)
@@ -50,10 +47,14 @@ set(CMAKE_MODULE_PATH ${CMAKE_MODULE_PATH} ${RCMAKE_SOURCE_DIR})
include(rcmake)
include(rcmake_runtime)
-include_directories(${Embree_INCLUDE_DIR} ${RSys_INCLUDE_DIR})
+include_directories(${EMBREE_INCLUDE_DIRS} ${RSys_INCLUDE_DIR})
rcmake_append_runtime_dirs(_runtime_dirs RSys Embree)
+if(CMAKE_COMPILER_IS_GNUCC)
+ set(CMAKE_C_FLAGS "${CMAKE_C_FLAGS} -std=c99")
+endif()
+
################################################################################
# Configure and define targets
################################################################################
@@ -92,27 +93,25 @@ rcmake_prepend_path(S3D_FILES_SRC ${S3D_SOURCE_DIR})
rcmake_prepend_path(S3D_FILES_INC ${S3D_SOURCE_DIR})
rcmake_prepend_path(S3D_FILES_INC_API ${S3D_SOURCE_DIR})
rcmake_prepend_path(S3D_FILES_DOC ${PROJECT_SOURCE_DIR}/../)
-set_source_files_properties(${S3D_FILES_SRC} PROPERTIES LANGUAGE CXX)
add_library(s3d SHARED ${S3D_FILES_SRC} ${S3D_FILES_INC} ${S3D_FILES_INC_API})
-target_link_libraries(s3d RSys Embree)
+target_link_libraries(s3d RSys ${EMBREE_LIBRARIES})
-if(CMAKE_COMPILER_IS_GNUCXX)
+if(CMAKE_COMPILER_IS_GNUCC)
target_link_libraries(s3d m)
endif()
set_target_properties(s3d PROPERTIES
DEFINE_SYMBOL S3D_SHARED_BUILD
- LINKER_LANGUAGE CXX
VERSION ${VERSION}
SOVERSION ${VERSION_MAJOR})
-if(CMAKE_COMPILER_IS_GNUCXX)
- # Shut up the use of variadic macros by Embree and the use of long long
- # constants by RSys
- set_target_properties(s3d PROPERTIES
- COMPILE_FLAGS "-Wno-variadic-macros -Wno-long-long")
-endif()
+#if(CMAKE_COMPILER_IS_GNUCXX) FIXME remove this ?
+# # Shut up the use of variadic macros by Embree and the use of long long
+# # constants by RSys
+# set_target_properties(s3d PROPERTIES
+# COMPILE_FLAGS "-Wno-variadic-macros -Wno-long-long")
+#endif()
rcmake_setup_devel(s3d Star3D ${VERSION} star/s3d_version.h)
diff --git a/cmake/EmbreeConfig.cmake b/cmake/EmbreeConfig.cmake
@@ -1,60 +0,0 @@
-# Copyright (C) |Meso|Star> 2015 (contact@meso-star.com)
-#
-# This software is a computer program whose purpose is to generate files
-# used to build the Star-3D library.
-#
-# This software is governed by the CeCILL license under French law and
-# abiding by the rules of distribution of free software. You can use,
-# modify and/or redistribute the software under the terms of the CeCILL
-# license as circulated by CEA, CNRS and INRIA at the following URL
-# "http://www.cecill.info".
-#
-# As a counterpart to the access to the source code and rights to copy,
-# modify and redistribute granted by the license, users are provided only
-# with a limited warranty and the software's author, the holder of the
-# economic rights, and the successive licensors have only limited
-# liability.
-#
-# In this respect, the user's attention is drawn to the risks associated
-# with loading, using, modifying and/or developing or reproducing the
-# software by the user in light of its specific status of free software,
-# that may mean that it is complicated to manipulate, and that also
-# therefore means that it is reserved for developers and experienced
-# professionals having in-depth computer knowledge. Users are therefore
-# encouraged to load and test the software's suitability as regards their
-# requirements in conditions enabling the security of their systems and/or
-# data to be ensured and, more generally, to use and operate it in the
-# same conditions as regards security.
-#
-# The fact that you are presently reading this means that you have had
-# knowledge of the CeCILL license and that you accept its terms.
-
-cmake_minimum_required(VERSION 2.6)
-
-# Try to find the Embree devel. Once done this will define:
-# - Embree_FOUND: system has Embree
-# - Embree_INCLUDE_DIR: the include directory
-# - Embree: Link this to use Embree
-
-find_path(Embree_INCLUDE_DIR embree2/rtcore.h)
-unset(Embree_LIBRARY CACHE)
-unset(Embree_LIBRARY_DEBUG CACHE)
-unset(Embree_LIBRARY_RELWITHDEBINFO CACHE)
-unset(Embree_LIBRARY_MINSIZEREL CACHE)
-find_library(Embree_LIBRARY embree DOC "Path to the library embree.")
-
-# Create the imported library target
-if(CMAKE_HOST_WIN32)
- set(_property IMPORTED_IMPLIB)
-else(CMAKE_HOST_WIN32)
- set(_property IMPORTED_LOCATION)
-endif(CMAKE_HOST_WIN32)
-add_library(Embree SHARED IMPORTED)
-set_target_properties(Embree PROPERTIES ${_property} ${Embree_LIBRARY})
-
-# Check the package
-include(FindPackageHandleStandardArgs)
-FIND_PACKAGE_HANDLE_STANDARD_ARGS(Embree DEFAULT_MSG
- Embree_INCLUDE_DIR
- Embree_LIBRARY)
-
diff --git a/src/s3d_backend.h b/src/s3d_backend.h
@@ -37,15 +37,14 @@
#ifdef COMPILER_GCC
#pragma GCC diagnostic push
- #pragma GCC diagnostic ignored "-Wpedantic"
- #pragma GCC diagnostic ignored "-Wsign-compare"
+/* #pragma GCC diagnostic ignored "-Wpedantic"
+ #pragma GCC diagnostic ignored "-Wsign-compare" */
#elif defined(COMPILER_CL)
#pragma warning(push)
- #pragma warning(disable:4324) /* Structure was padded due to alignment */
+/* #pragma warning(disable:4324)*/ /* Structure was padded due to alignment */
#endif
-#include <embree2/rtcore.h>
-#include <embree2/rtcore_ray.h>
+#include <embree3/rtcore.h>
#ifdef COMPILER_GCC
#pragma GCC diagnostic pop
diff --git a/src/s3d_buffer.h b/src/s3d_buffer.h
@@ -53,7 +53,7 @@
#define BUFFER_DARRAY_FUNC__(Func) CONCAT(CONCAT(BUFFER_DARRAY, _), Func)
struct BUFFER_NAME {
- BUFFER_DARRAY data;
+ struct BUFFER_DARRAY data;
struct mem_allocator* allocator;
ref_T ref;
};
diff --git a/src/s3d_c.h b/src/s3d_c.h
@@ -44,20 +44,39 @@ struct hit_filter {
void* data;
};
+
static FINLINE res_T
rtc_error_to_res_T(const enum RTCError err)
{
switch(err) {
- case RTC_NO_ERROR: return RES_OK;
- case RTC_UNKNOWN_ERROR: return RES_UNKNOWN_ERR;
- case RTC_INVALID_ARGUMENT: return RES_BAD_ARG;
- case RTC_INVALID_OPERATION: return RES_BAD_ARG;
- case RTC_OUT_OF_MEMORY: return RES_MEM_ERR;
- case RTC_UNSUPPORTED_CPU: return RES_BAD_ARG;
+ case RTC_ERROR_NONE: return RES_OK;
+ case RTC_ERROR_UNKNOWN: return RES_UNKNOWN_ERR;
+ case RTC_ERROR_INVALID_ARGUMENT: return RES_BAD_ARG;
+ case RTC_ERROR_INVALID_OPERATION: return RES_BAD_ARG;
+ case RTC_ERROR_OUT_OF_MEMORY: return RES_MEM_ERR;
+ case RTC_ERROR_UNSUPPORTED_CPU: return RES_BAD_ARG;
+ case RTC_ERROR_CANCELLED: return RES_UNKNOWN_ERR;
default: FATAL("Unreachable code\n"); break;
}
}
+static INLINE const char*
+rtc_error_string(const enum RTCError err)
+{
+ const char* str = NULL;
+ switch(err) {
+ case RTC_ERROR_NONE: str = "No error"; break;
+ case RTC_ERROR_UNKNOWN: str = "Unknown error"; break;
+ case RTC_ERROR_INVALID_ARGUMENT: str = "Invalid argument"; break;
+ case RTC_ERROR_INVALID_OPERATION: str = "Invalid operation"; break;
+ case RTC_ERROR_OUT_OF_MEMORY: str = "Out of memory"; break;
+ case RTC_ERROR_UNSUPPORTED_CPU: str = "Unsupported CPU"; break;
+ case RTC_ERROR_CANCELLED: str = "Cancelled operation"; break;
+ default: FATAL("Unreachable code\n"); break;
+ }
+ return str;
+}
+
static INLINE unsigned
s3d_type_get_dimension(const enum s3d_type type)
{
@@ -70,5 +89,100 @@ s3d_type_get_dimension(const enum s3d_type type)
}
}
+#define RAYN_GRAB(RayN, N, i, Type, Attr) \
+ (((Type*)((char*)(RayN)+(offsetof(struct RTCRay, Attr)*N)))[i])
+#define HITN_GRAB(HitN, N, i, Type, Attr) \
+ (((Type*)((char*)(HitN)+(offsetof(struct RTCHit, Attr)*N)))[i])
+#define RAYHITN_GET_RAYN(RayHitN, N) \
+ ((struct RTCRayN*)((char*)RayHitN+offsetof(struct RTCRayHit, ray)*N))
+#define RAYHITN_GET_HITN(RayHitN, N) \
+ ((struct RTCHitN*)((char*)RayHitN+offsetof(struct RTCRayHit, hit)*N))
+
+static FINLINE void
+rtc_rayN_get_ray
+ (const struct RTCRayN* rayN, /* SoA layout */
+ const size_t N, /* SoA width */
+ const size_t i, /* Id of the ray */
+ struct RTCRay* ray)
+{
+ ASSERT(rayN && ray && i < N);
+ ray->org_x = RAYN_GRAB(rayN, N, i, float, org_x);
+ ray->org_y = RAYN_GRAB(rayN, N, i, float, org_y);
+ ray->org_z = RAYN_GRAB(rayN, N, i, float, org_z);
+ ray->tnear = RAYN_GRAB(rayN, N, i, float, tnear);
+ ray->dir_x = RAYN_GRAB(rayN, N, i, float, dir_x);
+ ray->dir_y = RAYN_GRAB(rayN, N, i, float, dir_y);
+ ray->dir_z = RAYN_GRAB(rayN, N, i, float, dir_z);
+ ray->time = RAYN_GRAB(rayN, N, i, float, time);
+ ray->tfar = RAYN_GRAB(rayN, N, i, float, tfar);
+ ray->mask = RAYN_GRAB(rayN, N, i, unsigned, mask);
+ ray->id = RAYN_GRAB(rayN, N, i, unsigned, id);
+ ray->flags = RAYN_GRAB(rayN, N, i, unsigned, flags);
+}
+
+static FINLINE void
+rtc_hitN_get_hit
+ (const struct RTCHitN* hitN,
+ const size_t N,
+ const size_t i,
+ struct RTCHit* hit)
+{
+ size_t id;
+ ASSERT(hitN && hit && i < N);
+ hit->Ng_x = HITN_GRAB(hitN, N, i, float, Ng_x);
+ hit->Ng_y = HITN_GRAB(hitN, N, i, float, Ng_y);
+ hit->Ng_z = HITN_GRAB(hitN, N, i, float, Ng_z);
+ hit->u = HITN_GRAB(hitN, N, i, float, u);
+ hit->v = HITN_GRAB(hitN, N, i, float, v);
+ hit->primID = HITN_GRAB(hitN, N, i, unsigned, primID);
+ hit->geomID = HITN_GRAB(hitN, N, i, unsigned, geomID);
+ FOR_EACH(id, 0, RTC_MAX_INSTANCE_LEVEL_COUNT) {
+ hit->instID[id] = HITN_GRAB(hitN, N, i, unsigned, instID[id]);
+ }
+}
+
+static FINLINE void
+rtc_rayN_set_ray
+ (struct RTCRayN* rayN,
+ const size_t N,
+ const size_t i,
+ const struct RTCRay* ray)
+{
+ ASSERT(rayN && ray && i < N);
+ RAYN_GRAB(rayN, N, i, float, org_x) = ray->org_x;
+ RAYN_GRAB(rayN, N, i, float, org_y) = ray->org_y;
+ RAYN_GRAB(rayN, N, i, float, org_z) = ray->org_z;
+ RAYN_GRAB(rayN, N, i, float, tnear) = ray->tnear;
+ RAYN_GRAB(rayN, N, i, float, dir_x) = ray->dir_x;
+ RAYN_GRAB(rayN, N, i, float, dir_y) = ray->dir_y;
+ RAYN_GRAB(rayN, N, i, float, dir_z) = ray->dir_z;
+ RAYN_GRAB(rayN, N, i, float, time) = ray->time;
+ RAYN_GRAB(rayN, N, i, float, tfar) = ray->tfar;
+ RAYN_GRAB(rayN, N, i, unsigned, mask) = ray->mask;
+ RAYN_GRAB(rayN, N, i, unsigned, id) = ray->id;
+ RAYN_GRAB(rayN, N, i, unsigned, flags) = ray->flags;
+}
+
+static FINLINE void
+rtc_hitN_set_hit
+ (const struct RTCHitN* hitN,
+ const size_t N,
+ const size_t i,
+ struct RTCHit* hit)
+{
+ size_t id;
+ ASSERT(hitN && hit && i < N);
+ HITN_GRAB(hitN, N, i, float, Ng_x) = hit->Ng_x;
+ HITN_GRAB(hitN, N, i, float, Ng_y) = hit->Ng_y;
+ HITN_GRAB(hitN, N, i, float, Ng_z) = hit->Ng_z;
+ HITN_GRAB(hitN, N, i, float, u) = hit->u;
+ HITN_GRAB(hitN, N, i, float, v) = hit->v;
+ HITN_GRAB(hitN, N, i, unsigned, primID) = hit->primID;
+ HITN_GRAB(hitN, N, i, unsigned, geomID) = hit->geomID;
+ FOR_EACH(id, 0, RTC_MAX_INSTANCE_LEVEL_COUNT) {
+ HITN_GRAB(hitN, N, i, unsigned, instID[id]) = hit->instID[id];
+ }
+}
+
#endif /* S3D_C_H */
diff --git a/src/s3d_device.c b/src/s3d_device.c
@@ -31,6 +31,8 @@
* knowledge of the CeCILL license and that you accept its terms. */
#include "s3d.h"
+#include "s3d_c.h"
+
#include "s3d_backend.h"
#include "s3d_device_c.h"
@@ -40,6 +42,14 @@
/*******************************************************************************
* Helper functions
******************************************************************************/
+static void
+rtc_error_func(void* context, enum RTCError err, const char* str)
+{
+ struct s3d_device* dev = (struct s3d_device*)context;
+ ASSERT(dev);
+ log_error(dev, "Embree:error: %s -- %s", str, rtc_error_string(err));
+}
+
static INLINE void
log_msg
(struct s3d_device* dev,
@@ -63,7 +73,7 @@ device_release(ref_T* ref)
dev = CONTAINER_OF(ref, struct s3d_device, ref);
ASSERT(flist_name_is_empty(&dev->names) == 1);
flist_name_release(&dev->names);
- rtcDeleteDevice(dev->rtc);
+ rtcReleaseDevice(dev->rtc);
MEM_RM(dev->allocator, dev);
}
@@ -97,7 +107,17 @@ s3d_device_create
dev->verbose = verbose;
flist_name_init(allocator, &dev->names);
ref_init(&dev->ref);
+
dev->rtc = rtcNewDevice(verbose ? "verbose=1" : NULL);
+ if(dev->rtc == NULL) {
+ const enum RTCError err = rtcGetDeviceError(NULL);
+ log_error(dev, "Could not create the embree device -- %s.\n",
+ rtc_error_string(err));
+ res = rtc_error_to_res_T(err);
+ goto error;
+ }
+
+ rtcSetDeviceErrorFunction(dev->rtc, rtc_error_func, dev);
exit:
if(out_dev) *out_dev = dev;
diff --git a/src/s3d_geometry.c b/src/s3d_geometry.c
@@ -43,41 +43,76 @@
* Helper functions
******************************************************************************/
static FINLINE void
-sphere_ray_setup
- (RTCRay& ray,
- const float tfar,
- struct geometry* geom)
+sphere_ray_hit_setup
+ (const struct RTCIntersectFunctionNArguments* args, const float tfar)
{
- const RTCRay r = ray;
+ struct geometry* geom = args->geometryUserPtr;
+ struct RTCRayN* rayN;
+ struct RTCHitN* hitN;
+ struct RTCHit hit;
+ struct RTCRay ray;
float cos_theta;
- ASSERT(tfar >= 0 && geom);
+ float Ng[3];
+ float u, v;
+ size_t i;
+ ASSERT(args && args->primID == 0 && args->N == 1 && args->valid[0] != 0);
+ geom = args->geometryUserPtr;
+ ASSERT(geom && geom->type == GEOM_SPHERE);
+
+ rayN = RAYHITN_GET_RAYN(args->rayhit, args->N);
+ hitN = RAYHITN_GET_HITN(args->rayhit, args->N);
+
+ rtc_rayN_get_ray(rayN, args->N, 0, &ray);
ray.tfar = tfar;
- ray.geomID = geom->irtc;
- ray.primID = 0;
- f3_mulf(ray.Ng, ray.dir, tfar);
- f3_add(ray.Ng, ray.Ng, ray.org);
- f3_sub(ray.Ng, ray.Ng, geom->data.sphere->pos);
-
- /* Compute the parametric coordinate */
- f3_normalize(ray.Ng, ray.Ng);
- cos_theta = ray.Ng[2];
- ray.v = (1.f - cos_theta) * 0.5f;
+
+ Ng[0] = ray.dir_x*tfar + ray.org_x - geom->data.sphere->pos[0];
+ Ng[1] = ray.dir_y*tfar + ray.org_y - geom->data.sphere->pos[1];
+ Ng[2] = ray.dir_z*tfar + ray.org_z - geom->data.sphere->pos[2];
+
+ f3_normalize(Ng, Ng);
+ cos_theta = Ng[2];
+
+ v = (1.f - cos_theta) * 0.5f;
if(absf(cos_theta) == 1) {
- ray.u = 0;
- } else if(eq_epsf(ray.Ng[0], 0.f, 1.e-6f)) {
- ray.u = ray.Ng[1] > 0 ? 0.25f : 0.75f;
+ u = 0;
+ } else if(eq_epsf(Ng[0], 0.f, 1.e-6f)) {
+ u = Ng[1] > 0 ? 0.25f : 0.75f;
} else {
- double phi = atan2f(ray.Ng[1], ray.Ng[0]); /* phi in [-PI, PI] */
+ double phi = atan2f(Ng[1], Ng[0]); /* phi in [-PI, PI] */
if(phi < 0) phi = 2*PI + phi; /* phi in [0, 2PI] */
- ray.u = (float)(phi / (2*PI));
+ u = (float)(phi / (2*PI));
+ }
+
+ hit.Ng_x = Ng[0];
+ hit.Ng_y = Ng[1];
+ hit.Ng_z = Ng[2];
+ hit.u = u;
+ hit.v = v;
+ hit.primID = 0;
+ hit.geomID = geom->rtc_id;
+ FOR_EACH(i, 0, RTC_MAX_INSTANCE_LEVEL_COUNT) {
+ hit.instID[i] = args->context->instID[i];
}
/* Filter the intersection if required */
if(geom->data.sphere->filter.func) {
- rtc_hit_filter_wrapper(&geom->data.sphere->filter, ray);
- if(ray.geomID == RTC_INVALID_GEOMETRY_ID) ray = r; /* Restore the ray */
+ struct RTCFilterFunctionNArguments filter_args;
+ int valid;
+
+ filter_args.valid = &valid;
+ filter_args.geometryUserPtr = args->geometryUserPtr;
+ filter_args.context = args->context;
+ filter_args.ray = (struct RTCRayN*)&ray;
+ filter_args.hit = (struct RTCHitN*)&hit;
+ filter_args.N = args->N;
+
+ rtc_hit_filter_wrapper(&filter_args);
+ if(!filter_args.valid[0]) return;
}
+
+ RAYN_GRAB(rayN, args->N, 0, float, tfar) = tfar;
+ rtc_hitN_set_hit(hitN, args->N, 0, &hit);
}
static void
@@ -126,11 +161,12 @@ geometry_create
S3D(device_ref_get(dev));
geom->dev = dev;
geom->name = S3D_INVALID_ID;
- geom->irtc = RTC_INVALID_GEOMETRY_ID;
geom->flip_surface = 0;
geom->is_enabled = 1;
geom->type = GEOM_NONE;
geom->data.mesh = NULL;
+ geom->rtc = NULL;
+ geom->rtc_id = RTC_INVALID_GEOMETRY_ID;
exit:
*out_geom = geom;
@@ -158,62 +194,51 @@ geometry_ref_put(struct geometry* geom)
}
void
-geometry_rtc_sphere_bounds(void* data, size_t item, RTCBounds& bounds)
+geometry_rtc_sphere_bounds(const struct RTCBoundsFunctionArguments* args)
{
- struct geometry* geom = (struct geometry*)data;
+ struct geometry* geom;
struct sphere sphere;
- ASSERT(geom && item == 0 && geom->type == GEOM_SPHERE);
- (void)item;
- sphere = *geom->data.sphere;
- bounds.lower_x = sphere.pos[0] - sphere.radius;
- bounds.lower_y = sphere.pos[1] - sphere.radius;
- bounds.lower_z = sphere.pos[2] - sphere.radius;
- bounds.upper_x = sphere.pos[0] + sphere.radius;
- bounds.upper_y = sphere.pos[1] + sphere.radius;
- bounds.upper_z = sphere.pos[2] + sphere.radius;
-}
+ ASSERT(args && args->primID == 0 && args->timeStep == 0);
-void
-geometry_rtc_sphere_intersect(void* data, RTCRay& ray, size_t item)
-{
- float v[3];
- float A, B, C, D, Q, rcpA, t0, t1;
- struct geometry* geom = (struct geometry*)data;
- struct sphere sphere;
- ASSERT(geom && item == 0 && geom->type == GEOM_SPHERE);
- (void)item;
+ geom = args->geometryUserPtr;
+ ASSERT(geom && geom->type == GEOM_SPHERE);
sphere = *geom->data.sphere;
- f3_sub(v, ray.org, sphere.pos);
- A = f3_dot(ray.dir, ray.dir);
- B = 2*f3_dot(v, ray.dir);
- C = f3_dot(v, v) - sphere.radius*sphere.radius;
- D = B*B - 4*A*C;
-
- if(D < 0.0f) return;
- Q = (float)sqrt(D);
- rcpA = 1.f / A;
- t0 = 0.5f * rcpA * (-B - Q);
- t1 = 0.5f * rcpA * (-B + Q);
-
- if(ray.tnear < t0 && t0 < ray.tfar) sphere_ray_setup(ray, t0, geom);
- if(ray.tnear < t1 && t1 < ray.tfar) sphere_ray_setup(ray, t1, geom);
+ args->bounds_o->lower_x = sphere.pos[0] - sphere.radius;
+ args->bounds_o->lower_y = sphere.pos[1] - sphere.radius;
+ args->bounds_o->lower_z = sphere.pos[2] - sphere.radius;
+ args->bounds_o->upper_x = sphere.pos[0] + sphere.radius;
+ args->bounds_o->upper_y = sphere.pos[1] + sphere.radius;
+ args->bounds_o->upper_z = sphere.pos[2] + sphere.radius;
}
void
-geometry_rtc_sphere_occluded(void* data, RTCRay& ray, size_t item)
+geometry_rtc_sphere_intersect(const struct RTCIntersectFunctionNArguments* args)
{
float v[3];
+ float ray_org[3];
+ float ray_dir[3];
float A, B, C, D, Q, rcpA, t0, t1;
- struct geometry* geom = (struct geometry*)data;
+ struct geometry* geom;
struct sphere sphere;
- ASSERT(geom && item == 0 && geom->type == GEOM_SPHERE);
- (void)item;
+ struct RTCRayN* rayN;
+ ASSERT(args && args->primID == 0 && args->N == 1 && args->valid[0] != 0);
+
+ geom = args->geometryUserPtr;
+ ASSERT(geom && geom->type == GEOM_SPHERE);
+
+ rayN = RAYHITN_GET_RAYN(args->rayhit, args->N);
+ ray_org[0] = RAYN_GRAB(rayN, args->N, 0, float, org_x);
+ ray_org[1] = RAYN_GRAB(rayN, args->N, 0, float, org_y);
+ ray_org[2] = RAYN_GRAB(rayN, args->N, 0, float, org_z);
+ ray_dir[0] = RAYN_GRAB(rayN, args->N, 0, float, dir_x);
+ ray_dir[1] = RAYN_GRAB(rayN, args->N, 0, float, dir_y);
+ ray_dir[2] = RAYN_GRAB(rayN, args->N, 0, float, dir_z);
sphere = *geom->data.sphere;
- f3_sub(v, ray.org, sphere.pos);
- A = f3_dot(ray.dir, ray.dir);
- B = 2*f3_dot(v, ray.dir);
+ f3_sub(v, ray_org, sphere.pos);
+ A = f3_dot(ray_dir, ray_dir);
+ B = 2*f3_dot(v, ray_dir);
C = f3_dot(v, v) - sphere.radius*sphere.radius;
D = B*B - 4*A*C;
@@ -223,7 +248,13 @@ geometry_rtc_sphere_occluded(void* data, RTCRay& ray, size_t item)
t0 = 0.5f * rcpA * (-B - Q);
t1 = 0.5f * rcpA * (-B + Q);
- if(ray.tnear < t0 && t0 < ray.tfar) ray.geomID = 0;
- if(ray.tnear < t1 && t1 < ray.tfar) ray.geomID = 0;
+ if(RAYN_GRAB(rayN, args->N, 0, float, tnear) < t0
+ && RAYN_GRAB(rayN, args->N, 0, float, tfar) > t0) {
+ sphere_ray_hit_setup(args, t0);
+ }
+ if(RAYN_GRAB(rayN, args->N, 0, float, tnear) < t1
+ && RAYN_GRAB(rayN, args->N, 0, float, tfar) > t1) {
+ sphere_ray_hit_setup(args, t1);
+ }
}
diff --git a/src/s3d_geometry.h b/src/s3d_geometry.h
@@ -49,15 +49,16 @@ enum embree_attrib {
EMBREE_ENABLE = BIT(0),
EMBREE_FILTER_FUNCTION = BIT(1),
EMBREE_INDICES = BIT(2),
- EMBREE_TRANSFORM = BIT(4),
- EMBREE_VERTICES = BIT(5),
- EMBREE_USER_GEOMETRY = BIT(6)
+ EMBREE_TRANSFORM = BIT(3),
+ EMBREE_VERTICES = BIT(4),
+ EMBREE_USER_GEOMETRY = BIT(5)
};
/* Backend geometry */
struct geometry {
unsigned name; /* Client side identifier */
- unsigned irtc; /* Embree identifier */
+ RTCGeometry rtc; /* Embree geometry */
+ unsigned rtc_id; /* Embree geometry identifier */
unsigned scene_prim_id_offset; /* Offset from local to scene prim_id */
int embree_outdated_mask; /* Combination of embree_attrib */
@@ -91,21 +92,11 @@ geometry_ref_put
extern LOCAL_SYM void
geometry_rtc_sphere_bounds
- (void* data,
- size_t item,
- RTCBounds& bounds);
+ (const struct RTCBoundsFunctionArguments* args);
extern LOCAL_SYM void
geometry_rtc_sphere_intersect
- (void* data,
- RTCRay& ray,
- size_t item);
-
-extern LOCAL_SYM void
-geometry_rtc_sphere_occluded
- (void* data,
- RTCRay& ray,
- size_t item);
+ (const struct RTCIntersectFunctionNArguments* args);
#endif /* S3D_GEOMETRY_H */
diff --git a/src/s3d_scene_view.c b/src/s3d_scene_view.c
@@ -36,17 +36,17 @@
#include "s3d_scene_view_c.h"
#include "s3d_shape_c.h"
+#include <rsys/algorithm.h>
#include <rsys/float3.h>
#include <rsys/float33.h>
#include <rsys/mem_allocator.h>
-#include <algorithm>
-
-struct ray_extended : public RTCRay {
+struct intersect_context {
+ struct RTCIntersectContext rtc;
struct s3d_scene_view* scnview;
+ void* data; /* Per ray user defined data */
float ws_org[3]; /* World space ray origin */
float ws_dir[3]; /* World space ray direction */
- void* data; /* User defined data */
};
/*******************************************************************************
@@ -59,33 +59,47 @@ aabb_is_degenerated(const float lower[3], const float upper[3])
return lower[0] > upper[0] || lower[1] > upper[1] || lower[2] > upper[2];
}
-static FINLINE bool
-operator < (const struct fltui& it, const float val)
+static INLINE int
+cmp_float(const void* a, const void* b)
{
- /* This operator is used by the std::lower_bound algorithm that returns an
- * iterator to the first element that is not less than val while one expect
- * an iterator on the first element that is not less *or equal* than val.
- * That's why we use <= rather than < */
- return it.flt <= val;
+ const float key = *(const float*)a;
+ const float val = *(const float*)b;
+ if(key < val) return -1;
+ if(key > val) return +1;
+ return 0;
}
-static FINLINE bool
-operator < (const struct nprims_cdf& it, const size_t iprim)
+static INLINE int
+cmp_float_to_fltui(const void* a, const void* b)
{
- /* This operator is used by the std::lower_bound algorithm that returns an
- * iterator to the first element that is not less than iprim while one expect
- * an iterator on the first element that is not less *or equal* than iprim.
- * That's why we use <= rather than < */
- return it.nprims <= iprim;
+ const float key = *(const float*)a;
+ const struct fltui* fltui = (const struct fltui*)b;
+ if(key < fltui->flt) return -1;
+ if(key > fltui->flt) return +1;
+ return 0;
+}
+
+static INLINE int
+cmp_size_t_to_nprims_cdf(const void* a, const void* b)
+{
+ const size_t key = *(const size_t*)a;
+ const struct nprims_cdf* nprims_cdf = (const struct nprims_cdf*)b;
+ if(key < nprims_cdf->nprims) return -1;
+ if(key > nprims_cdf->nprims) return +1;
+ return 0;
}
static INLINE void
scene_view_destroy_geometry(struct s3d_scene_view* scnview, struct geometry* geom)
{
ASSERT(geom);
- if(geom->irtc != RTC_INVALID_GEOMETRY_ID) {
- rtcDeleteGeometry(scnview->rtc_scn, geom->irtc);
- geom->irtc = RTC_INVALID_GEOMETRY_ID;
+ if(geom->rtc) {
+ if(geom->rtc_id != RTC_INVALID_GEOMETRY_ID) {
+ rtcDetachGeometry(scnview->rtc_scn, geom->rtc_id);
+ geom->rtc_id = RTC_INVALID_GEOMETRY_ID;
+ }
+ rtcReleaseGeometry(geom->rtc);
+ geom->rtc = NULL;
scnview->rtc_scn_update = 1; /* Notify the scene upd */
}
geometry_ref_put(geom);
@@ -125,28 +139,32 @@ on_shape_detach
}
static INLINE void
-hit_setup(struct s3d_scene_view* scnview, const RTCRay* ray, struct s3d_hit* hit)
+hit_setup
+ (struct s3d_scene_view* scnview,
+ const struct RTCRayHit* ray_hit,
+ struct s3d_hit* hit)
{
float w;
char flip_surface = 0;
- ASSERT(scnview && hit && ray);
+ ASSERT(scnview && hit && ray_hit);
- if((unsigned)ray->geomID == RTC_INVALID_GEOMETRY_ID) { /* No hit */
+ if(ray_hit->hit.geomID == RTC_INVALID_GEOMETRY_ID) { /* No hit */
*hit = S3D_HIT_NULL;
return;
}
- f3_set(hit->normal, ray->Ng);
- hit->distance = ray->tfar;
+ hit->normal[0] = ray_hit->hit.Ng_x;
+ hit->normal[1] = ray_hit->hit.Ng_y;
+ hit->normal[2] = ray_hit->hit.Ng_z;
+ hit->distance = ray_hit->ray.tfar;
- if((unsigned)ray->instID == RTC_INVALID_GEOMETRY_ID) {
+ if(ray_hit->hit.instID[0] == RTC_INVALID_GEOMETRY_ID) {
struct geometry* geom_shape;
- ASSERT((unsigned)ray->geomID < darray_geom_size_get(&scnview->embree2geoms));
- geom_shape = scene_view_geometry_from_embree_id(scnview, ray->geomID);
+ geom_shape = scene_view_geometry_from_embree_id(scnview, ray_hit->hit.geomID);
hit->prim.shape__ = geom_shape;
hit->prim.inst__ = NULL;
- hit->prim.prim_id = ray->primID;
+ hit->prim.prim_id = ray_hit->hit.primID;
hit->prim.geom_id = geom_shape->name;
hit->prim.inst_id = S3D_INVALID_ID;
hit->prim.scene_prim_id = /* Compute the "scene space" primitive id */
@@ -158,16 +176,16 @@ hit_setup(struct s3d_scene_view* scnview, const RTCRay* ray, struct s3d_hit* hit
struct geometry* geom_inst;
struct geometry* geom_shape;
float transform[9];
- ASSERT((unsigned)ray->instID < darray_geom_size_get(&scnview->embree2geoms));
- geom_inst = scene_view_geometry_from_embree_id(scnview, ray->instID);
+ geom_inst = scene_view_geometry_from_embree_id
+ (scnview, ray_hit->hit.instID[0]);
geom_shape = scene_view_geometry_from_embree_id
- (geom_inst->data.instance->scnview, ray->geomID);
+ (geom_inst->data.instance->scnview, ray_hit->hit.geomID);
hit->prim.shape__ = geom_shape;
hit->prim.inst__ = geom_inst;
- hit->prim.prim_id = ray->primID;
+ hit->prim.prim_id = ray_hit->hit.primID;
hit->prim.geom_id = geom_shape->name;
hit->prim.inst_id = geom_inst->name;
- hit->prim.scene_prim_id = /* Compute the "scene space" */
+ hit->prim.scene_prim_id = /* Compute the "scene space" primitive id */
hit->prim.prim_id /* Shape space */
+ geom_shape->scene_prim_id_offset /* Inst space */
+ geom_inst->scene_prim_id_offset; /* Scene space */
@@ -184,9 +202,8 @@ hit_setup(struct s3d_scene_view* scnview, const RTCRay* ray, struct s3d_hit* hit
ASSERT(((struct geometry*)hit->prim.shape__)->type == GEOM_MESH
||((struct geometry*)hit->prim.shape__)->type == GEOM_SPHERE);
-
- hit->uv[0] = ray->u;
- hit->uv[1] = ray->v;
+ hit->uv[0] = ray_hit->hit.u;
+ hit->uv[1] = ray_hit->hit.v;
if(((struct geometry*)hit->prim.shape__)->type == GEOM_MESH) {
w = 1.f - hit->uv[0] - hit->uv[1];
@@ -218,75 +235,108 @@ embree_geometry_register
ASSERT(scnview && geom);
/* Create the Embree geometry if it is not valid */
- if(geom->irtc != RTC_INVALID_GEOMETRY_ID) {
+ if(geom->rtc != NULL) {
if(geom->type == GEOM_INSTANCE) {
/* If the geometry is an instance one have to update it if the
* instantiated geometry was updated. Currently, we have no simple way to
* know if the geometry was upd or not so we simply force the update. */
- rtcUpdate(scnview->rtc_scn, geom->irtc);
+ rtcCommitGeometry(geom->rtc);
}
} else {
switch(geom->type) {
case GEOM_MESH:
- geom->irtc = rtcNewTriangleMesh(scnview->rtc_scn, RTC_GEOMETRY_DYNAMIC,
- mesh_get_ntris(geom->data.mesh), mesh_get_nverts(geom->data.mesh));
+ geom->rtc = rtcNewGeometry
+ (scnview->scn->dev->rtc, RTC_GEOMETRY_TYPE_TRIANGLE);
break;
case GEOM_INSTANCE:
- geom->irtc = rtcNewInstance2
- (scnview->rtc_scn, geom->data.instance->scnview->rtc_scn);
+ geom->rtc = rtcNewGeometry
+ (scnview->scn->dev->rtc, RTC_GEOMETRY_TYPE_INSTANCE);
break;
case GEOM_SPHERE:
- geom->irtc = rtcNewUserGeometry3
- (scnview->rtc_scn, RTC_GEOMETRY_DYNAMIC, 1);
- rtcSetUserData(scnview->rtc_scn, geom->irtc, geom);
- rtcSetBoundsFunction
- (scnview->rtc_scn, geom->irtc, geometry_rtc_sphere_bounds);
- rtcSetIntersectFunction
- (scnview->rtc_scn, geom->irtc, geometry_rtc_sphere_intersect);
- rtcSetOccludedFunction
- (scnview->rtc_scn, geom->irtc, geometry_rtc_sphere_occluded);
+ geom->rtc = rtcNewGeometry
+ (scnview->scn->dev->rtc, RTC_GEOMETRY_TYPE_USER);
+ /* Setup geometry callbacks. Note that the "occluded" is not set since
+ * rtcOccluded is not used */
+ rtcSetGeometryUserPrimitiveCount(geom->rtc, 1);
+ rtcSetGeometryBoundsFunction(geom->rtc, geometry_rtc_sphere_bounds, NULL);
+ rtcSetGeometryIntersectFunction(geom->rtc, geometry_rtc_sphere_intersect);
break;
default: FATAL("Unreachable code\n"); break;
}
- if(geom->irtc == RTC_INVALID_GEOMETRY_ID)
+ if(geom->rtc == NULL)
return RES_UNKNOWN_ERR;
- scnview->rtc_scn_update = 1;
- }
+ /* Set the Star-3D representation of the geometry to the Embree geometry */
+ rtcSetGeometryUserData(geom->rtc, geom);
- /* Register the embree geometry in the embree2geoms associative array */
- if(geom->irtc >= darray_geom_size_get(&scnview->embree2geoms)) {
- const res_T res = darray_geom_resize(&scnview->embree2geoms, geom->irtc+1);
- if(res != RES_OK) {
- rtcDeleteGeometry(scnview->rtc_scn, geom->irtc);
- geom->irtc = RTC_INVALID_GEOMETRY_ID;
- return res;
- }
+ /* Attach the Embree geometry to the Embree scene of the scene view */
+ geom->rtc_id = rtcAttachGeometry(scnview->rtc_scn, geom->rtc);
+
+ scnview->rtc_scn_update = 1;
}
- darray_geom_data_get(&scnview->embree2geoms)[geom->irtc] = geom;
return RES_OK;
}
-static INLINE void
+static INLINE res_T
embree_geometry_setup_positions
(struct s3d_scene_view* scnview, struct geometry* geom)
{
- ASSERT(scnview && geom && geom->type == GEOM_MESH);
- ASSERT(geom->irtc != RTC_INVALID_GEOMETRY_ID);
- rtcSetBuffer(scnview->rtc_scn, geom->irtc, RTC_VERTEX_BUFFER,
- mesh_get_pos(geom->data.mesh), 0, sizeof(float[3]));
- rtcUpdateBuffer(scnview->rtc_scn, geom->irtc, RTC_VERTEX_BUFFER);
+ RTCBuffer buf = NULL;
+ float* verts;
+ size_t nverts;
+ res_T res = RES_OK;
+ ASSERT(scnview && geom && geom->type == GEOM_MESH && geom->rtc);
+
+ verts = mesh_get_pos(geom->data.mesh);
+ nverts = mesh_get_nverts(geom->data.mesh);
+
+ buf = rtcNewSharedBuffer
+ (scnview->scn->dev->rtc, verts, sizeof(float[3])*nverts);
+ if(!buf) {
+ res = rtc_error_to_res_T(rtcGetDeviceError(scnview->scn->dev->rtc));
+ goto error;
+ }
+
+ rtcSetGeometryBuffer(geom->rtc, RTC_BUFFER_TYPE_VERTEX, 0/*slot*/,
+ RTC_FORMAT_FLOAT3, buf, 0/*offset*/, sizeof(float[3])/*stride*/, nverts);
+ rtcUpdateGeometryBuffer(geom->rtc, RTC_BUFFER_TYPE_VERTEX, 0);
+
+exit:
+ if(buf) rtcReleaseBuffer(buf);
+ return res;
+error:
+ goto exit;
}
-static INLINE void
+static INLINE res_T
embree_geometry_setup_indices
(struct s3d_scene_view* scnview, struct geometry* geom)
{
- ASSERT(scnview && geom && geom->type == GEOM_MESH);
- ASSERT(geom->irtc != RTC_INVALID_GEOMETRY_ID);
- rtcSetBuffer(scnview->rtc_scn, geom->irtc, RTC_INDEX_BUFFER,
- mesh_get_ids(geom->data.mesh), 0, sizeof(uint32_t[3]));
- rtcUpdateBuffer(scnview->rtc_scn, geom->irtc, RTC_INDEX_BUFFER);
+ RTCBuffer buf = NULL;
+ size_t ntris;
+ uint32_t* ids;
+ res_T res = RES_OK;
+ ASSERT(scnview && geom && geom->type == GEOM_MESH && geom->rtc);
+
+ ids = mesh_get_ids(geom->data.mesh);
+ ntris = mesh_get_ntris(geom->data.mesh);
+
+ buf = rtcNewSharedBuffer
+ (scnview->scn->dev->rtc, ids, sizeof(uint32_t[3])*ntris*3);
+ if(!buf) {
+ res = rtc_error_to_res_T(rtcGetDeviceError(scnview->scn->dev->rtc));
+ goto error;
+ }
+
+ rtcSetGeometryBuffer(geom->rtc, RTC_BUFFER_TYPE_INDEX, 0/*slot*/,
+ RTC_FORMAT_UINT3, buf, 0/*offset*/, sizeof(uint32_t[3])/*stride*/, ntris*3);
+ rtcUpdateGeometryBuffer(geom->rtc, RTC_BUFFER_TYPE_INDEX, 0);
+
+exit:
+ if(buf) rtcReleaseBuffer(buf);
+ return res;
+error:
+ goto exit;
}
static INLINE void
@@ -295,9 +345,9 @@ embree_geometry_setup_enable_state
{
ASSERT(scnview && geom);
if(geom->is_enabled) {
- rtcEnable(scnview->rtc_scn, geom->irtc);
+ rtcEnableGeometry(geom->rtc);
} else {
- rtcDisable(scnview->rtc_scn, geom->irtc);
+ rtcDisableGeometry(geom->rtc);
}
}
@@ -305,15 +355,13 @@ static INLINE void
embree_geometry_setup_filter_function
(struct s3d_scene_view* scnview, struct geometry* geom)
{
- ASSERT(scnview && geom && geom->irtc != RTC_INVALID_GEOMETRY_ID);
- ASSERT(geom->type == GEOM_MESH);
+ ASSERT(scnview && geom && geom->rtc &&geom->type == GEOM_MESH);
+ (void)scnview;
if(!geom->data.mesh->filter.func) {
- rtcSetIntersectionFilterFunction(scnview->rtc_scn, geom->irtc, NULL);
+ rtcSetGeometryIntersectFilterFunction(geom->rtc, NULL);
} else {
- rtcSetIntersectionFilterFunction
- (scnview->rtc_scn, geom->irtc, rtc_hit_filter_wrapper);
- rtcSetUserData(scnview->rtc_scn, geom->irtc, &geom->data.mesh->filter);
+ rtcSetGeometryIntersectFilterFunction(geom->rtc, rtc_hit_filter_wrapper);
}
}
@@ -321,37 +369,29 @@ static INLINE void
embree_geometry_setup_transform
(struct s3d_scene_view* scnview, struct geometry* geom)
{
- ASSERT(scnview && geom && geom->irtc != RTC_INVALID_GEOMETRY_ID);
+ ASSERT(scnview && geom && geom->rtc != NULL);
ASSERT(geom->type == GEOM_INSTANCE);
- rtcSetTransform2
- (scnview->rtc_scn,
- geom->irtc,
- RTC_MATRIX_COLUMN_MAJOR,
- geom->data.instance->transform);
- rtcUpdate(scnview->rtc_scn, geom->irtc);
+ rtcSetGeometryTransform(geom->rtc, 0, RTC_FORMAT_FLOAT3X4_COLUMN_MAJOR,
+ geom->data.instance->transform);
}
static INLINE res_T
scene_view_setup_embree(struct s3d_scene_view* scnview)
{
struct htable_geom_iterator it, end;
- const RTCSceneFlags rtc_scene_mask =
- RTC_SCENE_DYNAMIC
- | RTC_SCENE_INCOHERENT
- | RTC_SCENE_ROBUST;
- const RTCAlgorithmFlags rtc_algorithm_mask = RTC_INTERSECT1;
int rtc_outdated = 0;
res_T res = RES_OK;
ASSERT(scnview);
/* The rtc_scn could be already allocated since the scene views are cached */
if(!scnview->rtc_scn) {
- scnview->rtc_scn = rtcDeviceNewScene
- (scnview->scn->dev->rtc, rtc_scene_mask, rtc_algorithm_mask);
+ scnview->rtc_scn = rtcNewScene(scnview->scn->dev->rtc);
if(!scnview->rtc_scn) {
- res = RES_MEM_ERR;
+ res = rtc_error_to_res_T(rtcGetDeviceError(scnview->scn->dev->rtc));
goto error;
}
+ rtcSetSceneFlags
+ (scnview->rtc_scn, RTC_SCENE_FLAG_ROBUST | RTC_SCENE_FLAG_DYNAMIC);
rtc_outdated = 1;
}
@@ -384,8 +424,10 @@ scene_view_setup_embree(struct s3d_scene_view* scnview)
embree_geometry_setup_filter_function(scnview, geom);
if((geom->embree_outdated_mask & EMBREE_TRANSFORM) != 0)
embree_geometry_setup_transform(scnview, geom);
- if((geom->embree_outdated_mask & EMBREE_USER_GEOMETRY) != 0)
- rtcUpdate(scnview->rtc_scn, geom->irtc);
+
+ /* Commit the updated geometry */
+ if(geom->embree_outdated_mask)
+ rtcCommitGeometry(geom->rtc);
geom->embree_outdated_mask = 0;
}
@@ -394,8 +436,8 @@ scene_view_setup_embree(struct s3d_scene_view* scnview)
/* Commit the embree changes */
if(rtc_outdated) {
- rtcCommit(scnview->rtc_scn);
- scnview->rtc_commit = 1; /* Notify that the scene view was committed */
+ rtcCommitScene(scnview->rtc_scn);
+ scnview->rtc_commit = 1; /* Notify that the RTC scene was committed */
scnview->rtc_scn_update = 0;
}
@@ -403,10 +445,9 @@ exit:
return res;
error:
if(scnview->rtc_scn) {
+ rtcReleaseScene(scnview->rtc_scn);
scnview->rtc_scn = NULL;
- rtcDeleteScene(scnview->rtc_scn);
}
- darray_geom_clear(&scnview->embree2geoms);
goto exit;
}
@@ -441,9 +482,13 @@ scene_view_register_mesh
/* Discard the geometry that is not geometrically valid */
if(!shape->data.mesh->indices || !shape->data.mesh->attribs[S3D_POSITION]) {
- if(geom->irtc != RTC_INVALID_GEOMETRY_ID) {
- rtcDeleteGeometry(scnview->rtc_scn, geom->irtc);
- geom->irtc = RTC_INVALID_GEOMETRY_ID;
+ if(geom->rtc != NULL) {
+ if(geom->rtc_id != RTC_INVALID_GEOMETRY_ID) {
+ rtcDetachGeometry(scnview->rtc_scn, geom->rtc_id);
+ geom->rtc_id = RTC_INVALID_GEOMETRY_ID;
+ }
+ rtcReleaseGeometry(geom->rtc);
+ geom->rtc = NULL;
}
mesh_clear(geom->data.mesh);
goto exit;
@@ -463,17 +508,21 @@ scene_view_register_mesh
/* Get a reference onto the shape mesh attribs */
FOR_EACH(iattr, 0, S3D_ATTRIBS_COUNT__) {
- geom->embree_outdated_mask |= EMBREE_VERTICES;
if(geom->data.mesh->attribs[iattr] == shape->data.mesh->attribs[iattr])
continue;
+ geom->embree_outdated_mask |= EMBREE_VERTICES;
+
if(geom->data.mesh->attribs[iattr]) { /* Release the previous buffer */
vertex_buffer_ref_put(geom->data.mesh->attribs[iattr]);
geom->data.mesh->attribs[iattr] = NULL;
}
+
+ /* This attrib does not exist anymore */
if(!shape->data.mesh->attribs[iattr])
continue;
+ /* Get the new buffer */
vertex_buffer_ref_get(shape->data.mesh->attribs[iattr]);
geom->data.mesh->attribs[iattr] = shape->data.mesh->attribs[iattr];
geom->data.mesh->attribs_type[iattr] = shape->data.mesh->attribs_type[iattr];
@@ -492,23 +541,6 @@ scene_view_register_mesh
geom->embree_outdated_mask |= EMBREE_FILTER_FUNCTION;
}
- if(geom->irtc != RTC_INVALID_GEOMETRY_ID) {
- struct index_buffer* shape_ids = shape->data.mesh->indices;
- struct index_buffer* geom_ids = geom->data.mesh->indices;
- struct vertex_buffer* shape_verts = shape->data.mesh->attribs[S3D_POSITION];
- struct vertex_buffer* geom_verts = geom->data.mesh->attribs[S3D_POSITION];
- const size_t shape_nids = darray_u32_size_get(&shape_ids->data);
- const size_t geom_nids = darray_u32_size_get(&geom_ids->data);
- const size_t shape_nverts = darray_float_size_get(&shape_verts->data);
- const size_t geom_nverts = darray_float_size_get(&geom_verts->data);
-
- /* The shape mesh was resize => the Embree geometry is no more valid */
- if(shape_nids != geom_nids || shape_nverts != geom_nverts) {
- rtcDeleteGeometry(scnview->rtc_scn, geom->irtc);
- geom->irtc = RTC_INVALID_GEOMETRY_ID;
- }
- }
-
geom->flip_surface = shape->flip_surface;
exit:
@@ -647,7 +679,7 @@ scene_view_register_instance
/* Update the enable flag */
if(geom->is_enabled != shape->is_enabled) {
geom->is_enabled = shape->is_enabled;
- geom->embree_outdated_mask |= EMBREE_TRANSFORM;
+ geom->embree_outdated_mask |= EMBREE_ENABLE;
}
geom->flip_surface = shape->flip_surface;
@@ -833,7 +865,7 @@ scene_view_compute_scene_aabb(struct s3d_scene_view* scnview)
}
}
-float
+static float
scene_view_compute_volume
(struct s3d_scene_view* scnview,
const char flip_surface)
@@ -960,7 +992,6 @@ scene_view_create(struct s3d_scene* scn, struct s3d_scene_view** out_scnview)
}
list_init(&scnview->node);
htable_geom_init(scn->dev->allocator, &scnview->cached_geoms);
- darray_geom_init(scn->dev->allocator, &scnview->embree2geoms);
darray_fltui_init(scn->dev->allocator, &scnview->cdf);
darray_nprims_cdf_init(scn->dev->allocator, &scnview->nprims_cdf);
htable_instview_init(scn->dev->allocator, &scnview->instviews);
@@ -1031,7 +1062,6 @@ scene_view_release(ref_T* ref)
/* Clear the scnview data structures excepted the cache of geometries that
* will be used to speed up the future scnview creation */
- darray_geom_clear(&scnview->embree2geoms);
darray_fltui_clear(&scnview->cdf);
darray_nprims_cdf_clear(&scnview->nprims_cdf);
f3_splat(scnview->lower, FLT_MAX);
@@ -1120,7 +1150,10 @@ s3d_scene_view_trace_ray
void* ray_data,
struct s3d_hit* hit)
{
- struct ray_extended ray_ex;
+ struct RTCRayHit ray_hit;
+ struct intersect_context intersect_ctx;
+ size_t i;
+
if(!scnview || !org || !dir || !range || !hit)
return RES_BAD_ARG;
if(!f3_is_normalized(dir)) {
@@ -1129,6 +1162,12 @@ s3d_scene_view_trace_ray
FUNC_NAME, SPLIT3(dir));
return RES_BAD_ARG;
}
+ if(range[0] < 0) {
+ log_error(scnview->scn->dev,
+ "%s: invalid ray range [%g, %g] - it must be in [0, INF).\n",
+ FUNC_NAME, range[0], range[1]);
+ return RES_BAD_ARG;
+ }
if((scnview->mask & S3D_TRACE) == 0) {
log_error(scnview->scn->dev,
"%s: the S3D_TRACE flag is not active onto the submitted scene view.\n",
@@ -1140,22 +1179,41 @@ s3d_scene_view_trace_ray
return RES_OK;
}
- f3_set(ray_ex.org, f3_set(ray_ex.ws_org, org));
- f3_set(ray_ex.dir, f3_set(ray_ex.ws_dir, dir));
- ray_ex.tnear = range[0];
- ray_ex.tfar = range[1];
- ray_ex.geomID = RTC_INVALID_GEOMETRY_ID;
- ray_ex.primID = RTC_INVALID_GEOMETRY_ID;
- ray_ex.instID = RTC_INVALID_GEOMETRY_ID;
- ray_ex.mask = 0xFFFFFFFF;
- ray_ex.time = 0.f;
- ray_ex.scnview = scnview;
- ray_ex.data = ray_data;
-
- rtcIntersect(scnview->rtc_scn, ray_ex);
-
- hit_setup(scnview, &ray_ex, hit);
-
+ /* Initialise the ray */
+ ray_hit.ray.org_x = org[0];
+ ray_hit.ray.org_y = org[1];
+ ray_hit.ray.org_z = org[2];
+ ray_hit.ray.dir_x = dir[0];
+ ray_hit.ray.dir_y = dir[1];
+ ray_hit.ray.dir_z = dir[2];
+ ray_hit.ray.tnear = range[0];
+ ray_hit.ray.tfar = range[1];
+ ray_hit.ray.time = FLT_MAX; /* Invalid fields */
+ ray_hit.ray.mask = UINT_MAX; /* Invalid fields */
+ ray_hit.ray.id = UINT_MAX; /* Invalid fields */
+ ray_hit.ray.flags = UINT_MAX; /* Invalid fields */
+
+ /* Initialise the hit */
+ ray_hit.hit.geomID = RTC_INVALID_GEOMETRY_ID;
+ FOR_EACH(i, 0, RTC_MAX_INSTANCE_LEVEL_COUNT) {
+ ray_hit.hit.instID[i] = RTC_INVALID_GEOMETRY_ID;
+ }
+
+ /* Initialise the intersect context */
+ rtcInitIntersectContext(&intersect_ctx.rtc);
+ intersect_ctx.ws_org[0] = org[0];
+ intersect_ctx.ws_org[1] = org[1];
+ intersect_ctx.ws_org[2] = org[2];
+ intersect_ctx.ws_dir[0] = dir[0];
+ intersect_ctx.ws_dir[1] = dir[1];
+ intersect_ctx.ws_dir[2] = dir[2];
+ intersect_ctx.scnview = scnview;
+ intersect_ctx.data = ray_data;
+
+ /* Here we go! */
+ rtcIntersect1(scnview->rtc_scn, &intersect_ctx.rtc, &ray_hit);
+
+ hit_setup(scnview, &ray_hit, hit);
return RES_OK;
}
@@ -1208,8 +1266,9 @@ s3d_scene_view_sample
{
struct geometry** pgeom;
struct geometry* geom;
- const struct fltui* fltui_begin, *fltui_end, *fltui_found;
- const float* flt_begin, *flt_end, *flt_found;
+ size_t sz;
+ const struct fltui* fltui, *fltui_found;
+ const float* flt, *flt_found;
unsigned ishape;
float f;
res_T res = RES_OK;
@@ -1244,15 +1303,16 @@ s3d_scene_view_sample
/* Map u to the CDF bounds */
f = u * darray_fltui_cdata_get(&scnview->cdf)[0].flt;
} else {
- fltui_begin = darray_fltui_cdata_get(&scnview->cdf);
- fltui_end = fltui_begin + darray_fltui_size_get(&scnview->cdf);
- f = u * fltui_end[-1].flt; /* Map u to the CDF bounds */
- fltui_found = std::lower_bound(fltui_begin, fltui_end, f);
- ASSERT(fltui_found != fltui_end);
+ fltui = darray_fltui_cdata_get(&scnview->cdf);
+ sz = darray_fltui_size_get(&scnview->cdf);
+ f = u * fltui[sz-1].flt; /* Map u to the CDF bounds */
+ fltui_found = search_lower_bound
+ (&f, fltui, sz, sizeof(*fltui), cmp_float_to_fltui);
+ ASSERT(fltui_found);
ishape = fltui_found->ui;
/* Transform u to the geometry CDF bounds */
- if(fltui_found != fltui_begin)
+ if(fltui_found != fltui)
f -= fltui_found[-1].flt;
}
pgeom = htable_geom_find(&scnview->cached_geoms, &ishape);
@@ -1272,14 +1332,15 @@ s3d_scene_view_sample
if(darray_fltui_size_get(&geom->data.instance->scnview->cdf) == 1) {
ishape = darray_fltui_cdata_get(&geom->data.instance->scnview->cdf)[0].ui;
} else {
- fltui_begin = darray_fltui_cdata_get(&geom->data.instance->scnview->cdf);
- fltui_end = fltui_begin + darray_fltui_size_get(&geom->data.instance->scnview->cdf);
- fltui_found = std::lower_bound(fltui_begin, fltui_end, f);
- ASSERT(fltui_found != fltui_end);
+ fltui = darray_fltui_cdata_get(&geom->data.instance->scnview->cdf);
+ sz = darray_fltui_size_get(&geom->data.instance->scnview->cdf);
+ fltui_found = search_lower_bound
+ (&f, fltui, sz, sizeof(*fltui), cmp_float_to_fltui);
+ ASSERT(fltui_found != NULL);
ishape = fltui_found->ui;
/* Transform u to the geometry CDF bounds */
- if(fltui_found != fltui_begin)
+ if(fltui_found != fltui)
f -= fltui_found[-1].flt;
}
pgeom = htable_geom_find(&geom->data.instance->scnview->cached_geoms, &ishape);
@@ -1294,11 +1355,11 @@ s3d_scene_view_sample
if(geom->type == GEOM_SPHERE) {
primitive->prim_id = 0;
} else if(geom->type == GEOM_MESH) {
- flt_begin = darray_float_cdata_get(&geom->data.mesh->cdf);
- flt_end = flt_begin + darray_float_size_get(&geom->data.mesh->cdf);
- flt_found = std::lower_bound(flt_begin, flt_end, f);
- ASSERT(flt_found != flt_end);
- primitive->prim_id = (unsigned)(flt_found - flt_begin);
+ flt = darray_float_cdata_get(&geom->data.mesh->cdf);
+ sz = darray_float_size_get(&geom->data.mesh->cdf);
+ flt_found = search_lower_bound(&f, flt, sz, sizeof(*flt), cmp_float);
+ ASSERT(flt_found != NULL);
+ primitive->prim_id = (unsigned)(flt_found - flt);
} else {
FATAL("Unreachable code\n");
}
@@ -1319,7 +1380,8 @@ s3d_scene_view_get_primitive
{
struct geometry** pgeom;
struct geometry* geom;
- const struct nprims_cdf* begin, *end, *found;
+ const struct nprims_cdf* begin, *found;
+ size_t sz;
size_t nprims;
unsigned ishape;
size_t i;
@@ -1350,9 +1412,10 @@ s3d_scene_view_get_primitive
ishape = darray_nprims_cdf_cdata_get(&scnview->nprims_cdf)[0].ishape;
} else {
begin = darray_nprims_cdf_cdata_get(&scnview->nprims_cdf);
- end = begin + darray_nprims_cdf_size_get(&scnview->nprims_cdf);
- found = std::lower_bound(begin, end, i);
- ASSERT(found != end);
+ sz = darray_nprims_cdf_size_get(&scnview->nprims_cdf);
+ found = search_lower_bound
+ (&i, begin, sz, sizeof(*begin), cmp_size_t_to_nprims_cdf);
+ ASSERT(found != NULL);
ishape = found->ishape;
if(found != begin) {
ASSERT(i >= found[-1].nprims);
@@ -1378,10 +1441,11 @@ s3d_scene_view_get_primitive
} else {
begin = darray_nprims_cdf_cdata_get
(&geom->data.instance->scnview->nprims_cdf);
- end = begin + darray_nprims_cdf_size_get
+ sz = darray_nprims_cdf_size_get
(&geom->data.instance->scnview->nprims_cdf);
- found = std::lower_bound(begin, end, i);
- ASSERT(found != end);
+ found = search_lower_bound
+ (&i, begin, sz, sizeof(*begin), cmp_size_t_to_nprims_cdf);
+ ASSERT(found != NULL);
ishape = found->ishape;
if(found != begin) {
ASSERT(i >= found[-1].nprims);
@@ -1542,7 +1606,7 @@ s3d_scene_view_get_aabb(struct s3d_scene_view* scnview, float lower[3], float up
void
scene_view_destroy(struct s3d_scene_view* scnview)
{
- htable_geom_iterator it, end;
+ struct htable_geom_iterator it, end;
ASSERT(scnview && !is_list_empty(&scnview->node)/*Not in use*/);
ASSERT(scnview->mask == 0);
@@ -1557,11 +1621,10 @@ scene_view_destroy(struct s3d_scene_view* scnview)
}
/* Delete the back-end scene */
- if(scnview->rtc_scn) rtcDeleteScene(scnview->rtc_scn);
+ if(scnview->rtc_scn) rtcReleaseScene(scnview->rtc_scn);
/* Release internal data structure */
htable_geom_release(&scnview->cached_geoms);
- darray_geom_release(&scnview->embree2geoms);
darray_fltui_release(&scnview->cdf);
darray_nprims_cdf_release(&scnview->nprims_cdf);
htable_instview_release(&scnview->instviews);
@@ -1577,16 +1640,35 @@ scene_view_destroy(struct s3d_scene_view* scnview)
/* Wrapper between an Embree and a Star-3D filter function */
void
-rtc_hit_filter_wrapper(void* user_ptr, RTCRay& ray)
+rtc_hit_filter_wrapper(const struct RTCFilterFunctionNArguments* args)
{
struct s3d_hit hit;
- struct hit_filter* filter = (struct hit_filter*)user_ptr;
- struct ray_extended* ray_ex = static_cast<struct ray_extended*>(&ray);
+ struct RTCRayHit ray_hit;
+ struct intersect_context* ctx;
+ struct geometry* geom;
+ struct hit_filter* filter;
+ ASSERT(args && args->N == 1 && args->context && args->valid[0] != 0);
+
+ rtc_rayN_get_ray(args->ray, args->N, 0, &ray_hit.ray);
+ rtc_hitN_get_hit(args->hit, args->N, 0, &ray_hit.hit);
+
+ ctx = CONTAINER_OF(args->context, struct intersect_context, rtc);
+
+ geom = args->geometryUserPtr;
+ switch(geom->type) {
+ case GEOM_MESH:
+ filter = &geom->data.mesh->filter;
+ break;
+ case GEOM_SPHERE:
+ filter = &geom->data.sphere->filter;
+ break;
+ default: FATAL("Unreachable code\n"); break;
+ }
+ ASSERT(filter->func);
- hit_setup(ray_ex->scnview, &ray, &hit);
- if(filter->func(&hit, ray_ex->ws_org, ray_ex->ws_dir, ray_ex->data, filter->data)) {
- /* Discard the intersection */
- ray.geomID = RTC_INVALID_GEOMETRY_ID;
+ hit_setup(ctx->scnview, &ray_hit, &hit);
+ if(filter->func(&hit, ctx->ws_org, ctx->ws_dir, ctx->data, filter->data)) {
+ args->valid[0] = 0;
}
}
diff --git a/src/s3d_scene_view_c.h b/src/s3d_scene_view_c.h
@@ -44,6 +44,7 @@
/* Forward declarations */
struct s3d_scene_view;
+struct geometry;
/*
* The geometry pointers must be initialized to NULL in order to define
@@ -89,7 +90,6 @@ struct s3d_scene_view {
struct list_node node; /* Attachment point to the scene scene_views pool */
struct htable_geom cached_geoms; /* Cached shape geometries */
- struct darray_geom embree2geoms; /* Embree index to geometry */
struct darray_fltui cdf; /* Unormalized CDF */
struct darray_nprims_cdf nprims_cdf;
@@ -115,8 +115,7 @@ struct s3d_scene_view {
extern LOCAL_SYM void
rtc_hit_filter_wrapper
- (void* user_ptr, /* struct hit_filter* */
- RTCRay& ray);
+ (const struct RTCFilterFunctionNArguments* args);
extern LOCAL_SYM void
scene_view_destroy
@@ -128,9 +127,11 @@ scene_view_geometry_from_embree_id
const unsigned irtc)
{
struct geometry* geom;
+ RTCGeometry rtc_geom;
ASSERT(scnview && irtc != RTC_INVALID_GEOMETRY_ID);
- ASSERT(irtc < darray_geom_size_get(&scnview->embree2geoms));
- geom = darray_geom_data_get(&scnview->embree2geoms)[irtc];
+ rtc_geom = rtcGetGeometry(scnview->rtc_scn, irtc);
+ ASSERT(rtc_geom);
+ geom = rtcGetGeometryUserData(rtc_geom);
ASSERT(geom);
return geom;
}