spongeLayer_box2d_sharpCornerPolyn6_for_coord Function

private function spongeLayer_box2d_sharpCornerPolyn6_for_coord(me, coord, n) result(res)

This function calculates the sigma for the 2d box shape sponge layer from coord for polynomial n6. Sponge profile: where, \sigma_A - sponge strength, L - thickness, x0 - start of sponge.

Arguments

Type IntentOptional Attributes Name
type(tem_spongeLayer_box_type) :: me

Spatial sponge layer to evaluate

real(kind=rk), intent(in) :: coord(n,3)

barycentric Ids of an elements. 1st index goes over number of elements and 2nd index goes over x,y,z coordinates

integer, intent(in) :: n

Number of arrays to return

Return Value real(kind=rk), (n)

return value


Called by

proc~~spongelayer_box2d_sharpcornerpolyn6_for_coord~~CalledByGraph proc~spongelayer_box2d_sharpcornerpolyn6_for_coord spongeLayer_box2d_sharpCornerPolyn6_for_coord proc~spongelayer_box2d_scalar_for_coord spongelayer_box2d_scalar_for_coord proc~spongelayer_box2d_scalar_for_coord->proc~spongelayer_box2d_sharpcornerpolyn6_for_coord proc~spongelayer_box2d_sharpcornerpolyn6_for_treeids spongelayer_box2d_sharpCornerPolyn6_for_treeIDs proc~spongelayer_box2d_sharpcornerpolyn6_for_treeids->proc~spongelayer_box2d_sharpcornerpolyn6_for_coord proc~spongelayer_box2d_vector_for_coord spongelayer_box2d_vector_for_coord proc~spongelayer_box2d_vector_for_coord->proc~spongelayer_box2d_scalar_for_coord proc~spongelayer_box2d_scalar_for_treeids spongelayer_box2d_scalar_for_treeIDs proc~spongelayer_box2d_scalar_for_treeids->proc~spongelayer_box2d_sharpcornerpolyn6_for_treeids interface~tem_spongelayer_box2d_for tem_spongeLayer_box2d_for interface~tem_spongelayer_box2d_for->proc~spongelayer_box2d_scalar_for_coord interface~tem_spongelayer_box2d_for->proc~spongelayer_box2d_vector_for_coord interface~tem_spongelayer_box2d_for->proc~spongelayer_box2d_scalar_for_treeids proc~spongelayer_box2d_vector_for_treeids spongelayer_box2d_vector_for_treeIDs interface~tem_spongelayer_box2d_for->proc~spongelayer_box2d_vector_for_treeids proc~tem_spatial_for_treeids tem_spatial_for_treeIDs proc~tem_spatial_for_treeids->interface~tem_spongelayer_box2d_for proc~tem_spatial_vector_for_coord tem_spatial_vector_for_coord proc~tem_spatial_vector_for_coord->interface~tem_spongelayer_box2d_for proc~spongelayer_box2d_vector_for_treeids->proc~spongelayer_box2d_scalar_for_treeids proc~tem_spatial_vector_for_treeids tem_spatial_vector_for_treeIDs proc~tem_spatial_vector_for_treeids->interface~tem_spongelayer_box2d_for proc~tem_spatial_for_coord tem_spatial_for_coord proc~tem_spatial_for_coord->interface~tem_spongelayer_box2d_for interface~tem_spatial_for tem_spatial_for interface~tem_spatial_for->proc~tem_spatial_for_treeids interface~tem_spatial_for->proc~tem_spatial_vector_for_coord interface~tem_spatial_for->proc~tem_spatial_vector_for_treeids interface~tem_spatial_for->proc~tem_spatial_for_coord proc~tem_spatial_scalar_for_index tem_spatial_scalar_for_index proc~tem_spatial_scalar_for_index->proc~tem_spatial_for_coord proc~tem_spatial_vector_for_index tem_spatial_vector_for_index proc~tem_spatial_vector_for_index->proc~tem_spatial_vector_for_coord

Contents


Source Code

  function spongeLayer_box2d_sharpCornerPolyn6_for_coord(me, coord, n) &
    & result(res)
    ! --------------------------------------------------------------------------
    !> Spatial sponge layer to evaluate
    type(tem_spongeLayer_box_type) :: me
    !> Number of arrays to return
    integer, intent(in) :: n
    !> barycentric Ids of an elements.
    !! 1st index goes over number of elements and
    !! 2nd index goes over x,y,z coordinates
    real(kind=rk), intent( in ) :: coord(n,3)
    !> return value
    real(kind=rk) :: res(n)
    ! --------------------------------------------------------------------------
    integer :: i
    real(kind=rk) :: sigma, origin(2), extent(2), coordLoc(2), normal
    real(kind=rk) :: proj_len1, proj_len2
    real(kind=rk) :: vec_min(2), vec_max(2), vec_minSqr(2), vec_maxSqr(2) 
    real(kind=rk) :: rad, const_fac, box_max(2)
    ! --------------------------------------------------------------------------
    origin(:) = me%origin(1:2)
    extent(:) = me%extent(1:2)
    box_max(:) = origin(:) + extent(:)

    const_fac = 729_rk/(16_rk*me%thickness**6)

    do i = 1,n
      coordLoc = coord(i,1:2)
      vec_min(:) = coordLoc(:) - origin(:)
      vec_max(:) = coordLoc(:) - box_max(:)
      vec_minSqr(:) = vec_min(:)**2
      vec_maxSqr(:) = vec_max(:)**2
      proj_len1 = 0_rk
      proj_len2 = 0_rk
      normal = 1_rk
      rad = 0_rk

      if (coordLoc(1) < origin(1) .and. coordLoc(2) < origin(2)) then
        ! South-West: -x,-y
        rad = sqrt( max(vec_minSqr(1), vec_minSqr(2)) )

      else if (coordLoc(1) < origin(1) .and. coordLoc(2) > box_max(2)) then
        ! North-West: -x, y
        rad = sqrt( max(vec_minSqr(1), vec_maxSqr(2)) )

      else if (coordLoc(1) > box_max(1) .and. coordLoc(2) < origin(2)) then
        ! South-East: x, -y
        rad = sqrt( max(vec_maxSqr(1), vec_minSqr(2)) )

      else if (coordLoc(1) > box_max(1) .and. coordLoc(2) > box_max(2)) then
        ! North-East: x, y
        rad = sqrt( max(vec_maxSqr(1), vec_maxSqr(2)) )

      else if (coordLoc(1) < origin(1)) then ! West: -x
        normal = -1_rk
        rad = vec_min(1)

      else if (coordLoc(2) < origin(2)) then ! South: -y
        normal = -1_rk
        rad = vec_min(2)

      else if (coordLoc(1) > box_max(1)) then ! East: x
        normal = 1_rk
        rad = vec_max(1)

      else if (coordLoc(2) > box_max(2)) then ! North: y
        normal = 1_rk
        rad = vec_max(2)
      end if

      proj_len1 = rad*normal
      proj_len2 = (me%thickness*normal - rad)*normal

      if (proj_len1 > 0) then
        sigma = const_fac * proj_len2**2 * (proj_len1**4)
        res(i) = min(1.0_rk, sigma) * me%dampFactor
      else
        res(i) = 0.0_rk
      end if
    end do

  end function spongeLayer_box2d_sharpCornerPolyn6_for_coord