m3d.F90 15.7 KB
Newer Older
1
!$Id: m3d.F90,v 1.23 2004-08-09 07:48:07 kbk Exp $
gotm's avatar
gotm committed
2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21
#include "cppdefs.h"
!-----------------------------------------------------------------------
!BOP
!
! !MODULE: m3d - 3D model component
!
! !INTERFACE:
   module m3d
!
! !DESCRIPTION:
!  This modules contains declarations for all variables related to 3D
!  hydrodynamical calculations. Information about the calculation domain
!  is included from the \emph{domain.F90} module.
!  The module contains public subroutines for initialisation, integration
!  and clean up of the 3D model component.
!  The actual calculation routines are called in integrate\_3d and is linked
!  in from the library lib3d.a.
!
! !USES:
   use parameters, only: avmmol
kbk's avatar
kbk committed
22
   use domain, only: maxdepth,vert_cord
kbk's avatar
kbk committed
23
   use m2d, only: Am
kbk's avatar
kbk committed
24
   use variables_2d, only: D,z,UEx,VEx
kbk's avatar
kbk committed
25
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
26 27
   use temperature,only: init_temperature, do_temperature
   use salinity,   only: init_salinity, do_salinity
kbk's avatar
kbk committed
28
   use eqstate,    only: init_eqstate, do_eqstate
29 30
   use internal_pressure, only: init_internal_pressure, do_internal_pressure
   use internal_pressure, only: ip_method
kbk's avatar
kbk committed
31
#endif
kbk's avatar
kbk committed
32
#ifdef SPM
kbk's avatar
kbk committed
33
   use suspended_matter, only: init_spm, do_spm
kbk's avatar
kbk committed
34
#endif
35
   use variables_3d
gotm's avatar
gotm committed
36
   use advection_3d, only: init_advection_3d
kbk's avatar
kbk committed
37
   use bdy_3d, only: init_bdy_3d, do_bdy_3d
kbk's avatar
kbk committed
38

gotm's avatar
gotm committed
39 40 41
   IMPLICIT NONE
!
! !PUBLIC DATA MEMBERS:
kbk's avatar
kbk committed
42
   integer                             :: M=1
43
   REALTYPE                            :: cord_relax=_ZERO_
kbk's avatar
kbk committed
44 45 46
   logical                             :: calc_temp=.true.
   logical                             :: calc_salt=.true.
   logical                             :: calc_spm=.false.
kbk's avatar
kbk committed
47
   logical                             :: hotstart_spm=.false.
kbk's avatar
kbk committed
48 49 50
   logical                             :: bdy3d=.false.
   integer                             :: bdyfmt_3d,bdyramp_3d
   character(len=PATH_MAX)             :: bdyfile_3d
