m2d.F90 7.5 KB
Newer Older
kbk's avatar
kbk committed
1
!$Id: m2d.F90,v 1.15 2006-03-01 15:54:07 kbk Exp $
gotm's avatar
gotm committed
2 3 4 5 6 7 8 9 10 11
#include "cppdefs.h"
!-----------------------------------------------------------------------
!BOP
!
! !MODULE: m2d - depth integrated hydrodynamical model (2D)
!
! !INTERFACE:
   module m2d
!
! !DESCRIPTION:
hb's avatar
hb committed
12
!  This module contains declarations for all variables related to 2D
gotm's avatar
gotm committed
13
!  hydrodynamical calculations. Information about the calculation domain
hb's avatar
hb committed
14
!  is included from the {\tt domain} module.
gotm's avatar
gotm committed
15 16
!  The module contains public subroutines for initialisation, integration
!  and clean up of the 2D model component.
hb's avatar
hb committed
17 18 19
!  The actual calculation routines are called in {\tt integrate\_2d}
!  and are linked
!  in from the library {\tt lib2d.a}.
gotm's avatar
gotm committed
20 21
!
! !USES:
kbk's avatar
kbk committed
22
   use time, only: julianday,secondsofday
gotm's avatar
gotm committed
23
   use parameters, only: avmmol
kbk's avatar
kbk committed
24 25
   use domain, only: imin,imax,jmin,jmax,az,au,av,H,HU,HV,min_depth
   use variables_2d
gotm's avatar
gotm committed
26 27 28
   IMPLICIT NONE
!
! !PUBLIC DATA MEMBERS:
kbk's avatar
kbk committed
29
   logical                   :: have_boundaries
30
   REALTYPE                  :: dtm, z0_const=0.010,Am=-_ONE_,An=-_ONE_
kbk's avatar
kbk committed
31 32 33 34 35 36 37 38 39
   integer                   :: MM=1,residual=-1
   logical                   :: bdy2d=.false.
   integer                   :: bdyfmt_2d,bdytype,bdyramp_2d=-1
   character(len=PATH_MAX)   :: bdyfile_2d
   REAL_4B                   :: bdy_old(1000)
   REAL_4B                   :: bdy_new(1000)
   REAL_4B                   :: bdy_data(1000)
   REAL_4B, allocatable      :: bdy_times(:)
   integer, parameter        :: comm_method=-1
gotm's avatar
gotm committed
40
!
kbk's avatar
kbk committed
41 42 43
! !REVISION HISTORY:
!  Original author(s): Karsten Bolding & Hans Burchard
!
gotm's avatar
gotm committed
44 45 46 47 48 49 50 51 52 53
! !LOCAL VARIABLES:
!
!EOP
!-----------------------------------------------------------------------

   contains

!-----------------------------------------------------------------------
!BOP
!
hb's avatar
hb committed
54
! !IROUTINE: init_2d - initialise 2D related stuff.
gotm's avatar
gotm committed
55 56 57 58 59 60
!
! !INTERFACE:
   subroutine init_2d(runtype,timestep,hotstart)
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
kbk's avatar
kbk committed
61 62 63
   integer, intent(in)                 :: runtype
   REALTYPE, intent(in)                :: timestep
   logical, intent(in)                 :: hotstart
gotm's avatar
gotm committed
64 65 66 67 68 69
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
hb's avatar
hb committed
70 71 72 73 74 75
!  Here, the {\tt m2d} namelist is read from {\tt getm.inp}, and the check
!  for the fulfilment of the CFL criterium for shallow water theory
!  {\tt cfl\_check} is called. A major part of this subroutine deals
!  then with the setting of local bathymetry values and initial surface
!  elevations in $u$- and $v$-points, also by calls to the subroutines 
!  {\tt uv\_depths} and {\tt depth\_update}.
gotm's avatar
gotm committed
76 77
!
! !LOCAL VARIABLES:
kbk's avatar
kbk committed
78
   integer                   :: rc
79
   integer                   :: vel_depth_method=0
kbk's avatar
kbk committed
80
   namelist /m2d/ &
kbk's avatar
kbk committed
81
          MM,z0_const,vel_depth_method,Am,An,residual, &
kbk's avatar
kbk committed
82
          bdy2d,bdyfmt_2d,bdyramp_2d,bdyfile_2d
gotm's avatar
gotm committed
83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101
!EOP
!-------------------------------------------------------------------------
!BOC
#ifdef DEBUG
   integer, save :: Ncall = 0
   Ncall = Ncall+1
   write(debug,*) 'init_2d() # ',Ncall
#endif

   LEVEL1 'init_2d'

!  Read 2D-model specific things from the namelist.
   read(NAMLST,m2d)

   dtm = timestep

!  Allocates memory for the public data members - if not static
   call init_variables_2d(runtype)

kbk's avatar
kbk committed
102
#ifdef PARALLEL
103 104
!   STDERR 'Not calling cfl_check() - PARALLEL'
   call cfl_check()
kbk's avatar
kbk committed
105
#else
106
   call cfl_check()
kbk's avatar
kbk committed
107
#endif
gotm's avatar
gotm committed
108 109

   if (Am .lt. _ZERO_) then
110 111 112 113
      LEVEL2 'Am < 0 --> horizontal momentum diffusion not included'
   end if
   if (An .lt. _ZERO_) then
      LEVEL2 'An < 0 --> numerical momentum diffusion not included'
