htrdr

Solving radiative transfer in heterogeneous media
git clone git://git.meso-star.fr/htrdr.git
Log | Files | Refs | README | LICENSE

commit ba238c823df8b7dfefa3eb9c0a8fcdbad1114106
parent 1ac4a8bd590835683bd6287edf708fcd9589b2d6
Author: Vincent Forest <vincent.forest@meso-star.com>
Date:   Mon, 22 Oct 2018 16:20:01 +0200

Fix the progress status message with MPI

Only the progress of the process 0 was reported. With this commit, the
process 0 now reports the overall progress status wrt all participating
processes

Diffstat:
Msrc/htrdr_draw_radiance_sw.c | 86+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++----------
1 file changed, 76 insertions(+), 10 deletions(-)

diff --git a/src/htrdr_draw_radiance_sw.c b/src/htrdr_draw_radiance_sw.c @@ -26,6 +26,7 @@ #include <omp.h> #include <mpi.h> +#include <unistd.h> #define RNG_SEQUENCE_SIZE 1000000 @@ -207,6 +208,36 @@ draw_tile return RES_OK; } +/* Return the overall percentage */ +static size_t +fetch_process_progress(struct htrdr* htrdr, size_t* progress) +{ + size_t overall_pcent; + int iproc; + ASSERT(htrdr && progress && htrdr->mpi_rank == 0); + + overall_pcent = progress[0]; + FOR_EACH(iproc, 1, htrdr->mpi_nprocs) { + size_t proc_pcent = progress[iproc]; + + /* Flush the last sent percentage of the process `iproc' */ + for(;;) { + int flag; + + CHK(MPI_Iprobe(iproc, 0/*tag*/, MPI_COMM_WORLD, &flag, + MPI_STATUS_IGNORE) == MPI_SUCCESS); + if(flag == 0) break; /* No more message */ + + CHK(MPI_Recv(&proc_pcent, sizeof(size_t), MPI_CHAR, iproc, 0/*tag*/, + MPI_COMM_WORLD, MPI_STATUS_IGNORE) == MPI_SUCCESS); + } + + progress[iproc] = proc_pcent; + overall_pcent += progress[iproc]; + } + return overall_pcent / (size_t)htrdr->mpi_nprocs; +} + /******************************************************************************* * Local functions ******************************************************************************/ @@ -220,12 +251,13 @@ htrdr_draw_radiance_sw struct ssp_rng_proxy* rng_proxy = NULL; struct ssp_rng** rngs = NULL; size_t ntiles_x, ntiles_y, ntiles, ntiles_adjusted; - size_t progress = 0; + size_t* progress = NULL; size_t i; int32_t mcode; /* Morton code of the tile */ struct htrdr_buffer_layout layout; double pix_sz[2]; /* Pixel size in the normalized image plane */ size_t spp; + size_t overall_pcent = 0; ATOMIC nsolved_tiles = 0; ATOMIC res = RES_OK; ASSERT(htrdr && cam && buf); @@ -258,7 +290,8 @@ htrdr_draw_radiance_sw RNG_SEQUENCE_SIZE * (size_t)htrdr->mpi_nprocs, /* Pitch */ htrdr->nthreads, &rng_proxy); if(res != RES_OK) { - htrdr_log_err(htrdr, "%s: could not create the RNG proxy.\n", FUNC_NAME); + htrdr_log_err(htrdr, "%s: could not create the RNG proxy -- %s.\n", + FUNC_NAME, res_to_cstr((res_T)res)); goto error; } @@ -272,11 +305,26 @@ htrdr_draw_radiance_sw FOR_EACH(i, 0, htrdr->nthreads) { res = ssp_rng_proxy_create_rng(rng_proxy, i, &rngs[i]); if(res != RES_OK) { - htrdr_log_err(htrdr,"%s: could not create the per thread RNG.\n", FUNC_NAME); + htrdr_log_err(htrdr,"%s: could not create the per thread RNG -- %s.\n", + FUNC_NAME, res_to_cstr((res_T)res)); goto error; } } + if(htrdr->mpi_rank == 0) { + progress = MEM_CALLOC + (htrdr->allocator, (size_t)htrdr->mpi_nprocs, sizeof(size_t)); + } else { + progress = MEM_CALLOC + (htrdr->allocator, 1, sizeof(size_t)); + } + if(!progress) { + htrdr_log_err(htrdr, + "%s: could not allocate the process progress counter.\n", FUNC_NAME); + res = RES_MEM_ERR; + goto error; + } + ntiles_x = (layout.width + (TILE_SIZE-1)/*ceil*/)/TILE_SIZE; ntiles_y = (layout.height+ (TILE_SIZE-1)/*ceil*/)/TILE_SIZE; ntiles_adjusted = round_up_pow2(MMAX(ntiles_x, ntiles_y)); @@ -290,7 +338,6 @@ htrdr_draw_radiance_sw htrdr_fflush(htrdr, stderr); omp_set_num_threads((int)htrdr->nthreads); - #pragma omp parallel for schedule(static, 1/*chunck size*/) for(mcode=0; mcode<(int64_t)ntiles_adjusted; ++mcode) { const int ithread = omp_get_thread_num(); @@ -326,16 +373,34 @@ htrdr_draw_radiance_sw pcent = n * 100 / ntiles; #pragma omp critical - if(pcent > progress) { - progress = pcent; - htrdr_fprintf(htrdr, stderr, "%c[2K\rRendering: %3lu%%", - 27, (unsigned long)pcent); - htrdr_fflush(htrdr, stderr); + if((size_t)pcent > progress[0]) { + progress[0] = pcent; + + if(htrdr->mpi_rank != 0) { + /* Send the progress percentage of the process to the master process */ + CHK(MPI_Send(&pcent, sizeof(size_t), MPI_CHAR, 0/*dst*/, 0/*tag*/, + MPI_COMM_WORLD) == MPI_SUCCESS); + } else { + overall_pcent = fetch_process_progress(htrdr, progress); + htrdr_fprintf(htrdr, stderr, "%c[2K\rRendering: %3lu%%", + 27, (unsigned long)overall_pcent); + htrdr_fflush(htrdr, stderr); + } } if(ATOMIC_GET(&res) != RES_OK) continue; } - htrdr_fprintf(htrdr, stderr, "\n"); + + if(htrdr->mpi_rank == 0) { + while(overall_pcent != 100) { + overall_pcent = fetch_process_progress(htrdr, progress); + htrdr_fprintf(htrdr, stderr, "%c[2K\rRendering: %3lu%%", + 27, (unsigned long)overall_pcent); + htrdr_fflush(htrdr, stderr); + sleep(1); + } + htrdr_fprintf(htrdr, stderr, "\n"); + } /* Gather accum buffers from the group of processes */ res = gather_buffer(htrdr, buf); @@ -343,6 +408,7 @@ htrdr_draw_radiance_sw exit: if(rng_proxy) SSP(rng_proxy_ref_put(rng_proxy)); + if(progress) MEM_RM(htrdr->allocator, progress); if(rngs) { FOR_EACH(i, 0, htrdr->nthreads) { if(rngs[i]) SSP(rng_ref_put(rngs[i]));