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/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: trunk/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftmx.F90 @ 2715

Last change on this file since 2715 was 2715, checked in by rblod, 13 years ago

First attempt to put dynamic allocation on the trunk

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