gotm's avatar
gotm committed
51 52 53 54 55
!
! !REVISION HISTORY:
!  Original author(s): Karsten Bolding & Hans Burchard
!
!  $Log: m3d.F90,v $
56 57 58 59
!  Revision 1.23  2004-08-09 07:48:07  kbk
!  checking for negative avmback and avhback
!
!  Revision 1.22  2004/08/06 15:14:35  hb
60 61 62
!  num and nuh now properly initialised and no gotm call for CONSTANT_VISCOSITY
!
!  Revision 1.21  2004/07/29 19:48:44  hb
hb's avatar
hb committed
63 64 65
!  num, nuh now initialised with _ZERO_
!
!  Revision 1.20  2004/07/28 14:58:18  hb
66 67 68
!  Changing subroutine calling order via MUDFLAT
!
!  Revision 1.19  2004/06/15 08:25:57  kbk
kbk's avatar
kbk committed
69 70 71
!  added supoort for spm - Ruiz
!
!  Revision 1.18  2004/05/04 09:23:51  kbk
72 73 74
!  hydrostatic consistency criteria stored in .3d.nc file
!
!  Revision 1.17  2004/04/23 09:03:59  kbk
kbk's avatar
kbk committed
75
!  reverted to pre-adaptive grid version
hb's avatar
hb committed
76 77
!
!  Revision 1.15  2004/04/20 16:49:37  hb
78 79 80
!  call to coordinates moved for better consistency (see JMB)
!
!  Revision 1.14  2004/04/06 12:42:50  kbk
81 82 83
!  internal pressure calculations now uses wrapper
!
!  Revision 1.13  2004/01/08 10:23:20  kbk
84 85 86
!  NN not needed for barotropic runs, NO_SUSP_MATTER works
!
!  Revision 1.12  2004/01/06 15:04:00  kbk
87 88 89
!  FCT advection + split of advection_3d.F90 + extra adv. input checks
!
!  Revision 1.11  2004/01/05 13:23:27  kbk
kbk's avatar
kbk committed
90 91 92
!  Poor Man's Z-coordinates
!
!  Revision 1.10  2004/01/02 13:54:24  kbk
93 94 95
!  read equation of state info from namelist - Ruiz
!
!  Revision 1.9  2003/12/16 17:02:44  kbk
kbk's avatar
kbk committed
96 97 98
!  removed TABS - 0. -> _ZERO_
!
!  Revision 1.8  2003/12/16 15:58:54  kbk
99 100 101
!  back ground viscosity and diffusivity (manuel)
!
!  Revision 1.7  2003/09/12 16:23:38  kbk
102 103 104
!  fixed order of use statement - now compiles on with IFC 7.1
!
!  Revision 1.6  2003/08/28 15:13:57  kbk
kbk's avatar
kbk committed
105 106 107
!  explicit set UEx and VEx to 0
!
!  Revision 1.5  2003/04/23 12:16:34  kbk
kbk's avatar
kbk committed
108 109 110
!  cleaned code + TABS to spaces
!
!  Revision 1.4  2003/04/07 16:28:34  kbk
kbk's avatar
kbk committed
111
!  parallel support
gotm's avatar
gotm committed
112 113 114
!
!  Revision 1.1.1.1  2002/05/02 14:00:51  gotm
!  recovering after CVS crash
gotm's avatar
gotm committed
115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210
!
!  Revision 1.30  2001/10/26 09:13:24  bbh
!  Only call slow_diffusion() if Am > 0.
!
!  Revision 1.29  2001/10/23 12:38:58  bbh
!  Forgot cord_relax in one call to coordinates()
!
!  Revision 1.28  2001/10/23 07:53:10  bbh
!  spm -> suspended_matter
!
!  Revision 1.27  2001/10/23 07:05:11  bbh
!  Parabolic viscosity - -DPARABOLIC_VISCOSITY
!
!  Revision 1.26  2001/10/22 11:16:09  bbh
!  Added support for particulate suspended matter - no IO yet
!
!  Revision 1.25  2001/10/22 09:26:41  bbh
!  Added cord_relax
!
!  Revision 1.24  2001/10/22 08:48:30  bbh
!  Am moved from paramters.F90 to m2d.F90
!
!  Revision 1.23  2001/10/12 11:39:20  bbh
!  TVD moved out of ??_momentum_3d.F90 and into uv_advect_3d.F90
!
!  Revision 1.22  2001/10/12 08:49:56  bbh
!  Support for horizontal diffusion
!
!  Revision 1.21  2001/09/19 13:07:00  bbh
!  Moved advection related 3D fields to global allocation
!
!  Revision 1.20  2001/09/03 20:14:03  bbh
!  Fixed reading of vel_{hor,vel,strang
!
!  Revision 1.19  2001/09/03 20:04:21  bbh
!  Allow individual advection settings for momentum, salinity and temperature
!
!  Revision 1.18  2001/09/03 13:02:52  bbh
!  Small changes due to NOMADS
!
!  Revision 1.17  2001/08/31 15:39:49  bbh
!  ifdef for CONSTANCE added
!
!  Revision 1.16  2001/08/30 08:56:26  bbh
!  Preparing for 3D boundary conditions
!
!  Revision 1.15  2001/08/29 12:08:05  bbh
!  temp_adv and salt_adv needs to be public
!
!  Revision 1.14  2001/08/29 11:21:46  bbh
!  namelists read in salinity and temperature + initialisation
!
!  Revision 1.13  2001/08/27 11:51:45  bbh
!  TVD-advection for momentum added, some bugs removed
!
!  Revision 1.12  2001/08/01 08:31:22  bbh
!  CURVILINEAR now implemented
!
!  Revision 1.11  2001/07/26 13:46:06  bbh
!  Included info for salinity and temperature in namelist
!
!  Revision 1.10  2001/06/22 08:19:10  bbh
!  Compiler options such as USE_MASK and OLD_DRY deleted.
!  Open and passive boundary for z created.
!  Various inconsistencies removed.
!  wait_halo added.
!  Checked loop boundaries
!
!  Revision 1.9  2001/05/25 19:36:36  bbh
!  Added call to eqstate
!
!  Revision 1.8  2001/05/21 13:07:19  bbh
!  dt and cnpar is in variables_3d.F90
!
!  Revision 1.7  2001/05/20 07:51:40  bbh
!  Internal pressure included
!
!  Revision 1.6  2001/05/18 08:17:49  bbh
!  Tidying up initialization of salinity and temperature
!
!  Revision 1.5  2001/05/15 11:47:57  bbh
!  Added update_3d_halo for initial S and T
!
!  Revision 1.4  2001/05/10 11:30:16  bbh
!  Added further support for baroclinicity
!
!  Revision 1.3  2001/05/03 20:12:31  bbh
!  Use of variables_3d
!
!  Revision 1.2  2001/04/24 08:39:21  bbh
!  Included runtype as argument to integrate_3d
!
!  Revision 1.1.1.1  2001/04/17 08:43:08  bbh
!  initial import into CVS
!
! !LOCAL VARIABLES:
211
   integer         :: vel_hor_adv=1,vel_ver_adv=1,vel_adv_split=0
gotm's avatar
gotm committed
212 213 214 215 216 217 218 219 220 221 222 223 224 225 226
!EOP
!-----------------------------------------------------------------------

   contains

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: init_3d - initialise 3D relatedstuff.
!
! !INTERFACE:
   subroutine init_3d(runtype,timestep,hotstart)
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
kbk's avatar
kbk committed
227 228 229
   integer, intent(in)                 :: runtype
   REALTYPE, intent(in)                :: timestep
   logical, intent(in)                 :: hotstart
gotm's avatar
gotm committed
230 231 232 233 234 235 236 237 238 239 240 241 242
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
!  Allocates memiory for 3D related fields.
!
! !REVISION HISTORY:
!
!  22Apr99   Karsten Bolding & Hans Burchard  Initial code.
!
! !LOCAL VARIABLES:
243
   integer         :: rc
kbk's avatar
kbk committed
244 245 246
   NAMELIST /m3d/ &
             M,cnpar,cord_relax,                        &
             bdy3d,bdyfmt_3d,bdyramp_3d,bdyfile_3d,     &
247
             vel_hor_adv,vel_ver_adv,vel_adv_split,     &
248
             calc_temp,calc_salt,calc_spm,              &
249
             avmback,avhback,ip_method
gotm's avatar
gotm committed
250 251 252 253 254 255 256 257 258 259 260 261 262 263 264
!
!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)

