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 branches/UKMO/r8395_mix-lyr_diag/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/UKMO/r8395_mix-lyr_diag/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfmxl.F90 @ 11290

Last change on this file since 11290 was 11290, checked in by jcastill, 5 years ago

Remove svn keywords

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   !!----------------------------------------------------------------------
13   USE oce             ! ocean dynamics and tracers variables
14   USE dom_oce         ! ocean space and time domain variables
[7646]15   USE trc_oce, ONLY: l_offline         ! ocean space and time domain variables
[1585]16   USE zdf_oce         ! ocean vertical physics
[3]17   USE in_out_manager  ! I/O manager
[258]18   USE prtctl          ! Print control
[4990]19   USE phycst          ! physical constants
[2715]20   USE iom             ! I/O library
21   USE lib_mpp         ! MPP library
[3294]22   USE wrk_nemo        ! work arrays
23   USE timing          ! Timing
[3]24
25   IMPLICIT NONE
26   PRIVATE
27
[2715]28   PUBLIC   zdf_mxl       ! called by step.F90
[3]29
[2715]30   INTEGER , PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   nmln    !: number of level in the mixed layer (used by TOP)
31   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hmld    !: mixing layer depth (turbocline)      [m]
32   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hmlp    !: mixed layer depth  (rho=rho0+zdcrit) [m]
[6352]33   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hmlpt   !: depth of the last T-point inside the mixed layer [m]
[3]34
[4990]35   REAL(wp), PUBLIC ::   rho_c = 0.01_wp    !: density criterion for mixed layer depth
36   REAL(wp)         ::   avt_c = 5.e-4_wp   ! Kz criterion for the turbocline depth
37
[3]38   !!----------------------------------------------------------------------
[2715]39   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
[11290]40   !! $Id$
[2715]41   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[3]42   !!----------------------------------------------------------------------
43CONTAINS
44
[2715]45   INTEGER FUNCTION zdf_mxl_alloc()
46      !!----------------------------------------------------------------------
47      !!               ***  FUNCTION zdf_mxl_alloc  ***
48      !!----------------------------------------------------------------------
[2787]49      zdf_mxl_alloc = 0      ! set to zero if no array to be allocated
[2758]50      IF( .NOT. ALLOCATED( nmln ) ) THEN
51         ALLOCATE( nmln(jpi,jpj), hmld(jpi,jpj), hmlp(jpi,jpj), hmlpt(jpi,jpj), STAT= zdf_mxl_alloc )
52         !
53         IF( lk_mpp             )   CALL mpp_sum ( zdf_mxl_alloc )
54         IF( zdf_mxl_alloc /= 0 )   CALL ctl_warn('zdf_mxl_alloc: failed to allocate arrays.')
55         !
56      ENDIF
[2715]57   END FUNCTION zdf_mxl_alloc
58
59
[3]60   SUBROUTINE zdf_mxl( kt )
61      !!----------------------------------------------------------------------
62      !!                  ***  ROUTINE zdfmxl  ***
63      !!                   
[1585]64      !! ** Purpose :   Compute the turbocline depth and the mixed layer depth
65      !!              with density criteria.
[3]66      !!
[1577]67      !! ** Method  :   The mixed layer depth is the shallowest W depth with
68      !!      the density of the corresponding T point (just bellow) bellow a
[4245]69      !!      given value defined locally as rho(10m) + rho_c
[1585]70      !!               The turbocline depth is the depth at which the vertical
71      !!      eddy diffusivity coefficient (resulting from the vertical physics
72      !!      alone, not the isopycnal part, see trazdf.F) fall below a given
[4990]73      !!      value defined locally (avt_c here taken equal to 5 cm/s2 by default)
[3]74      !!
[1585]75      !! ** Action  :   nmln, hmld, hmlp, hmlpt
[1559]76      !!----------------------------------------------------------------------
[2715]77      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
[4990]78      !
[6352]79      INTEGER  ::   ji, jj, jk      ! dummy loop indices
80      INTEGER  ::   iikn, iiki, ikt ! local integer
81      REAL(wp) ::   zN2_c           ! local scalar
[4990]82      INTEGER, POINTER, DIMENSION(:,:) ::   imld   ! 2D workspace
[3]83      !!----------------------------------------------------------------------
[3294]84      !
85      IF( nn_timing == 1 )  CALL timing_start('zdf_mxl')
86      !
87      CALL wrk_alloc( jpi,jpj, imld )
[3]88
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
96
[1559]97      ! w-level of the mixing and mixed layers
[7753]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_rau0   ! convert density criteria into N^2 criteria
[4990]101      DO jk = nlb10, jpkm1
102         DO jj = 1, jpj                ! Mixed layer level: w-level
103            DO ji = 1, jpi
104               ikt = mbkt(ji,jj)
[6140]105               hmlp(ji,jj) = hmlp(ji,jj) + MAX( rn2b(ji,jj,jk) , 0._wp ) * e3w_n(ji,jj,jk)
[4990]106               IF( hmlp(ji,jj) < zN2_c )   nmln(ji,jj) = MIN( jk , ikt ) + 1   ! Mixed layer level
107            END DO
108         END DO
109      END DO
110      !
[6140]111      ! w-level of the turbocline and mixing layer (iom_use)
[7753]112      imld(:,:) = mbkt(:,:) + 1        ! Initialization to the number of w ocean point
[4990]113      DO jk = jpkm1, nlb10, -1         ! from the bottom to nlb10
[3]114         DO jj = 1, jpj
115            DO ji = 1, jpi
[6140]116               IF( avt (ji,jj,jk) < avt_c * wmask(ji,jj,jk) )   imld(ji,jj) = jk      ! Turbocline
[3]117            END DO
118         END DO
119      END DO
[1559]120      ! depth of the mixing and mixed layers
[3]121      DO jj = 1, jpj
122         DO ji = 1, jpi
[1585]123            iiki = imld(ji,jj)
[1577]124            iikn = nmln(ji,jj)
[6140]125            hmld (ji,jj) = gdepw_n(ji,jj,iiki  ) * ssmask(ji,jj)    ! Turbocline depth
126            hmlp (ji,jj) = gdepw_n(ji,jj,iikn  ) * ssmask(ji,jj)    ! Mixed layer depth
127            hmlpt(ji,jj) = gdept_n(ji,jj,iikn-1) * ssmask(ji,jj)    ! depth of the last T-point inside the mixed layer
[3]128         END DO
129      END DO
[7646]130      !
131      IF( .NOT.l_offline ) THEN
132         IF( iom_use("mldr10_1") ) THEN
133            IF( ln_isfcav ) THEN  ;  CALL iom_put( "mldr10_1", hmlp - risfdep)   ! mixed layer thickness
134            ELSE                  ;  CALL iom_put( "mldr10_1", hmlp )            ! mixed layer depth
[6352]135            END IF
[6140]136         END IF
[7646]137         IF( iom_use("mldkz5") ) THEN
138            IF( ln_isfcav ) THEN  ;  CALL iom_put( "mldkz5"  , hmld - risfdep )   ! turbocline thickness
139            ELSE                  ;  CALL iom_put( "mldkz5"  , hmld )             ! turbocline depth
[6352]140            END IF
[7646]141         ENDIF
[2758]142      ENDIF
[1482]143     
[1577]144      IF(ln_ctl)   CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hmlp, clinfo2=' hmlp : ', ovlap=1 )
[1559]145      !
[3294]146      CALL wrk_dealloc( jpi,jpj, imld )
[2715]147      !
[3294]148      IF( nn_timing == 1 )  CALL timing_stop('zdf_mxl')
149      !
[3]150   END SUBROUTINE zdf_mxl
151
152   !!======================================================================
153END MODULE zdfmxl
Note: See TracBrowser for help on using the repository browser.