htrdr

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

commit 572eabf7ade0c9d1927a6eb5c4bc838a276569f9
parent d1befaf9fd97cb5756d697b742eb4ac7735b828d
Author: Vincent Forest <vincent.forest@meso-star.com>
Date:   Mon, 22 Feb 2021 11:30:13 +0100

Update the way progress messages are logged

Use the htrdr logger functions rather than explicit calls to fprintf
functions and thus remove the now useless htrdr_<fprintf|fflush>
functions. Take advantage of code cleaning to transform the
htrdr_mpi_error_string public function to the mpi_error_string helper
function.

Diffstat:
Msrc/core/htrdr.c | 68+++++++++++++++++++++++++-------------------------------------------
Msrc/core/htrdr.h | 25-------------------------
2 files changed, 25 insertions(+), 68 deletions(-)

diff --git a/src/core/htrdr.c b/src/core/htrdr.c @@ -60,6 +60,19 @@ mpi_thread_support_string(const int val) } } +static const char* +mpi_error_string(struct htrdr* htrdr, const int mpi_err) +{ + const int ithread = omp_get_thread_num(); + char* str; + int strlen_err; + int err; + ASSERT(htrdr && (size_t)ithread < htrdr->nthreads); + str = htrdr->mpi_err_str + ithread*MPI_MAX_ERROR_STRING; + err = MPI_Error_string(mpi_err, str, &strlen_err); + return err == MPI_SUCCESS ? str : "Invalid MPI error"; +} + static void release_mpi(struct htrdr* htrdr) { @@ -167,7 +180,7 @@ init_mpi(struct htrdr* htrdr) if(err != MPI_SUCCESS) { htrdr_log_err(htrdr, "could not determine the MPI rank of the calling process -- %s.\n", - htrdr_mpi_error_string(htrdr, err)); + mpi_error_string(htrdr, err)); res = RES_UNKNOWN_ERR; goto error; } @@ -176,7 +189,7 @@ init_mpi(struct htrdr* htrdr) if(err != MPI_SUCCESS) { htrdr_log_err(htrdr, "could retrieve the size of the MPI group -- %s.\n", - htrdr_mpi_error_string(htrdr, err)); + mpi_error_string(htrdr, err)); res = RES_UNKNOWN_ERR; goto error; } @@ -443,19 +456,6 @@ htrdr_get_s3d(struct htrdr* htrdr) return htrdr->s3d; } -const char* -htrdr_mpi_error_string(struct htrdr* htrdr, const int mpi_err) -{ - const int ithread = omp_get_thread_num(); - char* str; - int strlen_err; - int err; - ASSERT(htrdr && (size_t)ithread < htrdr->nthreads); - str = htrdr->mpi_err_str + ithread*MPI_MAX_ERROR_STRING; - err = MPI_Error_string(mpi_err, str, &strlen_err); - return err == MPI_SUCCESS ? str : "Invalid MPI error"; -} - res_T htrdr_open_output_stream (struct htrdr* htrdr, @@ -511,28 +511,6 @@ error: goto exit; } - -void -htrdr_fprintf(struct htrdr* htrdr, FILE* stream, const char* msg, ...) -{ - ASSERT(htrdr && msg); - if(htrdr->mpi_rank == 0) { - va_list vargs_list; - va_start(vargs_list, msg); - vfprintf(stream, msg, vargs_list); - va_end(vargs_list); - } -} - -void -htrdr_fflush(struct htrdr* htrdr, FILE* stream) -{ - ASSERT(htrdr); - if(htrdr->mpi_rank == 0) { - fflush(stream); - } -} - /******************************************************************************* * Local functions ******************************************************************************/ @@ -601,19 +579,17 @@ print_mpi_progress(struct htrdr* htrdr, const enum htrdr_mpi_message msg) if(htrdr->mpi_nprocs == 1) { switch(msg) { case HTRDR_MPI_PROGRESS_RENDERING: - htrdr_fprintf(htrdr, stderr, "\033[2K\rRendering: %3d%%", - htrdr->mpi_progress_render[0]); + htrdr_log(htrdr, "\33[2K\r"); /* Erase the line */ + htrdr_log(htrdr, "Rendering: %3d%%\r", htrdr->mpi_progress_render[0]); break; default: FATAL("Unreachable code.\n"); break; } - htrdr_fflush(htrdr, stderr); } else { int iproc; FOR_EACH(iproc, 0, htrdr->mpi_nprocs) { switch(msg) { case HTRDR_MPI_PROGRESS_RENDERING: - htrdr_fprintf(htrdr, stderr, - "\033[2K\rProcess %d -- rendering: %3d%%%c", + htrdr_log(htrdr, "Process %d -- rendering: %3d%%%c", iproc, htrdr->mpi_progress_render[iproc], iproc == htrdr->mpi_nprocs - 1 ? '\r' : '\n'); break; @@ -629,7 +605,13 @@ clear_mpi_progress(struct htrdr* htrdr, const enum htrdr_mpi_message msg) ASSERT(htrdr); (void)msg; if(htrdr->mpi_nprocs > 1) { - htrdr_fprintf(htrdr, stderr, "\033[%dA", htrdr->mpi_nprocs-1); + int iproc; + FOR_EACH(iproc, 0, htrdr->mpi_nprocs) { + htrdr_log(htrdr, "\33[2K\r"); /* Erase the line */ + if(iproc != htrdr->mpi_nprocs-1) { + htrdr_log(htrdr, "\033[1A\r"); /* Move up */ + } + } } } diff --git a/src/core/htrdr.h b/src/core/htrdr.h @@ -156,31 +156,6 @@ htrdr_open_output_stream int force_overwrite, FILE** out_fp); -/* TODO do not expose publicly this function(?) */ -HTRDR_CORE_API const char* -htrdr_mpi_error_string - (struct htrdr* htrdr, - const int mpi_err); - -/* TODO replace them by regular log message */ -HTRDR_CORE_API void -htrdr_fprintf - (struct htrdr* htrdr, - FILE* stream, - const char* msg, - ...) -#ifdef COMPILER_GCC - __attribute((format(printf, 3, 4))) -#endif - ; - -/* TODO remove this */ -HTRDR_CORE_API void -htrdr_fflush - (struct htrdr* htrdr, - FILE* stream); - - END_DECLS #endif /* HTRDR_H */