ywallfunmin Subroutine

public subroutine ywallfunmin()

Uses

  • proc~~ywallfunmin~~UsesGraph proc~ywallfunmin modibm::ywallfunmin module~initfac initfac proc~ywallfunmin->module~initfac module~modfields modfields proc~ywallfunmin->module~modfields module~modglobal modglobal proc~ywallfunmin->module~modglobal module~initfac->module~modglobal module~modmpi modmpi module~initfac->module~modmpi netcdf netcdf module~initfac->netcdf mpi mpi module~modmpi->mpi

Arguments

None

Calls

proc~~ywallfunmin~~CallsGraph proc~ywallfunmin modibm::ywallfunmin proc~wfgr wf_gr.f90::wfGR proc~ywallfunmin->proc~wfgr proc~wfmneutral wfmneutral.f90::wfmneutral proc~ywallfunmin->proc~wfmneutral proc~wfuno wf_uno.f90::wfuno proc~ywallfunmin->proc~wfuno proc~ywallscalarmin modibm::ywallscalarmin proc~ywallfunmin->proc~ywallscalarmin proc~unoh wf_uno.f90::unoh proc~wfuno->proc~unoh

Called by

proc~~ywallfunmin~~CalledByGraph proc~ywallfunmin modibm::ywallfunmin proc~ibmwallfun modibm::ibmwallfun proc~ibmwallfun->proc~ywallfunmin program~dalesurban DALESURBAN program~dalesurban->proc~ibmwallfun

Contents

Source Code


Source Code

   subroutine ywallfunmin
      use modglobal, only:dxf, dxhi, dy, dyi, dzhiq, dzf, dzhi, lles, nsv, numol, ltempeq, lmoist, &
         ih, jh, kh, ihc, jhc, khc, iwallmom, iwalltemp, iwallmoist, iwallscal, nblocks
      use modfields, only:u0, w0, up, wp, shear, thl0, thlp, qt0, qtp, sv0, svp, tfluxb, momfluxb, exnf, cth, qfluxb
      use initfac, only:fachf, block, faclGR, facqsat, facef, fachurel, facf, facT, facz0, facz0h
      !      use modsurfdata,     only : wtsurf
      integer i, j, k, n, nc, il, iu, kl, ku, im, jp, km, m

      if (iwallmom == 1) then
      !fixed flux, not implemented
      else if (iwallmom == 2) then
      do n = 1, nyminwall
      k = block(iyminwall(n, 1), 11)
      call wfuno(ih, jh, kh, up, wp, thlp, momfluxb, tfluxb, cth, bcTfluxA, u0, w0, thl0, facT(k,1), facz0(k), facz0h(k), iyminwall(n, 1), iyminwall(n, 2), 41)
      end do
      else if (iwallmom == 3) then
      do n = 1, nyminwall
      k = block(iyminwall(n, 1), 11)
      call wfmneutral(ih, jh, kh, up, wp, momfluxb, u0, w0, facz0(k), iyminwall(n, 1), iyminwall(n, 2), 41)
      end do
      end if !

      if (ltempeq) then
      if (iwalltemp == 1) then
         do n = 1, nyminwall !
            call ywallscalarmin(ih, jh, kh, thl0, thlp, bctfym, n)
         end do
      else if (iwalltemp == 2) then
         do n = 1, nyminwall
            k = block(iyminwall(n, 1), 11)
            call wfuno(ih, jh, kh, up, wp, thlp, momfluxb, tfluxb, cth, bcTfluxA, u0, w0, thl0, facT(k,1), facz0(k), facz0h(k), iyminwall(n, 1), iyminwall(n, 2), 42)
            fachf(k) = fachf(k) + bcTfluxA
         end do
      end if
      end if

      if (lmoist) then
      if (iwallmoist == 1) then
         do n = 1, nyminwall !
            call ywallscalarmin(ih, jh, kh, qt0, qtp, bcqfym, n)
         end do
      end if
      if ((ltempeq) .and. (iwallmoist == 2)) then
         do n = 1, nyminwall
            k = block(iyminwall(n, 1), 11)
            if (faclGR(k)) then
            call wfGR(ih, jh, kh, qtp, qfluxb, cth, bcqfluxA, qt0, facqsat(k),fachurel(k),facf(k,4),facf(k,5),iyminwall(n, 1), iyminwall(n, 2), 42)
               facef(k) = facef(k) + bcqfluxA
            end if
         end do
      end if
      end if

      if (nsv>0) then
      if (iwallscal == 1) then
         do n = 1, nyminwall !
            do m = 1, nsv
               call ywallscalarmin(ihc, jhc, khc, sv0(:,:,:,m), svp(:,:,:,m), 0., n)
            end do
         end do
      end if
      end if

   end subroutine ywallfunmin