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