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

source: trunk/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfmxl.F90 @ 7698

Last change on this file since 7698 was 7698, checked in by mocavero, 7 years ago

update trunk with OpenMP parallelization

  • Property svn:keywords set to Id
File size: 7.7 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 trc_oce, ONLY: l_offline         ! ocean space and time domain variables
16   USE zdf_oce         ! ocean vertical physics
17   USE in_out_manager  ! I/O manager
18   USE prtctl          ! Print control
19   USE phycst          ! physical constants
20   USE iom             ! I/O library
21   USE lib_mpp         ! MPP library
22   USE wrk_nemo        ! work arrays
23   USE timing          ! Timing
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      ! w-level of the mixing and mixed layers
98      zN2_c = grav * rho_c * r1_rau0           ! convert density criteria into N^2 criteria
99!$OMP PARALLEL
100!$OMP DO schedule(static) private(jj, ji)
101      DO jj = 1, jpj
102         DO ji = 1, jpi
103            nmln(ji,jj)  = nlb10               ! Initialization to the number of w ocean point
104            hmlp(ji,jj)  = 0._wp               ! here hmlp used as a dummy variable, integrating vertically N^2
105         END DO
106      END DO
107      DO jk = nlb10, jpkm1
108!$OMP DO schedule(static) private(jj, ji, ikt)
109         DO jj = 1, jpj                ! Mixed layer level: w-level
110            DO ji = 1, jpi
111               ikt = mbkt(ji,jj)
112               hmlp(ji,jj) = hmlp(ji,jj) + MAX( rn2b(ji,jj,jk) , 0._wp ) * e3w_n(ji,jj,jk)
113               IF( hmlp(ji,jj) < zN2_c )   nmln(ji,jj) = MIN( jk , ikt ) + 1   ! Mixed layer level
114            END DO
115         END DO
116      END DO
117      !
118      ! w-level of the turbocline and mixing layer (iom_use)
119!$OMP DO schedule(static) private(jj, ji)
120      DO jj = 1, jpj
121         DO ji = 1, jpi
122            imld(ji,jj) = mbkt(ji,jj) + 1        ! Initialization to the number of w ocean point
123         END DO
124      END DO
125      DO jk = jpkm1, nlb10, -1         ! from the bottom to nlb10
126!$OMP DO schedule(static) private(jj, ji)
127         DO jj = 1, jpj
128            DO ji = 1, jpi
129               IF( avt (ji,jj,jk) < avt_c * wmask(ji,jj,jk) )   imld(ji,jj) = jk      ! Turbocline
130            END DO
131         END DO
132      END DO
133      ! depth of the mixing and mixed layers
134!$OMP DO schedule(static) private(jj, ji, iiki, iikn)
135      DO jj = 1, jpj
136         DO ji = 1, jpi
137            iiki = imld(ji,jj)
138            iikn = nmln(ji,jj)
139            hmld (ji,jj) = gdepw_n(ji,jj,iiki  ) * ssmask(ji,jj)    ! Turbocline depth
140            hmlp (ji,jj) = gdepw_n(ji,jj,iikn  ) * ssmask(ji,jj)    ! Mixed layer depth
141            hmlpt(ji,jj) = gdept_n(ji,jj,iikn-1) * ssmask(ji,jj)    ! depth of the last T-point inside the mixed layer
142         END DO
143      END DO
144!$OMP END PARALLEL
145      !
146      IF( .NOT.l_offline ) THEN
147         IF( iom_use("mldr10_1") ) THEN
148            IF( ln_isfcav ) THEN  ;  CALL iom_put( "mldr10_1", hmlp - risfdep)   ! mixed layer thickness
149            ELSE                  ;  CALL iom_put( "mldr10_1", hmlp )            ! mixed layer depth
150            END IF
151         END IF
152         IF( iom_use("mldkz5") ) THEN
153            IF( ln_isfcav ) THEN  ;  CALL iom_put( "mldkz5"  , hmld - risfdep )   ! turbocline thickness
154            ELSE                  ;  CALL iom_put( "mldkz5"  , hmld )             ! turbocline depth
155            END IF
156         ENDIF
157      ENDIF
158     
159      IF(ln_ctl)   CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hmlp, clinfo2=' hmlp : ', ovlap=1 )
160      !
161      CALL wrk_dealloc( jpi,jpj, imld )
162      !
163      IF( nn_timing == 1 )  CALL timing_stop('zdf_mxl')
164      !
165   END SUBROUTINE zdf_mxl
166
167   !!======================================================================
168END MODULE zdfmxl
Note: See TracBrowser for help on using the repository browser.