advecu_2nd Subroutine

subroutine advecu_2nd(putin, putout)

Uses

  • proc~~advecu_2nd~~UsesGraph proc~advecu_2nd advec_2nd.f90::advecu_2nd module~modfields modfields proc~advecu_2nd->module~modfields module~modglobal modglobal proc~advecu_2nd->module~modglobal module~modibm modibm proc~advecu_2nd->module~modibm module~modmpi modmpi proc~advecu_2nd->module~modmpi module~modibmdata modibmdata module~modibm->module~modibmdata mpi mpi module~modmpi->mpi

Arguments

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

Called by

proc~~advecu_2nd~~CalledByGraph proc~advecu_2nd advec_2nd.f90::advecu_2nd proc~advection advection.f90::advection proc~advection->proc~advecu_2nd program~dalesurban DALESURBAN program~dalesurban->proc~advection

Contents

Source Code


Source Code

subroutine advecu_2nd(putin, putout)

   use modglobal, only:ih, ib, ie, jb, je, jh, kb, ke, kh, dxhiq, dyiq, dzf, dzfi5, dzhi, dxhi, libm
   use modfields, only:u0, v0, w0, pres0
   use modibm, only:nxwallsnorm, nzwallsnorm, nywallsm, nywallsp, ywallsm, ywallsp, &
      xwallsnorm, zwallsnorm
   use modmpi, only:myid
   implicit none

   real, dimension(ib - ih:ie + ih, jb - jh:je + jh, kb - kh:ke + kh), intent(in)  :: putin !< Input: the u-field
   real, dimension(ib - ih:ie + ih, jb - jh:je + jh, kb:ke + kh), intent(inout) :: putout !< Output: the tendency

   integer :: i, j, k, ip, im, jp, jm, kp, km, il, iu, jl, ju, kl, ku, n

   do k = kb, ke
      km = k - 1
      kp = k + 1
      do j = jb, je
         jm = j - 1
         jp = j + 1
         do i = ib, ie
            im = i - 1
            ip = i + 1
            putout(i, j, k) = putout(i, j, k) - ( &
                              ( &
                              (putin(i, j, k) + putin(ip, j, k))*(u0(i, j, k) + u0(ip, j, k)) &
                              - (putin(i, j, k) + putin(im, j, k))*(u0(i, j, k) + u0(im, j, k)) & ! d(uu)/dx
                              )*dxhiq(i) &
                              + ( &
                              (putin(i, j, k) + putin(i, jp, k))*(v0(i, jp, k) + v0(im, jp, k)) &
                              - (putin(i, j, k) + putin(i, jm, k))*(v0(i, j, k) + v0(im, j, k)) & ! d(vu)/dy
                              )*dyiq) &
                              - ((pres0(i, j, k) - pres0(i - 1, j, k))*dxhi(i)) ! - dp/dx

         end do
      end do
   end do

   do j = jb, je
      jm = j - 1
      jp = j + 1
      do i = ib, ie
         im = i - 1
         ip = i + 1
         do k = kb, ke
            km = k - 1
            kp = k + 1
            putout(i, j, k) = putout(i, j, k) - ( &
                              (putin(i, j, kp)*dzf(k) + putin(i, j, k)*dzf(kp))*dzhi(kp) &
                              *(w0(i, j, kp) + w0(im, j, kp)) &
                              - (putin(i, j, k)*dzf(km) + putin(i, j, km)*dzf(k))*dzhi(k) &
                              *(w0(i, j, k) + w0(im, j, k)) &
                              )*0.5*dzfi5(k)
         end do
      end do
   end do

end subroutine advecu_2nd