1 | MODULE 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 | ! |
---|
18 | USE in_out_manager ! I/O manager |
---|
19 | USE prtctl ! Print control |
---|
20 | USE phycst ! physical constants |
---|
21 | USE iom ! I/O library |
---|
22 | USE lib_mpp ! MPP library |
---|
23 | |
---|
24 | IMPLICIT NONE |
---|
25 | PRIVATE |
---|
26 | |
---|
27 | PUBLIC zdf_mxl ! called by zdfphy.F90 |
---|
28 | |
---|
29 | INTEGER , PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: nmln !: number of level in the mixed layer (used by LDF, ZDF, TRD, TOP) |
---|
30 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hmld !: mixing layer depth (turbocline) [m] (used by TOP) |
---|
31 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hmlp !: mixed layer depth (rho=rho0+zdcrit) [m] (used by LDF) |
---|
32 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hmlpt !: depth of the last T-point inside the mixed layer [m] (used by LDF) |
---|
33 | |
---|
34 | REAL(wp), PUBLIC :: rho_c = 0.01_wp !: density criterion for mixed layer depth |
---|
35 | REAL(wp), PUBLIC :: avt_c = 5.e-4_wp ! Kz criterion for the turbocline depth |
---|
36 | |
---|
37 | !!---------------------------------------------------------------------- |
---|
38 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
39 | !! $Id$ |
---|
40 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
41 | !!---------------------------------------------------------------------- |
---|
42 | CONTAINS |
---|
43 | |
---|
44 | INTEGER FUNCTION zdf_mxl_alloc() |
---|
45 | !!---------------------------------------------------------------------- |
---|
46 | !! *** FUNCTION zdf_mxl_alloc *** |
---|
47 | !!---------------------------------------------------------------------- |
---|
48 | zdf_mxl_alloc = 0 ! set to zero if no array to be allocated |
---|
49 | IF( .NOT. ALLOCATED( nmln ) ) THEN |
---|
50 | ALLOCATE( nmln(jpi,jpj), hmld(jpi,jpj), hmlp(jpi,jpj), hmlpt(jpi,jpj), STAT= zdf_mxl_alloc ) |
---|
51 | ! |
---|
52 | IF( lk_mpp ) CALL mpp_sum ( zdf_mxl_alloc ) |
---|
53 | IF( zdf_mxl_alloc /= 0 ) CALL ctl_warn('zdf_mxl_alloc: failed to allocate arrays.') |
---|
54 | ! |
---|
55 | ENDIF |
---|
56 | END FUNCTION zdf_mxl_alloc |
---|
57 | |
---|
58 | |
---|
59 | SUBROUTINE zdf_mxl( kt ) |
---|
60 | !!---------------------------------------------------------------------- |
---|
61 | !! *** ROUTINE zdfmxl *** |
---|
62 | !! |
---|
63 | !! ** Purpose : Compute the turbocline depth and the mixed layer depth |
---|
64 | !! with density criteria. |
---|
65 | !! |
---|
66 | !! ** Method : The mixed layer depth is the shallowest W depth with |
---|
67 | !! the density of the corresponding T point (just bellow) bellow a |
---|
68 | !! given value defined locally as rho(10m) + rho_c |
---|
69 | !! The turbocline depth is the depth at which the vertical |
---|
70 | !! eddy diffusivity coefficient (resulting from the vertical physics |
---|
71 | !! alone, not the isopycnal part, see trazdf.F) fall below a given |
---|
72 | !! value defined locally (avt_c here taken equal to 5 cm/s2 by default) |
---|
73 | !! |
---|
74 | !! ** Action : nmln, hmld, hmlp, hmlpt |
---|
75 | !!---------------------------------------------------------------------- |
---|
76 | INTEGER, INTENT(in) :: kt ! ocean time-step index |
---|
77 | ! |
---|
78 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
79 | INTEGER :: iikn, iiki, ikt ! local integer |
---|
80 | REAL(wp) :: zN2_c ! local scalar |
---|
81 | INTEGER, DIMENSION(jpi,jpj) :: imld ! 2D workspace |
---|
82 | !!---------------------------------------------------------------------- |
---|
83 | ! |
---|
84 | IF( kt == nit000 ) THEN |
---|
85 | IF(lwp) WRITE(numout,*) |
---|
86 | IF(lwp) WRITE(numout,*) 'zdf_mxl : mixed layer depth' |
---|
87 | IF(lwp) WRITE(numout,*) '~~~~~~~ ' |
---|
88 | ! ! allocate zdfmxl arrays |
---|
89 | IF( zdf_mxl_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'zdf_mxl : unable to allocate arrays' ) |
---|
90 | ENDIF |
---|
91 | ! |
---|
92 | ! w-level of the mixing and mixed layers |
---|
93 | nmln(:,:) = nlb10 ! Initialization to the number of w ocean point |
---|
94 | hmlp(:,:) = 0._wp ! here hmlp used as a dummy variable, integrating vertically N^2 |
---|
95 | zN2_c = grav * rho_c * r1_rau0 ! convert density criteria into N^2 criteria |
---|
96 | DO jk = nlb10, jpkm1 |
---|
97 | DO jj = 1, jpj ! Mixed layer level: w-level |
---|
98 | DO ji = 1, jpi |
---|
99 | ikt = mbkt(ji,jj) |
---|
100 | hmlp(ji,jj) = hmlp(ji,jj) + MAX( rn2b(ji,jj,jk) , 0._wp ) * e3w_n(ji,jj,jk) |
---|
101 | IF( hmlp(ji,jj) < zN2_c ) nmln(ji,jj) = MIN( jk , ikt ) + 1 ! Mixed layer level |
---|
102 | END DO |
---|
103 | END DO |
---|
104 | END DO |
---|
105 | ! |
---|
106 | ! w-level of the turbocline and mixing layer (iom_use) |
---|
107 | imld(:,:) = mbkt(:,:) + 1 ! Initialization to the number of w ocean point |
---|
108 | DO jk = jpkm1, nlb10, -1 ! from the bottom to nlb10 |
---|
109 | DO jj = 1, jpj |
---|
110 | DO ji = 1, jpi |
---|
111 | IF( avt (ji,jj,jk) < avt_c * wmask(ji,jj,jk) ) imld(ji,jj) = jk ! Turbocline |
---|
112 | END DO |
---|
113 | END DO |
---|
114 | END DO |
---|
115 | ! depth of the mixing and mixed layers |
---|
116 | DO jj = 1, jpj |
---|
117 | DO ji = 1, jpi |
---|
118 | iiki = imld(ji,jj) |
---|
119 | iikn = nmln(ji,jj) |
---|
120 | hmld (ji,jj) = gdepw_n(ji,jj,iiki ) * ssmask(ji,jj) ! Turbocline depth |
---|
121 | hmlp (ji,jj) = gdepw_n(ji,jj,iikn ) * ssmask(ji,jj) ! Mixed layer depth |
---|
122 | hmlpt(ji,jj) = gdept_n(ji,jj,iikn-1) * ssmask(ji,jj) ! depth of the last T-point inside the mixed layer |
---|
123 | END DO |
---|
124 | END DO |
---|
125 | ! |
---|
126 | IF( .NOT.l_offline ) THEN |
---|
127 | IF( iom_use("mldr10_1") ) THEN |
---|
128 | IF( ln_isfcav ) THEN ; CALL iom_put( "mldr10_1", hmlp - risfdep) ! mixed layer thickness |
---|
129 | ELSE ; CALL iom_put( "mldr10_1", hmlp ) ! mixed layer depth |
---|
130 | END IF |
---|
131 | END IF |
---|
132 | IF( iom_use("mldkz5") ) THEN |
---|
133 | IF( ln_isfcav ) THEN ; CALL iom_put( "mldkz5" , hmld - risfdep ) ! turbocline thickness |
---|
134 | ELSE ; CALL iom_put( "mldkz5" , hmld ) ! turbocline depth |
---|
135 | END IF |
---|
136 | ENDIF |
---|
137 | ENDIF |
---|
138 | ! |
---|
139 | IF(ln_ctl) CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hmlp, clinfo2=' hmlp : ' ) |
---|
140 | ! |
---|
141 | END SUBROUTINE zdf_mxl |
---|
142 | |
---|
143 | !!====================================================================== |
---|
144 | END MODULE zdfmxl |
---|