m3d.F90 20.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 30
   use m2d, only: uv_advect,uv_diffusion
   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
kbk's avatar
kbk committed
56 57 58
   logical                             :: calc_temp=.true.
   logical                             :: calc_salt=.true.
   logical                             :: bdy3d=.false.
59
   integer                             :: bdyfmt_3d,bdy3d_ramp
kbk's avatar
kbk committed
60
   character(len=PATH_MAX)             :: bdyfile_3d
kbk's avatar
kbk committed
61
   REALTYPE                            :: ip_fac=_ONE_
62
   integer                             :: vel_check=0
63
   REALTYPE                            :: min_vel=-4*_ONE_,max_vel=4*_ONE_
gotm's avatar
gotm committed
64
!
kbk's avatar
kbk committed
65 66 67
! !REVISION HISTORY:
!  Original author(s): Karsten Bolding & Hans Burchard
!
gotm's avatar
gotm committed
68
! !LOCAL VARIABLES:
69
   logical         :: advect_turbulence=.false.
kbk's avatar
kbk committed
70 71 72
#ifdef NO_BAROCLINIC
   integer         :: ip_method
#endif
kbk's avatar
kbk committed
73
   integer         :: ip_ramp=-1
gotm's avatar
gotm committed
74 75 76 77 78 79 80 81
!EOP
!-----------------------------------------------------------------------

   contains

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

142 143 144 145 146 147 148 149 150
   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
151 152
! Allocates memory for the public data members - if not static
   call init_variables_3d(runtype)
153
   call init_advection_3d()
gotm's avatar
gotm committed
154

155
!  Sanity checks for advection specifications
156 157 158 159 160 161 162 163 164 165 166 167 168
   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
169
   end if
170 171
#endif
   call print_adv_settings_3d(vel3d_adv_split,vel3d_adv_hor,vel3d_adv_ver,_ZERO_)
172

bjb's avatar
bjb committed
173
!  Sanity checks for bdy 3d
Knut's avatar
Knut committed
174
   if (.not.openbdy .or. runtype.eq.2) bdy3d=.false.
175 176 177 178 179 180 181 182 183 184 185 186 187 188
   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
189 190
   end if

kbk's avatar
kbk committed
191 192
   LEVEL2 'ip_ramp=',ip_ramp

193 194
   LEVEL2 'vel_check=',vel_check
   if (vel_check .ne. 0) then
kbk's avatar
kbk committed
195 196 197
      LEVEL3 'doing sanity checks on velocities'
      LEVEL3 'min_vel=',min_vel
      LEVEL3 'max_vel=',max_vel
198 199 200 201 202 203 204 205
      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
206
   dt = M*timestep
207 208 209 210 211

#ifdef CONSTANT_VISCOSITY
   num=avmback
   nuh=avhback
#else
212 213
   num=1.d-15
   nuh=1.d-15
214 215 216 217 218 219 220 221 222 223 224 225

#ifdef TURB_ADV
   if (.not. advect_turbulence) then
      LEVEL2 "reenabled advection of TKE and eps due to"
      LEVEL2 "obsolete TURB_ADV macro. Note that this"
      LEVEL2 "behaviour will be removed in the future!"
      advect_turbulence = .true.
   end if
#endif

   LEVEL2 "advect_turbulence = ",advect_turbulence

226
#endif
gotm's avatar
gotm committed
227 228

!  Needed for interpolation of temperature and salinity
kbk's avatar
kbk committed
229
   if (.not. hotstart) then
230
      ssen = z
gotm's avatar
gotm committed
231
      call start_macro()
232
      call coordinates(hotstart)
233
      call hcc_check()
kbk's avatar
kbk committed
234
   end if
gotm's avatar
gotm committed
235

236 237 238 239
   if (runtype .eq. 2) then
      calc_temp = .false.
      calc_salt = .false.
   else
kbk's avatar
kbk committed
240
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
241
      T = _ZERO_ ; S = _ZERO_ ; rho = _ZERO_
242 243
      if(calc_temp) call init_temperature()
      if(calc_salt) call init_salinity()
244
#endif
245
   end if
kbk's avatar
kbk committed
246

247
#ifndef NO_BAROCLINIC
kbk's avatar
kbk committed
248
    if (runtype .eq. 3 .or. runtype .eq. 4) then
249
      call init_eqstate()
gotm's avatar
gotm committed
250 251
#ifndef PECS
      call do_eqstate()
252
      call ss_nn()
gotm's avatar
gotm committed
253
#endif
254 255

      if (bdy3d) call init_bdy_3d()
256 257
      if (runtype .ge. 3) call init_internal_pressure()
      if (runtype .eq. 3) call do_internal_pressure()
gotm's avatar
gotm committed
258
   end if
259
#endif
260

261
   if (vert_cord .eq. _ADAPTIVE_COORDS_) call preadapt_coordinates(preadapt)
