m3d.F90 15.3 KB
Newer Older
kbk's avatar
kbk committed
1
!$Id: m3d.F90,v 1.31 2006-03-17 11:06:33 kbk Exp $
gotm's avatar
gotm committed
2 3 4 5 6 7 8 9 10 11
#include "cppdefs.h"
!-----------------------------------------------------------------------
!BOP
!
! !MODULE: m3d - 3D model component
!
! !INTERFACE:
   module m3d
!
! !DESCRIPTION:
hb's avatar
hb committed
12
!  This module contains declarations for all variables related to 3D
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 3D model component.
hb's avatar
hb committed
17 18 19 20 21 22 23 24
!  The {\tt m3d} module is initialised in the routine {\tt init\_3d}, see
!  section \ref{sec-init-3d} described on page
!  \pageref{sec-init-3d}.
!  The actual calculation routines are called in {\tt integrate\_3d}
!  (see section \ref{sec-integrate-3d} on page \pageref{sec-integrate-3d}).
!  and are linked in from the library {\tt lib3d.a}.
!  After the simulation, the module is closed in {\tt clean\_3d}, see
!  section \ref{sec-clean-3d} on page \pageref{sec-clean-3d}.
gotm's avatar
gotm committed
25 26
!
! !USES:
27
   use exceptions
gotm's avatar
gotm committed
28
   use parameters, only: avmmol
kbk's avatar
kbk committed
29
   use domain, only: maxdepth,vert_cord
kbk's avatar
kbk committed
30
   use m2d, only: Am
kbk's avatar
kbk committed
31
   use variables_2d, only: D,z,UEx,VEx
kbk's avatar
kbk committed
32
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
33 34
   use temperature,only: init_temperature, do_temperature
   use salinity,   only: init_salinity, do_salinity
kbk's avatar
kbk committed
35
   use eqstate,    only: init_eqstate, do_eqstate
36 37
   use internal_pressure, only: init_internal_pressure, do_internal_pressure
   use internal_pressure, only: ip_method
kbk's avatar
kbk committed
38
#endif
39 40 41
#ifdef GETM_BIO
   use bio, only: bio_calc
   use getm_bio, only: init_getm_bio, do_getm_bio
kbk's avatar
kbk committed
42
#endif
43
   use variables_3d
gotm's avatar
gotm committed
44
   use advection_3d, only: init_advection_3d
kbk's avatar
kbk committed
45
   use bdy_3d, only: init_bdy_3d, do_bdy_3d
kbk's avatar
kbk committed
46

gotm's avatar
gotm committed
47 48 49
   IMPLICIT NONE
!
! !PUBLIC DATA MEMBERS:
kbk's avatar
kbk committed
50
   integer                             :: M=1
51
   REALTYPE                            :: cord_relax=_ZERO_
kbk's avatar
kbk committed
52 53 54 55 56
   logical                             :: calc_temp=.true.
   logical                             :: calc_salt=.true.
   logical                             :: bdy3d=.false.
   integer                             :: bdyfmt_3d,bdyramp_3d
   character(len=PATH_MAX)             :: bdyfile_3d
gotm's avatar
gotm committed
57
!
kbk's avatar
kbk committed
58 59 60
! !REVISION HISTORY:
!  Original author(s): Karsten Bolding & Hans Burchard
!
gotm's avatar
gotm committed
61
! !LOCAL VARIABLES:
62
   integer         :: vel_hor_adv=1,vel_ver_adv=1,vel_adv_split=0
kbk's avatar
kbk committed
63 64 65
#ifdef NO_BAROCLINIC
   integer         :: ip_method
#endif
gotm's avatar
gotm committed
66 67 68 69 70 71 72 73
!EOP
!-----------------------------------------------------------------------

   contains

!-----------------------------------------------------------------------
!BOP
!
hb's avatar
hb committed
74
! !IROUTINE: init_3d - initialise 3D related stuff \label{sec-init-3d}
gotm's avatar
gotm committed
75 76 77 78 79 80
!
! !INTERFACE:
   subroutine init_3d(runtype,timestep,hotstart)
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
kbk's avatar
kbk committed
81 82 83
   integer, intent(in)                 :: runtype
   REALTYPE, intent(in)                :: timestep
   logical, intent(in)                 :: hotstart
