advecc_upw Subroutine

subroutine advecc_upw(hi, hj, hk, putin, putout)

Uses

  • proc~~advecc_upw~~UsesGraph proc~advecc_upw advec_upw.f90::advecc_upw module~modfields modfields proc~advecc_upw->module~modfields module~modglobal modglobal proc~advecc_upw->module~modglobal

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: hi
integer, intent(in) :: hj
integer, intent(in) :: hk
real, intent(in), dimension(ib - hi:ie + hi, jb - hj:je + hj, kb - hk:ke + hk) :: putin
real, intent(inout), dimension(ib - hi:ie + hi, jb - hj:je + hj, kb:ke + hk) :: putout

Called by

proc~~advecc_upw~~CalledByGraph proc~advecc_upw advec_upw.f90::advecc_upw proc~advection advection.f90::advection proc~advection->proc~advecc_upw program~dalesurban DALESURBAN program~dalesurban->proc~advection

Contents

Source Code


Source Code

subroutine advecc_upw(hi, hj, hk, putin, putout)

   use modglobal, only:ib, ie, ih, jb, je, jh, kb, ke, dyi, dxfci, dzfci
   use modfields, only:u0, v0, w0
   implicit none

   integer, intent(in) :: hi !< size of halo in i
   integer, intent(in) :: hj !< size of halo in j
   integer, intent(in) :: hk !< size of halo in k
   real, dimension(ib - hi:ie + hi, jb - hj:je + hj, kb - hk:ke + hk), intent(in)  :: putin !< Input: the cell centered field
   real, dimension(ib - hi:ie + hi, jb - hj:je + hj, kb:ke + hk), intent(inout) :: putout !< Output: the tendency

   real, allocatable, dimension(:, :, :) :: put
   integer :: i, j, k

   allocate (put(ib - hi:ie + hi, jb - hj:je + hj, kb - hk:ke + hk))

   do k = kb, ke
      do j = jb, je
         do i = ib, ie + 1
            if (u0(i, j, k) > 0) then
               put(i, j, k) = putin(i - 1, j, k)
            else
               put(i, j, k) = putin(i, j, k)
            endif
         enddo
      enddo
   enddo

   do k = kb, ke
      do j = jb, je
         do i = ib, ie
            putout(i, j, k) = putout(i, j, k) - &
                              (u0(i + 1, j, k)*put(i + 1, j, k) - u0(i, j, k)*put(i, j, k))*dxfci(i)
         enddo
      enddo
   enddo

   do k = kb, ke
      do j = jb, je + 1
         do i = ib, ie
            if (v0(i, j, k) > 0) then
               put(i, j, k) = putin(i, j - 1, k)
            else
               put(i, j, k) = putin(i, j, k)
            endif
         enddo
      enddo
   enddo
   do k = kb, ke
      do j = jb, je
         do i = ib, ie
            putout(i, j, k) = putout(i, j, k) - &
                              (v0(i, j + 1, k)*put(i, j + 1, k) - v0(i, j, k)*put(i, j, k))*dyi
         enddo
      enddo
   enddo

   do k = kb, ke + 1
      do j = jb, je
         do i = ib, ie
            if (w0(i, j, k) > 0) then
               put(i, j, k) = putin(i, j, k - 1)
            else
               put(i, j, k) = putin(i, j, k)
            endif
         enddo
      enddo
   enddo
   do k = kb, ke
      do j = jb, je
         do i = ib, ie
            putout(i, j, k) = putout(i, j, k) - &
                              (w0(i, j, k + 1)*put(i, j, k + 1) - w0(i, j, k)*put(i, j, k))*dzfci(k)
         enddo
      enddo
   enddo

   deallocate (put)

end subroutine advecc_upw