kb's avatar
kb committed
262

gotm's avatar
gotm committed
263 264 265 266 267 268 269 270
#ifdef DEBUG
   write(debug,*) 'Leaving init_3d()'
   write(debug,*)
#endif
   return
   end subroutine init_3d
!EOC

271 272 273 274 275 276
!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: postinit_3d - re-initialise some 3D after hotstart read.
!
! !INTERFACE:
277
   subroutine postinit_3d(runtype,timestep,hotstart)
278 279 280 281 282 283
! !USES:
   use domain, only: imin,imax,jmin,jmax, az,au,av
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
   integer, intent(in)                 :: runtype
284 285
   REALTYPE, intent(in)                :: timestep
   logical, intent(in)                 :: hotstart
286 287 288 289 290 291
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
292 293
!  This routine provides possibility to reset/initialize 3D variables to
!  ensure that velocities are correctly set on land cells after read
294 295 296
!  of a hotstart file.
!
! !LOCAL VARIABLES:
297
   integer                   :: i,j,rc
298 299 300 301 302 303 304 305 306 307 308
!EOP
!-------------------------------------------------------------------------
!BOC
#ifdef DEBUG
   integer, save :: Ncall = 0
   Ncall = Ncall+1
   write(debug,*) 'postinit_3d() # ',Ncall
#endif

   LEVEL1 'postinit_3d'

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
   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

350 351
! Hotstart fix - see postinit_2d
   if (hotstart) then
352

353 354 355 356 357 358
      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
359
      end do
360 361 362 363 364 365
      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
366
      end do
367 368 369 370 371 372 373 374 375 376 377 378 379
!     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
380
      end do
381 382 383

      call coordinates(hotstart)

384
   end if
385 386 387 388 389

   return
   end subroutine postinit_3d
!EOC

gotm's avatar
gotm committed
390 391 392
!-----------------------------------------------------------------------
!BOP
!
hb's avatar
hb committed
393 394
! !IROUTINE: integrate_3d - calls to do 3D model integration
! \label{sec-integrate-3d}
gotm's avatar
gotm committed
395 396 397
!
! !INTERFACE:
   subroutine integrate_3d(runtype,n)
bjb's avatar
bjb committed
398
   use getm_timers, only: tic, toc, TIM_INTEGR3D
399 400 401
#ifndef NO_BAROCLINIC
   use getm_timers, only: TIM_TEMPH, TIM_SALTH
#endif
gotm's avatar
gotm committed
402 403 404
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
kbk's avatar
kbk committed
405
   integer, intent(in)                 :: runtype,n
gotm's avatar
gotm committed
406 407 408 409 410 411
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
hb's avatar
hb committed
412
! This is a wrapper routine to call all 3D related subroutines.
413 414 415
! 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
416 417
! 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
418
! Lagrangian component which are currently under development are supported.
hb's avatar
hb committed
419 420 421 422 423
! Both, drying and flooding and
! Lagrangian coordinates does not go together yet.
! The call sequence is as follows:
!
! \vspace{0.5cm}
424
!
hb's avatar
hb committed
425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459
! \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} \\
460
! {\tt do\_spm}                & suspended matter equation & see page
hb's avatar
hb committed
461 462 463 464 465 466 467 468 469 470 471 472 473
! \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}
474
!
hb's avatar
hb committed
475 476 477 478 479 480
! \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
481 482
!
! !LOCAL VARIABLES:
kbk's avatar
kbk committed
483
  logical, save              :: ufirst=.true.
gotm's avatar
gotm committed
484 485 486 487 488 489 490 491 492 493
!
!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
494 495
#ifndef NO_BAROCLINIC
#endif
496
#ifdef MUDFLAT
497
   call coordinates(.false.)
498
#endif
gotm's avatar
gotm committed
499
#ifndef NO_BOTTFRIC
kbk's avatar
kbk committed
500
   if (kmax .gt. 1) then
gotm's avatar
gotm committed
501 502 503
      call bottom_friction_3d()
   end if
#endif
bjb's avatar
bjb committed
504
   call tic(TIM_INTEGR3D)
gotm's avatar
gotm committed
505
   SS = _ZERO_
bjb's avatar
bjb committed
506
   call toc(TIM_INTEGR3D)
kbk's avatar
kbk committed
507
#ifndef NO_BAROCLINIC
bjb's avatar
bjb committed
508
   call tic(TIM_INTEGR3D)
gotm's avatar
gotm committed
509
   NN = _ZERO_
bjb's avatar
bjb committed
510
   call toc(TIM_INTEGR3D)
511
   if (runtype .eq. 4) call do_internal_pressure()
kbk's avatar
kbk committed
512
#endif
bjb's avatar
bjb committed
513
   call tic(TIM_INTEGR3D)
514 515
   huo=hun
   hvo=hvn
kbk's avatar
kbk committed
516
   ip_fac=_ONE_
