eulerian_obj_AddUpGhostCells_y Subroutine

private impure subroutine eulerian_obj_AddUpGhostCells_y(this)

Adds up ghostcells in the y direction with non-blocking mpi directives.

Type Bound

eulerian_obj_base

Arguments

Type IntentOptional Attributes Name
class(eulerian_obj_base), intent(inout) :: this

An Eulerian object


Calls

proc~~eulerian_obj_addupghostcells_y~~CallsGraph proc~eulerian_obj_addupghostcells_y eulerian_obj_base%eulerian_obj_AddUpGhostCells_y cell cell proc~eulerian_obj_addupghostcells_y->cell mpi_irecv mpi_irecv proc~eulerian_obj_addupghostcells_y->mpi_irecv mpi_isend mpi_isend proc~eulerian_obj_addupghostcells_y->mpi_isend mpi_waitall mpi_waitall proc~eulerian_obj_addupghostcells_y->mpi_waitall

Called by

proc~~eulerian_obj_addupghostcells_y~~CalledByGraph proc~eulerian_obj_addupghostcells_y eulerian_obj_base%eulerian_obj_AddUpGhostCells_y proc~eulerian_obj_addupghostcells eulerian_obj_base%eulerian_obj_AddUpGhostCells proc~eulerian_obj_addupghostcells->proc~eulerian_obj_addupghostcells_y proc~marker_set_filter marker_set%marker_set_Filter proc~marker_set_filter->proc~eulerian_obj_addupghostcells proc~particle_set_filter particle_set%particle_set_Filter proc~particle_set_filter->proc~eulerian_obj_addupghostcells proc~marker_set_computesolidvolfrac marker_set%marker_set_ComputeSolidVolFrac proc~marker_set_computesolidvolfrac->proc~marker_set_filter proc~marker_set_getibforcing marker_set%marker_set_GetIBForcing proc~marker_set_getibforcing->proc~marker_set_filter proc~marker_set_updatenormals marker_set%marker_set_UpdateNormals proc~marker_set_updatenormals->proc~marker_set_filter proc~marker_set_updatesdf marker_set%marker_set_UpdateSDF proc~marker_set_updatesdf->proc~marker_set_filter proc~respart_set_filter ResPart_set%ResPart_set_Filter proc~respart_set_filter->proc~marker_set_filter proc~solid_set_filter solid_set%solid_set_Filter proc~solid_set_filter->proc~marker_set_filter proc~respart_set_getibforcing ResPart_set%ResPart_set_GetIBForcing proc~respart_set_getibforcing->proc~marker_set_getibforcing proc~respart_set_updatenormals ResPart_set%ResPart_set_UpdateNormals proc~respart_set_updatenormals->proc~marker_set_updatenormals proc~respart_set_updatesdf ResPart_set%ResPart_set_UpdateSDF proc~respart_set_updatesdf->proc~marker_set_updatesdf proc~cdifs_obj_advancesolutionib cdifs_obj_AdvanceSolutionIB proc~cdifs_obj_advancesolutionib->proc~respart_set_getibforcing proc~cdifs_obj_advancesolutionib->proc~respart_set_updatesdf proc~cdifs_obj_computesolidvf cdifs_obj%cdifs_obj_ComputeSolidVF proc~cdifs_obj_computesolidvf->proc~respart_set_updatenormals proc~cdifs_obj_preparesolver cdifs_obj_PrepareSolver proc~cdifs_obj_preparesolver->proc~respart_set_updatesdf proc~cdifs_obj_preparesolver->proc~cdifs_obj_computesolidvf proc~grans_obj_computesolidvf grans_obj%grans_obj_ComputeSolidVF proc~grans_obj_computesolidvf->proc~respart_set_updatenormals proc~grans_obj_preparesolver grans_obj_PrepareSolver proc~grans_obj_preparesolver->proc~respart_set_updatesdf proc~grans_obj_preparesolver->proc~grans_obj_computesolidvf interface~grans_obj_writeoutputdata grans_obj%grans_obj_WriteOutputData proc~grans_obj_preparesolver->interface~grans_obj_writeoutputdata proc~grans_obj_writeoutputdata grans_obj_WriteOutputData proc~grans_obj_writeoutputdata->proc~respart_set_updatesdf proc~grans_obj_writeoutputdata->proc~grans_obj_computesolidvf interface~cdifs_obj_preparesolver cdifs_obj%cdifs_obj_PrepareSolver interface~cdifs_obj_preparesolver->proc~cdifs_obj_preparesolver interface~grans_obj_preparesolver grans_obj%grans_obj_PrepareSolver interface~grans_obj_preparesolver->proc~grans_obj_preparesolver interface~grans_obj_writeoutputdata->proc~grans_obj_writeoutputdata proc~cdifs_obj_advancesolution cdifs_obj_AdvanceSolution proc~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolutionib proc~cdifs_obj_advancesolutionrp cdifs_obj_AdvanceSolutionRP proc~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolutionrp proc~cdifs_obj_advancesolutionrp->proc~cdifs_obj_computesolidvf interface~cdifs_obj_advancesolution cdifs_obj%cdifs_obj_AdvanceSolution interface~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolution

