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.
zdftmx.F90 in trunk/NEMO/OPA_SRC/ZDF – NEMO

source: trunk/NEMO/OPA_SRC/ZDF/zdftmx.F90 @ 1527

Last change on this file since 1527 was 1527, checked in by ctlod, 15 years ago

tke2: minor bug correction on the TKE vertical diffusion term, see ticket: #484

File size: 23.8 KB
Line 
1MODULE zdftmx
2   !!========================================================================
3   !!                       ***  MODULE  zdftmx  ***
4   !! Ocean physics: vertical tidal mixing coefficient
5   !!========================================================================
6   !! History :  1.0  !  2004-04  (L. Bessieres, G. Madec)  Original code
7   !!             -   !  2006-08  (A. Koch-Larrouy) Indonesian strait
8   !!----------------------------------------------------------------------
9#if defined key_zdftmx
10   !!----------------------------------------------------------------------
11   !!   'key_zdftmx'                                  Tidal vertical mixing
12   !!----------------------------------------------------------------------
13   !!   zdf_tmx      : global     momentum & tracer Kz with tidal induced Kz
14   !!   tmx_itf      : Indonesian momentum & tracer Kz with tidal induced Kz
15   !!----------------------------------------------------------------------
16   USE oce             ! ocean dynamics and tracers variables
17   USE dom_oce         ! ocean space and time domain variables
18   USE zdf_oce         ! ocean vertical physics variables
19   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
20   USE eosbn2          ! ocean equation of state
21   USE phycst          ! physical constants
22   USE prtctl          ! Print control
23   USE in_out_manager  ! I/O manager
24   USE iom             ! I/O Manager
25
26   IMPLICIT NONE
27   PRIVATE
28
29   PUBLIC   zdf_tmx    ! called in step module
30
31   LOGICAL, PUBLIC, PARAMETER ::   lk_zdftmx = .TRUE.    !: tidal mixing flag
32
33   !                                  !!* Namelist  namtmx : tidal mixing *
34   REAL(wp) ::  rn_htmx    = 500.      ! vertical decay scale for turbulence (meters)
35   REAL(wp) ::  rn_n2min   = 1.e-8     ! threshold of the Brunt-Vaisala frequency (s-1)
36   REAL(wp) ::  rn_tfe     = 1./3.     ! tidal dissipation efficiency (St Laurent et al. 2002)
37   REAL(wp) ::  rn_me      = 0.2       ! mixing efficiency (Osborn 1980)
38   LOGICAL  ::  ln_tmx_itf = .TRUE.    ! Indonesian Through Flow (ITF): Koch-Larrouy et al. (2007) parameterization
39   REAL(wp) ::  rn_tfe_itf = 1.        ! ITF tidal dissipation efficiency (St Laurent et al. 2002)
40
41   REAL(wp), DIMENSION(jpi,jpj)     ::   en_tmx     ! energy available for tidal mixing (W/m2)
42   REAL(wp), DIMENSION(jpi,jpj)     ::   mask_itf   ! mask to use over Indonesian area
43   REAL(wp), DIMENSION(jpi,jpj,jpk) ::   az_tmx     ! coefficient used to evaluate the tidal induced Kz
44
45   !! * Substitutions
46#  include "domzgr_substitute.h90"
47#  include "vectopt_loop_substitute.h90"
48   !!----------------------------------------------------------------------
49   !! NEMO/OPA 3.2 , LOCEAN-IPSL (2009)
50   !! $Id$
51   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
52   !!----------------------------------------------------------------------
53
54CONTAINS
55
56   SUBROUTINE zdf_tmx( kt )
57      !!----------------------------------------------------------------------
58      !!                  ***  ROUTINE zdf_tmx  ***
59      !!                   
60      !! ** Purpose :   add to the vertical mixing coefficients the effect of
61      !!              tidal mixing (Simmons et al 2004).
62      !!
63      !! ** Method  : - tidal-induced vertical mixing is given by:
64      !!                  Kz_tides = az_tmx / max( rn_n2min, N^2 )
65      !!              where az_tmx is a coefficient that specified the 3D space
66      !!              distribution of the faction of tidal energy taht is used
67      !!              for mixing. Its expression is set in zdf_tmx_init routine,
68      !!              following Simmons et al. 2004.
69      !!                NB: a specific bounding procedure is performed on av_tide
70      !!              so that the input tidal energy is actually almost used. The
71      !!              basic maximum value is 60 cm2/s, but values of 300 cm2/s
72      !!              can be reached in area where bottom stratification is too
73      !!              weak.
74      !!
75      !!              - update av_tide in the Indonesian Through Flow area
76      !!              following Koch-Larrouy et al. (2007) parameterisation
77      !!              (see tmx_itf routine).
78      !!
79      !!              - update the model vertical eddy viscosity and diffusivity:
80      !!                     avt  = avt  +    av_tides
81      !!                     avm  = avm  +    av_tides
82      !!                     avmu = avmu + mi(av_tides)
83      !!                     avmv = avmv + mj(av_tides)
84      !!
85      !! ** Action  :   avt, avm, avmu, avmv   increased by tidal mixing
86      !!
87      !! References : Simmons et al. 2004, Ocean Modelling, 6, 3-4, 245-263.
88      !!              Koch-Larrouy et al. 2007, GRL.
89      !!----------------------------------------------------------------------
90      INTEGER, INTENT(in) ::   kt   ! ocean time-step
91      !!
92      INTEGER  ::   ji, jj, jk   ! dummy loop indices
93      REAL(wp) ::   ztpc         ! scalar workspace
94      REAL(wp), DIMENSION(jpi,jpj) ::   zkz   ! temporary 2D workspace
95      !!----------------------------------------------------------------------
96
97      !                       
98      IF( kt == nit000  )   CALL zdf_tmx_init      ! Initialization (first time-step only)
99
100      !                                        ! ----------------------- !
101      !                                        !  Standard tidal mixing  !  (compute av_tide)
102      !                                        ! ----------------------- !
103      !                             !* First estimation (with n2 bound by rn_n2min) bounded by 60 cm2/s
104      av_tide(:,:,:) = MIN(  60.e-4, az_tmx(:,:,:) / MAX( rn_n2min, rn2(:,:,:) )  )
105
106      zkz(:,:) = 0.e0               !* Associated potential energy consummed over the whole water column
107      DO jk = 2, jpkm1
108         zkz(:,:) = zkz(:,:) + fse3w(:,:,jk) * MAX( 0.e0, rn2(:,:,jk) ) * rau0 * av_tide(:,:,jk)* tmask(:,:,jk)
109      END DO
110
111      DO jj = 1, jpj                !* Here zkz should be equal to en_tmx ==> multiply by en_tmx/zkz to recover en_tmx
112         DO ji = 1, jpi
113            IF( zkz(ji,jj) /= 0.e0 )   zkz(ji,jj) = en_tmx(ji,jj) / zkz(ji,jj)
114         END DO
115      END DO
116
117      DO jk = 2, jpkm1              !* Mutiply by zkz to recover en_tmx, BUT bound by 30/6 ==> av_tide bound by 300 cm2/s
118         av_tide(:,:,jk) = av_tide(:,:,jk) * MIN( zkz(:,:), 30./6. )   !kz max = 300 cm2/s
119      END DO
120
121      IF( kt == nit000 ) THEN       !* check at first time-step: diagnose the energy consumed by av_tide
122         ztpc = 0.e0
123         DO jk= 1, jpk
124            DO jj= 1, jpj
125               DO ji= 1, jpi
126                  ztpc = ztpc + fse3w(ji,jj,jk) * e1t(ji,jj) * e2t(ji,jj)   &
127                     &         * MAX( 0.e0, rn2(ji,jj,jk) ) * av_tide(ji,jj,jk) * tmask(ji,jj,jk) * tmask_i(ji,jj)
128               END DO
129            END DO
130         END DO
131         ztpc= rau0 / ( rn_tfe * rn_me ) * ztpc
132         IF(lwp) WRITE(numout,*) 
133         IF(lwp) WRITE(numout,*) '          N Total power consumption by av_tide    : ztpc = ', ztpc * 1.e-12 ,'TW'
134      ENDIF
135       
136      !                                        ! ----------------------- !
137      IF( ln_tmx_itf )   CALL tmx_itf( kt )    !    ITF  tidal mixing    !  (update av_tide)
138      !                                        ! ----------------------- !
139
140      !                                        ! ----------------------- !
141      !                                        !   Update  mixing coefs  !                         
142      !                                        ! ----------------------- !
143      DO jk = 2, jpkm1              !* update momentum & tracer diffusivity with tidal mixing
144         avt(:,:,jk) = avt(:,:,jk) + av_tide(:,:,jk)
145         avm(:,:,jk) = avm(:,:,jk) + av_tide(:,:,jk)
146         DO jj = 2, jpjm1
147            DO ji = fs_2, fs_jpim1  ! vector opt.
148               avmu(ji,jj,jk) = avmu(ji,jj,jk) + 0.5 * ( av_tide(ji,jj,jk) + av_tide(ji+1,jj  ,jk) ) * umask(ji,jj,jk)
149               avmv(ji,jj,jk) = avmv(ji,jj,jk) + 0.5 * ( av_tide(ji,jj,jk) + av_tide(ji  ,jj+1,jk) ) * vmask(ji,jj,jk)
150            END DO
151         END DO
152      END DO
153      CALL lbc_lnk( avmu, 'U', 1. )   ;   CALL lbc_lnk( avmv, 'V', 1. )      ! lateral boundary condition
154
155      IF(ln_ctl)   CALL prt_ctl(tab3d_1=av_tide ,clinfo1=' tmx - av_tide: ',tab3d_2=avt,clinfo2=' avt: ',ovlap=1,kdim=jpk)
156      !
157   END SUBROUTINE zdf_tmx
158
159
160   SUBROUTINE tmx_itf( kt )
161      !!----------------------------------------------------------------------
162      !!                  ***  ROUTINE tmx_itf  ***
163      !!                   
164      !! ** Purpose :   modify the vertical eddy diffusivity coefficients
165      !!              (av_tide) in the Indonesian Through Flow area (ITF).
166      !!
167      !! ** Method  : - Following Koch-Larrouy et al. (2007), in the ITF defined
168      !!                by msk_itf (read in a file, see tmx_init), the tidal
169      !!                mixing coefficient is computed with :
170      !!                  * q=1 (i.e. all the tidal energy remains trapped in
171      !!                         the area and thus is used for mixing)
172      !!                  * the vertical distribution of the tifal energy is a
173      !!                    proportional to N above the thermocline (d(N^2)/dz > 0)
174      !!                    and to N^2 below the thermocline (d(N^2)/dz < 0)
175      !!
176      !! ** Action  :   av_tide   updated in the ITF area (msk_itf)
177      !!
178      !! References :  Koch-Larrouy et al. 2007, GRL
179      !!----------------------------------------------------------------------
180      INTEGER, INTENT(in) ::   kt   ! ocean time-step
181      !!
182      INTEGER  ::   ji, jj, jk    ! dummy loop indices
183      REAL(wp) ::   zcoef, ztpc   ! temporary scalar
184      REAL(wp), DIMENSION(jpi,jpj)     ::   zkz                        ! 2D workspace
185      REAL(wp), DIMENSION(jpi,jpj)     ::   zsum1 , zsum2 , zsum       !  -      -
186      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zempba_3d_1, zempba_3d_2   ! 3D workspace
187      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zempba_3d  , zdn2dz        !  -      -
188      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zavt_itf                   !  -      -
189      !!----------------------------------------------------------------------
190
191      !                             ! compute the form function using N2 at each time step
192      zempba_3d_1(:,:,jpk) = 0.e0
193      zempba_3d_2(:,:,jpk) = 0.e0
194      DO jk = 1, jpkm1             
195         zdn2dz     (:,:,jk) = rn2(:,:,jk) - rn2(:,:,jk+1)           ! Vertical profile of dN2/dz
196!CDIR NOVERRCHK
197         zempba_3d_1(:,:,jk) = SQRT(  MAX( 0.e0, rn2(:,:,jk) )  )    !    -        -    of N
198         zempba_3d_2(:,:,jk) =        MAX( 0.e0, rn2(:,:,jk) )       !    -        -    of N^2
199      END DO
200      !
201      zsum (:,:) = 0.e0
202      zsum1(:,:) = 0.e0
203      zsum2(:,:) = 0.e0
204      DO jk= 2, jpk
205         zsum1(:,:) = zsum1(:,:) + zempba_3d_1(:,:,jk) * fse3w(:,:,jk)
206         zsum2(:,:) = zsum2(:,:) + zempba_3d_2(:,:,jk) * fse3w(:,:,jk)               
207      END DO
208      DO jj = 1, jpj
209         DO ji = 1, jpi
210            IF( zsum1(ji,jj) /= 0.e0 )   zsum1(ji,jj) = 1.e0 / zsum1(ji,jj)
211            IF( zsum2(ji,jj) /= 0.e0 )   zsum2(ji,jj) = 1.e0 / zsum2(ji,jj)               
212         END DO
213      END DO
214
215      DO jk= 1, jpk
216         DO jj = 1, jpj
217            DO ji = 1, jpi
218               zcoef = 0.5 - SIGN( 0.5, zdn2dz(ji,jj,jk) )       ! =0 if dN2/dz > 0, =1 otherwise
219               ztpc  = zempba_3d_1(ji,jj,jk) * zsum1(ji,jj) *        zcoef     &
220                  &  + zempba_3d_2(ji,jj,jk) * zsum2(ji,jj) * ( 1. - zcoef )
221               !
222               zempba_3d(ji,jj,jk) =               ztpc 
223               zsum     (ji,jj)    = zsum(ji,jj) + ztpc * fse3w(ji,jj,jk)
224            END DO
225         END DO
226       END DO
227       DO jj = 1, jpj
228          DO ji = 1, jpi
229             IF( zsum(ji,jj) > 0.e0 )   zsum(ji,jj) = 1.e0 / zsum(ji,jj)               
230          END DO
231       END DO
232
233      !                             ! first estimation bounded by 10 cm2/s (with n2 bounded by rn_n2min)
234      zcoef = rn_tfe_itf / ( rn_tfe * rau0 )
235      DO jk = 1, jpk
236         zavt_itf(:,:,jk) = MIN(  10.e-4, zcoef * en_tmx(:,:) * zsum(:,:) * zempba_3d(:,:,jk)   &
237            &                                      / MAX( rn_n2min, rn2(:,:,jk) ) * tmask(:,:,jk)  )
238      END DO           
239
240      zkz(:,:) = 0.e0               ! Associated potential energy consummed over the whole water column
241      DO jk = 2, jpkm1
242         zkz(:,:) = zkz(:,:) + fse3w(:,:,jk) * MAX( 0.e0, rn2(:,:,jk) ) * rau0 * zavt_itf(:,:,jk) * tmask(:,:,jk)
243      END DO
244
245      DO jj = 1, jpj                ! Here zkz should be equal to en_tmx ==> multiply by en_tmx/zkz to recover en_tmx
246         DO ji = 1, jpi
247            IF( zkz(ji,jj) /= 0.e0 )   zkz(ji,jj) = en_tmx(ji,jj) * rn_tfe_itf / rn_tfe / zkz(ji,jj)
248         END DO
249      END DO
250
251      DO jk = 2, jpkm1              ! Mutiply by zkz to recover en_tmx, BUT bound by 30/6 ==> zavt_itf bound by 300 cm2/s
252         zavt_itf(:,:,jk) = zavt_itf(:,:,jk) * MIN( zkz(:,:), 120./10. )   ! kz max = 120 cm2/s
253      END DO
254
255      IF( kt == nit000 ) THEN       ! diagnose the nergy consumed by zavt_itf
256         ztpc = 0.e0
257         DO jk= 1, jpk
258            DO jj= 1, jpj
259               DO ji= 1, jpi
260                  ztpc = ztpc + e1t(ji,jj) * e2t(ji,jj) * fse3w(ji,jj,jk) * MAX( 0.e0, rn2(ji,jj,jk) )   &
261                     &                     * zavt_itf(ji,jj,jk) * tmask(ji,jj,jk) * tmask_i(ji,jj)
262               END DO
263            END DO
264         END DO
265         ztpc= rau0 * ztpc / ( rn_me * rn_tfe_itf )
266         IF(lwp) WRITE(numout,*) '          N Total power consumption by zavt_itf: ztpc = ', ztpc * 1.e-12 ,'TW'
267      ENDIF
268
269      !                             ! Update av_tide with the ITF mixing coefficient
270      DO jk = 2, jpkm1
271         av_tide(:,:,jk) = av_tide (:,:,jk) * ( 1.e0 - mask_itf(:,:) )   &
272            &            + zavt_itf(:,:,jk) *          mask_itf(:,:) 
273      END DO
274      !
275   END SUBROUTINE tmx_itf
276
277
278   SUBROUTINE zdf_tmx_init
279      !!----------------------------------------------------------------------
280      !!                  ***  ROUTINE zdf_tmx_init  ***
281      !!                     
282      !! ** Purpose :   Initialization of the vertical tidal mixing, Reading
283      !!              of M2 and K1 tidal energy in nc files
284      !!
285      !! ** Method  : - Read the namtmx namelist and check the parameters
286      !!
287      !!              - Read the input data in NetCDF files :
288      !!              M2 and K1 tidal energy. The total tidal energy, en_tmx,
289      !!              is the sum of M2, K1 and S2 energy where S2 is assumed
290      !!              to be: S2=(1/2)^2 * M2
291      !!              mask_itf, a mask array that determine where substituing
292      !!              the standard Simmons et al. (2005) formulation with the
293      !!              one of Koch_Larrouy et al. (2007).
294      !!
295      !!              - Compute az_tmx, a 3D coefficient that allows to compute
296      !!             the standard tidal-induced vertical mixing as follows:
297      !!                  Kz_tides = az_tmx / max( rn_n2min, N^2 )
298      !!             with az_tmx a bottom intensified coefficient is given by:
299      !!                 az_tmx(z) = en_tmx / ( rau0 * rn_htmx ) * EXP( -(H-z)/rn_htmx )
300      !!                                                  / ( 1. - EXP( - H   /rn_htmx ) )
301      !!             where rn_htmx the characteristic length scale of the bottom
302      !!             intensification, en_tmx the tidal energy, and H the ocean depth
303      !!
304      !! ** input   :   - Namlist namtmx
305      !!                - NetCDF file : M2_ORCA2.nc, K1_ORCA2.nc, and mask_itf.nc
306      !!
307      !! ** Action  : - Increase by 1 the nstop flag is setting problem encounter
308      !!              - defined az_tmx used to compute tidal-induced mixing
309      !!
310      !! References : Simmons et al. 2004, Ocean Modelling, 6, 3-4, 245-263.
311      !!              Koch-Larrouy et al. 2007, GRL.
312      !!----------------------------------------------------------------------
313      INTEGER ::   ji, jj, jk    ! dummy loop indices
314      INTEGER ::   inum          ! temporary logical unit
315      REAL(wp) ::   ztpc, ze_z   ! total power consumption
316      REAL(wp), DIMENSION(jpi,jpj) ::  zem2, zek1   ! read M2 and K1 tidal energy
317      REAL(wp), DIMENSION(jpi,jpj) ::  zkz          ! total M2, K1 and S2 tidal energy
318      REAL(wp), DIMENSION(jpi,jpj) ::  zfact        ! used for vertical structure function
319      REAL(wp), DIMENSION(jpi,jpj) ::  zhdep        ! Ocean depth
320      REAL(wp), DIMENSION(jpi,jpj,jpk) ::  zpc      ! power consumption
321      !!
322      NAMELIST/namtmx/ rn_htmx, rn_n2min, rn_tfe, rn_me, ln_tmx_itf, rn_tfe_itf
323      !!----------------------------------------------------------------------
324
325      IF(lwp) THEN
326         WRITE(numout,*)
327         WRITE(numout,*) 'zdf_tmx_init : tidal mixing'
328         WRITE(numout,*) '~~~~~~~~~~~~'
329      ENDIF
330
331      ! Read Namelist namtmx : Tidal Mixing
332      ! --------------------
333      REWIND ( numnam )
334      READ   ( numnam, namtmx )
335
336      ! Control print
337      IF(lwp) THEN
338         WRITE(numout,*) '          Namelist namtmx : set tidal mixing parameters'
339         WRITE(numout,*) '             Vertical decay scale for turbulence   = ', rn_htmx 
340         WRITE(numout,*) '             Brunt-Vaisala frequency threshold     = ', rn_n2min
341         WRITE(numout,*) '             Tidal dissipation efficiency          = ', rn_tfe
342         WRITE(numout,*) '             Mixing efficiency                     = ', rn_me
343         WRITE(numout,*) '             ITF specific parameterisation         = ', ln_tmx_itf
344         WRITE(numout,*) '             ITF tidal dissipation efficiency      = ', rn_tfe_itf
345         WRITE(numout,*)
346      ENDIF
347
348      IF( ln_tmx_itf ) THEN   ! read the Indonesian Through Flow mask
349         CALL iom_open('mask_itf',inum)
350         CALL iom_get (inum, jpdom_data, 'tmaskitf',mask_itf,1) !
351         CALL iom_close(inum)
352      ENDIF
353
354      ! read M2 tidal energy flux : W/m2  ( zem2 < 0 )
355      CALL iom_open('M2rowdrg',inum)
356      CALL iom_get (inum, jpdom_data, 'field',zem2,1) !
357      CALL iom_close(inum)
358
359      ! read K1 tidal energy flux : W/m2  ( zek1 < 0 )
360      CALL iom_open('K1rowdrg',inum)
361      CALL iom_get (inum, jpdom_data, 'field',zek1,1) !
362      CALL iom_close(inum)
363 
364      ! Total tidal energy ( M2, S2 and K1  with S2=(1/2)^2 * M2 )
365      ! only the energy available for mixing is taken into account,
366      ! (mixing efficiency tidal dissipation efficiency)
367      en_tmx(:,:) = - rn_tfe * rn_me * ( zem2(:,:) * 1.25 + zek1(:,:) ) * tmask(:,:,1)
368
369      ! Vertical structure (az_tmx)
370      DO jj = 1, jpj                ! part independent of the level
371         DO ji = 1, jpi
372            zhdep(ji,jj) = fsdepw(ji,jj,mbathy(ji,jj))       ! depth of the ocean
373            zfact(ji,jj) = rau0 * rn_htmx * ( 1. - EXP( -zhdep(ji,jj) / rn_htmx ) )
374            IF( zfact(ji,jj) /= 0 )   zfact(ji,jj) = en_tmx(ji,jj) / zfact(ji,jj)
375         END DO
376      END DO
377      DO jk= 1, jpk                 ! complete with the level-dependent part
378         DO jj = 1, jpj
379            DO ji = 1, jpi
380               az_tmx(ji,jj,jk) = zfact(ji,jj) * EXP( -( zhdep(ji,jj)-fsdepw(ji,jj,jk) ) / rn_htmx ) * tmask(ji,jj,jk)
381            END DO
382         END DO
383      END DO
384
385      IF( nprint == 1 .AND. lwp ) THEN
386         ! Control print
387         ! Total power consumption due to vertical mixing
388         ! zpc = rau0 * 1/rn_me * rn2 * av_tide
389         av_tide(:,:,:) = 0.e0
390         DO jk = 2, jpkm1
391            av_tide(:,:,jk) = az_tmx(:,:,jk) / MAX( rn_n2min, rn2(:,:,jk) )
392         END DO
393
394         ztpc = 0.e0
395         zpc(:,:,:) = MAX(rn_n2min,rn2(:,:,:)) * av_tide(:,:,:)
396         DO jk= 2, jpkm1
397            DO jj = 1, jpj
398               DO ji = 1, jpi
399                  ztpc = ztpc + fse3w(ji,jj,jk) * e1t(ji,jj) * e2t(ji,jj) * zpc(ji,jj,jk) * tmask(ji,jj,jk) * tmask_i(ji,jj)
400               END DO
401            END DO
402         END DO
403         ztpc= rau0 * 1/(rn_tfe * rn_me) * ztpc
404
405         WRITE(numout,*) 
406         WRITE(numout,*) '          Total power consumption of the tidally driven part of Kz : ztpc = ', ztpc * 1.e-12 ,'TW'
407
408
409         ! control print 2
410         av_tide(:,:,:) = MIN( av_tide(:,:,:), 60.e-4 )   
411         zkz(:,:) = 0.e0
412         DO jk = 2, jpkm1
413         DO jj = 1, jpj
414            DO ji = 1, jpi
415               zkz(ji,jj) = zkz(ji,jj) + fse3w(ji,jj,jk) * MAX( 0.e0, rn2(ji,jj,jk) ) * rau0 * av_tide(ji,jj,jk)* tmask(ji,jj,jk)
416            END DO
417         END DO
418         END DO
419         ! Here zkz should be equal to en_tmx ==> multiply by en_tmx/zkz
420         DO jj = 1, jpj
421            DO ji = 1, jpi
422               IF( zkz(ji,jj) /= 0.e0 )   THEN
423                   zkz(ji,jj) = en_tmx(ji,jj) / zkz(ji,jj)
424               ENDIF
425            END DO
426         END DO
427         ztpc = 1.e50
428         DO jj = 1, jpj
429            DO ji = 1, jpi
430               IF( zkz(ji,jj) /= 0.e0 )   THEN
431                   ztpc = Min( zkz(ji,jj), ztpc)
432               ENDIF
433            END DO
434         END DO
435         WRITE(numout,*) '          Min de zkz ', ztpc, ' Max = ', maxval(zkz(:,:) )
436
437         DO jk = 2, jpkm1
438            av_tide(:,:,jk) = av_tide(:,:,jk) * MIN( zkz(:,:), 30./6. )   !kz max = 300 cm2/s
439         END DO
440         ztpc = 0.e0
441         zpc(:,:,:) = Max(0.e0,rn2(:,:,:)) * av_tide(:,:,:)
442         DO jk= 1, jpk
443            DO jj = 1, jpj
444               DO ji = 1, jpi
445                  ztpc = ztpc + fse3w(ji,jj,jk) * e1t(ji,jj) * e2t(ji,jj) * zpc(ji,jj,jk) * tmask(ji,jj,jk) * tmask_i(ji,jj)
446               END DO
447            END DO
448         END DO
449         ztpc= rau0 * 1/(rn_tfe * rn_me) * ztpc
450         WRITE(numout,*) '          2 Total power consumption of the tidally driven part of Kz : ztpc = ', ztpc * 1.e-12 ,'TW'
451
452         DO jk = 1, jpk
453            ze_z =                  SUM( e1t(:,:) * e2t(:,:) * av_tide(:,:,jk)     * tmask_i(:,:) )   &
454               &     / MAX( 1.e-20, SUM( e1t(:,:) * e2t(:,:) * tmask (:,:,jk) * tmask_i(:,:) ) )
455            ztpc = 1.E50
456            DO jj = 1, jpj
457               DO ji = 1, jpi
458                  IF( av_tide(ji,jj,jk) /= 0.e0 )   ztpc =Min( ztpc, av_tide(ji,jj,jk) )
459               END DO
460            END DO
461            WRITE(numout,*) '            N2 min - jk= ', jk,'   ', ze_z * 1.e4,' cm2/s min= ',ztpc*1.e4,   &
462               &       'max= ', MAXVAL(av_tide(:,:,jk) )*1.e4, ' cm2/s'
463         END DO
464
465         WRITE(numout,*) '          e_tide : ', SUM( e1t*e2t*en_tmx ) / ( rn_tfe * rn_me ) * 1.e-12, 'TW'
466         WRITE(numout,*) 
467         WRITE(numout,*) '          Initial profile of tidal vertical mixing'
468         DO jk = 1, jpk
469            DO jj = 1,jpj
470               DO ji = 1,jpi
471                  zkz(ji,jj) = az_tmx(ji,jj,jk) /MAX( rn_n2min, rn2(ji,jj,jk) )
472               END DO
473            END DO
474            ze_z =                  SUM( e1t(:,:) * e2t(:,:) * zkz(:,:)     * tmask_i(:,:) )   &
475               &     / MAX( 1.e-20, SUM( e1t(:,:) * e2t(:,:) * tmask (:,:,jk) * tmask_i(:,:) ) )
476            WRITE(numout,*) '                jk= ', jk,'   ', ze_z * 1.e4,' cm2/s'
477         END DO
478         DO jk = 1, jpk
479            zkz(:,:) = az_tmx(:,:,jk) /rn_n2min
480            ze_z =                  SUM( e1t(:,:) * e2t(:,:) * zkz(:,:)     * tmask_i(:,:) )   &
481               &     / MAX( 1.e-20, SUM( e1t(:,:) * e2t(:,:) * tmask (:,:,jk) * tmask_i(:,:) ) )
482            WRITE(numout,*) 
483            WRITE(numout,*) '          N2 min - jk= ', jk,'   ', ze_z * 1.e4,' cm2/s min= ',MINVAL(zkz)*1.e4,   &
484               &       'max= ', MAXVAL(zkz)*1.e4, ' cm2/s'
485         END DO
486         !
487      ENDIF
488
489   END SUBROUTINE zdf_tmx_init
490
491#else
492   !!----------------------------------------------------------------------
493   !!   Default option          Dummy module                NO Tidal MiXing
494   !!----------------------------------------------------------------------
495   LOGICAL, PUBLIC, PARAMETER ::   lk_zdftmx = .FALSE.   !: tidal mixing flag
496CONTAINS
497   SUBROUTINE zdf_tmx( kt )          ! Empty routine
498      WRITE(*,*) 'zdf_tmx: You should not have seen this print! error?', kt
499   END SUBROUTINE zdf_tmx
500#endif
501
502   !!======================================================================
503END MODULE zdftmx
Note: See TracBrowser for help on using the repository browser.