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/UKMO/AddC3SDiags/NEMOGCM/NEMO/OPA_SRC/DIA – NEMO

source: branches/UKMO/AddC3SDiags/NEMOGCM/NEMO/OPA_SRC/DIA/diahth.F90 @ 13063

Last change on this file since 13063 was 13063, checked in by pdavis, 4 years ago

Ported across calculation of salinity and temp over the first 300m

File size: 20.6 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   USE timing          ! preformance summary
26
27   IMPLICIT NONE
28   PRIVATE
29
30   PUBLIC   dia_hth       ! routine called by step.F90
31   PUBLIC   dia_hth_alloc ! routine called by nemogcm.F90
32
33   LOGICAL , PUBLIC, PARAMETER          ::   lk_diahth = .TRUE.    !: thermocline-20d depths flag
34   ! note: following variables should move to local variables once iom_put is always used
35   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hth    !: depth of the max vertical temperature gradient [m]
36   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hd20   !: depth of 20 C isotherm                         [m]
37   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hd28   !: depth of 28 C isotherm                         [m]
38   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   htc3   !: heat content of first 300 m                    [W]
39   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   tav3   !: Average temp of first 300 m                    [degC]
40   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sav3   !: Average salinity of first 300 m                [PSU]
41
42   !! * Substitutions
43#  include "domzgr_substitute.h90"
44   !!----------------------------------------------------------------------
45   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
46   !! $Id$
47   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
48   !!----------------------------------------------------------------------
49CONTAINS
50
51   FUNCTION dia_hth_alloc()
52      !!---------------------------------------------------------------------
53      INTEGER :: dia_hth_alloc
54      !!---------------------------------------------------------------------
55      !
56      ALLOCATE(hth(jpi,jpj), hd20(jpi,jpj), hd28(jpi,jpj), htc3(jpi,jpj), STAT=dia_hth_alloc)
57      !
58      ALLOCATE(tav3(jpi,jpj), sav3(jpi,jpj), STAT=dia_hth_alloc)
59      !
60      IF( lk_mpp           )   CALL mpp_sum ( dia_hth_alloc )
61      IF(dia_hth_alloc /= 0)   CALL ctl_warn('dia_hth_alloc: failed to allocate arrays.')
62      !
63   END FUNCTION dia_hth_alloc
64
65
66   SUBROUTINE dia_hth( kt )
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      !!
88      INTEGER                          ::   ji, jj, jk            ! dummy loop arguments
89      INTEGER                          ::   iid, ilevel           ! temporary integers
90      INTEGER, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   ik20, ik28  ! levels
91      REAL(wp)                         ::   zavt5 = 5.e-4_wp      ! Kz criterion for the turbocline depth
92      REAL(wp)                         ::   zrho3 = 0.03_wp       ! density     criterion for mixed layer depth
93      REAL(wp)                         ::   zrho1 = 0.01_wp       ! density     criterion for mixed layer depth
94      REAL(wp)                         ::   ztem2 = 0.2_wp        ! temperature criterion for mixed layer depth
95      REAL(wp)                         ::   zthick_0, zcoef       ! temporary scalars
96      REAL(wp)                         ::   zztmp, zzdep          ! temporary scalars inside do loop
97      REAL(wp)                         ::   zu, zv, zw, zut, zvt  ! temporary workspace
98      REAL(wp)                         ::   depth_remain          ! temporary float
99      REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   zabs2      ! MLD: abs( tn - tn(10m) ) = ztem2
100      REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   ztm2       ! Top of thermocline: tn = tn(10m) - ztem2     
101      REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   zrho10_3   ! MLD: rho = rho10m + zrho3     
102      REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   zpycn      ! pycnocline: rho = rho10m + (dr/dT)(T,S,10m)*(-0.2 degC)
103      REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   ztinv      ! max of temperature inversion
104      REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   zdepinv    ! depth of temperature inversion
105      REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   zrho0_3    ! MLD rho = rho(surf) = 0.03
106      REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   zrho0_1    ! MLD rho = rho(surf) = 0.01
107      REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   zmaxdzT    ! max of dT/dz
108      REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   zthick     ! vertical integration thickness
109      REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   zdelr      ! delta rho equivalent to deltaT = 0.2
110      !!----------------------------------------------------------------------
111      IF( nn_timing == 1 )   CALL timing_start('dia_hth')
112
113      IF( kt == nit000 ) THEN
114         !                                      ! allocate dia_hth array
115         IF( dia_hth_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'lim_sbc_init : unable to allocate standard arrays' )
116
117         IF(.not. ALLOCATED(ik20))THEN
118            ALLOCATE(ik20(jpi,jpj), ik28(jpi,jpj), &
119               &      zabs2(jpi,jpj),   &
120               &      ztm2(jpi,jpj),    &
121               &      zrho10_3(jpi,jpj),&
122               &      zpycn(jpi,jpj),   &
123               &      ztinv(jpi,jpj),   &
124               &      zdepinv(jpi,jpj), &
125               &      zrho0_3(jpi,jpj), &
126               &      zrho0_1(jpi,jpj), &
127               &      zmaxdzT(jpi,jpj), &
128               &      zthick(jpi,jpj),  &
129               &      zdelr(jpi,jpj), STAT=ji)
130            IF( lk_mpp  )   CALL mpp_sum(ji)
131            IF( ji /= 0 )   CALL ctl_stop( 'STOP', 'dia_hth : unable to allocate standard ocean arrays' )
132         END IF
133
134         IF(lwp) WRITE(numout,*)
135         IF(lwp) WRITE(numout,*) 'dia_hth : diagnostics of the thermocline depth'
136         IF(lwp) WRITE(numout,*) '~~~~~~~ '
137         IF(lwp) WRITE(numout,*)
138         IF(lwp .AND. lflush) CALL flush(numout)
139      ENDIF
140
141      ! initialization
142      ztinv  (:,:) = 0._wp 
143      zdepinv(:,:) = 0._wp 
144      zmaxdzT(:,:) = 0._wp 
145      DO jj = 1, jpj
146         DO ji = 1, jpi
147            zztmp = bathy(ji,jj)
148            hth     (ji,jj) = zztmp
149            zabs2   (ji,jj) = zztmp
150            ztm2    (ji,jj) = zztmp
151            zrho10_3(ji,jj) = zztmp
152            zpycn   (ji,jj) = zztmp
153        END DO
154      END DO
155      IF( nla10 > 1 ) THEN
156         DO jj = 1, jpj
157            DO ji = 1, jpi
158               zztmp = bathy(ji,jj)
159               zrho0_3(ji,jj) = zztmp
160               zrho0_1(ji,jj) = zztmp
161            END DO
162         END DO
163      ENDIF
164     
165      ! Preliminary computation
166      ! computation of zdelr = (dr/dT)(T,S,10m)*(-0.2 degC)
167      DO jj = 1, jpj
168         DO ji = 1, jpi
169            IF( tmask(ji,jj,nla10) == 1. ) THEN
170               zu  =  1779.50 + 11.250 * tsn(ji,jj,nla10,jp_tem) - 3.80   * tsn(ji,jj,nla10,jp_sal)                             &
171                  &                                              - 0.0745 * tsn(ji,jj,nla10,jp_tem) * tsn(ji,jj,nla10,jp_tem)   &
172                  &                                              - 0.0100 * tsn(ji,jj,nla10,jp_tem) * tsn(ji,jj,nla10,jp_sal)
173               zv  =  5891.00 + 38.000 * tsn(ji,jj,nla10,jp_tem) + 3.00   * tsn(ji,jj,nla10,jp_sal)                             &
174                  &                                              - 0.3750 * tsn(ji,jj,nla10,jp_tem) * tsn(ji,jj,nla10,jp_tem)
175               zut =    11.25 -  0.149 * tsn(ji,jj,nla10,jp_tem) - 0.01   * tsn(ji,jj,nla10,jp_sal)
176               zvt =    38.00 -  0.750 * tsn(ji,jj,nla10,jp_tem)
177               zw  = (zu + 0.698*zv) * (zu + 0.698*zv)
178               zdelr(ji,jj) = ztem2 * (1000.*(zut*zv - zvt*zu)/zw)
179            ELSE
180               zdelr(ji,jj) = 0._wp
181            ENDIF
182         END DO
183      END DO
184
185      ! ------------------------------------------------------------- !
186      ! thermocline depth: strongest vertical gradient of temperature !
187      ! turbocline depth (mixing layer depth): avt = zavt5            !
188      ! MLD: rho = rho(1) + zrho3                                     !
189      ! MLD: rho = rho(1) + zrho1                                     !
190      ! ------------------------------------------------------------- !
191      DO jk = jpkm1, 2, -1   ! loop from bottom to 2
192         DO jj = 1, jpj
193            DO ji = 1, jpi
194               !
195               zzdep = fsdepw(ji,jj,jk)
196               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)
197               zzdep = zzdep * tmask(ji,jj,1)
198
199               IF( zztmp > zmaxdzT(ji,jj) ) THEN                       
200                  zmaxdzT(ji,jj) = zztmp   ;   hth    (ji,jj) = zzdep                ! max and depth of dT/dz
201               ENDIF
202               
203               IF( nla10 > 1 ) THEN
204                  zztmp = rhop(ji,jj,jk) - rhop(ji,jj,1)                             ! delta rho(1)
205                  IF( zztmp > zrho3 )          zrho0_3(ji,jj) = zzdep                ! > 0.03
206                  IF( zztmp > zrho1 )          zrho0_1(ji,jj) = zzdep                ! > 0.01
207               ENDIF
208
209            END DO
210         END DO
211      END DO
212     
213      CALL iom_put( "mlddzt", hth )            ! depth of the thermocline
214      IF( nla10 > 1 ) THEN
215         CALL iom_put( "mldr0_3", zrho0_3 )   ! MLD delta rho(surf) = 0.03
216         CALL iom_put( "mldr0_1", zrho0_1 )   ! MLD delta rho(surf) = 0.01
217      ENDIF
218
219      ! ------------------------------------------------------------- !
220      ! MLD: abs( tn - tn(10m) ) = ztem2                              !
221      ! Top of thermocline: tn = tn(10m) - ztem2                      !
222      ! MLD: rho = rho10m + zrho3                                     !
223      ! pycnocline: rho = rho10m + (dr/dT)(T,S,10m)*(-0.2 degC)       !
224      ! temperature inversion: max( 0, max of tn - tn(10m) )          !
225      ! depth of temperature inversion                                !
226      ! ------------------------------------------------------------- !
227      DO jk = jpkm1, nlb10, -1   ! loop from bottom to nlb10
228         DO jj = 1, jpj
229            DO ji = 1, jpi
230               !
231               zzdep = fsdepw(ji,jj,jk) * tmask(ji,jj,1)
232               !
233               zztmp = tsn(ji,jj,nla10,jp_tem) - tsn(ji,jj,jk,jp_tem)  ! - delta T(10m)
234               IF( ABS(zztmp) > ztem2 )      zabs2   (ji,jj) = zzdep   ! abs > 0.2
235               IF(     zztmp  > ztem2 )      ztm2    (ji,jj) = zzdep   ! > 0.2
236               zztmp = -zztmp                                          ! delta T(10m)
237               IF( zztmp >  ztinv(ji,jj) ) THEN                        ! temperature inversion
238                  ztinv(ji,jj) = zztmp   ;   zdepinv (ji,jj) = zzdep   ! max value and depth
239               ENDIF
240
241               zztmp = rhop(ji,jj,jk) - rhop(ji,jj,nla10)              ! delta rho(10m)
242               IF( zztmp > zrho3        )    zrho10_3(ji,jj) = zzdep   ! > 0.03
243               IF( zztmp > zdelr(ji,jj) )    zpycn   (ji,jj) = zzdep   ! > equi. delta T(10m) - 0.2
244               !
245            END DO
246         END DO
247      END DO
248
249      CALL iom_put( "mld_dt02", zabs2        )   ! MLD abs(delta t) - 0.2
250      CALL iom_put( "topthdep", ztm2         )   ! T(10) - 0.2
251      CALL iom_put( "mldr10_3", zrho10_3     )   ! MLD delta rho(10m) = 0.03
252      CALL iom_put( "pycndep" , zpycn        )   ! MLD delta rho equi. delta T(10m) = 0.2
253      CALL iom_put( "tinv"    , ztinv        )   ! max. temp. inv. (t10 ref)
254      CALL iom_put( "depti"   , zdepinv      )   ! depth of max. temp. inv. (t10 ref)
255
256
257      ! ----------------------------------- !
258      ! search deepest level above 20C/28C  !
259      ! ----------------------------------- !
260      ik20(:,:) = 1
261      ik28(:,:) = 1
262      DO jk = 1, jpkm1   ! beware temperature is not always decreasing with depth => loop from top to bottom
263         DO jj = 1, jpj
264            DO ji = 1, jpi
265               zztmp = tsn(ji,jj,jk,jp_tem)
266               IF( zztmp >= 20. )   ik20(ji,jj) = jk
267               IF( zztmp >= 28. )   ik28(ji,jj) = jk
268            END DO
269         END DO
270      END DO
271
272      ! --------------------------- !
273      !  Depth of 20C/28C isotherm  !
274      ! --------------------------- !
275      DO jj = 1, jpj
276         DO ji = 1, jpi
277            !
278            zzdep = fsdepw(ji,jj,mbkt(ji,jj)+1)       ! depth of the oean bottom
279            !
280            iid = ik20(ji,jj)
281            IF( iid /= 1 ) THEN
282               zztmp =      fsdept(ji,jj,iid  )   &                     ! linear interpolation
283                  &  + (    fsdept(ji,jj,iid+1) - fsdept(ji,jj,iid)                       )   &
284                  &  * ( 20.*tmask(ji,jj,iid+1) - tsn(ji,jj,iid,jp_tem)                       )   &
285                  &  / ( tsn(ji,jj,iid+1,jp_tem) - tsn(ji,jj,iid,jp_tem) + (1.-tmask(ji,jj,1)) )
286               hd20(ji,jj) = MIN( zztmp , zzdep) * tmask(ji,jj,1)       ! bound by the ocean depth
287            ELSE
288               hd20(ji,jj) = 0._wp
289            ENDIF
290            !
291            iid = ik28(ji,jj)
292            IF( iid /= 1 ) THEN
293               zztmp =      fsdept(ji,jj,iid  )   &                     ! linear interpolation
294                  &  + (    fsdept(ji,jj,iid+1) - fsdept(ji,jj,iid)                       )   &
295                  &  * ( 28.*tmask(ji,jj,iid+1) -    tsn(ji,jj,iid,jp_tem)                       )   &
296                  &  / (  tsn(ji,jj,iid+1,jp_tem) -    tsn(ji,jj,iid,jp_tem) + (1.-tmask(ji,jj,1)) )
297               hd28(ji,jj) = MIN( zztmp , zzdep ) * tmask(ji,jj,1)      ! bound by the ocean depth
298            ELSE
299               hd28(ji,jj) = 0._wp
300            ENDIF
301
302         END DO
303      END DO
304      CALL iom_put( "20d", hd20 )   ! depth of the 20 isotherm
305      CALL iom_put( "28d", hd28 )   ! depth of the 28 isotherm
306
307      ! ---------------------------------------- !
308      !  Average temp & salinity of first 300 m  !
309      !  Copied from the hc300 section below     !
310      ! ---------------------------------------- !
311      ! find ilevel with (ilevel+1) the deepest W-level above 300m (we assume we can use e3t_0 to do this search...)
312      ilevel   = 0
313      zthick_0 = 0._wp
314      DO jk = 1, jpkm1                     
315         zthick_0 = zthick_0 + e3t_1d(jk)
316         IF( zthick_0 < 300. )   ilevel = jk
317      END DO
318      ! surface boundary condition
319      IF( lk_vvl ) THEN   ;   zthick(:,:) = 0._wp       ;   tav3(:,:) = 0._wp  ; sav3(:,:) = 0._wp 
320      ELSE
321        ! For the surface layer calculate the integrated temp and salinity
322        ! The SSH "layer" is considered to be at the same temp/salinity as the
323        ! first ocean level.
324        zthick(:,:) = sshn(:,:) * tmask(:,:,1)
325        tav3(:,:) = tsn(:,:,1,jp_tem) * tmask(:,:,1) * zthick(:,:) 
326        sav3(:,:) = tsn(:,:,1,jp_sal) * tmask(:,:,1) * zthick(:,:) 
327      ENDIF
328
329      ! integration down to ilevel
330      ! For each level above 300m (1 to ilevel) calculate the integrated
331      ! temp and salinty (T|S * layer thickness), and sum the thickness
332      DO jk = 1, ilevel
333         zthick(:,:) = zthick(:,:) + fse3t(:,:,jk) * tmask(:,:,jk)
334         tav3  (:,:) = tav3  (:,:) + tsn(:,:,jk,jp_tem) * tmask(:,:,jk) * fse3t(:,:,jk)
335         sav3  (:,:) = sav3  (:,:) + tsn(:,:,jk,jp_sal) * tmask(:,:,jk) * fse3t(:,:,jk)
336      END DO
337
338      DO jj = 1, jpj
339         DO ji = 1, jpi
340            ! What distance is there remaining to get to 300m
341            depth_remain = 300.0 - zthick(ji,jj)
342           
343            ! add the value of reamining thickness * T|S of the level beyond
344            ! 300m.
345            tav3(ji,jj) = tav3(ji,jj) + tsn(ji,jj,ilevel+1,jp_tem) * MIN( fse3t(ji,jj,ilevel+1), depth_remain )  &
346                                                                   * tmask(ji,jj,ilevel+1)
347            sav3(ji,jj) = sav3(ji,jj) + tsn(ji,jj,ilevel+1,jp_sal) * MIN( fse3t(ji,jj,ilevel+1), depth_remain )  &
348                                                                   * tmask(ji,jj,ilevel+1)
349            ! Add the remaining depth, this should take us to 300m (in most
350            ! places)
351            zthick(ji,jj) = zthick(ji,jj) + MIN( fse3t(ji,jj,ilevel+1), depth_remain ) * tmask(ji,jj,ilevel+1)
352            ! Divide by the depth to get the integrated temp|salinity
353            tav3(ji,jj) = tav3(ji,jj) / zthick(ji,jj)
354            sav3(ji,jj) = sav3(ji,jj) / zthick(ji,jj)
355
356         END DO
357      END DO
358
359      CALL iom_put( "tav300", tav3 )      ! first 300m average temperature
360      CALL iom_put( "sav300", sav3 )      ! first 300m average salinity
361
362      ! ----------------------------- !
363      !  Heat content of first 300 m  !
364      ! ----------------------------- !
365
366      ! find ilevel with (ilevel+1) the deepest W-level above 300m (we assume we can use e3t_1d to do this search...)
367      ilevel   = 0
368      zthick_0 = 0._wp
369      DO jk = 1, jpkm1                     
370         zthick_0 = zthick_0 + e3t_1d(jk)
371         IF( zthick_0 < 300. )   ilevel = jk
372      END DO
373      ! surface boundary condition
374      IF( lk_vvl ) THEN   ;   zthick(:,:) = 0._wp       ;   htc3(:,:) = 0._wp                                   
375      ELSE                ;   zthick(:,:) = sshn(:,:)   ;   htc3(:,:) = tsn(:,:,1,jp_tem) * sshn(:,:) * tmask(:,:,1)   
376      ENDIF
377      ! integration down to ilevel
378      DO jk = 1, ilevel
379         zthick(:,:) = zthick(:,:) + fse3t(:,:,jk)
380         htc3  (:,:) = htc3  (:,:) + fse3t(:,:,jk) * tsn(:,:,jk,jp_tem) * tmask(:,:,jk)
381      END DO
382      ! deepest layer
383      zthick(:,:) = 300. - zthick(:,:)   !   remaining thickness to reach 300m
384      DO jj = 1, jpj
385         DO ji = 1, jpi
386            htc3(ji,jj) = htc3(ji,jj) + tsn(ji,jj,ilevel+1,jp_tem) * MIN( fse3t(ji,jj,ilevel+1), zthick(ji,jj) )  &
387                                                                   * tmask(ji,jj,ilevel+1)
388         END DO
389      END DO
390      ! from temperature to heat contain
391      zcoef = rau0 * rcp
392      htc3(:,:) = zcoef * htc3(:,:)
393      CALL iom_put( "hc300", htc3 )      ! first 300m heat content
394      !
395      IF( nn_timing == 1 )   CALL timing_stop('dia_hth')
396      !
397   END SUBROUTINE dia_hth
398
399#else
400   !!----------------------------------------------------------------------
401   !!   Default option :                                       Empty module
402   !!----------------------------------------------------------------------
403   LOGICAL , PUBLIC, PARAMETER ::   lk_diahth = .FALSE.  !: thermocline-20d depths flag
404CONTAINS
405   SUBROUTINE dia_hth( kt )         ! Empty routine
406   IMPLICIT NONE
407    INTEGER, INTENT( in ) ::   kt      ! ocean time-step index
408    WRITE(*,*) 'dia_hth: You should not have seen this print! error?', kt
409   END SUBROUTINE dia_hth
410#endif
411
412   !!======================================================================
413END MODULE diahth
Note: See TracBrowser for help on using the repository browser.