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/dev_r11943_MERGE_2019/src/OCE/DIA – NEMO

source: NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/DIA/diahth.F90 @ 12242

Last change on this file since 12242 was 12193, checked in by davestorkey, 4 years ago

2019/dev_r11943_MERGE_2019: Merge in dev_r12072_TOP-01_ENHANCE-11_cethe

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