initsubgrid Subroutine

public subroutine initsubgrid()

Uses

  • proc~~initsubgrid~~UsesGraph proc~initsubgrid modsubgrid::initsubgrid module~modglobal modglobal proc~initsubgrid->module~modglobal module~modmpi modmpi proc~initsubgrid->module~modmpi mpi mpi module~modmpi->mpi

Arguments

None

Calls

proc~~initsubgrid~~CallsGraph proc~initsubgrid modsubgrid::initsubgrid proc~subgridnamelist modsubgrid::subgridnamelist proc~initsubgrid->proc~subgridnamelist mpi_bcast mpi_bcast proc~subgridnamelist->mpi_bcast

Called by

proc~~initsubgrid~~CalledByGraph proc~initsubgrid modsubgrid::initsubgrid proc~startup modstartup::startup proc~startup->proc~initsubgrid program~dalesurban DALESURBAN program~dalesurban->proc~startup

Contents

Source Code


Source Code

  subroutine initsubgrid
    use modglobal, only : ih,ib,ie,jh,jb,je,kb,ke,kh,delta,zf,fkar, &
         pi,ifnamopt,fname_options
    use modmpi, only : myid

    implicit none

    integer   :: i, k

    real :: ceps, ch
    real :: mlen

    call subgridnamelist

    allocate(ekm(ib-ih:ie+ih,jb-jh:je+jh,kb-kh:ke+kh))
    allocate(ekh(ib-ih:ie+ih,jb-jh:je+jh,kb-kh:ke+kh))
    allocate(zlt(ib-ih:ie+ih,jb-jh:je+jh,kb:ke+kh))
    allocate(sbdiss(ib-ih:ie+ih,jb-jh:je+jh,kb:ke+kh))
    allocate(sbshr(ib-ih:ie+ih,jb-jh:je+jh,kb:ke+kh))
    allocate(sbbuo(ib-ih:ie+ih,jb-jh:je+jh,kb:ke+kh))
    allocate(csz(ib-ih:ie+ih,kb:ke+kh))
    allocate(damp(ib:ie,jb:je,kb:ke))


    damp = 1. 
    cm = cf / (2. * pi) * (1.5*alpha_kolm)**(-1.5)

    !     ch   = 2. * alpha_kolm / beta_kolm
    ch   = prandtl
    ch2  = ch-ch1

    ceps = 2. * pi / cf * (1.5*alpha_kolm)**(-1.5)
    ce1  = (cn**2)* (cm/Rigc - ch1*cm )
    ce2  = ceps - ce1

    if(cs == -1.) then
       csz(:,:)  = (cm**3/ceps)**0.25   !< Smagorinsky constant
    else
       csz(:,:)  = cs
    end if

    !    if(lmason) then
    !      do k = kb,ke+kh
    !        do i=ib-ih,ie+ih
    !          mlen   = (1. / (csz(i,k) * delta(i,k))**nmason + 1. / (fkar * zf(k))**nmason)**(-1./nmason)
    !          csz(i,k) = mlen / delta(i,k)
    !        end do
    !      end do
    !    end if

    if (myid==0) then
       write (6,*) 'cf    = ',cf
       write (6,*) 'cm    = ',cm
       write (6,*) 'ch    = ',ch
       write (6,*) 'ch1   = ',ch1
       write (6,*) 'ch2   = ',ch2
       write (6,*) 'ceps  = ',ceps
       write (6,*) 'ceps1 = ',ce1
       write (6,*) 'ceps2 = ',ce2
       write (6,*) 'cs    = ',cs
       write (6,*) 'Rigc  = ',Rigc
    endif

  end subroutine initsubgrid