Adds up ghostcells in the z direction with non-blocking mpi directives.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(eulerian_obj_base), | intent(inout) | :: | this |
An Eulerian object |
impure subroutine eulerian_obj_AddUpGhostCells_z(this) !> Adds up ghostcells in the z direction ! with non-blocking mpi directives. implicit none class(eulerian_obj_base), intent(inout) :: this !! An Eulerian object ! Work variables integer :: Ng(3) integer :: sendL(3),sendR(3) type(MPI_STATUS) :: statuses(4) type(MPI_REQUEST) :: requests(4) integer :: ierr real(wp),allocatable :: buffR_r(:,:,:),buffL_r(:,:,:) integer, allocatable :: buffR_i(:,:,:),buffL_i(:,:,:) integer :: i,j,k,ii,jj,kk associate (lo => this%block%lo, hi => this%block%hi, & ngc=> this%block%ngc, mpi => this%parallel, & gc_slab_r => this%block%gc_slab_r, & gc_slab_i => this%block%gc_slab_i) ! Number of grid points (including ghostcells) Ng=hi-lo+1+2*ngc ! 2D exceptions if (Ng(3)-2*ngc.le.ngc) then select type (this) type is (eulerian_obj_r) do i=1,ngc this%cell(:,:,hi(3)) = this%cell(:,:,hi(3)) + this%cell(:,:,hi(3)+i) this%cell(:,:,lo(3)) = this%cell(:,:,lo(3)) + this%cell(:,:,lo(3)-i) end do type is (eulerian_obj_i) do i=1,ngc this%cell(:,:,hi(3)) = this%cell(:,:,hi(3)) + this%cell(:,:,hi(3)+i) this%cell(:,:,lo(3)) = this%cell(:,:,lo(3)) + this%cell(:,:,lo(3)-i) end do end select return end if ! Address of first element in buffer sendL(1)=lo(1)-ngc ; sendL(2)=lo(2)-ngc ; sendL(3)=lo(3)-ngc sendR(1)=lo(1)-ngc ; sendR(2)=lo(2)-ngc ; sendR(3)=hi(3)+1 select type (this) type is (eulerian_obj_r) ! Allocate receive buffers allocate(buffR_r(Ng(1),Ng(2),ngc)) allocate(buffL_r(Ng(1),Ng(2),ngc)) ! Post receives from Left and right ranks call MPI_IRECV( buffL_r(1,1,1), ngc*Ng(2)*Ng(1), mpi%REAL_WP, & mpi%rank%L(3)-1, 0, mpi%comm%g, requests(1), ierr) call MPI_IRECV( buffR_r(1,1,1), ngc*Ng(2)*Ng(1), mpi%REAL_WP, & mpi%rank%R(3)-1, 0, mpi%comm%g, requests(2), ierr) ! Send buffers to left and right ranks call MPI_ISEND( this%cell(sendR(1),sendR(2),sendR(3)), 1, gc_slab_r(3), & mpi%rank%R(3)-1, 0, mpi%comm%g, requests(4), ierr) call MPI_ISEND( this%cell(sendL(1),sendL(2),sendL(3)), 1, gc_slab_r(3), & mpi%rank%L(3)-1, 0, mpi%comm%g, requests(3), ierr) type is (eulerian_obj_i) ! Allocate receive buffers allocate(buffR_i(Ng(1),Ng(2),ngc)) allocate(buffL_i(Ng(1),Ng(2),ngc)) ! Post receives from Left and right ranks call MPI_IRECV( buffL_i(1,1,1), ngc*Ng(2)*Ng(1), mpi%INTEGER, & mpi%rank%L(3)-1, 0, mpi%comm%g, requests(1), ierr) call MPI_IRECV( buffR_i(1,1,1), ngc*Ng(2)*Ng(1), mpi%INTEGER, & mpi%rank%R(3)-1, 0, mpi%comm%g, requests(2), ierr) ! Send buffers to left and right ranks call MPI_ISEND( this%cell(sendR(1),sendR(2),sendR(3)), 1, gc_slab_i(3), & mpi%rank%R(3)-1, 0, mpi%comm%g, requests(4), ierr) call MPI_ISEND( this%cell(sendL(1),sendL(2),sendL(3)), 1, gc_slab_i(3), & mpi%rank%L(3)-1, 0, mpi%comm%g, requests(3), ierr) end select ! Synchronize call MPI_WAITALL( 4, requests, statuses, ierr ) select type(this) type is (eulerian_obj_r) ! Add left buffer to left ghostcells if (mpi%rank%L(3)-1.ne.MPI_PROC_NULL) then do k=lo(3),lo(3)+ngc-1 do j=lo(2)-ngc,hi(2)+ngc do i=lo(1)-ngc,hi(1)+ngc ii=i-(lo(1)-ngc )+1 jj=j-(lo(2)-ngc )+1 kk=k-(lo(3) )+1 this%cell(i,j,k)=this%cell(i,j,k) + buffL_r(ii,jj,kk) end do end do end do end if ! Add right buffer to right ghostcells if (mpi%rank%R(3)-1.ne.MPI_PROC_NULL) then do k=hi(3)-ngc+1,hi(3) do j=lo(2)-ngc,hi(2)+ngc do i=lo(1)-ngc,hi(1)+ngc ii=i-(lo(1)-ngc )+1 jj=j-(lo(2)-ngc )+1 kk=k-(hi(3)-ngc+1)+1 this%cell(i,j,k)=this%cell(i,j,k) + buffR_r(ii,jj,kk) end do end do end do end if deallocate(buffL_r,buffR_r) type is (eulerian_obj_i) ! Add left buffer to left ghostcells if (mpi%rank%L(3)-1.ne.MPI_PROC_NULL) then do k=lo(3),lo(3)+ngc-1 do j=lo(2)-ngc,hi(2)+ngc do i=lo(1)-ngc,hi(1)+ngc ii=i-(lo(1)-ngc )+1 jj=j-(lo(2)-ngc )+1 kk=k-(lo(3) )+1 this%cell(i,j,k)=this%cell(i,j,k) + buffL_i(ii,jj,kk) end do end do end do end if ! Add right buffer to right ghostcells if (mpi%rank%R(3)-1.ne.MPI_PROC_NULL) then do k=hi(3)-ngc+1,hi(3) do j=lo(2)-ngc,hi(2)+ngc do i=lo(1)-ngc,hi(1)+ngc ii=i-(lo(1)-ngc )+1 jj=j-(lo(2)-ngc )+1 kk=k-(hi(3)-ngc+1)+1 this%cell(i,j,k)=this%cell(i,j,k) + buffR_i(ii,jj,kk) end do end do end do end if deallocate(buffL_i,buffR_i) end select end associate return end subroutine eulerian_obj_AddUpGhostCells_z