stardis-solver

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

commit c41b28861fc419fdb65e87e509a1fc87a1e6f325
parent e8fef26eeb9058cf9c043a4811e6446e688ec9f7
Author: Vincent Forest <vincent.forest@meso-star.com>
Date:   Thu,  9 Dec 2021 16:10:59 +0100

Add MPI support to the function sdis_solve_boundary

Diffstat:
Msrc/sdis_solve_boundary_Xd.h | 217+++++++++++++++++++++++++++++++++++++++++++++----------------------------------
1 file changed, 125 insertions(+), 92 deletions(-)

diff --git a/src/sdis_solve_boundary_Xd.h b/src/sdis_solve_boundary_Xd.h @@ -117,24 +117,40 @@ XD(solve_boundary) struct sdis_green_function** out_green, struct sdis_estimator** out_estimator) { + /* 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** per_thread_green = NULL; + + /* Random number generator */ + struct ssp_rng_proxy* rng_proxy = NULL; + struct ssp_rng** per_thread_rng = NULL; + + /* XD variables */ struct XD(boundary_context) ctx = XD(BOUNDARY_CONTEXT_NULL); struct sXd(vertex_data) vdata = SXD_VERTEX_DATA_NULL; struct sXd(scene)* scene = NULL; struct sXd(shape)* shape = NULL; struct sXd(scene_view)* view = NULL; - struct sdis_estimator* estimator = NULL; - struct sdis_green_function* green = NULL; - struct sdis_green_function** greens = NULL; - struct ssp_rng_proxy* rng_proxy = NULL; - struct ssp_rng** rngs = NULL; - struct accum* acc_temps = NULL; - struct accum* acc_times = NULL; + + /* Miscellaneous */ + struct accum* per_thread_acc_temp = NULL; + struct accum* per_thread_acc_time = NULL; size_t nrealisations = 0; int64_t irealisation = 0; - size_t i; size_t view_nprims; - int progress = 0; + size_t i; + int32_t* progress = NULL; /* Per process progress bar */ int register_paths = SDIS_HEAT_PATH_NONE; + int is_master_process = 1; ATOMIC nsolved_realisations = 0; ATOMIC res = RES_OK; @@ -150,7 +166,7 @@ XD(solve_boundary) res = RES_BAD_ARG; goto error; } - + if(out_green && args->picard_order != 1) { log_err(scn->dev, "%s: the evaluation of the green function does not make " "sense when dealing with the non-linearities of the system; i.e. picard " @@ -186,6 +202,10 @@ XD(solve_boundary) } } +#ifdef SDIS_ENABLE_MPI + is_master_process = !scn->dev->use_mpi || scn->dev->mpi_rank == 0; +#endif + /* Create the Star-XD shape of the boundary */ #if SDIS_XD_DIMENSION == 2 res = s2d_shape_create_line_segments(scn->dev->sXd(dev), &shape); @@ -221,60 +241,59 @@ XD(solve_boundary) res = sXd(scene_view_create)(scene, SXD_SAMPLE, &view); if(res != RES_OK) goto error; - /* 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; } /* Create the per thread green function */ if(out_green) { - 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 surface 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 boundary_realisation_args realis_args = BOUNDARY_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; struct sdis_heat_path* pheat_path = NULL; @@ -298,7 +317,8 @@ XD(solve_boundary) time = sample_time(rng, args->time_range); if(out_green) { - res_local = green_function_create_path(greens[ithread], &green_path); + 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; @@ -392,9 +412,9 @@ XD(solve_boundary) 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 boundary 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); @@ -402,58 +422,77 @@ XD(solve_boundary) 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 boundary 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, "Surface probe 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; + struct accum acc_temp, 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); + time_current(&time0); - 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); + res = gather_accumulators + (scn->dev, MPI_SDIS_MSG_ACCUM_TEMP, per_thread_acc_temp, &acc_temp); + if(res != RES_OK) goto error; + res = gather_accumulators + (scn->dev, MPI_SDIS_MSG_ACCUM_TIME, per_thread_acc_time, &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, "Accumulators gathered in %s.\n", buf); + + /* Return an estimator only on master process */ + if(is_master_process) { + res = setup_estimator + (estimator, rng_proxy, &acc_temp, &acc_time, args->nrealisations); + if(res != RES_OK) goto error; + } } /* Setup the green function */ 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])); - } - 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])); + /* 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, 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(scene) SXD(scene_ref_put(scene)); if(shape) SXD(shape_ref_put(shape)); if(view) SXD(scene_view_ref_put(view)); @@ -462,14 +501,8 @@ exit: if(out_estimator) *out_estimator = estimator; return (res_T)res; error: - 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; } + if(green) { SDIS(green_function_ref_put(green)); green = NULL; } goto exit; }