gotm's avatar
gotm committed
84 85 86 87 88 89
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
hb's avatar
hb committed
90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107
!  Here, the {\tt m3d} namelist is read from {\tt getm.inp}, and the
!  initialisation of variables is called (see routine {\tt init\_variables}
!  described on page \pageref{sec-init-variables}).
!  Furthermore, a number of consistency checks are made for the choices
!  of the momentum advection schemes. When higher-order advection schemes
!  are chosen for the momentum advection, the compiler option {\tt UV\_TVD}
!  has to be set. Here, the macro time step $\Delta t$ is calculated 
!  from the micro time step $\Delta t_m$ and the split factor {\tt M}.
!  Then, in order to have the vertical coordinate system present already here,
!  {\tt coordinates} (see page \pageref{sec-coordinates}) needs to be called,
!  in order to enable proper interpolation of initial values for 
!  potential temperature $\theta$ and salinity $S$ for cold starts. 
!  Those initial values are afterwards read in via the routines
!  {\tt init\_temperature} (page \pageref{sec-init-temperature}) and
!  {\tt init\_salinity} (page \pageref{sec-init-salinity}).
!  Finally, in order to prepare for the first time step, the momentum advection
!  and internal pressure gradient routines are initialised and the
!  internal pressure gradient routine is called.
gotm's avatar
gotm committed
108 109
!
! !LOCAL VARIABLES:
110
   integer         :: rc
kbk's avatar
kbk committed
111 112 113
   NAMELIST /m3d/ &
             M,cnpar,cord_relax,                        &
             bdy3d,bdyfmt_3d,bdyramp_3d,bdyfile_3d,     &
114
             vel_hor_adv,vel_ver_adv,vel_adv_split,     &
kbk's avatar
kbk committed
115
             calc_temp,calc_salt,                       &
116
             avmback,avhback,ip_method
gotm's avatar
gotm committed
117 118 119 120 121 122 123 124 125 126 127 128 129 130 131
!
!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)

132 133 134 135 136 137 138 139 140
   if (avmback .lt. _ZERO_) then
      LEVEL2 "setting avmback to 0."
      avmback = _ZERO_
   end if
   if (avhback .lt. _ZERO_) then
      LEVEL2 "setting avhback to 0."
      avhback = _ZERO_
   end if

gotm's avatar
gotm committed
141 142 143
! Allocates memory for the public data members - if not static
   call init_variables_3d(runtype)

144 145 146 147 148
!  Sanity checks for advection specifications
   LEVEL2 'vel_hor_adv=   ',vel_hor_adv
   LEVEL2 'vel_ver_adv=   ',vel_ver_adv
   LEVEL2 'vel_adv_split= ',vel_adv_split
   if(vel_hor_adv .gt. 1) then
gotm's avatar
gotm committed
149 150 151
#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.'
152
      call getm_error("init_3d()","GETM needs recompilation")
gotm's avatar
gotm committed
153
#endif
154 155 156 157 158 159
   else
      vel_adv_split=-1
      if(vel_ver_adv .ne. 1) then
         LEVEL2 "setting vel_ver_adv to 1 - since vel_hor_adv is 1"
         vel_ver_adv=1
      end if
gotm's avatar
gotm committed
160
   end if
