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.
traadv_mle.F90 in NEMO/branches/UKMO/r8395_coupling_sequence/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: NEMO/branches/UKMO/r8395_coupling_sequence/NEMOGCM/NEMO/OPA_SRC/TRA/traadv_mle.F90 @ 10763

Last change on this file since 10763 was 10763, checked in by jcastill, 5 years ago

Remove svn keywords properly

File size: 19.0 KB
Line 
1MODULE traadv_mle
2   !!======================================================================
3   !!                    ***  MODULE  traadv_mle  ***
4   !! Ocean tracers: Mixed Layer Eddy induced transport
5   !!======================================================================
6   !! History :  3.3  !  2010-08  (G. Madec)  Original code
7   !!----------------------------------------------------------------------
8
9   !!----------------------------------------------------------------------
10   !!   tra_adv_mle      : update the effective transport with the Mixed Layer Eddy induced transport
11   !!   tra_adv_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   USE lbclnk         ! lateral boundary condition / mpp link
18   USE in_out_manager ! I/O manager
19   USE iom            ! IOM library
20   USE lib_mpp        ! MPP library
21   USE wrk_nemo       ! work arrays
22   USE timing         ! Timing
23
24   IMPLICIT NONE
25   PRIVATE
26
27   PUBLIC   tra_adv_mle        ! routine called in traadv.F90
28   PUBLIC   tra_adv_mle_init   ! routine called in traadv.F90
29
30   !                                       !!* namelist namtra_adv_mle *
31   LOGICAL, PUBLIC ::   ln_mle              ! flag to activate the Mixed Layer Eddy (MLE) parameterisation
32   INTEGER         ::   nn_mle              ! MLE type: =0 standard Fox-Kemper ; =1 new formulation
33   INTEGER         ::   nn_mld_uv           ! space interpolation of MLD at u- & v-pts (0=min,1=averaged,2=max)
34   INTEGER         ::   nn_conv             ! =1 no MLE in case of convection ; =0 always MLE
35   REAL(wp)        ::   rn_ce               ! MLE coefficient
36   !                                        ! parameters used in nn_mle = 0 case
37   REAL(wp)        ::   rn_lf                  ! typical scale of mixed layer front
38   REAL(wp)        ::   rn_time                ! time scale for mixing momentum across the mixed layer
39   !                                        ! parameters used in nn_mle = 1 case
40   REAL(wp)        ::   rn_lat                 ! reference latitude for a 5 km scale of ML front
41   REAL(wp)        ::   rn_rho_c_mle           ! Density criterion for definition of MLD used by FK
42
43   REAL(wp) ::   r5_21 = 5.e0 / 21.e0   ! factor used in mle streamfunction computation
44   REAL(wp) ::   rb_c                   ! ML buoyancy criteria = g rho_c /rau0 where rho_c is defined in zdfmld
45   REAL(wp) ::   rc_f                   ! MLE coefficient (= rn_ce / (5 km * fo) ) in nn_mle=1 case
46
47   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   rfu, rfv   ! modified Coriolis parameter (f+tau) at u- & v-pts
48   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   r1_ft      ! inverse of the modified Coriolis parameter at t-pts
49
50   !! * Substitutions
51#  include "vectopt_loop_substitute.h90"
52   !!----------------------------------------------------------------------
53   !! NEMO/OPA 4.0 , NEMO Consortium (2015)
54   !! $Id$
55   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
56   !!----------------------------------------------------------------------
57CONTAINS
58
59   SUBROUTINE tra_adv_mle( kt, kit000, pu, pv, pw, cdtype )
60      !!----------------------------------------------------------------------
61      !!                  ***  ROUTINE adv_mle  ***
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  : - (pun,pvn,pwn) 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      CHARACTER(len=3)                , INTENT(in   ) ::   cdtype     ! =TRA or TRC (tracer indicator)
84      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pu         ! in : 3 ocean transport components
85      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pv         ! out: same 3  transport components
86      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pw         !   increased by the MLE induced transport
87      !
88      INTEGER  ::   ji, jj, jk   ! dummy loop indices
89      INTEGER  ::   ikmax        ! temporary integer
90      REAL(wp) ::   zcuw, zmuw   ! local scalar
91      REAL(wp) ::   zcvw, zmvw   !   -      -
92      REAL(wp) ::   zc                                     !   -      -
93      !
94      INTEGER  ::   ii, ij, ik              ! local integers
95      INTEGER, DIMENSION(3) ::   ilocu      !
96      INTEGER, DIMENSION(2) ::   ilocs      !
97      REAL(wp), POINTER, DIMENSION(:,:  ) :: zpsim_u, zpsim_v, zmld, zbm, zhu, zhv, zn2, zLf_NH, zLf_MH
98      REAL(wp), POINTER, DIMENSION(:,:,:) :: zpsi_uw, zpsi_vw
99      INTEGER, POINTER, DIMENSION(:,:) :: inml_mle
100      !!----------------------------------------------------------------------
101      !
102      IF( nn_timing == 1 )  CALL timing_start('tra_adv_mle')
103      CALL wrk_alloc( jpi, jpj, zpsim_u, zpsim_v, zmld, zbm, zhu, zhv, zn2, zLf_NH, zLf_MH)
104      CALL wrk_alloc( jpi, jpj, jpk, zpsi_uw, zpsi_vw)
105      CALL wrk_alloc( jpi, jpj, inml_mle)
106      !
107      !                                      !==  MLD used for MLE  ==!
108      !                                                ! compute from the 10m density to deal with the diurnal cycle
109      inml_mle(:,:) = mbkt(:,:) + 1                    ! init. to number of ocean w-level (T-level + 1)
110      DO jk = jpkm1, nlb10, -1                         ! from the bottom to nlb10 (10m)
111         DO jj = 1, jpj
112            DO ji = 1, jpi                             ! index of the w-level at the ML based
113               IF( rhop(ji,jj,jk) > rhop(ji,jj,nla10) + rn_rho_c_mle )   inml_mle(ji,jj) = jk      ! Mixed layer
114            END DO
115         END DO
116      END DO
117      ikmax = MIN( MAXVAL( inml_mle(:,:) ), jpkm1 )                  ! max level of the computation
118      !
119      !
120      zmld(:,:) = 0._wp                      !==   Horizontal shape of the MLE  ==!
121      zbm (:,:) = 0._wp
122      zn2 (:,:) = 0._wp
123      DO jk = 1, ikmax                                 ! MLD and mean buoyancy and N2 over the mixed layer
124         DO jj = 1, jpj
125            DO ji = 1, jpi
126               zc = e3t_n(ji,jj,jk) * REAL( MIN( MAX( 0, inml_mle(ji,jj)-jk ) , 1  )  )    ! zc being 0 outside the ML t-points
127               zmld(ji,jj) = zmld(ji,jj) + zc
128               zbm (ji,jj) = zbm (ji,jj) + zc * (rau0 - rhop(ji,jj,jk) ) * r1_rau0
129               zn2 (ji,jj) = zn2 (ji,jj) + zc * (rn2(ji,jj,jk)+rn2(ji,jj,jk+1))*0.5_wp
130            END DO
131         END DO
132      END DO
133
134      SELECT CASE( nn_mld_uv )                         ! MLD at u- & v-pts
135      CASE ( 0 )                                               != min of the 2 neighbour MLDs
136         DO jj = 1, jpjm1
137            DO ji = 1, fs_jpim1   ! vector opt.
138               zhu(ji,jj) = MIN( zmld(ji+1,jj), zmld(ji,jj) )
139               zhv(ji,jj) = MIN( zmld(ji,jj+1), zmld(ji,jj) )
140            END DO
141         END DO
142      CASE ( 1 )                                               != average of the 2 neighbour MLDs
143         DO jj = 1, jpjm1
144            DO ji = 1, fs_jpim1   ! vector opt.
145               zhu(ji,jj) = ( zmld(ji+1,jj) + zmld(ji,jj) ) * 0.5_wp
146               zhv(ji,jj) = ( zmld(ji,jj+1) + zmld(ji,jj) ) * 0.5_wp
147            END DO
148         END DO
149      CASE ( 2 )                                               != max of the 2 neighbour MLDs
150         DO jj = 1, jpjm1
151            DO ji = 1, fs_jpim1   ! vector opt.
152               zhu(ji,jj) = MAX( zmld(ji+1,jj), zmld(ji,jj) )
153               zhv(ji,jj) = MAX( zmld(ji,jj+1), zmld(ji,jj) )
154            END DO
155         END DO
156      END SELECT
157      !                                                ! convert density into buoyancy
158      zbm(:,:) = + grav * zbm(:,:) / MAX( e3t_n(:,:,1), zmld(:,:) )
159      !
160      !
161      !                                      !==  Magnitude of the MLE stream function  ==!
162      !
163      !                 di[bm]  Ds
164      ! Psi = Ce  H^2 ---------------- e2u  mu(z)   where fu Lf = MAX( fu*rn_fl , (Db H)^1/2 )
165      !                  e1u   Lf fu                      and the e2u for the "transport"
166      !                                                      (not *e3u as divided by e3u at the end)
167      !
168      IF( nn_mle == 0 ) THEN           ! Fox-Kemper et al. 2010 formulation
169         DO jj = 1, jpjm1
170            DO ji = 1, fs_jpim1   ! vector opt.
171               zpsim_u(ji,jj) = rn_ce * zhu(ji,jj) * zhu(ji,jj)  * e2_e1u(ji,jj)                                            &
172                  &           * ( zbm(ji+1,jj) - zbm(ji,jj) ) * MIN( 111.e3_wp , e1u(ji,jj) )   &
173                  &           / (  MAX( rn_lf * rfu(ji,jj) , SQRT( rb_c * zhu(ji,jj) ) )   )
174                  !
175               zpsim_v(ji,jj) = rn_ce * zhv(ji,jj) * zhv(ji,jj)  * e1_e2v(ji,jj)                                            &
176                  &           * ( zbm(ji,jj+1) - zbm(ji,jj) ) * MIN( 111.e3_wp , e2v(ji,jj) )   &
177                  &           / (  MAX( rn_lf * rfv(ji,jj) , SQRT( rb_c * zhv(ji,jj) ) )   )
178            END DO
179         END DO
180         !
181      ELSEIF( nn_mle == 1 ) THEN       ! New formulation (Lf = 5km fo/ff with fo=Coriolis parameter at latitude rn_lat)
182         DO jj = 1, jpjm1
183            DO ji = 1, fs_jpim1   ! vector opt.
184               zpsim_u(ji,jj) = rc_f *   zhu(ji,jj)   * zhu(ji,jj)   * e2_e1u(ji,jj)               &
185                  &                  * ( zbm(ji+1,jj) - zbm(ji,jj) ) * MIN( 111.e3_wp , e1u(ji,jj) )
186                  !
187               zpsim_v(ji,jj) = rc_f *   zhv(ji,jj)   * zhv(ji,jj)   * e1_e2v(ji,jj)               &
188                  &                  * ( zbm(ji,jj+1) - zbm(ji,jj) ) * MIN( 111.e3_wp , e2v(ji,jj) )
189            END DO
190         END DO
191      ENDIF
192      !
193      IF( nn_conv == 1 ) THEN              ! No MLE in case of convection
194         DO jj = 1, jpjm1
195            DO ji = 1, fs_jpim1   ! vector opt.
196               IF( MIN( zn2(ji,jj) , zn2(ji+1,jj) ) < 0._wp )   zpsim_u(ji,jj) = 0._wp
197               IF( MIN( zn2(ji,jj) , zn2(ji,jj+1) ) < 0._wp )   zpsim_v(ji,jj) = 0._wp
198            END DO
199         END DO
200      ENDIF
201      !
202      !                                      !==  structure function value at uw- and vw-points  ==!
203      DO jj = 1, jpjm1
204         DO ji = 1, fs_jpim1   ! vector opt.
205            zhu(ji,jj) = 1._wp / zhu(ji,jj)                   ! hu --> 1/hu
206            zhv(ji,jj) = 1._wp / zhv(ji,jj)
207         END DO
208      END DO
209      !
210      zpsi_uw(:,:,:) = 0._wp
211      zpsi_vw(:,:,:) = 0._wp
212      !
213      DO jk = 2, ikmax                                ! start from 2 : surface value = 0
214         DO jj = 1, jpjm1
215            DO ji = 1, fs_jpim1   ! vector opt.
216               zcuw = 1._wp - ( gdepw_n(ji+1,jj,jk) + gdepw_n(ji,jj,jk) ) * zhu(ji,jj)
217               zcvw = 1._wp - ( gdepw_n(ji,jj+1,jk) + gdepw_n(ji,jj,jk) ) * zhv(ji,jj)
218               zcuw = zcuw * zcuw
219               zcvw = zcvw * zcvw
220               zmuw = MAX(  0._wp , ( 1._wp - zcuw ) * ( 1._wp + r5_21 * zcuw )  )
221               zmvw = MAX(  0._wp , ( 1._wp - zcvw ) * ( 1._wp + r5_21 * zcvw )  )
222               !
223               zpsi_uw(ji,jj,jk) = zpsim_u(ji,jj) * zmuw * umask(ji,jj,jk)
224               zpsi_vw(ji,jj,jk) = zpsim_v(ji,jj) * zmvw * vmask(ji,jj,jk)
225            END DO
226         END DO
227      END DO
228      !
229      !                                      !==  transport increased by the MLE induced transport ==!
230      DO jk = 1, ikmax
231         DO jj = 1, jpjm1                          ! CAUTION pu,pv must be defined at row/column i=1 / j=1
232            DO ji = 1, fs_jpim1   ! vector opt.
233               pu(ji,jj,jk) = pu(ji,jj,jk) + ( zpsi_uw(ji,jj,jk) - zpsi_uw(ji,jj,jk+1) )
234               pv(ji,jj,jk) = pv(ji,jj,jk) + ( zpsi_vw(ji,jj,jk) - zpsi_vw(ji,jj,jk+1) )
235            END DO
236         END DO
237         DO jj = 2, jpjm1
238            DO ji = fs_2, fs_jpim1   ! vector opt.
239               pw(ji,jj,jk) = pw(ji,jj,jk) - ( zpsi_uw(ji,jj,jk) - zpsi_uw(ji-1,jj,jk)   &
240                  &                          + zpsi_vw(ji,jj,jk) - zpsi_vw(ji,jj-1,jk) )
241            END DO
242         END DO
243      END DO
244
245      IF( cdtype == 'TRA') THEN              !==  outputs  ==!
246         !
247         zLf_NH(:,:) = SQRT( rb_c * zmld(:,:) ) * r1_ft(:,:)      ! Lf = N H / f
248         CALL iom_put( "Lf_NHpf" , zLf_NH  )    ! Lf = N H / f
249         !
250         ! divide by cross distance to give streamfunction with dimensions m^2/s
251         DO jk = 1, ikmax+1
252            zpsi_uw(:,:,jk) = zpsi_uw(:,:,jk) * r1_e2u(:,:)
253            zpsi_vw(:,:,jk) = zpsi_vw(:,:,jk) * r1_e1v(:,:)
254         END DO
255         CALL iom_put( "psiu_mle", zpsi_uw )    ! i-mle streamfunction
256         CALL iom_put( "psiv_mle", zpsi_vw )    ! j-mle streamfunction
257      ENDIF
258      CALL wrk_dealloc( jpi, jpj, zpsim_u, zpsim_v, zmld, zbm, zhu, zhv, zn2, zLf_NH, zLf_MH)
259      CALL wrk_dealloc( jpi, jpj, jpk, zpsi_uw, zpsi_vw)
260      CALL wrk_dealloc( jpi, jpj, inml_mle)
261
262      IF( nn_timing == 1 )  CALL timing_stop('tra_adv_mle')
263      !
264   END SUBROUTINE tra_adv_mle
265
266
267   SUBROUTINE tra_adv_mle_init
268      !!---------------------------------------------------------------------
269      !!                  ***  ROUTINE tra_adv_mle_init  ***
270      !!
271      !! ** Purpose :   Control the consistency between namelist options for
272      !!              tracer advection schemes and set nadv
273      !!----------------------------------------------------------------------
274      INTEGER  ::   ji, jj, jk   ! dummy loop indices
275      INTEGER  ::   ierr
276      INTEGER ::    ios                 ! Local integer output status for namelist read
277      REAL(wp) ::   z1_t2, zfu, zfv                                !    -         -
278      !
279      NAMELIST/namtra_adv_mle/ ln_mle , nn_mle, rn_ce, rn_lf, rn_time, rn_lat, nn_mld_uv, nn_conv, rn_rho_c_mle
280      !!----------------------------------------------------------------------
281
282      REWIND( numnam_ref )              ! Namelist namtra_adv_mle in reference namelist : Tracer advection scheme
283      READ  ( numnam_ref, namtra_adv_mle, IOSTAT = ios, ERR = 901)
284901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtra_adv_mle in reference namelist', lwp )
285
286      REWIND( numnam_cfg )              ! Namelist namtra_adv_mle in configuration namelist : Tracer advection scheme
287      READ  ( numnam_cfg, namtra_adv_mle, IOSTAT = ios, ERR = 902 )
288902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtra_adv_mle in configuration namelist', lwp )
289      IF(lwm) WRITE ( numond, namtra_adv_mle )
290
291      IF(lwp) THEN                     ! Namelist print
292         WRITE(numout,*)
293         WRITE(numout,*) 'tra_adv_mle_init : mixed layer eddy (MLE) advection acting on tracers'
294         WRITE(numout,*) '~~~~~~~~~~~~~~~~'
295         WRITE(numout,*) '   Namelist namtra_adv_mle : mixed layer eddy advection on tracers'
296         WRITE(numout,*) '      use mixed layer eddy (MLE, i.e. Fox-Kemper param) (T/F)      ln_mle    = ', ln_mle
297         WRITE(numout,*) '      MLE type: =0 standard Fox-Kemper ; =1 new formulation        nn_mle    = ', nn_mle
298         WRITE(numout,*) '      magnitude of the MLE (typical value: 0.06 to 0.08)           rn_ce     = ', rn_ce
299         WRITE(numout,*) '      scale of ML front (ML radius of deformation) (rn_mle=0)      rn_lf     = ', rn_lf, 'm'
300         WRITE(numout,*) '      maximum time scale of MLE                    (rn_mle=0)      rn_time   = ', rn_time, 's'
301         WRITE(numout,*) '      reference latitude (degrees) of MLE coef.    (rn_mle=1)      rn_lat    = ', rn_lat, 'deg'
302         WRITE(numout,*) '      space interp. of MLD at u-(v-)pts (0=min,1=averaged,2=max)   nn_mld_uv = ', nn_mld_uv
303         WRITE(numout,*) '      =1 no MLE in case of convection ; =0 always MLE              nn_conv   = ', nn_conv
304         WRITE(numout,*) '      Density difference used to define ML for FK              rn_rho_c_mle  = ', rn_rho_c_mle
305      ENDIF
306      !
307      IF(lwp) THEN
308         WRITE(numout,*)
309         IF( ln_mle ) THEN
310            WRITE(numout,*) '      ===>>   Mixed Layer Eddy induced transport added to tracer advection'
311            IF( nn_mle == 0 )   WRITE(numout,*) '              Fox-Kemper et al 2010 formulation'
312            IF( nn_mle == 1 )   WRITE(numout,*) '              New formulation'
313         ELSE
314            WRITE(numout,*) '      ===>>   Mixed Layer Eddy parametrisation NOT used'
315         ENDIF
316      ENDIF
317      !
318      IF( ln_mle ) THEN                ! MLE initialisation
319         !
320         rb_c = grav * rn_rho_c_mle /rau0        ! Mixed Layer buoyancy criteria
321         IF(lwp) WRITE(numout,*)
322         IF(lwp) WRITE(numout,*) '      ML buoyancy criteria = ', rb_c, ' m/s2 '
323         IF(lwp) WRITE(numout,*) '      associated ML density criteria defined in zdfmxl = ', rho_c, 'kg/m3'
324         !
325         IF( nn_mle == 0 ) THEN           ! MLE array allocation & initialisation
326            ALLOCATE( rfu(jpi,jpj) , rfv(jpi,jpj) , STAT= ierr )
327            IF( ierr /= 0 )   CALL ctl_stop( 'tra_adv_mle_init: failed to allocate arrays' )
328            z1_t2 = 1._wp / ( rn_time * rn_time )
329            DO jj = 2, jpj                           ! "coriolis+ time^-1" at u- & v-points
330               DO ji = fs_2, jpi   ! vector opt.
331                  zfu = ( ff_f(ji,jj) + ff_f(ji,jj-1) ) * 0.5_wp
332                  zfv = ( ff_f(ji,jj) + ff_f(ji-1,jj) ) * 0.5_wp
333                  rfu(ji,jj) = SQRT(  zfu * zfu + z1_t2 )
334                  rfv(ji,jj) = SQRT(  zfv * zfv + z1_t2 )
335               END DO
336            END DO
337            CALL lbc_lnk( rfu, 'U', 1. )   ;   CALL lbc_lnk( rfv, 'V', 1. )
338            !
339         ELSEIF( nn_mle == 1 ) THEN           ! MLE array allocation & initialisation
340            rc_f = rn_ce / (  5.e3_wp * 2._wp * omega * SIN( rad * rn_lat )  )
341            !
342         ENDIF
343         !
344         !                                ! 1/(f^2+tau^2)^1/2 at t-point (needed in both nn_mle case)
345         ALLOCATE( r1_ft(jpi,jpj) , STAT= ierr )
346         IF( ierr /= 0 )   CALL ctl_stop( 'tra_adv_mle_init: failed to allocate r1_ft array' )
347         !
348         z1_t2 = 1._wp / ( rn_time * rn_time )
349         r1_ft(:,:) = 1._wp / SQRT(  ff_t(:,:) * ff_t(:,:) + z1_t2  )
350         !
351      ENDIF
352      !
353   END SUBROUTINE tra_adv_mle_init
354
355   !!==============================================================================
356END MODULE traadv_mle
Note: See TracBrowser for help on using the repository browser.