265 266 267 268 269 270 271 272 273
   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

kbk's avatar
kbk committed
274 275 276 277
#ifndef SPM
   if(calc_spm) stop 'To use SPM you have to recompile with -DSPM'
#endif

gotm's avatar
gotm committed
278 279 280
! Allocates memory for the public data members - if not static
   call init_variables_3d(runtype)

281 282 283 284 285
!  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
286 287 288
#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.'
289
      call getm_error("init_3d()","GETM needs recompilation")
gotm's avatar
gotm committed
290
#endif
291 292 293 294 295 296
   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
297
   end if
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
   LEVEL2 "horizontal: ",trim(adv_schemes(vel_hor_adv))," of momentum"
   LEVEL2 "vertical:   ",trim(adv_schemes(vel_ver_adv))," of momentum"

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

gotm's avatar
gotm committed
348
   dt = M*timestep
349 350 351 352 353 354 355 356

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

!  Needed for interpolation of temperature and salinity
kbk's avatar
kbk committed
359
   if (.not. hotstart) then
gotm's avatar
gotm committed
360
      call start_macro()
kbk's avatar
kbk committed
361
      call coordinates(vert_cord,cord_relax,maxdepth)
362
      call hcc_check()
kbk's avatar
kbk committed
363
   end if
gotm's avatar
gotm committed
364

kbk's avatar
kbk committed
365
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
366
   if (runtype .eq. 3 .or. runtype .eq. 4) then
gotm's avatar
gotm committed
367
      T = _ZERO_ ; S = _ZERO_ ; rho = _ZERO_
gotm's avatar
gotm committed
368 369
      if(calc_temp) call init_temperature(1)
      if(calc_salt) call init_salinity(1)
kbk's avatar
kbk committed
370
   end if
371
#endif
kbk's avatar
kbk committed
372

373
#ifndef NO_BAROCLINIC
kbk's avatar
kbk committed
374
    if (runtype .eq. 3 .or. runtype .eq. 4) then
375
      call init_eqstate()
gotm's avatar
gotm committed
376 377 378
#ifndef PECS
      call do_eqstate()
#endif
379 380
      if (runtype .ge. 3) call init_internal_pressure()
      if (runtype .eq. 3) call do_internal_pressure()
gotm's avatar
gotm committed
381 382 383 384
      if (runtype .eq. 4) then
         call init_advection_3d(2)
      end if
   end if
kbk's avatar
kbk committed
385
#endif
kbk's avatar
kbk committed
386 387 388 389
#ifdef SPM
      if(calc_spm)  call init_spm(hotstart_spm,runtype)
      if(runtype.ne.4) call init_advection_3d(2)
#endif
gotm's avatar
gotm committed
390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410

   if (bdy3d) call init_bdy_3d()

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

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: integrate_3d - sequence of calls to do 3D model integration
!
! !INTERFACE:
   subroutine integrate_3d(runtype,n)
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
kbk's avatar
kbk committed
411
   integer, intent(in)                 :: runtype,n
