New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
zdfmxl.F90 in NEMO/branches/2021/dev_r14393_HPC-03_Mele_Comm_Cleanup/src/OCE/ZDF – NEMO

source: NEMO/branches/2021/dev_r14393_HPC-03_Mele_Comm_Cleanup/src/OCE/ZDF/zdfmxl.F90 @ 14601

Last change on this file since 14601 was 14601, checked in by francesca, 3 years ago

[comm_cleanup: ZDF files] - ticket #2607

  • Property svn:keywords set to Id
File size: 7.3 KB
RevLine 
[3]1MODULE zdfmxl
2   !!======================================================================
3   !!                       ***  MODULE  zdfmxl  ***
4   !! Ocean physics: mixed layer depth
5   !!======================================================================
[1559]6   !! History :  1.0  ! 2003-08  (G. Madec)  original code
[1585]7   !!            3.2  ! 2009-07  (S. Masson, G. Madec)  IOM + merge of DO-loop
[4990]8   !!            3.7  ! 2012-03  (G. Madec)  make public the density criteria for trdmxl
9   !!             -   ! 2014-02  (F. Roquet)  mixed layer depth calculated using N2 instead of rhop
[3]10   !!----------------------------------------------------------------------
[1585]11   !!   zdf_mxl      : Compute the turbocline and mixed layer depths.
[3]12   !!----------------------------------------------------------------------
[9019]13   USE oce            ! ocean dynamics and tracers variables
[12377]14   USE isf_oce        ! ice shelf
[9019]15   USE dom_oce        ! ocean space and time domain variables
16   USE trc_oce  , ONLY: l_offline         ! ocean space and time domain variables
17   USE zdf_oce        ! ocean vertical physics
[9104]18   !
[9019]19   USE in_out_manager ! I/O manager
20   USE prtctl         ! Print control
21   USE phycst         ! physical constants
22   USE iom            ! I/O library
23   USE lib_mpp        ! MPP library
[3]24
25   IMPLICIT NONE
26   PRIVATE
27
[9019]28   PUBLIC   zdf_mxl   ! called by zdfphy.F90
[3]29
[9019]30   INTEGER , PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   nmln    !: number of level in the mixed layer (used by LDF, ZDF, TRD, TOP)
31   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hmld    !: mixing layer depth (turbocline)      [m]   (used by TOP)
32   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hmlp    !: mixed layer depth  (rho=rho0+zdcrit) [m]   (used by LDF)
33   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hmlpt   !: depth of the last T-point inside the mixed layer [m] (used by LDF)
[3]34
[4990]35   REAL(wp), PUBLIC ::   rho_c = 0.01_wp    !: density criterion for mixed layer depth
[10351]36   REAL(wp), PUBLIC ::   avt_c = 5.e-4_wp   ! Kz criterion for the turbocline depth
[4990]37
[12377]38   !! * Substitutions
39#  include "do_loop_substitute.h90"
[13237]40#  include "domzgr_substitute.h90"
[3]41   !!----------------------------------------------------------------------
[9598]42   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
[1152]43   !! $Id$
[10068]44   !! Software governed by the CeCILL license (see ./LICENSE)
[3]45   !!----------------------------------------------------------------------
46CONTAINS
47
[2715]48   INTEGER FUNCTION zdf_mxl_alloc()
49      !!----------------------------------------------------------------------
50      !!               ***  FUNCTION zdf_mxl_alloc  ***
51      !!----------------------------------------------------------------------
[2787]52      zdf_mxl_alloc = 0      ! set to zero if no array to be allocated
[2758]53      IF( .NOT. ALLOCATED( nmln ) ) THEN
54         ALLOCATE( nmln(jpi,jpj), hmld(jpi,jpj), hmlp(jpi,jpj), hmlpt(jpi,jpj), STAT= zdf_mxl_alloc )
55         !
[10425]56         CALL mpp_sum ( 'zdfmxl', zdf_mxl_alloc )
57         IF( zdf_mxl_alloc /= 0 )   CALL ctl_stop( 'STOP', 'zdf_mxl_alloc: failed to allocate arrays.' )
[2758]58         !
59      ENDIF
[2715]60   END FUNCTION zdf_mxl_alloc
61
62
[12377]63   SUBROUTINE zdf_mxl( kt, Kmm )
[3]64      !!----------------------------------------------------------------------
65      !!                  ***  ROUTINE zdfmxl  ***
66      !!                   
[1585]67      !! ** Purpose :   Compute the turbocline depth and the mixed layer depth
68      !!              with density criteria.
[3]69      !!
[1577]70      !! ** Method  :   The mixed layer depth is the shallowest W depth with
71      !!      the density of the corresponding T point (just bellow) bellow a
[4245]72      !!      given value defined locally as rho(10m) + rho_c
[1585]73      !!               The turbocline depth is the depth at which the vertical
74      !!      eddy diffusivity coefficient (resulting from the vertical physics
75      !!      alone, not the isopycnal part, see trazdf.F) fall below a given
[4990]76      !!      value defined locally (avt_c here taken equal to 5 cm/s2 by default)
[3]77      !!
[1585]78      !! ** Action  :   nmln, hmld, hmlp, hmlpt
[1559]79      !!----------------------------------------------------------------------
[2715]80      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
[12377]81      INTEGER, INTENT(in) ::   Kmm  ! ocean time level index
[4990]82      !
[6352]83      INTEGER  ::   ji, jj, jk      ! dummy loop indices
84      INTEGER  ::   iikn, iiki, ikt ! local integer
85      REAL(wp) ::   zN2_c           ! local scalar
[9019]86      INTEGER, DIMENSION(jpi,jpj) ::   imld   ! 2D workspace
[3]87      !!----------------------------------------------------------------------
[3294]88      !
[3]89      IF( kt == nit000 ) THEN
90         IF(lwp) WRITE(numout,*)
91         IF(lwp) WRITE(numout,*) 'zdf_mxl : mixed layer depth'
92         IF(lwp) WRITE(numout,*) '~~~~~~~ '
[2715]93         !                             ! allocate zdfmxl arrays
94         IF( zdf_mxl_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_mxl : unable to allocate arrays' )
[3]95      ENDIF
[9104]96      !
[1559]97      ! w-level of the mixing and mixed layers
[13497]98      nmln(:,:)  = nlb10                  ! Initialization to the number of w ocean point
99      hmlp(:,:)  = 0._wp                  ! here hmlp used as a dummy variable, integrating vertically N^2
100      zN2_c = grav * rho_c * r1_rho0      ! convert density criteria into N^2 criteria
[14601]101      ! [comm_cleanup] ! DO_3D( 1, 1, 1, 1, nlb10, jpkm1 )   ! Mixed layer level: w-level
102      DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, nlb10, jpkm1 )   ! Mixed layer level: w-level
[12377]103         ikt = mbkt(ji,jj)
[13237]104         hmlp(ji,jj) =   &
105            & hmlp(ji,jj) + MAX( rn2b(ji,jj,jk) , 0._wp ) * e3w(ji,jj,jk,Kmm)
[12377]106         IF( hmlp(ji,jj) < zN2_c )   nmln(ji,jj) = MIN( jk , ikt ) + 1   ! Mixed layer level
107      END_3D
[4990]108      !
[6140]109      ! w-level of the turbocline and mixing layer (iom_use)
[13497]110      imld(:,:) = mbkt(:,:) + 1                ! Initialization to the number of w ocean point
[14601]111      ! [comm_cleanup] ! DO_3DS( 1, 1, 1, 1, jpkm1, nlb10, -1 )   ! from the bottom to nlb10
112      DO_3DS( nn_hls, nn_hls, nn_hls, nn_hls, jpkm1, nlb10, -1 )   ! from the bottom to nlb10
[12377]113         IF( avt (ji,jj,jk) < avt_c * wmask(ji,jj,jk) )   imld(ji,jj) = jk      ! Turbocline
114      END_3D
[1559]115      ! depth of the mixing and mixed layers
[14601]116      ! [comm_cleanup] ! DO_2D( 1, 1, 1, 1 )
117      DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
[12377]118         iiki = imld(ji,jj)
119         iikn = nmln(ji,jj)
120         hmld (ji,jj) = gdepw(ji,jj,iiki  ,Kmm) * ssmask(ji,jj)    ! Turbocline depth
121         hmlp (ji,jj) = gdepw(ji,jj,iikn  ,Kmm) * ssmask(ji,jj)    ! Mixed layer depth
122         hmlpt(ji,jj) = gdept(ji,jj,iikn-1,Kmm) * ssmask(ji,jj)    ! depth of the last T-point inside the mixed layer
123      END_2D
[7646]124      !
125      IF( .NOT.l_offline ) THEN
126         IF( iom_use("mldr10_1") ) THEN
127            IF( ln_isfcav ) THEN  ;  CALL iom_put( "mldr10_1", hmlp - risfdep)   ! mixed layer thickness
128            ELSE                  ;  CALL iom_put( "mldr10_1", hmlp )            ! mixed layer depth
[6352]129            END IF
[6140]130         END IF
[7646]131         IF( iom_use("mldkz5") ) THEN
132            IF( ln_isfcav ) THEN  ;  CALL iom_put( "mldkz5"  , hmld - risfdep )   ! turbocline thickness
133            ELSE                  ;  CALL iom_put( "mldkz5"  , hmld )             ! turbocline depth
[6352]134            END IF
[7646]135         ENDIF
[2758]136      ENDIF
[9104]137      !
[12377]138      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hmlp, clinfo2=' hmlp : ' )
[1559]139      !
[3]140   END SUBROUTINE zdf_mxl
141
142   !!======================================================================
143END MODULE zdfmxl
Note: See TracBrowser for help on using the repository browser.