source: CONFIG/UNIFORM/v6/IPSLCM6.5.1/SOURCES/NEMO/diahth.F90 @ 6429

Last change on this file since 6429 was 6429, checked in by cetlod, 14 months ago

CM6.5.2 : adding new diagnostics for NEMO

File size: 18.4 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: diahth.F90 15234 2021-09-08 14:07:02Z clem $
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'    )    .OR.  &
115            1  iom_use( 'maxdzT'   )  ) l_hth = .TRUE.
116         !                                      ! allocate dia_hth array
117         IF( l_hth ) THEN
118            IF( dia_hth_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_hth : unable to allocate standard arrays' )
119            IF(lwp) WRITE(numout,*)
120            IF(lwp) WRITE(numout,*) 'dia_hth : diagnostics of the thermocline depth'
121            IF(lwp) WRITE(numout,*) '~~~~~~~ '
122            IF(lwp) WRITE(numout,*)
123         ENDIF
124      ENDIF
125
126      IF( l_hth ) THEN
127         !
128         IF(  iom_use( 'mlddzt' ) .OR. iom_use( 'mldr0_3' ) .OR. iom_use( 'mldr0_1' ) .OR.  &
129            & iom_use( 'maxdzT' ) ) THEN
130
131            ! initialization
132            ztinv  (:,:) = 0._wp 
133            zdepinv(:,:) = 0._wp 
134            zmaxdzT(:,:) = 0._wp 
135            DO_2D( 1, 1, 1, 1 )
136               zztmp = gdepw(ji,jj,mbkt(ji,jj)+1,Kmm) 
137               hth     (ji,jj) = zztmp
138               zabs2   (ji,jj) = zztmp
139               ztm2    (ji,jj) = zztmp
140               zrho10_3(ji,jj) = zztmp
141               zpycn   (ji,jj) = zztmp
142            END_2D
143            IF( nla10 > 1 ) THEN
144               DO_2D( 1, 1, 1, 1 )
145                  zztmp = gdepw(ji,jj,mbkt(ji,jj)+1,Kmm) 
146                  zrho0_3(ji,jj) = zztmp
147                  zrho0_1(ji,jj) = zztmp
148               END_2D
149            ENDIF
150     
151            ! Preliminary computation
152            ! computation of zdelr = (dr/dT)(T,S,10m)*(-0.2 degC)
153            DO_2D( 1, 1, 1, 1 )
154               IF( tmask(ji,jj,nla10) == 1. ) THEN
155                  zu  =  1779.50 + 11.250 * ts(ji,jj,nla10,jp_tem,Kmm) - 3.80   * ts(ji,jj,nla10,jp_sal,Kmm)  &
156                     &           - 0.0745 * ts(ji,jj,nla10,jp_tem,Kmm) * ts(ji,jj,nla10,jp_tem,Kmm)   &
157                     &           - 0.0100 * ts(ji,jj,nla10,jp_tem,Kmm) * ts(ji,jj,nla10,jp_sal,Kmm)
158                  zv  =  5891.00 + 38.000 * ts(ji,jj,nla10,jp_tem,Kmm) + 3.00   * ts(ji,jj,nla10,jp_sal,Kmm)  &
159                     &           - 0.3750 * ts(ji,jj,nla10,jp_tem,Kmm) * ts(ji,jj,nla10,jp_tem,Kmm)
160                  zut =    11.25 -  0.149 * ts(ji,jj,nla10,jp_tem,Kmm) - 0.01   * ts(ji,jj,nla10,jp_sal,Kmm)
161                  zvt =    38.00 -  0.750 * ts(ji,jj,nla10,jp_tem,Kmm)
162                  zw  = (zu + 0.698*zv) * (zu + 0.698*zv)
163                  zdelr(ji,jj) = ztem2 * (1000.*(zut*zv - zvt*zu)/zw)
164               ELSE
165                  zdelr(ji,jj) = 0._wp
166               ENDIF
167            END_2D
168
169            ! ------------------------------------------------------------- !
170            ! thermocline depth: strongest vertical gradient of temperature !
171            ! turbocline depth (mixing layer depth): avt = zavt5            !
172            ! MLD: rho = rho(1) + zrho3                                     !
173            ! MLD: rho = rho(1) + zrho1                                     !
174            ! ------------------------------------------------------------- !
175            DO_3DS( 1, 1, 1, 1, jpkm1, 2, -1 )   ! loop from bottom to 2
176               !
177               zzdep = gdepw(ji,jj,jk,Kmm)
178               zztmp = ( ts(ji,jj,jk-1,jp_tem,Kmm) - ts(ji,jj,jk,jp_tem,Kmm) ) &
179                      & / zzdep * tmask(ji,jj,jk)   ! vertical gradient of temperature (dT/dz)
180               zzdep = zzdep * tmask(ji,jj,1)
181
182               IF( zztmp > zmaxdzT(ji,jj) ) THEN                       
183                   zmaxdzT(ji,jj) = zztmp   
184                   hth    (ji,jj) = zzdep                ! max and depth of dT/dz
185               ENDIF
186         
187               IF( nla10 > 1 ) THEN
188                  zztmp = rhop(ji,jj,jk) - rhop(ji,jj,1)                       ! delta rho(1)
189                  IF( zztmp > zrho3 )   zrho0_3(ji,jj) = zzdep                ! > 0.03
190                  IF( zztmp > zrho1 )   zrho0_1(ji,jj) = zzdep                ! > 0.01
191               ENDIF
192            END_3D
193         
194            CALL iom_put( 'maxdzT', zmaxdzT )        ! max of dT/dz
195            CALL iom_put( 'mlddzt', hth )            ! depth of the thermocline
196
197            IF( nla10 > 1 ) THEN
198               CALL iom_put( 'mldr0_3', zrho0_3 )   ! MLD delta rho(surf) = 0.03
199               CALL iom_put( 'mldr0_1', zrho0_1 )   ! MLD delta rho(surf) = 0.01
200            ENDIF
201            !
202         ENDIF
203         !
204         IF(  iom_use( 'mld_dt02' ) .OR. iom_use( 'topthdep' ) .OR. iom_use( 'mldr10_3' ) .OR.  &   
205            &  iom_use( 'pycndep' ) .OR. iom_use( 'tinv'     ) .OR. iom_use( 'depti'    )  ) THEN
206            ! ------------------------------------------------------------- !
207            ! MLD: abs( tn - tn(10m) ) = ztem2                              !
208            ! Top of thermocline: tn = tn(10m) - ztem2                      !
209            ! MLD: rho = rho10m + zrho3                                     !
210            ! pycnocline: rho = rho10m + (dr/dT)(T,S,10m)*(-0.2 degC)       !
211            ! temperature inversion: max( 0, max of tn - tn(10m) )          !
212            ! depth of temperature inversion                                !
213            ! ------------------------------------------------------------- !
214            DO_3DS( 1, 1, 1, 1, jpkm1, nlb10, -1 )   ! loop from bottom to nlb10
215               !
216               zzdep = gdepw(ji,jj,jk,Kmm) * tmask(ji,jj,1)
217               !
218               zztmp = ts(ji,jj,nla10,jp_tem,Kmm) - ts(ji,jj,jk,jp_tem,Kmm)  ! - delta T(10m)
219               IF( ABS(zztmp) > ztem2 )      zabs2   (ji,jj) = zzdep   ! abs > 0.2
220               IF(     zztmp  > ztem2 )      ztm2    (ji,jj) = zzdep   ! > 0.2
221               zztmp = -zztmp                                          ! delta T(10m)
222               IF( zztmp >  ztinv(ji,jj) ) THEN                        ! temperature inversion
223                  ztinv(ji,jj) = zztmp   
224                  zdepinv (ji,jj) = zzdep   ! max value and depth
225               ENDIF
226
227               zztmp = rhop(ji,jj,jk) - rhop(ji,jj,nla10)              ! delta rho(10m)
228               IF( zztmp > zrho3        )    zrho10_3(ji,jj) = zzdep   ! > 0.03
229               IF( zztmp > zdelr(ji,jj) )    zpycn   (ji,jj) = zzdep   ! > equi. delta T(10m) - 0.2
230               !
231            END_3D
232
233            CALL iom_put( 'mld_dt02', zabs2    )   ! MLD abs(delta t) - 0.2
234            CALL iom_put( 'topthdep', ztm2     )   ! T(10) - 0.2
235            CALL iom_put( 'mldr10_3', zrho10_3 )   ! MLD delta rho(10m) = 0.03
236            CALL iom_put( 'pycndep' , zpycn    )   ! MLD delta rho equi. delta T(10m) = 0.2
237            CALL iom_put( 'tinv'    , ztinv    )   ! max. temp. inv. (t10 ref)
238            CALL iom_put( 'depti'   , zdepinv  )   ! depth of max. temp. inv. (t10 ref)
239            !
240         ENDIF
241 
242         ! ------------------------------- !
243         !  Depth of 20C/26C/28C isotherm  !
244         ! ------------------------------- !
245         IF( iom_use ('20d') ) THEN  ! depth of the 20 isotherm
246            ztem2 = 20.
247            CALL dia_hth_dep( Kmm, ztem2, hd20 ) 
248            CALL iom_put( '20d', hd20 )   
249         ENDIF
250         !
251         IF( iom_use ('26d') ) THEN  ! depth of the 26 isotherm
252            ztem2 = 26.
253            CALL dia_hth_dep( Kmm, ztem2, hd26 ) 
254            CALL iom_put( '26d', hd26 )   
255         ENDIF
256         !
257         IF( iom_use ('28d') ) THEN  ! depth of the 28 isotherm
258            ztem2 = 28.
259            CALL dia_hth_dep( Kmm, ztem2, hd28 ) 
260            CALL iom_put( '28d', hd28 )   
261         ENDIF
262       
263         ! ----------------------------- !
264         !  Heat content of first 300 m  !
265         ! ----------------------------- !
266         IF( iom_use ('hc300') ) THEN 
267            zzdep = 300.
268            CALL  dia_hth_htc( Kmm, zzdep, ts(:,:,:,jp_tem,Kmm), htc3 )
269            CALL iom_put( 'hc300', rho0_rcp * htc3 )  ! vertically integrated heat content (J/m2)
270         ENDIF
271         !
272         ! ----------------------------- !
273         !  Heat content of first 700 m  !
274         ! ----------------------------- !
275         IF( iom_use ('hc700') ) THEN 
276            zzdep = 700.
277            CALL  dia_hth_htc( Kmm, zzdep, ts(:,:,:,jp_tem,Kmm), htc7 )
278            CALL iom_put( 'hc700', rho0_rcp * htc7 )  ! vertically integrated heat content (J/m2)
279 
280         ENDIF
281         !
282         ! ----------------------------- !
283         !  Heat content of first 2000 m  !
284         ! ----------------------------- !
285         IF( iom_use ('hc2000') ) THEN 
286            zzdep = 2000.
287            CALL  dia_hth_htc( Kmm, zzdep, ts(:,:,:,jp_tem,Kmm), htc20 )
288            CALL iom_put( 'hc2000', rho0_rcp * htc20 )  ! vertically integrated heat content (J/m2) 
289         ENDIF
290         !
291      ENDIF
292
293      !
294      IF( ln_timing )   CALL timing_stop('dia_hth')
295      !
296   END SUBROUTINE dia_hth
297
298   SUBROUTINE dia_hth_dep( Kmm, ptem, pdept )
299      !
300      INTEGER , INTENT(in) :: Kmm      ! ocean time level index
301      REAL(wp), INTENT(in) :: ptem
302      REAL(wp), DIMENSION(jpi,jpj), INTENT(out) :: pdept     
303      !
304      INTEGER  :: ji, jj, jk, iid
305      REAL(wp) :: zztmp, zzdep
306      INTEGER, DIMENSION(jpi,jpj) :: iktem
307     
308      ! --------------------------------------- !
309      ! search deepest level above ptem         !
310      ! --------------------------------------- !
311      iktem(:,:) = 1
312      DO_3D( 1, 1, 1, 1, 1, jpkm1 )   ! beware temperature is not always decreasing with depth => loop from top to bottom
313         zztmp = ts(ji,jj,jk,jp_tem,Kmm)
314         IF( zztmp >= ptem )   iktem(ji,jj) = jk
315      END_3D
316
317      ! ------------------------------- !
318      !  Depth of ptem isotherm         !
319      ! ------------------------------- !
320      DO_2D( 1, 1, 1, 1 )
321         !
322         zzdep = gdepw(ji,jj,mbkt(ji,jj)+1,Kmm)       ! depth of the ocean bottom
323         !
324         iid = iktem(ji,jj)
325         IF( iid /= 1 ) THEN
326             zztmp =     gdept(ji,jj,iid  ,Kmm)   &                     ! linear interpolation
327               &  + (    gdept(ji,jj,iid+1,Kmm) - gdept(ji,jj,iid,Kmm)                       )   &
328               &  * ( 20.*tmask(ji,jj,iid+1) - ts(ji,jj,iid,jp_tem,Kmm)                       )   &
329               &  / ( ts(ji,jj,iid+1,jp_tem,Kmm) - ts(ji,jj,iid,jp_tem,Kmm) + (1.-tmask(ji,jj,1)) )
330            pdept(ji,jj) = MIN( zztmp , zzdep) * tmask(ji,jj,1)       ! bound by the ocean depth
331         ELSE
332            pdept(ji,jj) = 0._wp
333         ENDIF
334      END_2D
335      !
336   END SUBROUTINE dia_hth_dep
337
338
339   SUBROUTINE dia_hth_htc( Kmm, pdep, pt, phtc )
340      !
341      INTEGER , INTENT(in) ::   Kmm      ! ocean time level index
342      REAL(wp), INTENT(in) ::   pdep     ! depth over the heat content
343      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in)    ::   pt   
344      REAL(wp), DIMENSION(jpi,jpj),     INTENT(inout) ::   phtc 
345      !
346      INTEGER  ::   ji, jj, jk, ik
347      REAL(wp), DIMENSION(jpi,jpj) ::   zthick
348      INTEGER , DIMENSION(jpi,jpj) ::   ilevel
349
350
351      ! surface boundary condition
352     
353      IF( .NOT. ln_linssh ) THEN   ;   zthick(:,:) = 0._wp          ;   phtc(:,:) = 0._wp                                   
354      ELSE                         ;   zthick(:,:) = ssh(:,:,Kmm)   ;   phtc(:,:) = pt(:,:,1) * ssh(:,:,Kmm) * tmask(:,:,1)   
355      ENDIF
356      !
357      ilevel(:,:) = 1
358      DO_3D( 1, 1, 1, 1, 1, jpkm1 )
359         IF( ( gdepw(ji,jj,jk+1,Kmm) < pdep ) .AND. ( tmask(ji,jj,jk) == 1 ) ) THEN
360             ilevel(ji,jj) = jk+1
361             zthick(ji,jj) = zthick(ji,jj) + e3t(ji,jj,jk,Kmm)
362             phtc  (ji,jj) = phtc  (ji,jj) + e3t(ji,jj,jk,Kmm) * pt(ji,jj,jk)
363         ENDIF
364      END_3D
365      !
366      DO_2D( 1, 1, 1, 1 )
367         ik = ilevel(ji,jj)
368         IF( tmask(ji,jj,ik) == 1 ) THEN
369            zthick(ji,jj) = MIN ( gdepw(ji,jj,ik+1,Kmm), pdep ) - zthick(ji,jj)   ! remaining thickness to reach dephw pdep
370            phtc(ji,jj)   = phtc(ji,jj) + pt(ji,jj,ik) * zthick(ji,jj)
371         ENDIF
372      END_2D
373      !
374   END SUBROUTINE dia_hth_htc
375
376   !!======================================================================
377END MODULE diahth
Note: See TracBrowser for help on using the repository browser.