eulerian_obj_UpdateGhostCells_z Subroutine

private impure subroutine eulerian_obj_UpdateGhostCells_z(this)

Updates the ghostcells in the z 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_updateghostcells_z~~CallsGraph proc~eulerian_obj_updateghostcells_z eulerian_obj_base%eulerian_obj_UpdateGhostCells_z cell cell proc~eulerian_obj_updateghostcells_z->cell mpi_irecv mpi_irecv proc~eulerian_obj_updateghostcells_z->mpi_irecv mpi_isend mpi_isend proc~eulerian_obj_updateghostcells_z->mpi_isend mpi_waitall mpi_waitall proc~eulerian_obj_updateghostcells_z->mpi_waitall

Called by

proc~~eulerian_obj_updateghostcells_z~~CalledByGraph proc~eulerian_obj_updateghostcells_z eulerian_obj_base%eulerian_obj_UpdateGhostCells_z proc~eulerian_obj_updateghostcells eulerian_obj_base%eulerian_obj_UpdateGhostCells proc~eulerian_obj_updateghostcells->proc~eulerian_obj_updateghostcells_z proc~bc_set_buildmask bc_set%bc_set_BuildMask proc~bc_set_buildmask->proc~eulerian_obj_updateghostcells proc~bc_set_updateboundaryscalar bc_set%bc_set_UpdateBoundaryScalar proc~bc_set_updateboundaryscalar->proc~eulerian_obj_updateghostcells proc~bc_set_updateboundaryvector bc_set%bc_set_UpdateBoundaryVector proc~bc_set_updateboundaryvector->proc~eulerian_obj_updateghostcells proc~cdifs_obj_advancesolutioncorrector cdifs_obj_AdvanceSolutionCorrector proc~cdifs_obj_advancesolutioncorrector->proc~eulerian_obj_updateghostcells none~updateboundary bc_set%UpdateBoundary proc~cdifs_obj_advancesolutioncorrector->none~updateboundary proc~hypre_obj_solve hypre_obj%hypre_obj_Solve proc~cdifs_obj_advancesolutioncorrector->proc~hypre_obj_solve proc~cdifs_obj_advancesolutionib cdifs_obj_AdvanceSolutionIB proc~cdifs_obj_advancesolutionib->proc~eulerian_obj_updateghostcells proc~cdifs_obj_advancesolutionib->none~updateboundary proc~respart_set_getibforcing ResPart_set%ResPart_set_GetIBForcing proc~cdifs_obj_advancesolutionib->proc~respart_set_getibforcing proc~respart_set_updatesdf ResPart_set%ResPart_set_UpdateSDF proc~cdifs_obj_advancesolutionib->proc~respart_set_updatesdf proc~cdifs_obj_writeoutputdata cdifs_obj_WriteOutputData proc~cdifs_obj_writeoutputdata->proc~eulerian_obj_updateghostcells proc~eulerian_obj_addupghostcells eulerian_obj_base%eulerian_obj_AddUpGhostCells proc~eulerian_obj_addupghostcells->proc~eulerian_obj_updateghostcells proc~eulerian_set_readsingleh5hut eulerian_set%eulerian_set_ReadSingleH5HUT proc~eulerian_set_readsingleh5hut->proc~eulerian_obj_updateghostcells proc~eulerian_set_readsinglehdf5 eulerian_set%eulerian_set_ReadSingleHDF5 proc~eulerian_set_readsinglehdf5->proc~eulerian_obj_updateghostcells proc~hypre_obj_setuprowsij hypre_obj%hypre_obj_SetupRowsIJ proc~hypre_obj_setuprowsij->proc~eulerian_obj_updateghostcells proc~hypre_obj_solveij hypre_obj%hypre_obj_SolveIJ proc~hypre_obj_solveij->proc~eulerian_obj_updateghostcells proc~marker_set_computesolidvolfrac marker_set%marker_set_ComputeSolidVolFrac proc~marker_set_computesolidvolfrac->proc~eulerian_obj_updateghostcells proc~hypre_obj_setup hypre_obj%hypre_obj_Setup proc~marker_set_computesolidvolfrac->proc~hypre_obj_setup proc~marker_set_computesolidvolfrac->proc~hypre_obj_solve proc~marker_set_filter marker_set%marker_set_Filter proc~marker_set_computesolidvolfrac->proc~marker_set_filter proc~marker_set_updatenormals marker_set%marker_set_UpdateNormals proc~marker_set_updatenormals->proc~eulerian_obj_updateghostcells proc~marker_set_updatenormals->none~updateboundary proc~marker_set_updatenormals->proc~marker_set_filter proc~op_obj_conv11 op_obj%op_obj_conv11 proc~op_obj_conv11->proc~eulerian_obj_updateghostcells proc~op_obj_conv12 op_obj%op_obj_conv12 proc~op_obj_conv12->proc~eulerian_obj_updateghostcells proc~op_obj_conv13 op_obj%op_obj_conv13 proc~op_obj_conv13->proc~eulerian_obj_updateghostcells proc~op_obj_conv21 op_obj%op_obj_conv21 proc~op_obj_conv21->proc~eulerian_obj_updateghostcells proc~op_obj_conv22 op_obj%op_obj_conv22 proc~op_obj_conv22->proc~eulerian_obj_updateghostcells proc~op_obj_conv23 op_obj%op_obj_conv23 proc~op_obj_conv23->proc~eulerian_obj_updateghostcells proc~op_obj_conv31 op_obj%op_obj_conv31 proc~op_obj_conv31->proc~eulerian_obj_updateghostcells proc~op_obj_conv32 op_obj%op_obj_conv32 proc~op_obj_conv32->proc~eulerian_obj_updateghostcells proc~op_obj_conv33 op_obj%op_obj_conv33 proc~op_obj_conv33->proc~eulerian_obj_updateghostcells proc~op_obj_strainrate op_obj%op_obj_StrainRate proc~op_obj_strainrate->proc~eulerian_obj_updateghostcells interface~cdifs_obj_writeoutputdata cdifs_obj%cdifs_obj_WriteOutputData interface~cdifs_obj_writeoutputdata->proc~cdifs_obj_writeoutputdata none~readsingle eulerian_set%ReadSingle none~readsingle->proc~eulerian_set_readsingleh5hut none~readsingle->proc~eulerian_set_readsinglehdf5 none~updateboundary->proc~bc_set_updateboundaryscalar none~updateboundary->proc~bc_set_updateboundaryvector proc~cdifs_obj_advancesolution cdifs_obj_AdvanceSolution proc~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolutioncorrector proc~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolutionib proc~cdifs_obj_advancesolutionpredictor cdifs_obj_AdvanceSolutionPredictor proc~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolutionpredictor proc~cdifs_obj_advancesolutionrp cdifs_obj_AdvanceSolutionRP proc~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolutionrp proc~cdifs_obj_preparesolverbcs cdifs_obj_PrepareSolverBCS proc~cdifs_obj_preparesolverbcs->proc~bc_set_buildmask proc~cdifs_obj_preparesolverbcs->none~updateboundary proc~hypre_obj_setup->proc~hypre_obj_setuprowsij proc~hypre_obj_solve->proc~hypre_obj_solveij 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~respart_set_getsurfacestresses ResPart_set%ResPart_set_GetSurfaceStresses proc~respart_set_getsurfacestresses->proc~op_obj_strainrate proc~respart_set_updatenormals ResPart_set%ResPart_set_UpdateNormals proc~respart_set_updatenormals->proc~marker_set_updatenormals interface~cdifs_obj_advancesolution cdifs_obj%cdifs_obj_AdvanceSolution interface~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolution proc~cdifs_obj_advancesolutionpredictor->none~updateboundary proc~cdifs_obj_computesolidvf cdifs_obj%cdifs_obj_ComputeSolidVF proc~cdifs_obj_computesolidvf->none~updateboundary proc~cdifs_obj_computesolidvf->proc~hypre_obj_solve proc~cdifs_obj_computesolidvf->proc~respart_set_updatenormals proc~cdifs_obj_preparesolver cdifs_obj_PrepareSolver proc~cdifs_obj_preparesolver->interface~cdifs_obj_writeoutputdata proc~cdifs_obj_preparesolver->none~updateboundary proc~cdifs_obj_preparesolver->proc~cdifs_obj_preparesolverbcs proc~cdifs_obj_preparesolver->proc~cdifs_obj_computesolidvf proc~cdifs_obj_preparesolveroperators cdifs_obj_PrepareSolverOperators proc~cdifs_obj_preparesolver->proc~cdifs_obj_preparesolveroperators proc~cdifs_obj_preparesolver->proc~respart_set_updatesdf proc~cdifs_obj_preparesolveroperatorsplap cdifs_obj_PrepareSolverOperatorsPLAP proc~cdifs_obj_preparesolveroperatorsplap->proc~hypre_obj_setup proc~cdifs_obj_preparesolveroperatorsvflap cdifs_obj_PrepareSolverOperatorsVFLAP proc~cdifs_obj_preparesolveroperatorsvflap->proc~hypre_obj_setup proc~eulerian_set_readh5hut eulerian_set%eulerian_set_ReadH5HUT proc~eulerian_set_readh5hut->none~readsingle proc~eulerian_set_readhdf5 eulerian_set%eulerian_set_ReadHDF5 proc~eulerian_set_readhdf5->none~readsingle proc~grans_obj_computesolidvf grans_obj%grans_obj_ComputeSolidVF proc~grans_obj_computesolidvf->none~updateboundary proc~grans_obj_computesolidvf->proc~hypre_obj_solve proc~grans_obj_computesolidvf->proc~respart_set_updatenormals proc~grans_obj_preparesolver grans_obj_PrepareSolver proc~grans_obj_preparesolver->none~updateboundary proc~grans_obj_preparesolver->proc~grans_obj_computesolidvf proc~grans_obj_preparesolveroperators grans_obj_PrepareSolverOperators proc~grans_obj_preparesolver->proc~grans_obj_preparesolveroperators interface~grans_obj_writeoutputdata grans_obj%grans_obj_WriteOutputData proc~grans_obj_preparesolver->interface~grans_obj_writeoutputdata proc~grans_obj_preparesolver->proc~respart_set_updatesdf proc~grans_obj_preparesolveroperators->proc~hypre_obj_setup proc~grans_obj_writeoutputdata grans_obj_WriteOutputData proc~grans_obj_writeoutputdata->none~updateboundary proc~grans_obj_writeoutputdata->proc~grans_obj_computesolidvf proc~grans_obj_writeoutputdata->proc~respart_set_updatesdf proc~marker_set_getibforcing marker_set%marker_set_GetIBForcing proc~marker_set_getibforcing->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~respart_set_gethydroforces ResPart_set%ResPart_set_GetHydroForces proc~respart_set_gethydroforces->proc~respart_set_getsurfacestresses proc~solid_set_filter solid_set%solid_set_Filter proc~solid_set_filter->proc~marker_set_filter 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_advancesolutionrp->proc~cdifs_obj_computesolidvf proc~cdifs_obj_advancesolutionrp->proc~respart_set_gethydroforces proc~cdifs_obj_preparesolveroperators->proc~cdifs_obj_preparesolveroperatorsplap proc~cdifs_obj_preparesolveroperators->proc~cdifs_obj_preparesolveroperatorsvflap proc~respart_set_getibforcing->proc~marker_set_getibforcing proc~respart_set_updatesdf->proc~marker_set_updatesdf

