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.
zdfiwm.F90 in branches/2017/wrk_OMP_test_for_Silvia/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/2017/wrk_OMP_test_for_Silvia/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfiwm.F90 @ 8279

Last change on this file since 8279 was 8279, checked in by mocavero, 7 years ago

Implementation of OMP coarse-grained parallelization on ZDF new package

  • Property svn:keywords set to Id
File size: 26.6 KB
Line 
1MODULE zdfiwm
2   !!========================================================================
3   !!                       ***  MODULE  zdfiwm  ***
4   !! Ocean physics: Internal gravity wave-driven vertical mixing
5   !!========================================================================
6   !! History :  1.0  !  2004-04  (L. Bessieres, G. Madec)  Original code
7   !!             -   !  2006-08  (A. Koch-Larrouy)  Indonesian strait
8   !!            3.3  !  2010-10  (C. Ethe, G. Madec)  reorganisation of initialisation phase
9   !!            3.6  !  2016-03  (C. de Lavergne)  New param: internal wave-driven mixing
10   !!            4.0  !  2017-04  (G. Madec)  renamed module, remove the old param. and the CPP keys
11   !!----------------------------------------------------------------------
12
13   !!----------------------------------------------------------------------
14   !!   zdf_iwm       : global     momentum & tracer Kz with wave induced Kz
15   !!   zdf_iwm_init  : global     momentum & tracer Kz with wave induced Kz
16   !!----------------------------------------------------------------------
17   USE oce            ! ocean dynamics and tracers variables
18   USE dom_oce        ! ocean space and time domain variables
19   USE zdf_oce        ! ocean vertical physics variables
20   USE zdfddm         ! ocean vertical physics: double diffusive mixing
21   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
22   USE eosbn2         ! ocean equation of state
23   USE phycst         ! physical constants
24   USE prtctl         ! Print control
25   USE in_out_manager ! I/O manager
26   USE iom            ! I/O Manager
27   USE lib_mpp        ! MPP library
28   USE wrk_nemo       ! work arrays
29   USE timing         ! Timing
30   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
31
32   IMPLICIT NONE
33   PRIVATE
34
35   PUBLIC   zdf_iwm         ! called in step module
36   PUBLIC   zdf_iwm_init    ! called in nemogcm module
37
38   !                       !!* Namelist  namzdf_iwm : internal wave-driven mixing *
39   INTEGER  ::  nn_zpyc     ! pycnocline-intensified mixing energy proportional to N (=1) or N^2 (=2)
40   LOGICAL  ::  ln_mevar    ! variable (=T) or constant (=F) mixing efficiency
41   LOGICAL  ::  ln_tsdiff   ! account for differential T/S wave-driven mixing (=T) or not (=F)
42
43   REAL(wp) ::  r1_6 = 1._wp / 6._wp
44
45   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   ebot_iwm     ! power available from high-mode wave breaking (W/m2)
46   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   epyc_iwm     ! power available from low-mode, pycnocline-intensified wave breaking (W/m2)
47   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   ecri_iwm     ! power available from low-mode, critical slope wave breaking (W/m2)
48   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   hbot_iwm     ! WKB decay scale for high-mode energy dissipation (m)
49   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   hcri_iwm     ! decay scale for low-mode critical slope dissipation (m)
50   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   emix_iwm     ! local energy density available for mixing (W/kg)
51   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   bflx_iwm     ! buoyancy flux Kz * N^2 (W/kg)
52   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   pcmap_iwm    ! vertically integrated buoyancy flux (W/m2)
53   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   zav_ratio    ! S/T diffusivity ratio (only for ln_tsdiff=T)
54   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   zav_wave     ! Internal wave-induced diffusivity
55
56   !! * Substitutions
57#  include "domain_substitute.h90"   
58   !!----------------------------------------------------------------------
59   !! NEMO/OPA 4.0 , NEMO Consortium (2016)
60   !! $Id$
61   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
62   !!----------------------------------------------------------------------
63CONTAINS
64
65   INTEGER FUNCTION zdf_iwm_alloc()
66      !!----------------------------------------------------------------------
67      !!                ***  FUNCTION zdf_iwm_alloc  ***
68      !!----------------------------------------------------------------------
69      ALLOCATE(     ebot_iwm(jpi,jpj),  epyc_iwm(jpi,jpj),  ecri_iwm(jpi,jpj)    ,   &
70      &             hbot_iwm(jpi,jpj),  hcri_iwm(jpi,jpj),  emix_iwm(jpi,jpj,jpk),   &
71      &         bflx_iwm(jpi,jpj,jpk), pcmap_iwm(jpi,jpj), zav_ratio(jpi,jpj,jpk),   & 
72      &         zav_wave(jpi,jpj,jpk), STAT=zdf_iwm_alloc     )
73      !
74      IF( lk_mpp             )   CALL mpp_sum ( zdf_iwm_alloc )
75      IF( zdf_iwm_alloc /= 0 )   CALL ctl_warn('zdf_iwm_alloc: failed to allocate arrays')
76   END FUNCTION zdf_iwm_alloc
77
78
79   SUBROUTINE zdf_iwm( ARG_2D, kt, p_avm, p_avt, p_avs )
80      !!----------------------------------------------------------------------
81      !!                  ***  ROUTINE zdf_iwm  ***
82      !!                   
83      !! ** Purpose :   add to the vertical mixing coefficients the effect of
84      !!              breaking internal waves.
85      !!
86      !! ** Method  : - internal wave-driven vertical mixing is given by:
87      !!                  Kz_wave = min(  100 cm2/s, f(  Reb = emix_iwm /( Nu * N^2 )  )
88      !!              where emix_iwm is the 3D space distribution of the wave-breaking
89      !!              energy and Nu the molecular kinematic viscosity.
90      !!              The function f(Reb) is linear (constant mixing efficiency)
91      !!              if the namelist parameter ln_mevar = F and nonlinear if ln_mevar = T.
92      !!
93      !!              - Compute emix_iwm, the 3D power density that allows to compute
94      !!              Reb and therefrom the wave-induced vertical diffusivity.
95      !!              This is divided into three components:
96      !!                 1. Bottom-intensified low-mode dissipation at critical slopes
97      !!                     emix_iwm(z) = ( ecri_iwm / rau0 ) * EXP( -(H-z)/hcri_iwm )
98      !!                                   / ( 1. - EXP( - H/hcri_iwm ) ) * hcri_iwm
99      !!              where hcri_iwm is the characteristic length scale of the bottom
100      !!              intensification, ecri_iwm a map of available power, and H the ocean depth.
101      !!                 2. Pycnocline-intensified low-mode dissipation
102      !!                     emix_iwm(z) = ( epyc_iwm / rau0 ) * ( sqrt(rn2(z))^nn_zpyc )
103      !!                                   / SUM( sqrt(rn2(z))^nn_zpyc * e3w(z) )
104      !!              where epyc_iwm is a map of available power, and nn_zpyc
105      !!              is the chosen stratification-dependence of the internal wave
106      !!              energy dissipation.
107      !!                 3. WKB-height dependent high mode dissipation
108      !!                     emix_iwm(z) = ( ebot_iwm / rau0 ) * rn2(z) * EXP(-z_wkb(z)/hbot_iwm)
109      !!                                   / SUM( rn2(z) * EXP(-z_wkb(z)/hbot_iwm) * e3w(z) )
110      !!              where hbot_iwm is the characteristic length scale of the WKB bottom
111      !!              intensification, ebot_iwm is a map of available power, and z_wkb is the
112      !!              WKB-stretched height above bottom defined as
113      !!                    z_wkb(z) = H * SUM( sqrt(rn2(z'>=z)) * e3w(z'>=z) )
114      !!                                 / SUM( sqrt(rn2(z'))    * e3w(z')    )
115      !!
116      !!              - update the model vertical eddy viscosity and diffusivity:
117      !!                     avt  = avt  +    av_wave
118      !!                     avm  = avm  +    av_wave
119      !!
120      !!              - if namelist parameter ln_tsdiff = T, account for differential mixing:
121      !!                     avs  = avt  +    av_wave * diffusivity_ratio(Reb)
122      !!
123      !! ** Action  : - Define emix_iwm used to compute internal wave-induced mixing
124      !!              - avt, avs, avm, increased by internal wave-driven mixing   
125      !!
126      !! References :  de Lavergne et al. 2015, JPO; 2016, in prep.
127      !!----------------------------------------------------------------------
128      INTEGER                    , INTENT(in   ) ::   ARG_2D         ! inner domain start-end i-indices
129      INTEGER                    , INTENT(in   ) ::   kt             ! ocean time step
130      REAL(wp), DIMENSION(:,:,:) , INTENT(inout) ::   p_avm          ! momentum Kz (w-points)
131      REAL(wp), DIMENSION(:,:,:) , INTENT(inout) ::   p_avt, p_avs   ! tracer   Kz (w-points)
132      !
133      INTEGER  ::   ji, jj, jk   ! dummy loop indices
134      REAL(wp) ::   zztemp       ! scalar workspace
135      REAL(wp), DIMENSION(WRK_2D) ::  zfact     ! Used for vertical structure
136      REAL(wp), DIMENSION(WRK_2D) ::  zhdep     ! Ocean depth
137      REAL(wp), DIMENSION(WRK_3D) ::  zwkb      ! WKB-stretched height above bottom
138      REAL(wp), DIMENSION(WRK_3D) ::  zweight   ! Weight for high mode vertical distribution
139      REAL(wp), DIMENSION(WRK_3D) ::  znu_t     ! Molecular kinematic viscosity (T grid)
140      REAL(wp), DIMENSION(WRK_3D) ::  znu_w     ! Molecular kinematic viscosity (W grid)
141      REAL(wp), DIMENSION(WRK_3D) ::  zReb      ! Turbulence intensity parameter
142      !!----------------------------------------------------------------------
143      !
144      IF( nn_timing == 1 )   CALL timing_start('zdf_iwm')
145      !
146      !                          ! ----------------------------- !
147      !                          !  Internal wave-driven mixing  !  (compute zav_wave)
148      !                          ! ----------------------------- !
149      !                             
150      !                        !* Critical slope mixing: distribute energy over the time-varying ocean depth,
151      !                                                 using an exponential decay from the seafloor.
152      DO jj = tnldj, tnlej
153         DO ji = tnldi, tnlei     ! part independent of the level
154            zhdep(ji,jj) = gdepw_0(ji,jj,mbkt(ji,jj)+1)       ! depth of the ocean
155            zfact(ji,jj) = rau0 * (  1._wp - EXP( -zhdep(ji,jj) / hcri_iwm(ji,jj) )  )
156            IF( zfact(ji,jj) /= 0._wp )   zfact(ji,jj) = ecri_iwm(ji,jj) / zfact(ji,jj)
157         END DO
158      END DO
159      DO jk = 2, jpkm1              ! complete with the level-dependent part
160         DO jj = tnldj, tnlej
161            DO ji = tnldi, tnlei
162               emix_iwm(ji,jj,jk) = zfact(ji,jj) * wmask(ji,jj,jk)                              & 
163                  &   * (  EXP( ( gde3w_n(ji,jj,jk  ) - zhdep(ji,jj) ) / hcri_iwm(ji,jj) )      &
164                  &      - EXP( ( gde3w_n(ji,jj,jk-1) - zhdep(ji,jj) ) / hcri_iwm(ji,jj) )  )   &
165                  &   / ( gde3w_n(ji,jj,jk) - gde3w_n(ji,jj,jk-1) )
166
167!!gm delta(gde3w_n) = e3t_n  !!  Please verify the grid-point position w versus t-point
168!!gm it seems to me that only 1/hcri_iwm  is used ==>  compute it one for all
169
170            END DO
171         END DO
172      END DO
173
174      !                        !* Pycnocline-intensified mixing: distribute energy over the time-varying
175      !                        !* ocean depth as proportional to sqrt(rn2)^nn_zpyc
176      !                                          ! (NB: N2 is masked, so no use of wmask here)
177      SELECT CASE ( nn_zpyc )
178      !
179      CASE ( 1 )               ! Dissipation scales as N (recommended)
180         !
181         zfact(WRK_2D) = 0._wp         ! part independent of the level
182         DO jk = 2, jpkm1
183            zfact(WRK_2D) = zfact(WRK_2D) + e3w_n(WRK_2D,jk) * SQRT(  MAX( 0._wp, rn2(WRK_2D,jk) )  )
184         END DO
185         !
186         WHERE( zfact(WRK_2D) /= 0 )   zfact(WRK_2D) = epyc_iwm(WRK_2D) / ( rau0 * zfact(WRK_2D) )
187!         DO jj = tnldj, tnlej
188!            DO ji = tnldi, tnlei
189!               IF( zfact(ji,jj) /= 0 )   zfact(ji,jj) = epyc_iwm(ji,jj) / ( rau0 * zfact(ji,jj) )
190!            END DO
191!         END DO
192         !                             ! complete with the level-dependent part
193         DO jk = 2, jpkm1
194            emix_iwm(WRK_2D,jk) = emix_iwm(WRK_2D,jk) + zfact(WRK_2D) * SQRT(  MAX( 0._wp, rn2(WRK_2D,jk) )  )
195         END DO
196         !
197      CASE ( 2 )               ! Dissipation scales as N^2
198         !
199         zfact(WRK_2D) = 0._wp
200         DO jk = 2, jpkm1              ! part independent of the level
201            zfact(WRK_2D) = zfact(WRK_2D) + e3w_n(WRK_2D,jk) * MAX( 0._wp, rn2(WRK_2D,jk) )
202         END DO
203         !
204         WHERE( zfact(WRK_2D) /= 0 )   zfact(WRK_2D) = epyc_iwm(WRK_2D) / ( rau0 * zfact(WRK_2D) )
205!         DO jj = tnldj, tnlej
206!            DO ji = tnldi, tnlei
207!               IF( zfact(ji,jj) /= 0 )   zfact(ji,jj) = epyc_iwm(ji,jj) / ( rau0 * zfact(ji,jj) )
208!            END DO
209!         END DO
210         !
211         DO jk = 2, jpkm1              ! complete with the level-dependent part
212            emix_iwm(WRK_2D,jk) = emix_iwm(WRK_2D,jk) + zfact(WRK_2D) * MAX( 0._wp, rn2(WRK_2D,jk) )
213         END DO
214         !
215      END SELECT
216
217      !                        !* WKB-height dependent mixing: distribute energy over the time-varying
218      !                        !* ocean depth as proportional to rn2 * exp(-z_wkb/rn_hbot)
219      !
220      zwkb (WRK_3D) = 0._wp
221      zfact(WRK_2D) = 0._wp
222      DO jk = 2, jpkm1
223!!gm  this is potentially dangerous
224!         zfact(WRK_2D)    = zfact(WRK_2D) + e3w_n(WRK_2D,jk) * SQRT(  MAX( 0._wp, rn2(WRK_2D,jk) )  )
225!         zwkb (WRK_2D,jk) = zfact(WRK_2D)
226         DO jj = tnldj, tnlej
227            DO ji = tnldi, tnlei
228               zztemp = zfact(ji,jj) + e3w_n(ji,jj,jk) * SQRT(  MAX( 0._wp, rn2(ji,jj,jk) )  )
229               zfact(ji,jj)    = zztemp
230               zwkb (ji,jj,jk) = zztemp
231            END DO
232         END DO         
233      END DO
234      !
235      DO jk = 2, jpkm1
236         DO jj = tnldj, tnlej
237            DO ji = tnldi, tnlei
238               IF( zfact(ji,jj) /= 0 )   zwkb(ji,jj,jk) = zhdep(ji,jj) * ( zfact(ji,jj) - zwkb(ji,jj,jk) )   &
239                  &                                     * wmask(ji,jj,jk) / zfact(ji,jj)
240            END DO
241         END DO
242      END DO
243      zwkb(WRK_2D,1) = zhdep(WRK_2D) * wmask(WRK_2D,1)
244      !
245      zweight(WRK_3D) = 0._wp
246      DO jk = 2, jpkm1
247         zweight(WRK_2D,jk) = MAX( 0._wp, rn2(WRK_2D,jk) ) * hbot_iwm(WRK_2D)  &
248            &               * (  EXP( -zwkb(WRK_2D,jk  ) / hbot_iwm(WRK_2D) )  &
249            &                  - EXP( -zwkb(WRK_2D,jk-1) / hbot_iwm(WRK_2D) )  )
250      END DO
251      !
252      zfact(WRK_2D) = 0._wp
253      DO jk = 2, jpkm1              ! part independent of the level
254         zfact(WRK_2D) = zfact(WRK_2D) + zweight(WRK_2D,jk)
255      END DO
256      !
257      WHERE( zfact(WRK_2D) /= 0 )   zfact(WRK_2D) = ebot_iwm(WRK_2D) / ( rau0 * zfact(WRK_2D) )
258!      DO jj = tnldj, tnlej
259!         DO ji = tnldi, tnlei
260!            IF( zfact(ji,jj) /= 0 )   zfact(ji,jj) = ebot_iwm(ji,jj) / ( rau0 * zfact(ji,jj) )
261!         END DO
262!      END DO
263      !
264      DO jk = 2, jpkm1              ! complete with the level-dependent part
265         emix_iwm(WRK_2D,jk) = emix_iwm(WRK_2D,jk) + zweight(WRK_2D,jk) * zfact(WRK_2D) * wmask(WRK_2D,jk)   &
266            &                                / ( gde3w_n(WRK_2D,jk) - gde3w_n(WRK_2D,jk-1) )
267!!gm  use of e3t_n just above?
268      END DO
269      !
270!!gm  this is to be replaced by just a constant value znu=1.e-6 m2/s
271      ! Calculate molecular kinematic viscosity
272      znu_t(WRK_3D) = 1.e-4_wp * (  17.91_wp - 0.53810_wp * tsn(WRK_3D,jp_tem)             &
273         &                        + 0.00694_wp * tsn(WRK_3D,jp_tem) * tsn(WRK_3D,jp_tem)   &
274         &                        + 0.02305_wp * tsn(WRK_3D,jp_sal)  ) * tmask(WRK_3D) * r1_rau0
275      DO jk = 2, jpkm1
276         znu_w(WRK_2D,jk) = 0.5_wp * ( znu_t(WRK_2D,jk-1) + znu_t(WRK_2D,jk) ) * wmask(WRK_2D,jk)
277      END DO
278!!gm end
279      !
280      ! Calculate turbulence intensity parameter Reb
281      DO jk = 2, jpkm1
282         zReb(WRK_2D,jk) = emix_iwm(WRK_2D,jk) / MAX( 1.e-20_wp, znu_w(WRK_2D,jk) * rn2(WRK_2D,jk) )
283      END DO
284      !
285      ! Define internal wave-induced diffusivity
286      DO jk = 2, jpkm1
287         zav_wave(WRK_2D,jk) = znu_w(WRK_2D,jk) * zReb(WRK_2D,jk) * r1_6   ! This corresponds to a constant mixing efficiency of 1/6
288      END DO
289      !
290      IF( ln_mevar ) THEN              ! Variable mixing efficiency case : modify zav_wave in the
291         DO jk = 2, jpkm1              ! energetic (Reb > 480) and buoyancy-controlled (Reb <10.224 ) regimes
292            DO jj = tnldj, tnlej
293               DO ji = tnldi, tnlei
294                  IF( zReb(ji,jj,jk) > 480.00_wp ) THEN
295                     zav_wave(ji,jj,jk) = 3.6515_wp * znu_w(ji,jj,jk) * SQRT( zReb(ji,jj,jk) )
296                  ELSEIF( zReb(ji,jj,jk) < 10.224_wp ) THEN
297                     zav_wave(ji,jj,jk) = 0.052125_wp * znu_w(ji,jj,jk) * zReb(ji,jj,jk) * SQRT( zReb(ji,jj,jk) )
298                  ENDIF
299               END DO
300            END DO
301         END DO
302      ENDIF
303      !
304      DO jk = 2, jpkm1                 ! Bound diffusivity by molecular value and 100 cm2/s
305         zav_wave(WRK_2D,jk) = MIN(  MAX( 1.4e-7_wp, zav_wave(WRK_2D,jk) ), 1.e-2_wp  ) * wmask(WRK_2D,jk)
306      END DO
307      !
308      IF( kt == nit000 ) THEN        !* Control print at first time-step: diagnose the energy consumed by zav_wave
309         zztemp = 0._wp
310!!gm used of glosum 3D....
311         DO jk = 2, jpkm1
312            DO jj = tnldj, tnlej
313               DO ji = tnldi, tnlei
314                  zztemp = zztemp + e3w_n(ji,jj,jk) * e1e2t(ji,jj)   &
315                     &         * MAX( 0._wp, rn2(ji,jj,jk) ) * zav_wave(ji,jj,jk) * wmask(ji,jj,jk) * tmask_i(ji,jj)
316               END DO
317            END DO
318         END DO
319!         IF( lk_mpp )   CALL mpp_sum( zztemp )
320         zztemp = rau0 * zztemp ! Global integral of rauo * Kz * N^2 = power contributing to mixing
321         !
322         IF(lwp) THEN
323            WRITE(numout,*)
324            WRITE(numout,*) 'zdf_iwm : Internal wave-driven mixing (iwm)'
325            WRITE(numout,*) '~~~~~~~ '
326            WRITE(numout,*)
327            WRITE(numout,*) '      Total power consumption by av_wave =  ', zztemp * 1.e-12_wp, 'TW'
328         ENDIF
329      ENDIF
330
331      !                          ! ----------------------- !
332      !                          !   Update  mixing coefs  !                         
333      !                          ! ----------------------- !
334      !     
335      IF( ln_tsdiff ) THEN          !* Option for differential mixing of salinity and temperature
336         DO jk = 2, jpkm1              ! Calculate S/T diffusivity ratio as a function of Reb
337            DO jj = tnldj, tnlej
338               DO ji = tnldi, tnlei
339                  zav_ratio(ji,jj,jk) = ( 0.505_wp + 0.495_wp *                                                                  &
340                      &   TANH(    0.92_wp * (   LOG10(  MAX( 1.e-20_wp, zReb(ji,jj,jk) * 5._wp * r1_6 )  ) - 0.60_wp   )    )   &
341                      &                 ) * wmask(ji,jj,jk)
342               END DO
343            END DO
344         END DO
345         CALL iom_put( "av_ratio", zav_ratio )
346         DO jk = 2, jpkm1           !* update momentum & tracer diffusivity with wave-driven mixing
347            p_avs(WRK_2D,jk) = p_avs(WRK_2D,jk) + zav_wave(WRK_2D,jk) * zav_ratio(WRK_2D,jk)
348            p_avt(WRK_2D,jk) = p_avt(WRK_2D,jk) + zav_wave(WRK_2D,jk)
349            p_avm(WRK_2D,jk) = p_avm(WRK_2D,jk) + zav_wave(WRK_2D,jk)
350         END DO
351         !
352      ELSE                          !* update momentum & tracer diffusivity with wave-driven mixing
353         DO jk = 2, jpkm1
354            p_avs(WRK_2D,jk) = p_avs(WRK_2D,jk) + zav_wave(WRK_2D,jk)
355            p_avt(WRK_2D,jk) = p_avt(WRK_2D,jk) + zav_wave(WRK_2D,jk)
356            p_avm(WRK_2D,jk) = p_avm(WRK_2D,jk) + zav_wave(WRK_2D,jk)
357         END DO
358      ENDIF
359
360      !                             !* output internal wave-driven mixing coefficient
361      CALL iom_put( "av_wave", zav_wave )
362                                    !* output useful diagnostics: N^2, Kz * N^2 (bflx_iwm),
363                                    !  vertical integral of rau0 * Kz * N^2 (pcmap_iwm), energy density (emix_iwm)
364      IF( iom_use("bflx_iwm") .OR. iom_use("pcmap_iwm") ) THEN
365         DO jk = 2, jpkm1
366            bflx_iwm(WRK_2D,jk) = MAX( 0._wp, rn2(WRK_2D,jk) ) * zav_wave(WRK_2D,jk)
367         END DO
368         pcmap_iwm(WRK_2D) = 0._wp
369         DO jk = 2, jpkm1
370            pcmap_iwm(WRK_2D) = pcmap_iwm(WRK_2D) + e3w_n(WRK_2D,jk) * bflx_iwm(WRK_2D,jk)
371         END DO
372         pcmap_iwm(WRK_2D) = rau0 * pcmap_iwm(WRK_2D)
373         CALL iom_put( "bflx_iwm" , bflx_iwm  )
374         CALL iom_put( "pcmap_iwm", pcmap_iwm )
375      ENDIF
376      CALL iom_put( "bn2", rn2 )
377      CALL iom_put( "emix_iwm", emix_iwm )
378     
379      IF(ln_ctl)   CALL prt_ctl(tab3d_1=zav_wave , clinfo1=' iwm - av_wave: ', tab3d_2=avt, clinfo2=' avt: ', ovlap=1, kdim=jpk)
380      !
381      IF( nn_timing == 1 )   CALL timing_stop('zdf_iwm')
382      !
383   END SUBROUTINE zdf_iwm
384
385
386   SUBROUTINE zdf_iwm_init
387      !!----------------------------------------------------------------------
388      !!                  ***  ROUTINE zdf_iwm_init  ***
389      !!                     
390      !! ** Purpose :   Initialization of the wave-driven vertical mixing, reading
391      !!              of input power maps and decay length scales in netcdf files.
392      !!
393      !! ** Method  : - Read the namzdf_iwm namelist and check the parameters
394      !!
395      !!              - Read the input data in NetCDF files :
396      !!              power available from high-mode wave breaking (mixing_power_bot.nc)
397      !!              power available from pycnocline-intensified wave-breaking (mixing_power_pyc.nc)
398      !!              power available from critical slope wave-breaking (mixing_power_cri.nc)
399      !!              WKB decay scale for high-mode wave-breaking (decay_scale_bot.nc)
400      !!              decay scale for critical slope wave-breaking (decay_scale_cri.nc)
401      !!
402      !! ** input   : - Namlist namzdf_iwm
403      !!              - NetCDF files : mixing_power_bot.nc, mixing_power_pyc.nc, mixing_power_cri.nc,
404      !!              decay_scale_bot.nc decay_scale_cri.nc
405      !!
406      !! ** Action  : - Increase by 1 the nstop flag is setting problem encounter
407      !!              - Define ebot_iwm, epyc_iwm, ecri_iwm, hbot_iwm, hcri_iwm
408      !!
409      !! References : de Lavergne et al. JPO, 2015 ; de Lavergne PhD 2016
410      !!              de Lavergne et al. in prep., 2017
411      !!----------------------------------------------------------------------
412      INTEGER  ::   ji, jj, jk   ! dummy loop indices
413      INTEGER  ::   inum         ! local integer
414      INTEGER  ::   ios
415      REAL(wp) ::   zbot, zpyc, zcri   ! local scalars
416      !!
417      NAMELIST/namzdf_iwm_new/ nn_zpyc, ln_mevar, ln_tsdiff
418      !!----------------------------------------------------------------------
419      !
420      IF( nn_timing == 1 )  CALL timing_start('zdf_iwm_init')
421      !
422      REWIND( numnam_ref )              ! Namelist namzdf_iwm in reference namelist : Wave-driven mixing
423      READ  ( numnam_ref, namzdf_iwm_new, IOSTAT = ios, ERR = 901)
424901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_iwm in reference namelist', lwp )
425      !
426      REWIND( numnam_cfg )              ! Namelist namzdf_iwm in configuration namelist : Wave-driven mixing
427      READ  ( numnam_cfg, namzdf_iwm_new, IOSTAT = ios, ERR = 902 )
428902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_iwm in configuration namelist', lwp )
429      IF(lwm) WRITE ( numond, namzdf_iwm_new )
430      !
431      IF(lwp) THEN                  ! Control print
432         WRITE(numout,*)
433         WRITE(numout,*) 'zdf_iwm_init : internal wave-driven mixing'
434         WRITE(numout,*) '~~~~~~~~~~~~'
435         WRITE(numout,*) '   Namelist namzdf_iwm_new : set wave-driven mixing parameters'
436         WRITE(numout,*) '      Pycnocline-intensified diss. scales as N (=1) or N^2 (=2) = ', nn_zpyc
437         WRITE(numout,*) '      Variable (T) or constant (F) mixing efficiency            = ', ln_mevar
438         WRITE(numout,*) '      Differential internal wave-driven mixing (T) or not (F)   = ', ln_tsdiff
439      ENDIF
440     
441      ! The new wave-driven mixing parameterization elevates avt and avm in the interior, and
442      ! ensures that avt remains larger than its molecular value (=1.4e-7). Therefore, avtb should
443      ! be set here to a very small value, and avmb to its (uniform) molecular value (=1.4e-6).
444      avmb(:) = 1.4e-6_wp        ! viscous molecular value
445      avtb(:) = 1.e-10_wp        ! very small diffusive minimum (background avt is specified in zdf_iwm)   
446      avtb_2d(:,:) = 1.e0_wp     ! uniform
447      IF(lwp) THEN                  ! Control print
448         WRITE(numout,*)
449         WRITE(numout,*) '   Force the background value applied to avm & avt in TKE to be everywhere ',   &
450            &               'the viscous molecular value & a very small diffusive value, resp.'
451      ENDIF
452           
453      !                             ! allocate iwm arrays
454      IF( zdf_iwm_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_iwm_init : unable to allocate iwm arrays' )
455      !
456      !                             ! read necessary fields
457      CALL iom_open('mixing_power_bot',inum)       ! energy flux for high-mode wave breaking [W/m2]
458      CALL iom_get  (inum, jpdom_data, 'field', ebot_iwm, 1 ) 
459      CALL iom_close(inum)
460      !
461      CALL iom_open('mixing_power_pyc',inum)       ! energy flux for pynocline-intensified wave breaking [W/m2]
462      CALL iom_get  (inum, jpdom_data, 'field', epyc_iwm, 1 )
463      CALL iom_close(inum)
464      !
465      CALL iom_open('mixing_power_cri',inum)       ! energy flux for critical slope wave breaking [W/m2]
466      CALL iom_get  (inum, jpdom_data, 'field', ecri_iwm, 1 )
467      CALL iom_close(inum)
468      !
469      CALL iom_open('decay_scale_bot',inum)        ! spatially variable decay scale for high-mode wave breaking [m]
470      CALL iom_get  (inum, jpdom_data, 'field', hbot_iwm, 1 )
471      CALL iom_close(inum)
472      !
473      CALL iom_open('decay_scale_cri',inum)        ! spatially variable decay scale for critical slope wave breaking [m]
474      CALL iom_get  (inum, jpdom_data, 'field', hcri_iwm, 1 )
475      CALL iom_close(inum)
476
477      ebot_iwm(:,:) = ebot_iwm(:,:) * ssmask(:,:)
478      epyc_iwm(:,:) = epyc_iwm(:,:) * ssmask(:,:)
479      ecri_iwm(:,:) = ecri_iwm(:,:) * ssmask(:,:)
480
481      ! Set once for all to zero the first and last vertical levels of appropriate variables
482      emix_iwm (:,:, 1 ) = 0._wp
483      emix_iwm (:,:,jpk) = 0._wp
484      zav_ratio(:,:, 1 ) = 0._wp
485      zav_ratio(:,:,jpk) = 0._wp
486      zav_wave (:,:, 1 ) = 0._wp
487      zav_wave (:,:,jpk) = 0._wp
488
489      zbot = glob_sum( e1e2t(:,:) * ebot_iwm(:,:) )
490      zpyc = glob_sum( e1e2t(:,:) * epyc_iwm(:,:) )
491      zcri = glob_sum( e1e2t(:,:) * ecri_iwm(:,:) )
492      IF(lwp) THEN
493         WRITE(numout,*) '      High-mode wave-breaking energy:             ', zbot * 1.e-12_wp, 'TW'
494         WRITE(numout,*) '      Pycnocline-intensifed wave-breaking energy: ', zpyc * 1.e-12_wp, 'TW'
495         WRITE(numout,*) '      Critical slope wave-breaking energy:        ', zcri * 1.e-12_wp, 'TW'
496      ENDIF
497      !
498      IF( nn_timing == 1 )  CALL timing_stop('zdf_iwm_init')
499      !
500   END SUBROUTINE zdf_iwm_init
501
502   !!======================================================================
503END MODULE zdfiwm
Note: See TracBrowser for help on using the repository browser.