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/UKMO/NEMO_4.0_mirror_text_diagnostics/src/OCE/DIA – NEMO

source: NEMO/branches/UKMO/NEMO_4.0_mirror_text_diagnostics/src/OCE/DIA/diahth.F90 @ 10986

Last change on this file since 10986 was 10986, checked in by andmirek, 5 years ago

GMED 462 add flush

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