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 | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(block_obj), | intent(inout) | :: | this |
A block object |
||
| integer, | intent(in) | :: | Nb(3) |
Number of blocks in each direction |
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