star-gf

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

commit 174a238ae6151575f0590dc4fdf2320652112226
parent 260323f4d715bf7e811ff94dd38d2bbd7e1601f9
Author: Vincent Forest <vincent.forest@meso-star.com>
Date:   Tue, 18 Oct 2016 13:24:07 +0200

Test the medium attenuation in 2D

Diffstat:
Msrc/test_sgf_cube.c | 54++++++++++++++++++++++++++++++++++--------------------
Msrc/test_sgf_square.c | 85++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++-----------------
Msrc/test_sgf_utils.h | 10+++++-----
3 files changed, 106 insertions(+), 43 deletions(-)

diff --git a/src/test_sgf_cube.c b/src/test_sgf_cube.c @@ -70,7 +70,7 @@ static const double emissivity[] = { }; /* Emissivity used to simulate 2 infinite planes */ -static const double emissivity_inf_plane[] = { +static const double emissivity_inf_bottom_top[] = { 0, 0, /* Front */ 0, 0, /* Left */ 0, 0, /* Back */ @@ -89,7 +89,7 @@ static const double specularity[] = { }; /* Specularity used to simulate 2 infinite planes */ -static const double specularity_inf_plane[] = { +static const double specularity_inf_bottom_top[] = { 1.0, 1.0, /* Front */ 1.0, 1.0, /* Left */ 1.0, 1.0, /* Back */ @@ -98,36 +98,48 @@ static const double specularity_inf_plane[] = { 0.0, 0.0 /* Bottom */ }; +/* Check the estimation of the bottom/top & bottom/medium Gebhart factors */ static void -check_inf_plane +check_bottom_top_medium_gf (struct sgf_device* sgf, struct ssp_rng* rng, struct sgf_scene_desc* desc, - const double flux_bottom_top, - const double flux_bottom_medium) + const double gf_bottom_top, /* Ref of the bottom/top Gebhart factor */ + const double gf_bottom_medium) /* Ref of the bottom/medium Gebhart Factor */ { + /* Indices of the top/bottom primitives */ + const size_t TOP0 = 8; + const size_t TOP1 = 9; + const size_t BOTTOM0 = 10; + const size_t BOTTOM1 = 11; + struct sgf_estimator* estimator; struct sgf_status status[6]; double E, SE; - CHECK(sgf_integrate(sgf, rng, 10, desc, NSTEPS, &estimator), RES_OK); - CHECK(sgf_estimator_get_status(estimator, 8, 0, status + 0), RES_OK); - CHECK(sgf_estimator_get_status(estimator, 9, 0, status + 1), 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_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); CHECK(sgf_estimator_get_status_medium(estimator, 0, status + 3), RES_OK); CHECK(sgf_estimator_ref_put(estimator), RES_OK); - CHECK(sgf_integrate(sgf, rng, 11, desc, NSTEPS, &estimator), RES_OK); - CHECK(sgf_estimator_get_status(estimator, 8, 0, status + 4), RES_OK); - CHECK(sgf_estimator_get_status(estimator, 9, 0, status + 5), 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_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); + /* Check the Gebhart factor between the bottom and the top plane. */ E = (status[0].E + status[1].E + status[4].E + status[5].E)/2; SE = (status[0].SE + status[1].SE + status[4].SE + status[5].SE)/2; - CHECK(eq_eps(E, flux_bottom_top, SE), 1); + CHECK(eq_eps(E, gf_bottom_top, SE), 1); + + /* Check the Gebhart factor between the bottom plane and the medium */ E = (status[2].E + status[3].E)/2; SE = (status[2].SE + status[3].SE)/2; - CHECK(eq_eps(E, flux_bottom_medium, SE), 1); + CHECK(eq_eps(E, gf_bottom_medium, SE), 1); } int @@ -234,19 +246,21 @@ 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 + * 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_plane; - mtr.specularity = specularity_inf_plane; + mtr.emissivity = emissivity_inf_bottom_top; + mtr.specularity = specularity_inf_bottom_top; mtr.absorption = ka; - FOR_EACH(i, 0, 12) ka[i] = 0; - check_inf_plane(sgf, rngs[0], &scn_desc, 1, 0); + check_bottom_top_medium_gf(sgf, rngs[0], &scn_desc, 1, 0); FOR_EACH(i, 0, 12) ka[i] = 0.1; - check_inf_plane(sgf, rngs[0], &scn_desc, 0.832583, 0.167417); + check_bottom_top_medium_gf(sgf, rngs[0], &scn_desc, 0.832583, 0.167417); FOR_EACH(i, 0, 12) ka[i] = 1; - check_inf_plane(sgf, rngs[0], &scn_desc, 0.219384, 0.780616); + check_bottom_top_medium_gf(sgf, rngs[0], &scn_desc, 0.219384, 0.780616); FOR_EACH(i, 0, 12) ka[i] = 10; - check_inf_plane(sgf, rngs[0], &scn_desc, 7.0975e-6, 0.999992902); + check_bottom_top_medium_gf(sgf, rngs[0], &scn_desc, 7.0975e-6, 0.999992902); CHECK(s3d_scene_end_session(scn), RES_OK); diff --git a/src/test_sgf_square.c b/src/test_sgf_square.c @@ -47,6 +47,14 @@ static const double emissivity[] = { 1.0 /* Right */ }; +/* Emissivity used to simulate 2 infinite planes */ +static const double emissivity_inf_bottom_top[] = { + 1, /* Bottom */ + 0, /* Left */ + 1, /* Top */ + 0 /* Right */ +}; + static const double specularity[] = { 0.0, /* Bottom */ 0.0, /* Left */ @@ -54,6 +62,15 @@ static const double specularity[] = { 0.0 /* Right */ }; +/* Specularity used to simulate 2 infinite segments */ +static const double specularity_inf_bottom_top[] = { + 0.0, /* Bottom */ + 1.0, /* Left */ + 0.0, /* Top */ + 1.0 /* Right */ +}; + + static void square_get_ids(const unsigned iseg, unsigned ids[2], void* data) { @@ -74,23 +91,33 @@ square_get_pos(const unsigned ivert, float pos[2], void* data) pos[1] = vertices[id + 1]; } -static double -get_property - (void* material, - const enum sgf_material_property prop, - const size_t iprim, - const size_t ispec_band) +/* Check the estimation of the bottom/top & bottom/medium Gebhart factors */ +static void +check_bottom_top_medium_gf + (struct sgf_device* sgf, + 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 */ { - CHECK(material, NULL); - CHECK(iprim < nsegs, 1); - CHECK(ispec_band, 0); - - switch(prop) { - case SGF_MATERIAL_EMISSIVITY: return emissivity[iprim]; break; - case SGF_MATERIAL_REFLECTIVITY: return 1.0 - emissivity[iprim]; break; - case SGF_MATERIAL_SPECULARITY: return specularity[iprim]; break; - default: FATAL("Unreachable code\n"); break; - } + /* Indices of the top/bottom primitives */ + const size_t TOP = 2; + const size_t BOTTOM = 0; + + struct sgf_estimator* estimator; + struct sgf_status status[2]; + + /* Estimate the Gebhart factors for the bottom segment */ + CHECK(sgf_integrate(sgf, rng, BOTTOM, desc, NSTEPS, &estimator), RES_OK); + CHECK(sgf_estimator_get_status(estimator, TOP, 0, status + 0), RES_OK); + CHECK(sgf_estimator_get_status_medium(estimator, 0, status + 1), RES_OK); + CHECK(sgf_estimator_ref_put(estimator), RES_OK); + + /* Check the Gebhart factor between the bottom and the top plane. */ + CHECK(eq_eps(status[0].E, gf_bottom_top, status[0].SE), 1); + + /* Check the Gebhart factor between the bottom plane and the medium */ + CHECK(eq_eps(status[1].E, gf_bottom_medium, status[1].SE), 1); } int @@ -106,6 +133,8 @@ main(int argc, char** argv) struct sgf_status status; struct sgf_estimator* estimator; struct ssp_rng* rng; + struct material mtr; + double ka[4]; size_t iprim, i; (void)argc, (void)argv; @@ -124,8 +153,11 @@ main(int argc, char** argv) CHECK(s2d_line_segments_setup_indexed_vertices (shape, nsegs, square_get_ids, nverts, &vdata, 1, &shape), RES_OK); - scn_desc.get_material_property = get_property; - scn_desc.material = NULL; + mtr.emissivity = emissivity; + mtr.specularity = specularity; + + scn_desc.get_material_property = get_material_property; + scn_desc.material = &mtr; scn_desc.spectral_bands_count = 1; scn_desc.dimensionality = SGF_2D; scn_desc.geometry.scn2d = scn; @@ -158,6 +190,23 @@ main(int argc, char** argv) scn_desc.dimensionality = INT_MAX; CHECK(sgf_integrate(sgf, rng, 0, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG); + /* Check medium attenuation with 2 parallel infinite segments. To simulate + * this configuration, the top and bottom segments of the square are fully + * emissive. The other ones are fully specular and have no emissivity */ + scn_desc.has_medium = 1; + scn_desc.dimensionality = SGF_2D; + mtr.emissivity = emissivity_inf_bottom_top; + mtr.specularity = specularity_inf_bottom_top; + mtr.absorption = ka; + FOR_EACH(i, 0, 4) ka[i] = 0; + check_bottom_top_medium_gf(sgf, rng, &scn_desc, 1, 0); + FOR_EACH(i, 0, 4) ka[i] = 0.1; + check_bottom_top_medium_gf(sgf, rng, &scn_desc, 0.832583, 0.167417); + FOR_EACH(i, 0, 4) ka[i] = 1; + check_bottom_top_medium_gf(sgf, rng, &scn_desc, 0.219384, 0.780616); + FOR_EACH(i, 0, 4) ka[i] = 10; + check_bottom_top_medium_gf(sgf, rng, &scn_desc, 7.0975e-6, 0.999992902); + CHECK(s2d_scene_end_session(scn), RES_OK); CHECK(s2d_device_ref_put(s2d), RES_OK); diff --git a/src/test_sgf_utils.h b/src/test_sgf_utils.h @@ -62,17 +62,17 @@ static INLINE double get_material_property (void* material, const enum sgf_material_property prop, - const size_t itri, + const size_t iprim, const size_t ispec_band) { struct material* mtr = material; NCHECK(material, NULL); (void)ispec_band; switch(prop) { - case SGF_MATERIAL_EMISSIVITY: return mtr->emissivity[itri]; break; - case SGF_MATERIAL_REFLECTIVITY: return 1.0 - mtr->emissivity[itri]; break; - case SGF_MATERIAL_SPECULARITY: return mtr->specularity[itri]; break; - case SGF_MEDIUM_ABSORPTION: return mtr->absorption[itri]; break; + case SGF_MATERIAL_EMISSIVITY: return mtr->emissivity[iprim]; break; + case SGF_MATERIAL_REFLECTIVITY: return 1.0 - mtr->emissivity[iprim]; break; + case SGF_MATERIAL_SPECULARITY: return mtr->specularity[iprim]; break; + case SGF_MEDIUM_ABSORPTION: return mtr->absorption[iprim]; break; default: FATAL("Unreachable code\n"); break; } }