zwallscalar Subroutine

public subroutine zwallscalar(hi, hj, hk, putin, putout, bcvalue, n)

Uses

  • proc~~zwallscalar~~UsesGraph proc~zwallscalar modibm::zwallscalar module~initfac initfac proc~zwallscalar->module~initfac module~modglobal modglobal proc~zwallscalar->module~modglobal module~modmpi modmpi proc~zwallscalar->module~modmpi module~modsubgriddata modsubgriddata proc~zwallscalar->module~modsubgriddata module~initfac->module~modglobal module~initfac->module~modmpi netcdf netcdf module~initfac->netcdf mpi mpi module~modmpi->mpi

Arguments

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

Called by

proc~~zwallscalar~~CalledByGraph proc~zwallscalar modibm::zwallscalar proc~zwallfun modibm::zwallfun proc~zwallfun->proc~zwallscalar proc~ibmwallfun modibm::ibmwallfun proc~ibmwallfun->proc~zwallfun program~dalesurban DALESURBAN program~dalesurban->proc~ibmwallfun

Contents

Source Code


Source Code

   subroutine zwallscalar(hi, hj, hk, putin, putout, bcvalue, n)
      use modglobal, only:jmax, dzf, dzfi, dzhi, dzh2i, ib, ie, jb, je, kb, ke, prandtlmoli, numol
      use modsubgriddata, only:ekh
      use modmpi, only:myid
      use initfac, only:block
      integer i, j, k, il, iu, jl, ju, km
      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, intent(in)    :: putin(ib - hi:ie + hi, jb - hj:je + hj, kb - hk:ke + hk)
      real, intent(inout) :: putout(ib - hi:ie + hi, jb - hj:je + hj, kb:ke + hk)
      real, intent(in)    :: bcvalue
      integer, intent(in)    :: n

      k = block(n, 6) + 1 !block location
      km = k - 1 !
      il = block(n, 1)
      iu = block(n, 2)
      jl = MAX(block(n, 3) - myid*jmax, 1)
      ju = MIN(block(n, 4) - myid*jmax, jmax)

      !  delta=putout(3,1,2)
      do j = jl, ju
         do i = il, iu
            putout(i, j, k) = putout(i, j, k) + ( &
                              0.5*(dzf(km)*ekh(i, j, k) + dzf(k)*ekh(i, j, km))* & ! zero flux
                              (putin(i, j, k) - putin(i, j, km))*dzh2i(k) - &
                              bcvalue)*dzfi(k)
         end do
      end do
   end subroutine zwallscalar