gotm's avatar
gotm committed
114 115 116 117 118 119 120 121 122
   end if

   LEVEL2 'Open boundary=',bdy2d
   if (bdy2d) then
      if(hotstart) bdyramp_2d = -1
      LEVEL2 TRIM(bdyfile_2d)
      LEVEL2 'Format=',bdyfmt_2d
   end if

123
   call uv_depths(vel_depth_method)
gotm's avatar
gotm committed
124

kbk's avatar
kbk committed
125 126 127 128
   where ( -H+min_depth .gt. _ZERO_ )
      z = -H+min_depth
   end where
   zo=z
gotm's avatar
gotm committed
129

kbk's avatar
kbk committed
130 131 132 133
   where (-HU+min_depth .gt. _ZERO_ )
      zu = -HU+min_depth
   end where
   zub = z0_const ; zub0 = z0_const
gotm's avatar
gotm committed
134

kbk's avatar
kbk committed
135 136 137 138
   where (-HV+min_depth .gt. _ZERO_ )
      zv = -HV+min_depth
   end where
   zvb = z0_const ; zvb0 = z0_const
gotm's avatar
gotm committed
139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159

   call depth_update()

#ifdef DEBUG
   write(debug,*) 'Leaving init_2d()'
   write(debug,*)
#endif
   return
   end subroutine init_2d
!EOC

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: integrate_2d - sequence of calls to do 2D model integration
!
! !INTERFACE:
   subroutine integrate_2d(runtype,loop,tausx,tausy,airp)
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
kbk's avatar
kbk committed
160 161 162 163
   integer, intent(in)                 :: runtype,loop
   REALTYPE, intent(in)                :: tausx(E2DFIELD)
   REALTYPE, intent(in)                :: tausy(E2DFIELD)
   REALTYPE, intent(in)                :: airp(E2DFIELD)
gotm's avatar
gotm committed
164 165 166 167 168 169
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
hb's avatar
hb committed
170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189
!  Here, all 2D related subroutines are called. The major calls and their
!  meaning are:
!
!  \vspace{0.5cm}
!
!  \begin{tabular}{ll}
!  {\tt call update\_2d\_bdy} & read in new lateral boundary conditions \\
!  {\tt call bottom\_friction} & update bottom friction\\
!  {\tt call uv\_advect} & calculate 2D advection terms\\
!  {\tt call uv\_diffusion} & calculate 2D  diffusion terms\\
!  {\tt call momentum} & iterate 2D momemtum equations\\
!  {\tt call sealevel} & update sea surface elevation\\
!  {\tt call depth\_update} & update water depths\\
!  {\tt call do\_residual} & calculate intermdediate values for residual currents
!  \end{tabular}
!
!  \vspace{0.5cm}
!
!  It should be noted that some of these calls may be excluded for certain
!  compiler options set in the {\tt Makefile} of the application.
gotm's avatar
gotm committed
190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207
!
! !LOCAL VARIABLES:
!EOP
!-------------------------------------------------------------------------
!BOC
#ifdef DEBUG
   integer, save :: Ncall = 0
   Ncall = Ncall+1
   write(debug,*) 'integrate_2d() # ',Ncall
#endif

   if (have_boundaries) call update_2d_bdy(loop,bdyramp_2d)

   if (mod(loop-1,MM) .eq. 0) then        ! MacroMicro time step
#ifndef NO_BOTTFRIC
      call bottom_friction(runtype)
#endif
   end if
208
   UEx=_ZERO_ ; VEx=_ZERO_
gotm's avatar
gotm committed
209 210 211 212 213
#ifdef NO_ADVECT
   STDERR 'NO_ADVECT 2D'
#else
#ifndef UV_ADV_DIRECT
   call uv_advect()
214 215
   if (Am .gt. _ZERO_ .or. An .gt. _ZERO_) then
      call uv_diffusion(Am,An) ! Has to be called after uv_advect.
gotm's avatar
gotm committed
216 217 218 219 220 221 222 223
   end if
#endif
#endif
   call momentum(loop,tausx,tausy,airp)
   if (runtype .gt. 1) then
      Uint=Uint+U
      Vint=Vint+V
   end if
224 225 226
   call sealevel()
   call depth_update()

gotm's avatar
gotm committed
227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254
   if(residual .gt. 0 .and. loop .ge. residual) then
      call do_residual(0)
   end if

#ifdef DEBUG
     write(debug,*) 'Leaving integrate_2d()'
     write(debug,*)
#endif
   return
   end subroutine integrate_2d
!EOC

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: clean_2d - cleanup after 2D run.
!
! !INTERFACE:
   subroutine clean_2d()
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
hb's avatar
hb committed
255 256
!  This routine executes a final call to {\tt do\_residual} where the residual
!  current calculations are finished.
gotm's avatar
gotm committed
257 258 259 260 261 262 263 264 265 266 267 268
!
! !LOCAL VARIABLES:
!
!EOP
!-----------------------------------------------------------------------
!BOC
#ifdef DEBUG
   integer, save :: Ncall = 0
   Ncall = Ncall+1
   write(debug,*) 'clean_2d() # ',Ncall
#endif

kbk's avatar
kbk committed
269 270 271
   if(residual .gt. 0) then
      call do_residual(1)
   end if
gotm's avatar
gotm committed
272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287

#ifdef DEBUG
   write(debug,*) 'Leaving clean_2d()'
   write(debug,*)
#endif
   return
   end subroutine clean_2d
!EOC

!-----------------------------------------------------------------------

   end module m2d

!-----------------------------------------------------------------------
! Copyright (C) 2001 - Hans Burchard and Karsten Bolding (BBH)         !
!-----------------------------------------------------------------------