star-gf

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

commit 389933d613bfee1b7e4e0db805c6e403d6d607cb
parent 10fbe01e1892bb1ad8c9f0f96caed8a11930ffee
Author: Vincent Forest <vincent.forest@meso-star.com>
Date:   Thu,  7 Jul 2016 10:55:46 +0200

Test the Gebhart Factor integration in 2D

Diffstat:
Mcmake/CMakeLists.txt | 1+
Msrc/sgf_estimator.c | 21+++++++++++++++++----
Asrc/test_sgf_square.c | 173+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 files changed, 191 insertions(+), 4 deletions(-)

diff --git a/cmake/CMakeLists.txt b/cmake/CMakeLists.txt @@ -82,6 +82,7 @@ if(NOT NO_TEST) endfunction(new_test) new_test(test_sgf_device) new_test(test_sgf_cube) + new_test(test_sgf_square) new_test(test_sgf_tetrahedron) rcmake_copy_runtime_libraries(test_sgf_tetrahedron) diff --git a/src/sgf_estimator.c b/src/sgf_estimator.c @@ -255,12 +255,25 @@ sgf_integrate SXD_FUNC(scene_get_aabb(scene, lower, upper)); epsilon = upper[0] - lower[0]; epsilon = MMIN(epsilon, upper[1] - lower[1]); - epsilon = MMIN(epsilon, upper[2] - lower[2]); + if(desc->dimensionality == SGF_3D) { + epsilon = MMIN(epsilon, upper[2] - lower[2]); + } if(epsilon <= 0.f) { if(dev->verbose) { - logger_print(dev->logger, LOG_ERROR, - "Degenerated scene geometry. lower: %g, %g, %g; upper: %g %g %g\n", - SPLIT3(lower), SPLIT3(upper)); + logger_print(dev->logger, LOG_ERROR,"Degenerated scene geometry. "); + switch(desc->dimensionality) { + case SGF_2D: + logger_print(dev->logger, LOG_ERROR, + "lower: %g, %g; upper: %g, %g\n", + SPLIT2(lower), SPLIT2(upper)); + break; + case SGF_3D: + logger_print(dev->logger, LOG_ERROR, + "lower: %g, %g, %g; upper: %g, %g, %g\n", + SPLIT3(lower), SPLIT3(upper)); + break; + default: FATAL("Unreachable code\n"); break; + } } res = RES_BAD_ARG; goto error; diff --git a/src/test_sgf_square.c b/src/test_sgf_square.c @@ -0,0 +1,173 @@ +/* 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/logger.h> +#include <star/s2d.h> +#include <star/ssp.h> + +#include <limits.h> + +#define NSTEPS 100000 + +static const float vertices[] = { + 1.f, 0.f, + 0.f, 0.f, + 0.f, 1.f, + 1.f, 1.f +}; +const unsigned nverts = sizeof(vertices)/sizeof(float[2]); + +const unsigned indices[] = { + 0, 1, /* Bottom */ + 1, 2, /* Left */ + 2, 3, /* Top */ + 3, 0 /* Right */ +}; +const unsigned nsegs = sizeof(indices)/sizeof(unsigned[2]); + +static const double emissivity[] = { + 1.0, /* Bottom */ + 1.0, /* Left */ + 1.0, /* Top */ + 1.0 /* Right */ +}; + +static const double specularity[] = { + 0.0, /* Bottom */ + 0.0, /* Left */ + 0.0, /* Top */ + 0.0 /* Right */ +}; + +static void +square_get_ids(const unsigned iseg, unsigned ids[2], void* data) +{ + const unsigned id = iseg * 2/*#coords per segment*/; + (void)data; + NCHECK(ids, NULL); + ids[0] = indices[id + 0]; + ids[1] = indices[id + 1]; +} + +static void +square_get_pos(const unsigned ivert, float pos[2], void* data) +{ + const unsigned id = ivert * 2/*#vertices per segment*/; + (void)data; + NCHECK(pos, NULL); + pos[0] = vertices[id + 0]; + 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(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; + } +} + +int +main(int argc, char** argv) +{ + struct mem_allocator allocator; + struct s2d_device* s2d; + struct s2d_scene* scn; + struct s2d_shape* shape; + struct s2d_vertex_data vdata; + struct sgf_device* sgf; + struct sgf_scene_desc scn_desc = SGF_SCENE_DESC_NULL; + struct sgf_status status; + struct sgf_estimator* estimator; + struct ssp_rng* rng; + size_t iprim, i; + (void)argc, (void)argv; + + mem_init_proxy_allocator(&allocator, &mem_default_allocator); + + CHECK(s2d_device_create(NULL, &allocator, 1, &s2d), RES_OK); + CHECK(s2d_shape_create_line_segments(s2d, &shape), RES_OK); + CHECK(s2d_scene_create(s2d, &scn), RES_OK); + CHECK(s2d_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); + + vdata.type = S2D_FLOAT2; + vdata.usage = S2D_POSITION; + vdata.get = square_get_pos; + 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; + scn_desc.spectral_bands_count = 1; + scn_desc.dimensionality = SGF_2D; + scn_desc.geometry.scn2d = scn; + + CHECK(sgf_integrate(sgf, rng, 0, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG); + + CHECK(s2d_scene_begin_session(scn, S2D_TRACE|S2D_GET_PRIMITIVE), RES_OK); + + FOR_EACH(iprim, 0, 3) { + CHECK(sgf_integrate(sgf, rng, iprim, &scn_desc, NSTEPS, &estimator), RES_OK); + FOR_EACH(i, 0, 3) { + CHECK(sgf_estimator_get_status(estimator, i, 0, &status), RES_OK); + + if(i == iprim) { + CHECK(eq_eps(status.E, 0, status.SE), 1); + } else if(i == iprim + 2 || i == iprim - 2) { /* parallel faces */ + CHECK(eq_eps(status.E, sqrt(2) - 1, 2*status.SE), 1); + } else { /* Orthogonal faces */ + CHECK(iprim == i+1 || iprim == i-1, 1); + CHECK(eq_eps(status.E, 1 - sqrt(2)/2, 2*status.SE), 1); + } + CHECK(sgf_estimator_get_status_infinity(estimator, i, 0, &status), RES_OK); + CHECK(eq_eps(status.E, 0, status.SE), 1); + } + CHECK(sgf_estimator_ref_put(estimator), RES_OK); + } + + CHECK(sgf_integrate(sgf, rng, 4, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG); + scn_desc.dimensionality = INT_MAX; + CHECK(sgf_integrate(sgf, rng, 0, &scn_desc, NSTEPS, &estimator), RES_BAD_ARG); + + CHECK(s2d_scene_end_session(scn), RES_OK); + + CHECK(s2d_device_ref_put(s2d), RES_OK); + CHECK(s2d_shape_ref_put(shape), RES_OK); + CHECK(s2d_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; +} +