161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210
   LEVEL2 "horizontal: ",trim(adv_schemes(vel_hor_adv))," of momentum"
   LEVEL2 "vertical:   ",trim(adv_schemes(vel_ver_adv))," of momentum"

   select case (vel_adv_split)
      case (-1)
      case (0)
         select case (vel_hor_adv)
            case (2,3,4,5,6)
            case default
               call getm_error("init_3d()", &
                    "vel_adv_split=0: vel_hor_adv not valid (2-6)")
         end select
         select case (vel_ver_adv)
            case (2,3,4,5,6)
            case default
               call getm_error("init_3d()", &
                    "vel_adv_split=0: vel_ver_adv not valid (2-6)")
         end select
         LEVEL2 "1D split --> full u, full v, full w"
      case (1)
         select case (vel_hor_adv)
            case (2,3,4,5,6)
            case default
               call getm_error("init_3d()", &
                    "vel_adv_split=1: vel_hor_adv not valid (2-6)")
         end select
         select case (vel_ver_adv)
            case (2,3,4,5,6)
            case default
               call getm_error("init_3d()", &
                    "vel_adv_split=1: vel_ver_adv not valid (2-6)")
         end select
         LEVEL2 "1D split --> half u, half v, full w, half v, half u"
      case (2)
         select case (vel_hor_adv)
            case (2,7)
            case default
               call getm_error("init_3d()", &
                    "vel_adv_split=2: vel_hor_adv not valid (2,7)")
         end select
         select case (vel_ver_adv)
            case (2,3,4,5,6)
            case default
               call getm_error("init_3d()", &
                    "vel_adv_split=2: vel_ver_adv not valid (2-6)")
         end select
         LEVEL2 "2D-hor, 1D-vert split --> full uv, full w"
      case default
   end select

gotm's avatar
gotm committed
211
   dt = M*timestep
212 213 214 215 216 217 218 219

#ifdef CONSTANT_VISCOSITY
   num=avmback
   nuh=avhback
#else
   num=1.e-15
   nuh=1.e-15
#endif
gotm's avatar
gotm committed
220 221

!  Needed for interpolation of temperature and salinity
kbk's avatar
kbk committed
222
   if (.not. hotstart) then
gotm's avatar
gotm committed
223
      call start_macro()
kbk's avatar
kbk committed
224
      call coordinates(vert_cord,cord_relax,maxdepth)
225
      call hcc_check()
kbk's avatar
kbk committed
226
   end if
gotm's avatar
gotm committed
227

kbk's avatar
kbk committed
228
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
229
   if (runtype .eq. 3 .or. runtype .eq. 4) then
gotm's avatar
gotm committed
230
      T = _ZERO_ ; S = _ZERO_ ; rho = _ZERO_
gotm's avatar
gotm committed
231 232
      if(calc_temp) call init_temperature(1)
      if(calc_salt) call init_salinity(1)
kbk's avatar
kbk committed
233
   end if
234
#endif
kbk's avatar
kbk committed
235

236
#ifndef NO_BAROCLINIC
kbk's avatar
kbk committed
237
    if (runtype .eq. 3 .or. runtype .eq. 4) then
238
      call init_eqstate()
gotm's avatar
gotm committed
239 240 241
#ifndef PECS
      call do_eqstate()
#endif
242 243
      if (runtype .ge. 3) call init_internal_pressure()
      if (runtype .eq. 3) call do_internal_pressure()
gotm's avatar
gotm committed
244 245 246 247
      if (runtype .eq. 4) then
         call init_advection_3d(2)
      end if
   end if
kbk's avatar
kbk committed
248
#endif
249 250 251 252

#ifdef GETM_BIO
!KBK      if(bio_calc) call init_getm_bio()
      call init_getm_bio()
kbk's avatar
kbk committed
253
#endif
gotm's avatar
gotm committed
254 255 256 257 258 259 260 261 262 263 264 265 266 267

   if (bdy3d) call init_bdy_3d()

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

!-----------------------------------------------------------------------
!BOP
!
hb's avatar
hb committed
268 269
! !IROUTINE: integrate_3d - calls to do 3D model integration
! \label{sec-integrate-3d}
gotm's avatar
gotm committed
270 271 272 273 274 275
!
! !INTERFACE:
   subroutine integrate_3d(runtype,n)
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
kbk's avatar
kbk committed
276
   integer, intent(in)                 :: runtype,n
