m2d.F90 8.58 KB
Newer Older
1
!$Id: m2d.F90,v 1.9 2004-01-05 08:59:38 kbk Exp $
gotm's avatar
gotm committed
2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20
#include "cppdefs.h"
!-----------------------------------------------------------------------
!BOP
!
! !MODULE: m2d - depth integrated hydrodynamical model (2D)
!
! !INTERFACE:
   module m2d
!
! !DESCRIPTION:
!  This modules contains declarations for all variables related to 2D
!  hydrodynamical calculations. Information about the calculation domain
!  is included from the \emph{domain.F90} module.
!  The module contains public subroutines for initialisation, integration
!  and clean up of the 2D model component.
!  The actual calculation routines are called in integrate\_2d and is linked
!  in from the library lib2d.a.
!
! !USES:
kbk's avatar
kbk committed
21
   use time, only: julianday,secondsofday
gotm's avatar
gotm committed
22
   use parameters, only: avmmol
kbk's avatar
kbk committed
23 24
   use domain, only: imin,imax,jmin,jmax,az,au,av,H,HU,HV,min_depth
   use variables_2d
gotm's avatar
gotm committed
25 26 27
   IMPLICIT NONE
!
! !PUBLIC DATA MEMBERS:
kbk's avatar
kbk committed
28
   logical                   :: have_boundaries
29
   REALTYPE                  :: dtm, z0_const=0.010,Am=-_ONE_,An=-_ONE_
kbk's avatar
kbk committed
30 31 32 33 34 35 36 37 38
   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
39 40 41 42 43
!
! !REVISION HISTORY:
!  Original author(s): Karsten Bolding & Hans Burchard
!
!  $Log: m2d.F90,v $
44 45 46 47
!  Revision 1.9  2004-01-05 08:59:38  kbk
!  different velocity point depth calculations using vel_depth_method
!
!  Revision 1.8  2003/09/13 10:00:51  kbk
48 49 50
!  added numerical diffusion - An - to namelist - passed to uv_diffusion()
!
!  Revision 1.7  2003/08/28 10:28:40  kbk
51 52 53
!  explict setting UEx and VEx to 0 every micro time step
!
!  Revision 1.6  2003/08/15 12:47:40  kbk
54 55 56
!  also calling cfl_check() when parallel run
!
!  Revision 1.5  2003/08/03 08:53:30  kbk
57 58 59
!  not calling cfl_check when PARALLEL - should be fixed
!
!  Revision 1.4  2003/05/12 09:22:28  kbk
kbk's avatar
kbk committed
60 61 62
!  removed use halo_zones - not used
!
!  Revision 1.3  2003/04/23 12:09:43  kbk
kbk's avatar
kbk committed
63 64 65
!  cleaned code + TABS to spaces
!
!  Revision 1.2  2003/04/07 12:17:08  kbk
kbk's avatar
kbk committed
66 67 68 69
!  parallel version
!
!  Revision 1.1.1.1  2002/05/02 14:00:41  gotm
!  recovering after CVS crash
gotm's avatar
gotm committed
70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124
!
!  Revision 1.11  2001/10/22 11:55:30  bbh
!  Only call uv_diffusion() if Am > zero
!
!  Revision 1.10  2001/10/22 08:48:30  bbh
!  Am moved from paramters.F90 to m2d.F90
!
!  Revision 1.9  2001/10/11 13:07:03  bbh
!  Support for horizontal diffusion
!
!  Revision 1.8  2001/09/01 17:07:10  bbh
!  Ramping of surface elevation boundaries - via namelist
!
!  Revision 1.7  2001/08/27 11:53:13  bbh
!  TVD-advection for momentum added, some bugs removed
!
!  Revision 1.6  2001/06/22 08:19:10  bbh
!  Compiler options such as USE_MASK and OLD_DRY deleted.
!  Open and passive boundary for z created.
!  Various inconsistencies removed.
!  wait_halo added.
!  Checked loop boundaries
!
!  Revision 1.5  2001/05/18 12:55:13  bbh
!  Included masks in calls to update_2d_halo()
!
!  Revision 1.4  2001/05/06 18:51:55  bbh
!  Towards proper implementation of specified 2D bdy.
!
!  Revision 1.3  2001/05/03 19:35:01  bbh
!  Use of variables_2d
!
!  Revision 1.2  2001/04/24 08:24:58  bbh
!  Use runtype instead of macro
!
!  Revision 1.1.1.1  2001/04/17 08:43:07  bbh
!  initial import into CVS
!
! !LOCAL VARIABLES:
!
!EOP
!-----------------------------------------------------------------------

   contains

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: init_2d - initialise 2D relatedstuff.
!
! !INTERFACE:
   subroutine init_2d(runtype,timestep,hotstart)
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
kbk's avatar
kbk committed
125 126 127
   integer, intent(in)                 :: runtype
   REALTYPE, intent(in)                :: timestep
   logical, intent(in)                 :: hotstart
