m3d.F90 22.7 KB
Newer Older
gotm's avatar
gotm committed
1 2 3 4 5 6 7 8 9 10
#include "cppdefs.h"
!-----------------------------------------------------------------------
!BOP
!
! !MODULE: m3d - 3D model component
!
! !INTERFACE:
   module m3d
!
! !DESCRIPTION:
hb's avatar
hb committed
11
!  This module contains declarations for all variables related to 3D
gotm's avatar
gotm committed
12
!  hydrodynamical calculations. Information about the calculation domain
hb's avatar
hb committed
13
!  is included from the {\tt domain} module.
gotm's avatar
gotm committed
14 15
!  The module contains public subroutines for initialisation, integration
!  and clean up of the 3D model component.
hb's avatar
hb committed
16 17 18 19 20 21 22 23
!  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
24 25
!
! !USES:
26
   use exceptions
gotm's avatar
gotm committed
27
   use parameters, only: avmmol
28
   use domain, only: openbdy,maxdepth,vert_cord,az
29
   use m2d, only: depth_update,uv_advect,uv_diffusion
30
   use variables_2d, only: z,Uint,Vint,UEx,VEx
kbk's avatar
kbk committed
31
#ifndef NO_BAROCLINIC
32 33 34
   use temperature,only: init_temperature, do_temperature, &
            init_temperature_field
   use salinity,   only: init_salinity, do_salinity, init_salinity_field
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
   use variables_3d
40 41
   use advection, only: NOADV
   use advection_3d, only: init_advection_3d,print_adv_settings_3d,adv_ver_iterations
kbk's avatar
kbk committed
42
   use bdy_3d, only: init_bdy_3d, do_bdy_3d
43
   use bdy_3d, only: bdy3d_tmrlx, bdy3d_tmrlx_ucut, bdy3d_tmrlx_max, bdy3d_tmrlx_min
44
!  Necessary to use halo_zones because update_3d_halos() have been moved out
kbk's avatar
kbk committed
45 46
!  temperature.F90 and salinity.F90 - should be changed at a later stage
   use halo_zones, only: update_3d_halo,wait_halo,D_TAG
kbk's avatar
kbk committed
47

gotm's avatar
gotm committed
48 49 50
   IMPLICIT NONE
!
! !PUBLIC DATA MEMBERS:
kbk's avatar
kbk committed
51
   integer                             :: M=1
52
   REALTYPE                            :: cord_relax=_ZERO_
53 54 55
   integer                             :: vel3d_adv_split=0
   integer                             :: vel3d_adv_hor=1
   integer                             :: vel3d_adv_ver=1
56 57 58
   integer                             :: turb_adv_split=0
   integer                             :: turb_adv_hor=0
   integer                             :: turb_adv_ver=0
Knut's avatar
Knut committed
59 60
   logical                             :: smooth_bvf_hor=.false.
   logical                             :: smooth_bvf_ver=.false.
61 62
   logical                             :: calc_temp=.true.
   logical                             :: calc_salt=.true.
kbk's avatar
kbk committed
63
   logical                             :: bdy3d=.false.
64
   integer                             :: bdyfmt_3d,bdy3d_ramp
kbk's avatar
kbk committed
65
   character(len=PATH_MAX)             :: bdyfile_3d
kbk's avatar
kbk committed
66
   REALTYPE                            :: ip_fac=_ONE_
67
   integer                             :: vel_check=0
68
   REALTYPE                            :: min_vel=-4*_ONE_,max_vel=4*_ONE_
gotm's avatar
gotm committed
69
!
kbk's avatar
kbk committed
70 71 72
! !REVISION HISTORY:
!  Original author(s): Karsten Bolding & Hans Burchard
!
gotm's avatar
gotm committed
73
! !LOCAL VARIABLES:
Knut's avatar
Knut committed
74
   logical,private :: ufirst=.true.
75
   logical         :: advect_turbulence=.false.
kbk's avatar
kbk committed
76 77 78
#ifdef NO_BAROCLINIC
   integer         :: ip_method
