m3d.F90 16.6 KB
Newer Older
1
!$Id: m3d.F90,v 1.33 2006-12-15 09:57:50 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,az
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 47 48
!  Necessary to use halo_zones because update_3d_halos() have been moved out 
!  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
49

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

   contains

!-----------------------------------------------------------------------
!BOP
!
hb's avatar
hb committed
79
! !IROUTINE: init_3d - initialise 3D related stuff \label{sec-init-3d}
gotm's avatar
gotm committed
80 81 82 83 84 85
!
! !INTERFACE:
   subroutine init_3d(runtype,timestep,hotstart)
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
kbk's avatar
kbk committed
86 87 88
   integer, intent(in)                 :: runtype
   REALTYPE, intent(in)                :: timestep
   logical, intent(in)                 :: hotstart
gotm's avatar
gotm committed
89 90 91 92 93 94
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
hb's avatar
hb committed
95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112
!  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
113 114
!
! !LOCAL VARIABLES:
115
   integer         :: rc
kbk's avatar
kbk committed
116 117 118
   NAMELIST /m3d/ &
             M,cnpar,cord_relax,                        &
             bdy3d,bdyfmt_3d,bdyramp_3d,bdyfile_3d,     &
119
             vel_hor_adv,vel_ver_adv,vel_adv_split,     &
kbk's avatar
kbk committed
120
             calc_temp,calc_salt,                       &
121 122
             avmback,avhback,ip_method,                 &
             vel_check,min_vel,max_vel
gotm's avatar
gotm committed
123 124 125 126 127 128 129 130 131 132 133 134 135 136 137
!
!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)

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

150 151 152 153 154
!  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
155 156 157
#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.'
158
      call getm_error("init_3d()","GETM needs recompilation")
gotm's avatar
gotm committed
159
#endif
160 161 162 163 164 165
   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
166
   end if
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 214 215 216
   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

217 218 219 220 221 222 223 224 225 226 227
   LEVEL2 'vel_check=',vel_check
   if (vel_check .ne. 0) then
      LEVEL3 'doing sanity checks on velocities every ',vel_check,' timesteps'
      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
228
   dt = M*timestep
229 230 231 232 233 234 235 236

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

!  Needed for interpolation of temperature and salinity
kbk's avatar
kbk committed
239
   if (.not. hotstart) then
gotm's avatar
gotm committed
240
      call start_macro()
kbk's avatar
kbk committed
241
      call coordinates(vert_cord,cord_relax,maxdepth)
242
      call hcc_check()
kbk's avatar
kbk committed
243
   end if
gotm's avatar
gotm committed
244

kbk's avatar
kbk committed
245
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
246
   if (runtype .eq. 3 .or. runtype .eq. 4) then
gotm's avatar
gotm committed
247
      T = _ZERO_ ; S = _ZERO_ ; rho = _ZERO_
gotm's avatar
gotm committed
248 249
      if(calc_temp) call init_temperature(1)
      if(calc_salt) call init_salinity(1)
kbk's avatar
kbk committed
250
   end if
251
#endif
kbk's avatar
kbk committed
252

253
#ifndef NO_BAROCLINIC
kbk's avatar
kbk committed
254
    if (runtype .eq. 3 .or. runtype .eq. 4) then
255
      call init_eqstate()
gotm's avatar
gotm committed
256 257 258
#ifndef PECS
      call do_eqstate()
#endif
259 260
      if (runtype .ge. 3) call init_internal_pressure()
      if (runtype .eq. 3) call do_internal_pressure()
gotm's avatar
gotm committed
261 262 263 264
      if (runtype .eq. 4) then
         call init_advection_3d(2)
      end if
   end if
kbk's avatar
kbk committed
265
#endif
266 267 268 269

#ifdef GETM_BIO
!KBK      if(bio_calc) call init_getm_bio()
      call init_getm_bio()
kbk's avatar
kbk committed
270
#endif
gotm's avatar
gotm committed
271 272 273 274 275 276 277 278 279 280 281 282 283 284

   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
285 286
! !IROUTINE: integrate_3d - calls to do 3D model integration
! \label{sec-integrate-3d}
gotm's avatar
gotm committed
287 288 289 290 291 292
!
! !INTERFACE:
   subroutine integrate_3d(runtype,n)
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
kbk's avatar
kbk committed
293
   integer, intent(in)                 :: runtype,n
gotm's avatar
gotm committed
294 295 296 297 298 299
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
hb's avatar
hb committed
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
! 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
369 370
!
! !LOCAL VARIABLES:
kbk's avatar
kbk committed
371
  logical, save              :: ufirst=.true.
gotm's avatar
gotm committed
372 373 374 375 376 377 378 379 380 381
!
!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
382 383
#ifndef NO_BAROCLINIC
#endif
384 385 386
#ifdef MUDFLAT
   call coordinates(vert_cord,cord_relax,maxdepth)
