collision_obj_SetupCollisionBlock Subroutine

private impure subroutine collision_obj_SetupCollisionBlock(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 complete and global view of the collision block.

Type Bound

collision_obj

Arguments

Type IntentOptional 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


Calls

proc~~collision_obj_setupcollisionblock~~CallsGraph proc~collision_obj_setupcollisionblock collision_obj%collision_obj_SetupCollisionBlock none~initialize~6 block_obj%Initialize proc~collision_obj_setupcollisionblock->none~initialize~6 proc~block_obj_setperiodicity block_obj%block_obj_SetPeriodicity proc~collision_obj_setupcollisionblock->proc~block_obj_setperiodicity proc~block_obj_setupuniformgrid block_obj%block_obj_SetupUniformGrid proc~collision_obj_setupcollisionblock->proc~block_obj_setupuniformgrid proc~block_obj_init block_obj%block_obj_Init none~initialize~6->proc~block_obj_init proc~block_obj_init2 block_obj%block_obj_Init2 none~initialize~6->proc~block_obj_init2 proc~axis_obj_init axis_obj%axis_obj_Init proc~block_obj_setupuniformgrid->proc~axis_obj_init proc~block_obj_setconveniencepointers block_obj%block_obj_SetConveniencePointers proc~block_obj_setupuniformgrid->proc~block_obj_setconveniencepointers proc~block_obj_setupmpitypes block_obj%block_obj_SetupMPITypes proc~block_obj_setupuniformgrid->proc~block_obj_setupmpitypes proc~block_obj_updategridghostcells block_obj%block_obj_UpdateGridGhostCells proc~block_obj_setupuniformgrid->proc~block_obj_updategridghostcells proc~block_obj_updatemidpoints block_obj%block_obj_UpdateMidPoints proc~block_obj_setupuniformgrid->proc~block_obj_updatemidpoints proc~block_obj_updatespacing block_obj%block_obj_UpdateSpacing proc~block_obj_setupuniformgrid->proc~block_obj_updatespacing proc~block_obj_init2->proc~block_obj_setupuniformgrid 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 proc~block_obj_updateextents->proc~axis_obj_init

Called by

proc~~collision_obj_setupcollisionblock~~CalledByGraph proc~collision_obj_setupcollisionblock collision_obj%collision_obj_SetupCollisionBlock proc~cdifs_obj_preparesolver cdifs_obj_PrepareSolver proc~cdifs_obj_preparesolver->proc~collision_obj_setupcollisionblock proc~grans_obj_preparesolvercollision grans_obj_PrepareSolverCollision proc~grans_obj_preparesolvercollision->proc~collision_obj_setupcollisionblock interface~cdifs_obj_preparesolver cdifs_obj%cdifs_obj_PrepareSolver interface~cdifs_obj_preparesolver->proc~cdifs_obj_preparesolver proc~grans_obj_preparesolver grans_obj_PrepareSolver proc~grans_obj_preparesolver->proc~grans_obj_preparesolvercollision interface~grans_obj_preparesolver grans_obj%grans_obj_PrepareSolver interface~grans_obj_preparesolver->proc~grans_obj_preparesolver

Source Code

    impure subroutine collision_obj_SetupCollisionBlock(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
      ! complete and global view of the collision block.
      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)

      if (present(block)) then
        xlo = block%pmin
        xhi = block%pmax
        periods = block%periods
      else if (associated(this%RP)) then
        xlo = this%RP%block%pmin
        xhi = this%RP%block%pmax
        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
      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)

      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_SetupCollisionBlock