[3] | 1 | MODULE ldfslp |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE ldfslp *** |
---|
| 4 | !! Ocean physics: slopes of neutral surfaces |
---|
| 5 | !!====================================================================== |
---|
[1515] | 6 | !! History : OPA ! 1994-12 (G. Madec, M. Imbard) Original code |
---|
| 7 | !! 8.0 ! 1997-06 (G. Madec) optimization, lbc |
---|
| 8 | !! 8.1 ! 1999-10 (A. Jouzeau) NEW profile in the mixed layer |
---|
[2528] | 9 | !! NEMO 1.0 ! 2002-10 (G. Madec) Free form, F90 |
---|
| 10 | !! - ! 2005-10 (A. Beckmann) correction for s-coordinates |
---|
| 11 | !! 3.3 ! 2010-10 (G. Nurser, C. Harris, G. Madec) add Griffies operator |
---|
| 12 | !! - ! 2010-11 (F. Dupond, G. Madec) bug correction in slopes just below the ML |
---|
[1515] | 13 | !!---------------------------------------------------------------------- |
---|
[3] | 14 | #if defined key_ldfslp || defined key_esopa |
---|
| 15 | !!---------------------------------------------------------------------- |
---|
| 16 | !! 'key_ldfslp' Rotation of lateral mixing tensor |
---|
| 17 | !!---------------------------------------------------------------------- |
---|
[3625] | 18 | !! ldf_slp_grif : calculates the triads of isoneutral slopes (Griffies operator) |
---|
| 19 | !! ldf_slp : calculates the slopes of neutral surface (Madec operator) |
---|
| 20 | !! ldf_slp_mxl : calculates the slopes at the base of the mixed layer (Madec operator) |
---|
| 21 | !! ldf_slp_init : initialization of the slopes computation |
---|
[3] | 22 | !!---------------------------------------------------------------------- |
---|
[3848] | 23 | USE oce ! ocean dynamics and tracers |
---|
| 24 | USE dom_oce ! ocean space and time domain |
---|
| 25 | USE ldftra_oce ! lateral diffusion: traceur |
---|
| 26 | USE ldfdyn_oce ! lateral diffusion: dynamics |
---|
| 27 | USE phycst ! physical constants |
---|
| 28 | USE zdfmxl ! mixed layer depth |
---|
| 29 | USE eosbn2 ! equation of states |
---|
[4990] | 30 | ! |
---|
| 31 | USE in_out_manager ! I/O manager |
---|
[3848] | 32 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
| 33 | USE prtctl ! Print control |
---|
| 34 | USE wrk_nemo ! work arrays |
---|
| 35 | USE timing ! Timing |
---|
| 36 | USE lib_fortran ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) |
---|
[3] | 37 | |
---|
| 38 | IMPLICIT NONE |
---|
| 39 | PRIVATE |
---|
| 40 | |
---|
[2528] | 41 | PUBLIC ldf_slp ! routine called by step.F90 |
---|
| 42 | PUBLIC ldf_slp_grif ! routine called by step.F90 |
---|
| 43 | PUBLIC ldf_slp_init ! routine called by opa.F90 |
---|
[3] | 44 | |
---|
[2528] | 45 | LOGICAL , PUBLIC, PARAMETER :: lk_ldfslp = .TRUE. !: slopes flag |
---|
| 46 | ! !! Madec operator |
---|
[2715] | 47 | ! Arrays allocated in ldf_slp_init() routine once we know whether we're using the Griffies or Madec operator |
---|
| 48 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: uslp, wslpi !: i_slope at U- and W-points |
---|
| 49 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: vslp, wslpj !: j-slope at V- and W-points |
---|
| 50 | ! !! Griffies operator |
---|
| 51 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: wslp2 !: wslp**2 from Griffies quarter cells |
---|
[3294] | 52 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:,:) :: triadi_g, triadj_g !: skew flux slopes relative to geopotentials |
---|
[2715] | 53 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:,:) :: triadi , triadj !: isoneutral slopes relative to model-coordinate |
---|
[3] | 54 | |
---|
[2528] | 55 | ! !! Madec operator |
---|
[2715] | 56 | ! Arrays allocated in ldf_slp_init() routine once we know whether we're using the Griffies or Madec operator |
---|
| 57 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: omlmask ! mask of the surface mixed layer at T-pt |
---|
| 58 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: uslpml, wslpiml ! i_slope at U- and W-points just below the mixed layer |
---|
| 59 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: vslpml, wslpjml ! j_slope at V- and W-points just below the mixed layer |
---|
[2528] | 60 | |
---|
| 61 | REAL(wp) :: repsln = 1.e-25_wp ! tiny value used as minium of di(rho), dj(rho) and dk(rho) |
---|
| 62 | |
---|
[3] | 63 | !! * Substitutions |
---|
| 64 | # include "domzgr_substitute.h90" |
---|
[2528] | 65 | # include "ldftra_substitute.h90" |
---|
| 66 | # include "ldfeiv_substitute.h90" |
---|
[3] | 67 | # include "vectopt_loop_substitute.h90" |
---|
| 68 | !!---------------------------------------------------------------------- |
---|
[2715] | 69 | !! NEMO/OPA 4.0 , NEMO Consortium (2011) |
---|
[1156] | 70 | !! $Id$ |
---|
[2528] | 71 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
[3] | 72 | !!---------------------------------------------------------------------- |
---|
| 73 | CONTAINS |
---|
| 74 | |
---|
| 75 | SUBROUTINE ldf_slp( kt, prd, pn2 ) |
---|
| 76 | !!---------------------------------------------------------------------- |
---|
| 77 | !! *** ROUTINE ldf_slp *** |
---|
[3294] | 78 | !! |
---|
[1515] | 79 | !! ** Purpose : Compute the slopes of neutral surface (slope of isopycnal |
---|
[2528] | 80 | !! surfaces referenced locally) (ln_traldf_iso=T). |
---|
[3294] | 81 | !! |
---|
| 82 | !! ** Method : The slope in the i-direction is computed at U- and |
---|
| 83 | !! W-points (uslp, wslpi) and the slope in the j-direction is |
---|
[3] | 84 | !! computed at V- and W-points (vslp, wslpj). |
---|
| 85 | !! They are bounded by 1/100 over the whole ocean, and within the |
---|
| 86 | !! surface layer they are bounded by the distance to the surface |
---|
| 87 | !! ( slope<= depth/l where l is the length scale of horizontal |
---|
| 88 | !! diffusion (here, aht=2000m2/s ==> l=20km with a typical velocity |
---|
| 89 | !! of 10cm/s) |
---|
| 90 | !! A horizontal shapiro filter is applied to the slopes |
---|
[461] | 91 | !! ln_sco=T, s-coordinate, add to the previously computed slopes |
---|
[3] | 92 | !! the slope of the model level surface. |
---|
| 93 | !! macro-tasked on horizontal slab (jk-loop) (2, jpk-1) |
---|
| 94 | !! [slopes already set to zero at level 1, and to zero or the ocean |
---|
[461] | 95 | !! bottom slope (ln_sco=T) at level jpk in inildf] |
---|
[3] | 96 | !! |
---|
[3294] | 97 | !! ** Action : - uslp, wslpi, and vslp, wslpj, the i- and j-slopes |
---|
[3] | 98 | !! of now neutral surfaces at u-, w- and v- w-points, resp. |
---|
[1515] | 99 | !!---------------------------------------------------------------------- |
---|
[2715] | 100 | INTEGER , INTENT(in) :: kt ! ocean time-step index |
---|
| 101 | REAL(wp), INTENT(in), DIMENSION(:,:,:) :: prd ! in situ density |
---|
| 102 | REAL(wp), INTENT(in), DIMENSION(:,:,:) :: pn2 ! Brunt-Vaisala frequency (locally ref.) |
---|
[1515] | 103 | !! |
---|
| 104 | INTEGER :: ji , jj , jk ! dummy loop indices |
---|
| 105 | INTEGER :: ii0, ii1, iku ! temporary integer |
---|
| 106 | INTEGER :: ij0, ij1, ikv ! temporary integer |
---|
[3294] | 107 | REAL(wp) :: zeps, zm1_g, zm1_2g, z1_16, zcofw ! local scalars |
---|
[2528] | 108 | REAL(wp) :: zci, zfi, zau, zbu, zai, zbi ! - - |
---|
| 109 | REAL(wp) :: zcj, zfj, zav, zbv, zaj, zbj ! - - |
---|
| 110 | REAL(wp) :: zck, zfk, zbw ! - - |
---|
[4990] | 111 | REAL(wp) :: zdepv, zdepu ! - - |
---|
[3294] | 112 | REAL(wp), POINTER, DIMENSION(:,:,:) :: zwz, zww |
---|
| 113 | REAL(wp), POINTER, DIMENSION(:,:,:) :: zdzr |
---|
| 114 | REAL(wp), POINTER, DIMENSION(:,:,:) :: zgru, zgrv |
---|
[4990] | 115 | REAL(wp), POINTER, DIMENSION(:,: ) :: zhmlpu, zhmlpv |
---|
[3] | 116 | !!---------------------------------------------------------------------- |
---|
[3294] | 117 | ! |
---|
| 118 | IF( nn_timing == 1 ) CALL timing_start('ldf_slp') |
---|
| 119 | ! |
---|
| 120 | CALL wrk_alloc( jpi,jpj,jpk, zwz, zww, zdzr, zgru, zgrv ) |
---|
[4990] | 121 | CALL wrk_alloc( jpi,jpj, zhmlpu, zhmlpv ) |
---|
[2715] | 122 | |
---|
[4488] | 123 | IF ( ln_traldf_iso .OR. ln_dynldf_iso ) THEN |
---|
| 124 | |
---|
| 125 | zeps = 1.e-20_wp !== Local constant initialization ==! |
---|
| 126 | z1_16 = 1.0_wp / 16._wp |
---|
| 127 | zm1_g = -1.0_wp / grav |
---|
| 128 | zm1_2g = -0.5_wp / grav |
---|
| 129 | ! |
---|
| 130 | zww(:,:,:) = 0._wp |
---|
| 131 | zwz(:,:,:) = 0._wp |
---|
| 132 | ! |
---|
| 133 | DO jk = 1, jpk !== i- & j-gradient of density ==! |
---|
| 134 | DO jj = 1, jpjm1 |
---|
| 135 | DO ji = 1, fs_jpim1 ! vector opt. |
---|
| 136 | zgru(ji,jj,jk) = umask(ji,jj,jk) * ( prd(ji+1,jj ,jk) - prd(ji,jj,jk) ) |
---|
| 137 | zgrv(ji,jj,jk) = vmask(ji,jj,jk) * ( prd(ji ,jj+1,jk) - prd(ji,jj,jk) ) |
---|
| 138 | END DO |
---|
[3] | 139 | END DO |
---|
| 140 | END DO |
---|
[4488] | 141 | IF( ln_zps ) THEN ! partial steps correction at the bottom ocean level |
---|
| 142 | DO jj = 1, jpjm1 |
---|
| 143 | DO ji = 1, jpim1 |
---|
[4990] | 144 | ! IF should be useless check zpshde (PM) |
---|
[5098] | 145 | zgru(ji,jj,mbku(ji,jj)) = gru(ji,jj) |
---|
| 146 | zgrv(ji,jj,mbkv(ji,jj)) = grv(ji,jj) |
---|
| 147 | END DO |
---|
| 148 | END DO |
---|
| 149 | ENDIF |
---|
| 150 | IF( ln_zps .AND. ln_isfcav ) THEN ! partial steps correction at the bottom ocean level |
---|
| 151 | DO jj = 1, jpjm1 |
---|
| 152 | DO ji = 1, jpim1 |
---|
[4990] | 153 | IF ( miku(ji,jj) > 1 ) zgru(ji,jj,miku(ji,jj)) = grui(ji,jj) |
---|
| 154 | IF ( mikv(ji,jj) > 1 ) zgrv(ji,jj,mikv(ji,jj)) = grvi(ji,jj) |
---|
[4488] | 155 | END DO |
---|
[461] | 156 | END DO |
---|
[4488] | 157 | ENDIF |
---|
| 158 | ! |
---|
[5098] | 159 | !== Local vertical density gradient at T-point == ! (evaluated from N^2) |
---|
| 160 | ! surface initialisation |
---|
| 161 | zdzr(:,:,1) = 0._wp |
---|
| 162 | ! interior value |
---|
| 163 | DO jk = 2, jpkm1 |
---|
[4488] | 164 | ! ! zdzr = d/dz(prd)= - ( prd ) / grav * mk(pn2) -- at t point |
---|
| 165 | ! ! trick: tmask(ik ) = 0 => all pn2 = 0 => zdzr = 0 |
---|
| 166 | ! ! else tmask(ik+1) = 0 => pn2(ik+1) = 0 => zdzr divides by 1 |
---|
| 167 | ! ! umask(ik+1) /= 0 => all pn2 /= 0 => zdzr divides by 2 |
---|
| 168 | ! ! NB: 1/(tmask+1) = (1-.5*tmask) substitute a / by a * ==> faster |
---|
| 169 | zdzr(:,:,jk) = zm1_g * ( prd(:,:,jk) + 1._wp ) & |
---|
| 170 | & * ( pn2(:,:,jk) + pn2(:,:,jk+1) ) * ( 1._wp - 0.5_wp * tmask(:,:,jk+1) ) |
---|
[3] | 171 | END DO |
---|
[5098] | 172 | IF ( ln_isfcav ) THEN |
---|
| 173 | ! if isf need to overwrite the interior value at at the first ocean point |
---|
| 174 | DO jj = 1, jpjm1 |
---|
| 175 | DO ji = 1, jpim1 |
---|
| 176 | zdzr(ji,jj,1:mikt(ji,jj)) = 0._wp |
---|
| 177 | END DO |
---|
[4990] | 178 | END DO |
---|
[5098] | 179 | END IF |
---|
[4488] | 180 | ! |
---|
| 181 | ! !== Slopes just below the mixed layer ==! |
---|
| 182 | CALL ldf_slp_mxl( prd, pn2, zgru, zgrv, zdzr ) ! output: uslpml, vslpml, wslpiml, wslpjml |
---|
[2389] | 183 | |
---|
[3294] | 184 | |
---|
[4488] | 185 | ! I. slopes at u and v point | uslp = d/di( prd ) / d/dz( prd ) |
---|
| 186 | ! =========================== | vslp = d/dj( prd ) / d/dz( prd ) |
---|
| 187 | ! |
---|
[4990] | 188 | DO jj = 2, jpjm1 |
---|
| 189 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 190 | IF (miku(ji,jj) .GT. miku(ji+1,jj)) zhmlpu(ji,jj) = hmlpt(ji ,jj) |
---|
| 191 | IF (miku(ji,jj) .LT. miku(ji+1,jj)) zhmlpu(ji,jj) = hmlpt(ji+1,jj) |
---|
| 192 | IF (miku(ji,jj) .EQ. miku(ji+1,jj)) zhmlpu(ji,jj) = MAX(hmlpt(ji ,jj), hmlpt(ji+1,jj)) |
---|
| 193 | IF (mikv(ji,jj) .GT. miku(ji,jj+1)) zhmlpv(ji,jj) = hmlpt(ji ,jj) |
---|
| 194 | IF (mikv(ji,jj) .LT. miku(ji,jj+1)) zhmlpv(ji,jj) = hmlpt(ji,jj+1) |
---|
| 195 | IF (mikv(ji,jj) .EQ. miku(ji,jj+1)) zhmlpv(ji,jj) = MAX(hmlpt(ji,jj), hmlpt(ji,jj+1)) |
---|
| 196 | ENDDO |
---|
| 197 | ENDDO |
---|
[4488] | 198 | DO jk = 2, jpkm1 !* Slopes at u and v points |
---|
| 199 | DO jj = 2, jpjm1 |
---|
| 200 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 201 | ! ! horizontal and vertical density gradient at u- and v-points |
---|
| 202 | zau = zgru(ji,jj,jk) / e1u(ji,jj) |
---|
| 203 | zav = zgrv(ji,jj,jk) / e2v(ji,jj) |
---|
| 204 | zbu = 0.5_wp * ( zdzr(ji,jj,jk) + zdzr(ji+1,jj ,jk) ) |
---|
| 205 | zbv = 0.5_wp * ( zdzr(ji,jj,jk) + zdzr(ji ,jj+1,jk) ) |
---|
| 206 | ! ! bound the slopes: abs(zw.)<= 1/100 and zb..<0 |
---|
| 207 | ! ! + kxz max= ah slope max =< e1 e3 /(pi**2 2 dt) |
---|
| 208 | zbu = MIN( zbu, -100._wp* ABS( zau ) , -7.e+3_wp/fse3u(ji,jj,jk)* ABS( zau ) ) |
---|
| 209 | zbv = MIN( zbv, -100._wp* ABS( zav ) , -7.e+3_wp/fse3v(ji,jj,jk)* ABS( zav ) ) |
---|
| 210 | ! ! uslp and vslp output in zwz and zww, resp. |
---|
[4990] | 211 | zfi = MAX( omlmask(ji,jj,jk), omlmask(ji+1,jj,jk) ) |
---|
| 212 | zfj = MAX( omlmask(ji,jj,jk), omlmask(ji,jj+1,jk) ) |
---|
| 213 | ! thickness of water column between surface and level k at u/v point |
---|
| 214 | zdepu = 0.5_wp * (( fsdept(ji,jj,jk) + fsdept(ji+1,jj ,jk) ) & |
---|
| 215 | - 2 * MAX( risfdep(ji,jj), risfdep(ji+1,jj ) ) & |
---|
| 216 | - fse3u(ji,jj,miku(ji,jj)) ) |
---|
| 217 | zdepv = 0.5_wp * (( fsdept(ji,jj,jk) + fsdept(ji ,jj+1,jk) ) & |
---|
| 218 | - 2 * MAX( risfdep(ji,jj), risfdep(ji,jj+1) ) & |
---|
| 219 | - fse3v(ji,jj,mikv(ji,jj)) ) |
---|
| 220 | zwz(ji,jj,jk) = ( 1. - zfi) * zau / ( zbu - zeps ) & |
---|
| 221 | & + zfi * uslpml(ji,jj) & |
---|
| 222 | & * zdepu / MAX( zhmlpu(ji,jj), 5._wp ) |
---|
[5098] | 223 | zwz(ji,jj,jk) = zwz(ji,jj,jk) * wumask(ji,jj,jk) |
---|
[4990] | 224 | zww(ji,jj,jk) = ( 1. - zfj) * zav / ( zbv - zeps ) & |
---|
| 225 | & + zfj * vslpml(ji,jj) & |
---|
| 226 | & * zdepv / MAX( zhmlpv(ji,jj), 5._wp ) |
---|
[5098] | 227 | zww(ji,jj,jk) = zww(ji,jj,jk) * wvmask(ji,jj,jk) |
---|
[4990] | 228 | |
---|
| 229 | |
---|
[2528] | 230 | !!gm modif to suppress omlmask.... (as in Griffies case) |
---|
[4488] | 231 | ! ! ! jk must be >= ML level for zf=1. otherwise zf=0. |
---|
| 232 | ! zfi = REAL( 1 - 1/(1 + jk / MAX( nmln(ji+1,jj), nmln(ji,jj) ) ), wp ) |
---|
| 233 | ! zfj = REAL( 1 - 1/(1 + jk / MAX( nmln(ji,jj+1), nmln(ji,jj) ) ), wp ) |
---|
| 234 | ! zci = 0.5 * ( fsdept(ji+1,jj,jk)+fsdept(ji,jj,jk) ) / MAX( hmlpt(ji,jj), hmlpt(ji+1,jj), 10. ) ) |
---|
| 235 | ! zcj = 0.5 * ( fsdept(ji,jj+1,jk)+fsdept(ji,jj,jk) ) / MAX( hmlpt(ji,jj), hmlpt(ji,jj+1), 10. ) ) |
---|
| 236 | ! zwz(ji,jj,jk) = ( zfi * zai / ( zbi - zeps ) + ( 1._wp - zfi ) * wslpiml(ji,jj) * zci ) * tmask(ji,jj,jk) |
---|
| 237 | ! zww(ji,jj,jk) = ( zfj * zaj / ( zbj - zeps ) + ( 1._wp - zfj ) * wslpjml(ji,jj) * zcj ) * tmask(ji,jj,jk) |
---|
[2528] | 238 | !!gm end modif |
---|
[4488] | 239 | END DO |
---|
[3] | 240 | END DO |
---|
| 241 | END DO |
---|
[4488] | 242 | CALL lbc_lnk( zwz, 'U', -1. ) ; CALL lbc_lnk( zww, 'V', -1. ) ! lateral boundary conditions |
---|
| 243 | ! |
---|
| 244 | ! !* horizontal Shapiro filter |
---|
| 245 | DO jk = 2, jpkm1 |
---|
| 246 | DO jj = 2, jpjm1, MAX(1, jpj-3) ! rows jj=2 and =jpjm1 only |
---|
| 247 | DO ji = 2, jpim1 |
---|
| 248 | uslp(ji,jj,jk) = z1_16 * ( zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk) & |
---|
| 249 | & + zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk) & |
---|
| 250 | & + 2.*( zwz(ji ,jj-1,jk) + zwz(ji-1,jj ,jk) & |
---|
| 251 | & + zwz(ji+1,jj ,jk) + zwz(ji ,jj+1,jk) ) & |
---|
| 252 | & + 4.* zwz(ji ,jj ,jk) ) |
---|
| 253 | vslp(ji,jj,jk) = z1_16 * ( zww(ji-1,jj-1,jk) + zww(ji+1,jj-1,jk) & |
---|
| 254 | & + zww(ji-1,jj+1,jk) + zww(ji+1,jj+1,jk) & |
---|
| 255 | & + 2.*( zww(ji ,jj-1,jk) + zww(ji-1,jj ,jk) & |
---|
| 256 | & + zww(ji+1,jj ,jk) + zww(ji ,jj+1,jk) ) & |
---|
| 257 | & + 4.* zww(ji,jj ,jk) ) |
---|
| 258 | END DO |
---|
[3] | 259 | END DO |
---|
[4488] | 260 | DO jj = 3, jpj-2 ! other rows |
---|
| 261 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 262 | uslp(ji,jj,jk) = z1_16 * ( zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk) & |
---|
| 263 | & + zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk) & |
---|
| 264 | & + 2.*( zwz(ji ,jj-1,jk) + zwz(ji-1,jj ,jk) & |
---|
| 265 | & + zwz(ji+1,jj ,jk) + zwz(ji ,jj+1,jk) ) & |
---|
| 266 | & + 4.* zwz(ji ,jj ,jk) ) |
---|
| 267 | vslp(ji,jj,jk) = z1_16 * ( zww(ji-1,jj-1,jk) + zww(ji+1,jj-1,jk) & |
---|
| 268 | & + zww(ji-1,jj+1,jk) + zww(ji+1,jj+1,jk) & |
---|
| 269 | & + 2.*( zww(ji ,jj-1,jk) + zww(ji-1,jj ,jk) & |
---|
| 270 | & + zww(ji+1,jj ,jk) + zww(ji ,jj+1,jk) ) & |
---|
| 271 | & + 4.* zww(ji,jj ,jk) ) |
---|
| 272 | END DO |
---|
[3] | 273 | END DO |
---|
[4488] | 274 | ! !* decrease along coastal boundaries |
---|
| 275 | DO jj = 2, jpjm1 |
---|
| 276 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 277 | uslp(ji,jj,jk) = uslp(ji,jj,jk) * ( umask(ji,jj+1,jk) + umask(ji,jj-1,jk ) ) * 0.5_wp & |
---|
[4990] | 278 | & * ( umask(ji,jj ,jk) + umask(ji,jj ,jk+1) ) * 0.5_wp & |
---|
| 279 | & * umask(ji,jj,jk-1) !* umask(ji,jj,jk) * umask(ji,jj,jk+1) |
---|
[4488] | 280 | vslp(ji,jj,jk) = vslp(ji,jj,jk) * ( vmask(ji+1,jj,jk) + vmask(ji-1,jj,jk ) ) * 0.5_wp & |
---|
[4990] | 281 | & * ( vmask(ji ,jj,jk) + vmask(ji ,jj,jk+1) ) * 0.5_wp & |
---|
| 282 | & * vmask(ji,jj,jk-1) !* vmask(ji,jj,jk) * vmask(ji,jj,jk+1) |
---|
[4488] | 283 | END DO |
---|
[3] | 284 | END DO |
---|
| 285 | END DO |
---|
| 286 | |
---|
| 287 | |
---|
[4488] | 288 | ! II. slopes at w point | wslpi = mij( d/di( prd ) / d/dz( prd ) |
---|
| 289 | ! =========================== | wslpj = mij( d/dj( prd ) / d/dz( prd ) |
---|
| 290 | ! |
---|
| 291 | DO jk = 2, jpkm1 |
---|
| 292 | DO jj = 2, jpjm1 |
---|
| 293 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 294 | ! !* Local vertical density gradient evaluated from N^2 |
---|
[4990] | 295 | zbw = zm1_2g * pn2 (ji,jj,jk) * ( prd (ji,jj,jk) + prd (ji,jj,jk-1) + 2. ) * tmask(ji,jj,jk) * tmask(ji,jj,jk-1) |
---|
[4488] | 296 | ! !* Slopes at w point |
---|
| 297 | ! ! i- & j-gradient of density at w-points |
---|
| 298 | zci = MAX( umask(ji-1,jj,jk ) + umask(ji,jj,jk ) & |
---|
| 299 | & + umask(ji-1,jj,jk-1) + umask(ji,jj,jk-1) , zeps ) * e1t(ji,jj) |
---|
| 300 | zcj = MAX( vmask(ji,jj-1,jk ) + vmask(ji,jj,jk-1) & |
---|
| 301 | & + vmask(ji,jj-1,jk-1) + vmask(ji,jj,jk ) , zeps ) * e2t(ji,jj) |
---|
| 302 | zai = ( zgru (ji-1,jj,jk ) + zgru (ji,jj,jk-1) & |
---|
[4990] | 303 | & + zgru (ji-1,jj,jk-1) + zgru (ji,jj,jk ) ) / zci |
---|
[4488] | 304 | zaj = ( zgrv (ji,jj-1,jk ) + zgrv (ji,jj,jk-1) & |
---|
[4990] | 305 | & + zgrv (ji,jj-1,jk-1) + zgrv (ji,jj,jk ) ) / zcj |
---|
[4488] | 306 | ! ! bound the slopes: abs(zw.)<= 1/100 and zb..<0. |
---|
| 307 | ! ! + kxz max= ah slope max =< e1 e3 /(pi**2 2 dt) |
---|
| 308 | zbi = MIN( zbw ,- 100._wp* ABS( zai ) , -7.e+3_wp/fse3w(ji,jj,jk)* ABS( zai ) ) |
---|
| 309 | zbj = MIN( zbw , -100._wp* ABS( zaj ) , -7.e+3_wp/fse3w(ji,jj,jk)* ABS( zaj ) ) |
---|
| 310 | ! ! wslpi and wslpj with ML flattening (output in zwz and zww, resp.) |
---|
[4990] | 311 | zfk = MAX( omlmask(ji,jj,jk), omlmask(ji,jj,jk-1) ) ! zfk=1 in the ML otherwise zfk=0 |
---|
| 312 | zck = ( fsdepw(ji,jj,jk) - fsdepw(ji,jj,mikt(ji,jj) ) ) / MAX( hmlp(ji,jj), 10._wp ) |
---|
| 313 | zwz(ji,jj,jk) = ( zai / ( zbi - zeps ) * ( 1._wp - zfk ) & |
---|
[5098] | 314 | & + zck * wslpiml(ji,jj) * zfk ) * wmask(ji,jj,jk) |
---|
[4990] | 315 | zww(ji,jj,jk) = ( zaj / ( zbj - zeps ) * ( 1._wp - zfk ) & |
---|
[5098] | 316 | & + zck * wslpjml(ji,jj) * zfk ) * wmask(ji,jj,jk) |
---|
[2528] | 317 | |
---|
| 318 | !!gm modif to suppress omlmask.... (as in Griffies operator) |
---|
[4488] | 319 | ! ! ! jk must be >= ML level for zfk=1. otherwise zfk=0. |
---|
| 320 | ! zfk = REAL( 1 - 1/(1 + jk / nmln(ji+1,jj)), wp ) |
---|
| 321 | ! zck = fsdepw(ji,jj,jk) / MAX( hmlp(ji,jj), 10. ) |
---|
| 322 | ! zwz(ji,jj,jk) = ( zfk * zai / ( zbi - zeps ) + ( 1._wp - zfk ) * wslpiml(ji,jj) * zck ) * tmask(ji,jj,jk) |
---|
| 323 | ! zww(ji,jj,jk) = ( zfk * zaj / ( zbj - zeps ) + ( 1._wp - zfk ) * wslpjml(ji,jj) * zck ) * tmask(ji,jj,jk) |
---|
[2528] | 324 | !!gm end modif |
---|
[4488] | 325 | END DO |
---|
[3] | 326 | END DO |
---|
| 327 | END DO |
---|
[4488] | 328 | CALL lbc_lnk( zwz, 'T', -1. ) ; CALL lbc_lnk( zww, 'T', -1. ) ! lateral boundary conditions |
---|
| 329 | ! |
---|
| 330 | ! !* horizontal Shapiro filter |
---|
| 331 | DO jk = 2, jpkm1 |
---|
| 332 | DO jj = 2, jpjm1, MAX(1, jpj-3) ! rows jj=2 and =jpjm1 only |
---|
| 333 | DO ji = 2, jpim1 |
---|
| 334 | zcofw = tmask(ji,jj,jk) * z1_16 |
---|
| 335 | wslpi(ji,jj,jk) = ( zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk) & |
---|
| 336 | & + zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk) & |
---|
| 337 | & + 2.*( zwz(ji ,jj-1,jk) + zwz(ji-1,jj ,jk) & |
---|
| 338 | & + zwz(ji+1,jj ,jk) + zwz(ji ,jj+1,jk) ) & |
---|
| 339 | & + 4.* zwz(ji ,jj ,jk) ) * zcofw |
---|
[3] | 340 | |
---|
[4488] | 341 | wslpj(ji,jj,jk) = ( zww(ji-1,jj-1,jk) + zww(ji+1,jj-1,jk) & |
---|
| 342 | & + zww(ji-1,jj+1,jk) + zww(ji+1,jj+1,jk) & |
---|
| 343 | & + 2.*( zww(ji ,jj-1,jk) + zww(ji-1,jj ,jk) & |
---|
| 344 | & + zww(ji+1,jj ,jk) + zww(ji ,jj+1,jk) ) & |
---|
| 345 | & + 4.* zww(ji ,jj ,jk) ) * zcofw |
---|
| 346 | END DO |
---|
[3] | 347 | END DO |
---|
[4488] | 348 | DO jj = 3, jpj-2 ! other rows |
---|
| 349 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 350 | zcofw = tmask(ji,jj,jk) * z1_16 |
---|
| 351 | wslpi(ji,jj,jk) = ( zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk) & |
---|
| 352 | & + zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk) & |
---|
| 353 | & + 2.*( zwz(ji ,jj-1,jk) + zwz(ji-1,jj ,jk) & |
---|
| 354 | & + zwz(ji+1,jj ,jk) + zwz(ji ,jj+1,jk) ) & |
---|
| 355 | & + 4.* zwz(ji ,jj ,jk) ) * zcofw |
---|
[3] | 356 | |
---|
[4488] | 357 | wslpj(ji,jj,jk) = ( zww(ji-1,jj-1,jk) + zww(ji+1,jj-1,jk) & |
---|
| 358 | & + zww(ji-1,jj+1,jk) + zww(ji+1,jj+1,jk) & |
---|
| 359 | & + 2.*( zww(ji ,jj-1,jk) + zww(ji-1,jj ,jk) & |
---|
| 360 | & + zww(ji+1,jj ,jk) + zww(ji ,jj+1,jk) ) & |
---|
| 361 | & + 4.* zww(ji ,jj ,jk) ) * zcofw |
---|
| 362 | END DO |
---|
[3] | 363 | END DO |
---|
[4488] | 364 | ! !* decrease along coastal boundaries |
---|
| 365 | DO jj = 2, jpjm1 |
---|
| 366 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 367 | zck = ( umask(ji,jj,jk) + umask(ji-1,jj,jk) ) & |
---|
| 368 | & * ( vmask(ji,jj,jk) + vmask(ji,jj-1,jk) ) * 0.25 |
---|
[5098] | 369 | wslpi(ji,jj,jk) = wslpi(ji,jj,jk) * zck * wmask(ji,jj,jk) |
---|
| 370 | wslpj(ji,jj,jk) = wslpj(ji,jj,jk) * zck * wmask(ji,jj,jk) |
---|
[4488] | 371 | END DO |
---|
[3] | 372 | END DO |
---|
| 373 | END DO |
---|
[3294] | 374 | |
---|
[4488] | 375 | ! III. Specific grid points |
---|
| 376 | ! =========================== |
---|
[1515] | 377 | ! |
---|
[4488] | 378 | IF( cp_cfg == "orca" .AND. jp_cfg == 4 ) THEN ! ORCA_R4 configuration: horizontal diffusion in specific area |
---|
| 379 | ! ! Gibraltar Strait |
---|
| 380 | ij0 = 50 ; ij1 = 53 |
---|
| 381 | ii0 = 69 ; ii1 = 71 ; uslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0._wp |
---|
| 382 | ij0 = 51 ; ij1 = 53 |
---|
| 383 | ii0 = 68 ; ii1 = 71 ; vslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0._wp |
---|
| 384 | ii0 = 69 ; ii1 = 71 ; wslpi( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0._wp |
---|
| 385 | ii0 = 69 ; ii1 = 71 ; wslpj( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0._wp |
---|
| 386 | ! |
---|
| 387 | ! ! Mediterrannean Sea |
---|
| 388 | ij0 = 49 ; ij1 = 56 |
---|
| 389 | ii0 = 71 ; ii1 = 90 ; uslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0._wp |
---|
| 390 | ij0 = 50 ; ij1 = 56 |
---|
| 391 | ii0 = 70 ; ii1 = 90 ; vslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0._wp |
---|
| 392 | ii0 = 71 ; ii1 = 90 ; wslpi( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0._wp |
---|
| 393 | ii0 = 71 ; ii1 = 90 ; wslpj( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0._wp |
---|
| 394 | ENDIF |
---|
[3] | 395 | |
---|
[3294] | 396 | |
---|
[4488] | 397 | ! IV. Lateral boundary conditions |
---|
| 398 | ! =============================== |
---|
| 399 | CALL lbc_lnk( uslp , 'U', -1. ) ; CALL lbc_lnk( vslp , 'V', -1. ) |
---|
| 400 | CALL lbc_lnk( wslpi, 'W', -1. ) ; CALL lbc_lnk( wslpj, 'W', -1. ) |
---|
[3] | 401 | |
---|
[1515] | 402 | |
---|
[4488] | 403 | IF(ln_ctl) THEN |
---|
| 404 | CALL prt_ctl(tab3d_1=uslp , clinfo1=' slp - u : ', tab3d_2=vslp, clinfo2=' v : ', kdim=jpk) |
---|
| 405 | CALL prt_ctl(tab3d_1=wslpi, clinfo1=' slp - wi: ', tab3d_2=wslpj, clinfo2=' wj: ', kdim=jpk) |
---|
| 406 | ENDIF |
---|
| 407 | ! |
---|
| 408 | |
---|
| 409 | ELSEIF ( lk_vvl ) THEN |
---|
| 410 | |
---|
| 411 | IF(lwp) THEN |
---|
| 412 | WRITE(numout,*) ' Horizontal mixing in s-coordinate: slope = slope of s-surfaces' |
---|
| 413 | ENDIF |
---|
| 414 | |
---|
| 415 | ! geopotential diffusion in s-coordinates on tracers and/or momentum |
---|
| 416 | ! The slopes of s-surfaces are computed at each time step due to vvl |
---|
| 417 | ! The slopes for momentum diffusion are i- or j- averaged of those on tracers |
---|
| 418 | |
---|
| 419 | ! set the slope of diffusion to the slope of s-surfaces |
---|
| 420 | ! ( c a u t i o n : minus sign as fsdep has positive value ) |
---|
[4990] | 421 | DO jj = 2, jpjm1 |
---|
| 422 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 423 | uslp(ji,jj,1) = -1./e1u(ji,jj) * ( fsdept_b(ji+1,jj,1) - fsdept_b(ji ,jj ,1) ) * umask(ji,jj,1) |
---|
| 424 | vslp(ji,jj,1) = -1./e2v(ji,jj) * ( fsdept_b(ji,jj+1,1) - fsdept_b(ji ,jj ,1) ) * vmask(ji,jj,1) |
---|
| 425 | wslpi(ji,jj,1) = -1./e1t(ji,jj) * ( fsdepw_b(ji+1,jj,1) - fsdepw_b(ji-1,jj,1) ) * tmask(ji,jj,1) * 0.5 |
---|
| 426 | wslpj(ji,jj,1) = -1./e2t(ji,jj) * ( fsdepw_b(ji,jj+1,1) - fsdepw_b(ji,jj-1,1) ) * tmask(ji,jj,1) * 0.5 |
---|
| 427 | END DO |
---|
| 428 | END DO |
---|
| 429 | |
---|
| 430 | DO jk = 2, jpk |
---|
[4488] | 431 | DO jj = 2, jpjm1 |
---|
| 432 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 433 | uslp(ji,jj,jk) = -1./e1u(ji,jj) * ( fsdept_b(ji+1,jj,jk) - fsdept_b(ji ,jj ,jk) ) * umask(ji,jj,jk) |
---|
| 434 | vslp(ji,jj,jk) = -1./e2v(ji,jj) * ( fsdept_b(ji,jj+1,jk) - fsdept_b(ji ,jj ,jk) ) * vmask(ji,jj,jk) |
---|
[4990] | 435 | wslpi(ji,jj,jk) = -1./e1t(ji,jj) * ( fsdepw_b(ji+1,jj,jk) - fsdepw_b(ji-1,jj,jk) ) & |
---|
[5098] | 436 | & * wmask(ji,jj,jk) * 0.5 |
---|
[4990] | 437 | wslpj(ji,jj,jk) = -1./e2t(ji,jj) * ( fsdepw_b(ji,jj+1,jk) - fsdepw_b(ji,jj-1,jk) ) & |
---|
[5098] | 438 | & * wmask(ji,jj,jk) * 0.5 |
---|
[4488] | 439 | END DO |
---|
| 440 | END DO |
---|
| 441 | END DO |
---|
| 442 | |
---|
| 443 | ! Lateral boundary conditions on the slopes |
---|
| 444 | CALL lbc_lnk( uslp , 'U', -1. ) ; CALL lbc_lnk( vslp , 'V', -1. ) |
---|
| 445 | CALL lbc_lnk( wslpi, 'W', -1. ) ; CALL lbc_lnk( wslpj, 'W', -1. ) |
---|
| 446 | |
---|
| 447 | if( kt == nit000 ) then |
---|
| 448 | IF(lwp) WRITE(numout,*) ' max slop: u',SQRT( MAXVAL(uslp*uslp)), ' v ', SQRT(MAXVAL(vslp)), & |
---|
| 449 | & ' wi', sqrt(MAXVAL(wslpi)), ' wj', sqrt(MAXVAL(wslpj)) |
---|
| 450 | endif |
---|
| 451 | |
---|
| 452 | IF(ln_ctl) THEN |
---|
| 453 | CALL prt_ctl(tab3d_1=uslp , clinfo1=' slp - u : ', tab3d_2=vslp, clinfo2=' v : ', kdim=jpk) |
---|
| 454 | CALL prt_ctl(tab3d_1=wslpi, clinfo1=' slp - wi: ', tab3d_2=wslpj, clinfo2=' wj: ', kdim=jpk) |
---|
| 455 | ENDIF |
---|
| 456 | |
---|
[49] | 457 | ENDIF |
---|
[4488] | 458 | |
---|
[3294] | 459 | CALL wrk_dealloc( jpi,jpj,jpk, zwz, zww, zdzr, zgru, zgrv ) |
---|
[4990] | 460 | CALL wrk_dealloc( jpi,jpj, zhmlpu, zhmlpv) |
---|
[2715] | 461 | ! |
---|
[3294] | 462 | IF( nn_timing == 1 ) CALL timing_stop('ldf_slp') |
---|
| 463 | ! |
---|
[3] | 464 | END SUBROUTINE ldf_slp |
---|
| 465 | |
---|
[3294] | 466 | |
---|
[2528] | 467 | SUBROUTINE ldf_slp_grif ( kt ) |
---|
| 468 | !!---------------------------------------------------------------------- |
---|
| 469 | !! *** ROUTINE ldf_slp_grif *** |
---|
| 470 | !! |
---|
| 471 | !! ** Purpose : Compute the squared slopes of neutral surfaces (slope |
---|
| 472 | !! of iso-pycnal surfaces referenced locally) (ln_traldf_grif=T) |
---|
[3294] | 473 | !! at W-points using the Griffies quarter-cells. |
---|
[2528] | 474 | !! |
---|
[3294] | 475 | !! ** Method : calculates alpha and beta at T-points |
---|
[2528] | 476 | !! |
---|
| 477 | !! ** Action : - triadi_g, triadj_g T-pts i- and j-slope triads relative to geopot. (used for eiv) |
---|
| 478 | !! - triadi , triadj T-pts i- and j-slope triads relative to model-coordinate |
---|
| 479 | !! - wslp2 squared slope of neutral surfaces at w-points. |
---|
| 480 | !!---------------------------------------------------------------------- |
---|
[3294] | 481 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
| 482 | !! |
---|
[2528] | 483 | INTEGER :: ji, jj, jk, jl, ip, jp, kp ! dummy loop indices |
---|
[3294] | 484 | INTEGER :: iku, ikv ! local integer |
---|
| 485 | REAL(wp) :: zfacti, zfactj ! local scalars |
---|
| 486 | REAL(wp) :: znot_thru_surface ! local scalars |
---|
| 487 | REAL(wp) :: zdit, zdis, zdjt, zdjs, zdkt, zdks, zbu, zbv, zbti, zbtj |
---|
| 488 | REAL(wp) :: zdxrho_raw, zti_coord, zti_raw, zti_lim, zti_g_raw, zti_g_lim |
---|
| 489 | REAL(wp) :: zdyrho_raw, ztj_coord, ztj_raw, ztj_lim, ztj_g_raw, ztj_g_lim |
---|
[2528] | 490 | REAL(wp) :: zdzrho_raw |
---|
[3294] | 491 | REAL(wp), POINTER, DIMENSION(:,:) :: z1_mlbw |
---|
| 492 | REAL(wp), POINTER, DIMENSION(:,:,:,:) :: zdxrho , zdyrho, zdzrho ! Horizontal and vertical density gradients |
---|
| 493 | REAL(wp), POINTER, DIMENSION(:,:,:,:) :: zti_mlb, ztj_mlb ! for Griffies operator only |
---|
[2528] | 494 | !!---------------------------------------------------------------------- |
---|
[3294] | 495 | ! |
---|
| 496 | IF( nn_timing == 1 ) CALL timing_start('ldf_slp_grif') |
---|
| 497 | ! |
---|
| 498 | CALL wrk_alloc( jpi,jpj, z1_mlbw ) |
---|
| 499 | CALL wrk_alloc( jpi,jpj,jpk,2, zdxrho , zdyrho, zdzrho, klstart = 0 ) |
---|
| 500 | CALL wrk_alloc( jpi,jpj, 2,2, zti_mlb, ztj_mlb, kkstart = 0, klstart = 0 ) |
---|
| 501 | ! |
---|
[2528] | 502 | !--------------------------------! |
---|
| 503 | ! Some preliminary calculation ! |
---|
| 504 | !--------------------------------! |
---|
| 505 | ! |
---|
[3294] | 506 | DO jl = 0, 1 !== unmasked before density i- j-, k-gradients ==! |
---|
| 507 | ! |
---|
| 508 | ip = jl ; jp = jl ! guaranteed nonzero gradients ( absolute value larger than repsln) |
---|
| 509 | DO jk = 1, jpkm1 ! done each pair of triad |
---|
| 510 | DO jj = 1, jpjm1 ! NB: not masked ==> a minimum value is set |
---|
| 511 | DO ji = 1, fs_jpim1 ! vector opt. |
---|
| 512 | zdit = ( tsb(ji+1,jj,jk,jp_tem) - tsb(ji,jj,jk,jp_tem) ) ! i-gradient of T & S at u-point |
---|
| 513 | zdis = ( tsb(ji+1,jj,jk,jp_sal) - tsb(ji,jj,jk,jp_sal) ) |
---|
| 514 | zdjt = ( tsb(ji,jj+1,jk,jp_tem) - tsb(ji,jj,jk,jp_tem) ) ! j-gradient of T & S at v-point |
---|
| 515 | zdjs = ( tsb(ji,jj+1,jk,jp_sal) - tsb(ji,jj,jk,jp_sal) ) |
---|
[4990] | 516 | zdxrho_raw = ( - rab_b(ji+ip,jj ,jk,jp_tem) * zdit + rab_b(ji+ip,jj ,jk,jp_sal) * zdis ) / e1u(ji,jj) |
---|
| 517 | zdyrho_raw = ( - rab_b(ji ,jj+jp,jk,jp_tem) * zdjt + rab_b(ji ,jj+jp,jk,jp_sal) * zdjs ) / e2v(ji,jj) |
---|
| 518 | zdxrho(ji+ip,jj ,jk,1-ip) = SIGN( MAX( repsln, ABS( zdxrho_raw ) ), zdxrho_raw ) ! keep the sign |
---|
[3294] | 519 | zdyrho(ji ,jj+jp,jk,1-jp) = SIGN( MAX( repsln, ABS( zdyrho_raw ) ), zdyrho_raw ) |
---|
| 520 | END DO |
---|
[2528] | 521 | END DO |
---|
| 522 | END DO |
---|
[3294] | 523 | ! |
---|
[4990] | 524 | IF( ln_zps .AND. l_grad_zps ) THEN ! partial steps: correction of i- & j-grad on bottom |
---|
[3294] | 525 | DO jj = 1, jpjm1 |
---|
| 526 | DO ji = 1, jpim1 |
---|
| 527 | iku = mbku(ji,jj) ; ikv = mbkv(ji,jj) ! last ocean level (u- & v-points) |
---|
| 528 | zdit = gtsu(ji,jj,jp_tem) ; zdjt = gtsv(ji,jj,jp_tem) ! i- & j-gradient of Temperature |
---|
| 529 | zdis = gtsu(ji,jj,jp_sal) ; zdjs = gtsv(ji,jj,jp_sal) ! i- & j-gradient of Salinity |
---|
[4990] | 530 | zdxrho_raw = ( - rab_b(ji+ip,jj ,iku,jp_tem) * zdit + rab_b(ji+ip,jj ,iku,jp_sal) * zdis ) / e1u(ji,jj) |
---|
| 531 | zdyrho_raw = ( - rab_b(ji ,jj+jp,ikv,jp_tem) * zdjt + rab_b(ji ,jj+jp,ikv,jp_sal) * zdjs ) / e2v(ji,jj) |
---|
[3294] | 532 | zdxrho(ji+ip,jj ,iku,1-ip) = SIGN( MAX( repsln, ABS( zdxrho_raw ) ), zdxrho_raw ) ! keep the sign |
---|
| 533 | zdyrho(ji ,jj+jp,ikv,1-jp) = SIGN( MAX( repsln, ABS( zdyrho_raw ) ), zdyrho_raw ) |
---|
[2528] | 534 | END DO |
---|
| 535 | END DO |
---|
[3294] | 536 | ENDIF |
---|
| 537 | ! |
---|
[2528] | 538 | END DO |
---|
[3294] | 539 | |
---|
| 540 | DO kp = 0, 1 !== unmasked before density i- j-, k-gradients ==! |
---|
| 541 | DO jk = 1, jpkm1 ! done each pair of triad |
---|
| 542 | DO jj = 1, jpj ! NB: not masked ==> a minimum value is set |
---|
| 543 | DO ji = 1, jpi ! vector opt. |
---|
| 544 | IF( jk+kp > 1 ) THEN ! k-gradient of T & S a jk+kp |
---|
| 545 | zdkt = ( tsb(ji,jj,jk+kp-1,jp_tem) - tsb(ji,jj,jk+kp,jp_tem) ) |
---|
| 546 | zdks = ( tsb(ji,jj,jk+kp-1,jp_sal) - tsb(ji,jj,jk+kp,jp_sal) ) |
---|
| 547 | ELSE |
---|
| 548 | zdkt = 0._wp ! 1st level gradient set to zero |
---|
| 549 | zdks = 0._wp |
---|
| 550 | ENDIF |
---|
[4990] | 551 | zdzrho_raw = ( - rab_b(ji,jj,jk,jp_tem) * zdkt + rab_b(ji,jj,jk,jp_sal) * zdks ) / fse3w(ji,jj,jk+kp) |
---|
| 552 | zdzrho(ji,jj,jk,kp) = - MIN( - repsln, zdzrho_raw ) ! force zdzrho >= repsln |
---|
[3294] | 553 | END DO |
---|
[2528] | 554 | END DO |
---|
| 555 | END DO |
---|
| 556 | END DO |
---|
| 557 | ! |
---|
[3294] | 558 | DO jj = 1, jpj !== Reciprocal depth of the w-point below ML base ==! |
---|
[2528] | 559 | DO ji = 1, jpi |
---|
| 560 | jk = MIN( nmln(ji,jj), mbkt(ji,jj) ) + 1 ! MIN in case ML depth is the ocean depth |
---|
| 561 | z1_mlbw(ji,jj) = 1._wp / fsdepw(ji,jj,jk) |
---|
| 562 | END DO |
---|
| 563 | END DO |
---|
| 564 | ! |
---|
[3294] | 565 | ! !== intialisations to zero ==! |
---|
[2528] | 566 | ! |
---|
[3294] | 567 | wslp2 (:,:,:) = 0._wp ! wslp2 will be cumulated 3D field set to zero |
---|
| 568 | triadi_g(:,:,1,:,:) = 0._wp ; triadi_g(:,:,jpk,:,:) = 0._wp ! set surface and bottom slope to zero |
---|
[2528] | 569 | triadj_g(:,:,1,:,:) = 0._wp ; triadj_g(:,:,jpk,:,:) = 0._wp |
---|
[3294] | 570 | !!gm _iso set to zero missing |
---|
| 571 | triadi (:,:,1,:,:) = 0._wp ; triadj (:,:,jpk,:,:) = 0._wp ! set surface and bottom slope to zero |
---|
| 572 | triadj (:,:,1,:,:) = 0._wp ; triadj (:,:,jpk,:,:) = 0._wp |
---|
| 573 | |
---|
[2528] | 574 | !-------------------------------------! |
---|
| 575 | ! Triads just below the Mixed Layer ! |
---|
| 576 | !-------------------------------------! |
---|
| 577 | ! |
---|
[3294] | 578 | DO jl = 0, 1 ! calculate slope of the 4 triads immediately ONE level below mixed-layer base |
---|
| 579 | DO kp = 0, 1 ! with only the slope-max limit and MASKED |
---|
[2528] | 580 | DO jj = 1, jpjm1 |
---|
| 581 | DO ji = 1, fs_jpim1 |
---|
| 582 | ip = jl ; jp = jl |
---|
[3845] | 583 | ! |
---|
| 584 | jk = nmln(ji+ip,jj) + 1 |
---|
| 585 | IF( jk .GT. mbkt(ji+ip,jj) ) THEN !ML reaches bottom |
---|
| 586 | zti_mlb(ji+ip,jj ,1-ip,kp) = 0.0_wp |
---|
| 587 | ELSE |
---|
| 588 | ! Add s-coordinate slope at t-points (do this by *subtracting* gradient of depth) |
---|
| 589 | zti_g_raw = ( zdxrho(ji+ip,jj,jk-kp,1-ip) / zdzrho(ji+ip,jj,jk-kp,kp) & |
---|
| 590 | & - ( fsdept(ji+1,jj,jk-kp) - fsdept(ji,jj,jk-kp) ) / e1u(ji,jj) ) * umask(ji,jj,jk) |
---|
| 591 | zti_mlb(ji+ip,jj ,1-ip,kp) = SIGN( MIN( rn_slpmax, ABS( zti_g_raw ) ), zti_g_raw ) |
---|
| 592 | ENDIF |
---|
| 593 | ! |
---|
| 594 | jk = nmln(ji,jj+jp) + 1 |
---|
| 595 | IF( jk .GT. mbkt(ji,jj+jp) ) THEN !ML reaches bottom |
---|
| 596 | ztj_mlb(ji ,jj+jp,1-jp,kp) = 0.0_wp |
---|
| 597 | ELSE |
---|
| 598 | ztj_g_raw = ( zdyrho(ji,jj+jp,jk-kp,1-jp) / zdzrho(ji,jj+jp,jk-kp,kp) & |
---|
| 599 | & - ( fsdept(ji,jj+1,jk-kp) - fsdept(ji,jj,jk-kp) ) / e2v(ji,jj) ) * vmask(ji,jj,jk) |
---|
| 600 | ztj_mlb(ji ,jj+jp,1-jp,kp) = SIGN( MIN( rn_slpmax, ABS( ztj_g_raw ) ), ztj_g_raw ) |
---|
| 601 | ENDIF |
---|
[2528] | 602 | END DO |
---|
| 603 | END DO |
---|
| 604 | END DO |
---|
| 605 | END DO |
---|
| 606 | |
---|
| 607 | !-------------------------------------! |
---|
| 608 | ! Triads with surface limits ! |
---|
| 609 | !-------------------------------------! |
---|
| 610 | ! |
---|
[3294] | 611 | DO kp = 0, 1 ! k-index of triads |
---|
[2528] | 612 | DO jl = 0, 1 |
---|
[3294] | 613 | ip = jl ; jp = jl ! i- and j-indices of triads (i-k and j-k planes) |
---|
[2528] | 614 | DO jk = 1, jpkm1 |
---|
[3294] | 615 | ! Must mask contribution to slope from dz/dx at constant s for triads jk=1,kp=0 that poke up though ocean surface |
---|
| 616 | znot_thru_surface = REAL( 1-1/(jk+kp), wp ) !jk+kp=1,=0.; otherwise=1.0 |
---|
[2528] | 617 | DO jj = 1, jpjm1 |
---|
[3294] | 618 | DO ji = 1, fs_jpim1 ! vector opt. |
---|
[2528] | 619 | ! |
---|
| 620 | ! Calculate slope relative to geopotentials used for GM skew fluxes |
---|
[3294] | 621 | ! Add s-coordinate slope at t-points (do this by *subtracting* gradient of depth) |
---|
[2528] | 622 | ! Limit by slope *relative to geopotentials* by rn_slpmax, and mask by psi-point |
---|
| 623 | ! masked by umask taken at the level of dz(rho) |
---|
| 624 | ! |
---|
| 625 | ! raw slopes: unmasked unbounded slopes (relative to geopotential (zti_g) and model surface (zti) |
---|
| 626 | ! |
---|
| 627 | zti_raw = zdxrho(ji+ip,jj ,jk,1-ip) / zdzrho(ji+ip,jj ,jk,kp) ! unmasked |
---|
| 628 | ztj_raw = zdyrho(ji ,jj+jp,jk,1-jp) / zdzrho(ji ,jj+jp,jk,kp) |
---|
[3294] | 629 | |
---|
| 630 | ! Must mask contribution to slope for triad jk=1,kp=0 that poke up though ocean surface |
---|
| 631 | zti_coord = znot_thru_surface * ( fsdept(ji+1,jj ,jk) - fsdept(ji,jj,jk) ) / e1u(ji,jj) |
---|
| 632 | ztj_coord = znot_thru_surface * ( fsdept(ji ,jj+1,jk) - fsdept(ji,jj,jk) ) / e2v(ji,jj) ! unmasked |
---|
| 633 | zti_g_raw = zti_raw - zti_coord ! ref to geopot surfaces |
---|
| 634 | ztj_g_raw = ztj_raw - ztj_coord |
---|
[2528] | 635 | zti_g_lim = SIGN( MIN( rn_slpmax, ABS( zti_g_raw ) ), zti_g_raw ) |
---|
| 636 | ztj_g_lim = SIGN( MIN( rn_slpmax, ABS( ztj_g_raw ) ), ztj_g_raw ) |
---|
| 637 | ! |
---|
[3294] | 638 | ! Below ML use limited zti_g as is & mask |
---|
| 639 | ! Inside ML replace by linearly reducing sx_mlb towards surface & mask |
---|
[2528] | 640 | ! |
---|
| 641 | zfacti = REAL( 1 - 1/(1 + (jk+kp-1)/nmln(ji+ip,jj)), wp ) ! k index of uppermost point(s) of triad is jk+kp-1 |
---|
| 642 | zfactj = REAL( 1 - 1/(1 + (jk+kp-1)/nmln(ji,jj+jp)), wp ) ! must be .ge. nmln(ji,jj) for zfact=1 |
---|
| 643 | ! ! otherwise zfact=0 |
---|
[3294] | 644 | zti_g_lim = ( zfacti * zti_g_lim & |
---|
[2528] | 645 | & + ( 1._wp - zfacti ) * zti_mlb(ji+ip,jj,1-ip,kp) & |
---|
[3294] | 646 | & * fsdepw(ji+ip,jj,jk+kp) * z1_mlbw(ji+ip,jj) ) * umask(ji,jj,jk+kp) |
---|
| 647 | ztj_g_lim = ( zfactj * ztj_g_lim & |
---|
[2528] | 648 | & + ( 1._wp - zfactj ) * ztj_mlb(ji,jj+jp,1-jp,kp) & |
---|
[3294] | 649 | & * fsdepw(ji,jj+jp,jk+kp) * z1_mlbw(ji,jj+jp) ) * vmask(ji,jj,jk+kp) |
---|
[2528] | 650 | ! |
---|
[3294] | 651 | triadi_g(ji+ip,jj ,jk,1-ip,kp) = zti_g_lim |
---|
| 652 | triadj_g(ji ,jj+jp,jk,1-jp,kp) = ztj_g_lim |
---|
[2528] | 653 | ! |
---|
| 654 | ! Get coefficients of isoneutral diffusion tensor |
---|
| 655 | ! 1. Utilise gradients *relative* to s-coordinate, so add t-point slopes (*subtract* depth gradients) |
---|
| 656 | ! 2. We require that isoneutral diffusion gives no vertical buoyancy flux |
---|
| 657 | ! i.e. 33 term = (real slope* 31, 13 terms) |
---|
| 658 | ! To do this, retain limited sx**2 in vertical flux, but divide by real slope for 13/31 terms |
---|
| 659 | ! Equivalent to tapering A_iso = sx_limited**2/(real slope)**2 |
---|
| 660 | ! |
---|
[3294] | 661 | zti_lim = ( zti_g_lim + zti_coord ) * umask(ji,jj,jk+kp) ! remove coordinate slope => relative to coordinate surfaces |
---|
| 662 | ztj_lim = ( ztj_g_lim + ztj_coord ) * vmask(ji,jj,jk+kp) |
---|
[2528] | 663 | ! |
---|
[3294] | 664 | IF( ln_triad_iso ) THEN |
---|
| 665 | zti_raw = zti_lim**2 / zti_raw |
---|
| 666 | ztj_raw = ztj_lim**2 / ztj_raw |
---|
| 667 | zti_raw = SIGN( MIN( ABS(zti_lim), ABS( zti_raw ) ), zti_raw ) |
---|
| 668 | ztj_raw = SIGN( MIN( ABS(ztj_lim), ABS( ztj_raw ) ), ztj_raw ) |
---|
| 669 | zti_lim = zfacti * zti_lim & |
---|
| 670 | & + ( 1._wp - zfacti ) * zti_raw |
---|
| 671 | ztj_lim = zfactj * ztj_lim & |
---|
| 672 | & + ( 1._wp - zfactj ) * ztj_raw |
---|
| 673 | ENDIF |
---|
| 674 | triadi(ji+ip,jj ,jk,1-ip,kp) = zti_lim |
---|
| 675 | triadj(ji ,jj+jp,jk,1-jp,kp) = ztj_lim |
---|
| 676 | ! |
---|
[2528] | 677 | zbu = e1u(ji ,jj) * e2u(ji ,jj) * fse3u(ji ,jj,jk ) |
---|
| 678 | zbv = e1v(ji ,jj) * e2v(ji ,jj) * fse3v(ji ,jj,jk ) |
---|
| 679 | zbti = e1t(ji+ip,jj) * e2t(ji+ip,jj) * fse3w(ji+ip,jj,jk+kp) |
---|
| 680 | zbtj = e1t(ji,jj+jp) * e2t(ji,jj+jp) * fse3w(ji,jj+jp,jk+kp) |
---|
| 681 | ! |
---|
[3294] | 682 | !!gm this may inhibit vectorization on Vect Computers, and even on scalar computers.... ==> to be checked |
---|
| 683 | wslp2 (ji+ip,jj,jk+kp) = wslp2(ji+ip,jj,jk+kp) + 0.25_wp * zbu / zbti * zti_g_lim**2 ! masked |
---|
| 684 | wslp2 (ji,jj+jp,jk+kp) = wslp2(ji,jj+jp,jk+kp) + 0.25_wp * zbv / zbtj * ztj_g_lim**2 |
---|
[2528] | 685 | END DO |
---|
| 686 | END DO |
---|
| 687 | END DO |
---|
| 688 | END DO |
---|
| 689 | END DO |
---|
| 690 | ! |
---|
| 691 | wslp2(:,:,1) = 0._wp ! force the surface wslp to zero |
---|
[3294] | 692 | |
---|
[2528] | 693 | CALL lbc_lnk( wslp2, 'W', 1. ) ! lateral boundary confition on wslp2 only ==>>> gm : necessary ? to be checked |
---|
| 694 | ! |
---|
[3294] | 695 | CALL wrk_dealloc( jpi,jpj, z1_mlbw ) |
---|
| 696 | CALL wrk_dealloc( jpi,jpj,jpk,2, zdxrho , zdyrho, zdzrho, klstart = 0 ) |
---|
| 697 | CALL wrk_dealloc( jpi,jpj, 2,2, zti_mlb, ztj_mlb, kkstart = 0, klstart = 0 ) |
---|
[2715] | 698 | ! |
---|
[3294] | 699 | IF( nn_timing == 1 ) CALL timing_stop('ldf_slp_grif') |
---|
| 700 | ! |
---|
[2528] | 701 | END SUBROUTINE ldf_slp_grif |
---|
| 702 | |
---|
| 703 | |
---|
| 704 | SUBROUTINE ldf_slp_mxl( prd, pn2, p_gru, p_grv, p_dzr ) |
---|
[3] | 705 | !!---------------------------------------------------------------------- |
---|
| 706 | !! *** ROUTINE ldf_slp_mxl *** |
---|
| 707 | !! |
---|
[3294] | 708 | !! ** Purpose : Compute the slopes of iso-neutral surface just below |
---|
[1515] | 709 | !! the mixed layer. |
---|
| 710 | !! |
---|
[2528] | 711 | !! ** Method : The slope in the i-direction is computed at u- & w-points |
---|
| 712 | !! (uslpml, wslpiml) and the slope in the j-direction is computed |
---|
| 713 | !! at v- and w-points (vslpml, wslpjml) with the same bounds as |
---|
| 714 | !! in ldf_slp. |
---|
[3] | 715 | !! |
---|
[2389] | 716 | !! ** Action : uslpml, wslpiml : i- & j-slopes of neutral surfaces |
---|
[3294] | 717 | !! vslpml, wslpjml just below the mixed layer |
---|
[2389] | 718 | !! omlmask : mixed layer mask |
---|
[1515] | 719 | !!---------------------------------------------------------------------- |
---|
[2715] | 720 | REAL(wp), DIMENSION(:,:,:), INTENT(in) :: prd ! in situ density |
---|
| 721 | REAL(wp), DIMENSION(:,:,:), INTENT(in) :: pn2 ! Brunt-Vaisala frequency (locally ref.) |
---|
| 722 | REAL(wp), DIMENSION(:,:,:), INTENT(in) :: p_gru, p_grv ! i- & j-gradient of density (u- & v-pts) |
---|
| 723 | REAL(wp), DIMENSION(:,:,:), INTENT(in) :: p_dzr ! z-gradient of density (T-point) |
---|
[3] | 724 | !! |
---|
[3294] | 725 | INTEGER :: ji , jj , jk ! dummy loop indices |
---|
| 726 | INTEGER :: iku, ikv, ik, ikm1 ! local integers |
---|
[2528] | 727 | REAL(wp) :: zeps, zm1_g, zm1_2g ! local scalars |
---|
| 728 | REAL(wp) :: zci, zfi, zau, zbu, zai, zbi ! - - |
---|
| 729 | REAL(wp) :: zcj, zfj, zav, zbv, zaj, zbj ! - - |
---|
| 730 | REAL(wp) :: zck, zfk, zbw ! - - |
---|
[3] | 731 | !!---------------------------------------------------------------------- |
---|
[3294] | 732 | ! |
---|
| 733 | IF( nn_timing == 1 ) CALL timing_start('ldf_slp_mxl') |
---|
| 734 | ! |
---|
[2528] | 735 | zeps = 1.e-20_wp !== Local constant initialization ==! |
---|
| 736 | zm1_g = -1.0_wp / grav |
---|
| 737 | zm1_2g = -0.5_wp / grav |
---|
[1515] | 738 | ! |
---|
[2528] | 739 | uslpml (1,:) = 0._wp ; uslpml (jpi,:) = 0._wp |
---|
| 740 | vslpml (1,:) = 0._wp ; vslpml (jpi,:) = 0._wp |
---|
| 741 | wslpiml(1,:) = 0._wp ; wslpiml(jpi,:) = 0._wp |
---|
| 742 | wslpjml(1,:) = 0._wp ; wslpjml(jpi,:) = 0._wp |
---|
| 743 | ! |
---|
[3294] | 744 | ! !== surface mixed layer mask ! |
---|
| 745 | DO jk = 1, jpk ! =1 inside the mixed layer, =0 otherwise |
---|
[3] | 746 | DO jj = 1, jpj |
---|
| 747 | DO ji = 1, jpi |
---|
| 748 | ik = nmln(ji,jj) - 1 |
---|
[4990] | 749 | IF( jk <= ik .AND. jk >= mikt(ji,jj) ) THEN ; omlmask(ji,jj,jk) = 1._wp |
---|
[2528] | 750 | ELSE ; omlmask(ji,jj,jk) = 0._wp |
---|
[3] | 751 | ENDIF |
---|
| 752 | END DO |
---|
| 753 | END DO |
---|
| 754 | END DO |
---|
| 755 | |
---|
| 756 | |
---|
| 757 | ! Slopes of isopycnal surfaces just before bottom of mixed layer |
---|
| 758 | ! -------------------------------------------------------------- |
---|
[1515] | 759 | ! The slope are computed as in the 3D case. |
---|
| 760 | ! A key point here is the definition of the mixed layer at u- and v-points. |
---|
| 761 | ! It is assumed to be the maximum of the two neighbouring T-point mixed layer depth. |
---|
| 762 | ! Otherwise, a n2 value inside the mixed layer can be involved in the computation |
---|
| 763 | ! of the slope, resulting in a too steep diagnosed slope and thus a spurious eddy |
---|
| 764 | ! induce velocity field near the base of the mixed layer. |
---|
[3] | 765 | !----------------------------------------------------------------------- |
---|
[1515] | 766 | ! |
---|
[3] | 767 | DO jj = 2, jpjm1 |
---|
| 768 | DO ji = 2, jpim1 |
---|
[3294] | 769 | ! !== Slope at u- & v-points just below the Mixed Layer ==! |
---|
[2528] | 770 | ! |
---|
[3294] | 771 | ! !- vertical density gradient for u- and v-slopes (from dzr at T-point) |
---|
[4990] | 772 | iku = MIN( MAX( miku(ji,jj)+1, nmln(ji,jj) , nmln(ji+1,jj) ) , jpkm1 ) ! ML (MAX of T-pts, bound by jpkm1) |
---|
| 773 | ikv = MIN( MAX( mikv(ji,jj)+1, nmln(ji,jj) , nmln(ji,jj+1) ) , jpkm1 ) ! |
---|
[2528] | 774 | zbu = 0.5_wp * ( p_dzr(ji,jj,iku) + p_dzr(ji+1,jj ,iku) ) |
---|
| 775 | zbv = 0.5_wp * ( p_dzr(ji,jj,ikv) + p_dzr(ji ,jj+1,ikv) ) |
---|
[3294] | 776 | ! !- horizontal density gradient at u- & v-points |
---|
[2528] | 777 | zau = p_gru(ji,jj,iku) / e1u(ji,jj) |
---|
| 778 | zav = p_grv(ji,jj,ikv) / e2v(ji,jj) |
---|
[3294] | 779 | ! !- bound the slopes: abs(zw.)<= 1/100 and zb..<0 |
---|
| 780 | ! kxz max= ah slope max =< e1 e3 /(pi**2 2 dt) |
---|
[2772] | 781 | zbu = MIN( zbu , -100._wp* ABS( zau ) , -7.e+3_wp/fse3u(ji,jj,iku)* ABS( zau ) ) |
---|
| 782 | zbv = MIN( zbv , -100._wp* ABS( zav ) , -7.e+3_wp/fse3v(ji,jj,ikv)* ABS( zav ) ) |
---|
[3294] | 783 | ! !- Slope at u- & v-points (uslpml, vslpml) |
---|
[2772] | 784 | uslpml(ji,jj) = zau / ( zbu - zeps ) * umask(ji,jj,iku) |
---|
| 785 | vslpml(ji,jj) = zav / ( zbv - zeps ) * vmask(ji,jj,ikv) |
---|
[2528] | 786 | ! |
---|
[3294] | 787 | ! !== i- & j-slopes at w-points just below the Mixed Layer ==! |
---|
[2528] | 788 | ! |
---|
| 789 | ik = MIN( nmln(ji,jj) + 1, jpk ) |
---|
| 790 | ikm1 = MAX( 1, ik-1 ) |
---|
[3294] | 791 | ! !- vertical density gradient for w-slope (from N^2) |
---|
[2528] | 792 | zbw = zm1_2g * pn2 (ji,jj,ik) * ( prd (ji,jj,ik) + prd (ji,jj,ikm1) + 2. ) |
---|
[3294] | 793 | ! !- horizontal density i- & j-gradient at w-points |
---|
[2528] | 794 | zci = MAX( umask(ji-1,jj,ik ) + umask(ji,jj,ik ) & |
---|
[3294] | 795 | & + umask(ji-1,jj,ikm1) + umask(ji,jj,ikm1) , zeps ) * e1t(ji,jj) |
---|
[2528] | 796 | zcj = MAX( vmask(ji,jj-1,ik ) + vmask(ji,jj,ik ) & |
---|
| 797 | & + vmask(ji,jj-1,ikm1) + vmask(ji,jj,ikm1) , zeps ) * e2t(ji,jj) |
---|
| 798 | zai = ( p_gru(ji-1,jj,ik ) + p_gru(ji,jj,ik) & |
---|
| 799 | & + p_gru(ji-1,jj,ikm1) + p_gru(ji,jj,ikm1 ) ) / zci * tmask(ji,jj,ik) |
---|
| 800 | zaj = ( p_grv(ji,jj-1,ik ) + p_grv(ji,jj,ik ) & |
---|
| 801 | & + p_grv(ji,jj-1,ikm1) + p_grv(ji,jj,ikm1) ) / zcj * tmask(ji,jj,ik) |
---|
[3294] | 802 | ! !- bound the slopes: abs(zw.)<= 1/100 and zb..<0. |
---|
| 803 | ! kxz max= ah slope max =< e1 e3 /(pi**2 2 dt) |
---|
[2528] | 804 | zbi = MIN( zbw , -100._wp* ABS( zai ) , -7.e+3_wp/fse3w(ji,jj,ik)* ABS( zai ) ) |
---|
| 805 | zbj = MIN( zbw , -100._wp* ABS( zaj ) , -7.e+3_wp/fse3w(ji,jj,ik)* ABS( zaj ) ) |
---|
[3294] | 806 | ! !- i- & j-slope at w-points (wslpiml, wslpjml) |
---|
[5098] | 807 | wslpiml(ji,jj) = zai / ( zbi - zeps ) * wmask (ji,jj,ik) |
---|
| 808 | wslpjml(ji,jj) = zaj / ( zbj - zeps ) * wmask (ji,jj,ik) |
---|
[3] | 809 | END DO |
---|
| 810 | END DO |
---|
[3294] | 811 | !!gm this lbc_lnk should be useless.... |
---|
[2528] | 812 | CALL lbc_lnk( uslpml , 'U', -1. ) ; CALL lbc_lnk( vslpml , 'V', -1. ) ! lateral boundary cond. (sign change) |
---|
| 813 | CALL lbc_lnk( wslpiml, 'W', -1. ) ; CALL lbc_lnk( wslpjml, 'W', -1. ) ! lateral boundary conditions |
---|
[1515] | 814 | ! |
---|
[3294] | 815 | IF( nn_timing == 1 ) CALL timing_stop('ldf_slp_mxl') |
---|
| 816 | ! |
---|
[3] | 817 | END SUBROUTINE ldf_slp_mxl |
---|
| 818 | |
---|
| 819 | |
---|
| 820 | SUBROUTINE ldf_slp_init |
---|
| 821 | !!---------------------------------------------------------------------- |
---|
| 822 | !! *** ROUTINE ldf_slp_init *** |
---|
| 823 | !! |
---|
| 824 | !! ** Purpose : Initialization for the isopycnal slopes computation |
---|
| 825 | !! |
---|
[3294] | 826 | !! ** Method : read the nammbf namelist and check the parameter |
---|
| 827 | !! values called by tra_dmp at the first timestep (nit000) |
---|
[3] | 828 | !!---------------------------------------------------------------------- |
---|
| 829 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
[2528] | 830 | INTEGER :: ierr ! local integer |
---|
[3] | 831 | !!---------------------------------------------------------------------- |
---|
[3294] | 832 | ! |
---|
| 833 | IF( nn_timing == 1 ) CALL timing_start('ldf_slp_init') |
---|
| 834 | ! |
---|
| 835 | IF(lwp) THEN |
---|
[3] | 836 | WRITE(numout,*) |
---|
[2528] | 837 | WRITE(numout,*) 'ldf_slp_init : direction of lateral mixing' |
---|
| 838 | WRITE(numout,*) '~~~~~~~~~~~~' |
---|
[3] | 839 | ENDIF |
---|
[3294] | 840 | |
---|
[2528] | 841 | IF( ln_traldf_grif ) THEN ! Griffies operator : triad of slopes |
---|
| 842 | ALLOCATE( triadi_g(jpi,jpj,jpk,0:1,0:1) , triadj_g(jpi,jpj,jpk,0:1,0:1) , wslp2(jpi,jpj,jpk) , STAT=ierr ) |
---|
| 843 | ALLOCATE( triadi (jpi,jpj,jpk,0:1,0:1) , triadj (jpi,jpj,jpk,0:1,0:1) , STAT=ierr ) |
---|
[2715] | 844 | IF( ierr > 0 ) CALL ctl_stop( 'STOP', 'ldf_slp_init : unable to allocate Griffies operator slope' ) |
---|
[2528] | 845 | ! |
---|
| 846 | IF( ln_dynldf_iso ) CALL ctl_stop( 'ldf_slp_init: Griffies operator on momentum not supported' ) |
---|
| 847 | ! |
---|
| 848 | ELSE ! Madec operator : slopes at u-, v-, and w-points |
---|
| 849 | ALLOCATE( uslp(jpi,jpj,jpk) , vslp(jpi,jpj,jpk) , wslpi(jpi,jpj,jpk) , wslpj(jpi,jpj,jpk) , & |
---|
| 850 | & omlmask(jpi,jpj,jpk) , uslpml(jpi,jpj) , vslpml(jpi,jpj) , wslpiml(jpi,jpj) , wslpjml(jpi,jpj) , STAT=ierr ) |
---|
[2715] | 851 | IF( ierr > 0 ) CALL ctl_stop( 'STOP', 'ldf_slp_init : unable to allocate Madec operator slope ' ) |
---|
[3] | 852 | |
---|
[2528] | 853 | ! Direction of lateral diffusion (tracers and/or momentum) |
---|
| 854 | ! ------------------------------ |
---|
| 855 | uslp (:,:,:) = 0._wp ; uslpml (:,:) = 0._wp ! set the slope to zero (even in s-coordinates) |
---|
| 856 | vslp (:,:,:) = 0._wp ; vslpml (:,:) = 0._wp |
---|
| 857 | wslpi(:,:,:) = 0._wp ; wslpiml(:,:) = 0._wp |
---|
| 858 | wslpj(:,:,:) = 0._wp ; wslpjml(:,:) = 0._wp |
---|
[3] | 859 | |
---|
[5016] | 860 | IF(ln_sco .AND. (ln_traldf_hor .OR. ln_dynldf_hor )) THEN |
---|
[2715] | 861 | IF(lwp) WRITE(numout,*) ' Horizontal mixing in s-coordinate: slope = slope of s-surfaces' |
---|
[3] | 862 | |
---|
[2528] | 863 | ! geopotential diffusion in s-coordinates on tracers and/or momentum |
---|
| 864 | ! The slopes of s-surfaces are computed once (no call to ldfslp in step) |
---|
| 865 | ! The slopes for momentum diffusion are i- or j- averaged of those on tracers |
---|
[3] | 866 | |
---|
[2528] | 867 | ! set the slope of diffusion to the slope of s-surfaces |
---|
| 868 | ! ( c a u t i o n : minus sign as fsdep has positive value ) |
---|
| 869 | DO jk = 1, jpk |
---|
| 870 | DO jj = 2, jpjm1 |
---|
| 871 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[4488] | 872 | uslp (ji,jj,jk) = -1./e1u(ji,jj) * ( fsdept_b(ji+1,jj,jk) - fsdept_b(ji ,jj ,jk) ) * umask(ji,jj,jk) |
---|
| 873 | vslp (ji,jj,jk) = -1./e2v(ji,jj) * ( fsdept_b(ji,jj+1,jk) - fsdept_b(ji ,jj ,jk) ) * vmask(ji,jj,jk) |
---|
| 874 | wslpi(ji,jj,jk) = -1./e1t(ji,jj) * ( fsdepw_b(ji+1,jj,jk) - fsdepw_b(ji-1,jj,jk) ) * tmask(ji,jj,jk) * 0.5 |
---|
| 875 | wslpj(ji,jj,jk) = -1./e2t(ji,jj) * ( fsdepw_b(ji,jj+1,jk) - fsdepw_b(ji,jj-1,jk) ) * tmask(ji,jj,jk) * 0.5 |
---|
[2528] | 876 | END DO |
---|
[3] | 877 | END DO |
---|
| 878 | END DO |
---|
[2715] | 879 | CALL lbc_lnk( uslp , 'U', -1. ) ; CALL lbc_lnk( vslp , 'V', -1. ) ! Lateral boundary conditions |
---|
| 880 | CALL lbc_lnk( wslpi, 'W', -1. ) ; CALL lbc_lnk( wslpj, 'W', -1. ) |
---|
[2528] | 881 | ENDIF |
---|
[2715] | 882 | ENDIF |
---|
| 883 | ! |
---|
[3294] | 884 | IF( nn_timing == 1 ) CALL timing_stop('ldf_slp_init') |
---|
| 885 | ! |
---|
[3] | 886 | END SUBROUTINE ldf_slp_init |
---|
| 887 | |
---|
| 888 | #else |
---|
| 889 | !!------------------------------------------------------------------------ |
---|
| 890 | !! Dummy module : NO Rotation of lateral mixing tensor |
---|
| 891 | !!------------------------------------------------------------------------ |
---|
[32] | 892 | LOGICAL, PUBLIC, PARAMETER :: lk_ldfslp = .FALSE. !: slopes flag |
---|
[3] | 893 | CONTAINS |
---|
[3294] | 894 | SUBROUTINE ldf_slp( kt, prd, pn2 ) ! Dummy routine |
---|
| 895 | INTEGER, INTENT(in) :: kt |
---|
[1515] | 896 | REAL, DIMENSION(:,:,:), INTENT(in) :: prd, pn2 |
---|
[32] | 897 | WRITE(*,*) 'ldf_slp: You should not have seen this print! error?', kt, prd(1,1,1), pn2(1,1,1) |
---|
[3] | 898 | END SUBROUTINE ldf_slp |
---|
[2605] | 899 | SUBROUTINE ldf_slp_grif( kt ) ! Dummy routine |
---|
| 900 | INTEGER, INTENT(in) :: kt |
---|
| 901 | WRITE(*,*) 'ldf_slp_grif: You should not have seen this print! error?', kt |
---|
| 902 | END SUBROUTINE ldf_slp_grif |
---|
[3294] | 903 | SUBROUTINE ldf_slp_init ! Dummy routine |
---|
[2528] | 904 | END SUBROUTINE ldf_slp_init |
---|
[3] | 905 | #endif |
---|
| 906 | |
---|
| 907 | !!====================================================================== |
---|
| 908 | END MODULE ldfslp |
---|