Initializes cblock to handle collisions. This extra block is expected to be coarser than the simulation block, but larger than the maximum object size. It is used to expedite neighbor searches.
With this subroutine, each MPI rank owns a LOCAL view of the collision block.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(collision_obj), | intent(inout) | :: | this |
Collision utility |
||
| real(kind=wp), | intent(in) | :: | ds |
Target grid spacing |
||
| integer, | intent(in) | :: | ngc |
Number of ghost cells for collision block |
||
| type(block_obj), | intent(in), | optional | :: | block |
Optional block to conform to |
impure subroutine collision_obj_SetupCollisionBlock2(this,ds,ngc,block) !> Initializes cblock to handle collisions. This ! extra block is expected to be coarser than the ! simulation block, but larger than the maximum ! object size. It is used to expedite neighbor ! searches. ! ! With this subroutine, each MPI rank owns a ! LOCAL view of the collision block. use leapBlock, only : axis_obj implicit none class(collision_obj), intent(inout) :: this !! Collision utility real(wp), intent(in) :: ds !! Target grid spacing integer, intent(in) :: ngc !! Number of ghost cells for collision block type(block_obj), intent(in), & optional :: block !! Optional block to conform to ! work variables logical :: periods(3) integer :: lo(3) integer :: hi(3) real(wp):: xhi(3) real(wp):: xlo(3) real(wp):: local_xhi(3) real(wp):: local_xlo(3) integer :: Nb(3) integer :: dir,i integer :: minlo(3) integer :: maxhi(3) type(axis_obj):: axis_global(3) ! Figure out domain extents and periodicity if (present(block)) then xlo = block%pmin xhi = block%pmax periods = block%periods local_xlo = block%xlo local_xhi = block%xhi else if (associated(this%RP)) then xlo = this%RP%block%pmin xhi = this%RP%block%pmax local_xlo = this%RP%block%xlo local_xhi = this%RP%block%xhi periods = this%RP%block%periods else if (associated(this%PP)) then xlo = this%PP%block%pmin xhi = this%PP%block%pmax periods = this%PP%block%periods local_xlo = this%PP%block%xlo local_xhi = this%PP%block%xhi else call this%parallel%Stop("Cannot setup collision block without associated data") end if ! Number of grid points lo = [1,1,1] hi = floor((xhi-xlo)/ds) ! Create a uniform grid within each block call this%cblock%Initialize(ngc,this%parallel) ! Impose same periodicity as main simulation block call this%cblock%SetPeriodicity(periods) ! Create a uniform grid call this%cblock%SetupUniformGrid(xlo,xhi,lo,hi) ! Partition to match simulation block as closely as possible ! Stoore global indices minlo = lo maxhi = hi ! Store global axis prior to subdividing do dir=1,3 call axis_global(dir)%Initialize(minlo(dir),maxhi(dir),ngc) axis_global(dir)%x = this%cblock%axis(dir)%x end do ! Switch flag on this%cblock%is_partitioned=.true. ! Decomposition Nb = this%parallel%np ! Get bounds of sub-block ownded by this MPI rank this%cblock%hi = min(max(this%cblock%GetOwningCellWGC(local_xhi), minlo),maxhi) this%cblock%lo = min(max(this%cblock%GetOwningCellWGC(local_xlo) + 1, minlo),maxhi) do dir=1,3 if (this%parallel%rank%dir(dir).eq.1) then this%cblock%lo(dir) = lo(dir) end if end do ! Redifine hi/lo to be the local ones lo = this%cblock%lo hi = this%cblock%hi ! Reinitialize axes do dir=1,3 call this%cblock%axis(dir)%Finalize() call this%cblock%axis(dir)%Initialize(lo(dir),hi(dir),ngc) end do ! Associate pointers call this%cblock%SetConveniencePointers ! Get local axes coordinates from global axes do dir=1,3 do i=lo(dir),hi(dir)+1 this%cblock%axis(dir)%x(i) = axis_global(dir)%x(i) end do end do ! Remove global arrays do dir=1,3 call axis_global(dir)%Finalize end do ! Update values within ghostcells call this%cblock%UpdateGridGhostCells() ! Update mid points (xm) call this%cblock%UpdateMidPoints() ! Update spacing (dxm) call this%cblock%UpdateSpacing() ! Create MPI type for ghostcell communication call this%cblock%SetupMPITypes() ! Allocate neighbor list if (associated(this%RP)) then ! Allocate neighbor list allocate(this%RPneighbors(lo(1)-ngc:hi(1)+ngc,lo(2)-ngc:hi(2)+ngc,lo(3)-ngc:hi(3)+ngc)) allocate(this%RPobjincell(lo(1)-ngc:hi(1)+ngc,lo(2)-ngc:hi(2)+ngc,lo(3)-ngc:hi(3)+ngc)) end if if (associated(this%IB)) then ! Allocate neighbor list allocate(this%IBneighbors(lo(1)-ngc:hi(1)+ngc,lo(2)-ngc:hi(2)+ngc,lo(3)-ngc:hi(3)+ngc)) allocate(this%IBobjincell(lo(1)-ngc:hi(1)+ngc,lo(2)-ngc:hi(2)+ngc,lo(3)-ngc:hi(3)+ngc)) end if if (associated(this%PP)) then ! Allocate neighbor list allocate(this%PPneighbors(lo(1)-ngc:hi(1)+ngc,lo(2)-ngc:hi(2)+ngc,lo(3)-ngc:hi(3)+ngc)) allocate(this%PPobjincell(lo(1)-ngc:hi(1)+ngc,lo(2)-ngc:hi(2)+ngc,lo(3)-ngc:hi(3)+ngc)) end if return end subroutine collision_obj_SetupCollisionBlock2