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/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfmxl.F90 @ 7037

Last change on this file since 7037 was 7037, checked in by mocavero, 8 years ago

ORCA2_LIM_PISCES hybrid version update

  • Property svn:keywords set to Id
File size: 7.6 KB
Line 
1MODULE zdfmxl
2   !!======================================================================
3   !!                       ***  MODULE  zdfmxl  ***
4   !! Ocean physics: mixed layer depth
5   !!======================================================================
6   !! History :  1.0  ! 2003-08  (G. Madec)  original code
7   !!            3.2  ! 2009-07  (S. Masson, G. Madec)  IOM + merge of DO-loop
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
10   !!----------------------------------------------------------------------
11   !!   zdf_mxl      : Compute the turbocline and mixed layer depths.
12   !!----------------------------------------------------------------------
13   USE oce             ! ocean dynamics and tracers variables
14   USE dom_oce         ! ocean space and time domain variables
15   USE zdf_oce         ! ocean vertical physics
16   USE in_out_manager  ! I/O manager
17   USE prtctl          ! Print control
18   USE phycst          ! physical constants
19   USE iom             ! I/O library
20   USE lib_mpp         ! MPP library
21   USE wrk_nemo        ! work arrays
22   USE timing          ! Timing
23   USE trc_oce, ONLY : lk_offline ! offline flag
24
25   IMPLICIT NONE
26   PRIVATE
27
28   PUBLIC   zdf_mxl       ! called by step.F90
29
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]
33   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hmlpt   !: depth of the last T-point inside the mixed layer [m]
34
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
38   !!----------------------------------------------------------------------
39   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
40   !! $Id$
41   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
42   !!----------------------------------------------------------------------
43CONTAINS
44
45   INTEGER FUNCTION zdf_mxl_alloc()
46      !!----------------------------------------------------------------------
47      !!               ***  FUNCTION zdf_mxl_alloc  ***
48      !!----------------------------------------------------------------------
49      zdf_mxl_alloc = 0      ! set to zero if no array to be allocated
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
57   END FUNCTION zdf_mxl_alloc
58
59
60   SUBROUTINE zdf_mxl( kt )
61      !!----------------------------------------------------------------------
62      !!                  ***  ROUTINE zdfmxl  ***
63      !!                   
64      !! ** Purpose :   Compute the turbocline depth and the mixed layer depth
65      !!              with density criteria.
66      !!
67      !! ** Method  :   The mixed layer depth is the shallowest W depth with
68      !!      the density of the corresponding T point (just bellow) bellow a
69      !!      given value defined locally as rho(10m) + rho_c
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
73      !!      value defined locally (avt_c here taken equal to 5 cm/s2 by default)
74      !!
75      !! ** Action  :   nmln, hmld, hmlp, hmlpt
76      !!----------------------------------------------------------------------
77      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
78      !
79      INTEGER  ::   ji, jj, jk      ! dummy loop indices
80      INTEGER  ::   iikn, iiki, ikt ! local integer
81      REAL(wp) ::   zN2_c           ! local scalar
82      INTEGER, POINTER, DIMENSION(:,:) ::   imld   ! 2D workspace
83      !!----------------------------------------------------------------------
84      !
85      IF( nn_timing == 1 )  CALL timing_start('zdf_mxl')
86      !
87      CALL wrk_alloc( jpi,jpj, imld )
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,*) '~~~~~~~ '
93         !                             ! allocate zdfmxl arrays
94         IF( zdf_mxl_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_mxl : unable to allocate arrays' )
95      ENDIF
96
97      zN2_c = grav * rho_c * r1_rau0   ! convert density criteria into N^2 criteria
98
99      ! w-level of the mixing and mixed layers
100!$OMP PARALLEL
101!$OMP WORKSHARE
102      nmln(:,:)  = nlb10               ! Initialization to the number of w ocean point
103      hmlp(:,:)  = 0._wp               ! here hmlp used as a dummy variable, integrating vertically N^2
104!$OMP END WORKSHARE
105      DO jk = nlb10, jpkm1
106!$OMP DO schedule(static) private(jj, ji, ikt)
107         DO jj = 1, jpj                ! Mixed layer level: w-level
108            DO ji = 1, jpi
109               ikt = mbkt(ji,jj)
110               hmlp(ji,jj) = hmlp(ji,jj) + MAX( rn2b(ji,jj,jk) , 0._wp ) * e3w_n(ji,jj,jk)
111               IF( hmlp(ji,jj) < zN2_c )   nmln(ji,jj) = MIN( jk , ikt ) + 1   ! Mixed layer level
112            END DO
113         END DO
114      END DO
115      !
116      ! w-level of the turbocline and mixing layer (iom_use)
117!$OMP WORKSHARE
118      imld(:,:) = mbkt(:,:) + 1        ! Initialization to the number of w ocean point
119!$OMP END WORKSHARE
120
121      DO jk = jpkm1, nlb10, -1         ! from the bottom to nlb10
122!$OMP DO schedule(static) private(jj, ji)
123         DO jj = 1, jpj
124            DO ji = 1, jpi
125               IF( avt (ji,jj,jk) < avt_c * wmask(ji,jj,jk) )   imld(ji,jj) = jk      ! Turbocline
126            END DO
127         END DO
128      END DO
129      ! depth of the mixing and mixed layers
130!$OMP DO schedule(static) private(jj, ji, iiki, iikn)
131      DO jj = 1, jpj
132         DO ji = 1, jpi
133            iiki = imld(ji,jj)
134            iikn = nmln(ji,jj)
135            hmld (ji,jj) = gdepw_n(ji,jj,iiki  ) * ssmask(ji,jj)    ! Turbocline depth
136            hmlp (ji,jj) = gdepw_n(ji,jj,iikn  ) * ssmask(ji,jj)    ! Mixed layer depth
137            hmlpt(ji,jj) = gdept_n(ji,jj,iikn-1) * ssmask(ji,jj)    ! depth of the last T-point inside the mixed layer
138         END DO
139      END DO
140!$OMP END DO NOWAIT
141!$OMP END PARALLEL
142      ! no need to output in offline mode
143      IF( .NOT.lk_offline ) THEN   
144         IF ( iom_use("mldr10_1") ) THEN
145            IF( ln_isfcav ) THEN
146               CALL iom_put( "mldr10_1", hmlp - risfdep)   ! mixed layer thickness
147            ELSE
148               CALL iom_put( "mldr10_1", hmlp )            ! mixed layer depth
149            END IF
150         END IF
151         IF ( iom_use("mldkz5") ) THEN
152            IF( ln_isfcav ) THEN
153               CALL iom_put( "mldkz5"  , hmld - risfdep )   ! turbocline thickness
154            ELSE
155               CALL iom_put( "mldkz5"  , hmld )             ! turbocline depth
156            END IF
157         END IF
158      ENDIF
159     
160      IF(ln_ctl)   CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hmlp, clinfo2=' hmlp : ', ovlap=1 )
161      !
162      CALL wrk_dealloc( jpi,jpj, imld )
163      !
164      IF( nn_timing == 1 )  CALL timing_stop('zdf_mxl')
165      !
166   END SUBROUTINE zdf_mxl
167
168   !!======================================================================
169END MODULE zdfmxl
Note: See TracBrowser for help on using the repository browser.