subgridnamelist Subroutine

public subroutine subgridnamelist()

Uses

  • proc~~subgridnamelist~~UsesGraph proc~subgridnamelist subgridnamelist module~modglobal modglobal proc~subgridnamelist->module~modglobal module~modmpi modmpi proc~subgridnamelist->module~modmpi mpi mpi module~modmpi->mpi

Arguments

None

Calls

proc~~subgridnamelist~~CallsGraph proc~subgridnamelist subgridnamelist mpi_bcast mpi_bcast proc~subgridnamelist->mpi_bcast

Called by

proc~~subgridnamelist~~CalledByGraph proc~subgridnamelist subgridnamelist proc~initsubgrid initsubgrid proc~initsubgrid->proc~subgridnamelist program~dalesurban DALESURBAN program~dalesurban->proc~initsubgrid

Namelists

Namelist NAMSUBGRID


Variables

Name Type Default Description
ldelta logical .false.
lmason logical .false.
cf real 2.5
cn real 0.76
Rigc real 0.25
Prandtl real 0.333
lsmagorinsky logical .false.
lvreman logical .false.
loneeqn logical .false.
c_vreman real 0.07
cs real -1.
nmason real 2.
lbuoycorr logical .false.

Source Code

  subroutine subgridnamelist
    use modglobal, only : pi,ifnamopt,fname_options,lles
    use modmpi,    only : myid, nprocs, comm3d, mpierr, my_real, mpi_logical, mpi_integer

    implicit none

    integer :: ierr

    namelist/NAMSUBGRID/ &
         ldelta,lmason, cf,cn,Rigc,Prandtl,lsmagorinsky,lvreman,loneeqn,c_vreman,cs,nmason,lbuoycorr

    if(myid==0)then
       open(ifnamopt,file=fname_options,status='old',iostat=ierr)
       read (ifnamopt,NAMSUBGRID,iostat=ierr)
       if (ierr > 0) then
          write(0, *) 'ERROR: Problem in namoptions NAMSUBGRID'
          write(0, *) 'iostat error: ', ierr
          stop 1
       endif
       !write(6 ,NAMSUBGRID)
       close(ifnamopt)
    end if

    call MPI_BCAST(ldelta     ,1,MPI_LOGICAL,0,comm3d,mpierr)
    call MPI_BCAST(lmason     ,1,MPI_LOGICAL,0,comm3d,mpierr)
    call MPI_BCAST(nmason     ,1,MY_REAL    ,0,comm3d,mpierr)
    call MPI_BCAST(lsmagorinsky,1,MPI_LOGICAL,0,comm3d,mpierr)
    call MPI_BCAST(lvreman    ,1,MPI_LOGICAL,0,comm3d,mpierr)
    call MPI_BCAST(lbuoycorr  ,1,MPI_LOGICAL,0,comm3d,mpierr)
    call MPI_BCAST(loneeqn    ,1,MPI_LOGICAL,0,comm3d,mpierr)
    call MPI_BCAST(c_vreman   ,1,MY_REAL   ,0,comm3d,mpierr)
    call MPI_BCAST(cs         ,1,MY_REAL   ,0,comm3d,mpierr)
    call MPI_BCAST(cf         ,1,MY_REAL   ,0,comm3d,mpierr)
    call MPI_BCAST(cn         ,1,MY_REAL   ,0,comm3d,mpierr)
    call MPI_BCAST(Rigc       ,1,MY_REAL   ,0,comm3d,mpierr)
    call MPI_BCAST(Prandtl    ,1,MY_REAL   ,0,comm3d,mpierr)
    prandtli = 1./Prandtl

    if ((lsmagorinsky) .or. (lvreman) .or. (loneeqn)) then
       lles =.true.
    endif


  end subroutine subgridnamelist