m3d.F90 20.2 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
kbk's avatar
kbk committed
29
   use m2d, only: Am
kbk's avatar
kbk committed
30
   use variables_2d, only: D,z,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
gotm's avatar
gotm committed
40
   use advection_3d, only: init_advection_3d
kbk's avatar
kbk committed
41
   use bdy_3d, only: init_bdy_3d, do_bdy_3d
42
   use bdy_3d, only: bdy3d_tmrlx, bdy3d_tmrlx_ucut, bdy3d_tmrlx_max, bdy3d_tmrlx_min
43
!  Necessary to use halo_zones because update_3d_halos() have been moved out
kbk's avatar
kbk committed
44 45
!  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
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
kbk's avatar
kbk committed
57
   REALTYPE                            :: ip_fac=_ONE_
58
   integer                             :: vel_check=0
59
   REALTYPE                            :: min_vel=-4*_ONE_,max_vel=4*_ONE_
gotm's avatar
gotm committed
60
!
kbk's avatar
kbk committed
61 62 63
! !REVISION HISTORY:
!  Original author(s): Karsten Bolding & Hans Burchard
!
gotm's avatar
gotm committed
64
! !LOCAL VARIABLES:
65
   integer         :: vel_hor_adv=1,vel_ver_adv=1,vel_adv_split=0
kbk's avatar
kbk committed
66 67 68
#ifdef NO_BAROCLINIC
   integer         :: ip_method
#endif
kbk's avatar
kbk committed
69
   integer         :: ip_ramp=-1
gotm's avatar
gotm committed
70 71 72 73 74 75 76 77
!EOP
!-----------------------------------------------------------------------

   contains

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

135 136 137 138 139 140 141 142 143
   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
144 145 146
! Allocates memory for the public data members - if not static
   call init_variables_3d(runtype)

147 148 149 150 151
!  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
152 153 154
#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.'
155
      call getm_error("init_3d()","GETM needs recompilation")
gotm's avatar
gotm committed
156
#endif
157 158 159 160 161 162
   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
163
   end if
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 211 212 213
   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

bjb's avatar
bjb committed
214
!  Sanity checks for bdy 3d
215 216 217 218 219 220 221 222 223 224 225 226 227 228
   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
229 230
   end if

kbk's avatar
kbk committed
231 232
   LEVEL2 'ip_ramp=',ip_ramp

233 234
   LEVEL2 'vel_check=',vel_check
   if (vel_check .ne. 0) then
kbk's avatar
kbk committed
235 236 237
      LEVEL3 'doing sanity checks on velocities'
      LEVEL3 'min_vel=',min_vel
      LEVEL3 'max_vel=',max_vel
238 239 240 241 242 243 244 245
      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
246
   dt = M*timestep
247 248 249 250 251

#ifdef CONSTANT_VISCOSITY
   num=avmback
   nuh=avhback
#else
252 253
   num=1.d-15
   nuh=1.d-15
254
#endif
gotm's avatar
gotm committed
255 256

!  Needed for interpolation of temperature and salinity
kbk's avatar
kbk committed
257
   if (.not. hotstart) then
gotm's avatar
gotm committed
258
      call start_macro()
259
      call coordinates(hotstart)
260
      call hcc_check()
kbk's avatar
kbk committed
261
   end if
gotm's avatar
gotm committed
262

kbk's avatar
kbk committed
263
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
264
   if (runtype .eq. 3 .or. runtype .eq. 4) then
gotm's avatar
gotm committed
265
      T = _ZERO_ ; S = _ZERO_ ; rho = _ZERO_
gotm's avatar
gotm committed
266 267
      if(calc_temp) call init_temperature(1)
      if(calc_salt) call init_salinity(1)
kbk's avatar
kbk committed
268
   end if
269
#endif
kbk's avatar
kbk committed
270

271 272
    call init_advection_3d(2)

