m3d.F90 10.9 KB
Newer Older
kbk's avatar
kbk committed
1
!$Id: m3d.F90,v 1.5 2003-04-23 12:16:34 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 21
#include "cppdefs.h"
!-----------------------------------------------------------------------
!BOP
!
! !MODULE: m3d - 3D model component
!
! !INTERFACE:
   module m3d
!
! !DESCRIPTION:
!  This modules contains declarations for all variables related to 3D
!  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 3D model component.
!  The actual calculation routines are called in integrate\_3d and is linked
!  in from the library lib3d.a.
!
! !USES:
   use parameters, only: avmmol
kbk's avatar
kbk committed
22 23
   use domain, only: vert_cord
   use m2d, only: Am
kbk's avatar
kbk committed
24 25
   use variables_2d, only: D,z
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
26 27
   use temperature,only: init_temperature, do_temperature
   use salinity,   only: init_salinity, do_salinity
kbk's avatar
kbk committed
28 29 30
   use eqstate,    only: init_eqstate, do_eqstate
#endif
#ifndef NO_BAROCLINIC
kbk's avatar
kbk committed
31
   use suspended_matter, only: init_spm, do_spm
kbk's avatar
kbk committed
32
#endif
gotm's avatar
gotm committed
33
   use advection_3d, only: init_advection_3d
kbk's avatar
kbk committed
34
   use bdy_3d, only: init_bdy_3d, do_bdy_3d
gotm's avatar
gotm committed
35
   use variables_3d
kbk's avatar
kbk committed
36 37 38 39
#ifdef PARALLEL
use halo_mpi
#endif

gotm's avatar
gotm committed
40 41 42
   IMPLICIT NONE
!
! !PUBLIC DATA MEMBERS:
kbk's avatar
kbk committed
43 44 45 46 47 48 49 50
   integer                             :: M=1
   REALTYPE                            :: cord_relax=_ZERO_
   logical                             :: calc_temp=.true.
   logical                             :: calc_salt=.true.
   logical                             :: calc_spm=.false.
   logical                             :: bdy3d=.false.
   integer                             :: bdyfmt_3d,bdyramp_3d
   character(len=PATH_MAX)             :: bdyfile_3d
