source: CONFIG/UNIFORM/v6/NEMO_v6.5/SOURCES/zdfiwm.F90 @ 6601

Last change on this file since 6601 was 6601, checked in by cetlod, 9 months ago

NEMOv6.5 : Update config to switch to NEMOv4.2.1 and PISCES gas

File size: 25.1 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   !!            4.0  !  2020-12  (C. de Lavergne)  Update param to match published one
12   !!            4.0  !  2021-09  (C. de Lavergne)  Add energy from trapped and shallow internal tides
13   !!----------------------------------------------------------------------
14
15   !!----------------------------------------------------------------------
16   !!   zdf_iwm       : global     momentum & tracer Kz with wave induced Kz
17   !!   zdf_iwm_init  : global     momentum & tracer Kz with wave induced Kz
18   !!----------------------------------------------------------------------
19   USE oce            ! ocean dynamics and tracers variables
20   USE dom_oce        ! ocean space and time domain variables
21   USE zdf_oce        ! ocean vertical physics variables
22   USE zdfddm         ! ocean vertical physics: double diffusive mixing
23   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
24   USE eosbn2         ! ocean equation of state
25   USE phycst         ! physical constants
26   !
27   USE fldread        ! field read
28   USE prtctl         ! Print control
29   USE in_out_manager ! I/O manager
30   USE iom            ! I/O Manager
31   USE lib_mpp        ! MPP library
32   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
33
34   IMPLICIT NONE
35   PRIVATE
36
37   PUBLIC   zdf_iwm        ! called in step module
38   PUBLIC   zdf_iwm_init   ! called in nemogcm module
39
40   !                      !!* Namelist  namzdf_iwm : internal wave-driven mixing *
41   LOGICAL ::  ln_mevar    ! variable (=T) or constant (=F) mixing efficiency
42   LOGICAL ::  ln_tsdiff   ! account for differential T/S wave-driven mixing (=T) or not (=F)
43
44   REAL(wp)::  r1_6 = 1._wp / 6._wp
45   REAL(wp)::  rnu  = 1.4e-6_wp   ! molecular kinematic viscosity
46
47   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   ebot_iwm   ! bottom-intensified dissipation above abyssal hills (W/m2)
48   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   ecri_iwm   ! bottom-intensified dissipation at topographic slopes (W/m2)
49   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   ensq_iwm   ! dissipation scaling with squared buoyancy frequency (W/m2)
50   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   esho_iwm   ! dissipation due to shoaling internal tides (W/m2)
51   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hbot_iwm   ! decay scale for abyssal hill dissipation (m)
52   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hcri_iwm   ! inverse decay scale for topographic slope dissipation (m-1)
53
54   !! * Substitutions
55#  include "do_loop_substitute.h90"
56#  include "single_precision_substitute.h90"
57#  include "domzgr_substitute.h90"
58   !!----------------------------------------------------------------------
59   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
60   !! $Id: zdfiwm.F90 15533 2021-11-24 12:07:20Z cdllod $
61   !! Software governed by the CeCILL license (see ./LICENSE)
62   !!----------------------------------------------------------------------
63CONTAINS
64
65   INTEGER FUNCTION zdf_iwm_alloc()
66      !!----------------------------------------------------------------------
67      !!                ***  FUNCTION zdf_iwm_alloc  ***
68      !!----------------------------------------------------------------------
69      ALLOCATE( ebot_iwm(jpi,jpj),  ecri_iwm(jpi,jpj),  ensq_iwm(jpi,jpj) ,     &
70      &         esho_iwm(jpi,jpj),  hbot_iwm(jpi,jpj),  hcri_iwm(jpi,jpj) , STAT=zdf_iwm_alloc )
71      !
72      CALL mpp_sum ( 'zdfiwm', zdf_iwm_alloc )
73      IF( zdf_iwm_alloc /= 0 )   CALL ctl_stop( 'STOP', 'zdf_iwm_alloc: failed to allocate arrays' )
74   END FUNCTION zdf_iwm_alloc
75
76
77   SUBROUTINE zdf_iwm( kt, Kmm, p_avm, p_avt, p_avs )
78      !!----------------------------------------------------------------------
79      !!                  ***  ROUTINE zdf_iwm  ***
80      !!                   
81      !! ** Purpose :   add to the vertical mixing coefficients the effect of
82      !!              breaking internal waves.
83      !!
84      !! ** Method  : - internal wave-driven vertical mixing is given by:
85      !!                  Kz_wave = min( f( Reb = zemx_iwm / (Nu * N^2) ), 100 cm2/s  )
86      !!              where zemx_iwm is the 3D space distribution of the wave-breaking
87      !!              energy and Nu the molecular kinematic viscosity.
88      !!              The function f(Reb) is linear (constant mixing efficiency)
89      !!              if the namelist parameter ln_mevar = F and nonlinear if ln_mevar = T.
90      !!
91      !!              - Compute zemx_iwm, the 3D power density that allows to compute
92      !!              Reb and therefrom the wave-induced vertical diffusivity.
93      !!              This is divided into four components:
94      !!                 1. Bottom-intensified dissipation at topographic slopes, expressed
95      !!              as an exponential decay above the bottom.
96      !!                     zemx_iwm(z) = ( ecri_iwm / rho0 ) * EXP( -(H-z)/hcri_iwm )
97      !!                                   / ( 1. - EXP( - H/hcri_iwm ) ) * hcri_iwm
98      !!              where hcri_iwm is the characteristic length scale of the bottom
99      !!              intensification, ecri_iwm a static 2D map of available power, and
100      !!              H the ocean depth.
101      !!                 2. Bottom-intensified dissipation above abyssal hills, expressed
102      !!              as an algebraic decay above bottom.
103      !!                     zemx_iwm(z) = ( ebot_iwm / rho0 ) * ( 1 + hbot_iwm/H )
104      !!                                   / ( 1 + (H-z)/hbot_iwm )^2                               
105      !!              where hbot_iwm is the characteristic length scale of the bottom
106      !!              intensification and ebot_iwm is a static 2D map of available power.
107      !!                 3. Dissipation scaling in the vertical with the squared buoyancy
108      !!              frequency (N^2).
109      !!                     zemx_iwm(z) = ( ensq_iwm / rho0 ) * rn2(z)
110      !!                                   / ZSUM( rn2 * e3w )
111      !!              where ensq_iwm is a static 2D map of available power.
112      !!                 4. Dissipation due to shoaling internal tides, scaling in the
113      !!              vertical with the buoyancy frequency (N).
114      !!                     zemx_iwm(z) = ( esho_iwm / rho0 ) * sqrt(rn2(z))
115      !!                                   / ZSUM( sqrt(rn2) * e3w )
116      !!              where esho_iwm is a static 2D map of available power.
117      !!
118      !!              - update the model vertical eddy viscosity and diffusivity:
119      !!                     avt  = avt  +    av_wave
120      !!                     avs  = avs  +    av_wave
121      !!                     avm  = avm  +    av_wave
122      !!
123      !!              - if namelist parameter ln_tsdiff = T, account for differential mixing:
124      !!                     avs  = avs  +    av_wave * diffusivity_ratio(Reb)
125      !!
126      !! ** Action  : - avt, avs, avm, increased by internal wave-driven mixing   
127      !!
128      !! References :  de Lavergne et al. JAMES 2020, https://doi.org/10.1029/2020MS002065
129      !!               de Lavergne et al. JPO 2016, https://doi.org/10.1175/JPO-D-14-0259.1
130      !!----------------------------------------------------------------------
131      INTEGER                    , INTENT(in   ) ::   kt             ! ocean time step
132      INTEGER                    , INTENT(in   ) ::   Kmm            ! time level index     
133      REAL(wp), DIMENSION(:,:,:) , INTENT(inout) ::   p_avm          ! momentum Kz (w-points)
134      REAL(wp), DIMENSION(:,:,:) , INTENT(inout) ::   p_avt, p_avs   ! tracer   Kz (w-points)
135      !
136      INTEGER  ::   ji, jj, jk   ! dummy loop indices
137      REAL(wp), SAVE :: zztmp
138      !
139      REAL(wp), DIMENSION(A2D(nn_hls))     ::   zfact       ! Used for vertical structure
140      REAL(wp), DIMENSION(A2D(nn_hls),jpk) ::   zReb        ! Turbulence intensity parameter
141      REAL(wp), DIMENSION(A2D(nn_hls),jpk) ::   zemx_iwm    ! local energy density available for mixing (W/kg)
142      REAL(wp), DIMENSION(A2D(nn_hls),jpk) ::   zav_ratio   ! S/T diffusivity ratio (only for ln_tsdiff=T)
143      REAL(wp), DIMENSION(A2D(nn_hls),jpk) ::   zav_wave    ! Internal wave-induced diffusivity
144      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::   z3d  ! 3D workspace used for iom_put
145      REAL(wp), ALLOCATABLE, DIMENSION(:,:)   ::   z2d  ! 2D     -      -    -     -
146      !!----------------------------------------------------------------------
147      !
148      !                       !* Initialize appropriately certain variables
149      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 1, jpk )
150         zav_ratio(ji,jj,jk) = 1._wp * wmask(ji,jj,jk)  ! important to set it to 1 here
151      END_3D
152      IF( iom_use("emix_iwm") )                         zemx_iwm (:,:,:) = 0._wp
153      IF( iom_use("av_wave") .OR. sn_cfctl%l_prtctl )   zav_wave (:,:,:) = 0._wp
154      !
155      !                       ! ----------------------------- !
156      !                       !  Internal wave-driven mixing  !  (compute zav_wave)
157      !                       ! ----------------------------- !
158      !                             
159      !                       !* 'cri' component: distribute energy over the time-varying
160      !                       !* ocean depth using an exponential decay from the seafloor.
161      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )                ! part independent of the level
162         IF( ht(ji,jj) /= 0._wp ) THEN ; zfact(ji,jj) = ecri_iwm(ji,jj) * r1_rho0 / ( 1._wp - EXP( -ht(ji,jj) * hcri_iwm(ji,jj) ) )
163         ELSE                          ; zfact(ji,jj) = 0._wp
164         ENDIF
165      END_2D
166
167      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )      ! complete with the level-dependent part
168         zemx_iwm(ji,jj,jk) = zfact(ji,jj) * (  EXP( ( gdept(ji,jj,jk  ,Kmm) - ht(ji,jj) ) * hcri_iwm(ji,jj) )   &
169            &                                 - EXP( ( gdept(ji,jj,jk-1,Kmm) - ht(ji,jj) ) * hcri_iwm(ji,jj) )   &
170            &                                ) * wmask(ji,jj,jk) / e3w(ji,jj,jk,Kmm)
171      END_3D
172
173                               !* 'bot' component: distribute energy over the time-varying
174                               !* ocean depth using an algebraic decay above the seafloor.
175      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )               ! part independent of the level
176         IF( ht(ji,jj) /= 0._wp ) THEN ; zfact(ji,jj) = ebot_iwm(ji,jj) * ( 1._wp +  hbot_iwm(ji,jj) / ht(ji,jj)  ) * r1_rho0
177         ELSE                          ; zfact(ji,jj) = 0._wp
178         ENDIF
179      END_2D
180
181      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )     ! complete with the level-dependent part
182         zemx_iwm(ji,jj,jk) = zemx_iwm(ji,jj,jk) +                                                                           & 
183            &                 zfact(ji,jj) * (  1._wp / ( 1._wp + ( ht(ji,jj) - gdept(ji,jj,jk  ,Kmm) ) / hbot_iwm(ji,jj) )  &
184            &                                 - 1._wp / ( 1._wp + ( ht(ji,jj) - gdept(ji,jj,jk-1,Kmm) ) / hbot_iwm(ji,jj) )  &
185            &                                ) * wmask(ji,jj,jk) / e3w(ji,jj,jk,Kmm)
186      END_3D
187
188                               !* 'nsq' component: distribute energy over the time-varying
189                               !* ocean depth as proportional to rn2
190      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
191         zfact(ji,jj) = 0._wp
192      END_2D
193      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )     ! part independent of the level
194         zfact(ji,jj) = zfact(ji,jj) + e3w(ji,jj,jk,Kmm) * MAX( 0._wp, rn2(ji,jj,jk) )
195      END_3D
196      !
197      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
198         IF( zfact(ji,jj) /= 0._wp )   zfact(ji,jj) = ensq_iwm(ji,jj) * r1_rho0 / zfact(ji,jj)
199      END_2D
200      !
201      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )     ! complete with the level-dependent part
202         zemx_iwm(ji,jj,jk) = zemx_iwm(ji,jj,jk) + zfact(ji,jj) * MAX( 0._wp, rn2(ji,jj,jk) )
203      END_3D
204
205                               !* 'sho' component: distribute energy over the time-varying
206                               !* ocean depth as proportional to sqrt(rn2)
207      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
208         zfact(ji,jj) = 0._wp
209      END_2D
210      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )     ! part independent of the level
211         zfact(ji,jj) = zfact(ji,jj) + e3w(ji,jj,jk,Kmm) * SQRT( MAX( 0._wp, rn2(ji,jj,jk) ) )
212      END_3D
213      !
214      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
215         IF( zfact(ji,jj) /= 0._wp )   zfact(ji,jj) = esho_iwm(ji,jj) * r1_rho0 / zfact(ji,jj)
216      END_2D
217      !
218      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )     ! complete with the level-dependent part
219         zemx_iwm(ji,jj,jk) = zemx_iwm(ji,jj,jk) + zfact(ji,jj) * SQRT( MAX( 0._wp, rn2(ji,jj,jk) ) )
220      END_3D
221
222      ! Calculate turbulence intensity parameter Reb
223      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
224         zReb(ji,jj,jk) = zemx_iwm(ji,jj,jk) / MAX( 1.e-20_wp, rnu * rn2(ji,jj,jk) )
225      END_3D
226      !
227      ! Define internal wave-induced diffusivity
228      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
229         zav_wave(ji,jj,jk) = zReb(ji,jj,jk) * r1_6 * rnu  ! This corresponds to a constant mixing efficiency of 1/6
230      END_3D
231      !
232      IF( ln_mevar ) THEN                                          ! Variable mixing efficiency case : modify zav_wave in the
233         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 ) ! energetic (Reb > 480) and buoyancy-controlled (Reb <10.224) regimes
234            IF( zReb(ji,jj,jk) > 480.00_wp ) THEN
235               zav_wave(ji,jj,jk) = 3.6515_wp * rnu * SQRT( zReb(ji,jj,jk) )
236            ELSEIF( zReb(ji,jj,jk) < 10.224_wp ) THEN
237               zav_wave(ji,jj,jk) = 0.052125_wp * rnu * zReb(ji,jj,jk) * SQRT( zReb(ji,jj,jk) )
238            ENDIF
239         END_3D
240      ENDIF
241      !
242      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )    ! Bound diffusivity by molecular value and 100 cm2/s
243         zav_wave(ji,jj,jk) = MIN( MAX( 1.4e-7_wp, zav_wave(ji,jj,jk) ), 1.e-2_wp ) * wmask(ji,jj,jk)
244      END_3D
245      !
246      !                          ! ----------------------- !
247      !                          !   Update  mixing coefs  !                         
248      !                          ! ----------------------- !
249      !
250      IF( ln_tsdiff ) THEN                !* Option for differential mixing of salinity and temperature
251         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 ) ! Calculate S/T diffusivity ratio as a function of Reb (else it is set to 1)
252            zav_ratio(ji,jj,jk) = ( 0.505_wp + &
253               &                    0.495_wp * TANH( 0.92_wp * ( LOG10( MAX( 1.e-20, zReb(ji,jj,jk) * 5._wp * r1_6 ) ) - 0.60_wp ) ) &
254               &                  ) * wmask(ji,jj,jk)
255         END_3D
256      ENDIF
257      CALL iom_put( "av_ratio", zav_ratio )
258      !
259      DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 ) !* update momentum & tracer diffusivity with wave-driven mixing
260         p_avs(ji,jj,jk) = p_avs(ji,jj,jk) + zav_wave(ji,jj,jk) * zav_ratio(ji,jj,jk)
261         p_avt(ji,jj,jk) = p_avt(ji,jj,jk) + zav_wave(ji,jj,jk)
262         p_avm(ji,jj,jk) = p_avm(ji,jj,jk) + zav_wave(ji,jj,jk)
263      END_3D
264      !                             !* output internal wave-driven mixing coefficient
265      CALL iom_put( "av_wave", zav_wave )
266                                    !* output useful diagnostics: Kz*N^2 ,
267                                    !  vertical integral of rho0 * Kz * N^2 , energy density (zemx_iwm)
268      IF( iom_use("bflx_iwm") .OR. iom_use("pcmap_iwm") ) THEN
269         ALLOCATE( z2d(A2D(nn_hls)) , z3d(A2D(nn_hls),jpk) )
270         z2d(:,:) = 0._wp ; z3d(:,:,:) = 0._wp   ! Initialisation for iom_put
271         DO_3D( 0, 0, 0, 0, 2, jpkm1 )
272            z3d(ji,jj,jk) = MAX( 0._wp, rn2(ji,jj,jk) ) * zav_wave(ji,jj,jk)
273            z2d(ji,jj) = z2d(ji,jj) + rho0 * e3w(ji,jj,jk,Kmm) * z3d(ji,jj,jk) * wmask(ji,jj,jk)
274         END_3D
275         CALL iom_put(  "bflx_iwm", z3d )
276         CALL iom_put( "pcmap_iwm", z2d )
277         DEALLOCATE( z2d , z3d )
278      ENDIF
279      CALL iom_put( "emix_iwm", zemx_iwm )
280
281      !
282      IF( kt == nit000 ) THEN        !* Control print at first time-step: diagnose the energy consumed by zav_wave
283         IF( .NOT. l_istiled .OR. ntile == 1 ) zztmp = 0._wp                    ! Do only on the first tile
284         DO_3D( 0, 0, 0, 0, 2, jpkm1 )
285            zztmp = zztmp + e3w(ji,jj,jk,Kmm) * e1e2t(ji,jj)   &
286               &          * MAX( 0._wp, rn2(ji,jj,jk) ) * zav_wave(ji,jj,jk) * wmask(ji,jj,jk) * tmask_i(ji,jj)
287         END_3D
288
289         IF( .NOT. l_istiled .OR. ntile == nijtile ) THEN                       ! Do only on the last tile
290            CALL mpp_sum( 'zdfiwm', zztmp )
291            zztmp = rho0 * zztmp ! Global integral of rho0 * Kz * N^2 = power contributing to mixing
292            !
293            IF(lwp) THEN
294               WRITE(numout,*)
295               WRITE(numout,*) 'zdf_iwm : Internal wave-driven mixing (iwm)'
296               WRITE(numout,*) '~~~~~~~ '
297               WRITE(numout,*)
298               WRITE(numout,*) '      Total power consumption by av_wave =  ', zztmp * 1.e-12_wp, 'TW'
299            ENDIF
300         ENDIF
301      ENDIF
302
303      IF(sn_cfctl%l_prtctl)   CALL prt_ctl(tab3d_1=CASTDP(zav_wave) , clinfo1=' iwm - av_wave: ', tab3d_2=CASTDP(avt), clinfo2=' avt: ')
304      !
305   END SUBROUTINE zdf_iwm
306
307
308   SUBROUTINE zdf_iwm_init
309      !!----------------------------------------------------------------------
310      !!                  ***  ROUTINE zdf_iwm_init  ***
311      !!                     
312      !! ** Purpose :   Initialization of the internal wave-driven vertical mixing, reading
313      !!              of input power maps and decay length scales in a netcdf file.
314      !!
315      !! ** Method  : - Read the namzdf_iwm namelist and check the parameters
316      !!
317      !!              - Read the input data in a NetCDF file (zdfiwm_forcing.nc) with variables:
318      !!              'power_bot' bottom-intensified dissipation above abyssal hills
319      !!              'power_cri' bottom-intensified dissipation at topographic slopes
320      !!              'power_nsq' dissipation scaling with squared buoyancy frequency
321      !!              'power_sho' dissipation due to shoaling internal tides
322      !!              'scale_bot' decay scale for abyssal hill dissipation
323      !!              'scale_cri' decay scale for topographic-slope dissipation
324      !!
325      !! ** input   : - Namlist namzdf_iwm
326      !!              - NetCDF file : zdfiwm_forcing.nc
327      !!
328      !! ** Action  : - Increase by 1 the nstop flag is setting problem encounter
329      !!              - Define ebot_iwm, ecri_iwm, ensq_iwm, esho_iwm, hbot_iwm, hcri_iwm
330      !!
331      !! References : de Lavergne et al. JAMES 2020, https://doi.org/10.1029/2020MS002065
332      !!----------------------------------------------------------------------
333      INTEGER  ::   ifpr               ! dummy loop indices
334      INTEGER  ::   inum               ! local integer
335      INTEGER  ::   ios
336      !
337      CHARACTER(len=256)            ::   cn_dir                 ! Root directory for location of ssr files
338      INTEGER, PARAMETER            ::   jpiwm  = 6             ! maximum number of variables to read
339      INTEGER, PARAMETER            ::   jp_mpb = 1
340      INTEGER, PARAMETER            ::   jp_mpc = 2
341      INTEGER, PARAMETER            ::   jp_mpn = 3
342      INTEGER, PARAMETER            ::   jp_mps = 4
343      INTEGER, PARAMETER            ::   jp_dsb = 5
344      INTEGER, PARAMETER            ::   jp_dsc = 6
345      !
346      TYPE(FLD_N), DIMENSION(jpiwm) ::   slf_iwm                        ! array of namelist informations
347      TYPE(FLD_N)                   ::   sn_mpb, sn_mpc, sn_mpn, sn_mps ! information about Mixing Power field to be read
348      TYPE(FLD_N)                   ::   sn_dsb, sn_dsc                 ! information about Decay Scale field to be read
349      TYPE(FLD  ), DIMENSION(jpiwm) ::   sf_iwm                         ! structure of input fields (file informations, fields read)
350      !
351      REAL(wp), DIMENSION(jpi,jpj,4) ::   ztmp
352      REAL(wp), DIMENSION(4)         ::   zdia
353      !
354      NAMELIST/namzdf_iwm/ ln_mevar, ln_tsdiff, &
355          &                cn_dir, sn_mpb, sn_mpc, sn_mpn, sn_mps, sn_dsb, sn_dsc
356      !!----------------------------------------------------------------------
357      !
358      READ  ( numnam_ref, namzdf_iwm, IOSTAT = ios, ERR = 901)
359901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namzdf_iwm in reference namelist' )
360      !
361      READ  ( numnam_cfg, namzdf_iwm, IOSTAT = ios, ERR = 902 )
362902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namzdf_iwm in configuration namelist' )
363      IF(lwm) WRITE ( numond, namzdf_iwm )
364      !
365      IF(lwp) THEN                  ! Control print
366         WRITE(numout,*)
367         WRITE(numout,*) 'zdf_iwm_init : internal wave-driven mixing'
368         WRITE(numout,*) '~~~~~~~~~~~~'
369         WRITE(numout,*) '   Namelist namzdf_iwm : set wave-driven mixing parameters'
370         WRITE(numout,*) '      Variable (T) or constant (F) mixing efficiency            = ', ln_mevar
371         WRITE(numout,*) '      Differential internal wave-driven mixing (T) or not (F)   = ', ln_tsdiff
372      ENDIF
373     
374      ! This internal-wave-driven mixing parameterization elevates avt and avm in the interior, and
375      ! ensures that avt remains larger than its molecular value (=1.4e-7). Therefore, avtb should
376      ! be set here to a very small value, and avmb to its (uniform) molecular value (=1.4e-6).
377      avmb(:) = 1.2e-4_wp        ! molecular value
378      avtb(:) = 1.e-10_wp        ! very small diffusive minimum (background avt is specified in zdf_iwm)   
379      avtb_2d(:,:) = 1._wp       ! uniform
380      IF(lwp) THEN                  ! Control print
381         WRITE(numout,*)
382         WRITE(numout,*) '   Force the background value applied to avm & avt in TKE to be everywhere ',   &
383            &               'the viscous molecular value & a very small diffusive value, resp.'
384      ENDIF
385           
386      !                             ! allocate iwm arrays
387      IF( zdf_iwm_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_iwm_init : unable to allocate iwm arrays' )
388      !
389      ! store namelist information in an array
390      slf_iwm(jp_mpb) = sn_mpb ; slf_iwm(jp_mpc) = sn_mpc ; slf_iwm(jp_mpn) = sn_mpn ; slf_iwm(jp_mps) = sn_mps
391      slf_iwm(jp_dsb) = sn_dsb ; slf_iwm(jp_dsc) = sn_dsc
392      !
393      DO ifpr= 1, jpiwm
394         ALLOCATE( sf_iwm(ifpr)%fnow(jpi,jpj,1)   )
395         IF( slf_iwm(ifpr)%ln_tint ) ALLOCATE( sf_iwm(ifpr)%fdta(jpi,jpj,1,2) )
396      END DO
397
398      ! fill sf_iwm with sf_iwm and control print
399      CALL fld_fill( sf_iwm, slf_iwm , cn_dir, 'zdfiwm_init', 'iwm input file', 'namiwm' )
400
401      !                             ! hard-coded default values
402      sf_iwm(jp_mpb)%fnow(:,:,1) = 1.e-10_wp
403      sf_iwm(jp_mpc)%fnow(:,:,1) = 1.e-10_wp
404      sf_iwm(jp_mpn)%fnow(:,:,1) = 1.e-5_wp
405      sf_iwm(jp_mps)%fnow(:,:,1) = 1.e-10_wp
406      sf_iwm(jp_dsb)%fnow(:,:,1) = 100._wp
407      sf_iwm(jp_dsc)%fnow(:,:,1) = 100._wp
408
409      !                             ! read necessary fields
410      CALL fld_read( nit000, 1, sf_iwm )
411
412      ebot_iwm(:,:) = sf_iwm(1)%fnow(:,:,1) * ssmask(:,:) ! energy flux for dissipation above abyssal hills [W/m2]
413      ecri_iwm(:,:) = sf_iwm(2)%fnow(:,:,1) * ssmask(:,:) ! energy flux for dissipation at topographic slopes [W/m2]
414      ensq_iwm(:,:) = sf_iwm(3)%fnow(:,:,1) * ssmask(:,:) ! energy flux for dissipation scaling with N^2 [W/m2]
415      esho_iwm(:,:) = sf_iwm(4)%fnow(:,:,1) * ssmask(:,:) ! energy flux for dissipation due to shoaling [W/m2]
416      hbot_iwm(:,:) = sf_iwm(5)%fnow(:,:,1)               ! spatially variable decay scale for abyssal hill dissipation [m]
417      hcri_iwm(:,:) = sf_iwm(6)%fnow(:,:,1)               ! spatially variable decay scale for topographic-slope [m]
418
419      hcri_iwm(:,:) = 1._wp / hcri_iwm(:,:) ! only the inverse height is used, hence calculated here once for all
420
421      ! diags
422      ztmp(:,:,1) = e1e2t(:,:) * ebot_iwm(:,:)
423      ztmp(:,:,2) = e1e2t(:,:) * ecri_iwm(:,:)
424      ztmp(:,:,3) = e1e2t(:,:) * ensq_iwm(:,:)
425      ztmp(:,:,4) = e1e2t(:,:) * esho_iwm(:,:)
426
427      zdia(1:4) =glob_sum_vec( 'zdfiwm', CASTDP(ztmp(:,:,1:4)) )
428
429      IF(lwp) THEN
430         WRITE(numout,*) '      Dissipation above abyssal hills:        ', zdia(1) * 1.e-12_wp, 'TW'
431         WRITE(numout,*) '      Dissipation along topographic slopes:   ', zdia(2) * 1.e-12_wp, 'TW'
432         WRITE(numout,*) '      Dissipation scaling with N^2:           ', zdia(3) * 1.e-12_wp, 'TW'
433         WRITE(numout,*) '      Dissipation due to shoaling:            ', zdia(4) * 1.e-12_wp, 'TW'
434      ENDIF
435      !
436   END SUBROUTINE zdf_iwm_init
437
438   !!======================================================================
439END MODULE zdfiwm
Note: See TracBrowser for help on using the repository browser.