Source Code

    impure subroutine eulerian_obj_AddUpGhostCells_y(this)
      !> Adds up ghostcells in the y 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(2)-2*ngc.le.ngc) then
          select type (this)
          type is (eulerian_obj_r)
            do i=1,ngc
              this%cell(:,hi(2),:) = this%cell(:,hi(2),:) + this%cell(:,hi(2)+i,:)
              this%cell(:,lo(2),:) = this%cell(:,lo(2),:) + this%cell(:,lo(2)-i,:)
            end do
          type is (eulerian_obj_i)
            do i=1,ngc
              this%cell(:,hi(2),:) = this%cell(:,hi(2),:) + this%cell(:,hi(2)+i,:)
              this%cell(:,lo(2),:) = this%cell(:,lo(2),:) + this%cell(:,lo(2)-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)=hi(2)+1    ; sendR(3)=lo(3)-ngc

        select type (this)
        type is (eulerian_obj_r)
          ! Allocate receive buffers
          allocate(buffR_r(Ng(1),ngc,Ng(3)))
          allocate(buffL_r(Ng(1),ngc,Ng(3)))

          ! Post receives from Left and right ranks
          call MPI_IRECV( buffL_r(1,1,1), ngc*Ng(1)*Ng(3), mpi%REAL_WP, &
            mpi%rank%L(2)-1, 0, mpi%comm%g, requests(1), ierr)
          call MPI_IRECV( buffR_r(1,1,1), ngc*Ng(1)*Ng(3), mpi%REAL_WP, &
            mpi%rank%R(2)-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(2), &
            mpi%rank%R(2)-1, 0, mpi%comm%g, requests(4), ierr)
          call MPI_ISEND( this%cell(sendL(1),sendL(2),sendL(3)), 1, gc_slab_r(2), &
            mpi%rank%L(2)-1, 0, mpi%comm%g, requests(3), ierr)
        type is (eulerian_obj_i)
          ! Allocate receive buffers
          allocate(buffR_i(Ng(1),ngc,Ng(3)))
          allocate(buffL_i(Ng(1),ngc,Ng(3)))

          ! Post receives from Left and right ranks
          call MPI_IRECV( buffL_i(1,1,1), ngc*Ng(1)*Ng(3), mpi%INTEGER, &
            mpi%rank%L(2)-1, 0, mpi%comm%g, requests(1), ierr)
          call MPI_IRECV( buffR_i(1,1,1), ngc*Ng(1)*Ng(3), mpi%INTEGER, &
            mpi%rank%R(2)-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(2), &
            mpi%rank%R(2)-1, 0, mpi%comm%g, requests(4), ierr)
          call MPI_ISEND( this%cell(sendL(1),sendL(2),sendL(3)), 1, gc_slab_i(2), &
            mpi%rank%L(2)-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(2)-1.ne.MPI_PROC_NULL) then
            do k=lo(3)-ngc,hi(3)+ngc
              do j=lo(2),lo(2)+ngc-1
                do i=lo(1)-ngc,hi(1)+ngc
                  ii=i-(lo(1)-ngc  )+1
                  jj=j-(lo(2)      )+1
                  kk=k-(lo(3)-ngc  )+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(2)-1.ne.MPI_PROC_NULL) then
            do k=lo(3)-ngc,hi(3)+ngc
              do j=hi(2)-ngc+1,hi(2)
                do i=lo(1)-ngc,hi(1)+ngc
                  ii=i-(lo(1)-ngc  )+1
                  jj=j-(hi(2)-ngc+1)+1
                  kk=k-(lo(3)-ngc  )+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(2)-1.ne.MPI_PROC_NULL) then
            do k=lo(3)-ngc,hi(3)+ngc
              do j=lo(2),lo(2)+ngc-1
                do i=lo(1)-ngc,hi(1)+ngc
                  ii=i-(lo(1)-ngc  )+1
                  jj=j-(lo(2)      )+1
                  kk=k-(lo(3)-ngc  )+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(2)-1.ne.MPI_PROC_NULL) then
            do k=lo(3)-ngc,hi(3)+ngc
              do j=hi(2)-ngc+1,hi(2)
                do i=lo(1)-ngc,hi(1)+ngc
                  ii=i-(lo(1)-ngc  )+1
                  jj=j-(hi(2)-ngc+1)+1
                  kk=k-(lo(3)-ngc  )+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_y