gotm's avatar
gotm committed
128 129 130 131 132 133 134 135 136 137 138 139 140
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
!  Allocates memiory for 2D related fields.
!
! !REVISION HISTORY:
!
!  22Apr99   Karsten Bolding & Hans Burchard  Initial code.
!
! !LOCAL VARIABLES:
kbk's avatar
kbk committed
141
   integer                   :: rc
142
   integer                   :: vel_depth_method=0
kbk's avatar
kbk committed
143
   namelist /m2d/ &
144
          MM,z0_const,Am,An,residual,bdy2d,bdyfmt_2d,bdyramp_2d,bdyfile_2d
gotm's avatar
gotm committed
145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163
!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
164
#ifdef PARALLEL
165 166
!   STDERR 'Not calling cfl_check() - PARALLEL'
   call cfl_check()
kbk's avatar
kbk committed
167
#else
168
   call cfl_check()
kbk's avatar
kbk committed
169
#endif
gotm's avatar
gotm committed
170 171

   if (Am .lt. _ZERO_) then
172 173 174 175
      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
176 177 178 179 180 181 182 183 184 185 186
   end if

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

!  Boundary related information
   if (bdy2d) then
kbk's avatar
kbk committed
187 188
!     call have_bdy()
!     call print_bdy('Local Boundary Information')
gotm's avatar
gotm committed
189 190 191
!kbk     if (have_boundaries) call init_2d_bdy(bdyfmt_2d,bdyfile_2d)
   end if

192
   call uv_depths(vel_depth_method)
gotm's avatar
gotm committed
193

kbk's avatar
kbk committed
194 195 196 197
   where ( -H+min_depth .gt. _ZERO_ )
      z = -H+min_depth
   end where
   zo=z
gotm's avatar
gotm committed
198

kbk's avatar
kbk committed
199 200 201 202
   where (-HU+min_depth .gt. _ZERO_ )
      zu = -HU+min_depth
   end where
   zub = z0_const ; zub0 = z0_const
gotm's avatar
gotm committed
203

kbk's avatar
kbk committed
204 205 206 207
   where (-HV+min_depth .gt. _ZERO_ )
      zv = -HV+min_depth
   end where
   zvb = z0_const ; zvb0 = z0_const
gotm's avatar
gotm committed
208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228

   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
229 230 231 232
   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
233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
!  A wrapper to call all 2D related subroutines in one subroutine.
!
! !REVISION HISTORY:
!  22Nov Author name Initial code
!
! !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
261
   UEx=_ZERO_ ; VEx=_ZERO_
gotm's avatar
gotm committed
262 263 264 265 266
#ifdef NO_ADVECT
   STDERR 'NO_ADVECT 2D'
#else
#ifndef UV_ADV_DIRECT
   call uv_advect()
267 268
   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
269 270 271 272 273 274 275 276
   end if
#endif
#endif
   call momentum(loop,tausx,tausy,airp)
   if (runtype .gt. 1) then
      Uint=Uint+U
      Vint=Vint+V
   end if
277 278 279
   call sealevel()
   call depth_update()

gotm's avatar
gotm committed
280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340
   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:
!  This routine cleans up after a 2D integration. Close open files etc.
!
! !REVISION HISTORY:
!  22Nov Author name Initial code
!
! !LOCAL VARIABLES:
!
!EOP
!-----------------------------------------------------------------------
!BOC
#ifdef DEBUG
   integer, save :: Ncall = 0
   Ncall = Ncall+1
   write(debug,*) 'clean_2d() # ',Ncall
#endif

   call do_residual(1)

#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)         !
!-----------------------------------------------------------------------