From d8e0f2954de985d9dca44a1948129016ce8e004e Mon Sep 17 00:00:00 2001 From: Grzegorz Kowal Date: Thu, 11 Nov 2021 17:31:50 -0300 Subject: [PATCH] MESH: Rewrite redistribute_blocks(). This rewrite removes large buffers used by the MPI exchange of data block arrays between processes. Signed-off-by: Grzegorz Kowal --- sources/mesh.F90 | 75 ++++++++---------------------------------------- 1 file changed, 12 insertions(+), 63 deletions(-) diff --git a/sources/mesh.F90 b/sources/mesh.F90 index b2e30c5..fb4d570 100644 --- a/sources/mesh.F90 +++ b/sources/mesh.F90 @@ -834,8 +834,6 @@ module mesh subroutine redistribute_blocks() #ifdef MPI -! import external procedures and variables -! use blocks , only : block_meta, block_data, list_meta use blocks , only : get_nleafs, nregs use blocks , only : append_datablock, remove_datablock, link_blocks @@ -845,46 +843,26 @@ module mesh use mpitools , only : send_array, receive_array #endif /* MPI */ -! local variables are not implicit by default -! implicit none #ifdef MPI -! local variables -! logical :: flag integer :: status integer(kind=4) :: np, nl, nc, nr -! local pointers -! type(block_meta), pointer :: pmeta type(block_data), pointer :: pdata -! tag for the MPI data exchange -! - integer(kind=4) :: itag + integer(kind=4) :: tag1, tag2 -! array for number of data block for autobalancing -! integer(kind=4), dimension(nodes) :: nb integer(kind=4), dimension(lprocs,nodes) :: lb - -! local buffer for data block exchange -! -#if NDIMS == 3 - real(kind=8) , dimension(nv,nn,nn,nn,nregs+1) :: rbuf -#else /* NDIMS == 3 */ - real(kind=8) , dimension(nv,nn,nn, 1,nregs+1) :: rbuf -#endif /* NDIMS == 3 */ #endif /* MPI */ !------------------------------------------------------------------------------- ! #ifdef MPI #ifdef PROFILE -! start accounting time for the block redistribution -! call start_timer(ima) #endif /* PROFILE */ @@ -906,12 +884,7 @@ module mesh np = 0 nl = 0 -! set the pointer to the first block on the meta block list -! pmeta => list_meta - -! iterate over all meta blocks and reassign their process numbers -! do while (associated(pmeta)) ! consider only meta blocks which belong to active processes @@ -922,8 +895,6 @@ module mesh ! if (pmeta%process /= np) then -! check if the block is the leaf -! if (pmeta%leaf) then ! indicate that the block structure will change @@ -932,49 +903,33 @@ module mesh ! generate a tag for communication ! - itag = 16 * (np * nprocs + pmeta%process) + tag1 = np * nprocs + pmeta%process + tag2 = tag1 + 1 -! sends the block to the right process +! send the data block arrays to the destination process, and +! remove the data block from the current process ! if (nproc == pmeta%process) then -! copy data to buffer -! - rbuf(:,:,:,:,1:nregs) = pmeta%data%uu(:,:,:,:,:) - rbuf(:,:,:,:,nregs+1) = pmeta%data%q (:,:,:,:) + call send_array(np, tag1, pmeta%data%uu(:,:,:,:,:)) + call send_array(np, tag2, pmeta%data%q (:,:,:,:)) -! send data -! - call send_array(np, itag, rbuf) - -! remove data block from the current process -! call remove_datablock(pmeta%data, status) -! send data block -! end if ! nproc == pmeta%process -! receive the block from another process +! on the destination process, allocate new data block, associate it with +! the corresponding meta block, and fill it with the received arrays ! if (nproc == np) then -! allocate a new data block and link it with the current meta block -! call append_datablock(pdata, status) call link_blocks(pmeta, pdata) -! receive the data -! - call receive_array(pmeta%process, itag, rbuf) - -! coppy the buffer to data block -! - pmeta%data%uu(:,:,:,:,:) = rbuf(:,:,:,:,1:nregs) - pmeta%data%q (:,:,:,:) = rbuf(:,:,:,:,nregs+1) + call receive_array(pmeta%process, tag1, pmeta%data%uu(:,:,:,:,:)) + call receive_array(pmeta%process, tag2, pmeta%data%q (:,:,:,:)) end if ! nproc == n - end if ! leaf ! set new processor number @@ -984,7 +939,7 @@ module mesh end if ! pmeta%process /= np ! increase the number of blocks on the current process; if it exceeds the -! allowed number reset the counter and increase the processor number +! allowed number, reset the counter and increase the processor number ! if (pmeta%leaf) then @@ -1020,18 +975,12 @@ module mesh end if ! nl >= lb(nr,np) end if ! leaf - end if ! pmeta%process < nprocs -! assign the pointer to the next meta block -! pmeta => pmeta%next - end do ! pmeta #ifdef PROFILE -! stop accounting time for the block redistribution -! call stop_timer(ima) #endif /* PROFILE */ #endif /* MPI */