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/2020/dev_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final/src/OCE/DIA – NEMO

source: NEMO/branches/2020/dev_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final/src/OCE/DIA/diahth.F90 @ 14219

Last change on this file since 14219 was 14219, checked in by mcastril, 4 years ago

Add Mixed Precision support by Oriol Tintó

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