Source Code

    impure subroutine eulerian_obj_UpdateGhostCells_z(this)
      !> Updates the 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 :: recvL(3),recvR(3)
      integer :: sendL(3),sendR(3)
      integer :: i
      type(MPI_STATUS)  :: statuses(4)
      type(MPI_REQUEST) :: requests(4)
      integer :: ierr

      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/1D and pseudo-2D/1D exceptions:
        ! more ghost cells than internal cells
        if (Ng(3)-2*ngc.le.ngc) then
          select type (this)
          type is (eulerian_obj_r)
            do i=1,ngc
              this%cell(:,:,hi(3)+i) = this%cell(:,:,lo(3)+i-1)
              this%cell(:,:,lo(3)-i) = this%cell(:,:,hi(3)-i+1)
            end do
          type is (eulerian_obj_i)
            do i=1,ngc
              this%cell(:,:,hi(3)+i) = this%cell(:,:,lo(3)+i-1)
              this%cell(:,:,lo(3)-i) = this%cell(:,:,hi(3)-i+1)
            end do
          end select

          return
        end if

        ! Address of first element in buffer
        recvL(1)=lo(1)-ngc  ; recvL(2)=lo(2)-ngc  ; recvL(3)=lo(3)-ngc
        recvR(1)=lo(1)-ngc  ; recvR(2)=lo(2)-ngc  ; recvR(3)=hi(3)+1
        sendL(1)=lo(1)-ngc  ; sendL(2)=lo(2)-ngc  ; sendL(3)=lo(3)
        sendR(1)=lo(1)-ngc  ; sendR(2)=lo(2)-ngc  ; sendR(3)=hi(3)-ngc+1

        select type (this)
        type is (eulerian_obj_r)
          ! Post receives from Left and right ranks
          call MPI_IRECV( this%cell(recvL(1),recvL(2),recvL(3)), 1, gc_slab_r(3), &
            mpi%rank%L(3)-1, 0, mpi%comm%g, requests(1), ierr)
          call MPI_IRECV( this%cell(recvR(1),recvR(2),recvR(3)), 1, gc_slab_r(3), &
            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)
          ! Post receives from Left and right ranks
          call MPI_IRECV( this%cell(recvL(1),recvL(2),recvL(3)), 1, gc_slab_i(3), &
            mpi%rank%L(3)-1, 0, mpi%comm%g, requests(1), ierr)
          call MPI_IRECV( this%cell(recvR(1),recvR(2),recvR(3)), 1, gc_slab_i(3), &
            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 )

      end associate

      return
    end subroutine eulerian_obj_UpdateGhostCells_z