cfl_check - check for explicit barotropic time step. (Source File: cfl_check.F90)

INTERFACE:

    subroutine cfl_check()
DESCRIPTION:

This routine loops over all horizontal grid points and calculated the maximum time step according to the CFL criterium by Beckers and Deleersnijder (1993):

$\displaystyle \Delta t_{\max} = \min_{i,j} \left\{\frac{\Delta x_{i,j} \Delta y_{i,j}}
{\sqrt{2} c_{i,j} \sqrt{\Delta x_{i,j}^2+ \Delta y_{i,j}^2}}\right\}$ (95)

with the local shallow water wave speed

$\displaystyle c_{i,j}=\sqrt{g H_{i,j}},$ (96)

where $g$ is the gravitational acceleration and $H_{i,j}$ is the local bathymetry value. In case that the chosen micro time step $\Delta t_m$ is larger than $\Delta t_{\max}$, the program will be aborted. In any case the CFL diagnostics will be written to standard output.

USES:

    use parameters, only: g
    use domain, only: imin,imax,jmin,jmax,H,az
 #if defined(SPHERICAL) || defined(CURVILINEAR)
    use domain, only: dyc,dxc
 #else
    use domain, only: dy,dx
 #endif
    use m2d, only: dtm
    IMPLICIT NONE
REVISION HISTORY:
    Original author(s): Karsten Bolding & Hans Burchard
LOCAL VARIABLES:
    integer                   :: pos(2),max_pos(2),rc,i,j
    REALTYPE                  :: h_max=-99.,c,max_dt,dtt
    logical, dimension(:,:), allocatable :: lmask