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.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(grans_obj), | intent(inout) | :: | this |
GRANS solver |
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