star-gf

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

commit f1763a37e41c06568935650bb4083ca47a89f5ef
parent a9633e482b551d456a5dfa959991d93b02edffb6
Author: Vincent Forest <vincent.forest@meso-star.com>
Date:   Mon, 17 Oct 2016 16:36:06 +0200

Test the medium attenuation in 3D

Diffstat:
Mcmake/CMakeLists.txt | 3++-
Msrc/test_sgf_cube.c | 174++++++++++++++++++++++++++++++++-----------------------------------------------
Asrc/test_sgf_estimator.c | 276+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Msrc/test_sgf_utils.h | 2++
4 files changed, 350 insertions(+), 105 deletions(-)

diff --git a/cmake/CMakeLists.txt b/cmake/CMakeLists.txt @@ -80,8 +80,9 @@ if(NOT NO_TEST) target_link_libraries(${_name} sgf RSys) add_test(${_name} ${_name}) endfunction(new_test) - new_test(test_sgf_device) new_test(test_sgf_cube) + new_test(test_sgf_device) + new_test(test_sgf_estimator) new_test(test_sgf_square) new_test(test_sgf_tetrahedron) rcmake_copy_runtime_libraries(test_sgf_tetrahedron) diff --git a/src/test_sgf_cube.c b/src/test_sgf_cube.c @@ -69,6 +69,16 @@ static const double emissivity[] = { 0.4, 0.4 /* Bottom */ }; +/* Emissivity used to simulate 2 infinite planes */ +static const double emissivity_inf_plane[] = { + 0, 0, /* Front */ + 0, 0, /* Left */ + 0, 0, /* Back */ + 0, 0, /* Right */ + 1, 1, /* Top */ + 1, 1 /* Bottom */ +}; + static const double specularity[] = { 0.0, 0.0, /* Front */ 0.0, 0.0, /* Left */ @@ -78,6 +88,48 @@ static const double specularity[] = { 0.0, 0.0 /* Bottom */ }; +/* Specularity used to simulate 2 infinite planes */ +static const double specularity_inf_plane[] = { + 1.0, 1.0, /* Front */ + 1.0, 1.0, /* Left */ + 1.0, 1.0, /* Back */ + 1.0, 1.0, /* Right */ + 0.0, 0.0, /* Top */ + 0.0, 0.0 /* Bottom */ +}; + +static void +check_inf_plane + (struct sgf_device* sgf, + struct ssp_rng* rng, + struct sgf_scene_desc* desc, + const double flux_bottom_top, + const double flux_bottom_medium) +{ + 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); + 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); + CHECK(sgf_estimator_ref_put(estimator), RES_OK); + + 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); + 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); +} + int main(int argc, char** argv) { @@ -91,9 +143,9 @@ main(int argc, char** argv) struct sgf_device* sgf = NULL; struct sgf_scene_desc scn_desc = SGF_SCENE_DESC_NULL; struct sgf_status* status = NULL; - struct sgf_estimator* estimator; struct ssp_rng_proxy* proxy; struct ssp_rng** rngs = NULL; + double ka[12]; int iprim; unsigned i; unsigned nbuckets; @@ -138,114 +190,14 @@ main(int argc, char** argv) 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, SIZE_MAX, NULL), RES_BAD_ARG); - CHECK(sgf_estimator_get_status(estimator, SIZE_MAX, SIZE_MAX, NULL), RES_BAD_ARG); - CHECK(sgf_estimator_get_status(NULL, 0, SIZE_MAX, NULL), RES_BAD_ARG); - CHECK(sgf_estimator_get_status(estimator, 0, SIZE_MAX, NULL), RES_BAD_ARG); - CHECK(sgf_estimator_get_status(NULL, SIZE_MAX, SIZE_MAX, status), RES_BAD_ARG); - CHECK(sgf_estimator_get_status(estimator, SIZE_MAX, SIZE_MAX, status), RES_BAD_ARG); - CHECK(sgf_estimator_get_status(NULL, 0, SIZE_MAX, status), RES_BAD_ARG); - CHECK(sgf_estimator_get_status(estimator, 0, SIZE_MAX, status), RES_BAD_ARG); - CHECK(sgf_estimator_get_status(NULL, SIZE_MAX, 0, NULL), RES_BAD_ARG); - CHECK(sgf_estimator_get_status(estimator, SIZE_MAX, 0, NULL), RES_BAD_ARG); - CHECK(sgf_estimator_get_status(NULL, 0, 0, NULL), RES_BAD_ARG); - CHECK(sgf_estimator_get_status(estimator, 0, 0, NULL), RES_BAD_ARG); - CHECK(sgf_estimator_get_status(NULL, SIZE_MAX, 0, status), RES_BAD_ARG); - CHECK(sgf_estimator_get_status(estimator, SIZE_MAX, 0, status), RES_BAD_ARG); - CHECK(sgf_estimator_get_status(NULL, 0, 0, status), RES_BAD_ARG); - 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, status), 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(eq_eps(status->E, 0, status->SE), 1); - 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 */ #pragma omp parallel for for(iprim = 0; iprim < (int)nprims; ++iprim) { size_t iprim2; struct sgf_status* row = status + (size_t)iprim * nprims; + struct sgf_estimator* est = NULL; struct sgf_status infinity; const int ithread = omp_get_thread_num(); - struct sgf_estimator* est; CHECK(sgf_integrate (sgf, rngs[ithread], (size_t)iprim, &scn_desc, NSTEPS, &est), RES_OK); @@ -261,8 +213,6 @@ main(int argc, char** argv) CHECK(sgf_estimator_ref_put(est), RES_OK); } - CHECK(s3d_scene_end_session(scn), RES_OK); - /* Merge the radiative flux of coplanar primitives */ for(iprim=0; iprim < (int)(nprims/2); ++iprim) { const struct sgf_status* row_src0 = status + (size_t)iprim * 2 * nprims; @@ -284,6 +234,22 @@ main(int argc, char** argv) CHECK(eq_eps(sum, 1.0, 1.e-3), 1); /* Ensure the energy conservation */ } + scn_desc.has_medium = 1; + mtr.emissivity = emissivity_inf_plane; + mtr.specularity = specularity_inf_plane; + mtr.absorption = ka; + + FOR_EACH(i, 0, 12) ka[i] = 0; + check_inf_plane(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); + FOR_EACH(i, 0, 12) ka[i] = 1; + check_inf_plane(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(s3d_scene_end_session(scn), RES_OK); + CHECK(s3d_shape_ref_put(shape), RES_OK); CHECK(s3d_scene_ref_put(scn), RES_OK); CHECK(s3d_device_ref_put(s3d), RES_OK); diff --git a/src/test_sgf_estimator.c b/src/test_sgf_estimator.c @@ -0,0 +1,276 @@ +/* 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/>. */ + +#include "sgf.h" +#include "test_sgf_utils.h" + +#include <rsys/stretchy_array.h> + +#include <star/s3d.h> +#include <star/ssp.h> + +#define NSTEPS 10000L + +static const float vertices[] = { + 0.f, 0.f, 0.f, + 1.f, 0.f, 0.f, + 0.f, 1.f, 0.f, + 1.f, 1.f, 0.f, + 0.f, 0.f, 1.f, + 1.f, 0.f, 1.f, + 0.f, 1.f, 1.f, + 1.f, 1.f, 1.f +}; +static const size_t nvertices = sizeof(vertices) / sizeof(float[3]); + +/* Front faces are CW. The normals point into the cube */ +static const unsigned indices[] = { + 0, 2, 1, 1, 2, 3, /* Front */ + 0, 4, 2, 2, 4, 6, /* Left */ + 4, 5, 6, 6, 5, 7, /* Back */ + 3, 7, 1, 1, 7, 5, /* Right */ + 2, 6, 3, 3, 6, 7, /* Top */ + 0, 1, 4, 4, 1, 5 /* Bottom */ +}; +static const size_t nprims = (int)(sizeof(indices) / sizeof(unsigned[3])); + +static const double emissivity[] = { + 0.6, 0.6, /* Front */ + 0.8, 0.8, /* Left */ + 0.9, 0.9, /* Back */ + 0.7, 0.7, /* Right */ + 0.3, 0.3, /* Top */ + 0.4, 0.4 /* Bottom */ +}; + +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 */ +}; + +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 sgf_device* sgf = NULL; + struct sgf_estimator* estimator = NULL; + struct sgf_scene_desc scn_desc = SGF_SCENE_DESC_NULL; + 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_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); + CHECK(sgf_estimator_get_status(NULL, 0, SIZE_MAX, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status(estimator, 0, SIZE_MAX, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status(NULL, SIZE_MAX, SIZE_MAX, &status), RES_BAD_ARG); + CHECK(sgf_estimator_get_status(estimator, SIZE_MAX, SIZE_MAX, &status), RES_BAD_ARG); + CHECK(sgf_estimator_get_status(NULL, 0, SIZE_MAX, &status), RES_BAD_ARG); + CHECK(sgf_estimator_get_status(estimator, 0, SIZE_MAX, &status), RES_BAD_ARG); + CHECK(sgf_estimator_get_status(NULL, SIZE_MAX, 0, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status(estimator, SIZE_MAX, 0, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status(NULL, 0, 0, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status(estimator, 0, 0, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status(NULL, SIZE_MAX, 0, &status), RES_BAD_ARG); + CHECK(sgf_estimator_get_status(estimator, SIZE_MAX, 0, &status), RES_BAD_ARG); + CHECK(sgf_estimator_get_status(NULL, 0, 0, &status), RES_BAD_ARG); + 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, &status), 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(eq_eps(status.E, 0, status.SE), 1); + CHECK(status.nsteps, NSTEPS); + + CHECK(sgf_estimator_get_status_medium(NULL, SIZE_MAX, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_medium(estimator, SIZE_MAX, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_medium(NULL, 0, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_medium(estimator, 0, NULL), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_medium(NULL, SIZE_MAX, &status), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_medium(estimator, SIZE_MAX,&status), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_medium(NULL, 0, &status), RES_BAD_ARG); + CHECK(sgf_estimator_get_status_medium(estimator, 0, &status), RES_OK); + CHECK(status.E, 0); + CHECK(status.V, 0); + CHECK(status.SE, 0); + 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); + + 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(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_memory_allocator(&allocator); + mem_shutdown_proxy_allocator(&allocator); + CHECK(mem_allocated_size(), 0); + return 0; +} + diff --git a/src/test_sgf_utils.h b/src/test_sgf_utils.h @@ -29,6 +29,7 @@ struct triangle_mesh { struct material { const double* emissivity; const double* specularity; + const double* absorption; }; static INLINE void @@ -71,6 +72,7 @@ get_material_property 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; default: FATAL("Unreachable code\n"); break; } }