#endif
gotm's avatar
gotm committed
387
#ifndef NO_BOTTFRIC
kbk's avatar
kbk committed
388
   if (kmax .gt. 1) then
gotm's avatar
gotm committed
389 390 391 392
      call bottom_friction_3d()
   end if
#endif
   SS = _ZERO_
kbk's avatar
kbk committed
393
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
394
   NN = _ZERO_
395
   if (runtype .eq. 4) call do_internal_pressure()
kbk's avatar
kbk committed
396
#endif
397 398
   huo=hun
   hvo=hvn
gotm's avatar
gotm committed
399
   if (ufirst) then
400 401
      call uu_momentum_3d(n,bdy3d)
      call vv_momentum_3d(n,bdy3d)
gotm's avatar
gotm committed
402 403
      ufirst=.false.
   else
404 405
      call vv_momentum_3d(n,bdy3d)
      call uu_momentum_3d(n,bdy3d)
gotm's avatar
gotm committed
406 407
      ufirst=.true.
   end if
408
#ifndef MUDFLAT
409
   call coordinates(vert_cord,cord_relax,maxdepth)
410
#endif
gotm's avatar
gotm committed
411 412 413 414
   if (kmax .gt. 1) then
      call ww_momentum_3d()
   end if
#ifndef NO_ADVECT
kbk's avatar
kbk committed
415
   if (kmax .gt. 1) then
416
      call uv_advect_3d(vel_hor_adv,vel_ver_adv,vel_adv_split)
gotm's avatar
gotm committed
417 418 419 420 421 422 423
      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
424 425

   if (kmax .gt. 1) then
426
#ifndef NO_BOTTFRIC
gotm's avatar
gotm committed
427
      call stresses_3d()
428 429
#endif
#ifndef CONSTANT_VISCOSITY
gotm's avatar
gotm committed
430 431 432 433 434
#ifndef PARABOLIC_VISCOSITY
      call ss_nn()
#endif
      call gotm()
#endif
435
   end if
kbk's avatar
kbk committed
436
#ifndef NO_BAROCLINIC
kbk's avatar
kbk committed
437
   if(runtype .eq. 4) then        ! prognostic T and S
gotm's avatar
gotm committed
438 439
      if (calc_temp) call do_temperature(n)
      if (calc_salt) call do_salinity(n)
kbk's avatar
kbk committed
440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456

!     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.
#ifndef NO_BAROCLINIC
      if (bdy3d) call do_bdy_3d(0,T)
      if (calc_temp) then
         call update_3d_halo(T,T,az,iimin,jjmin,iimax,jjmax,kmax,D_TAG)
         call wait_halo(D_TAG)
         call mirror_bdy_3d(T,D_TAG)
      end if
      if (calc_salt) then
         call update_3d_halo(S,S,az,iimin,jjmin,iimax,jjmax,kmax,D_TAG)
         call wait_halo(D_TAG)
         call mirror_bdy_3d(S,D_TAG)
      end if
#endif
kbk's avatar
kbk committed
457
   end if
458
#endif
kbk's avatar
kbk committed
459

460 461
#ifndef NO_BAROCLINIC
   if(runtype .eq. 4) then
gotm's avatar
gotm committed
462 463 464 465
#ifndef PECS
      call do_eqstate()
#endif
   end if
kbk's avatar
kbk committed
466
#endif
gotm's avatar
gotm committed
467

468 469
#ifdef GETM_BIO
   if (bio_calc) call do_getm_bio(dt)
kbk's avatar
kbk committed
470 471
#endif

kbk's avatar
kbk committed
472
   UEx=_ZERO_ ; VEx=_ZERO_
kbk's avatar
kbk committed
473
   if (kmax .gt. 1) then
gotm's avatar
gotm committed
474
#ifndef NO_BOTTFRIC
kbk's avatar
kbk committed
475
      call slow_bottom_friction()
gotm's avatar
gotm committed
476 477 478
#endif
#ifndef NO_ADVECT
#ifndef UV_ADV_DIRECT
kbk's avatar
kbk committed
479 480 481 482
      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
483 484 485
#endif
#endif
   end if
gotm's avatar
gotm committed
486

gotm's avatar
gotm committed
487 488 489 490 491 492 493 494 495 496 497 498 499 500
   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
501
! !IROUTINE: clean_3d - cleanup after 3D run \label{sec-clean-3d}
gotm's avatar
gotm committed
502 503 504 505 506 507 508 509 510 511 512 513
!
! !INTERFACE:
   subroutine clean_3d()
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
hb's avatar
hb committed
514 515
! Here, a call to the routine {\tt clean\_variables\_3d} which howewer
! does not do anything yet.
gotm's avatar
gotm committed
516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544
!
! !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)         !
!-----------------------------------------------------------------------