star-gf

Compute Gebhart factors
git clone git://git.meso-star.fr/star-gf.git
Log | Files | Refs | README | LICENSE

commit e1009f719d4b84131e95b5b2ac0e07cb0ac8b09f
parent 73d5158a39b2ea230665aa786ee8f35826720cd5
Author: Vincent Forest <vincent.forest@meso-star.com>
Date:   Wed, 19 Oct 2016 15:13:13 +0200

Major refactoring of the integration process

Rely on the new sgf_scene abstraction rather than the user defined
sgf_scene_desc to retrieve the scene data.

TODO: 2D integration is currently not supported.

Diffstat:
Mcmake/CMakeLists.txt | 2+-
Msrc/sgf.h | 9++++-----
Msrc/sgf_estimator.c | 90++++++++++++++++++++++++++++++-------------------------------------------------
Msrc/sgf_realisation.h | 71++++++++++++++++++++++++++---------------------------------------------
Msrc/sgf_scene.c | 29+++++------------------------
Asrc/sgf_scene_c.h | 97+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Msrc/test_sgf_cube.c | 110++++++++++++++++++++++++++++++++++++++++---------------------------------------
Msrc/test_sgf_estimator.c | 209+++++++++++++++++++++++++------------------------------------------------------
Msrc/test_sgf_scene3d.c | 12++++++------
Msrc/test_sgf_tetrahedron.c | 67+++++++++++++++++++++++++++----------------------------------------
Msrc/test_sgf_utils.h | 2+-
11 files changed, 323 insertions(+), 375 deletions(-)