273
#ifndef NO_BAROCLINIC
kbk's avatar
kbk committed
274
    if (runtype .eq. 3 .or. runtype .eq. 4) then
275
      call init_eqstate()
gotm's avatar
gotm committed
276 277
#ifndef PECS
      call do_eqstate()
278
      call ss_nn()
gotm's avatar
gotm committed
279
#endif
280 281 282

      if (.not. openbdy) bdy3d=.false.
      if (bdy3d) call init_bdy_3d()
283 284
      if (runtype .ge. 3) call init_internal_pressure()
      if (runtype .eq. 3) call do_internal_pressure()
gotm's avatar
gotm committed
285
   end if
286
#endif
287

288
   if (vert_cord .eq. _ADAPTIVE_COORDS_) call preadapt_coordinates(preadapt)
kb's avatar
kb committed
289

gotm's avatar
gotm committed
290 291 292 293 294 295 296 297
#ifdef DEBUG
   write(debug,*) 'Leaving init_3d()'
   write(debug,*)
#endif
   return
   end subroutine init_3d
!EOC

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 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372
!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: postinit_3d - re-initialise some 3D after hotstart read.
!
! !INTERFACE:
   subroutine postinit_3d(runtype)
! !USES:
   use domain, only: imin,imax,jmin,jmax, az,au,av
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
   integer, intent(in)                 :: runtype
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
!  This routine provides possibility to reset/initialize 3D variables to 
!  ensure that velocities are correctly set on land cells after read 
!  of a hotstart file.
!
! !LOCAL VARIABLES:
   integer                   :: i,j
!EOP
!-------------------------------------------------------------------------
!BOC
#ifdef DEBUG
   integer, save :: Ncall = 0
   Ncall = Ncall+1
   write(debug,*) 'postinit_3d() # ',Ncall
#endif

   LEVEL1 'postinit_3d'

! Go over mask and make sure that velocities are zero on land.
   do j=jmin-HALO,jmax+HALO
      do i=imin-HALO,imax+HALO
         if(az(i,j) .eq. 0) then
            ho(i,j,:)  = _ZERO_
            hn(i,j,:)  = _ZERO_
            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
   end do
   do j=jmin-HALO,jmax+HALO
      do i=imin-HALO,imax+HALO
         if (au(i,j) .eq. 0) then
            uu(i,j,:)  = _ZERO_
            huo(i,j,:) = _ZERO_
            hun(i,j,:) = _ZERO_
         end if
      end do
   end do
   do j=jmin-HALO,jmax+HALO
      do i=imin-HALO,imax+HALO
         if (av(i,j) .eq. 0) then
            vv(i,j,:)  = _ZERO_
            hvo(i,j,:) = _ZERO_
            hvn(i,j,:) = _ZERO_
         end if
      end do
   end do

   return
   end subroutine postinit_3d
!EOC

gotm's avatar
gotm committed
373 374 375
!-----------------------------------------------------------------------
!BOP
!
hb's avatar
hb committed
376 377
! !IROUTINE: integrate_3d - calls to do 3D model integration
! \label{sec-integrate-3d}
gotm's avatar
gotm committed
378 379 380
!
! !INTERFACE:
   subroutine integrate_3d(runtype,n)
bjb's avatar
bjb committed
381
   use getm_timers, only: tic, toc, TIM_INTEGR3D
382 383 384
#ifndef NO_BAROCLINIC
   use getm_timers, only: TIM_TEMPH, TIM_SALTH
#endif
gotm's avatar
gotm committed
385 386 387
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
kbk's avatar
kbk committed
388
   integer, intent(in)                 :: runtype,n
