PUBLIC INTERFACE / ROUTINES / NAMELIST / DIAGNOSTICS / ERRORS / REFERENCES / NOTES


Module shallow_physics_mod

     Contact:   Isaac Held, Bruce Wyman
     Reviewers:
     Change history: WebCVS Log for shallow_physics.f90


OVERVIEW


    Provides simple forcing for a shallow water model.


    Forcing is generated by relaxing geopotential height to an
    equilibrium state. The velocities are also relaxed linearly
    to zero.

    Two mass sources may be specified. The ITCZ source results in a 
    zonally symmetric Hadley cell response. The monsoonal source
    results in a subtropical anticyclone. When both the ITCZ and
    monsoonal sources are present (the default), one also generates a
    quasi-stationary Rossby wave emanating from the monsoonal source.


OTHER MODULES USED


           fms_mod
  time_manager_mod


PUBLIC INTERFACE


  use shallow_physics_mod [,only: shallow_physics_init,
                             shallow_physics,
                             shallow_physics_end     ]

shallow_physics_init           : initializes module
shallow_physics                : adds tendencies to existing tendencies
shallow_physics_end            : releases memory and terminates the module

Notes:

 Optional namelist input can be read from file input.nml.
 No restart files are needed or generated by this module.
 No data files are needed.


PUBLIC ROUTINES


call shallow_physics_init ( axes, Time, lon, lat )

 INPUT:

    axes -- integer, dimension(4)
        axis indices (x,y,zfull,zhalf) returned by module diag_manager_mod

    Time -- type(time_type)
        the current time

    lon  -- longitude in radians at the center of a grid box
              [real, dimension(:) or dimension(:,:)]

    lat  -- latitude in radians at the center of a grid box
              [real, dimension(:) or dimension(:,:)]


call shallow_physics ( is, ie, js, je, timelev, dt, Time, um, vm, hm, u, v, h, u_dt, v_dt, h_dt ) INPUT: is, ie, js, je, timelev, dt, Time, u, v, h, um, vm, hm is, ie, js, je -- integer, scalar indices of the data window being passed timelev -- integer, scalar time level used to compute the forcing if timelev = -1, values at previous time level (tau-1) are used (um,vm,hm) if timelev = 0, values at current time level (tau) are used (u,v,h) if timelev = +1, values at next time level (tau+1) are used (e.g., um+dt*u_dt) dt -- real, scalar the time step (seconds) for leap-frog schemes use 2*dt Time -- type(time_type) the current time um, vm, hm -- real, dimension(:.:) zonal wind (m/s), meridional wind(m/s), geopotential height (m2/s2) at the previous time level u, v, h -- real, dimension(:.:) zonal wind (m/s), meridional wind(m/s), geopotential height (m2/s2) at the current time level IN/OUT: u_dt, v_dt, h_dt -- real, dimension(:,:) u_dt = du/dt = tendency of u (m/(s^2)), etc output has heating and friction added to input tendencies

NAMELIST


(all values listed are default values)   

logical :: no_forcing = .false.  Turns off all forcing.
                                 (Interfaces become dummy routines)

real    :: h_0             = 3.e04 (m2/s2)   mean height of flow
real    :: h_monsoon       = 2.e04 (m2/s2)   amplitude of monsoon source
real    :: lon_monsoon     =  90.0 (degrees) longitude of monsoon source origin
real    :: lat_monsoon     =  25.0 (degrees) latitude of monsoon source origin
real    :: width_monsoon   =  15.0 (degrees) width of monsoon source
real    :: h_itcz          = 1.e05 (m2/s2)   amplitude of ITCZ source
real    :: width_itcz      =  4.0  (degrees) latitude width of ITCZ source

real    :: fric_damp_time  = -20.0  frictional damping time
real    :: therm_damp_time = -10.0  thermal damping time

damping times units are (1/s) 
  The namelist value is the inverse of damping time,
  if positive the units are in (seconds), and
  if negative the units are in (days).


DIAGNOSTIC FIELDS


None.


ERROR MESSAGES


FATAL errors in shallow_physics_mod

    module has already been initialized
        You have attempted to call shallow_physics_init more than once.
        If this is what you intend to do, then call shallow_physics_end
        before recalling shallow_physics_init.

    module has not been initialized
        A FATAL error will occur in shallow_physics if you 
        have not previously called shallow_physics_init.
        Also, a WARNING will occur in shallow_physics_end if you
        have not previously called shallow_physics_init.

    invalid value for timelev argument
        The value of timelev must be set as -1, 0, or +1.
        Refer to the namelist documentation for details.


REFERENCES


None.


NOTES


None.


FUTURE PLANS


None.