#endif
kbk's avatar
kbk committed
79
   integer         :: ip_ramp=-1
gotm's avatar
gotm committed
80 81 82 83 84 85 86 87
!EOP
!-----------------------------------------------------------------------

   contains

!-----------------------------------------------------------------------
!BOP
!
hb's avatar
hb committed
88
! !IROUTINE: init_3d - initialise 3D related stuff \label{sec-init-3d}
gotm's avatar
gotm committed
89 90 91 92 93 94
!
! !INTERFACE:
   subroutine init_3d(runtype,timestep,hotstart)
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
95
   integer, intent(in)                 :: runtype
kbk's avatar
kbk committed
96 97
   REALTYPE, intent(in)                :: timestep
   logical, intent(in)                 :: hotstart
gotm's avatar
gotm committed
98 99 100
!
!
! !DESCRIPTION:
hb's avatar
hb committed
101 102 103 104 105 106
!  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}
107
!  has to be set. Here, the macro time step $\Delta t$ is calculated
hb's avatar
hb committed
108 109 110
!  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,
111 112
!  in order to enable proper interpolation of initial values for
!  potential temperature $\theta$ and salinity $S$ for cold starts.
hb's avatar
hb committed
113 114 115 116 117 118
!  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
119 120
!
! !LOCAL VARIABLES:
121
   integer         :: rc
kbk's avatar
kbk committed
122
   NAMELIST /m3d/ &
123
             M,cnpar,cord_relax,adv_ver_iterations,       &
124
             bdy3d,bdyfmt_3d,bdy3d_ramp,bdyfile_3d,       &
125 126 127
             bdy3d_tmrlx,bdy3d_tmrlx_ucut,                &
             bdy3d_tmrlx_max,bdy3d_tmrlx_min,             &
             vel3d_adv_split,vel3d_adv_hor,vel3d_adv_ver, &
128
             turb_adv_split,turb_adv_hor,turb_adv_ver,    &
129
             calc_temp,calc_salt,                         &
Knut's avatar
Knut committed
130
             avmback,avhback,smooth_bvf_hor,smooth_bvf_ver, &
131
             ip_method,ip_ramp,                           &
132
             vel_check,min_vel,max_vel
gotm's avatar
gotm committed
133 134 135 136 137 138 139 140 141 142 143 144 145 146
!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)

Knut's avatar
Knut committed
147 148
   LEVEL2 "splitting factor M: ",M

149 150 151 152 153 154 155 156 157
   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
158 159
! Allocates memory for the public data members - if not static
   call init_variables_3d(runtype)
160
   call init_advection_3d()
gotm's avatar
gotm committed
161

162
!  Sanity checks for advection specifications
163 164 165 166 167 168 169 170 171 172 173 174 175
   LEVEL2 'Advection of horizontal 3D velocities'
#ifdef NO_ADVECT
   if (vel3d_adv_hor .ne. NOADV) then
      LEVEL2 "reset vel3d_adv_hor= ",NOADV," because of"
      LEVEL2 "obsolete NO_ADVECT macro. Note that this"
      LEVEL2 "behaviour will be removed in the future."
      vel3d_adv_hor = NOADV
   end if
   if (vel3d_adv_ver .ne. NOADV) then
      LEVEL2 "reset vel3d_adv_ver= ",NOADV," because of"
      LEVEL2 "obsolete NO_ADVECT macro. Note that this"
      LEVEL2 "behaviour will be removed in the future."
      vel3d_adv_ver = NOADV
gotm's avatar
gotm committed
176
   end if
177 178
#endif
   call print_adv_settings_3d(vel3d_adv_split,vel3d_adv_hor,vel3d_adv_ver,_ZERO_)
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

#ifndef CONSTANT_VISCOSITY

   LEVEL2 'Advection of TKE and eps'

