advecc2nd_corr_liberal Subroutine

public subroutine advecc2nd_corr_liberal(var, rhs)

Uses

  • proc~~advecc2nd_corr_liberal~~UsesGraph proc~advecc2nd_corr_liberal advecc2nd_corr_liberal decomp_2d decomp_2d proc~advecc2nd_corr_liberal->decomp_2d module~modfields modfields proc~advecc2nd_corr_liberal->module~modfields module~modglobal modglobal proc~advecc2nd_corr_liberal->module~modglobal module~modsubgriddata modsubgriddata proc~advecc2nd_corr_liberal->module~modsubgriddata module~modfields->decomp_2d

Arguments

Type IntentOptional Attributes Name
real, intent(in) :: var(ib-ih:ie+ih,jb-jh:je+jh,kb-kh:ke+kh)
real, intent(inout) :: rhs(ib-ih:ie+ih,jb-jh:je+jh,kb:ke+kh)

Calls

proc~~advecc2nd_corr_liberal~~CallsGraph proc~advecc2nd_corr_liberal advecc2nd_corr_liberal zstart zstart proc~advecc2nd_corr_liberal->zstart

Called by

proc~~advecc2nd_corr_liberal~~CalledByGraph proc~advecc2nd_corr_liberal advecc2nd_corr_liberal proc~ibmnorm ibmnorm proc~ibmnorm->proc~advecc2nd_corr_liberal program~dalesurban DALESURBAN program~dalesurban->proc~ibmnorm

Source Code

   subroutine advecc2nd_corr_liberal(var, rhs)
     ! Removes the advection contribution from solid scalar points as calculated
     ! by the 2nd order scheme, and replaces it with a contribution in which the
     ! value inside the solid is equal to the value outside, thereby modelling
     ! a zero (advective) flux condition.
     ! Due to potentially nonzero solid velocities due to the pressure correction,
     ! the IBM will not be conservative.
     use modglobal,      only : eps1, ib, ie, ih, jb, je, jh, kb, ke, kh, &
                                dx2i, dxi5, dy2i, dyi5, dzf, dzh2i, dzfi, dzhi, dzfi5
     use modfields,      only : u0, v0, w0
     use modsubgriddata, only : ekh
     use decomp_2d,      only : zstart

     real, intent(in)    :: var(ib-ih:ie+ih,jb-jh:je+jh,kb-kh:ke+kh)
     real, intent(inout) :: rhs(ib-ih:ie+ih,jb-jh:je+jh,kb   :ke+kh)
     integer :: i, j, k, n, m

     do n = 1,bound_info_c%nbndptsrank
      !n = bound_info_c%bndptsrank(m)
         i = bound_info_c%bndpts_loc(n,1) - zstart(1) + 1
         j = bound_info_c%bndpts_loc(n,2) - zstart(2) + 1
         k = bound_info_c%bndpts_loc(n,3) - zstart(3) + 1

         if (abs(mask_c(i+1,j,k)) < eps1) then ! var(i+1) is solid
           rhs(i,j,k) = rhs(i,j,k) + u0(i+1,j,k)*(var(i+1,j,k) + var(i,j,k))*dxi5 & ! negate contribution added in advection using var(i+1)
                                   - u0(i+1,j,k)*(var(i  ,j,k) + var(i,j,k))*dxi5   ! add corresponding contribution with var(i+1) = var(i)
         end if

         if (abs(mask_c(i-1,j,k)) < eps1) then ! var(i-1) is solid
           rhs(i,j,k) = rhs(i,j,k) - u0(i,j,k)*(var(i-1,j,k) + var(i,j,k))*dxi5 & ! negate contribution added in advection using var(i-1)
                                   + u0(i,j,k)*(var(i  ,j,k) + var(i,j,k))*dxi5   ! add corresponding contribution with var(i-1) = var(i)
         end if

         if (abs(mask_c(i,j+1,k)) < eps1) then ! var(j+1) is solid
           rhs(i,j,k) = rhs(i,j,k) + v0(i,j+1,k)*(var(i,j+1,k) + var(i,j,k))*dyi5 & ! negate contribution added in advection using var(j+1)
                                   - v0(i,j+1,k)*(var(i,j  ,k) + var(i,j,k))*dyi5   ! add corresponding contribution with var(j+1) = var(j)
         end if

         if (abs(mask_c(i,j-1,k)) < eps1) then ! var(j-1) is solid
           rhs(i,j,k) = rhs(i,j,k) - v0(i,j,k)*(var(i,j-1,k) + var(i,j,k))*dyi5 & ! negate contribution added in advection using var(j-1)
                                   + v0(i,j,k)*(var(i,j  ,k) + var(i,j,k))*dyi5   ! add corresponding contribution with var(j-1) = var(j)
         end if

         if (abs(mask_c(i,j,k+1)) < eps1) then ! var(k+1) is solid
           rhs(i,j,k) = rhs(i,j,k) + w0(i,j,k+1)*(var(i,j,k+1)*dzf(k) + var(i,j,k)*dzf(k+1))*dzhi(k+1)*dzfi5(k) & ! negate contribution added in advection using var(k+1)
                                   - w0(i,j,k+1)*(var(i,j,k  )*dzf(k) + var(i,j,k)*dzf(k+1))*dzhi(k+1)*dzfi5(k)   ! add corresponding contribution with var(k+1) = var(k)
         end if

         if (abs(mask_c(i,j,k-1)) < eps1) then ! var(k-1) is solid
           rhs(i,j,k) = rhs(i,j,k) - w0(i,j,k)*(var(i,j,k-1)*dzf(k) + var(i,j,k)*dzf(k-1))*dzhi(k)*dzfi5(k) & ! negate contribution added in advection using var(k-1)
                                   + w0(i,j,k)*(var(i,j,k  )*dzf(k) + var(i,j,k)*dzf(k-1))*dzhi(k)*dzfi5(k)   ! add corresponding contribution with var(k-1) = var(k)
         end if

     end do
   end subroutine advecc2nd_corr_liberal