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 branches/2011/dev_LOCEAN_CMCC_INGV_MERCATOR_2011/NEMOGCM/NEMO/OPA_SRC/DIA – NEMO

source: branches/2011/dev_LOCEAN_CMCC_INGV_MERCATOR_2011/NEMOGCM/NEMO/OPA_SRC/DIA/diahth.F90 @ 3103

Last change on this file since 3103 was 2977, checked in by cetlod, 13 years ago

Add in branch 2011/dev_LOCEAN_2011 changes from 2011/dev_r2787_PISCES_improvment, 2011/dev_r2787_LOCEAN_offline_fldread and 2011/dev_r2787_LOCEAN3_TRA_TRP branches, see ticket #877

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