#ifdef NO_ADVECT
   if (turb_adv_hor .ne. NOADV) then
      LEVEL2 "reset turb_adv_hor= ",NOADV," because of"
      LEVEL2 "obsolete NO_ADVECT macro. Note that this"
      LEVEL2 "behaviour will be removed in the future."
      turb_adv_hor = NOADV
   end if
   if (turb_adv_ver .ne. NOADV) then
      LEVEL2 "reset turb_adv_ver= ",NOADV," because of"
      LEVEL2 "obsolete NO_ADVECT macro. Note that this"
      LEVEL2 "behaviour will be removed in the future."
      turb_adv_ver = NOADV
   end if
#endif
   call print_adv_settings_3d(turb_adv_split,turb_adv_hor,turb_adv_ver,_ZERO_)

   advect_turbulence = (turb_adv_hor.ne.NOADV .or. turb_adv_ver.ne.NOADV)

#ifdef TURB_ADV
   if (.not. advect_turbulence) then
      LEVEL2 "WARNING: ignored obsolete TURB_ADV macro!"
   end if
#endif

Knut's avatar
Knut committed
209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228
#ifdef SMOOTH_BVF_HOR
   if (.not. smooth_bvf_hor) then
      LEVEL2 "reset smooth_bvf_hor=T because of obsolete"
      LEVEL2 "SMOOTH_BVF_HOR macro. Note that this"
      LEVEL2 "behaviour will be removed in the future."
      smooth_bvf_hor = .true.
   end if
#endif
   LEVEL2 "smooth_bvf_hor = ",smooth_bvf_hor

#ifdef _SMOOTH_BVF_VERT_
   if (.not. smooth_bvf_ver) then
      LEVEL2 "reset smooth_bvf_ver=T because of obsolete"
      LEVEL2 "_SMOOTH_BVF_VERT_ macro. Note that this"
      LEVEL2 "behaviour will be removed in the future."
      smooth_bvf_ver = .true.
   end if
#endif
   LEVEL2 "smooth_bvf_ver = ",smooth_bvf_ver

229 230 231
#endif


bjb's avatar
bjb committed
232
!  Sanity checks for bdy 3d
Knut's avatar
Knut committed
233
   if (.not.openbdy .or. runtype.eq.2) bdy3d=.false.
234 235 236 237 238 239 240 241 242 243 244 245 246 247
   if (bdy3d .and. bdy3d_tmrlx) then
      LEVEL2 'bdy3d_tmrlx=.true.'
      LEVEL2 'bdy3d_tmrlx_max=   ',bdy3d_tmrlx_max
      LEVEL2 'bdy3d_tmrlx_min=   ',bdy3d_tmrlx_min
      LEVEL2 'bdy3d_tmrlx_ucut=  ',bdy3d_tmrlx_ucut
      if (bdy3d_tmrlx_min<_ZERO_ .or. bdy3d_tmrlx_min>_ONE_)          &
           call getm_error("init_3d()",                               &
           "bdy3d_tmrlx_min is out of valid range [0:1]")
      if (bdy3d_tmrlx_max<bdy3d_tmrlx_min .or. bdy3d_tmrlx_min>_ONE_) &
           call getm_error("init_3d()",                               &
           "bdy3d_tmrlx_max is out of valid range [bdy3d_tmrlx_max:1]")
      if (bdy3d_tmrlx_ucut<_ZERO_)                                    &
           call getm_error("init_3d()",                               &
           "bdy3d_tmrlx_max is out of valid range [0:inf[")
bjb's avatar
bjb committed
248 249
   end if

kbk's avatar
kbk committed
250 251
   LEVEL2 'ip_ramp=',ip_ramp

252 253
   LEVEL2 'vel_check=',vel_check
   if (vel_check .ne. 0) then
kbk's avatar
kbk committed
254 255 256
      LEVEL3 'doing sanity checks on velocities'
      LEVEL3 'min_vel=',min_vel
      LEVEL3 'max_vel=',max_vel
257 258 259 260 261 262 263 264
      if (vel_check .gt. 0) then
         LEVEL3 'out-of-bound values result in termination of program'
      end if
      if (vel_check .lt. 0) then
         LEVEL3 'out-of-bound values result in warnings only'
      end if
   end if

gotm's avatar
gotm committed
265
   dt = M*timestep
266 267 268 269 270

