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/2020/dev_r12558_HPC-08_epico_Extra_Halo/src/OCE/DIA – NEMO

source: NEMO/branches/2020/dev_r12558_HPC-08_epico_Extra_Halo/src/OCE/DIA/diahth.F90 @ 13248

Last change on this file since 13248 was 13248, checked in by francesca, 4 years ago

dev_r12558_HPC-08_epico_Extra_Halo: merge with trunk@13237, see #2366

  • Property svn:keywords set to Id
File size: 18.0 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   USE timing          ! preformance summary
23
24   IMPLICIT NONE
25   PRIVATE
26
27   PUBLIC   dia_hth       ! routine called by step.F90
28   PUBLIC   dia_hth_alloc ! routine called by nemogcm.F90
29
30   LOGICAL, SAVE  ::   l_hth     !: thermocline-20d depths flag
31   
32   ! note: following variables should move to local variables once iom_put is always used
33   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hth    !: depth of the max vertical temperature gradient [m]
34   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hd20   !: depth of 20 C isotherm                         [m]
35   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hd26   !: depth of 26 C isotherm                         [m]
36   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hd28   !: depth of 28 C isotherm                         [m]
37   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   htc3   !: heat content of first 300 m                    [W]
38   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   htc7   !: heat content of first 700 m                    [W]
39   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   htc20  !: heat content of first 2000 m                   [W]
40
41
42   !! * Substitutions
43#  include "do_loop_substitute.h90"
44#  include "domzgr_substitute.h90"
45   !!----------------------------------------------------------------------
46   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
47   !! $Id$
48   !! Software governed by the CeCILL license (see ./LICENSE)
49   !!----------------------------------------------------------------------
50CONTAINS
51
52   FUNCTION dia_hth_alloc()
53      !!---------------------------------------------------------------------
54      INTEGER :: dia_hth_alloc
55      !!---------------------------------------------------------------------
56      !
57      ALLOCATE( hth(jpi,jpj), hd20(jpi,jpj), hd26(jpi,jpj), hd28(jpi,jpj), &
58         &      htc3(jpi,jpj), htc7(jpi,jpj), htc20(jpi,jpj), STAT=dia_hth_alloc )
59      !
60      CALL mpp_sum ( 'diahth', dia_hth_alloc )
61      IF(dia_hth_alloc /= 0)   CALL ctl_stop( 'STOP', 'dia_hth_alloc: failed to allocate arrays.' )
62      !
63   END FUNCTION dia_hth_alloc
64
65
66   SUBROUTINE dia_hth( kt, Kmm )
67      !!---------------------------------------------------------------------
68      !!                  ***  ROUTINE dia_hth  ***
69      !!
70      !! ** Purpose : Computes
71      !!      the mixing layer depth (turbocline): avt = 5.e-4
72      !!      the depth of strongest vertical temperature gradient
73      !!      the mixed layer depth with density     criteria: rho = rho(10m or surf) + 0.03(or 0.01)
74      !!      the mixed layer depth with temperature criteria: abs( tn - tn(10m) ) = 0.2       
75      !!      the top of the thermochine: tn = tn(10m) - ztem2
76      !!      the pycnocline depth with density criteria equivalent to a temperature variation
77      !!                rho = rho10m + (dr/dT)(T,S,10m)*(-0.2 degC)
78      !!      the barrier layer thickness
79      !!      the maximal verical inversion of temperature and its depth max( 0, max of tn - tn(10m) )
80      !!      the depth of the 20 degree isotherm (linear interpolation)
81      !!      the depth of the 28 degree isotherm (linear interpolation)
82      !!      the heat content of first 300 m
83      !!
84      !! ** Method :
85      !!-------------------------------------------------------------------
86      INTEGER, INTENT( in ) ::   kt      ! ocean time-step index
87      INTEGER, INTENT( in ) ::   Kmm     ! ocean time level index
88      !!
89      INTEGER                      ::   ji, jj, jk            ! dummy loop arguments
90      REAL(wp)                     ::   zrho3 = 0.03_wp       ! density     criterion for mixed layer depth
91      REAL(wp)                     ::   zrho1 = 0.01_wp       ! density     criterion for mixed layer depth
92      REAL(wp)                     ::   ztem2 = 0.2_wp        ! temperature criterion for mixed layer depth
93      REAL(wp)                     ::   zztmp, zzdep          ! temporary scalars inside do loop
94      REAL(wp)                     ::   zu, zv, zw, zut, zvt  ! temporary workspace
95      REAL(wp), DIMENSION(jpi,jpj) ::   zabs2      ! MLD: abs( tn - tn(10m) ) = ztem2
96      REAL(wp), DIMENSION(jpi,jpj) ::   ztm2       ! Top of thermocline: tn = tn(10m) - ztem2     
97      REAL(wp), DIMENSION(jpi,jpj) ::   zrho10_3   ! MLD: rho = rho10m + zrho3     
98      REAL(wp), DIMENSION(jpi,jpj) ::   zpycn      ! pycnocline: rho = rho10m + (dr/dT)(T,S,10m)*(-0.2 degC)
99      REAL(wp), DIMENSION(jpi,jpj) ::   ztinv      ! max of temperature inversion
100      REAL(wp), DIMENSION(jpi,jpj) ::   zdepinv    ! depth of temperature inversion
101      REAL(wp), DIMENSION(jpi,jpj) ::   zrho0_3    ! MLD rho = rho(surf) = 0.03
102      REAL(wp), DIMENSION(jpi,jpj) ::   zrho0_1    ! MLD rho = rho(surf) = 0.01
103      REAL(wp), DIMENSION(jpi,jpj) ::   zmaxdzT    ! max of dT/dz
104      REAL(wp), DIMENSION(jpi,jpj) ::   zdelr      ! delta rho equivalent to deltaT = 0.2
105      !!----------------------------------------------------------------------
106      IF( ln_timing )   CALL timing_start('dia_hth')
107
108      IF( kt == nit000 ) THEN
109         l_hth = .FALSE.
110         IF(   iom_use( 'mlddzt'   ) .OR. iom_use( 'mldr0_3'  ) .OR. iom_use( 'mldr0_1'  )    .OR.  & 
111            &  iom_use( 'mld_dt02' ) .OR. iom_use( 'topthdep' ) .OR. iom_use( 'mldr10_3' )    .OR.  &   
112            &  iom_use( '20d'      ) .OR. iom_use( '26d'      ) .OR. iom_use( '28d'      )    .OR.  &   
113            &  iom_use( 'hc300'    ) .OR. iom_use( 'hc700'    ) .OR. iom_use( 'hc2000'   )    .OR.  &   
114            &  iom_use( 'pycndep'  ) .OR. iom_use( 'tinv'     ) .OR. iom_use( 'depti'    )  ) l_hth = .TRUE.
115         !                                      ! allocate dia_hth array
116         IF( l_hth ) THEN
117            IF( dia_hth_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_hth : unable to allocate standard arrays' )
118            IF(lwp) WRITE(numout,*)
119            IF(lwp) WRITE(numout,*) 'dia_hth : diagnostics of the thermocline depth'
120            IF(lwp) WRITE(numout,*) '~~~~~~~ '
121            IF(lwp) WRITE(numout,*)
122         ENDIF
123      ENDIF
124
125      IF( l_hth ) THEN
126         !
127         IF( iom_use( 'mlddzt' ) .OR. iom_use( 'mldr0_3' ) .OR. iom_use( 'mldr0_1' ) ) THEN
128            ! initialization
129            ztinv  (:,:) = 0._wp 
130            zdepinv(:,:) = 0._wp 
131            zmaxdzT(:,:) = 0._wp 
132            DO_2D_11_11
133               zztmp = gdepw(ji,jj,mbkt(ji,jj)+1,Kmm) 
134               hth     (ji,jj) = zztmp
135               zabs2   (ji,jj) = zztmp
136               ztm2    (ji,jj) = zztmp
137               zrho10_3(ji,jj) = zztmp
138               zpycn   (ji,jj) = zztmp
139            END_2D
140            IF( nla10 > 1 ) THEN
141               DO_2D_11_11
142                  zztmp = gdepw(ji,jj,mbkt(ji,jj)+1,Kmm) 
143                  zrho0_3(ji,jj) = zztmp
144                  zrho0_1(ji,jj) = zztmp
145               END_2D
146            ENDIF
147     
148            ! Preliminary computation
149            ! computation of zdelr = (dr/dT)(T,S,10m)*(-0.2 degC)
150            DO_2D_11_11
151               IF( tmask(ji,jj,nla10) == 1. ) THEN
152                  zu  =  1779.50 + 11.250 * ts(ji,jj,nla10,jp_tem,Kmm) - 3.80   * ts(ji,jj,nla10,jp_sal,Kmm)  &
153                     &           - 0.0745 * ts(ji,jj,nla10,jp_tem,Kmm) * ts(ji,jj,nla10,jp_tem,Kmm)   &
154                     &           - 0.0100 * ts(ji,jj,nla10,jp_tem,Kmm) * ts(ji,jj,nla10,jp_sal,Kmm)
155                  zv  =  5891.00 + 38.000 * ts(ji,jj,nla10,jp_tem,Kmm) + 3.00   * ts(ji,jj,nla10,jp_sal,Kmm)  &
156                     &           - 0.3750 * ts(ji,jj,nla10,jp_tem,Kmm) * ts(ji,jj,nla10,jp_tem,Kmm)
157                  zut =    11.25 -  0.149 * ts(ji,jj,nla10,jp_tem,Kmm) - 0.01   * ts(ji,jj,nla10,jp_sal,Kmm)
158                  zvt =    38.00 -  0.750 * ts(ji,jj,nla10,jp_tem,Kmm)
159                  zw  = (zu + 0.698*zv) * (zu + 0.698*zv)
160                  zdelr(ji,jj) = ztem2 * (1000.*(zut*zv - zvt*zu)/zw)
161               ELSE
162                  zdelr(ji,jj) = 0._wp
163               ENDIF
164            END_2D
165
166            ! ------------------------------------------------------------- !
167            ! thermocline depth: strongest vertical gradient of temperature !
168            ! turbocline depth (mixing layer depth): avt = zavt5            !
169            ! MLD: rho = rho(1) + zrho3                                     !
170            ! MLD: rho = rho(1) + zrho1                                     !
171            ! ------------------------------------------------------------- !
172            DO_3DS_11_11( jpkm1, 2, -1 )
173               !
174               zzdep = gdepw(ji,jj,jk,Kmm)
175               zztmp = ( ts(ji,jj,jk-1,jp_tem,Kmm) - ts(ji,jj,jk,jp_tem,Kmm) ) &
176                      & / zzdep * tmask(ji,jj,jk)   ! vertical gradient of temperature (dT/dz)
177               zzdep = zzdep * tmask(ji,jj,1)
178
179               IF( zztmp > zmaxdzT(ji,jj) ) THEN                       
180                   zmaxdzT(ji,jj) = zztmp   
181                   hth    (ji,jj) = zzdep                ! max and depth of dT/dz
182               ENDIF
183         
184               IF( nla10 > 1 ) THEN
185                  zztmp = rhop(ji,jj,jk) - rhop(ji,jj,1)                       ! delta rho(1)
186                  IF( zztmp > zrho3 )   zrho0_3(ji,jj) = zzdep                ! > 0.03
187                  IF( zztmp > zrho1 )   zrho0_1(ji,jj) = zzdep                ! > 0.01
188               ENDIF
189            END_3D
190         
191            CALL iom_put( 'mlddzt', hth )            ! depth of the thermocline
192            IF( nla10 > 1 ) THEN
193               CALL iom_put( 'mldr0_3', zrho0_3 )   ! MLD delta rho(surf) = 0.03
194               CALL iom_put( 'mldr0_1', zrho0_1 )   ! MLD delta rho(surf) = 0.01
195            ENDIF
196            !
197         ENDIF
198         !
199         IF(  iom_use( 'mld_dt02' ) .OR. iom_use( 'topthdep' ) .OR. iom_use( 'mldr10_3' ) .OR.  &   
200            &  iom_use( 'pycndep' ) .OR. iom_use( 'tinv'     ) .OR. iom_use( 'depti'    )  ) THEN
201            ! ------------------------------------------------------------- !
202            ! MLD: abs( tn - tn(10m) ) = ztem2                              !
203            ! Top of thermocline: tn = tn(10m) - ztem2                      !
204            ! MLD: rho = rho10m + zrho3                                     !
205            ! pycnocline: rho = rho10m + (dr/dT)(T,S,10m)*(-0.2 degC)       !
206            ! temperature inversion: max( 0, max of tn - tn(10m) )          !
207            ! depth of temperature inversion                                !
208            ! ------------------------------------------------------------- !
209            DO_3DS_11_11( jpkm1, nlb10, -1 )
210               !
211               zzdep = gdepw(ji,jj,jk,Kmm) * tmask(ji,jj,1)
212               !
213               zztmp = ts(ji,jj,nla10,jp_tem,Kmm) - ts(ji,jj,jk,jp_tem,Kmm)  ! - delta T(10m)
214               IF( ABS(zztmp) > ztem2 )      zabs2   (ji,jj) = zzdep   ! abs > 0.2
215               IF(     zztmp  > ztem2 )      ztm2    (ji,jj) = zzdep   ! > 0.2
216               zztmp = -zztmp                                          ! delta T(10m)
217               IF( zztmp >  ztinv(ji,jj) ) THEN                        ! temperature inversion
218                  ztinv(ji,jj) = zztmp   
219                  zdepinv (ji,jj) = zzdep   ! max value and depth
220               ENDIF
221
222               zztmp = rhop(ji,jj,jk) - rhop(ji,jj,nla10)              ! delta rho(10m)
223               IF( zztmp > zrho3        )    zrho10_3(ji,jj) = zzdep   ! > 0.03
224               IF( zztmp > zdelr(ji,jj) )    zpycn   (ji,jj) = zzdep   ! > equi. delta T(10m) - 0.2
225               !
226            END_3D
227
228            CALL iom_put( 'mld_dt02', zabs2    )   ! MLD abs(delta t) - 0.2
229            CALL iom_put( 'topthdep', ztm2     )   ! T(10) - 0.2
230            CALL iom_put( 'mldr10_3', zrho10_3 )   ! MLD delta rho(10m) = 0.03
231            CALL iom_put( 'pycndep' , zpycn    )   ! MLD delta rho equi. delta T(10m) = 0.2
232            CALL iom_put( 'tinv'    , ztinv    )   ! max. temp. inv. (t10 ref)
233            CALL iom_put( 'depti'   , zdepinv  )   ! depth of max. temp. inv. (t10 ref)
234            !
235         ENDIF
236 
237         ! ------------------------------- !
238         !  Depth of 20C/26C/28C isotherm  !
239         ! ------------------------------- !
240         IF( iom_use ('20d') ) THEN  ! depth of the 20 isotherm
241            ztem2 = 20.
242            CALL dia_hth_dep( Kmm, ztem2, hd20 ) 
243            CALL iom_put( '20d', hd20 )   
244         ENDIF
245         !
246         IF( iom_use ('26d') ) THEN  ! depth of the 26 isotherm
247            ztem2 = 26.
248            CALL dia_hth_dep( Kmm, ztem2, hd26 ) 
249            CALL iom_put( '26d', hd26 )   
250         ENDIF
251         !
252         IF( iom_use ('28d') ) THEN  ! depth of the 28 isotherm
253            ztem2 = 28.
254            CALL dia_hth_dep( Kmm, ztem2, hd28 ) 
255            CALL iom_put( '28d', hd28 )   
256         ENDIF
257       
258         ! ----------------------------- !
259         !  Heat content of first 300 m  !
260         ! ----------------------------- !
261         IF( iom_use ('hc300') ) THEN 
262            zzdep = 300.
263            CALL  dia_hth_htc( Kmm, zzdep, ts(:,:,:,jp_tem,Kmm), htc3 )
264            CALL iom_put( 'hc300', rho0_rcp * htc3 )  ! vertically integrated heat content (J/m2)
265         ENDIF
266         !
267         ! ----------------------------- !
268         !  Heat content of first 700 m  !
269         ! ----------------------------- !
270         IF( iom_use ('hc700') ) THEN 
271            zzdep = 700.
272            CALL  dia_hth_htc( Kmm, zzdep, ts(:,:,:,jp_tem,Kmm), htc7 )
273            CALL iom_put( 'hc700', rho0_rcp * htc7 )  ! vertically integrated heat content (J/m2)
274 
275         ENDIF
276         !
277         ! ----------------------------- !
278         !  Heat content of first 2000 m  !
279         ! ----------------------------- !
280         IF( iom_use ('hc2000') ) THEN 
281            zzdep = 2000.
282            CALL  dia_hth_htc( Kmm, zzdep, ts(:,:,:,jp_tem,Kmm), htc20 )
283            CALL iom_put( 'hc2000', rho0_rcp * htc20 )  ! vertically integrated heat content (J/m2) 
284         ENDIF
285         !
286      ENDIF
287
288      !
289      IF( ln_timing )   CALL timing_stop('dia_hth')
290      !
291   END SUBROUTINE dia_hth
292
293   SUBROUTINE dia_hth_dep( Kmm, ptem, pdept )
294      !
295      INTEGER , INTENT(in) :: Kmm      ! ocean time level index
296      REAL(wp), INTENT(in) :: ptem
297      REAL(wp), DIMENSION(jpi,jpj), INTENT(out) :: pdept     
298      !
299      INTEGER  :: ji, jj, jk, iid
300      REAL(wp) :: zztmp, zzdep
301      INTEGER, DIMENSION(jpi,jpj) :: iktem
302     
303      ! --------------------------------------- !
304      ! search deepest level above ptem         !
305      ! --------------------------------------- !
306      iktem(:,:) = 1
307      DO_3D_11_11( 1, jpkm1 )
308         zztmp = ts(ji,jj,jk,jp_tem,Kmm)
309         IF( zztmp >= ptem )   iktem(ji,jj) = jk
310      END_3D
311
312      ! ------------------------------- !
313      !  Depth of ptem isotherm         !
314      ! ------------------------------- !
315      DO_2D_11_11
316         !
317         zzdep = gdepw(ji,jj,mbkt(ji,jj)+1,Kmm)       ! depth of the ocean bottom
318         !
319         iid = iktem(ji,jj)
320         IF( iid /= 1 ) THEN
321             zztmp =     gdept(ji,jj,iid  ,Kmm)   &                     ! linear interpolation
322               &  + (    gdept(ji,jj,iid+1,Kmm) - gdept(ji,jj,iid,Kmm)                       )   &
323               &  * ( 20.*tmask(ji,jj,iid+1) - ts(ji,jj,iid,jp_tem,Kmm)                       )   &
324               &  / ( ts(ji,jj,iid+1,jp_tem,Kmm) - ts(ji,jj,iid,jp_tem,Kmm) + (1.-tmask(ji,jj,1)) )
325            pdept(ji,jj) = MIN( zztmp , zzdep) * tmask(ji,jj,1)       ! bound by the ocean depth
326         ELSE
327            pdept(ji,jj) = 0._wp
328         ENDIF
329      END_2D
330      !
331   END SUBROUTINE dia_hth_dep
332
333
334   SUBROUTINE dia_hth_htc( Kmm, pdep, pt, phtc )
335      !
336      INTEGER , INTENT(in) :: Kmm      ! ocean time level index
337      REAL(wp), INTENT(in) :: pdep     ! depth over the heat content
338      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in) :: pt   
339      REAL(wp), DIMENSION(jpi,jpj), INTENT(inout) :: phtc 
340      !
341      INTEGER  :: ji, jj, jk, ik
342      REAL(wp), DIMENSION(jpi,jpj) :: zthick
343      INTEGER , DIMENSION(jpi,jpj) :: ilevel
344
345
346      ! surface boundary condition
347     
348      IF( .NOT. ln_linssh ) THEN   ;   zthick(:,:) = 0._wp       ;   phtc(:,:) = 0._wp                                   
349      ELSE                         ;   zthick(:,:) = ssh(:,:,Kmm)   ;   phtc(:,:) = pt(:,:,1) * ssh(:,:,Kmm) * tmask(:,:,1)   
350      ENDIF
351      !
352      ilevel(:,:) = 1
353      DO_3D_11_11( 2, jpkm1 )
354         IF( ( gdept(ji,jj,jk,Kmm) < pdep ) .AND. ( tmask(ji,jj,jk) == 1 ) ) THEN
355             ilevel(ji,jj) = jk
356             zthick(ji,jj) = zthick(ji,jj) + e3t(ji,jj,jk,Kmm)
357             phtc  (ji,jj) = phtc  (ji,jj) + e3t(ji,jj,jk,Kmm) * pt(ji,jj,jk)
358         ENDIF
359      END_3D
360      !
361      DO_2D_11_11
362         ik = ilevel(ji,jj)
363         zthick(ji,jj) = pdep - zthick(ji,jj)   !   remaining thickness to reach depht pdep
364         phtc(ji,jj)   = phtc(ji,jj)    &
365            &           + pt (ji,jj,ik+1) * MIN( e3t(ji,jj,ik+1,Kmm), zthick(ji,jj) ) &
366                                                       * tmask(ji,jj,ik+1)
367      END_2D
368      !
369      !
370   END SUBROUTINE dia_hth_htc
371
372   !!======================================================================
373END MODULE diahth
Note: See TracBrowser for help on using the repository browser.