diff --git a/cmake/CMakeLists.txt b/cmake/CMakeLists.txt @@ -89,7 +89,7 @@ if(NOT NO_TEST) new_test(test_sgf_device) new_test(test_sgf_estimator) new_test(test_sgf_scene3d) - new_test(test_sgf_square) + # new_test(test_sgf_square) new_test(test_sgf_tetrahedron) rcmake_copy_runtime_libraries(test_sgf_tetrahedron) diff --git a/src/sgf.h b/src/sgf.h @@ -162,11 +162,11 @@ sgf_scene_primitives_count unsigned* nprims); SGF_API res_T -sgf_scene_integration_begin +sgf_scene_begin_integration (struct sgf_scene* scn); SGF_API res_T -sgf_scene_integration_end +sgf_scene_end_integration (struct sgf_scene* scn); /******************************************************************************* @@ -178,10 +178,9 @@ sgf_scene_integration_end * by several threads. */ SGF_API res_T sgf_integrate - (struct sgf_device* dev, - struct ssp_rng* rng, + (struct sgf_scene* scn, const size_t primitive_id, - struct sgf_scene_desc* desc, + struct ssp_rng* rng, const size_t steps_count, /* #realisations */ struct sgf_estimator** out_estimator); diff --git a/src/sgf_estimator.c b/src/sgf_estimator.c @@ -33,8 +33,8 @@ #include <rsys/dynamic_array.h> /* Generate the 2D realisation function */ -#define SGF_DIMENSIONALITY 2 -#include "sgf_realisation.h" +/*#define SGF_DIMENSIONALITY 2 +#include "sgf_realisation.h"*/ /* Generate the 3D realisation function */ #define SGF_DIMENSIONALITY 3 @@ -141,38 +141,40 @@ estimator_release(ref_T* ref) ******************************************************************************/ res_T sgf_integrate - (struct sgf_device* dev, - struct ssp_rng* rng, + (struct sgf_scene* scn, const size_t iprim, - struct sgf_scene_desc* desc, + struct ssp_rng* rng, const size_t steps_count, struct sgf_estimator** out_estimator) { - #define SXD_FUNC(Func) (desc->dimensionality == SGF_2D ? S2D(Func) : S3D(Func)) - #define SXD_ENUM(Enum) (desc->dimensionality == SGF_2D ? S2D_##Enum : S3D_##Enum) + #define SXD_FUNC(Func) (dimensionality == SGF_2D ? S2D(Func) : S3D(Func)) + #define SXD_ENUM(Enum) (dimensionality == SGF_2D ? S2D_##Enum : S3D_##Enum) + const enum sgf_dimensionality dimensionality = SGF_3D; /* FIXME */ struct htable_bounce path; struct sgf_estimator* estimator = NULL; size_t istep; size_t iband; size_t nprims; - float epsilon; - float lower[3], upper[3]; void* scene; int mask; res_T res = RES_OK; gebhart_radiative_path_T gebhart_radiative_path; - if(!dev) return RES_BAD_ARG; - htable_bounce_init(dev->allocator, &path); + if(!scn) return RES_BAD_ARG; + htable_bounce_init(scn->dev->allocator, &path); - if(!steps_count || !rng || !desc || !check_scene_desc(desc) || !out_estimator) { + if(!steps_count || !rng || !scn || !out_estimator) { res = RES_BAD_ARG; goto error; } - res = estimator_create(dev, &estimator); + res = estimator_create(scn->dev, &estimator); if(res != RES_OK) goto error; +#if 1 + scene = scn->s3d_scn; + gebhart_radiative_path = gebhart_radiative_path_3d; +#else switch(desc->dimensionality) { case SGF_2D: scene = desc->geometry.scn2d; @@ -184,63 +186,41 @@ sgf_integrate break; default: FATAL("Unreachable code\n"); break; } +#endif /* Check scene active sessions */ SXD_FUNC(scene_get_session_mask(scene, &mask)); if(!(mask & SXD_ENUM(TRACE)) || !(mask & SXD_ENUM(GET_PRIMITIVE))) { - log_error(dev, - "%s: no active trace/get_primitive session on the Star-%dD scene.\n", - FUNC_NAME, desc->dimensionality == SGF_2D ? 2 : 3); - res = RES_BAD_ARG; + log_error(scn->dev, + "%s: no active integration on subimitted scene.\n", FUNC_NAME); + res = RES_BAD_OP; goto error; } /* Check submitted primitive_id */ SXD_FUNC(scene_primitives_count(scene, &nprims)); if(iprim >= nprims) { - log_error(dev, "%s: invalid primitive index `%lu'\n", + log_error(scn->dev, "%s: invalid primitive index `%lu'\n", FUNC_NAME, (unsigned long)iprim); res = RES_BAD_ARG; goto error; } - /* Compute an epsilon relatively to the scene size */ - SXD_FUNC(scene_get_aabb(scene, lower, upper)); - epsilon = upper[0] - lower[0]; - epsilon = MMIN(epsilon, upper[1] - lower[1]); - if(desc->dimensionality == SGF_3D) { - epsilon = MMIN(epsilon, upper[2] - lower[2]); - } - if(epsilon <= 0.f) { - log_error(dev, "%s: degenerated scene geometry. ", FUNC_NAME); - if(desc->dimensionality == SGF_2D) { - log_error(dev, "lower: %g, %g; upper: %g, %g\n", - SPLIT2(lower), SPLIT2(upper)); - } else { - log_error(dev, "lower: %g, %g, %g; upper: %g, %g, %g\n", - SPLIT3(lower), SPLIT3(upper)); - } - res = RES_BAD_ARG; - goto error; - } - epsilon = MMAX(epsilon*1.e-6f, FLT_MIN); - /* Allocate and init the accumulators of radiative flux */ - res = darray_accum_resize(&estimator->buf, nprims*desc->spectral_bands_count); + res = darray_accum_resize(&estimator->buf, nprims*scn->nbands); if(res != RES_OK) { - log_error(dev, "%s: couldn't allocate the Gebhart Factor result buffer.\n", + log_error(scn->dev, "%s: couldn't allocate the Gebhart Factor result buffer.\n", FUNC_NAME); goto error; } memset(darray_accum_data_get(&estimator->buf), 0, darray_accum_size_get(&estimator->buf)*sizeof(struct accum)); - if(desc->has_medium) { + if(scn->has_medium) { /* Allocate and init the accumulators of absorbed radiative flux */ - res = darray_accum_resize - (&estimator->buf_medium, desc->spectral_bands_count); + res = darray_accum_resize(&estimator->buf_medium, scn->nbands); if(res != RES_OK) { - log_error(dev, + log_error(scn->dev, "%s: couldn't allocate the accumulators of the per spectral radiative flux " "absorbed by the medium.\n", FUNC_NAME); goto error; @@ -249,10 +229,9 @@ sgf_integrate darray_accum_size_get(&estimator->buf_medium)*sizeof(struct accum)); } else { /* Allocate and init the accumulators of infinite radiative flux */ - res = darray_accum_resize - (&estimator->buf_infinity, desc->spectral_bands_count); + res = darray_accum_resize(&estimator->buf_infinity, scn->nbands); if(res != RES_OK) { - log_error(dev, + log_error(scn->dev, "%s: couldn't allocate the accumulators of the per spectral band radiative flux" " that goes to the infinite.\n", FUNC_NAME); goto error; @@ -262,21 +241,20 @@ sgf_integrate } /* Invoke the Gebhart factor integration. */ - FOR_EACH(iband, 0, desc->spectral_bands_count) { + FOR_EACH(iband, 0, scn->nbands) { struct accum* accums = NULL; struct accum* accum_infinity = NULL; struct accum* accum_medium = NULL; double ka = -1; /* Absorption coefficient of the medium */ accums = darray_accum_data_get(&estimator->buf) + iband * nprims; - if(!desc->has_medium) { + if(!scn->has_medium) { accum_infinity = darray_accum_data_get(&estimator->buf_infinity) + iband; } else { accum_medium = darray_accum_data_get(&estimator->buf_medium) + iband; - ka = desc->get_material_property - (desc->material, SGF_MEDIUM_ABSORPTION, iprim, iband); + ka = scene_get_absorption(scn, iprim, iband); if(ka < 0) { - log_error(dev, "%s: invalid absorption coefficient `%g'.\n", + log_error(scn->dev, "%s: invalid absorption coefficient `%g'.\n", FUNC_NAME, ka); res = RES_BAD_ARG; goto error; @@ -284,14 +262,14 @@ sgf_integrate } FOR_EACH(istep, 0, steps_count) { - res = gebhart_radiative_path(dev, accums, accum_infinity, accum_medium, - rng, &path, ka, iband, epsilon, iprim, desc); + res = gebhart_radiative_path(scn->dev, accums, accum_infinity, + accum_medium, rng, &path, ka, iband, iprim, scn); if(res != RES_OK) goto error; } } estimator->nsteps = steps_count; estimator->nprims = nprims; - estimator->nbands = desc->spectral_bands_count; + estimator->nbands = scn->nbands; exit: if(out_estimator) *out_estimator = estimator; diff --git a/src/sgf_realisation.h b/src/sgf_realisation.h @@ -23,6 +23,8 @@ #ifndef SGF_REALISATION_H #define SGF_REALISATION_H +#include "sgf_scene_c.h" + #include <rsys/float2.h> #include <rsys/float3.h> #include <rsys/hash_table.h> @@ -59,9 +61,8 @@ typedef res_T struct htable_bounce* path, /* Store the path */ const double absorption_coef, /* In m^-1 */ const size_t ispectral_band, - const float epsilon, const size_t primitive_id, - struct sgf_scene_desc* desc); + struct sgf_scene* scn); static FINLINE void accum_weight(struct accum* accum, const double weight) @@ -173,11 +174,11 @@ hit3d_get_normal(struct s3d_hit* hit, float normal[3]) primitive2d_get_normal(Prim, Normal) #define sXd_primitive_sample_position(Prim, Rng, Position) \ primitive2d_sample_position(Prim, Rng, Position) - #define sXd_scene_trace_ray(Scn, Org, Dir, Range, Hit) \ - S2D(scene_trace_ray_3d(Scn, Org, Dir, Range, NULL, Hit)) + #define sXd_scene_trace_ray(Scn, Org, Dir, Range, RayData, Hit) \ + S2D(scene_trace_ray_3d(Scn, Org, Dir, Range, RayData, Hit)) #define sXd_hit_get_normal(Hit, Normal) \ hit2d_get_normal(Hit, Normal) - #define sgf_scene_desc_get_sXd_scene(Desc) (Desc)->geometry.scn2d + #define sgf_scene_get_sXd_scene(Scn) (Scn)->s3d_scn #elif SGF_DIMENSIONALITY == 3 #define GEBHART_RADIATIVE_PATH gebhart_radiative_path_3d @@ -198,11 +199,11 @@ hit3d_get_normal(struct s3d_hit* hit, float normal[3]) primitive3d_get_normal(Prim, Normal) #define sXd_primitive_sample_position(Prim, Rng, Position) \ primitive3d_sample_position(Prim, Rng, Position) - #define sXd_scene_trace_ray(Scn, Org, Dir, Range, Hit) \ - S3D(scene_trace_ray(Scn, Org, Dir, Range, NULL, Hit)) + #define sXd_scene_trace_ray(Scn, Org, Dir, Range, RayData, Hit) \ + S3D(scene_trace_ray(Scn, Org, Dir, Range, RayData, Hit)) #define sXd_hit_get_normal(Hit, Normal) \ hit3d_get_normal(Hit, Normal) - #define sgf_scene_desc_get_sXd_scene(Desc) (Desc)->geometry.scn3d + #define sgf_scene_get_sXd_scene(Scn) (Scn)->s3d_scn #else #error Unexpected dimensionility #endif @@ -217,15 +218,13 @@ GEBHART_RADIATIVE_PATH struct htable_bounce* path, const double absorption_coef, /* m^-1 */ const size_t ispectral_band, - const float epsilon, const size_t primitive_id, - struct sgf_scene_desc* desc) + struct sgf_scene* scn) { - sXd_scene_T* scene; + sXd_scene_T* sXd_scn; sXd_hit_T hit; - sXd_primitive_T prim_from, prim; + sXd_primitive_T prim; struct htable_bounce_iterator it, end; - size_t nself_hits = 0; const double trans_min = 1.e-8; /* Minimum transmissivity threshold */ double proba_reflec_spec; double transmissivity, emissivity, reflectivity, specularity; @@ -241,21 +240,20 @@ GEBHART_RADIATIVE_PATH float dir[4]; /* Radiative path direction. dir[3] <=> sampled dir pdf */ float range[2]; /* Traced ray range */ res_T res = RES_OK; - ASSERT(accums && epsilon > 0.f && rng && desc); + ASSERT(accums && rng && scn); ASSERT(absorption_coef < 0 || accum_medium); ASSERT(absorption_coef >= 0 || accum_infinity); - scene = sgf_scene_desc_get_sXd_scene(desc); + sXd_scn = sgf_scene_get_sXd_scene(scn); htable_bounce_clear(path); /* Discard faces with no emissivity */ - emissivity = desc->get_material_property(desc->material, - SGF_MATERIAL_EMISSIVITY, primitive_id, 0); + emissivity = scene_get_emissivity(scn, primitive_id, ispectral_band); if(emissivity <= 0.0) return RES_OK; - /* Retrieve the scene primitive */ - sXd_scene_get_primitive(scene, (unsigned)primitive_id, &prim); + /* Retrieve the sXd_scn primitive */ + sXd_scene_get_primitive(sXd_scn, (unsigned)primitive_id, &prim); /* Get the geometric normal of the input primitive */ sXd_primitive_get_normal(&prim, normal); /* Uniformly sample prim to define the origin of the radiative path */ @@ -265,23 +263,11 @@ GEBHART_RADIATIVE_PATH transmissivity = 1.0; for(;;) { /* Here we go */ - prim_from = prim; - range[0] = epsilon, range[1] = FLT_MAX; - - nself_hits = 0; - do { /* Handle auto intersection */ - sXd_scene_trace_ray(scene, pos, dir, range, &hit); - if(!SXD_HIT_NONE(&hit)) range[0] = MMAX(range[0], hit.distance); - range[0] = nextafterf(range[0], range[1]); - } while(SXD_PRIMITIVE_EQ(&hit.prim, &prim_from) /* Self hit */ - && ++nself_hits < NSELF_HITS_MAX); /* Too many self hits */ - - if(nself_hits >= NSELF_HITS_MAX) { - log_error(dev, - "Too many self hit on a given ray. Ray origin: %g, %g, %g\n", - SPLIT3(pos)); - return RES_BAD_OP; - } + struct ray_data ray_data; + ray_data.prim_from = prim; + range[0] = FLT_MIN, range[1] = FLT_MAX; + + sXd_scene_trace_ray(sXd_scn, pos, dir, range, &ray_data, &hit); /* Handle medium absorption */ if(absorption_coef >= 0) { @@ -296,8 +282,7 @@ GEBHART_RADIATIVE_PATH /* Check the consistency of the medium description, i.e. the absorption * coefficient must be the same when it is fetched from any primitive * surrounding the current enclosure */ - ka = desc->get_material_property - (desc->material, SGF_MEDIUM_ABSORPTION, primitive_id, ispectral_band); + ka = scene_get_absorption(scn, primitive_id, ispectral_band); if(ka != absorption_coef) { log_error(dev, "Inconsistent medium description.\n"); return RES_BAD_ARG; @@ -321,13 +306,9 @@ GEBHART_RADIATIVE_PATH prim = hit.prim; /* Fetch material property */ - emissivity = desc->get_material_property(desc->material, - SGF_MATERIAL_EMISSIVITY, prim.scene_prim_id, ispectral_band); - specularity = desc->get_material_property(desc->material, - SGF_MATERIAL_SPECULARITY, prim.scene_prim_id, ispectral_band); - reflectivity = desc->get_material_property(desc->material, - SGF_MATERIAL_REFLECTIVITY, prim.scene_prim_id, ispectral_band); - + emissivity = scene_get_emissivity(scn, prim.scene_prim_id, ispectral_band); + specularity = scene_get_specularity(scn, prim.scene_prim_id, ispectral_band); + reflectivity = scene_get_reflectivity(scn, prim.scene_prim_id, ispectral_band); if(transmissivity > trans_min) { const double weight = transmissivity * emissivity; double* pweight = htable_bounce_find(path, &prim.scene_prim_id); diff --git a/src/sgf_scene.c b/src/sgf_scene.c @@ -15,31 +15,10 @@ #include "sgf.h" #include "sgf_device_c.h" - -#include <rsys/dynamic_array_double.h> -#include <rsys/ref_count.h> +#include "sgf_scene_c.h" #include <star/s3d.h> -struct ray_data { - struct s3d_primitive prim_from; -}; - -struct sgf_scene { - struct s3d_scene* s3d_scn; - struct s3d_shape* s3d_shape; - - struct darray_double abs; /* Per primitive absorption */ - struct darray_double emi; /* Per primitive emissivity */ - struct darray_double refl; /* Per primitive reflectivity */ - struct darray_double spec; /* Per primitive specularity */ - unsigned nprims; /* #primitives */ - unsigned nbands; /* #spectral bands */ - - ref_T ref; - struct sgf_device* dev; -}; - /******************************************************************************* * Helper functions ******************************************************************************/ @@ -227,6 +206,7 @@ sgf_scene3d_setup(struct sgf_scene* scn, const struct sgf_scene3d_desc* desc) scn->nbands = desc->nbands; scn->nprims = desc->nprims; + scn->has_medium = desc->get_absorption != NULL; exit: return res; @@ -238,6 +218,7 @@ error: darray_double_clear(&scn->spec); scn->nprims = 0; scn->nbands = 0; + scn->has_medium = 0; } goto exit; } @@ -251,7 +232,7 @@ sgf_scene_primitives_count(struct sgf_scene* scn, unsigned* nprims) } res_T -sgf_scene_integration_begin(struct sgf_scene* scn) +sgf_scene_begin_integration(struct sgf_scene* scn) { int session_mask; if(!scn) return RES_BAD_ARG; @@ -265,7 +246,7 @@ sgf_scene_integration_begin(struct sgf_scene* scn) } res_T -sgf_scene_integration_end(struct sgf_scene* scn) +sgf_scene_end_integration(struct sgf_scene* scn) { int session_mask; if(!scn) return RES_BAD_ARG; diff --git a/src/sgf_scene_c.h b/src/sgf_scene_c.h @@ -0,0 +1,97 @@ +/* Copyright (C) 2015-2016 EDF S.A., France (syrthes-support@edf.fr) + * + * 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/>. */ + +#ifndef SGF_SCENE_C_H +#define SGF_SCENE_C_H + +#include <rsys/dynamic_array_double.h> +#include <rsys/ref_count.h> +#include <rsys/rsys.h> + +#include <star/s3d.h> + +struct sgf_device; + +struct ray_data { + struct s3d_primitive prim_from; +}; + +struct sgf_scene { + struct s3d_scene* s3d_scn; + struct s3d_shape* s3d_shape; + + struct darray_double abs; /* Per primitive absorption */ + struct darray_double emi; /* Per primitive emissivity */ + struct darray_double refl; /* Per primitive reflectivity */ + struct darray_double spec; /* Per primitive specularity */ + unsigned nprims; /* #primitives */ + unsigned nbands; /* #spectral bands */ + int has_medium; + + ref_T ref; + struct sgf_device* dev; +}; + + +static FINLINE double +get_property__ + (const struct darray_double* property, + const size_t nprims, + const size_t nbands, + const size_t iprim, + const size_t iband) +{ + size_t i; + ASSERT(iband < nbands && iprim < nprims && property); + (void)nbands; + i = iband * nprims + iprim; + ASSERT(i < darray_double_size_get(property)); + return darray_double_cdata_get(property)[i]; +} + +static FINLINE double +scene_get_absorption + (const struct sgf_scene* scn, const size_t iprim, const size_t iband) +{ + ASSERT(scn && scn->has_medium); + return get_property__(&scn->abs, scn->nprims, scn->nbands, iprim, iband); +} + +static FINLINE double +scene_get_emissivity + (const struct sgf_scene* scn, const size_t iprim, const size_t iband) +{ + ASSERT(scn); + return get_property__(&scn->emi, scn->nprims, scn->nbands, iprim, iband); +} + +static FINLINE double +scene_get_reflectivity + (const struct sgf_scene* scn, const size_t iprim, const size_t iband) +{ + ASSERT(scn); + return get_property__(&scn->refl, scn->nprims, scn->nbands, iprim, iband); +} + +static FINLINE double +scene_get_specularity + (const struct sgf_scene* scn, const size_t iprim, const size_t iband) +{ + ASSERT(scn); + return get_property__(&scn->spec, scn->nprims, scn->nbands, iprim, iband); +} + +#endif /* SGF_SCENE_C_H */ + diff --git a/src/test_sgf_cube.c b/src/test_sgf_cube.c @@ -19,7 +19,6 @@ #include <rsys/logger.h> #include <rsys/stretchy_array.h> -#include <star/s3d.h> #include <star/ssp.h> #include <omp.h> @@ -101,9 +100,8 @@ static const double specularity_inf_bottom_top[] = { /* Check the estimation of the bottom/top & bottom/medium Gebhart factors */ static void check_bottom_top_medium_gf - (struct sgf_device* sgf, + (struct sgf_scene* scn, struct ssp_rng* rng, - struct sgf_scene_desc* desc, const double gf_bottom_top, /* Ref of the bottom/top Gebhart factor */ const double gf_bottom_medium) /* Ref of the bottom/medium Gebhart Factor */ { @@ -117,8 +115,10 @@ check_bottom_top_medium_gf struct sgf_status status[6]; double E, SE; + CHECK(sgf_scene_begin_integration(scn), RES_OK); + /* Estimate the Gebhart factors for the 1st triangle of the bottom face */ - CHECK(sgf_integrate(sgf, rng, BOTTOM0, desc, NSTEPS, &estimator), RES_OK); + CHECK(sgf_integrate(scn, BOTTOM0, rng, NSTEPS, &estimator), RES_OK); CHECK(sgf_estimator_get_status(estimator, TOP0, 0, status + 0), RES_OK); CHECK(sgf_estimator_get_status(estimator, TOP1, 0, status + 1), RES_OK); CHECK(sgf_estimator_get_status_medium(estimator, 0, status + 2), RES_OK); @@ -126,7 +126,7 @@ check_bottom_top_medium_gf CHECK(sgf_estimator_ref_put(estimator), RES_OK); /* Estimate the Gebhart factors for the 2nd triangle of the bottom face */ - CHECK(sgf_integrate(sgf, rng, BOTTOM1, desc, NSTEPS, &estimator), RES_OK); + CHECK(sgf_integrate(scn, BOTTOM1, rng, NSTEPS, &estimator), RES_OK); CHECK(sgf_estimator_get_status(estimator, TOP0, 0, status + 4), RES_OK); CHECK(sgf_estimator_get_status(estimator, TOP1, 0, status + 5), RES_OK); CHECK(sgf_estimator_ref_put(estimator), RES_OK); @@ -140,20 +140,18 @@ check_bottom_top_medium_gf E = (status[2].E + status[3].E)/2; SE = (status[2].SE + status[3].SE)/2; CHECK(eq_eps(E, gf_bottom_medium, SE), 1); + + CHECK(sgf_scene_end_integration(scn), RES_OK); } int main(int argc, char** argv) { struct mem_allocator allocator; - struct s3d_device* s3d; - struct s3d_scene* scn; - struct s3d_shape* shape; - struct s3d_vertex_data attrib; - struct triangle_mesh mesh; - struct material mtr; + struct shape_context shape; + struct sgf_scene* scn; + struct sgf_scene3d_desc desc = SGF_SCENE3D_DESC_NULL; struct sgf_device* sgf = NULL; - struct sgf_scene_desc scn_desc = SGF_SCENE_DESC_NULL; struct sgf_status* status = NULL; struct ssp_rng_proxy* proxy; struct ssp_rng** rngs = NULL; @@ -166,33 +164,28 @@ main(int argc, char** argv) mem_init_proxy_allocator(&allocator, &mem_default_allocator); nbuckets = (unsigned)omp_get_num_procs(); - CHECK(s3d_device_create(NULL, &allocator, 0, &s3d), RES_OK); - CHECK(s3d_shape_create_mesh(s3d, &shape), RES_OK); - 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, &allocator, 1, &sgf), RES_OK); - - attrib.type = S3D_FLOAT3; - attrib.usage = S3D_POSITION; - attrib.get = get_pos; - - mesh.vertices = vertices; - mesh.nvertices = nvertices; - mesh.indices = indices; - mesh.ntriangles = nprims; - - mtr.emissivity = emissivity; - mtr.specularity = specularity; - - CHECK(s3d_mesh_setup_indexed_vertices(shape, (unsigned)nprims, get_ids, - (unsigned)nvertices, &attrib, 1, &mesh), RES_OK); - - scn_desc.get_material_property = get_material_property; - scn_desc.material = &mtr; - scn_desc.spectral_bands_count = 1; - scn_desc.dimensionality = SGF_3D; - scn_desc.geometry.scn3d = scn; + CHECK(sgf_scene3d_create(sgf, &scn), RES_OK); + + shape.vertices = vertices; + shape.nvertices = nvertices; + shape.indices = indices; + shape.nprimitives = nprims; + shape.emissivity = emissivity; + shape.specularity = specularity; + + desc.get_position = get_position; + desc.get_indices = get_indices; + desc.get_emissivity = get_emissivity; + desc.get_specularity = get_specularity; + desc.get_reflectivity = get_reflectivity; + desc.nprims = (unsigned)nprims; + desc.nverts = (unsigned)nvertices; + desc.nbands = 1; + desc.context = &shape; + + CHECK(sgf_scene3d_setup(scn, &desc), RES_OK); status = sa_add(status, nprims*nprims); rngs = sa_add(rngs, nbuckets); @@ -200,7 +193,7 @@ main(int argc, char** argv) FOR_EACH(i, 0, nbuckets) CHECK(ssp_rng_proxy_create_rng(proxy, i, rngs + i), RES_OK); - CHECK(s3d_scene_begin_session(scn, S3D_TRACE|S3D_GET_PRIMITIVE), RES_OK); + CHECK(sgf_scene_begin_integration(scn), RES_OK); /* Integrate the Gebhart Factors */ #pragma omp parallel for @@ -212,7 +205,7 @@ main(int argc, char** argv) const int ithread = omp_get_thread_num(); CHECK(sgf_integrate - (sgf, rngs[ithread], (size_t)iprim, &scn_desc, NSTEPS, &est), RES_OK); + (scn, (size_t)iprim, rngs[ithread], NSTEPS, &est), RES_OK); FOR_EACH(iprim2, 0, nprims) { CHECK(sgf_estimator_get_status(est, iprim2, 0, row + iprim2), RES_OK); @@ -246,29 +239,38 @@ main(int argc, char** argv) CHECK(eq_eps(sum, 1.0, 1.e-3), 1); /* Ensure the energy conservation */ } - /* Check medium attenuation with 2 parallel infinite planes. To simulate this + CHECK(sgf_scene_end_integration(scn), RES_OK); + + /* + * Check medium attenuation with 2 parallel infinite planes. To simulate this * configuration, the top and bottom faces of the cube are fully emissive. - * The other ones are fully specular and have no emissivity */ - scn_desc.has_medium = 1; - mtr.emissivity = emissivity_inf_bottom_top; - mtr.specularity = specularity_inf_bottom_top; - mtr.absorption = ka; + * The other ones are fully specular and have no emissivity + */ + + shape.absorption = ka; + shape.emissivity = emissivity_inf_bottom_top; + shape.specularity = specularity_inf_bottom_top; + desc.get_absorption = get_absorption; + FOR_EACH(i, 0, 12) ka[i] = 0; - check_bottom_top_medium_gf(sgf, rngs[0], &scn_desc, 1, 0); + CHECK(sgf_scene3d_setup(scn, &desc), RES_OK); + check_bottom_top_medium_gf(scn, rngs[0], 1, 0); + FOR_EACH(i, 0, 12) ka[i] = 0.1; - check_bottom_top_medium_gf(sgf, rngs[0], &scn_desc, 0.832583, 0.167417); + CHECK(sgf_scene3d_setup(scn, &desc), RES_OK); + check_bottom_top_medium_gf(scn, rngs[0], 0.832583, 0.167417); + FOR_EACH(i, 0, 12) ka[i] = 1; - check_bottom_top_medium_gf(sgf, rngs[0], &scn_desc, 0.219384, 0.780616); - FOR_EACH(i, 0, 12) ka[i] = 10; - check_bottom_top_medium_gf(sgf, rngs[0], &scn_desc, 7.0975e-6, 0.999992902); + CHECK(sgf_scene3d_setup(scn, &desc), RES_OK); + check_bottom_top_medium_gf(scn, rngs[0], 0.219384, 0.780616); - CHECK(s3d_scene_end_session(scn), RES_OK); + FOR_EACH(i, 0, 12) ka[i] = 10; + CHECK(sgf_scene3d_setup(scn, &desc), RES_OK); + check_bottom_top_medium_gf(scn, rngs[0], 7.0975e-6, 0.999992902); - 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); + CHECK(sgf_scene_ref_put(scn), RES_OK); FOR_EACH(i, 0, nbuckets) CHECK(ssp_rng_ref_put(rngs[i]), RES_OK); sa_release(rngs); diff --git a/src/test_sgf_estimator.c b/src/test_sgf_estimator.c @@ -56,137 +56,85 @@ static const double emissivity[] = { }; static const double specularity[] = { - 1.0, 1.0, /* Front */ - 1.0, 1.0, /* Left */ - 1.0, 1.0, /* Back */ - 1.0, 1.0, /* Right */ - 1.0, 1.0, /* Top */ - 1.0, 1.0 /* Bottom */ -}; - -static const double ka[] = { /* Absorption coefficient of the medium */ - 1, 1, /* Front */ - 1, 1, /* Left */ - 1, 1, /* Back */ - 1, 1, /* Right */ - 1, 1, /* Top */ - 1, 1 /* Bottom */ + 0.0, 0.0, /* Front */ + 0.0, 0.0, /* Left */ + 0.0, 0.0, /* Back */ + 0.0, 0.0, /* Right */ + 0.0, 0.0, /* Top */ + 0.0, 0.0 /* Bottom */ }; int main(int argc, char** argv) { struct mem_allocator allocator; - struct s3d_device* s3d = NULL; - struct s3d_scene* scn = NULL; - struct s3d_shape* shape = NULL; - struct s3d_vertex_data attrib = S3D_VERTEX_DATA_NULL; + struct shape_context shape; struct sgf_device* sgf = NULL; struct sgf_estimator* estimator = NULL; - struct sgf_scene_desc scn_desc = SGF_SCENE_DESC_NULL; + struct sgf_scene3d_desc desc = SGF_SCENE3D_DESC_NULL; + struct sgf_scene* scn; struct sgf_status status; struct ssp_rng* rng = NULL; - struct triangle_mesh mesh; - struct material mtr; (void)argc, (void)argv; mem_init_proxy_allocator(&allocator, &mem_default_allocator); - CHECK(s3d_device_create(NULL, &allocator, 1, &s3d), RES_OK); - CHECK(s3d_shape_create_mesh(s3d, &shape), RES_OK); - CHECK(s3d_scene_create(s3d, &scn), RES_OK); - CHECK(s3d_scene_attach_shape(scn, shape), RES_OK); CHECK(ssp_rng_create(&allocator, &ssp_rng_threefry, &rng), RES_OK); CHECK(sgf_device_create(NULL, &allocator, 1, &sgf), RES_OK); - - attrib.type = S3D_FLOAT3; - attrib.usage = S3D_POSITION; - attrib.get = get_pos; - - mesh.vertices = vertices; - mesh.nvertices = nvertices; - mesh.indices = indices; - mesh.ntriangles = nprims; - - mtr.emissivity = emissivity; - mtr.specularity = specularity; - mtr.absorption = NULL; - - CHECK(s3d_mesh_setup_indexed_vertices(shape, (unsigned)nprims, get_ids, - (unsigned)nvertices, &attrib, 1, &mesh), RES_OK); - - scn_desc.get_material_property = get_material_property; - scn_desc.material = &mtr; - scn_desc.spectral_bands_count = 1; - scn_desc.dimensionality = SGF_3D; - scn_desc.geometry.scn3d = scn; - - 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, rng, SIZE_MAX, NULL, 0, NULL), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 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, rng, 0, NULL, 0, NULL), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 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, rng, SIZE_MAX, &scn_desc, 0, NULL), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 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, rng, 0, &scn_desc, 0, NULL), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 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, rng, SIZE_MAX, NULL, NSTEPS, NULL), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 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, rng, 0, NULL, NSTEPS, NULL), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 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, rng, SIZE_MAX, &scn_desc, NSTEPS, NULL), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 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, rng, 0, &scn_desc, NSTEPS, NULL), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 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, rng, SIZE_MAX, NULL, 0, &estimator), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 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, rng, 0, NULL, 0, &estimator), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 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, rng, SIZE_MAX, &scn_desc, 0, &estimator), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 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, rng, 0, &scn_desc, 0, &estimator), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 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, rng, SIZE_MAX, NULL, NSTEPS, &estimator), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 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, rng, 0, NULL, NSTEPS, &estimator), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 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, rng, SIZE_MAX, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 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, rng, 0, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 0, &scn_desc, NSTEPS, &estimator), RES_OK); + CHECK(sgf_scene3d_create(sgf, &scn), RES_OK); + + shape.vertices = vertices; + shape.nvertices = nvertices; + shape.indices = indices; + shape.nprimitives = nprims; + shape.emissivity = emissivity; + shape.specularity = specularity; + + desc.get_position = get_position; + desc.get_indices = get_indices; + desc.get_emissivity = get_emissivity; + desc.get_specularity = get_specularity; + desc.get_reflectivity = get_reflectivity; + desc.nprims = (unsigned)nprims; + desc.nverts = (unsigned)nvertices; + desc.nbands = 1; + desc.context = &shape; + + CHECK(sgf_scene3d_setup(scn, &desc), RES_OK); + CHECK(sgf_scene_begin_integration(scn), RES_OK); + + CHECK(sgf_integrate(NULL, SIZE_MAX, NULL, 0, NULL), RES_BAD_ARG); + CHECK(sgf_integrate(scn, SIZE_MAX, NULL, 0, NULL), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, 0, NULL, 0, NULL), RES_BAD_ARG); + CHECK(sgf_integrate(scn, 0, NULL, 0, NULL), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, SIZE_MAX, rng, 0, NULL), RES_BAD_ARG); + CHECK(sgf_integrate(scn, SIZE_MAX, rng, 0, NULL), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, 0, rng, 0, NULL), RES_BAD_ARG); + CHECK(sgf_integrate(scn, 0, rng, 0, NULL), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, SIZE_MAX, NULL, NSTEPS, NULL), RES_BAD_ARG); + CHECK(sgf_integrate(scn, SIZE_MAX, NULL, NSTEPS, NULL), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, 0, NULL, NSTEPS, NULL), RES_BAD_ARG); + CHECK(sgf_integrate(scn, 0, NULL, NSTEPS, NULL), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, SIZE_MAX, rng, NSTEPS, NULL), RES_BAD_ARG); + CHECK(sgf_integrate(scn, SIZE_MAX, rng, NSTEPS, NULL), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, 0, rng, NSTEPS, NULL), RES_BAD_ARG); + CHECK(sgf_integrate(scn, 0, rng, NSTEPS, NULL), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, SIZE_MAX, NULL, 0, &estimator), RES_BAD_ARG); + CHECK(sgf_integrate(scn, SIZE_MAX, NULL, 0, &estimator), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, 0, NULL, 0, &estimator), RES_BAD_ARG); + CHECK(sgf_integrate(scn, 0, NULL, 0, &estimator), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, SIZE_MAX, rng, 0, &estimator), RES_BAD_ARG); + CHECK(sgf_integrate(scn, SIZE_MAX, rng, 0, &estimator), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, 0, rng, 0, &estimator), RES_BAD_ARG); + CHECK(sgf_integrate(scn, 0, rng, 0, &estimator), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, SIZE_MAX, NULL, NSTEPS, &estimator), RES_BAD_ARG); + CHECK(sgf_integrate(scn, SIZE_MAX, NULL, NSTEPS, &estimator), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, 0, NULL, NSTEPS, &estimator), RES_BAD_ARG); + CHECK(sgf_integrate(scn, 0, NULL, NSTEPS, &estimator), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, SIZE_MAX, rng, NSTEPS, &estimator), RES_BAD_ARG); + CHECK(sgf_integrate(scn, SIZE_MAX, rng, NSTEPS, &estimator), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, 0, rng, NSTEPS, &estimator), RES_BAD_ARG); + CHECK(sgf_integrate(scn, 0, rng, NSTEPS, &estimator), RES_OK); CHECK(sgf_estimator_get_status(NULL, SIZE_MAX, SIZE_MAX, NULL), RES_BAD_ARG); CHECK(sgf_estimator_get_status(estimator, SIZE_MAX, SIZE_MAX, NULL), RES_BAD_ARG); @@ -236,37 +184,12 @@ main(int argc, char** argv) CHECK(sgf_estimator_ref_put(estimator), RES_OK); CHECK(sgf_estimator_ref_put(estimator), RES_OK); - scn_desc.get_material_property = NULL; - CHECK(sgf_integrate(sgf, rng, 0, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG); - scn_desc.get_material_property = get_material_property; - scn_desc.spectral_bands_count = 0; - CHECK(sgf_integrate(sgf, rng, 0, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG); - scn_desc.spectral_bands_count = 1; - scn_desc.geometry.scn3d = NULL; - CHECK(sgf_integrate(sgf, rng, 0, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG); - scn_desc.geometry.scn3d = scn; - - CHECK(s3d_scene_end_session(scn), RES_OK); - CHECK(sgf_integrate(sgf, rng, 0, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG); - CHECK(s3d_scene_begin_session(scn, S3D_TRACE|S3D_GET_PRIMITIVE), RES_OK); - - mtr.absorption = ka; - scn_desc.has_medium = 1; - CHECK(sgf_integrate(sgf, rng, 0, &scn_desc, NSTEPS, &estimator), RES_OK); - CHECK(sgf_estimator_get_status(estimator, 0, 0, &status), RES_OK); - CHECK(sgf_estimator_get_status_infinity(estimator, 0, &status), RES_OK); - CHECK(eq_eps(status.E, 0, status.SE), 1); - CHECK(sgf_estimator_get_status_medium(estimator, 0, &status), RES_OK); - CHECK(eq_eps(status.E, 0, status.SE), 0); - - CHECK(sgf_estimator_ref_put(estimator), RES_OK); - CHECK(s3d_scene_end_session(scn), RES_OK); + CHECK(sgf_scene_end_integration(scn), RES_OK); + CHECK(sgf_integrate(scn, 0, rng, NSTEPS, &estimator), RES_BAD_OP); - CHECK(s3d_device_ref_put(s3d), RES_OK); - CHECK(s3d_shape_ref_put(shape), RES_OK); - CHECK(s3d_scene_ref_put(scn), RES_OK); CHECK(ssp_rng_ref_put(rng), RES_OK); CHECK(sgf_device_ref_put(sgf), RES_OK); + CHECK(sgf_scene_ref_put(scn), RES_OK); check_memory_allocator(&allocator); mem_shutdown_proxy_allocator(&allocator); diff --git a/src/test_sgf_scene3d.c b/src/test_sgf_scene3d.c @@ -117,13 +117,13 @@ main(int argc, char** argv) ctx.absorption = absorption_bad; CHECK(sgf_scene3d_setup(scn, &desc), RES_BAD_ARG); - CHECK(sgf_scene_integration_begin(NULL), RES_BAD_ARG); - CHECK(sgf_scene_integration_begin(scn), RES_OK); - CHECK(sgf_scene_integration_begin(scn), RES_BAD_OP); + CHECK(sgf_scene_begin_integration(NULL), RES_BAD_ARG); + CHECK(sgf_scene_begin_integration(scn), RES_OK); + CHECK(sgf_scene_begin_integration(scn), RES_BAD_OP); - CHECK(sgf_scene_integration_end(NULL), RES_BAD_ARG); - CHECK(sgf_scene_integration_end(scn), RES_OK); - CHECK(sgf_scene_integration_end(scn), RES_BAD_OP); + CHECK(sgf_scene_end_integration(NULL), RES_BAD_ARG); + CHECK(sgf_scene_end_integration(scn), RES_OK); + CHECK(sgf_scene_end_integration(scn), RES_BAD_OP); CHECK(sgf_scene_ref_get(NULL), RES_BAD_ARG); CHECK(sgf_scene_ref_get(scn), RES_OK); diff --git a/src/test_sgf_tetrahedron.c b/src/test_sgf_tetrahedron.c @@ -19,7 +19,6 @@ #include <rsys/logger.h> #include <rsys/stretchy_array.h> -#include <star/s3d.h> #include <star/ssp.h> #include <omp.h> @@ -50,14 +49,10 @@ int main(int argc, char** argv) { struct mem_allocator allocator; - struct s3d_device* s3d; - struct s3d_scene* scn; - struct s3d_shape* shape; - struct s3d_vertex_data attrib; - struct triangle_mesh mesh; - struct material mtr; + struct shape_context shape; struct sgf_device* sgf = NULL; - struct sgf_scene_desc desc = SGF_SCENE_DESC_NULL; + struct sgf_scene3d_desc desc = SGF_SCENE3D_DESC_NULL; + struct sgf_scene* scn = NULL; struct sgf_status* status = NULL; struct ssp_rng_proxy* proxy = NULL; struct ssp_rng** rngs = NULL; @@ -68,34 +63,28 @@ main(int argc, char** argv) mem_init_proxy_allocator(&allocator, &mem_default_allocator); nbuckets = (unsigned)omp_get_num_procs(); - - CHECK(s3d_device_create(NULL, &allocator, 0, &s3d), RES_OK); - CHECK(s3d_shape_create_mesh(s3d, &shape), RES_OK); - 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, &allocator, 1, &sgf), RES_OK); - - attrib.type = S3D_FLOAT3; - attrib.usage = S3D_POSITION; - attrib.get = get_pos; - - mesh.vertices = vertices; - mesh.nvertices = nvertices; - mesh.indices = indices; - mesh.ntriangles = nprims; - - mtr.emissivity = emissivity; - mtr.specularity = specularity; - - CHECK(s3d_mesh_setup_indexed_vertices(shape, (unsigned)nprims, get_ids, - (unsigned)nvertices, &attrib, 1, &mesh), RES_OK); - - desc.get_material_property = get_material_property; - desc.material = &mtr; - desc.spectral_bands_count = 1; - desc.dimensionality = SGF_3D; - desc.geometry.scn3d = scn; + CHECK(sgf_scene3d_create(sgf, &scn), RES_OK); + + shape.vertices = vertices; + shape.nvertices = nvertices; + shape.indices = indices; + shape.nprimitives = nprims; + shape.emissivity = emissivity; + shape.specularity = specularity; + + desc.get_position = get_position; + desc.get_indices = get_indices; + desc.get_emissivity = get_emissivity; + desc.get_specularity = get_specularity; + desc.get_reflectivity = get_reflectivity; + desc.nprims = (unsigned)nprims; + desc.nverts = (unsigned)nvertices; + desc.nbands = 1; + desc.context = &shape; + + CHECK(sgf_scene3d_setup(scn, &desc), RES_OK); status = sa_add(status, nprims*nprims); rngs = sa_add(rngs, nbuckets); @@ -103,7 +92,7 @@ main(int argc, char** argv) FOR_EACH(i, 0, nbuckets) CHECK(ssp_rng_proxy_create_rng(proxy, i, rngs + i), RES_OK); - CHECK(s3d_scene_begin_session(scn, S3D_TRACE|S3D_GET_PRIMITIVE), RES_OK); + CHECK(sgf_scene_begin_integration(scn), RES_OK); /* Gebhart Factor integration */ #pragma omp parallel for @@ -115,7 +104,7 @@ main(int argc, char** argv) const int ithread = omp_get_thread_num(); CHECK(sgf_integrate - (sgf, rngs[ithread], (size_t)iprim, &desc, NSTEPS, &estimator), RES_OK); + (scn, (size_t)iprim, rngs[ithread], NSTEPS, &estimator), RES_OK); FOR_EACH(iprim2, 0, nprims) { CHECK(sgf_estimator_get_status(estimator, iprim2, 0, row + iprim2), RES_OK); @@ -127,7 +116,7 @@ main(int argc, char** argv) CHECK(sgf_estimator_ref_put(estimator), RES_OK); } - CHECK(s3d_scene_end_session(scn), RES_OK); + CHECK(sgf_scene_end_integration(scn), RES_OK); /* Expected value */ for(iprim=0; iprim < (int)nprims; ++iprim) { @@ -154,11 +143,9 @@ main(int argc, char** argv) printf("\n"); } - 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); + CHECK(sgf_scene_ref_put(scn), RES_OK); FOR_EACH(i, 0, nbuckets) CHECK(ssp_rng_ref_put(rngs[i]), RES_OK); diff --git a/src/test_sgf_utils.h b/src/test_sgf_utils.h @@ -143,7 +143,7 @@ get_specularity(const unsigned iprim, const unsigned iband, void* ctx) static INLINE double get_reflectivity(const unsigned iprim, const unsigned iband, void* ctx) { - return 1 - get_specularity(iprim, iband, ctx); + return 1 - get_emissivity(iprim, iband, ctx); } static INLINE void