grans_obj_PrepareSolverCollision Subroutine

subroutine grans_obj_PrepareSolverCollision(this)

Prepares the collision utility used by GRANS. If the collision grid spacing is not provided, this subroutine makes it equal to the largest particle diameter. Collisions with IBs are default to true, but can be turned off. Users can use collisions with walls instead (for simple geometries). Additional parameters are read in the Prepare subroutine of the collision object.

Arguments

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

GRANS solver


Calls

proc~~grans_obj_preparesolvercollision~~CallsGraph proc~grans_obj_preparesolvercollision grans_obj_PrepareSolverCollision center center proc~grans_obj_preparesolvercollision->center none~add~4 collision_obj%Add proc~grans_obj_preparesolvercollision->none~add~4 none~get~4 parser_obj%Get proc~grans_obj_preparesolvercollision->none~get~4 particle particle proc~grans_obj_preparesolvercollision->particle proc~collision_obj_init collision_obj%collision_obj_Init proc~grans_obj_preparesolvercollision->proc~collision_obj_init proc~collision_obj_prepare collision_obj%collision_obj_Prepare proc~grans_obj_preparesolvercollision->proc~collision_obj_prepare proc~collision_obj_setupcollisionblock collision_obj%collision_obj_SetupCollisionBlock proc~grans_obj_preparesolvercollision->proc~collision_obj_setupcollisionblock proc~parser_obj_isdefined parser_obj%parser_obj_IsDefined proc~grans_obj_preparesolvercollision->proc~parser_obj_isdefined proc~collision_obj_addimmersedboundaries collision_obj%collision_obj_AddImmersedBoundaries none~add~4->proc~collision_obj_addimmersedboundaries proc~collision_obj_addpointparticles collision_obj%collision_obj_AddPointParticles none~add~4->proc~collision_obj_addpointparticles proc~collision_obj_addresolvedparticles collision_obj%collision_obj_AddResolvedParticles none~add~4->proc~collision_obj_addresolvedparticles proc~parser_obj_read0d parser_obj%parser_obj_read0D none~get~4->proc~parser_obj_read0d proc~parser_obj_read1d parser_obj%parser_obj_read1D none~get~4->proc~parser_obj_read1d proc~collision_obj_prepare->none~get~4 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~parser_obj_fetchlabelid parser_obj%parser_obj_FetchLabelID proc~parser_obj_isdefined->proc~parser_obj_fetchlabelid 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~parser_obj_read0d->proc~parser_obj_fetchlabelid none~assigndefault parser_obj%AssignDefault proc~parser_obj_read0d->none~assigndefault proc~parser_obj_read1d->proc~parser_obj_fetchlabelid proc~parser_obj_read1d->none~assigndefault proc~parser_obj_assigndefault0d parser_obj%parser_obj_AssignDefault0D none~assigndefault->proc~parser_obj_assigndefault0d proc~parser_obj_assigndefault1d parser_obj%parser_obj_AssignDefault1D none~assigndefault->proc~parser_obj_assigndefault1d 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~~grans_obj_preparesolvercollision~~CalledByGraph proc~grans_obj_preparesolvercollision grans_obj_PrepareSolverCollision 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

    subroutine grans_obj_PrepareSolverCollision(this)
      !> Prepares the collision utility used by GRANS.
      ! If the collision grid spacing is not provided, this subroutine
      ! makes it equal to the largest particle diameter.
      ! Collisions with IBs are default to true, but can be turned off.
      ! Users can use collisions with walls instead (for simple geometries).
      ! Additional parameters are read in the Prepare subroutine of the
      ! collision object.
      implicit none
      class(grans_obj), intent(inout) :: this                                  !! GRANS solver
      ! Work variables
      real(wp):: ds
      real(wp):: buf
      integer :: n
      logical :: use_IB_col

      call this%parser%Get('Collision grid spacing', ds,          default = 0.0_wp )
      call this%parser%Get('Collision with IB',      use_IB_col,  default = .true. )

      ! Determine spacing for collision blocks
      ! - Default value
      if (.not.this%parser%IsDefined("Collision grid spacing")) then
        ds = minval(this%block%dx)
        if (this%use_RP) then
          select type(center => this%RP%p)
          type is (ResPart_obj)
            do n=1,this%RP%count_
              ds = max(ds,center(n)%d)
            end do
          end select
        end if
        if (this%use_PP) then
          select type(particle => this%PP%p)
          type is (particle_obj)
            do n=1,this%PP%count_
              ds = max(ds,particle(n)%d)
            end do
          end select
        end if
        call this%parallel%Max(ds, buf); ds=buf
      end if

      ! Prepare structure
      call this%collisions%Initialize(this%parallel)
      call this%collisions%Prepare(this%timer,this%parser,this%monitors)
      if (this%use_IB.and.use_IB_col) &
                       call this%collisions%Add(this%IB)
      if (this%use_RP) call this%collisions%Add(this%RP)
      if (this%use_PP) call this%collisions%Add(this%PP)

      call this%collisions%SetupCollisionBlock(ds,1)

      return
    end subroutine grans_obj_PrepareSolverCollision