gotm's avatar
gotm committed
277 278 279 280 281 282
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
hb's avatar
hb committed
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 341 342 343 344 345 346 347 348 349 350 351
! This is a wrapper routine to call all 3D related subroutines.
! The call position for the {\tt coordinates} routine depends on 
! the compiler option 
! {\tt MUDFLAT}: If it is defined, then the 
! call to {\tt coordinates} construction is made such that drying and flooding
! is stable. If {\tt MUDFLAT} is not defined, then the adaptive grids with
! Lagrangian component which are currently under development are supported. 
! Both, drying and flooding and
! Lagrangian coordinates does not go together yet.
! The call sequence is as follows:
!
! \vspace{0.5cm}
! 
! \begin{tabular}{lll}
! {\tt start\_macro}           & initialising a 3d step & see page
! \pageref{sec-start-macro} \\
! {\tt do\_bdy\_3d}            & boundary conditions for $\theta$ and $S$ & see
! page \pageref{sec-do-bdy-3d} \\
! {\tt coordinates}            & layer heights ({\tt MUTFLAT} defined) & see
! page \pageref{sec-coordinates} \\
! {\tt bottom\_friction\_3d}   & bottom friction & see page
! \pageref{sec-bottom-friction-3d} \\
! {\tt do\_internal\_pressure} & internal pressure gradient & see page
! \pageref{sec-do-internal-pressure} \\
! {\tt uu\_momentum\_3d}       & layer-integrated $u$-velocity & see page
! \pageref{sec-uu-momentum-3d} \\
! {\tt vv\_momentum\_3d}       & layer-integrated $v$-velocity & see page
! \pageref{sec-vv-momentum-3d} \\
! {\tt coordinates}            & layer heights ({\tt MUTFLAT} not defined) & see
! page \pageref{sec-coordinates} \\
! {\tt ww\_momentum\_3d}       & grid-related vertical velocity & see page
! \pageref{sec-ww-momentum-3d} \\
! {\tt uv\_advect\_3d}         & momentum advection & see page
! \pageref{sec-uv-advect-3d} \\
! {\tt uv\_diffusion\_3d}      & momentum diffusion & see page
! \pageref{sec-uv-diffusion-3d} \\
! {\tt stresses\_3d}           & stresses (for GOTM) & see page
! \pageref{sec-stresses-3d} \\
! {\tt ss\_nn}                 & shear and stratification (for GOTM) & see page
! \pageref{sec-ss-nn} \\
! {\tt gotm}                   & interface and call to GOTM & see page
! \pageref{sec-gotm} \\
! {\tt do\_temperature}        & potential temperature equation & see page
! \pageref{sec-do-temperature} \\
! {\tt do\_salinity}           & salinity equation & see page
! \pageref{sec-do-salinity} \\
! {\tt do\_eqstate}            & equation of state & see page
! \pageref{sec-do-eqstate} \\
! {\tt do\_spm}                & suspended matter equation & see page 
! \pageref{sec-do-spm} \\
! {\tt do\_getm\_bio}          & call to GOTM-BIO (not yet released) & \\
! {\tt slow\_bottom\_friction} & slow bottom friction & see page
! \pageref{sec-slow-bottom-friction} \\
! {\tt slow\_advection}        & slow advection terms & see page
! \pageref{sec-slow-advection} \\
! {\tt slow\_diffusion}        & slow diffusion terms & see page
! \pageref{sec-slow-diffusion} \\
! {\tt slow\_terms}            & sum of slow terms & see page
! \pageref{sec-slow-terms} \\
! {\tt stop\_macro}            & finishing a 3d step & see page
! \pageref{sec-stop-macro}
! \end{tabular}
! 
! \vspace{0.5cm}
!
! Several calls are only executed for certain compiler options. At each
! time step the call sequence for the horizontal momentum equations is
! changed in order to allow for higher order accuracy for the Coriolis
! rotation.
gotm's avatar
gotm committed
352 353
!
! !LOCAL VARIABLES:
kbk's avatar
kbk committed
354
  logical, save              :: ufirst=.true.
gotm's avatar
gotm committed
355 356 357 358 359 360 361 362 363 364
!
!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
365
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
366
   if (bdy3d) call do_bdy_3d(0,T)