#ifdef CONSTANT_VISCOSITY
   num=avmback
   nuh=avhback
#else
271 272
   num=1.d-15
   nuh=1.d-15
273
#endif
gotm's avatar
gotm committed
274 275

!  Needed for interpolation of temperature and salinity
kbk's avatar
kbk committed
276
   if (.not. hotstart) then
277
      ssen = z
gotm's avatar
gotm committed
278
      call start_macro()
279
      call coordinates(hotstart)
280
      call hcc_check()
kbk's avatar
kbk committed
281
   end if
gotm's avatar
gotm committed
282

283 284 285 286
   if (runtype .eq. 2) then
      calc_temp = .false.
      calc_salt = .false.
   else
kbk's avatar
kbk committed
287
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
288
      T = _ZERO_ ; S = _ZERO_ ; rho = _ZERO_
289 290
      if(calc_temp) call init_temperature()
      if(calc_salt) call init_salinity()
291
#endif
292
   end if
kbk's avatar
kbk committed
293

294
#ifndef NO_BAROCLINIC
kbk's avatar
kbk committed
295
    if (runtype .eq. 3 .or. runtype .eq. 4) then
296
      call init_eqstate()
gotm's avatar
gotm committed
297
      call do_eqstate()
298
      call buoyancy_frequency()
299 300

      if (bdy3d) call init_bdy_3d()
301 302
      if (runtype .ge. 3) call init_internal_pressure()
      if (runtype .eq. 3) call do_internal_pressure()
gotm's avatar
gotm committed
303
   end if
304
#endif
305

306
   if (vert_cord .eq. _ADAPTIVE_COORDS_) call preadapt_coordinates(preadapt)
kb's avatar
kb committed
307

gotm's avatar
gotm committed
308 309 310 311 312 313 314 315
#ifdef DEBUG
   write(debug,*) 'Leaving init_3d()'
   write(debug,*)
#endif
   return
   end subroutine init_3d
!EOC

316 317 318 319 320 321
!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: postinit_3d - re-initialise some 3D after hotstart read.
!
! !INTERFACE:
Knut's avatar
Knut committed
322
   subroutine postinit_3d(runtype,timestep,hotstart,MinN)
323 324 325 326 327
! !USES:
   use domain, only: imin,imax,jmin,jmax, az,au,av
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
Knut's avatar
Knut committed
328
   integer, intent(in)                 :: runtype,MinN
329 330
   REALTYPE, intent(in)                :: timestep
   logical, intent(in)                 :: hotstart
331 332 333 334 335 336
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
337 338
!  This routine provides possibility to reset/initialize 3D variables to
!  ensure that velocities are correctly set on land cells after read
339 340 341
!  of a hotstart file.
!
! !LOCAL VARIABLES:
342
   integer                   :: i,j,rc
343 344 345 346 347 348 349 350 351 352 353
!EOP
!-------------------------------------------------------------------------
!BOC
#ifdef DEBUG
   integer, save :: Ncall = 0
   Ncall = Ncall+1
   write(debug,*) 'postinit_3d() # ',Ncall
#endif

   LEVEL1 'postinit_3d'

Knut's avatar
Knut committed
354 355
   ufirst = ( mod(int(ceiling((_ONE_*MinN)/M)),2) .eq. 1 )

