stardis-solver

Solve coupled heat transfers
git clone git://git.meso-star.fr/stardis-solver.git
Log | Files | Refs | README | LICENSE

commit be785897501f6118aec8fe4547881cf1883f40e3
parent 4181cf76fe1b24f83d864084bccf6f9cc5774084
Author: Vincent Forest <vincent.forest@meso-star.com>
Date:   Fri, 17 Dec 2021 16:37:16 +0100

Add MPI support to the function sdis_solve_medium

Diffstat:
Msrc/sdis_solve_medium_Xd.h | 236++++++++++++++++++++++++++++++++++++++++++-------------------------------------
1 file changed, 127 insertions(+), 109 deletions(-)

diff --git a/src/sdis_solve_medium_Xd.h b/src/sdis_solve_medium_Xd.h @@ -234,35 +234,42 @@ XD(solve_medium) struct sdis_green_function** out_green, /* May be NULL <=> No green func */ struct sdis_estimator** out_estimator) /* May be NULL <=> No estimator */ { - struct darray_enclosure_cumul cumul; + /* Time registration */ + struct time time0, time1; + char buf[128]; /* Temporary buffer used to store formated time */ + + /* Device variables */ + struct mem_allocator* allocator = NULL; + size_t nthreads = 0; + + /* Stardis variables */ + struct sdis_estimator* estimator = NULL; struct sdis_green_function* green = NULL; - struct sdis_green_function** greens = NULL; + struct sdis_green_function** per_thread_green = NULL; + + /* Random number generator */ struct ssp_rng_proxy* rng_proxy = NULL; - struct ssp_rng** rngs = NULL; - struct sdis_estimator* estimator = NULL; - struct accum* acc_temps = NULL; - struct accum* acc_times = NULL; + struct ssp_rng** per_thread_rng = NULL; + + /* Miscellaneous */ + struct darray_enclosure_cumul cumul; + struct accum* per_thread_acc_temp = NULL; + struct accum* per_thread_acc_time = NULL; size_t nrealisations = 0; int64_t irealisation; + int32_t* progress = NULL; /* Per process progress bar */ + int is_master_process = 1; int cumul_is_init = 0; - size_t i; - int progress = 0; int register_paths = SDIS_HEAT_PATH_NONE; ATOMIC nsolved_realisations = 0; ATOMIC res = RES_OK; - if(!scn) { - res = RES_BAD_ARG; - goto error; - } - + if(!scn) { res = RES_BAD_ARG; goto error; } + if(!out_estimator && !out_green) { res = RES_BAD_ARG; goto error; } res = check_solve_medium_args(args); if(res != RES_OK) goto error; - - if(!out_estimator && !out_green) { - res = RES_BAD_ARG; - goto error; - } + res = XD(scene_check_dimensionality)(scn); + if(res != RES_OK) goto error; if(out_green && args->picard_order != 1) { log_err(scn->dev, "%s: the evaluation of the green function does not make " @@ -273,38 +280,27 @@ XD(solve_medium) goto error; } -#if SDIS_XD_DIMENSION == 2 - if(scene_is_2d(scn) == 0) { res = RES_BAD_ARG; goto error; } -#else - if(scene_is_2d(scn) != 0) { res = RES_BAD_ARG; goto error; } +#ifdef SDIS_ENABLE_MPI + is_master_process = !scn->dev->use_mpi || scn->dev->mpi_rank == 0; #endif - /* Create the proxy RNG */ - if(args->rng_state) { - res = ssp_rng_proxy_create_from_rng(scn->dev->allocator, args->rng_state, - scn->dev->nthreads, &rng_proxy); - if(res != RES_OK) goto error; - } else { - res = ssp_rng_proxy_create(scn->dev->allocator, SSP_RNG_MT19937_64, - scn->dev->nthreads, &rng_proxy); - if(res != RES_OK) goto error; - } + nthreads = scn->dev->nthreads; + allocator = scn->dev->allocator; - /* Create the per thread RNG */ - rngs = MEM_CALLOC(scn->dev->allocator, scn->dev->nthreads, sizeof(*rngs)); - if(!rngs) { res = RES_MEM_ERR; goto error; } - FOR_EACH(i, 0, scn->dev->nthreads) { - res = ssp_rng_proxy_create_rng(rng_proxy, i, rngs+i); - if(res != RES_OK) goto error; - } + /* Create the per thread RNGs */ + res = create_per_thread_rng + (scn->dev, args->rng_state, &rng_proxy, &per_thread_rng); + if(res != RES_OK) goto error; + + /* Allocate the per process progress status */ + res = alloc_process_progress(scn->dev, &progress); + if(res != RES_OK) goto error; /* Create the per thread accumulators */ - acc_temps = MEM_CALLOC - (scn->dev->allocator, scn->dev->nthreads, sizeof(*acc_temps)); - if(!acc_temps) { res = RES_MEM_ERR; goto error; } - acc_times = MEM_CALLOC - (scn->dev->allocator, scn->dev->nthreads, sizeof(*acc_times)); - if(!acc_times) { res = RES_MEM_ERR; goto error; } + per_thread_acc_temp = MEM_CALLOC(allocator, nthreads, sizeof(struct accum)); + per_thread_acc_time = MEM_CALLOC(allocator, nthreads, sizeof(struct accum)); + if(!per_thread_acc_temp) { res = RES_MEM_ERR; goto error; } + if(!per_thread_acc_time) { res = RES_MEM_ERR; goto error; } /* Compute the enclosure cumulative */ darray_enclosure_cumul_init(scn->dev->allocator, &cumul); @@ -313,32 +309,39 @@ XD(solve_medium) if(res != RES_OK) goto error; if(out_green) { - /* Create the per thread green function */ - greens = MEM_CALLOC(scn->dev->allocator, scn->dev->nthreads, sizeof(*greens)); - if(!greens) { res = RES_MEM_ERR; goto error; } - FOR_EACH(i, 0, scn->dev->nthreads) { - res = green_function_create(scn, &greens[i]); - if(res != RES_OK) goto error; - } + res = create_per_thread_green_function(scn, &per_thread_green); + if(res != RES_OK) goto error; } - /* Create the estimator */ - if(out_estimator) { + /* Create the estimator on the master process only. No estimator is needed + * for non master process */ + if(out_estimator && is_master_process) { res = estimator_create(scn->dev, SDIS_ESTIMATOR_TEMPERATURE, &estimator); if(res != RES_OK) goto error; } - nrealisations = args->nrealisations; - register_paths = out_estimator ? args->register_paths : SDIS_HEAT_PATH_NONE; + /* Synchronise the processes */ + process_barrier(scn->dev); + + #define PROGRESS_MSG "Solving medium temperature: " + print_progress(scn->dev, progress, PROGRESS_MSG); + + /* Begin time registration of the computation */ + time_current(&time0); + + /* Here we go! Launch the Monte Carlo estimation */ + nrealisations = compute_process_realisations_count(scn->dev, args->nrealisations); + register_paths = out_estimator && is_master_process + ? args->register_paths : SDIS_HEAT_PATH_NONE; omp_set_num_threads((int)scn->dev->nthreads); #pragma omp parallel for schedule(static) for(irealisation = 0; irealisation < (int64_t)nrealisations; ++irealisation) { struct probe_realisation_args realis_args = PROBE_REALISATION_ARGS_NULL; struct time t0, t1; const int ithread = omp_get_thread_num(); - struct ssp_rng* rng = rngs[ithread]; - struct accum* acc_temp = &acc_temps[ithread]; - struct accum* acc_time = &acc_times[ithread]; + struct ssp_rng* rng = per_thread_rng[ithread]; + struct accum* acc_temp = &per_thread_acc_temp[ithread]; + struct accum* acc_time = &per_thread_acc_time[ithread]; struct green_path_handle* pgreen_path = NULL; struct green_path_handle green_path = GREEN_PATH_HANDLE_NULL; const struct enclosure* enc = NULL; @@ -358,11 +361,8 @@ XD(solve_medium) time = sample_time(rng, args->time_range); if(out_green) { - res_local = green_function_create_path(greens[ithread], &green_path); - if(res_local != RES_OK) { - ATOMIC_SET(&res, res_local); - goto error_it; - } + res_local = green_function_create_path(per_thread_green[ithread], &green_path); + if(res_local != RES_OK) { ATOMIC_SET(&res, res_local); goto error_it; } pgreen_path = &green_path; } @@ -430,9 +430,9 @@ XD(solve_medium) n = (size_t)ATOMIC_INCR(&nsolved_realisations); pcent = (int)((double)n * 100.0 / (double)nrealisations + 0.5/*round*/); #pragma omp critical - if(pcent > progress) { - progress = pcent; - log_info(scn->dev, "Solving medium temperature: %3d%%\r", progress); + if(pcent > progress[0]) { + progress[0] = pcent; + print_progress_update(scn->dev, progress, PROGRESS_MSG); } exit_it: if(pheat_path) heat_path_release(pheat_path); @@ -440,71 +440,89 @@ XD(solve_medium) error_it: goto exit_it; } + /* Synchronise processes */ + process_barrier(scn->dev); + + res = gather_res_T(scn->dev, (res_T)res); if(res != RES_OK) goto error; - /* Add a new line after the progress status */ - log_info(scn->dev, "Solving medium temperature: %3d%%\n", progress); + print_progress_update(scn->dev, progress, PROGRESS_MSG); + log_info(scn->dev, "\n"); + #undef PROGRESS_MSG + + /* Report computation time */ + time_sub(&time0, time_current(&time1), &time0); + time_dump(&time0, TIME_ALL, NULL, buf, sizeof(buf)); + log_info(scn->dev, "Medium temperature solved in %s.\n", buf); + + /* Gather the RNG proxy sequence IDs and ensure that the RNG proxy state of + * the master process is greater than the RNG proxy state of all other + * processes */ + res = gather_rng_proxy_sequence_id(scn->dev, rng_proxy); + if(res != RES_OK) goto error; /* Setup the estimated temperature */ if(out_estimator) { struct accum acc_temp; struct accum acc_time; - sum_accums(acc_temps, scn->dev->nthreads, &acc_temp); - sum_accums(acc_times, scn->dev->nthreads, &acc_time); - ASSERT(acc_temp.count == acc_time.count); - - estimator_setup_realisations_count(estimator, nrealisations, acc_temp.count); - estimator_setup_temperature(estimator, acc_temp.sum, acc_temp.sum2); - estimator_setup_realisation_time(estimator, acc_time.sum, acc_time.sum2); - res = estimator_save_rng_state(estimator, rng_proxy); - if(res != RES_OK) goto error; + time_current(&time0); + + #define GATHER_ACCUMS(Msg, Acc) { \ + res = gather_accumulators(scn->dev, Msg, per_thread_##Acc, &Acc); \ + if(res != RES_OK) goto error; \ + } (void)0 + GATHER_ACCUMS(MPI_SDIS_MSG_ACCUM_TEMP, acc_temp); + GATHER_ACCUMS(MPI_SDIS_MSG_ACCUM_TIME, acc_time); + #undef GATHER_ACCUMS + + time_sub(&time0, time_current(&time1), &time0); + time_dump(&time0, TIME_ALL, NULL, buf, sizeof(buf)); + log_info(scn->dev, "Accumulators gathered in %s.\n", buf); + + /* Return an estimator only on master process */ + if(is_master_process) { + ASSERT(acc_temp.count == acc_time.count); + estimator_setup_realisations_count(estimator, args->nrealisations, acc_temp.count); + estimator_setup_temperature(estimator, acc_temp.sum, acc_temp.sum2); + estimator_setup_realisation_time(estimator, acc_time.sum, acc_time.sum2); + res = estimator_save_rng_state(estimator, rng_proxy); + if(res != RES_OK) goto error; + } } if(out_green) { - struct accum acc_time; + time_current(&time0); - /* Redux the per thread green function into the green of the 1st thread */ - green = greens[0]; /* Return the green of the 1st thread */ - greens[0] = NULL; /* Make invalid the 1st green for 'on exit' clean up*/ - res = green_function_redux_and_clear(green, greens+1, scn->dev->nthreads-1); + res = gather_green_functions + (scn, rng_proxy, per_thread_green, per_thread_acc_time, &green); if(res != RES_OK) goto error; - /* Finalize the estimated green */ - sum_accums(acc_times, scn->dev->nthreads, &acc_time); - res = green_function_finalize(green, rng_proxy, &acc_time); - if(res != RES_OK) goto error; - } + time_sub(&time0, time_current(&time1), &time0); + time_dump(&time0, TIME_ALL, NULL, buf, sizeof(buf)); + log_info(scn->dev, "Green functions gathered in %s.\n", buf); -exit: - if(rngs) { - FOR_EACH(i, 0, scn->dev->nthreads) { - if(rngs[i]) SSP(rng_ref_put(rngs[i])); + /* Return a green function only on master process */ + if(!is_master_process) { + SDIS(green_function_ref_put(green)); + green = NULL; } - MEM_RM(scn->dev->allocator, rngs); } - if(greens) { - FOR_EACH(i, 0, scn->dev->nthreads) { - if(greens[i]) SDIS(green_function_ref_put(greens[i])); - } - MEM_RM(scn->dev->allocator, greens); - } - if(acc_temps) MEM_RM(scn->dev->allocator, acc_temps); - if(acc_times) MEM_RM(scn->dev->allocator, acc_times); + +exit: + if(per_thread_rng) release_per_thread_rng(scn->dev, per_thread_rng); + if(per_thread_green) release_per_thread_green_function(scn, per_thread_green); + if(progress) free_process_progress(scn->dev, progress); + if(per_thread_acc_temp) MEM_RM(scn->dev->allocator, per_thread_acc_temp); + if(per_thread_acc_time) MEM_RM(scn->dev->allocator, per_thread_acc_time); if(cumul_is_init) darray_enclosure_cumul_release(&cumul); if(rng_proxy) SSP(rng_proxy_ref_put(rng_proxy)); if(out_estimator) *out_estimator = estimator; if(out_green) *out_green = green; return (res_T)res; error: - if(green) { - SDIS(green_function_ref_put(green)); - green = NULL; - } - if(estimator) { - SDIS(estimator_ref_put(estimator)); - estimator = NULL; - } + if(green) { SDIS(green_function_ref_put(green)); green = NULL; } + if(estimator) { SDIS(estimator_ref_put(estimator)); estimator = NULL; } goto exit; }