gotm's avatar
gotm committed
389 390 391 392 393 394
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
hb's avatar
hb committed
395
! This is a wrapper routine to call all 3D related subroutines.
396 397 398
! The call position for the {\tt coordinates} routine depends on
! the compiler option
! {\tt MUDFLAT}: If it is defined, then the
hb's avatar
hb committed
399 400
! 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
401
! Lagrangian component which are currently under development are supported.
hb's avatar
hb committed
402 403 404 405 406
! Both, drying and flooding and
! Lagrangian coordinates does not go together yet.
! The call sequence is as follows:
!
! \vspace{0.5cm}
407
!
hb's avatar
hb committed
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 435 436 437 438 439 440 441 442
! \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} \\
443
! {\tt do\_spm}                & suspended matter equation & see page
hb's avatar
hb committed
444 445 446 447 448 449 450 451 452 453 454 455 456
! \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}
457
!
hb's avatar
hb committed
458 459 460 461 462 463
! \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
464 465
!
! !LOCAL VARIABLES:
kbk's avatar
kbk committed
466
  logical, save              :: ufirst=.true.
gotm's avatar
gotm committed
467 468 469 470 471 472 473 474 475 476
!
!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
477 478
#ifndef NO_BAROCLINIC
#endif
479
#ifdef MUDFLAT
480
   call coordinates(.false.)
481
#endif
gotm's avatar
gotm committed
482
#ifndef NO_BOTTFRIC
kbk's avatar
kbk committed
483
   if (kmax .gt. 1) then
gotm's avatar
gotm committed
484 485 486
      call bottom_friction_3d()
   end if
#endif
bjb's avatar
bjb committed
487
   call tic(TIM_INTEGR3D)
gotm's avatar
gotm committed
488
   SS = _ZERO_
bjb's avatar
bjb committed
489
   call toc(TIM_INTEGR3D)
kbk's avatar
kbk committed
490
#ifndef NO_BAROCLINIC
bjb's avatar
bjb committed
491
   call tic(TIM_INTEGR3D)
gotm's avatar
gotm committed
492
   NN = _ZERO_
bjb's avatar
bjb committed
493
   call toc(TIM_INTEGR3D)
494
   if (runtype .eq. 4) call do_internal_pressure()
kbk's avatar
kbk committed
495
#endif
bjb's avatar
bjb committed
496
   call tic(TIM_INTEGR3D)
497 498
   huo=hun
   hvo=hvn
kbk's avatar
kbk committed
499
   ip_fac=_ONE_
frv-bjb's avatar
frv-bjb committed
500
   if (ip_ramp .gt. 0) ip_fac=min( _ONE_ , n*_ONE_/ip_ramp)
bjb's avatar
bjb committed
501
   call toc(TIM_INTEGR3D)
502 503 504
#ifdef STRUCTURE_FRICTION
   call structure_friction_3d
#endif
gotm's avatar
gotm committed
505
   if (ufirst) then
506 507
      call uu_momentum_3d(n,bdy3d)
      call vv_momentum_3d(n,bdy3d)
gotm's avatar
gotm committed
508 509
      ufirst=.false.
   else
510 511
      call vv_momentum_3d(n,bdy3d)
      call uu_momentum_3d(n,bdy3d)
gotm's avatar
gotm committed
512 513
      ufirst=.true.
   end if
514

515
#ifndef MUDFLAT
516 517 518 519
   if (kmax .gt. 1) then
      if (vert_cord .eq. _ADAPTIVE_COORDS_) call ss_nn()
   end if
   call coordinates(.false.)
520
#endif
gotm's avatar
gotm committed
521 522 523 524
   if (kmax .gt. 1) then
      call ww_momentum_3d()
   end if
#ifndef NO_ADVECT
kbk's avatar
kbk committed
525
   if (kmax .gt. 1) then
526
      call uv_advect_3d(vel_hor_adv,vel_ver_adv,vel_adv_split)
gotm's avatar
gotm committed
527 528 529 530 531 532 533
      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
534
   if (kmax .gt. 1) then
535
#ifndef NO_BOTTFRIC
gotm's avatar
gotm committed
536
      call stresses_3d()
537 538
#endif
#ifndef CONSTANT_VISCOSITY
gotm's avatar
gotm committed
539
#ifndef PARABOLIC_VISCOSITY
540
      if (vert_cord .ne. _ADAPTIVE_COORDS_) call ss_nn()