356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396
   if (do_numerical_analyses) then

      allocate(numdis3d(I3DFIELD),stat=rc)
      if (rc /= 0) stop 'postinit_3d: Error allocating memory (numdis3d)'
      numdis3d = _ZERO_
      allocate(numdis2d(I2DFIELD),stat=rc)
      if (rc /= 0) stop 'postinit_3d: Error allocating memory (numdis2d)'
      numdis2d = _ZERO_

      if (calc_temp) then
         allocate(phymix3d_T(I3DFIELD),stat=rc)
         if (rc /= 0) stop 'postinit_3d: Error allocating memory (phymix3d_T)'
         phymix3d_T = _ZERO_
         allocate(phymix2d_T(I2DFIELD),stat=rc)
         if (rc /= 0) stop 'postinit_3d: Error allocating memory (phymix2d_T)'
         phymix2d_T = _ZERO_
         allocate(nummix3d_T(I3DFIELD),stat=rc)
         if (rc /= 0) stop 'postinit_3d: Error allocating memory (nummix3d_T)'
         nummix3d_T = _ZERO_
         allocate(nummix2d_T(I2DFIELD),stat=rc)
         if (rc /= 0) stop 'postinit_3d: Error allocating memory (nummix2d_T)'
         nummix2d_T = _ZERO_
      end if

      if (calc_salt) then
         allocate(phymix3d_S(I3DFIELD),stat=rc)
         if (rc /= 0) stop 'postinit_3d: Error allocating memory (phymix3d_S)'
         phymix3d_S = _ZERO_
         allocate(phymix2d_S(I2DFIELD),stat=rc)
         if (rc /= 0) stop 'postinit_3d: Error allocating memory (phymix2d_S)'
         phymix2d_S = _ZERO_
         allocate(nummix3d_S(I3DFIELD),stat=rc)
         if (rc /= 0) stop 'postinit_3d: Error allocating memory (nummix3d_S)'
         nummix3d_S = _ZERO_
         allocate(nummix2d_S(I2DFIELD),stat=rc)
         if (rc /= 0) stop 'postinit_3d: Error allocating memory (nummix2d_S)'
         nummix2d_S = _ZERO_
      end if

   end if

397 398
! Hotstart fix - see postinit_2d
   if (hotstart) then
399

400 401 402 403 404 405
      do j=jmin-HALO,jmax+HALO
         do i=imin-HALO,imax+HALO
            if (au(i,j) .eq. 0) then
               uu(i,j,:)  = _ZERO_
            end if
         end do
406
      end do
407 408 409 410 411 412
      do j=jmin-HALO,jmax+HALO
         do i=imin-HALO,imax+HALO
            if (av(i,j) .eq. 0) then
               vv(i,j,:)  = _ZERO_
            end if
         end do
413
      end do
414 415 416 417 418 419 420 421 422 423 424 425 426
!     These may not be necessary, but we clean up anyway just in case.
      do j=jmin-HALO,jmax+HALO
         do i=imin-HALO,imax+HALO
            if(az(i,j) .eq. 0) then
               tke(i,j,:) = _ZERO_
               num(i,j,:) = 1.e-15
               nuh(i,j,:) = 1.e-15
#ifndef NO_BAROCLINIC
               S(i,j,:)   = _ZERO_
               T(i,j,:)   = _ZERO_
#endif
            end if
         end do
427
      end do
428

Knut's avatar
Knut committed
429
      call depth_update(sseo,ssen,Dn,Dveln,Dun,Dvn,from3d=.true.)
430
!     KK-TODO: do not store ss[u|v]n in hotstart file
Knut's avatar
Knut committed
431
!              can be calculated here (if needed at all... use of D[u|v]n)
432 433
!     ssun = Dun - HU
!     ssvn = Dvn - HV
434 435
      call coordinates(hotstart)

436 437
      if (vert_cord .eq. _ADAPTIVE_COORDS_) call shear_frequency()

438
   end if
439 440 441 442 443

   return
   end subroutine postinit_3d
!EOC

gotm's avatar
gotm committed
444 445 446
!-----------------------------------------------------------------------
!BOP
!
hb's avatar
hb committed
447 448
! !IROUTINE: integrate_3d - calls to do 3D model integration
! \label{sec-integrate-3d}
gotm's avatar
gotm committed
449 450 451
!
! !INTERFACE:
   subroutine integrate_3d(runtype,n)
bjb's avatar
bjb committed
452
   use getm_timers, only: tic, toc, TIM_INTEGR3D
453 454 455
#ifndef NO_BAROCLINIC
   use getm_timers, only: TIM_TEMPH, TIM_SALTH
#endif
gotm's avatar
gotm committed
456 457 458
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
kbk's avatar
kbk committed
459
   integer, intent(in)                 :: runtype,n
