commit 5301884f22d0bb7a8730624f6f399eeb28051748
parent f73c4fa7b729b19d304c973f22fb3bf9bfccce19
Author: Vincent Forest <vincent.forest@meso-star.com>
Date: Thu, 24 Sep 2015 18:59:31 +0200
Still another major API refactoring
The integration process does not compute the whole Gebhart Factor matrix
anymore. At each "integrate" invocation, an estimator of one line of the
matrix is returned.
Diffstat:
| M | src/sgf.c | | | 398 | ++++++++++++++++++++++++++++++------------------------------------------------- |
| M | src/sgf.h | | | 33 | ++++++++++++++++----------------- |
| M | src/test_sgf_cube.c | | | 174 | +++++++++++++++++++++++++++++++++++++++++++++++-------------------------------- |
3 files changed, 268 insertions(+), 337 deletions(-)
diff --git a/src/sgf.c b/src/sgf.c
@@ -31,11 +31,9 @@
#include <star/s3d.h>
#include <star/ssp.h>
-#include <rsys/dynamic_array_size_t.h>
-#include <rsys/hash_table.h>
+#include <rsys/dynamic_array.h>
#include <rsys/logger.h>
#include <rsys/mem_allocator.h>
-#include <rsys/mutex.h>
#include <rsys/ref_count.h>
struct accum {
@@ -43,35 +41,23 @@ struct accum {
double sqr_radiative_flux;
};
-#define HTABLE_NAME accum
-#define HTABLE_DATA struct accum
-#define HTABLE_KEY uint64_t
-#include <rsys/hash_table.h>
-
#define DARRAY_NAME accum
#define DARRAY_DATA struct accum
#include <rsys/dynamic_array.h>
-#define SPARSE
-
-#define CELL_KEY(From, To) (((uint64_t)From<<32) | (uint64_t)(To))
-
struct sgf_device {
-#ifdef SPARSE
- struct htable_accum matrix;
-#else
- struct darray_accum matrix;
- size_t nprims;
-#endif
- struct darray_size_t nsteps;
- struct mutex* mutex;
-
int verbose;
struct logger* logger;
+
struct mem_allocator* allocator;
+ ref_T ref;
+};
- struct sgf_scene_desc scn_desc;
- int is_integrating;
+struct sgf_estimator {
+ struct darray_accum buf;
+ size_t nsteps;
+
+ struct sgf_device* dev;
ref_T ref;
};
@@ -81,6 +67,39 @@ struct sgf_device {
/*******************************************************************************
* Helper function
******************************************************************************/
+static res_T
+estimator_create(struct sgf_device* dev, struct sgf_estimator** out_estimator)
+{
+ struct sgf_estimator* estimator = NULL;
+ res_T res = RES_OK;
+ ASSERT(dev && out_estimator);
+
+ estimator = MEM_CALLOC(dev->allocator, 1, sizeof(struct sgf_estimator));
+ if(!estimator) {
+ if(dev->verbose) {
+ logger_print(dev->logger, LOG_ERROR,
+"Not enough memory space: couldn't allocato the Gebhart Factor estimator.\n");
+ }
+ res = RES_MEM_ERR;
+ goto error;
+ }
+ ref_init(&estimator->ref);
+ SGF(device_ref_get(dev));
+ estimator->dev = dev;
+ darray_accum_init(dev->allocator, &estimator->buf);
+
+exit:
+ *out_estimator = estimator;
+ return res;
+
+error:
+ if(estimator) {
+ SGF(estimator_ref_put(estimator));
+ estimator = NULL;
+ }
+ goto exit;
+}
+
static INLINE char
check_scene_desc(struct sgf_scene_desc* desc)
{
@@ -88,48 +107,26 @@ check_scene_desc(struct sgf_scene_desc* desc)
return desc->get_material_property && desc->scene;
}
-static FINLINE res_T
+static FINLINE void
accum_weight
- (struct sgf_device* dev,
- const size_t ifrom,
- const size_t ito,
+ (struct sgf_estimator* estimator,
+ const size_t iface,
const double weight)
{
struct accum* acc;
- uint64_t key;
- res_T res = RES_OK;
- ASSERT(dev);
-#ifdef SPARSE
- key = CELL_KEY(ifrom, ito);
- mutex_lock(dev->mutex);
- acc = htable_accum_find(&dev->matrix, &key);
- if(acc) {
- acc->radiative_flux += weight;
- acc->sqr_radiative_flux += weight * weight;
- } else {
- struct accum tmp;
- tmp.radiative_flux = weight;
- tmp.sqr_radiative_flux = weight * weight;
- res = htable_accum_set(&dev->matrix, &key, &tmp);
- }
- mutex_unlock(dev->mutex);
-#else
- key = ifrom * dev->nprims + ito;
- acc = darray_accum_data_get(&dev->matrix) + key;
- mutex_lock(dev->mutex);
+ ASSERT(estimator);
+ acc = darray_accum_data_get(&estimator->buf) + iface;
acc->radiative_flux += weight;
acc->sqr_radiative_flux += weight * weight;
- mutex_unlock(dev->mutex);
-#endif
- return res;
}
/* Return RES_BAD_OP on numerical issue. i.e. the radiative path is skipped */
static res_T
gebhart_radiative_path
- (struct sgf_device* dev,
+ (struct sgf_estimator* estimator,
+ struct ssp_rng* rng,
const size_t primitive_id,
- struct ssp_rng* rng)
+ struct sgf_scene_desc* desc)
{
struct s3d_hit hit;
struct s3d_primitive prim_from, prim;
@@ -151,16 +148,15 @@ gebhart_radiative_path
float pos[3];
float dir[4];
float range[2];
- res_T res = RES_OK;
/* Discard faces with no emissivity */
- emissivity = dev->scn_desc.get_material_property(dev->scn_desc.material,
+ emissivity = desc->get_material_property(desc->material,
SGF_MATERIAL_EMISSIVITY, primitive_id, 0);
if(emissivity <= 0.0)
return RES_OK;
- S3D(scene_primitives_count(dev->scn_desc.scene, &nprims));
- S3D(scene_get_primitive(dev->scn_desc.scene, (unsigned)primitive_id, &prim));
+ S3D(scene_primitives_count(desc->scene, &nprims));
+ S3D(scene_get_primitive(desc->scene, (unsigned)primitive_id, &prim));
/* Compute the geometric normal of the primitive */
S3D(primitive_get_attrib(&prim, S3D_GEOMETRY_NORMAL, st, &attrib));
@@ -185,7 +181,7 @@ gebhart_radiative_path
range[0] = 1.e-6f, range[1] = FLT_MAX;
for(;;) {
- S3D(scene_trace_ray(dev->scn_desc.scene, pos, dir, range, &hit));
+ S3D(scene_trace_ray(desc->scene, pos, dir, range, &hit));
if(S3D_PRIMITIVE_EQ(&hit.prim, &prim_from)) range[0] += 1.e-6f; /* Self hit */
else break;
}
@@ -202,17 +198,16 @@ gebhart_radiative_path
prim = hit.prim;
/* Fetch material property */
- emissivity = dev->scn_desc.get_material_property(dev->scn_desc.material,
+ emissivity = desc->get_material_property(desc->material,
SGF_MATERIAL_EMISSIVITY, prim.scene_prim_id, 0);
- specularity = dev->scn_desc.get_material_property(dev->scn_desc.material,
+ specularity = desc->get_material_property(desc->material,
SGF_MATERIAL_SPECULARITY, prim.scene_prim_id, 0);
- reflectivity = dev->scn_desc.get_material_property(dev->scn_desc.material,
+ reflectivity = desc->get_material_property(desc->material,
SGF_MATERIAL_REFLECTIVITY, prim.scene_prim_id, 0);
if(transmissivity > trans_min) {
const double weight = transmissivity * emissivity;
- res = accum_weight(dev, primitive_id, prim.scene_prim_id, weight);
- if(res != RES_OK) goto error;
+ accum_weight(estimator, prim.scene_prim_id, weight);
transmissivity = transmissivity * (1.0 - emissivity);
#ifndef NDEBUG
sum_radiative_flux += weight;
@@ -221,8 +216,7 @@ gebhart_radiative_path
/* Russian roulette */
if(ssp_rng_canonical(rng) < emissivity) {
const double weight = transmissivity;
- res = accum_weight(dev, primitive_id, prim.scene_prim_id, weight);
- if(res != RES_OK) goto error;
+ accum_weight(estimator, prim.scene_prim_id, weight);
#ifndef NDEBUG
sum_radiative_flux += weight;
#endif
@@ -251,10 +245,7 @@ gebhart_radiative_path
ASSERT(eq_eps(sum_radiative_flux + sky_gebhart, 1.0, 1.e6));
#endif
-exit:
- return res;
-error:
- goto exit;
+ return RES_OK;
}
static void
@@ -263,16 +254,22 @@ device_release(ref_T* ref)
struct sgf_device* dev;
ASSERT(ref);
dev = CONTAINER_OF(ref, struct sgf_device, ref);
-#ifdef SPARSE
- htable_accum_release(&dev->matrix);
-#else
- darray_accum_release(&dev->matrix);
-#endif
- darray_size_t_release(&dev->nsteps);
- if(dev->mutex) mutex_destroy(dev->mutex);
MEM_RM(dev->allocator, dev);
}
+static void
+estimator_release(ref_T* ref)
+{
+ struct sgf_device* dev;
+ struct sgf_estimator* estimator;
+ ASSERT(ref);
+ estimator = CONTAINER_OF(ref, struct sgf_estimator, ref);
+ darray_accum_release(&estimator->buf);
+ dev = estimator->dev;
+ MEM_RM(dev->allocator, estimator);
+ SGF(device_ref_put(dev));
+}
+
/*******************************************************************************
* API functions
******************************************************************************/
@@ -309,21 +306,6 @@ sgf_device_create
dev->allocator = mem_allocator;
dev->logger = logger;
dev->verbose = verbose;
-#ifdef SPARSE
- htable_accum_init(dev->allocator, &dev->matrix);
-#else
- darray_accum_init(dev->allocator, &dev->matrix);
-#endif
- darray_size_t_init(dev->allocator, &dev->nsteps);
-
- dev->mutex = mutex_create();
- if(!dev->mutex) {
- if(verbose) {
- logger_print(logger, LOG_ERROR, "Couldn't create the SGF mutex.\n");
- }
- res = RES_MEM_ERR;
- goto error;
- }
exit:
if(out_dev) *out_dev = dev;
@@ -353,124 +335,42 @@ sgf_device_ref_put(struct sgf_device* dev)
}
res_T
-sgf_begin_integrate(struct sgf_device* dev, struct sgf_scene_desc* desc)
-{
- res_T res = RES_OK;
- size_t nprims;
-
- if(!dev || !desc || !check_scene_desc(desc)) return RES_BAD_ARG;
-
- if(dev->is_integrating) {
- if(dev->verbose) {
- logger_print(dev->logger, LOG_ERROR,
- "A Gebhart Factor integration is already running.\n");
- }
- return RES_BAD_OP;
- }
-
- res = s3d_scene_begin_session(desc->scene, S3D_TRACE|S3D_GET_PRIMITIVE);
- if(res != RES_OK) {
- if(dev->verbose) {
- logger_print(dev->logger, LOG_ERROR,
- "Couldn't begin a Star-3D sesion on the submitted scene\n");
- }
- return res;
- }
- dev->scn_desc = *desc;
- S3D(scene_primitives_count(dev->scn_desc.scene, &nprims));
-
- res = darray_size_t_resize(&dev->nsteps, nprims);
- if(res != RES_OK) {
- if(dev->verbose) {
- logger_print(dev->logger, LOG_ERROR,
-"Couldn't allocate the Gebhart Factor integration data structure.\n");
- }
- S3D(scene_end_session(desc->scene));
- return res;
- }
-
- memset(darray_size_t_data_get(&dev->nsteps), 0, nprims * sizeof(size_t));
-#ifdef SPARSE
- htable_accum_clear(&dev->matrix);
-#else
- res = darray_accum_resize(&dev->matrix, nprims*nprims);
- if(res != RES_OK) {
- if(dev->verbose) {
- logger_print(dev->logger, LOG_ERROR,
-"Couldn't allocate the Gebhart Factor integration data structure.\n");
- }
- S3D(scene_end_session(desc->scene));
- return res;
- }
- memset(darray_accum_data_get(&dev->matrix), 0, nprims*nprims*sizeof(struct accum));
- dev->nprims = nprims;
-#endif
- S3D(scene_ref_get(dev->scn_desc.scene));
- dev->is_integrating = 1;
-
- return RES_OK;
-}
-
-res_T
-sgf_end_integrate(struct sgf_device* dev)
-{
- res_T res = RES_OK;
-
- if(!dev) return RES_BAD_ARG;
- if(!dev->is_integrating) {
- if(dev->verbose) {
- logger_print(dev->logger, LOG_ERROR,
- "No Gebhart Factor integration is running\n");
- }
- return RES_BAD_OP;
- }
- res = s3d_scene_end_session(dev->scn_desc.scene);
- if(res != RES_OK) {
- if(dev->verbose) {
- logger_print(dev->logger, LOG_ERROR,
- "Couldn't stop the Star-3D session.\n");
- }
- return res;
- }
-
- S3D(scene_ref_put(dev->scn_desc.scene));
- dev->is_integrating = 0;
-
-#ifdef SPARSE
- htable_accum_clear(&dev->matrix);
-#else
- darray_accum_clear(&dev->matrix);
-#endif
- darray_size_t_clear(&dev->nsteps);
- return RES_OK;
-}
-
-res_T
sgf_integrate
(struct sgf_device* dev,
struct ssp_rng* rng,
const size_t primitive_id,
- const size_t steps_count)
+ struct sgf_scene_desc* desc,
+ const size_t steps_count,
+ struct sgf_estimator** out_estimator)
{
+ struct sgf_estimator* estimator = NULL;
+ size_t istep = 0;
size_t nprims;
- size_t i;
+ int mask;
res_T res = RES_OK;
- if(!dev || !steps_count || !rng) {
+ if(!dev || !steps_count || !rng || !desc || !check_scene_desc(desc)
+ || !out_estimator) {
res = RES_BAD_ARG;
goto error;
}
- if(!dev->is_integrating) {
+
+ res = estimator_create(dev, &estimator);
+ if(res != RES_OK) goto error;
+
+ /* Check scene active sessions */
+ S3D(scene_get_session_mask(desc->scene, &mask));
+ if(!(mask & S3D_TRACE) || !(mask & S3D_GET_PRIMITIVE)) {
if(dev->verbose) {
logger_print(dev->logger, LOG_ERROR,
-"The sgf_integrate function can be called only if sgf_begin_integrate is active"
-"on the SGF device.\n");
+ "No active trace/get_primitive session on the Star-3D scene\n");
}
- res = RES_BAD_OP;
+ res = RES_BAD_ARG;
goto error;
}
- S3D(scene_primitives_count(dev->scn_desc.scene, &nprims));
+ /* Check submitted primitive_id */
+ S3D(scene_primitives_count(desc->scene, &nprims));
if(primitive_id >= nprims) {
if(dev->verbose) {
logger_print(dev->logger, LOG_ERROR,
@@ -480,82 +380,80 @@ sgf_integrate
goto error;
}
- mutex_lock(dev->mutex);
- darray_size_t_data_get(&dev->nsteps)[primitive_id] += steps_count;
- mutex_unlock(dev->mutex);
-
- FOR_EACH(i, 0, steps_count) {
- res = gebhart_radiative_path(dev, primitive_id, rng);
- ASSERT(res == RES_OK);
+ /* Allocate internal memory space */
+ res = darray_accum_resize(&estimator->buf, nprims*desc->spectral_bands_count);
+ if(res != RES_OK) {
+ if(dev->verbose) {
+ logger_print(dev->logger, LOG_ERROR,
+"Couldn't allocate the Gebhart Factor integration data structure.\n");
+ }
+ goto error;
+ }
+ memset(darray_accum_data_get(&estimator->buf), 0,
+ darray_accum_size_get(&estimator->buf)*sizeof(struct accum));
+
+ /* Here we go ! Perform the Gebhart factor integration.
+ * TODO interate per spectral bands */
+ FOR_EACH(istep, 0, steps_count) {
+ res = gebhart_radiative_path(estimator, rng, primitive_id, desc);
+ if(res != RES_OK) goto error;
}
+ estimator->nsteps = steps_count;
exit:
+ if(out_estimator) *out_estimator = estimator;
return res;
error:
+ if(estimator) SGF(estimator_ref_put(estimator));
goto exit;
}
res_T
-sgf_get
- (struct sgf_device* dev,
- const size_t ifrom,
- const size_t ito,
- struct sgf_estimator* estimator)
+sgf_estimator_ref_get(struct sgf_estimator* estimator)
+{
+ if(!estimator) return RES_BAD_ARG;
+ ref_get(&estimator->ref);
+ return RES_OK;
+}
+
+res_T
+sgf_estimator_ref_put(struct sgf_estimator* estimator)
+{
+ if(!estimator) return RES_BAD_ARG;
+ ref_put(&estimator->ref, estimator_release);
+ return RES_OK;
+}
+
+res_T
+sgf_estimator_get_status
+ (struct sgf_estimator* estimator,
+ const size_t iprimitive,
+ struct sgf_status* status)
{
- uint64_t key;
-#ifdef SPARSE
struct accum* acc;
-#endif
- struct accum acc_cp;
- size_t nprims;
+ size_t nsteps;
- if(!dev || !estimator)
+ if(!estimator || !status)
return RES_BAD_ARG;
- if(!dev->is_integrating) {
- if(dev->verbose) {
- logger_print(dev->logger, LOG_ERROR,
-"The sgf_get function can be called only if sgf_begin_integrate is active"
-"on the SGF device.\n");
- }
- return RES_BAD_OP;
- }
-
- S3D(scene_primitives_count(dev->scn_desc.scene, &nprims));
- if(ifrom >= nprims || ito >= nprims) {
- if(dev->verbose) {
- logger_print(dev->logger, LOG_ERROR,
- "Invalid primitive index `%lu'\n",
- (unsigned long)(ifrom >= nprims ? ifrom : ito));
+ if(iprimitive >= darray_accum_size_get(&estimator->buf)) {
+ if(estimator->dev->verbose) {
+ logger_print(estimator->dev->logger, LOG_ERROR,
+ "Invalid primitive index `%lu'\n", iprimitive);
}
return RES_BAD_ARG;
}
-#ifdef SPARSE
- key = CELL_KEY(ifrom, ito);
- mutex_lock(dev->mutex);
- acc = htable_accum_find(&dev->matrix, &key);
- if(acc) acc_cp = *acc;
- mutex_unlock(dev->mutex);
- if(!acc) {
- estimator->E = estimator->V = estimator->SE = 0.0;
- estimator->nsteps = 0;
- } else
-#else
- key = ifrom * dev->nprims + ito;
- mutex_lock(dev->mutex);
- acc_cp = darray_accum_data_get(&dev->matrix)[key];
- mutex_unlock(dev->mutex);
-#endif
- {
- size_t nsteps = darray_size_t_data_get(&dev->nsteps)[ifrom];
- estimator->E = acc_cp.radiative_flux / (double)nsteps;
- estimator->V =
- acc_cp.sqr_radiative_flux / (double)nsteps
- - (acc_cp.radiative_flux * acc_cp.radiative_flux) / (double)(nsteps * nsteps);
- estimator->SE = sqrt(estimator->V / (double)nsteps);
- estimator->nsteps = nsteps;
- }
+ acc = darray_accum_data_get(&estimator->buf) + iprimitive;
+ nsteps = estimator->nsteps;
+
+ status->E = acc->radiative_flux / (double)nsteps;
+ status->V =
+ acc->sqr_radiative_flux / (double)nsteps
+ - (acc->radiative_flux * acc->radiative_flux) / (double)(nsteps * nsteps);
+ status->SE = sqrt(status->V / (double)nsteps);
+ status->nsteps = nsteps;
+
return RES_OK;
}
diff --git a/src/sgf.h b/src/sgf.h
@@ -78,7 +78,7 @@ static const struct sgf_scene_desc SGF_SCENE_DESC_NULL = {
NULL, NULL, 0, 0
};
-struct sgf_estimator {
+struct sgf_status {
double E; /* Expected value */
double V; /* Variance */
double SE; /* Standard error */
@@ -86,6 +86,7 @@ struct sgf_estimator {
};
struct sgf_device;
+struct sgf_estimator;
BEGIN_DECLS
@@ -105,29 +106,27 @@ sgf_device_ref_put
(struct sgf_device* dev);
SGF_API res_T
-sgf_begin_integrate
+sgf_integrate
(struct sgf_device* dev,
- struct sgf_scene_desc* scene);
+ struct ssp_rng* rng,
+ const size_t primitive_id,
+ struct sgf_scene_desc* desc,
+ const size_t steps_count,
+ struct sgf_estimator** out_estimator);
SGF_API res_T
-sgf_end_integrate
- (struct sgf_device* dev);
+sgf_estimator_ref_get
+ (struct sgf_estimator* estimator);
-/* Can be called only if sgf_begin_integrate is active on `dev' */
SGF_API res_T
-sgf_integrate
- (struct sgf_device* dev,
- struct ssp_rng* rng,
- const size_t primitive_id,
- const size_t steps_count);
+sgf_estimator_ref_put
+ (struct sgf_estimator* estimator);
-/* Can be called only if sgf_begin_integrate is active on `dev' */
SGF_API res_T
-sgf_get
- (struct sgf_device* dev,
- const size_t primitive_id_from,
- const size_t primitive_id_to,
- struct sgf_estimator* estimator);
+sgf_estimator_get_status
+ (struct sgf_estimator* estimator,
+ const size_t primitive_id,
+ struct sgf_status* status);
END_DECLS
diff --git a/src/test_sgf_cube.c b/src/test_sgf_cube.c
@@ -40,7 +40,7 @@
#define NSTEPS 100000L
/*
- * Matrix of Gebhart factors
+ * Analytic Gebhart Factor matrix
*
* 0.065871765 0.245053213 0.281090450 0.210375872 0.0838339939 0.113774706
* 0.183789910 0.100632183 0.291901621 0.218467251 0.0870583783 0.118150656
@@ -103,7 +103,8 @@ main(int argc, char** argv)
struct material mtr;
struct sgf_device* sgf = NULL;
struct sgf_scene_desc scn_desc = SGF_SCENE_DESC_NULL;
- struct sgf_estimator* estimators = NULL;
+ struct sgf_status* status = NULL;
+ struct sgf_estimator* estimator;
struct ssp_rng_proxy* proxy;
struct ssp_rng** rngs = NULL;
size_t iprim;
@@ -119,7 +120,7 @@ main(int argc, char** argv)
CHECK(s3d_scene_create(s3d, &scn), RES_OK);
CHECK(s3d_scene_attach_shape(scn, shape), RES_OK);
CHECK(ssp_rng_proxy_create(&allocator, &ssp_rng_threefry, nbuckets, &proxy), RES_OK);
- CHECK(sgf_device_create(NULL, NULL, 1, &sgf), RES_OK);
+ CHECK(sgf_device_create(NULL, &allocator, 1, &sgf), RES_OK);
attribs[0].type = S3D_FLOAT3;
attribs[0].usage = S3D_POSITION;
@@ -142,92 +143,127 @@ main(int argc, char** argv)
scn_desc.spectral_bands_count = 1;
scn_desc.scene = scn;
- estimators = sa_add(estimators, nprims * nprims);
+ status = sa_add(status, nprims*nprims);
rngs = sa_add(rngs, nbuckets);
- FOR_EACH(i, 0, nbuckets) {
+ FOR_EACH(i, 0, nbuckets)
CHECK(ssp_rng_proxy_create_rng(proxy, i, rngs + i), RES_OK);
- }
- CHECK(sgf_begin_integrate(NULL, NULL), RES_BAD_ARG);
- CHECK(sgf_begin_integrate(sgf, NULL), RES_BAD_ARG);
- CHECK(sgf_begin_integrate(NULL, &scn_desc), RES_BAD_ARG);
- CHECK(sgf_begin_integrate(sgf, &scn_desc), RES_OK);
- CHECK(sgf_begin_integrate(sgf, &scn_desc), RES_BAD_OP);
-
- CHECK(sgf_integrate(NULL, NULL, SIZE_MAX, 0), RES_BAD_ARG);
- CHECK(sgf_integrate(sgf, NULL, SIZE_MAX, 0), RES_BAD_ARG);
- CHECK(sgf_integrate(NULL, rngs[0], SIZE_MAX, 0), RES_BAD_ARG);
- CHECK(sgf_integrate(sgf, rngs[0], SIZE_MAX, 0), RES_BAD_ARG);
- CHECK(sgf_integrate(NULL, NULL, 0, 0), RES_BAD_ARG);
- CHECK(sgf_integrate(sgf, NULL, 0, 0), RES_BAD_ARG);
- CHECK(sgf_integrate(NULL, rngs[0], 0, 0), RES_BAD_ARG);
- CHECK(sgf_integrate(sgf, rngs[0], 0, 0), RES_BAD_ARG);
- CHECK(sgf_integrate(NULL, NULL, SIZE_MAX, NSTEPS), RES_BAD_ARG);
- CHECK(sgf_integrate(sgf, NULL, SIZE_MAX, NSTEPS), RES_BAD_ARG);
- CHECK(sgf_integrate(NULL, rngs[0], SIZE_MAX, NSTEPS), RES_BAD_ARG);
- CHECK(sgf_integrate(sgf, rngs[0], SIZE_MAX, NSTEPS), RES_BAD_ARG);
- CHECK(sgf_integrate(NULL, NULL, 0, NSTEPS), RES_BAD_ARG);
- CHECK(sgf_integrate(sgf, NULL, 0, NSTEPS), RES_BAD_ARG);
- CHECK(sgf_integrate(NULL, rngs[0], 0, NSTEPS), RES_BAD_ARG);
- CHECK(sgf_integrate(sgf, rngs[0], 0, NSTEPS), RES_OK);
-
- CHECK(sgf_get(NULL, SIZE_MAX, SIZE_MAX, NULL), RES_BAD_ARG);
- CHECK(sgf_get(sgf, SIZE_MAX, SIZE_MAX, NULL), RES_BAD_ARG);
- CHECK(sgf_get(NULL, 0, SIZE_MAX, NULL), RES_BAD_ARG);
- CHECK(sgf_get(sgf, 0, SIZE_MAX, NULL), RES_BAD_ARG);
- CHECK(sgf_get(NULL, SIZE_MAX, 0, NULL), RES_BAD_ARG);
- CHECK(sgf_get(sgf, SIZE_MAX, 0, NULL), RES_BAD_ARG);
- CHECK(sgf_get(NULL, 0, 0, NULL), RES_BAD_ARG);
- CHECK(sgf_get(sgf, 0, 0, NULL), RES_BAD_ARG);
- CHECK(sgf_get(NULL, SIZE_MAX, SIZE_MAX, estimators + 0), RES_BAD_ARG);
- CHECK(sgf_get(sgf, SIZE_MAX, SIZE_MAX, estimators + 0), RES_BAD_ARG);
- CHECK(sgf_get(NULL, 0, SIZE_MAX, estimators + 0), RES_BAD_ARG);
- CHECK(sgf_get(sgf, 0, SIZE_MAX, estimators + 0), RES_BAD_ARG);
- CHECK(sgf_get(NULL, SIZE_MAX, 0, estimators + 0), RES_BAD_ARG);
- CHECK(sgf_get(sgf, SIZE_MAX, 0, estimators + 0), RES_BAD_ARG);
- CHECK(sgf_get(NULL, 0, 0, estimators + 0), RES_BAD_ARG);
- CHECK(sgf_get(sgf, 0, 0, estimators + 0), RES_OK);
-
- CHECK(estimators[0].nsteps, NSTEPS);
-
- CHECK(sgf_end_integrate(NULL), RES_BAD_ARG);
- CHECK(sgf_end_integrate(sgf), RES_OK);
- CHECK(sgf_end_integrate(sgf), RES_BAD_OP);
-
- CHECK(sgf_integrate(sgf, rngs[0], 0, NSTEPS), RES_BAD_OP);
- CHECK(sgf_get(sgf, 0, 0, estimators + 0), RES_BAD_OP);
+ CHECK(s3d_scene_begin_session(scn, S3D_TRACE|S3D_GET_PRIMITIVE), RES_OK);
+
+ CHECK(sgf_integrate(NULL, NULL, SIZE_MAX, NULL, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, NULL, SIZE_MAX, NULL, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, rngs[0], SIZE_MAX, NULL, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, rngs[0], SIZE_MAX, NULL, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, NULL, 0, NULL, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, NULL, 0, NULL, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, rngs[0], 0, NULL, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, rngs[0], 0, NULL, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, NULL, SIZE_MAX, &scn_desc, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, NULL, SIZE_MAX, &scn_desc, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, rngs[0], SIZE_MAX, &scn_desc, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, rngs[0], SIZE_MAX, &scn_desc, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, NULL, 0, &scn_desc, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, NULL, 0, &scn_desc, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, rngs[0], 0, &scn_desc, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, rngs[0], 0, &scn_desc, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, NULL, SIZE_MAX, NULL, NSTEPS, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, NULL, SIZE_MAX, NULL, NSTEPS, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, rngs[0], SIZE_MAX, NULL, NSTEPS, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, rngs[0], SIZE_MAX, NULL, NSTEPS, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, NULL, 0, NULL, NSTEPS, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, NULL, 0, NULL, NSTEPS, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, rngs[0], 0, NULL, NSTEPS, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, rngs[0], 0, NULL, NSTEPS, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, NULL, SIZE_MAX, &scn_desc, NSTEPS, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, NULL, SIZE_MAX, &scn_desc, NSTEPS, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, rngs[0], SIZE_MAX, &scn_desc, NSTEPS, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, rngs[0], SIZE_MAX, &scn_desc, NSTEPS, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, NULL, 0, &scn_desc, NSTEPS, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, NULL, 0, &scn_desc, NSTEPS, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, rngs[0], 0, &scn_desc, NSTEPS, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, rngs[0], 0, &scn_desc, NSTEPS, NULL), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, NULL, SIZE_MAX, NULL, 0, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, NULL, SIZE_MAX, NULL, 0, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, rngs[0], SIZE_MAX, NULL, 0, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, rngs[0], SIZE_MAX, NULL, 0, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, NULL, 0, NULL, 0, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, NULL, 0, NULL, 0, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, rngs[0], 0, NULL, 0, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, rngs[0], 0, NULL, 0, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, NULL, SIZE_MAX, &scn_desc, 0, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, NULL, SIZE_MAX, &scn_desc, 0, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, rngs[0], SIZE_MAX, &scn_desc, 0, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, rngs[0], SIZE_MAX, &scn_desc, 0, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, NULL, 0, &scn_desc, 0, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, NULL, 0, &scn_desc, 0, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, rngs[0], 0, &scn_desc, 0, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, rngs[0], 0, &scn_desc, 0, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, NULL, SIZE_MAX, NULL, NSTEPS, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, NULL, SIZE_MAX, NULL, NSTEPS, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, rngs[0], SIZE_MAX, NULL, NSTEPS, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, rngs[0], SIZE_MAX, NULL, NSTEPS, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, NULL, 0, NULL, NSTEPS, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, NULL, 0, NULL, NSTEPS, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, rngs[0], 0, NULL, NSTEPS, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, rngs[0], 0, NULL, NSTEPS, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, NULL, SIZE_MAX, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, NULL, SIZE_MAX, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, rngs[0], SIZE_MAX, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, rngs[0], SIZE_MAX, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, NULL, 0, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, NULL, 0, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(NULL, rngs[0], 0, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG);
+ CHECK(sgf_integrate(sgf, rngs[0], 0, &scn_desc, NSTEPS, &estimator), RES_OK);
+
+ CHECK(sgf_estimator_get_status(NULL, SIZE_MAX, NULL), RES_BAD_ARG);
+ CHECK(sgf_estimator_get_status(estimator, SIZE_MAX, NULL), RES_BAD_ARG);
+ CHECK(sgf_estimator_get_status(NULL, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_estimator_get_status(estimator, 0, NULL), RES_BAD_ARG);
+ CHECK(sgf_estimator_get_status(NULL, SIZE_MAX, status), RES_BAD_ARG);
+ CHECK(sgf_estimator_get_status(estimator, SIZE_MAX, status), RES_BAD_ARG);
+ CHECK(sgf_estimator_get_status(NULL, 0, status), RES_BAD_ARG);
+ CHECK(sgf_estimator_get_status(estimator, 0, status), RES_OK);
+
+ CHECK(status->nsteps, NSTEPS);
+
+ CHECK(sgf_estimator_ref_get(NULL), RES_BAD_ARG);
+ CHECK(sgf_estimator_ref_get(estimator), RES_OK);
+ CHECK(sgf_estimator_ref_put(NULL), RES_BAD_ARG);
+ CHECK(sgf_estimator_ref_put(estimator), RES_OK);
+ CHECK(sgf_estimator_ref_put(estimator), RES_OK);
/* Integrate the Gebhart Factors */
- CHECK(sgf_begin_integrate(sgf, &scn_desc), RES_OK);
-
#pragma omp parallel for
for(iprim = 0; iprim < nprims; ++iprim) {
size_t iprim2;
- struct sgf_estimator* row = estimators + iprim * nprims;
+ struct sgf_status* row = status + iprim * nprims;
const int ithread = omp_get_thread_num();
+ struct sgf_estimator* est;
- CHECK(sgf_integrate(sgf, rngs[ithread], iprim, NSTEPS), RES_OK);
+ CHECK(sgf_integrate(sgf, rngs[ithread], iprim, &scn_desc, NSTEPS, &est), RES_OK);
FOR_EACH(iprim2, 0, nprims) {
- CHECK(sgf_get(sgf, iprim, iprim2, row + iprim2), RES_OK);
+ CHECK(sgf_estimator_get_status(est, iprim2, row + iprim2), RES_OK);
CHECK(row[iprim2].nsteps, NSTEPS);
}
+ CHECK(sgf_estimator_ref_put(est), RES_OK);
}
- CHECK(sgf_end_integrate(sgf), RES_OK);
+
+ CHECK(s3d_scene_end_session(scn), RES_OK);
/* Merge the radiative flux of coplanar primitives */
FOR_EACH(iprim, 0, nprims/2) {
- const struct sgf_estimator* row_src0 = estimators + iprim * 2 * nprims;
- const struct sgf_estimator* row_src1 = row_src0 + nprims;
+ const struct sgf_status* row_src0 = status + iprim * 2 * nprims;
+ const struct sgf_status* row_src1 = row_src0 + nprims;
size_t icol;
double sum = 0;
FOR_EACH(icol, 0, nprims/2) {
- const struct sgf_estimator* src0 = row_src0 + icol * 2;
- const struct sgf_estimator* src1 = src0 + 1;
- const struct sgf_estimator* src2 = row_src1 + icol * 2;
- const struct sgf_estimator* src3 = src2 + 1;
+ const struct sgf_status* src0 = row_src0 + icol * 2;
+ const struct sgf_status* src1 = src0 + 1;
+ const struct sgf_status* src2 = row_src1 + icol * 2;
+ const struct sgf_status* src3 = src2 + 1;
double E = (src0->E + src1->E + src2->E + src3->E) / 2;
sum += E;
@@ -237,17 +273,15 @@ main(int argc, char** argv)
CHECK(eq_eps(sum, 1.0, 1.e-3), 1); /* Ensure the energy conservation */
}
- sa_release(estimators);
-
CHECK(s3d_shape_ref_put(shape), RES_OK);
CHECK(s3d_scene_ref_put(scn), RES_OK);
CHECK(s3d_device_ref_put(s3d), RES_OK);
CHECK(ssp_rng_proxy_ref_put(proxy), RES_OK);
CHECK(sgf_device_ref_put(sgf), RES_OK);
- FOR_EACH(i, 0, nbuckets) {
+ FOR_EACH(i, 0, nbuckets)
CHECK(ssp_rng_ref_put(rngs[i]), RES_OK);
- }
sa_release(rngs);
+ sa_release(status);
check_memory_allocator(&allocator);
mem_shutdown_proxy_allocator(&allocator);