gotm's avatar
gotm committed
541 542
#endif
      call gotm()
hb's avatar
hb committed
543 544 545
#ifdef TURB_ADV
      call tke_eps_advect_3d(vel_hor_adv,vel_ver_adv,vel_adv_split)
#endif
gotm's avatar
gotm committed
546
#endif
547
   end if
kbk's avatar
kbk committed
548
#ifndef NO_BAROCLINIC
kbk's avatar
kbk committed
549
   if(runtype .eq. 4) then        ! prognostic T and S
gotm's avatar
gotm committed
550 551
      if (calc_temp) call do_temperature(n)
      if (calc_salt) call do_salinity(n)
kbk's avatar
kbk committed
552 553 554 555

!     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
556
      call tic(TIM_INTEGR3D)
kbk's avatar
kbk committed
557 558
      if (bdy3d) call do_bdy_3d(0,T)
      if (calc_temp) then
559
         call tic(TIM_TEMPH)
560
         call update_3d_halo(T,T,az,imin,jmin,imax,jmax,kmax,D_TAG)
kbk's avatar
kbk committed
561
         call wait_halo(D_TAG)
562
         call toc(TIM_TEMPH)
kbk's avatar
kbk committed
563 564 565
         call mirror_bdy_3d(T,D_TAG)
      end if
      if (calc_salt) then
566
         call tic(TIM_SALTH)
567
         call update_3d_halo(S,S,az,imin,jmin,imax,jmax,kmax,D_TAG)
kbk's avatar
kbk committed
568
         call wait_halo(D_TAG)
569
         call toc(TIM_SALTH)
kbk's avatar
kbk committed
570 571
         call mirror_bdy_3d(S,D_TAG)
      end if
bjb's avatar
bjb committed
572
      call toc(TIM_INTEGR3D)
kbk's avatar
kbk committed
573

gotm's avatar
gotm committed
574 575 576 577
#ifndef PECS
      call do_eqstate()
#endif
   end if
kbk's avatar
kbk committed
578
#endif
gotm's avatar
gotm committed
579

580
   call tic(TIM_INTEGR3D)
kbk's avatar
kbk committed
581
   UEx=_ZERO_ ; VEx=_ZERO_
582
   call toc(TIM_INTEGR3D)
583
#ifndef NO_BAROTROPIC
kbk's avatar
kbk committed
584
   if (kmax .gt. 1) then
gotm's avatar
gotm committed
585
#ifndef NO_BOTTFRIC
kbk's avatar
kbk committed
586
      call slow_bottom_friction()
gotm's avatar
gotm committed
587 588 589
#endif
#ifndef NO_ADVECT
#ifndef UV_ADV_DIRECT
kbk's avatar
kbk committed
590 591 592 593
      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
594 595 596
#endif
#endif
   end if
gotm's avatar
gotm committed
597

gotm's avatar
gotm committed
598
   call slow_terms()
599
#endif
bjb's avatar
bjb committed
600
   call tic(TIM_INTEGR3D)
gotm's avatar
gotm committed
601
   call stop_macro()
bjb's avatar
bjb committed
602
   call toc(TIM_INTEGR3D)
gotm's avatar
gotm committed
603 604 605 606 607 608 609 610 611 612 613 614

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

!-----------------------------------------------------------------------
!BOP
!
hb's avatar
hb committed
615
! !IROUTINE: clean_3d - cleanup after 3D run \label{sec-clean-3d}
gotm's avatar
gotm committed
616 617 618 619 620 621 622 623 624 625 626 627
!
! !INTERFACE:
   subroutine clean_3d()
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
hb's avatar
hb committed
628 629
! Here, a call to the routine {\tt clean\_variables\_3d} which howewer
! does not do anything yet.
gotm's avatar
gotm committed
630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658
!
! !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)         !
!-----------------------------------------------------------------------