gotm's avatar
gotm committed
460 461 462 463 464 465
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
hb's avatar
hb committed
466
! This is a wrapper routine to call all 3D related subroutines.
467 468
! The call position for the {\tt coordinates} routine depends on
! the compiler option
hb's avatar
hb committed
469 470 471
! The call sequence is as follows:
!
! \vspace{0.5cm}
472
!
hb's avatar
hb committed
473 474 475 476 477 478 479 480 481 482 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
! \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} \\
508
! {\tt do\_spm}                & suspended matter equation & see page
hb's avatar
hb committed
509 510 511 512 513 514 515 516 517 518 519 520 521
! \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}
522
!
hb's avatar
hb committed
523 524 525 526 527 528
! \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
529 530 531 532 533 534 535 536 537 538 539 540
!
! !LOCAL VARIABLES:
!
!EOP
!-------------------------------------------------------------------------
!BOC
#ifdef DEBUG
   integer, save :: Ncall = 0
   Ncall = Ncall+1
   write(debug,*) 'integrate_3d() # ',Ncall
#endif
   call start_macro()
541

542
   call coordinates(.false.)
Knut's avatar
Knut committed
543

gotm's avatar
gotm committed
544
#ifndef NO_BOTTFRIC
kbk's avatar
kbk committed
545
   if (kmax .gt. 1) then
gotm's avatar
gotm committed
546 547 548
      call bottom_friction_3d()
   end if
#endif
kbk's avatar
kbk committed
549
#ifndef NO_BAROCLINIC
550
   if (runtype .eq. 4) call do_internal_pressure()
kbk's avatar
kbk committed
551
#endif
Knut's avatar
Knut committed
552

bjb's avatar
bjb committed
553
   call tic(TIM_INTEGR3D)
kbk's avatar
kbk committed
554
   ip_fac=_ONE_
frv-bjb's avatar
frv-bjb committed
555
   if (ip_ramp .gt. 0) ip_fac=min( _ONE_ , n*_ONE_/ip_ramp)
bjb's avatar
bjb committed
556
   call toc(TIM_INTEGR3D)
Knut's avatar
Knut committed
557

558 559 560
#ifdef STRUCTURE_FRICTION
   call structure_friction_3d
#endif
gotm's avatar
gotm committed
561
   if (ufirst) then
562 563
      call uu_momentum_3d(n,bdy3d)
      call vv_momentum_3d(n,bdy3d)
gotm's avatar
gotm committed
564 565
      ufirst=.false.
   else
566 567
      call vv_momentum_3d(n,bdy3d)
      call uu_momentum_3d(n,bdy3d)
gotm's avatar
gotm committed
568 569
      ufirst=.true.
   end if
570

gotm's avatar
gotm committed
571
   if (kmax .gt. 1) then
Knut's avatar
Knut committed
572

Knut's avatar
Knut committed
573 574 575 576 577 578 579 580 581
!    KK-TODO: In realistic simulations (gotm) we need SS
!             in any case, therefore it is done here by default.
!             In the future one might check whether a very seldom case
!             is present, where it can be skipped.
!             We need SS: 1) #if (!defined(CONSTANT_VISCOSITY) && !defined(PARABOLIC_VISCOSITY))
!                         2) adpative coordinates
!                         3) if(do_numerical_analyses_3d) [physical dissipation analyses]
      call shear_frequency()

gotm's avatar
gotm committed
582
      call ww_momentum_3d()
583 584 585 586

      call uv_advect_3d()
      call uv_diffusion_3d()  ! Must be called after uv_advect_3d

587
#ifndef NO_BOTTFRIC
gotm's avatar
gotm committed
588
      call stresses_3d()
589 590
#endif
#ifndef CONSTANT_VISCOSITY
gotm's avatar
gotm committed
591
      call gotm()
592
      if (advect_turbulence) call tke_eps_advect_3d()
gotm's avatar
gotm committed
593
#endif
Knut's avatar
Knut committed
594

595
   end if
Knut's avatar
Knut committed
596

kbk's avatar
kbk committed
597
#ifndef NO_BAROCLINIC
kbk's avatar
kbk committed
598
   if(runtype .eq. 4) then        ! prognostic T and S
