diffc Subroutine

public subroutine diffc(hi, hj, hk, putin, putout)

Uses

  • proc~~diffc~~UsesGraph proc~diffc modsubgrid::diffc module~modglobal modglobal proc~diffc->module~modglobal module~modmpi modmpi proc~diffc->module~modmpi 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)

Called by

proc~~diffc~~CalledByGraph proc~diffc modsubgrid::diffc proc~subgrid modsubgrid::subgrid proc~subgrid->proc~diffc program~dalesurban DALESURBAN program~dalesurban->proc~subgrid

Contents

Source Code


Source Code

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

    use modglobal, only : ib,ie,ih,jb,je,jh,kb,ke,kh,dxh,dxhi,dxh2i,dxf,dxfi,dzf,dzfi,dyi,dy2i,&
         dzhi,dzh2i,jmax, numol, prandtlmoli,lles
    use modmpi, only : myid
    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, 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    cekh
    integer i,j,k,jm,jp,km,kp

    if (lles) then

       do k=kb,ke
          kp=k+1
          km=k-1

          do j=jb,je
             jp=j+1
             jm=j-1

             do i=ib,ie
                putout(i,j,k) = putout(i,j,k) &
                     +  0.5 *  ( &
                     ( (ekh(i+1,j,k)*dxf(i)+ekh(i,j,k)*dxf(i+1)) *(putin(i+1,j,k)-putin(i,j,k))*dxh2i(i+1) &
                     -(ekh(i,j,k)*dxf(i-1)+ekh(i-1,j,k)*dxf(i))*(putin(i,j,k)-putin(i-1,j,k))*dxh2i(i)) * dxfi(i) &
                     + &
                     ( (ekh(i,jp,k)+ekh(i,j,k)) *(putin(i,jp,k)-putin(i,j,k)) &
                     -(ekh(i,j,k)+ekh(i,jm,k)) *(putin(i,j,k)-putin(i,jm,k)) )*dy2i &
                     + &
                     ( (dzf(kp)*ekh(i,j,k) + dzf(k)*ekh(i,j,kp)) &
                     *  (putin(i,j,kp)-putin(i,j,k)) * dzh2i(kp) &
                     - &
                     (dzf(km)*ekh(i,j,k) + dzf(k)*ekh(i,j,km)) &
                     *  (putin(i,j,k)-putin(i,j,km)) * dzh2i(k)           )*dzfi(k) &
                     ) 
             end do
          end do
       end do

    else ! DNS
       cekh = numol* prandtlmoli
       do k=kb,ke
          kp=k+1
          km=k-1

          do j=jb,je
             jp=j+1
             jm=j-1

             do i=ib,ie
                putout(i,j,k) = putout(i,j,k) &
                     +( &
                     ( cekh *(putin(i+1,j,k)-putin(i,j,k))*dxhi(i+1) &
                     -cekh*(putin(i,j,k)-putin(i-1,j,k))*dxhi(i)) * dxfi(i) &
                     + &
                     ( cekh *(putin(i,jp,k)-putin(i,j,k)) &
                     -cekh *(putin(i,j,k)-putin(i,jm,k)) )*dy2i &
                     + &
                     ( cekh * (putin(i,j,kp)-putin(i,j,k)) * dzhi(kp) &
                     -cekh * (putin(i,j,k)-putin(i,j,km)) * dzhi(k)           )*dzfi(k) &
                     ) 
             end do
          end do
       end do

    end if ! lles=.true.

  end subroutine diffc