gotm's avatar
gotm committed
412 413 414 415 416 417 418
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
! A wrapper to call all 3D related subroutines in one subroutine.
419
! if #MUDFLAT is defined, then the sequence of velocity equations
420
! and coordinate construction is made such that drying and flooding
421 422
! is stable. If #MUDFLAT is not defined, then adaptive grids with
! Lagrangian component are supported. Both, drying and flooding and
423
! Lagrangian coordinates does not go together.
gotm's avatar
gotm committed
424 425 426 427 428
!
! !REVISION HISTORY:
!  See log for module
!
! !LOCAL VARIABLES:
kbk's avatar
kbk committed
429
  logical, save              :: ufirst=.true.
gotm's avatar
gotm committed
430 431 432 433 434 435 436 437 438 439
!
!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
440
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
441
   if (bdy3d) call do_bdy_3d(0,T)
kbk's avatar
kbk committed
442
#endif
443 444 445
#ifdef MUDFLAT
   call coordinates(vert_cord,cord_relax,maxdepth)
#endif
gotm's avatar
gotm committed
446
#ifndef NO_BOTTFRIC
kbk's avatar
kbk committed
447
   if (kmax .gt. 1) then
gotm's avatar
gotm committed
448 449 450 451
      call bottom_friction_3d()
   end if
#endif
   SS = _ZERO_
kbk's avatar
kbk committed
452
#ifndef NO_BAROCLINIC
gotm's avatar
gotm committed
453
   NN = _ZERO_
454
   if (runtype .eq. 4) call do_internal_pressure()
kbk's avatar
kbk committed
455
#endif
456 457
   huo=hun
   hvo=hvn
gotm's avatar
gotm committed
458 459 460 461 462 463 464 465 466
   if (ufirst) then
      call uu_momentum_3d(bdy3d)
      call vv_momentum_3d(bdy3d)
      ufirst=.false.
   else
      call vv_momentum_3d(bdy3d)
      call uu_momentum_3d(bdy3d)
      ufirst=.true.
   end if
467
#ifndef MUDFLAT
468
   call coordinates(vert_cord,cord_relax,maxdepth)
469
#endif
gotm's avatar
gotm committed
470 471 472 473
   if (kmax .gt. 1) then
      call ww_momentum_3d()
   end if
#ifndef NO_ADVECT
kbk's avatar
kbk committed
474
   if (kmax .gt. 1) then
475
      call uv_advect_3d(vel_hor_adv,vel_ver_adv,vel_adv_split)
gotm's avatar
gotm committed
476 477 478 479 480 481 482
      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
483 484

   if (kmax .gt. 1) then
485
#ifndef NO_BOTTFRIC
gotm's avatar
gotm committed
486
      call stresses_3d()
487 488
#endif
#ifndef CONSTANT_VISCOSITY
gotm's avatar
gotm committed
489 490 491 492 493
#ifndef PARABOLIC_VISCOSITY
      call ss_nn()
#endif
      call gotm()
#endif
494
   end if
kbk's avatar
kbk committed
495
#ifndef NO_BAROCLINIC
kbk's avatar
kbk committed
496
   if(runtype .eq. 4) then        ! prognostic T and S
gotm's avatar
gotm committed
497 498
      if (calc_temp) call do_temperature(n)
      if (calc_salt) call do_salinity(n)
kbk's avatar
kbk committed
499
   end if
500
#endif
kbk's avatar
kbk committed
501

502 503
#ifndef NO_BAROCLINIC
   if(runtype .eq. 4) then
gotm's avatar
gotm committed
504 505 506 507
#ifndef PECS
      call do_eqstate()
#endif
   end if
kbk's avatar
kbk committed
508
#endif
gotm's avatar
gotm committed
509

kbk's avatar
kbk committed
510 511 512 513
#ifdef SPM
      if (calc_spm) call do_spm()
#endif

kbk's avatar
kbk committed
514
   UEx=_ZERO_ ; VEx=_ZERO_
kbk's avatar
kbk committed
515
   if (kmax .gt. 1) then
gotm's avatar
gotm committed
516
#ifndef NO_BOTTFRIC
kbk's avatar
kbk committed
517
      call slow_bottom_friction()
gotm's avatar
gotm committed
518 519 520
#endif
#ifndef NO_ADVECT
#ifndef UV_ADV_DIRECT
kbk's avatar
kbk committed
521 522 523 524
      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
525 526 527
#endif
#endif
   end if
gotm's avatar
gotm committed
528

gotm's avatar
gotm committed
529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588
   call slow_terms()
   call stop_macro()

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

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: clean_3d - cleanup after 3D run.
!
! !INTERFACE:
   subroutine clean_3d()
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
!
! !INPUT/OUTPUT PARAMETERS:
!
! !OUTPUT PARAMETERS:
!
! !DESCRIPTION:
!  This routine cleans up after a 3D integration. Close open files etc.
!
! !REVISION HISTORY:
!  22Nov Author name Initial code
!
! !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)         !
!-----------------------------------------------------------------------