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.
tramle.F90 in NEMO/trunk/src/OCE/TRA – NEMO

source: NEMO/trunk/src/OCE/TRA/tramle.F90 @ 13982

Last change on this file since 13982 was 13982, checked in by smasson, 3 years ago

trunk: merge dev_r13923_Tiling_Cleanup_MPI3_LoopFusion into the trunk

  • Property svn:keywords set to Id
File size: 18.0 KB
Line 
1MODULE tramle
2   !!======================================================================
3   !!                    ***  MODULE  tramle  ***
4   !! Ocean tracers: Mixed Layer Eddy induced transport
5   !!======================================================================
6   !! History :  3.3  !  2010-08  (G. Madec)  Original code
7   !!----------------------------------------------------------------------
8
9   !!----------------------------------------------------------------------
10   !!   tra_mle_trp   : update the effective transport with the Mixed Layer Eddy induced transport
11   !!   tra_mle_init  : initialisation of the Mixed Layer Eddy induced transport computation
12   !!----------------------------------------------------------------------
13   USE oce            ! ocean dynamics and tracers variables
14   USE dom_oce        ! ocean space and time domain variables
15   USE phycst         ! physical constant
16   USE zdfmxl         ! mixed layer depth
17   !
18   USE in_out_manager ! I/O manager
19   USE iom            ! IOM library
20   USE lib_mpp        ! MPP library
21   USE lbclnk         ! lateral boundary condition / mpp link
22
23   IMPLICIT NONE
24   PRIVATE
25
26   PUBLIC   tra_mle_trp        ! routine called in traadv.F90
27   PUBLIC   tra_mle_init   ! routine called in traadv.F90
28
29   !                                    !!* namelist namtra_mle *
30   LOGICAL, PUBLIC ::   ln_mle           !: flag to activate the Mixed Layer Eddy (MLE) parameterisation
31   INTEGER         ::      nn_mle           ! MLE type: =0 standard Fox-Kemper ; =1 new formulation
32   INTEGER         ::      nn_mld_uv        ! space interpolation of MLD at u- & v-pts (0=min,1=averaged,2=max)
33   INTEGER         ::      nn_conv          ! =1 no MLE in case of convection ; =0 always MLE
34   REAL(wp)        ::      rn_ce            ! MLE coefficient
35   !                                        ! parameters used in nn_mle = 0 case
36   REAL(wp)        ::      rn_lf               ! typical scale of mixed layer front
37   REAL(wp)        ::      rn_time             ! time scale for mixing momentum across the mixed layer
38   !                                        ! parameters used in nn_mle = 1 case
39   REAL(wp)        ::      rn_lat              ! reference latitude for a 5 km scale of ML front
40   REAL(wp)        ::      rn_rho_c_mle        ! Density criterion for definition of MLD used by FK
41
42   REAL(wp) ::   r5_21 = 5.e0 / 21.e0   ! factor used in mle streamfunction computation
43   REAL(wp) ::   rb_c                   ! ML buoyancy criteria = g rho_c /rho0 where rho_c is defined in zdfmld
44   REAL(wp) ::   rc_f                   ! MLE coefficient (= rn_ce / (5 km * fo) ) in nn_mle=1 case
45
46   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   rfu, rfv   ! modified Coriolis parameter (f+tau) at u- & v-pts
47   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   r1_ft      ! inverse of the modified Coriolis parameter at t-pts
48
49   !! * Substitutions
50#  include "do_loop_substitute.h90"
51#  include "domzgr_substitute.h90"
52   !!----------------------------------------------------------------------
53   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
54   !! $Id$
55   !! Software governed by the CeCILL license (see ./LICENSE)
56   !!----------------------------------------------------------------------
57CONTAINS
58
59   SUBROUTINE tra_mle_trp( kt, kit000, pu, pv, pw, cdtype, Kmm )
60      !!----------------------------------------------------------------------
61      !!                  ***  ROUTINE tra_mle_trp  ***
62      !!
63      !! ** Purpose :   Add to the transport the Mixed Layer Eddy induced transport
64      !!
65      !! ** Method  :   The 3 components of the Mixed Layer Eddy (MLE) induced
66      !!              transport are computed as follows :
67      !!                zu_mle = dk[ zpsi_uw ]
68      !!                zv_mle = dk[ zpsi_vw ]
69      !!                zw_mle = - di[ zpsi_uw ] - dj[ zpsi_vw ]
70      !!                where zpsi is the MLE streamfunction at uw and vw points (see the doc)
71      !!              and added to the input velocity :
72      !!                p.n = p.n + z._mle
73      !!
74      !! ** Action  : - (pu,pv,pw) increased by the mle transport
75      !!                CAUTION, the transport is not updated at the last line/raw
76      !!                         this may be a problem for some advection schemes
77      !!
78      !! References: Fox-Kemper et al., JPO, 38, 1145-1165, 2008
79      !!             Fox-Kemper and Ferrari, JPO, 38, 1166-1179, 2008
80      !!----------------------------------------------------------------------
81      INTEGER                     , INTENT(in   ) ::   kt         ! ocean time-step index
82      INTEGER                     , INTENT(in   ) ::   kit000     ! first time step index
83      INTEGER                     , INTENT(in   ) ::   Kmm        ! ocean time level index
84      CHARACTER(len=3)            , INTENT(in   ) ::   cdtype     ! =TRA or TRC (tracer indicator)
85      REAL(wp), DIMENSION(A2D(nn_hls),jpk), INTENT(inout) ::   pu         ! in : 3 ocean transport components
86      REAL(wp), DIMENSION(A2D(nn_hls),jpk), INTENT(inout) ::   pv         ! out: same 3  transport components
87      REAL(wp), DIMENSION(A2D(nn_hls),jpk), INTENT(inout) ::   pw         !   increased by the MLE induced transport
88      !
89      INTEGER  ::   ji, jj, jk          ! dummy loop indices
90      INTEGER  ::   ii, ij, ik, ikmax   ! local integers
91      REAL(wp) ::   zcuw, zmuw, zc      ! local scalar
92      REAL(wp) ::   zcvw, zmvw          !   -      -
93      INTEGER , DIMENSION(A2D(nn_hls))     :: inml_mle
94      REAL(wp), DIMENSION(A2D(nn_hls))     :: zpsim_u, zpsim_v, zmld, zbm, zhu, zhv, zn2, zLf_MH
95      REAL(wp), DIMENSION(A2D(nn_hls),jpk) :: zpsi_uw, zpsi_vw
96      ! TEMP: [tiling] These changes not necessary if using XIOS (subdomain support)
97      REAL(wp), DIMENSION(:,:),   ALLOCATABLE, SAVE :: zLf_NH
98      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, SAVE :: zpsiu_mle, zpsiv_mle
99      !!----------------------------------------------------------------------
100      !
101      !                                      !==  MLD used for MLE  ==!
102      !                                                ! compute from the 10m density to deal with the diurnal cycle
103      DO_2D( 1, 1, 1, 1 )
104         inml_mle(ji,jj) = mbkt(ji,jj) + 1                    ! init. to number of ocean w-level (T-level + 1)
105      END_2D
106      IF ( nla10 > 0 ) THEN                            ! avoid case where first level is thicker than 10m
107         DO_3DS( 1, 1, 1, 1, jpkm1, nlb10, -1 )        ! from the bottom to nlb10 (10m)
108            IF( rhop(ji,jj,jk) > rhop(ji,jj,nla10) + rn_rho_c_mle )   inml_mle(ji,jj) = jk      ! Mixed layer
109         END_3D
110      ENDIF
111      ikmax = MIN( MAXVAL( inml_mle(:,:) ), jpkm1 )                  ! max level of the computation
112      !
113      !
114      zmld(:,:) = 0._wp                      !==   Horizontal shape of the MLE  ==!
115      zbm (:,:) = 0._wp
116      zn2 (:,:) = 0._wp
117      DO_3D( 1, 1, 1, 1, 1, ikmax )                    ! MLD and mean buoyancy and N2 over the mixed layer
118         zc = e3t(ji,jj,jk,Kmm) * REAL( MIN( MAX( 0, inml_mle(ji,jj)-jk ) , 1  )  )    ! zc being 0 outside the ML t-points
119         zmld(ji,jj) = zmld(ji,jj) + zc
120         zbm (ji,jj) = zbm (ji,jj) + zc * (rho0 - rhop(ji,jj,jk) ) * r1_rho0
121         zn2 (ji,jj) = zn2 (ji,jj) + zc * (rn2(ji,jj,jk)+rn2(ji,jj,jk+1))*0.5_wp
122      END_3D
123
124      SELECT CASE( nn_mld_uv )                         ! MLD at u- & v-pts
125      CASE ( 0 )                                               != min of the 2 neighbour MLDs
126         DO_2D( 1, 0, 1, 0 )
127            zhu(ji,jj) = MIN( zmld(ji+1,jj), zmld(ji,jj) )
128            zhv(ji,jj) = MIN( zmld(ji,jj+1), zmld(ji,jj) )
129         END_2D
130      CASE ( 1 )                                               != average of the 2 neighbour MLDs
131         DO_2D( 1, 0, 1, 0 )
132            zhu(ji,jj) = ( zmld(ji+1,jj) + zmld(ji,jj) ) * 0.5_wp
133            zhv(ji,jj) = ( zmld(ji,jj+1) + zmld(ji,jj) ) * 0.5_wp
134         END_2D
135      CASE ( 2 )                                               != max of the 2 neighbour MLDs
136         DO_2D( 1, 0, 1, 0 )
137            zhu(ji,jj) = MAX( zmld(ji+1,jj), zmld(ji,jj) )
138            zhv(ji,jj) = MAX( zmld(ji,jj+1), zmld(ji,jj) )
139         END_2D
140      END SELECT
141      !                                                ! convert density into buoyancy
142      DO_2D( 1, 1, 1, 1 )
143         zbm(ji,jj) = + grav * zbm(ji,jj) / MAX( e3t(ji,jj,1,Kmm), zmld(ji,jj) )
144      END_2D
145      !
146      !
147      !                                      !==  Magnitude of the MLE stream function  ==!
148      !
149      !                 di[bm]  Ds
150      ! Psi = Ce  H^2 ---------------- e2u  mu(z)   where fu Lf = MAX( fu*rn_fl , (Db H)^1/2 )
151      !                  e1u   Lf fu                      and the e2u for the "transport"
152      !                                                      (not *e3u as divided by e3u at the end)
153      !
154      IF( nn_mle == 0 ) THEN           ! Fox-Kemper et al. 2010 formulation
155         DO_2D( 1, 0, 1, 0 )
156            zpsim_u(ji,jj) = rn_ce * zhu(ji,jj) * zhu(ji,jj)  * e2_e1u(ji,jj)                                            &
157               &           * ( zbm(ji+1,jj) - zbm(ji,jj) ) * MIN( 111.e3_wp , e1u(ji,jj) )   &
158               &           / (  MAX( rn_lf * rfu(ji,jj) , SQRT( rb_c * zhu(ji,jj) ) )   )
159               !
160            zpsim_v(ji,jj) = rn_ce * zhv(ji,jj) * zhv(ji,jj)  * e1_e2v(ji,jj)                                            &
161               &           * ( zbm(ji,jj+1) - zbm(ji,jj) ) * MIN( 111.e3_wp , e2v(ji,jj) )   &
162               &           / (  MAX( rn_lf * rfv(ji,jj) , SQRT( rb_c * zhv(ji,jj) ) )   )
163         END_2D
164         !
165      ELSEIF( nn_mle == 1 ) THEN       ! New formulation (Lf = 5km fo/ff with fo=Coriolis parameter at latitude rn_lat)
166         DO_2D( 1, 0, 1, 0 )
167            zpsim_u(ji,jj) = rc_f *   zhu(ji,jj)   * zhu(ji,jj)   * e2_e1u(ji,jj)               &
168               &                  * ( zbm(ji+1,jj) - zbm(ji,jj) ) * MIN( 111.e3_wp , e1u(ji,jj) )
169               !
170            zpsim_v(ji,jj) = rc_f *   zhv(ji,jj)   * zhv(ji,jj)   * e1_e2v(ji,jj)               &
171               &                  * ( zbm(ji,jj+1) - zbm(ji,jj) ) * MIN( 111.e3_wp , e2v(ji,jj) )
172         END_2D
173      ENDIF
174      !
175      IF( nn_conv == 1 ) THEN              ! No MLE in case of convection
176         DO_2D( 1, 0, 1, 0 )
177            IF( MIN( zn2(ji,jj) , zn2(ji+1,jj) ) < 0._wp )   zpsim_u(ji,jj) = 0._wp
178            IF( MIN( zn2(ji,jj) , zn2(ji,jj+1) ) < 0._wp )   zpsim_v(ji,jj) = 0._wp
179         END_2D
180      ENDIF
181      !
182      !                                      !==  structure function value at uw- and vw-points  ==!
183      DO_2D( 1, 0, 1, 0 )
184         zhu(ji,jj) = 1._wp / zhu(ji,jj)                   ! hu --> 1/hu
185         zhv(ji,jj) = 1._wp / zhv(ji,jj)
186      END_2D
187      !
188      zpsi_uw(:,:,:) = 0._wp
189      zpsi_vw(:,:,:) = 0._wp
190      !
191      DO_3D( 1, 0, 1, 0, 2, ikmax )                ! start from 2 : surface value = 0
192         zcuw = 1._wp - ( gdepw(ji+1,jj,jk,Kmm) + gdepw(ji,jj,jk,Kmm) ) * zhu(ji,jj)
193         zcvw = 1._wp - ( gdepw(ji,jj+1,jk,Kmm) + gdepw(ji,jj,jk,Kmm) ) * zhv(ji,jj)
194         zcuw = zcuw * zcuw
195         zcvw = zcvw * zcvw
196         zmuw = MAX(  0._wp , ( 1._wp - zcuw ) * ( 1._wp + r5_21 * zcuw )  )
197         zmvw = MAX(  0._wp , ( 1._wp - zcvw ) * ( 1._wp + r5_21 * zcvw )  )
198         !
199         zpsi_uw(ji,jj,jk) = zpsim_u(ji,jj) * zmuw * umask(ji,jj,jk)
200         zpsi_vw(ji,jj,jk) = zpsim_v(ji,jj) * zmvw * vmask(ji,jj,jk)
201      END_3D
202      !
203      !                                      !==  transport increased by the MLE induced transport ==!
204      DO jk = 1, ikmax
205         DO_2D( 1, 0, 1, 0 )                      ! CAUTION pu,pv must be defined at row/column i=1 / j=1
206            pu(ji,jj,jk) = pu(ji,jj,jk) + ( zpsi_uw(ji,jj,jk) - zpsi_uw(ji,jj,jk+1) )
207            pv(ji,jj,jk) = pv(ji,jj,jk) + ( zpsi_vw(ji,jj,jk) - zpsi_vw(ji,jj,jk+1) )
208         END_2D
209         DO_2D( 0, 0, 0, 0 )
210            pw(ji,jj,jk) = pw(ji,jj,jk) - ( zpsi_uw(ji,jj,jk) - zpsi_uw(ji-1,jj,jk)   &
211               &                          + zpsi_vw(ji,jj,jk) - zpsi_vw(ji,jj-1,jk) )
212         END_2D
213      END DO
214
215      ! TEMP: [tiling] These changes not necessary if using XIOS (subdomain support)
216      IF( cdtype == 'TRA') THEN              !==  outputs  ==!
217         IF( ntile == 0 .OR. ntile == 1 ) THEN                             ! Do only on the first tile
218            ALLOCATE( zLf_NH(jpi,jpj), zpsiu_mle(jpi,jpj,jpk), zpsiv_mle(jpi,jpj,jpk) )
219            zpsiu_mle(:,:,:) = 0._wp ; zpsiv_mle(:,:,:) = 0._wp
220         ENDIF
221         !
222         DO_2D( 0, 0, 0, 0 )
223            zLf_NH(ji,jj) = SQRT( rb_c * zmld(ji,jj) ) * r1_ft(ji,jj)      ! Lf = N H / f
224         END_2D
225         !
226         ! divide by cross distance to give streamfunction with dimensions m^2/s
227         DO_3D( 0, 0, 0, 0, 1, ikmax+1 )
228            zpsiu_mle(ji,jj,jk) = zpsi_uw(ji,jj,jk) * r1_e2u(ji,jj)
229            zpsiv_mle(ji,jj,jk) = zpsi_vw(ji,jj,jk) * r1_e1v(ji,jj)
230         END_3D
231
232         IF( ntile == 0 .OR. ntile == nijtile ) THEN                       ! Do only on the last tile
233            CALL iom_put( "Lf_NHpf" , zLf_NH  )    ! Lf = N H / f
234            CALL iom_put( "psiu_mle", zpsiu_mle )    ! i-mle streamfunction
235            CALL iom_put( "psiv_mle", zpsiv_mle )    ! j-mle streamfunction
236            DEALLOCATE( zLf_NH, zpsiu_mle, zpsiv_mle )
237         ENDIF
238      ENDIF
239      !
240   END SUBROUTINE tra_mle_trp
241
242
243   SUBROUTINE tra_mle_init
244      !!---------------------------------------------------------------------
245      !!                  ***  ROUTINE tra_mle_init  ***
246      !!
247      !! ** Purpose :   Control the consistency between namelist options for
248      !!              tracer advection schemes and set nadv
249      !!----------------------------------------------------------------------
250      INTEGER  ::   ji, jj, jk   ! dummy loop indices
251      INTEGER  ::   ierr
252      INTEGER ::    ios                 ! Local integer output status for namelist read
253      REAL(wp) ::   z1_t2, zfu, zfv                                !    -         -
254      !
255      NAMELIST/namtra_mle/ ln_mle , nn_mle, rn_ce, rn_lf, rn_time, rn_lat, nn_mld_uv, nn_conv, rn_rho_c_mle
256      !!----------------------------------------------------------------------
257
258      READ  ( numnam_ref, namtra_mle, IOSTAT = ios, ERR = 901)
259901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namtra_mle in reference namelist' )
260
261      READ  ( numnam_cfg, namtra_mle, IOSTAT = ios, ERR = 902 )
262902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namtra_mle in configuration namelist' )
263      IF(lwm) WRITE ( numond, namtra_mle )
264
265      IF(lwp) THEN                     ! Namelist print
266         WRITE(numout,*)
267         WRITE(numout,*) 'tra_mle_init : mixed layer eddy (MLE) advection acting on tracers'
268         WRITE(numout,*) '~~~~~~~~~~~~~'
269         WRITE(numout,*) '   Namelist namtra_mle : mixed layer eddy advection applied on tracers'
270         WRITE(numout,*) '      use mixed layer eddy (MLE, i.e. Fox-Kemper param) (T/F)      ln_mle       = ', ln_mle
271         WRITE(numout,*) '         MLE type: =0 standard Fox-Kemper ; =1 new formulation        nn_mle    = ', nn_mle
272         WRITE(numout,*) '         magnitude of the MLE (typical value: 0.06 to 0.08)           rn_ce     = ', rn_ce
273         WRITE(numout,*) '         scale of ML front (ML radius of deformation) (rn_mle=0)      rn_lf     = ', rn_lf, 'm'
274         WRITE(numout,*) '         maximum time scale of MLE                    (rn_mle=0)      rn_time   = ', rn_time, 's'
275         WRITE(numout,*) '         reference latitude (degrees) of MLE coef.    (rn_mle=1)      rn_lat    = ', rn_lat, 'deg'
276         WRITE(numout,*) '         space interp. of MLD at u-(v-)pts (0=min,1=averaged,2=max)   nn_mld_uv = ', nn_mld_uv
277         WRITE(numout,*) '         =1 no MLE in case of convection ; =0 always MLE              nn_conv   = ', nn_conv
278         WRITE(numout,*) '         Density difference used to define ML for FK              rn_rho_c_mle  = ', rn_rho_c_mle
279      ENDIF
280      !
281      IF(lwp) THEN
282         WRITE(numout,*)
283         IF( ln_mle ) THEN
284            WRITE(numout,*) '   ==>>>   Mixed Layer Eddy induced transport added to tracer advection'
285            IF( nn_mle == 0 )   WRITE(numout,*) '              Fox-Kemper et al 2010 formulation'
286            IF( nn_mle == 1 )   WRITE(numout,*) '              New formulation'
287         ELSE
288            WRITE(numout,*) '   ==>>>   Mixed Layer Eddy parametrisation NOT used'
289         ENDIF
290      ENDIF
291      !
292      IF( ln_mle ) THEN                ! MLE initialisation
293         !
294         rb_c = grav * rn_rho_c_mle /rho0        ! Mixed Layer buoyancy criteria
295         IF(lwp) WRITE(numout,*)
296         IF(lwp) WRITE(numout,*) '      ML buoyancy criteria = ', rb_c, ' m/s2 '
297         IF(lwp) WRITE(numout,*) '      associated ML density criteria defined in zdfmxl = ', rho_c, 'kg/m3'
298         !
299         IF( nn_mle == 0 ) THEN           ! MLE array allocation & initialisation
300            ALLOCATE( rfu(jpi,jpj) , rfv(jpi,jpj) , STAT= ierr )
301            IF( ierr /= 0 )   CALL ctl_stop( 'tra_adv_mle_init: failed to allocate arrays' )
302            z1_t2 = 1._wp / ( rn_time * rn_time )
303            DO_2D( nn_hls-1, nn_hls, nn_hls-1, nn_hls )                      ! "coriolis+ time^-1" at u- & v-points
304               zfu = ( ff_f(ji,jj) + ff_f(ji,jj-1) ) * 0.5_wp
305               zfv = ( ff_f(ji,jj) + ff_f(ji-1,jj) ) * 0.5_wp
306               rfu(ji,jj) = SQRT(  zfu * zfu + z1_t2 )
307               rfv(ji,jj) = SQRT(  zfv * zfv + z1_t2 )
308            END_2D
309            IF (nn_hls.EQ.1) CALL lbc_lnk_multi( 'tramle', rfu, 'U', 1.0_wp , rfv, 'V', 1.0_wp )
310            !
311         ELSEIF( nn_mle == 1 ) THEN           ! MLE array allocation & initialisation
312            rc_f = rn_ce / (  5.e3_wp * 2._wp * omega * SIN( rad * rn_lat )  )
313            !
314         ENDIF
315         !
316         !                                ! 1/(f^2+tau^2)^1/2 at t-point (needed in both nn_mle case)
317         ALLOCATE( r1_ft(jpi,jpj) , STAT= ierr )
318         IF( ierr /= 0 )   CALL ctl_stop( 'tra_adv_mle_init: failed to allocate r1_ft array' )
319         !
320         z1_t2 = 1._wp / ( rn_time * rn_time )
321         r1_ft(:,:) = 1._wp / SQRT(  ff_t(:,:) * ff_t(:,:) + z1_t2  )
322         !
323      ENDIF
324      !
325   END SUBROUTINE tra_mle_init
326
327   !!==============================================================================
328END MODULE tramle
Note: See TracBrowser for help on using the repository browser.