source: CONFIG/UNIFORM/v7/IPSLCM7/SOURCES/NEMO/zdfiwm.F90

Last change on this file was 6346, checked in by cetlod, 16 months ago

IPSLCM7_WORK : 1st step of coupling DYNAMICO-LMDZ-ORCHIDEE-NEMO4.2 + Oasis3-mct.5.0

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