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