star-gf

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

commit 2816d41307f8ef530049c5a123dae2dc401809f1
parent 34368761dc359da86f418c7b2b52edefab6a9212
Author: Vincent Forest <vincent.forest@meso-star.com>
Date:   Thu, 24 Sep 2015 11:11:22 +0200

Add parallelism to the sgf_cube test

Diffstat:
Msrc/test_sgf_cube.c | 49+++++++++++++++++++++++++++++++++----------------
1 file changed, 33 insertions(+), 16 deletions(-)

diff --git a/src/test_sgf_cube.c b/src/test_sgf_cube.c @@ -35,7 +35,9 @@ #include <star/s3d.h> #include <star/ssp.h> -#define NSTEPS 10000L +#include <omp.h> + +#define NSTEPS 100000L /* * Matrix of Gebhart factors @@ -102,19 +104,24 @@ main(int argc, char** argv) struct sgf_device* sgf = NULL; struct sgf_scene_desc scn_desc = SGF_SCENE_DESC_NULL; struct sgf_estimator* estimators = NULL; - struct ssp_rng* rng; + struct ssp_rng_proxy* proxy; + struct ssp_rng** rngs = NULL; size_t iprim; + unsigned i; + unsigned nbuckets; (void)argc, (void)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_create(&allocator, &ssp_rng_threefry, &rng), RES_OK); + CHECK(ssp_rng_proxy_create(&allocator, &ssp_rng_threefry, nbuckets, &proxy), RES_OK); CHECK(sgf_device_create(NULL, NULL, 1, &sgf), RES_OK); + attribs[0].type = S3D_FLOAT3; attribs[0].usage = S3D_POSITION; attribs[0].get = get_pos; @@ -137,6 +144,11 @@ main(int argc, char** argv) scn_desc.scene = scn; estimators = sa_add(estimators, nprims * nprims); + rngs = sa_add(rngs, nbuckets); + + FOR_EACH(i, 0, nbuckets) { + CHECK(ssp_rng_proxy_create_rng(proxy, i, rngs + i), RES_OK); + } CHECK(sgf_begin_integrate(NULL, NULL), RES_BAD_ARG); CHECK(sgf_begin_integrate(sgf, NULL), RES_BAD_ARG); @@ -146,20 +158,20 @@ main(int argc, char** argv) CHECK(sgf_integrate(NULL, NULL, SIZE_MAX, 0), RES_BAD_ARG); CHECK(sgf_integrate(sgf, NULL, SIZE_MAX, 0), RES_BAD_ARG); - CHECK(sgf_integrate(NULL, rng, SIZE_MAX, 0), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, SIZE_MAX, 0), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, rngs[0], SIZE_MAX, 0), RES_BAD_ARG); + CHECK(sgf_integrate(sgf, rngs[0], SIZE_MAX, 0), RES_BAD_ARG); CHECK(sgf_integrate(NULL, NULL, 0, 0), RES_BAD_ARG); CHECK(sgf_integrate(sgf, NULL, 0, 0), RES_BAD_ARG); - CHECK(sgf_integrate(NULL, rng, 0, 0), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 0, 0), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, rngs[0], 0, 0), RES_BAD_ARG); + CHECK(sgf_integrate(sgf, rngs[0], 0, 0), RES_BAD_ARG); CHECK(sgf_integrate(NULL, NULL, SIZE_MAX, NSTEPS), RES_BAD_ARG); CHECK(sgf_integrate(sgf, NULL, SIZE_MAX, NSTEPS), RES_BAD_ARG); - CHECK(sgf_integrate(NULL, rng, SIZE_MAX, NSTEPS), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, SIZE_MAX, NSTEPS), RES_BAD_ARG); + CHECK(sgf_integrate(NULL, rngs[0], SIZE_MAX, NSTEPS), RES_BAD_ARG); + CHECK(sgf_integrate(sgf, rngs[0], SIZE_MAX, NSTEPS), RES_BAD_ARG); CHECK(sgf_integrate(NULL, NULL, 0, NSTEPS), RES_BAD_ARG); CHECK(sgf_integrate(sgf, NULL, 0, NSTEPS), RES_BAD_ARG); - CHECK(sgf_integrate(NULL, rng, 0, NSTEPS), RES_BAD_ARG); - CHECK(sgf_integrate(sgf, rng, 0, NSTEPS), RES_OK); + CHECK(sgf_integrate(NULL, rngs[0], 0, NSTEPS), RES_BAD_ARG); + CHECK(sgf_integrate(sgf, rngs[0], 0, NSTEPS), RES_OK); CHECK(sgf_get(NULL, SIZE_MAX, SIZE_MAX, NULL), RES_BAD_ARG); CHECK(sgf_get(sgf, SIZE_MAX, SIZE_MAX, NULL), RES_BAD_ARG); @@ -184,18 +196,19 @@ main(int argc, char** argv) CHECK(sgf_end_integrate(sgf), RES_OK); CHECK(sgf_end_integrate(sgf), RES_BAD_OP); - CHECK(sgf_integrate(sgf, rng, 0, NSTEPS), RES_BAD_OP); + CHECK(sgf_integrate(sgf, rngs[0], 0, NSTEPS), RES_BAD_OP); CHECK(sgf_get(sgf, 0, 0, estimators + 0), RES_BAD_OP); /* Integrate the Gebhart Factors */ CHECK(sgf_begin_integrate(sgf, &scn_desc), RES_OK); -/* #pragma omp parallel for */ + #pragma omp parallel for for(iprim = 0; iprim < nprims; ++iprim) { size_t iprim2; struct sgf_estimator* row = estimators + iprim * nprims; + const int ithread = omp_get_thread_num(); - CHECK(sgf_integrate(sgf, rng, iprim, NSTEPS), RES_OK); + CHECK(sgf_integrate(sgf, rngs[ithread], iprim, NSTEPS), RES_OK); FOR_EACH(iprim2, 0, nprims) { CHECK(sgf_get(sgf, iprim, iprim2, row + iprim2), RES_OK); @@ -222,7 +235,7 @@ main(int argc, char** argv) printf("%.6f ", E); } printf("\n"); - CHECK(eq_eps(sum, 1.0, 1.e-6), 1); /* Ensure the energy conservation */ + CHECK(eq_eps(sum, 1.0, 1.e-3), 1); /* Ensure the energy conservation */ } sa_release(estimators); @@ -230,8 +243,12 @@ main(int argc, char** argv) 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_ref_put(rng), RES_OK); + CHECK(ssp_rng_proxy_ref_put(proxy), RES_OK); CHECK(sgf_device_ref_put(sgf), RES_OK); + FOR_EACH(i, 0, nbuckets) { + CHECK(ssp_rng_ref_put(rngs[i]), RES_OK); + } + sa_release(rngs); check_memory_allocator(&allocator); mem_shutdown_proxy_allocator(&allocator);