gotm's avatar
gotm committed
599 600
      if (calc_temp) call do_temperature(n)
      if (calc_salt) call do_salinity(n)
kbk's avatar
kbk committed
601 602 603 604

!     The following is a bit clumpsy and should be changed when do_bdy_3d()
!     operates on individual fields and not as is the case now - on both
!     T and S.
bjb's avatar
bjb committed
605
      call tic(TIM_INTEGR3D)
kbk's avatar
kbk committed
606 607
      if (bdy3d) call do_bdy_3d(0,T)
      if (calc_temp) then
608
         call tic(TIM_TEMPH)
609
         call update_3d_halo(T,T,az,imin,jmin,imax,jmax,kmax,D_TAG)
kbk's avatar
kbk committed
610
         call wait_halo(D_TAG)
611
         call toc(TIM_TEMPH)
kbk's avatar
kbk committed
612 613 614
         call mirror_bdy_3d(T,D_TAG)
      end if
      if (calc_salt) then
615
         call tic(TIM_SALTH)
616
         call update_3d_halo(S,S,az,imin,jmin,imax,jmax,kmax,D_TAG)
kbk's avatar
kbk committed
617
         call wait_halo(D_TAG)
618
         call toc(TIM_SALTH)
kbk's avatar
kbk committed
619 620
         call mirror_bdy_3d(S,D_TAG)
      end if
bjb's avatar
bjb committed
621
      call toc(TIM_INTEGR3D)
kbk's avatar
kbk committed
622

gotm's avatar
gotm committed
623
      call do_eqstate()
624 625 626 627 628 629 630 631
!     KK-TODO: In realistic simulations (baroclinic+gotm) we need NN
!              in any case, therefore it is done here by default.
!              In the future one might check whether a very seldom case
!              is present, where it can be skipped.
!              We need NN (runtype .ge. 3):
!                          1) #if (!defined(CONSTANT_VISCOSITY) && !defined(PARABOLIC_VISCOSITY))
!                          2) adaptive coordinates
      call buoyancy_frequency()
gotm's avatar
gotm committed
632
   end if
kbk's avatar
kbk committed
633
#endif
gotm's avatar
gotm committed
634

635
#ifndef NO_BAROTROPIC
kbk's avatar
kbk committed
636
   if (kmax .gt. 1) then
gotm's avatar
gotm committed
637
#ifndef NO_BOTTFRIC
kbk's avatar
kbk committed
638
      call slow_bottom_friction()
gotm's avatar
gotm committed
639
#endif
640 641

      call tic(TIM_INTEGR3D)
Knut's avatar
Knut committed
642
      call uv_advect(Uint,Vint,Dn,Dveln,Dun,Dvn)
643 644 645
      call uv_diffusion(0,Uint,Vint,Dn,Dun,Dvn) ! Has to be called after uv_advect.
      call toc(TIM_INTEGR3D)

gotm's avatar
gotm committed
646
   end if
gotm's avatar
gotm committed
647

gotm's avatar
gotm committed
648
   call slow_terms()
649
#endif
bjb's avatar
bjb committed
650
   call tic(TIM_INTEGR3D)
gotm's avatar
gotm committed
651
   call stop_macro()
bjb's avatar
bjb committed
652
   call toc(TIM_INTEGR3D)
gotm's avatar
gotm committed
653 654 655 656 657 658 659 660 661 662 663 664

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

!-----------------------------------------------------------------------
!BOP
!
hb's avatar
hb committed
665
! !IROUTINE: clean_3d - cleanup after 3D run \label{sec-clean-3d}
gotm's avatar
gotm committed
666 667 668 669 670 671 672 673 674 675 676 677
!
! !INTERFACE:
   subroutine clean_3d()
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
hb's avatar
hb committed
678 679
! Here, a call to the routine {\tt clean\_variables\_3d} which howewer
! does not do anything yet.
gotm's avatar
gotm committed
680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708
!
! !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)         !
!-----------------------------------------------------------------------