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

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

trunk: merge dev_r14312_MPI_Interface into the trunk, #2598

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