kbk's avatar
kbk committed
367
#endif
368 369 370
#ifdef MUDFLAT
   call coordinates(vert_cord,cord_relax,maxdepth)
#endif
gotm's avatar
gotm committed
371
#ifndef NO_BOTTFRIC
kbk's avatar
kbk committed
372
   if (kmax .gt. 1) then
gotm's avatar
gotm committed
373 374 375 376
      call bottom_friction_3d()
   end if
#endif
   SS = _ZERO_
kbk's avatar
kbk committed
377
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
378
   NN = _ZERO_
379
   if (runtype .eq. 4) call do_internal_pressure()
kbk's avatar
kbk committed
380
#endif
381 382
   huo=hun
   hvo=hvn
gotm's avatar
gotm committed
383 384 385 386 387 388 389 390 391
   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
392
#ifndef MUDFLAT
393
   call coordinates(vert_cord,cord_relax,maxdepth)
394
#endif
gotm's avatar
gotm committed
395 396 397 398
   if (kmax .gt. 1) then
      call ww_momentum_3d()
   end if
#ifndef NO_ADVECT
kbk's avatar
kbk committed
399
   if (kmax .gt. 1) then
400
      call uv_advect_3d(vel_hor_adv,vel_ver_adv,vel_adv_split)
gotm's avatar
gotm committed
401 402 403 404 405 406 407
      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
408 409

   if (kmax .gt. 1) then
410
#ifndef NO_BOTTFRIC
gotm's avatar
gotm committed
411
      call stresses_3d()
412 413
#endif
#ifndef CONSTANT_VISCOSITY
gotm's avatar
gotm committed
414 415 416 417 418
#ifndef PARABOLIC_VISCOSITY
      call ss_nn()
#endif
      call gotm()
#endif
419
   end if
kbk's avatar
kbk committed
420
#ifndef NO_BAROCLINIC
kbk's avatar
kbk committed
421
   if(runtype .eq. 4) then        ! prognostic T and S
gotm's avatar
gotm committed
422 423
      if (calc_temp) call do_temperature(n)
      if (calc_salt) call do_salinity(n)
kbk's avatar
kbk committed
424
   end if
425
#endif
kbk's avatar
kbk committed
426

427 428
#ifndef NO_BAROCLINIC
   if(runtype .eq. 4) then
gotm's avatar
gotm committed
429 430 431 432
#ifndef PECS
      call do_eqstate()
#endif
   end if
kbk's avatar
kbk committed
433
#endif
gotm's avatar
gotm committed
434

435 436
#ifdef GETM_BIO
   if (bio_calc) call do_getm_bio(dt)
kbk's avatar
kbk committed
437 438
#endif

kbk's avatar
kbk committed
439
   UEx=_ZERO_ ; VEx=_ZERO_
kbk's avatar
kbk committed
440
   if (kmax .gt. 1) then
gotm's avatar
gotm committed
441
#ifndef NO_BOTTFRIC
kbk's avatar
kbk committed
442
      call slow_bottom_friction()
gotm's avatar
gotm committed
443 444 445
#endif
#ifndef NO_ADVECT
#ifndef UV_ADV_DIRECT
kbk's avatar
kbk committed
446 447 448 449
      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
450 451 452
#endif
#endif
   end if
gotm's avatar
gotm committed
453

gotm's avatar
gotm committed
454 455 456 457 458 459 460 461 462 463 464 465 466 467
   call slow_terms()
   call stop_macro()

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

!-----------------------------------------------------------------------
!BOP
!
hb's avatar
hb committed
468
! !IROUTINE: clean_3d - cleanup after 3D run \label{sec-clean-3d}
gotm's avatar
gotm committed
469 470 471 472 473 474 475 476 477 478 479 480
!
! !INTERFACE:
   subroutine clean_3d()
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
hb's avatar
hb committed
481 482
! Here, a call to the routine {\tt clean\_variables\_3d} which howewer
! does not do anything yet.
gotm's avatar
gotm committed
483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511
!
! !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)         !
!-----------------------------------------------------------------------