star-gf

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

commit 3944cc06922a700bb60ef1565315d5ad3f352a13
parent de44f0c3189f33c8115036c6f2d9fa51999e372c
Author: Vincent Forest <vincent.forest@meso-star.com>
Date:   Fri, 25 Sep 2015 16:53:51 +0200

Update the sgf_estimator_get_status_infinity API

The radiative flux that hit nothing is now estimated per primitive.

Diffstat:
Msrc/sgf.h | 1+
Msrc/sgf_estimator.c | 117+++++++++++++++++++++++++++++++++++++++++++------------------------------------
Msrc/test_sgf_cube.c | 30+++++++++++++++++++-----------
Msrc/test_sgf_tetrahedron.c | 7+++----
4 files changed, 87 insertions(+), 68 deletions(-)

diff --git a/src/sgf.h b/src/sgf.h @@ -144,6 +144,7 @@ sgf_estimator_get_status SGF_API res_T sgf_estimator_get_status_infinity (struct sgf_estimator* estimator, + const size_t primitive_id, const size_t spectral_band, struct sgf_status* status); diff --git a/src/sgf_estimator.c b/src/sgf_estimator.c @@ -37,6 +37,11 @@ #include <rsys/mem_allocator.h> #include <rsys/ref_count.h> +enum status_type { + GEBHART_FACTOR, + INFINITY_RADIATIVE_FLUX +}; + struct accum { double radiative_flux; double sqr_radiative_flux; @@ -97,7 +102,7 @@ setup_status static res_T gebhart_radiative_path (struct accum* accums, - struct accum* accum_infinity, + struct accum* accums_infinity, struct ssp_rng* rng, const size_t primitive_id, struct sgf_scene_desc* desc) @@ -122,7 +127,7 @@ gebhart_radiative_path float range[2]; /* Traced ray range */ float u, v; - ASSERT(accums && accum_infinity && rng && desc); + ASSERT(accums && accums_infinity && rng && desc); /* Discard faces with no emissivity */ emissivity = desc->get_material_property(desc->material, @@ -160,8 +165,7 @@ gebhart_radiative_path } if(S3D_HIT_NONE(&hit)) { /* The ray is outside the volume */ - accum_infinity->radiative_flux += transmissivity; - accum_infinity->sqr_radiative_flux += transmissivity * transmissivity; + accum_weight(accums_infinity, prim.scene_prim_id, transmissivity); #ifndef NDEBUG infinite_radiative_flux = transmissivity; #endif @@ -257,6 +261,52 @@ error: goto exit; } +static res_T +estimator_get_status + (struct sgf_estimator* estimator, + const enum status_type type, + const size_t iprimitive, + const size_t ispectral_band, + struct sgf_status* status) +{ + struct accum* acc; + size_t iacc; + + if(!estimator || !status) + return RES_BAD_ARG; + + if(iprimitive >= estimator->nprims) { + if(estimator->dev->verbose) { + logger_print(estimator->dev->logger, LOG_ERROR, + "Out of bound primitive index `%lu'\n", iprimitive); + } + return RES_BAD_ARG; + } + if(ispectral_band >= estimator->nbands) { + if(estimator->dev->verbose) { + logger_print(estimator->dev->logger, LOG_ERROR, + "Out of bound spectral band index `%lu'\n", ispectral_band); + } + return RES_BAD_ARG; + } + + iacc = ispectral_band * estimator->nprims + iprimitive; + switch(type) { + case GEBHART_FACTOR: + ASSERT(iacc < darray_accum_size_get(&estimator->buf)); + acc = darray_accum_data_get(&estimator->buf) + iacc; + break; + case INFINITY_RADIATIVE_FLUX: + ASSERT(iacc < darray_accum_size_get(&estimator->buf_infinity)); + acc = darray_accum_data_get(&estimator->buf_infinity) + iacc; + break; + default: FATAL("Unreachable code\n"); break; + } + setup_status(acc, estimator->nsteps, status); + return RES_OK; + +} + static void estimator_release(ref_T* ref) { @@ -334,7 +384,8 @@ sgf_integrate darray_accum_size_get(&estimator->buf)*sizeof(struct accum)); /* Allocate and init the infinite radiative flux buffer */ - res = darray_accum_resize(&estimator->buf_infinity, desc->spectral_bands_count); + res = darray_accum_resize + (&estimator->buf_infinity, nprims*desc->spectral_bands_count); if(res != RES_OK) { if(dev->verbose) { logger_print(dev->logger, LOG_ERROR, @@ -349,12 +400,12 @@ sgf_integrate FOR_EACH(ispectral_band, 0, desc->spectral_bands_count) { struct accum* accums = darray_accum_data_get(&estimator->buf) + ispectral_band * nprims; - struct accum* accum_infinity = - darray_accum_data_get(&estimator->buf_infinity) + ispectral_band; + struct accum* accums_infinity = + darray_accum_data_get(&estimator->buf_infinity) + ispectral_band * nprims; FOR_EACH(istep, 0, steps_count) { res = gebhart_radiative_path - (accums, accum_infinity, rng, primitive_id, desc); + (accums, accums_infinity, rng, primitive_id, desc); if(res != RES_OK) goto error; } } @@ -393,58 +444,18 @@ sgf_estimator_get_status const size_t ispectral_band, struct sgf_status* status) { - struct accum* acc; - size_t iacc; - - if(!estimator || !status) - return RES_BAD_ARG; - - if(iprimitive >= estimator->nprims) { - if(estimator->dev->verbose) { - logger_print(estimator->dev->logger, LOG_ERROR, - "Out of bound primitive index `%lu'\n", iprimitive); - } - return RES_BAD_ARG; - } - if(ispectral_band >= estimator->nbands) { - if(estimator->dev->verbose) { - logger_print(estimator->dev->logger, LOG_ERROR, - "Out of bound spectral band index `%lu'\n", ispectral_band); - } - return RES_BAD_ARG; - } - - iacc = ispectral_band * estimator->nprims + iprimitive; - ASSERT(iacc < darray_accum_size_get(&estimator->buf)); - acc = darray_accum_data_get(&estimator->buf) + iacc; - setup_status(acc, estimator->nsteps, status); - - return RES_OK; + return estimator_get_status + (estimator, GEBHART_FACTOR, iprimitive, ispectral_band, status); } res_T sgf_estimator_get_status_infinity (struct sgf_estimator* estimator, + const size_t iprimitive, const size_t ispectral_band, struct sgf_status* status) { - struct accum* acc; - - if(!estimator || !status) - return RES_BAD_ARG; - - if(ispectral_band >= estimator->nbands) { - if(estimator->dev->verbose) { - logger_print(estimator->dev->logger, LOG_ERROR, - "Out of bound spectral band index `%lu'\n", ispectral_band); - } - return RES_BAD_ARG; - } - - ASSERT(ispectral_band < darray_accum_size_get(&estimator->buf_infinity)); - acc = darray_accum_data_get(&estimator->buf_infinity) + ispectral_band; - setup_status(acc, estimator->nsteps, status); - - return RES_OK; + return estimator_get_status + (estimator, INFINITY_RADIATIVE_FLUX, iprimitive, ispectral_band, status); } diff --git a/src/test_sgf_cube.c b/src/test_sgf_cube.c @@ -234,14 +234,22 @@ main(int argc, char** argv) CHECK(sgf_estimator_get_status(estimator, 0, 0, status), RES_OK); CHECK(status->nsteps, NSTEPS); - CHECK(sgf_estimator_get_status_infinity(NULL, SIZE_MAX, NULL), RES_BAD_ARG); - CHECK(sgf_estimator_get_status_infinity(estimator, SIZE_MAX, NULL), RES_BAD_ARG); - CHECK(sgf_estimator_get_status_infinity(NULL, 0, NULL), RES_BAD_ARG); - CHECK(sgf_estimator_get_status_infinity(estimator, 0, NULL), RES_BAD_ARG); - CHECK(sgf_estimator_get_status_infinity(NULL, SIZE_MAX, NULL), RES_BAD_ARG); - CHECK(sgf_estimator_get_status_infinity(estimator, SIZE_MAX, status), RES_BAD_ARG); - CHECK(sgf_estimator_get_status_infinity(NULL, 0, status), RES_BAD_ARG); - CHECK(sgf_estimator_get_status_infinity(estimator, 0, status), RES_OK); + CHECK(sgf_estimator_get_status_infinity(NULL, SIZE_MAX, SIZE_MAX, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_infinity(estimator, SIZE_MAX, SIZE_MAX, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_infinity(NULL, 0, SIZE_MAX, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_infinity(estimator, 0, SIZE_MAX, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_infinity(NULL, SIZE_MAX, SIZE_MAX, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_infinity(estimator, SIZE_MAX, SIZE_MAX, status), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_infinity(NULL, 0, SIZE_MAX, status), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_infinity(estimator, 0, SIZE_MAX, status), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_infinity(NULL, SIZE_MAX, 0, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_infinity(estimator, SIZE_MAX, 0, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_infinity(NULL, 0, 0, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_infinity(estimator, 0, 0, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_infinity(NULL, SIZE_MAX, 0, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_infinity(estimator, SIZE_MAX, 0, status), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_infinity(NULL, 0, 0, status), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_infinity(estimator, 0, 0, status), RES_OK); CHECK(eq_eps(status->E, 0, status->SE), 1); CHECK(status->nsteps, NSTEPS); @@ -256,18 +264,18 @@ main(int argc, char** argv) for(iprim = 0; iprim < nprims; ++iprim) { size_t iprim2; struct sgf_status* row = status + iprim * nprims; - struct sgf_status inf; const int ithread = omp_get_thread_num(); struct sgf_estimator* est; CHECK(sgf_integrate(sgf, rngs[ithread], iprim, &scn_desc, NSTEPS, &est), RES_OK); - CHECK(sgf_estimator_get_status_infinity(est, 0, &inf), RES_OK); - CHECK(eq_eps(inf.E, 0, 3 * inf.SE), 1); FOR_EACH(iprim2, 0, nprims) { + struct sgf_status inf; CHECK(sgf_estimator_get_status(est, iprim2, 0, row + iprim2), RES_OK); CHECK(row[iprim2].nsteps, NSTEPS); + CHECK(sgf_estimator_get_status_infinity(est, iprim2, 0, &inf), RES_OK); + CHECK(eq_eps(inf.E, 0, 3 * inf.SE), 1); } CHECK(sgf_estimator_ref_put(est), RES_OK); diff --git a/src/test_sgf_tetrahedron.c b/src/test_sgf_tetrahedron.c @@ -121,19 +121,18 @@ main(int argc, char** argv) #pragma omp parallel for for(iprim = 0; iprim < nprims; ++iprim) { struct sgf_status* row = status + iprim*nprims; - struct sgf_status inf; struct sgf_estimator* estimator = NULL; size_t iprim2; const int ithread = omp_get_thread_num(); CHECK(sgf_integrate(sgf, rngs[ithread], iprim, &desc, NSTEPS, &estimator), RES_OK); - CHECK(sgf_estimator_get_status_infinity(estimator, 0, &inf), RES_OK); - CHECK(eq_eps(inf.E, 0, 3*inf.SE), 1); - FOR_EACH(iprim2, 0, nprims) { + struct sgf_status inf; CHECK(sgf_estimator_get_status(estimator, iprim2, 0, row + iprim2), RES_OK); CHECK(row[iprim2].nsteps, NSTEPS); + CHECK(sgf_estimator_get_status_infinity(estimator, iprim2, 0, &inf), RES_OK); + CHECK(eq_eps(inf.E, 0, 3*inf.SE), 1); } CHECK(sgf_estimator_ref_put(estimator), RES_OK); }