gotm's avatar
gotm committed
51 52 53 54 55
!
! !REVISION HISTORY:
!  Original author(s): Karsten Bolding & Hans Burchard
!
!  $Log: m3d.F90,v $
kbk's avatar
kbk committed
56 57 58 59
!  Revision 1.5  2003-04-23 12:16:34  kbk
!  cleaned code + TABS to spaces
!
!  Revision 1.4  2003/04/07 16:28:34  kbk
kbk's avatar
kbk committed
60
!  parallel support
gotm's avatar
gotm committed
61 62 63
!
!  Revision 1.1.1.1  2002/05/02 14:00:51  gotm
!  recovering after CVS crash
gotm's avatar
gotm committed
64 65 66 67 68 69 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 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159
!
!  Revision 1.30  2001/10/26 09:13:24  bbh
!  Only call slow_diffusion() if Am > 0.
!
!  Revision 1.29  2001/10/23 12:38:58  bbh
!  Forgot cord_relax in one call to coordinates()
!
!  Revision 1.28  2001/10/23 07:53:10  bbh
!  spm -> suspended_matter
!
!  Revision 1.27  2001/10/23 07:05:11  bbh
!  Parabolic viscosity - -DPARABOLIC_VISCOSITY
!
!  Revision 1.26  2001/10/22 11:16:09  bbh
!  Added support for particulate suspended matter - no IO yet
!
!  Revision 1.25  2001/10/22 09:26:41  bbh
!  Added cord_relax
!
!  Revision 1.24  2001/10/22 08:48:30  bbh
!  Am moved from paramters.F90 to m2d.F90
!
!  Revision 1.23  2001/10/12 11:39:20  bbh
!  TVD moved out of ??_momentum_3d.F90 and into uv_advect_3d.F90
!
!  Revision 1.22  2001/10/12 08:49:56  bbh
!  Support for horizontal diffusion
!
!  Revision 1.21  2001/09/19 13:07:00  bbh
!  Moved advection related 3D fields to global allocation
!
!  Revision 1.20  2001/09/03 20:14:03  bbh
!  Fixed reading of vel_{hor,vel,strang
!
!  Revision 1.19  2001/09/03 20:04:21  bbh
!  Allow individual advection settings for momentum, salinity and temperature
!
!  Revision 1.18  2001/09/03 13:02:52  bbh
!  Small changes due to NOMADS
!
!  Revision 1.17  2001/08/31 15:39:49  bbh
!  ifdef for CONSTANCE added
!
!  Revision 1.16  2001/08/30 08:56:26  bbh
!  Preparing for 3D boundary conditions
!
!  Revision 1.15  2001/08/29 12:08:05  bbh
!  temp_adv and salt_adv needs to be public
!
!  Revision 1.14  2001/08/29 11:21:46  bbh
!  namelists read in salinity and temperature + initialisation
!
!  Revision 1.13  2001/08/27 11:51:45  bbh
!  TVD-advection for momentum added, some bugs removed
!
!  Revision 1.12  2001/08/01 08:31:22  bbh
!  CURVILINEAR now implemented
!
!  Revision 1.11  2001/07/26 13:46:06  bbh
!  Included info for salinity and temperature in namelist
!
!  Revision 1.10  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.9  2001/05/25 19:36:36  bbh
!  Added call to eqstate
!
!  Revision 1.8  2001/05/21 13:07:19  bbh
!  dt and cnpar is in variables_3d.F90
!
!  Revision 1.7  2001/05/20 07:51:40  bbh
!  Internal pressure included
!
!  Revision 1.6  2001/05/18 08:17:49  bbh
!  Tidying up initialization of salinity and temperature
!
!  Revision 1.5  2001/05/15 11:47:57  bbh
!  Added update_3d_halo for initial S and T
!
!  Revision 1.4  2001/05/10 11:30:16  bbh
!  Added further support for baroclinicity
!
!  Revision 1.3  2001/05/03 20:12:31  bbh
!  Use of variables_3d
!
!  Revision 1.2  2001/04/24 08:39:21  bbh
!  Included runtype as argument to integrate_3d
!
!  Revision 1.1.1.1  2001/04/17 08:43:08  bbh
!  initial import into CVS
!
! !LOCAL VARIABLES:
kbk's avatar
kbk committed
160
   integer                   :: vel_hor_adv=1,vel_ver_adv=1,vel_strang=0
gotm's avatar
gotm committed
161 162 163 164 165 166 167 168 169 170 171 172 173 174 175
!EOP
!-----------------------------------------------------------------------

   contains

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: init_3d - initialise 3D relatedstuff.
!
! !INTERFACE:
   subroutine init_3d(runtype,timestep,hotstart)
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
kbk's avatar
kbk committed
176 177 178
   integer, intent(in)                 :: runtype
   REALTYPE, intent(in)                :: timestep
   logical, intent(in)                 :: hotstart
gotm's avatar
gotm committed
179 180 181 182 183 184 185 186 187 188 189 190 191
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
!  Allocates memiory for 3D related fields.
!
! !REVISION HISTORY:
!
!  22Apr99   Karsten Bolding & Hans Burchard  Initial code.
!
! !LOCAL VARIABLES:
kbk's avatar
kbk committed
192 193 194 195 196 197
   integer                   :: rc
   NAMELIST /m3d/ &
             M,cnpar,cord_relax,                        &
             bdy3d,bdyfmt_3d,bdyramp_3d,bdyfile_3d,     &
             vel_hor_adv,vel_ver_adv,vel_strang,        &
             calc_temp,calc_salt,calc_spm
gotm's avatar
gotm committed
198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227
!
!EOP
!-------------------------------------------------------------------------
!BOC
#ifdef DEBUG
   integer, save :: Ncall = 0
   Ncall = Ncall+1
   write(debug,*) 'init_3d() # ',Ncall
#endif

   LEVEL1 'init_3d'
!  Read 3D-model specific things from the namelist.
   read(NAMLST,m3d)
!   rewind(NAMLST)

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

   LEVEL2 'vel_hor_adv= ',vel_hor_adv
   LEVEL2 'vel_ver_adv= ',vel_ver_adv
   LEVEL2 'vel_strang=  ',vel_strang

   if(vel_hor_adv .gt. 1 .or. vel_ver_adv .gt. 1) then
#ifndef UV_TVD
      STDERR 'To run the model with higher order advection for momentum'
      STDERR 'you need to re-compile the model with the option -DUV_TVD.'
      stop 'init_3d'
#endif
   end if
   dt = M*timestep
gotm's avatar
gotm committed
228 229
   num=avmmol
   nuh=avmmol
gotm's avatar
gotm committed
230 231

!  Needed for interpolation of temperature and salinity
kbk's avatar
kbk committed
232
   if (.not. hotstart) then
gotm's avatar
gotm committed
233 234
      call start_macro()
      call coordinates(vert_cord,cord_relax)
kbk's avatar
kbk committed
235
   end if
gotm's avatar
gotm committed
236

kbk's avatar
kbk committed
237
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
238
   if (runtype .eq. 3 .or. runtype .eq. 4) then
gotm's avatar
gotm committed
239
      T = _ZERO_ ; S = _ZERO_ ; rho = _ZERO_
gotm's avatar
gotm committed
240 241 242 243 244 245 246 247 248 249 250 251
      if(calc_temp) call init_temperature(1)
      if(calc_salt) call init_salinity(1)
      if(calc_spm)  call init_spm(1)
      call init_eqstate(1)
#ifndef PECS
      call do_eqstate()
#endif
      if (runtype .eq. 3) call internal_pressure()
      if (runtype .eq. 4) then
         call init_advection_3d(2)
      end if
   end if
kbk's avatar
kbk committed
252
#endif
gotm's avatar
gotm committed
253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273

   if (bdy3d) call init_bdy_3d()

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

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: integrate_3d - sequence of calls to do 3D model integration
!
! !INTERFACE:
   subroutine integrate_3d(runtype,n)
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
kbk's avatar
kbk committed
274
   integer, intent(in)                 :: runtype,n
gotm's avatar
gotm committed
275 276 277 278 279 280 281 282 283 284 285 286
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
! A wrapper to call all 3D related subroutines in one subroutine.
!
! !REVISION HISTORY:
!  See log for module
!
! !LOCAL VARIABLES:
kbk's avatar
kbk committed
287
  logical, save              :: ufirst=.true.
gotm's avatar
gotm committed
288 289 290 291 292 293 294 295 296 297
!
!EOP
!-------------------------------------------------------------------------
!BOC
#ifdef DEBUG
   integer, save :: Ncall = 0
   Ncall = Ncall+1
   write(debug,*) 'integrate_3d() # ',Ncall
#endif
   call start_macro()
kbk's avatar
kbk committed
298
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
299
   if (bdy3d) call do_bdy_3d(0,T)
kbk's avatar
kbk committed
300
#endif
gotm's avatar
gotm committed
301 302 303

   call coordinates(vert_cord,cord_relax)
#ifndef NO_BOTTFRIC
kbk's avatar
kbk committed
304
   if (kmax .gt. 1) then
gotm's avatar
gotm committed
305 306 307 308
      call bottom_friction_3d()
   end if
#endif
   SS = _ZERO_
kbk's avatar
kbk committed
309
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
310 311
   NN = _ZERO_
   if (runtype .eq. 4) call internal_pressure()
kbk's avatar
kbk committed
312
#endif
gotm's avatar
gotm committed
313 314 315 316 317 318 319 320 321 322 323 324 325
   if (ufirst) then
      call uu_momentum_3d(bdy3d)
      call vv_momentum_3d(bdy3d)
      ufirst=.false.
   else
      call vv_momentum_3d(bdy3d)
      call uu_momentum_3d(bdy3d)
      ufirst=.true.
   end if
   if (kmax .gt. 1) then
      call ww_momentum_3d()
   end if
#ifndef NO_ADVECT
kbk's avatar
kbk committed
326
   if (kmax .gt. 1) then
gotm's avatar
gotm committed
327 328 329 330 331 332 333 334
      call uv_advect_3d(vel_hor_adv,vel_ver_adv,vel_strang)
      if (Am .gt. _ZERO_) then
         call uv_diffusion_3d(Am)  ! Must be called after uv_advect_3d
      end if
   end if
#else
   STDERR 'NO_ADVECT 3D'
#endif
kbk's avatar
kbk committed
335

gotm's avatar
gotm committed
336 337
#ifdef CONST_VISC
   num = 1.000e-4
kbk's avatar
kbk committed
338
   nuh = 1.000e-5
gotm's avatar
gotm committed
339 340
#else
#ifndef NO_BOTTFRIC
kbk's avatar
kbk committed
341
   if (kmax .gt. 1) then
gotm's avatar
gotm committed
342 343 344 345 346 347 348 349
      call stresses_3d()
#ifndef PARABOLIC_VISCOSITY
      call ss_nn()
#endif
      call gotm()
   end if
#endif
#endif
kbk's avatar
kbk committed
350
#ifndef NO_BAROCLINIC
kbk's avatar
kbk committed
351
   if(runtype .eq. 4) then        ! prognostic T and S
gotm's avatar
gotm committed
352 353 354 355 356 357 358
      if (calc_temp) call do_temperature(n)
      if (calc_salt) call do_salinity(n)
      if (calc_spm) call do_spm()
#ifndef PECS
      call do_eqstate()
#endif
   end if
kbk's avatar
kbk committed
359
#endif
gotm's avatar
gotm committed
360

kbk's avatar
kbk committed
361
   if (kmax .gt. 1) then
gotm's avatar
gotm committed
362
#ifndef NO_BOTTFRIC
kbk's avatar
kbk committed
363
      call slow_bottom_friction()
gotm's avatar
gotm committed
364 365 366
#endif
#ifndef NO_ADVECT
#ifndef UV_ADV_DIRECT
kbk's avatar
kbk committed
367 368 369 370
      call slow_advection()
      if (Am .gt. _ZERO_) then
         call slow_diffusion(Am) ! Has to be called after slow_advection.
      end if
gotm's avatar
gotm committed
371 372 373
#endif
#endif
   end if
gotm's avatar
gotm committed
374

gotm's avatar
gotm committed
375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434
   call slow_terms()
   call stop_macro()

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

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: clean_3d - cleanup after 3D run.
!
! !INTERFACE:
   subroutine clean_3d()
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
!  This routine cleans up after a 3D 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_3d() # ',Ncall
#endif

   call clean_variables_3d()

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

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

   end module m3d

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