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 @ 2528

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

Update NEMOGCM from branch nemo_v3_3_beta

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