1 | MODULE diahth |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE diahth *** |
---|
4 | !! Ocean diagnostics: thermocline and 20 degree depth |
---|
5 | !!====================================================================== |
---|
6 | !! History : OPA ! 1994-09 (J.-P. Boulanger) Original code |
---|
7 | !! ! 1996-11 (E. Guilyardi) OPA8 |
---|
8 | !! ! 1997-08 (G. Madec) optimization |
---|
9 | !! ! 1999-07 (E. Guilyardi) hd28 + heat content |
---|
10 | !! NEMO 1.0 ! 2002-06 (G. Madec) F90: Free form and module |
---|
11 | !! 3.2 ! 2009-07 (S. Masson) hc300 bugfix + cleaning + add new diag |
---|
12 | !!---------------------------------------------------------------------- |
---|
13 | !!---------------------------------------------------------------------- |
---|
14 | !! 'key_diahth' : thermocline depth diag. |
---|
15 | !!---------------------------------------------------------------------- |
---|
16 | !! dia_hth : Compute varius diagnostics associated with the mixed layer |
---|
17 | !!---------------------------------------------------------------------- |
---|
18 | USE oce ! ocean dynamics and tracers |
---|
19 | USE dom_oce ! ocean space and time domain |
---|
20 | USE phycst ! physical constants |
---|
21 | ! |
---|
22 | USE in_out_manager ! I/O manager |
---|
23 | USE lib_mpp ! MPP library |
---|
24 | USE iom ! I/O library |
---|
25 | |
---|
26 | IMPLICIT NONE |
---|
27 | PRIVATE |
---|
28 | |
---|
29 | PUBLIC dia_hth ! routine called by step.F90 |
---|
30 | PUBLIC dia_hth_init ! routine called by nemogcm.F90 |
---|
31 | |
---|
32 | LOGICAL, PUBLIC :: ln_diahth !: Compute further diagnostics of ML and thermocline depth |
---|
33 | |
---|
34 | !!---------------------------------------------------------------------- |
---|
35 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
36 | !! $Id$ |
---|
37 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
38 | !!---------------------------------------------------------------------- |
---|
39 | CONTAINS |
---|
40 | |
---|
41 | |
---|
42 | SUBROUTINE dia_hth( kt ) |
---|
43 | !!--------------------------------------------------------------------- |
---|
44 | !! *** ROUTINE dia_hth *** |
---|
45 | !! |
---|
46 | !! ** Purpose : Computes |
---|
47 | !! the mixing layer depth (turbocline): avt = 5.e-4 |
---|
48 | !! the depth of strongest vertical temperature gradient |
---|
49 | !! the mixed layer depth with density criteria: rho = rho(10m or surf) + 0.03(or 0.01) |
---|
50 | !! the mixed layer depth with temperature criteria: abs( tn - tn(10m) ) = 0.2 |
---|
51 | !! the top of the thermochine: tn = tn(10m) - ztem2 |
---|
52 | !! the pycnocline depth with density criteria equivalent to a temperature variation |
---|
53 | !! rho = rho10m + (dr/dT)(T,S,10m)*(-0.2 degC) |
---|
54 | !! the barrier layer thickness |
---|
55 | !! the maximal verical inversion of temperature and its depth max( 0, max of tn - tn(10m) ) |
---|
56 | !! the depth of the 20 degree isotherm (linear interpolation) |
---|
57 | !! the depth of the 28 degree isotherm (linear interpolation) |
---|
58 | !! the heat content of first 300 m |
---|
59 | !! |
---|
60 | !! ** Method : |
---|
61 | !!------------------------------------------------------------------- |
---|
62 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
63 | !! |
---|
64 | INTEGER :: ji, jj, jk ! dummy loop arguments |
---|
65 | INTEGER :: iid, ilevel ! temporary integers |
---|
66 | INTEGER, DIMENSION(jpi,jpj) :: ik20, ik28 ! levels |
---|
67 | REAL(wp) :: zavt5 = 5.e-4_wp ! Kz criterion for the turbocline depth |
---|
68 | REAL(wp) :: zrho3 = 0.03_wp ! density criterion for mixed layer depth |
---|
69 | REAL(wp) :: zrho1 = 0.01_wp ! density criterion for mixed layer depth |
---|
70 | REAL(wp) :: ztem2 = 0.2_wp ! temperature criterion for mixed layer depth |
---|
71 | REAL(wp) :: zthick_0, zcoef ! temporary scalars |
---|
72 | REAL(wp) :: zztmp, zzdep ! temporary scalars inside do loop |
---|
73 | REAL(wp) :: zu, zv, zw, zut, zvt ! temporary workspace |
---|
74 | REAL(wp), DIMENSION(jpi,jpj) :: zabs2 ! MLD: abs( tn - tn(10m) ) = ztem2 |
---|
75 | REAL(wp), DIMENSION(jpi,jpj) :: ztm2 ! Top of thermocline: tn = tn(10m) - ztem2 |
---|
76 | REAL(wp), DIMENSION(jpi,jpj) :: zrho10_3 ! MLD: rho = rho10m + zrho3 |
---|
77 | REAL(wp), DIMENSION(jpi,jpj) :: zpycn ! pycnocline: rho = rho10m + (dr/dT)(T,S,10m)*(-0.2 degC) |
---|
78 | REAL(wp), DIMENSION(jpi,jpj) :: ztinv ! max of temperature inversion |
---|
79 | REAL(wp), DIMENSION(jpi,jpj) :: zdepinv ! depth of temperature inversion |
---|
80 | REAL(wp), DIMENSION(jpi,jpj) :: zrho0_3 ! MLD rho = rho(surf) = 0.03 |
---|
81 | REAL(wp), DIMENSION(jpi,jpj) :: zrho0_1 ! MLD rho = rho(surf) = 0.01 |
---|
82 | REAL(wp), DIMENSION(jpi,jpj) :: zmaxdzT ! max of dT/dz |
---|
83 | REAL(wp), DIMENSION(jpi,jpj) :: zthick ! vertical integration thickness |
---|
84 | REAL(wp), DIMENSION(jpi,jpj) :: zdelr ! delta rho equivalent to deltaT = 0.2 |
---|
85 | ! note: following variables should move to local variables once iom_put is always used |
---|
86 | REAL(wp), DIMENSION(jpi,jpj) :: zhth !: depth of the max vertical temperature gradient [m] |
---|
87 | REAL(wp), DIMENSION(jpi,jpj) :: zhd20 !: depth of 20 C isotherm [m] |
---|
88 | REAL(wp), DIMENSION(jpi,jpj) :: zhd28 !: depth of 28 C isotherm [m] |
---|
89 | REAL(wp), DIMENSION(jpi,jpj) :: zhtc3 !: heat content of first 300 m [W] |
---|
90 | |
---|
91 | IF (iom_use("mlddzt").OR.iom_use("mldr0_3").OR.iom_use("mldr0_1")) THEN |
---|
92 | ! ------------------------------------------------------------- ! |
---|
93 | ! thermocline depth: strongest vertical gradient of temperature ! |
---|
94 | ! turbocline depth (mixing layer depth): avt = zavt5 ! |
---|
95 | ! MLD: rho = rho(1) + zrho3 ! |
---|
96 | ! MLD: rho = rho(1) + zrho1 ! |
---|
97 | ! ------------------------------------------------------------- ! |
---|
98 | zmaxdzT(:,:) = 0._wp |
---|
99 | IF( nla10 > 1 ) THEN |
---|
100 | DO jj = 1, jpj |
---|
101 | DO ji = 1, jpi |
---|
102 | zztmp = gdepw_n(ji,jj,mbkt(ji,jj)+1) |
---|
103 | zrho0_3(ji,jj) = zztmp |
---|
104 | zrho0_1(ji,jj) = zztmp |
---|
105 | zhth(ji,jj) = zztmp |
---|
106 | END DO |
---|
107 | END DO |
---|
108 | ELSE IF (iom_use("mlddzt")) THEN |
---|
109 | DO jj = 1, jpj |
---|
110 | DO ji = 1, jpi |
---|
111 | zztmp = gdepw_n(ji,jj,mbkt(ji,jj)+1) |
---|
112 | zhth(ji,jj) = zztmp |
---|
113 | END DO |
---|
114 | END DO |
---|
115 | ELSE |
---|
116 | zhth(:,:) = 0._wp |
---|
117 | ENDIF |
---|
118 | |
---|
119 | DO jk = jpkm1, 2, -1 ! loop from bottom to 2 |
---|
120 | DO jj = 1, jpj |
---|
121 | DO ji = 1, jpi |
---|
122 | ! |
---|
123 | zzdep = gdepw_n(ji,jj,jk) |
---|
124 | zztmp = ( tsn(ji,jj,jk-1,jp_tem) - tsn(ji,jj,jk,jp_tem) ) / zzdep * tmask(ji,jj,jk) ! vertical gradient of temperature (dT/dz) |
---|
125 | zzdep = zzdep * tmask(ji,jj,1) |
---|
126 | |
---|
127 | IF( zztmp > zmaxdzT(ji,jj) ) THEN |
---|
128 | zmaxdzT(ji,jj) = zztmp ; zhth (ji,jj) = zzdep ! max and depth of dT/dz |
---|
129 | ENDIF |
---|
130 | |
---|
131 | IF( nla10 > 1 ) THEN |
---|
132 | zztmp = rhop(ji,jj,jk) - rhop(ji,jj,1) ! delta rho(1) |
---|
133 | IF( zztmp > zrho3 ) zrho0_3(ji,jj) = zzdep ! > 0.03 |
---|
134 | IF( zztmp > zrho1 ) zrho0_1(ji,jj) = zzdep ! > 0.01 |
---|
135 | ENDIF |
---|
136 | |
---|
137 | END DO |
---|
138 | END DO |
---|
139 | END DO |
---|
140 | |
---|
141 | IF (iom_use("mlddzt")) CALL iom_put( "mlddzt", zhth*tmask(:,:,1) ) ! depth of the thermocline |
---|
142 | IF( nla10 > 1 ) THEN |
---|
143 | IF (iom_use("mldr0_3")) CALL iom_put( "mldr0_3", zrho0_3*tmask(:,:,1) ) ! MLD delta rho(surf) = 0.03 |
---|
144 | IF (iom_use("mldr0_1")) CALL iom_put( "mldr0_1", zrho0_1*tmask(:,:,1) ) ! MLD delta rho(surf) = 0.01 |
---|
145 | ENDIF |
---|
146 | ENDIF |
---|
147 | |
---|
148 | IF (iom_use("mld_dt02").OR.iom_use("topthdep").OR.iom_use("mldr10_3").OR. & |
---|
149 | & iom_use("pycndep").OR.iom_use("tinv").OR.iom_use("depti")) THEN |
---|
150 | DO jj = 1, jpj |
---|
151 | DO ji = 1, jpi |
---|
152 | zztmp = gdepw_n(ji,jj,mbkt(ji,jj)+1) |
---|
153 | zabs2 (ji,jj) = zztmp |
---|
154 | ztm2 (ji,jj) = zztmp |
---|
155 | zrho10_3(ji,jj) = zztmp |
---|
156 | zpycn (ji,jj) = zztmp |
---|
157 | END DO |
---|
158 | END DO |
---|
159 | ztinv (:,:) = 0._wp |
---|
160 | zdepinv(:,:) = 0._wp |
---|
161 | |
---|
162 | IF (iom_use("pycndep")) THEN |
---|
163 | ! Preliminary computation |
---|
164 | ! computation of zdelr = (dr/dT)(T,S,10m)*(-0.2 degC) |
---|
165 | DO jj = 1, jpj |
---|
166 | DO ji = 1, jpi |
---|
167 | IF( tmask(ji,jj,nla10) == 1. ) THEN |
---|
168 | zu = 1779.50 + 11.250 * tsn(ji,jj,nla10,jp_tem) - 3.80 * tsn(ji,jj,nla10,jp_sal) & |
---|
169 | & - 0.0745 * tsn(ji,jj,nla10,jp_tem) * tsn(ji,jj,nla10,jp_tem) & |
---|
170 | & - 0.0100 * tsn(ji,jj,nla10,jp_tem) * tsn(ji,jj,nla10,jp_sal) |
---|
171 | zv = 5891.00 + 38.000 * tsn(ji,jj,nla10,jp_tem) + 3.00 * tsn(ji,jj,nla10,jp_sal) & |
---|
172 | & - 0.3750 * tsn(ji,jj,nla10,jp_tem) * tsn(ji,jj,nla10,jp_tem) |
---|
173 | zut = 11.25 - 0.149 * tsn(ji,jj,nla10,jp_tem) - 0.01 * tsn(ji,jj,nla10,jp_sal) |
---|
174 | zvt = 38.00 - 0.750 * tsn(ji,jj,nla10,jp_tem) |
---|
175 | zw = (zu + 0.698*zv) * (zu + 0.698*zv) |
---|
176 | zdelr(ji,jj) = ztem2 * (1000.*(zut*zv - zvt*zu)/zw) |
---|
177 | ELSE |
---|
178 | zdelr(ji,jj) = 0._wp |
---|
179 | ENDIF |
---|
180 | END DO |
---|
181 | END DO |
---|
182 | ELSE |
---|
183 | zdelr(:,:) = 0._wp |
---|
184 | ENDIF |
---|
185 | |
---|
186 | ! ------------------------------------------------------------- ! |
---|
187 | ! MLD: abs( tn - tn(10m) ) = ztem2 ! |
---|
188 | ! Top of thermocline: tn = tn(10m) - ztem2 ! |
---|
189 | ! MLD: rho = rho10m + zrho3 ! |
---|
190 | ! pycnocline: rho = rho10m + (dr/dT)(T,S,10m)*(-0.2 degC) ! |
---|
191 | ! temperature inversion: max( 0, max of tn - tn(10m) ) ! |
---|
192 | ! depth of temperature inversion ! |
---|
193 | ! ------------------------------------------------------------- ! |
---|
194 | DO jk = jpkm1, nlb10, -1 ! loop from bottom to nlb10 |
---|
195 | DO jj = 1, jpj |
---|
196 | DO ji = 1, jpi |
---|
197 | ! |
---|
198 | zzdep = gdepw_n(ji,jj,jk) * tmask(ji,jj,1) |
---|
199 | ! |
---|
200 | zztmp = tsn(ji,jj,nla10,jp_tem) - tsn(ji,jj,jk,jp_tem) ! - delta T(10m) |
---|
201 | IF( ABS(zztmp) > ztem2 ) zabs2 (ji,jj) = zzdep ! abs > 0.2 |
---|
202 | IF( zztmp > ztem2 ) ztm2 (ji,jj) = zzdep ! > 0.2 |
---|
203 | zztmp = -zztmp ! delta T(10m) |
---|
204 | IF( zztmp > ztinv(ji,jj) ) THEN ! temperature inversion |
---|
205 | ztinv(ji,jj) = zztmp ; zdepinv (ji,jj) = zzdep ! max value and depth |
---|
206 | ENDIF |
---|
207 | |
---|
208 | zztmp = rhop(ji,jj,jk) - rhop(ji,jj,nla10) ! delta rho(10m) |
---|
209 | IF( zztmp > zrho3 ) zrho10_3(ji,jj) = zzdep ! > 0.03 |
---|
210 | IF( zztmp > zdelr(ji,jj) ) zpycn (ji,jj) = zzdep ! > equi. delta T(10m) - 0.2 |
---|
211 | ! |
---|
212 | END DO |
---|
213 | END DO |
---|
214 | END DO |
---|
215 | |
---|
216 | IF (iom_use("mld_dt02")) CALL iom_put( "mld_dt02", zabs2*tmask(:,:,1) ) ! MLD abs(delta t) - 0.2 |
---|
217 | IF (iom_use("topthdep")) CALL iom_put( "topthdep", ztm2*tmask(:,:,1) ) ! T(10) - 0.2 |
---|
218 | IF (iom_use("mldr10_3")) CALL iom_put( "mldr10_3", zrho10_3*tmask(:,:,1) ) ! MLD delta rho(10m) = 0.03 |
---|
219 | IF (iom_use("pycndep")) CALL iom_put( "pycndep" , zpycn*tmask(:,:,1) ) ! MLD delta rho equi. delta T(10m) = 0.2 |
---|
220 | IF (iom_use("tinv")) CALL iom_put( "tinv" , ztinv*tmask(:,:,1) ) ! max. temp. inv. (t10 ref) |
---|
221 | IF (iom_use("depti")) CALL iom_put( "depti" , zdepinv*tmask(:,:,1) ) ! depth of max. temp. inv. (t10 ref) |
---|
222 | ENDIF |
---|
223 | |
---|
224 | IF(iom_use("20d").OR.iom_use("28d")) THEN |
---|
225 | ! ----------------------------------- ! |
---|
226 | ! search deepest level above 20C/28C ! |
---|
227 | ! ----------------------------------- ! |
---|
228 | ik20(:,:) = 1 |
---|
229 | ik28(:,:) = 1 |
---|
230 | DO jk = 1, jpkm1 ! beware temperature is not always decreasing with depth => loop from top to bottom |
---|
231 | DO jj = 1, jpj |
---|
232 | DO ji = 1, jpi |
---|
233 | zztmp = tsn(ji,jj,jk,jp_tem) |
---|
234 | IF( zztmp >= 20. ) ik20(ji,jj) = jk |
---|
235 | IF( zztmp >= 28. ) ik28(ji,jj) = jk |
---|
236 | END DO |
---|
237 | END DO |
---|
238 | END DO |
---|
239 | |
---|
240 | ! --------------------------- ! |
---|
241 | ! Depth of 20C/28C isotherm ! |
---|
242 | ! --------------------------- ! |
---|
243 | DO jj = 1, jpj |
---|
244 | DO ji = 1, jpi |
---|
245 | ! |
---|
246 | zzdep = gdepw_n(ji,jj,mbkt(ji,jj)+1) ! depth of the oean bottom |
---|
247 | ! |
---|
248 | iid = ik20(ji,jj) |
---|
249 | IF( iid /= 1 ) THEN |
---|
250 | zztmp = gdept_n(ji,jj,iid ) & ! linear interpolation |
---|
251 | & + ( gdept_n(ji,jj,iid+1) - gdept_n(ji,jj,iid) ) & |
---|
252 | & * ( 20.*tmask(ji,jj,iid+1) - tsn(ji,jj,iid,jp_tem) ) & |
---|
253 | & / ( tsn(ji,jj,iid+1,jp_tem) - tsn(ji,jj,iid,jp_tem) + (1.-tmask(ji,jj,1)) ) |
---|
254 | zhd20(ji,jj) = MIN( zztmp , zzdep) * tmask(ji,jj,1) ! bound by the ocean depth |
---|
255 | ELSE |
---|
256 | zhd20(ji,jj) = 0._wp |
---|
257 | ENDIF |
---|
258 | ! |
---|
259 | iid = ik28(ji,jj) |
---|
260 | IF( iid /= 1 ) THEN |
---|
261 | zztmp = gdept_n(ji,jj,iid ) & ! linear interpolation |
---|
262 | & + ( gdept_n(ji,jj,iid+1) - gdept_n(ji,jj,iid) ) & |
---|
263 | & * ( 28.*tmask(ji,jj,iid+1) - tsn(ji,jj,iid,jp_tem) ) & |
---|
264 | & / ( tsn(ji,jj,iid+1,jp_tem) - tsn(ji,jj,iid,jp_tem) + (1.-tmask(ji,jj,1)) ) |
---|
265 | zhd28(ji,jj) = MIN( zztmp , zzdep ) * tmask(ji,jj,1) ! bound by the ocean depth |
---|
266 | ELSE |
---|
267 | zhd28(ji,jj) = 0._wp |
---|
268 | ENDIF |
---|
269 | |
---|
270 | END DO |
---|
271 | END DO |
---|
272 | CALL iom_put( "20d", zhd20 ) ! depth of the 20 isotherm |
---|
273 | CALL iom_put( "28d", zhd28 ) ! depth of the 28 isotherm |
---|
274 | ENDIF |
---|
275 | |
---|
276 | ! ----------------------------- ! |
---|
277 | ! Heat content of first 300 m ! |
---|
278 | ! ----------------------------- ! |
---|
279 | IF (iom_use("hc300")) THEN |
---|
280 | ! find ilevel with (ilevel+1) the deepest W-level above 300m (we assume we can use e3t_1d to do this search...) |
---|
281 | ilevel = 0 |
---|
282 | zthick_0 = 0._wp |
---|
283 | DO jk = 1, jpkm1 |
---|
284 | zthick_0 = zthick_0 + e3t_1d(jk) |
---|
285 | IF( zthick_0 < 300. ) ilevel = jk |
---|
286 | END DO |
---|
287 | ! surface boundary condition |
---|
288 | IF( ln_linssh ) THEN ; zthick(:,:) = sshn(:,:) ; zhtc3(:,:) = tsn(:,:,1,jp_tem) * sshn(:,:) * tmask(:,:,1) |
---|
289 | ELSE ; zthick(:,:) = 0._wp ; zhtc3(:,:) = 0._wp |
---|
290 | ENDIF |
---|
291 | ! integration down to ilevel |
---|
292 | DO jk = 1, ilevel |
---|
293 | zthick(:,:) = zthick(:,:) + e3t_n(:,:,jk) |
---|
294 | zhtc3 (:,:) = zhtc3 (:,:) + e3t_n(:,:,jk) * tsn(:,:,jk,jp_tem) * tmask(:,:,jk) |
---|
295 | END DO |
---|
296 | ! deepest layer |
---|
297 | zthick(:,:) = 300. - zthick(:,:) ! remaining thickness to reach 300m |
---|
298 | DO jj = 1, jpj |
---|
299 | DO ji = 1, jpi |
---|
300 | zhtc3(ji,jj) = zhtc3(ji,jj) + tsn(ji,jj,ilevel+1,jp_tem) & |
---|
301 | & * MIN( e3t_n(ji,jj,ilevel+1), zthick(ji,jj) ) * tmask(ji,jj,ilevel+1) |
---|
302 | END DO |
---|
303 | END DO |
---|
304 | ! from temperature to heat contain |
---|
305 | zcoef = rau0 * rcp |
---|
306 | zhtc3(:,:) = zcoef * zhtc3(:,:) |
---|
307 | CALL iom_put( "hc300", zhtc3*tmask(:,:,1) ) ! first 300m heat content |
---|
308 | ENDIF |
---|
309 | ! |
---|
310 | END SUBROUTINE dia_hth |
---|
311 | |
---|
312 | |
---|
313 | SUBROUTINE dia_hth_init |
---|
314 | !!--------------------------------------------------------------------------- |
---|
315 | !! *** ROUTINE dia_hth_init *** |
---|
316 | !! |
---|
317 | !! ** Purpose: Initialization for ML and thermocline depths |
---|
318 | !! |
---|
319 | !! ** Action : Read namelist flag to permit or not extra Ml depths and thermocline depths |
---|
320 | !!--------------------------------------------------------------------------- |
---|
321 | INTEGER :: ierror, ios ! local integer |
---|
322 | !! |
---|
323 | NAMELIST/namhth/ ln_diahth |
---|
324 | !!--------------------------------------------------------------------------- |
---|
325 | ! |
---|
326 | IF(lwp) THEN |
---|
327 | WRITE(numout,*) |
---|
328 | WRITE(numout,*) 'dia_hth_init : heat and salt budgets diagnostics' |
---|
329 | WRITE(numout,*) '~~~~~~~~~~~~ ' |
---|
330 | ENDIF |
---|
331 | REWIND( numnam_ref ) ! Namelist namhth in reference namelist |
---|
332 | READ ( numnam_ref, namhth, IOSTAT = ios, ERR = 901) |
---|
333 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namhth in reference namelist', lwp ) |
---|
334 | REWIND( numnam_cfg ) ! Namelist namhth in configuration namelist |
---|
335 | READ ( numnam_cfg, namhth, IOSTAT = ios, ERR = 902 ) |
---|
336 | 902 IF( ios > 0 ) CALL ctl_nam ( ios , 'namhth in configuration namelist', lwp ) |
---|
337 | IF(lwm) WRITE( numond, namhth ) |
---|
338 | |
---|
339 | IF(lwp) THEN |
---|
340 | WRITE(numout,*) ' Namelist namhth :' |
---|
341 | WRITE(numout,*) ' check the heat and salt budgets (T) or not (F) ln_diahth = ', ln_diahth |
---|
342 | ENDIF |
---|
343 | ! |
---|
344 | END SUBROUTINE dia_hth_init |
---|
345 | END MODULE diahth |
---|