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.
diahth.F90 in NEMO/branches/2019/fix_ticket2293/src/OCE/DIA – NEMO

source: NEMO/branches/2019/fix_ticket2293/src/OCE/DIA/diahth.F90 @ 11103

Last change on this file since 11103 was 11103, checked in by agn, 5 years ago

Added new code to branch

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