block_obj_Partition Subroutine

private impure subroutine block_obj_Partition(this, Nb)

Partitions a parent block into sub-blocks based on a given decomposition Nb(3) for parallel simulations. This will also define the partition axes, and update local bounds/extents, and global bounds/extents.

Type Bound

block_obj

Arguments

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

A block object

integer, intent(in) :: Nb(3)

Number of blocks in each direction


Calls

proc~~block_obj_partition~~CallsGraph proc~block_obj_partition block_obj%block_obj_Partition proc~axis_obj_final axis_obj%axis_obj_Final proc~block_obj_partition->proc~axis_obj_final proc~axis_obj_init axis_obj%axis_obj_Init proc~block_obj_partition->proc~axis_obj_init proc~block_obj_setconveniencepointers block_obj%block_obj_SetConveniencePointers proc~block_obj_partition->proc~block_obj_setconveniencepointers proc~block_obj_setupmpitypes block_obj%block_obj_SetupMPITypes proc~block_obj_partition->proc~block_obj_setupmpitypes proc~block_obj_subdivideblock block_obj%block_obj_SubDivideBlock proc~block_obj_partition->proc~block_obj_subdivideblock proc~block_obj_updategridghostcells block_obj%block_obj_UpdateGridGhostCells proc~block_obj_partition->proc~block_obj_updategridghostcells proc~block_obj_updatemidpoints block_obj%block_obj_UpdateMidPoints proc~block_obj_partition->proc~block_obj_updatemidpoints proc~block_obj_updatespacing block_obj%block_obj_UpdateSpacing proc~block_obj_partition->proc~block_obj_updatespacing proc~parallel_obj_topology parallel_obj%parallel_obj_Topology proc~block_obj_partition->proc~parallel_obj_topology mpi_type_commit mpi_type_commit proc~block_obj_setupmpitypes->mpi_type_commit mpi_type_free mpi_type_free proc~block_obj_setupmpitypes->mpi_type_free mpi_type_vector mpi_type_vector proc~block_obj_setupmpitypes->mpi_type_vector mpi_irecv mpi_irecv proc~block_obj_updategridghostcells->mpi_irecv mpi_isend mpi_isend proc~block_obj_updategridghostcells->mpi_isend mpi_wait mpi_wait proc~block_obj_updategridghostcells->mpi_wait mpi_waitall mpi_waitall proc~block_obj_updategridghostcells->mpi_waitall proc~block_obj_updateextents block_obj%block_obj_UpdateExtents proc~block_obj_updategridghostcells->proc~block_obj_updateextents mpi_cart_coords mpi_cart_coords proc~parallel_obj_topology->mpi_cart_coords mpi_cart_create mpi_cart_create proc~parallel_obj_topology->mpi_cart_create mpi_cart_rank mpi_cart_rank proc~parallel_obj_topology->mpi_cart_rank mpi_cart_shift mpi_cart_shift proc~parallel_obj_topology->mpi_cart_shift mpi_comm_rank mpi_comm_rank proc~parallel_obj_topology->mpi_comm_rank mpi_dims_create mpi_dims_create proc~parallel_obj_topology->mpi_dims_create proc~block_obj_updateextents->proc~axis_obj_init

Called by

proc~~block_obj_partition~~CalledByGraph proc~block_obj_partition block_obj%block_obj_Partition proc~cdifs_obj_preparesolverblock cdifs_obj_PrepareSolverBlock proc~cdifs_obj_preparesolverblock->proc~block_obj_partition proc~grans_obj_preparesolverblock grans_obj_PrepareSolverBlock proc~grans_obj_preparesolverblock->proc~block_obj_partition proc~cdifs_obj_preparesolver cdifs_obj_PrepareSolver proc~cdifs_obj_preparesolver->proc~cdifs_obj_preparesolverblock proc~grans_obj_preparesolver grans_obj_PrepareSolver proc~grans_obj_preparesolver->proc~grans_obj_preparesolverblock 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

Source Code

    impure subroutine block_obj_Partition(this,Nb)
      !> Partitions a parent block into sub-blocks
      ! based on a given decomposition Nb(3) for parallel simulations.
      ! This will also define the partition axes, and update
      ! local bounds/extents, and global bounds/extents.
      implicit none
      class(block_obj), intent(inout) :: this                                  !! A block object
      integer,          intent(in)    :: Nb(3)                                 !! Number of blocks in each direction
      ! Work variables
      integer :: i,dir
      integer :: minlo(3)
      integer :: maxhi(3)
      type(axis_obj):: axis_global(3)

      if (.not.associated(this%parallel)) return

      ! Store global axis prior to subdividing
      do dir=1,3
        call axis_global(dir)%Initialize(this%lo(dir),this%hi(dir),this%ngc)
        axis_global(dir)%x = this%axis(dir)%x
      end do
      minlo = this%lo
      maxhi = this%hi

      ! Switch flag on
      this%is_partitioned=.true.

      ! Total number of blocks must match total number of MPI ranks
      if (this%parallel%nproc.ne.Nb(1)*Nb(2)*Nb(3)) &
        call this%parallel%Stop('Unable to partition block due to mismatch in parallel decomposition')

      ! Generate an MPI cartesian topology
      call this%parallel%Topology(this%periods,maxhi-minlo+1,Nb)

      ! Get bounds of sub-block ownded by this MPI rank
      call this%SubDivideBlock(this%parallel%rank%dir,Nb,minlo,maxhi,this%lo,this%hi)

      ! Initialize axes
      do dir=1,3
        call this%axis(dir)%Finalize()
        call this%axis(dir)%Initialize(this%lo(dir),this%hi(dir),this%ngc)
      end do

      ! Associate pointers
      call this%SetConveniencePointers

      associate (lo => this%lo, hi => this%hi, axis=>this%axis)
        ! Get local axes coordinates from global axes
        do dir=1,3
          do i=lo(dir),hi(dir)+1
            axis(dir)%x(i) = axis_global(dir)%x(i)
          end do
        end do
      end associate

      ! Remove global arrays
      do dir=1,3
        call axis_global(dir)%Finalize
      end do

      ! Update values within ghostcells
      call this%UpdateGridGhostCells()

      ! Update mid points (xm)
      call this%UpdateMidPoints()

      ! Update spacing (dxm)
      call this%UpdateSpacing()

      ! Create MPI type for ghostcell communication
      call this%SetupMPITypes()

      return
    end subroutine block_obj_Partition