frv-bjb's avatar
frv-bjb committed
517
   if (ip_ramp .gt. 0) ip_fac=min( _ONE_ , n*_ONE_/ip_ramp)
bjb's avatar
bjb committed
518
   call toc(TIM_INTEGR3D)
519 520 521
#ifdef STRUCTURE_FRICTION
   call structure_friction_3d
#endif
gotm's avatar
gotm committed
522
   if (ufirst) then
523 524
      call uu_momentum_3d(n,bdy3d)
      call vv_momentum_3d(n,bdy3d)
gotm's avatar
gotm committed
525 526
      ufirst=.false.
   else
527 528
      call vv_momentum_3d(n,bdy3d)
      call uu_momentum_3d(n,bdy3d)
gotm's avatar
gotm committed
529 530
      ufirst=.true.
   end if
531

Knut's avatar
Knut committed
532 533
!  HERE WE HAVE TO UPDATE SS !!!

534
#ifndef MUDFLAT
535 536 537 538
   if (kmax .gt. 1) then
      if (vert_cord .eq. _ADAPTIVE_COORDS_) call ss_nn()
   end if
   call coordinates(.false.)
539
#endif
Knut's avatar
Knut committed
540

gotm's avatar
gotm committed
541
   if (kmax .gt. 1) then
Knut's avatar
Knut committed
542

gotm's avatar
gotm committed
543
      call ww_momentum_3d()
544 545 546 547

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

548
#ifndef NO_BOTTFRIC
gotm's avatar
gotm committed
549
      call stresses_3d()
550 551
#endif
#ifndef CONSTANT_VISCOSITY
gotm's avatar
gotm committed
552
#ifndef PARABOLIC_VISCOSITY
553
      if (vert_cord .ne. _ADAPTIVE_COORDS_) call ss_nn()
gotm's avatar
gotm committed
554 555
#endif
      call gotm()
556
      if (advect_turbulence) call tke_eps_advect_3d()
gotm's avatar
gotm committed
557
#endif
Knut's avatar
Knut committed
558

559
   end if
Knut's avatar
Knut committed
560

kbk's avatar
kbk committed
561
#ifndef NO_BAROCLINIC
kbk's avatar
kbk committed
562
   if(runtype .eq. 4) then        ! prognostic T and S
gotm's avatar
gotm committed
563 564
      if (calc_temp) call do_temperature(n)
      if (calc_salt) call do_salinity(n)
kbk's avatar
kbk committed
565 566 567 568

!     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
569
      call tic(TIM_INTEGR3D)
kbk's avatar
kbk committed
570 571
      if (bdy3d) call do_bdy_3d(0,T)
      if (calc_temp) then
572
         call tic(TIM_TEMPH)
573
         call update_3d_halo(T,T,az,imin,jmin,imax,jmax,kmax,D_TAG)
kbk's avatar
kbk committed
574
         call wait_halo(D_TAG)
575
         call toc(TIM_TEMPH)
kbk's avatar
kbk committed
576 577 578
         call mirror_bdy_3d(T,D_TAG)
      end if
      if (calc_salt) then
579
         call tic(TIM_SALTH)
580
         call update_3d_halo(S,S,az,imin,jmin,imax,jmax,kmax,D_TAG)
kbk's avatar
kbk committed
581
         call wait_halo(D_TAG)
582
         call toc(TIM_SALTH)
kbk's avatar
kbk committed
583 584
         call mirror_bdy_3d(S,D_TAG)
      end if
bjb's avatar
bjb committed
585
      call toc(TIM_INTEGR3D)
kbk's avatar
kbk committed
586

gotm's avatar
gotm committed
587 588 589
#ifndef PECS
      call do_eqstate()
#endif
Knut's avatar
Knut committed
590
!     HERE WE HAVE TO UPDATE NN !!!
gotm's avatar
gotm committed
591
   end if
kbk's avatar
kbk committed
592
#endif
gotm's avatar
gotm committed
593

594
#ifndef NO_BAROTROPIC
kbk's avatar
kbk committed
595
   if (kmax .gt. 1) then
gotm's avatar
gotm committed
596
#ifndef NO_BOTTFRIC
kbk's avatar
kbk committed
597
      call slow_bottom_friction()
gotm's avatar
gotm committed
598
#endif
599 600 601 602 603 604

      call tic(TIM_INTEGR3D)
      call uv_advect(Uint,Vint,Dun,Dvn)
      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
605
   end if
gotm's avatar
gotm committed
606

gotm's avatar
gotm committed
607
   call slow_terms()
608
#endif
bjb's avatar
bjb committed
609
   call tic(TIM_INTEGR3D)
gotm's avatar
gotm committed
610
   call stop_macro()
bjb's avatar
bjb committed
611
   call toc(TIM_INTEGR3D)
gotm's avatar
gotm committed
612 613 614 615 616 617 618 619 620 621 622 623

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

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