- Timestamp:
- 2015-06-19T18:07:11+02:00 (9 years ago)
- Location:
- branches/UKMO/dev_r5107_mld_zint/NEMOGCM/NEMO/OPA_SRC/LDF
- Files:
-
- 5 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/UKMO/dev_r5107_mld_zint/NEMOGCM/NEMO/OPA_SRC/LDF/ldfdyn_c2d.h90
r5247 r5447 140 140 !!---------------------------------------------------------------------- 141 141 USE ldftra_oce, ONLY: aht0 142 USE iom 142 143 ! 143 144 LOGICAL, INTENT (in) :: ld_print ! If true, output arrays on numout … … 150 151 CHARACTER (len=15) :: clexp 151 152 INTEGER, POINTER, DIMENSION(:,:) :: icof 152 INTEGER, ALLOCATABLE, DIMENSION(:,:) :: idata153 REAL(wp), DIMENSION(:,:), ALLOCATABLE :: ztemp2d ! temporary array to read ahmcoef file 153 154 !!---------------------------------------------------------------------- 154 155 ! … … 232 233 ! Read 2d integer array to specify western boundary increase in the 233 234 ! ===================== equatorial strip (20N-20S) defined at t-points 234 235 ALLOCATE( idata(jpidta,jpjdta), STAT=ierror )236 IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'ldf_dyn_c2d_orca: unable to allocate idata array' )237 235 ! 238 CALL ctl_opn( inum, 'ahmcoef', 'OLD', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp ) 239 READ(inum,9101) clexp, iim, ijm 240 READ(inum,'(/)') 241 ifreq = 40 242 il1 = 1 243 DO jn = 1, jpidta/ifreq+1 244 READ(inum,'(/)') 245 il2 = MIN( jpidta, il1+ifreq-1 ) 246 READ(inum,9201) ( ii, ji = il1, il2, 5 ) 247 READ(inum,'(/)') 248 DO jj = jpjdta, 1, -1 249 READ(inum,9202) ij, ( idata(ji,jj), ji = il1, il2 ) 250 END DO 251 il1 = il1 + ifreq 252 END DO 253 254 DO jj = 1, nlcj 255 DO ji = 1, nlci 256 icof(ji,jj) = idata( mig(ji), mjg(jj) ) 257 END DO 258 END DO 259 DO jj = nlcj+1, jpj 260 DO ji = 1, nlci 261 icof(ji,jj) = icof(ji,nlcj) 262 END DO 263 END DO 264 DO jj = 1, jpj 265 DO ji = nlci+1, jpi 266 icof(ji,jj) = icof(nlci,jj) 267 END DO 268 END DO 269 270 9101 FORMAT(1x,a15,2i8) 271 9201 FORMAT(3x,13(i3,12x)) 272 9202 FORMAT(i3,41i3) 273 274 DEALLOCATE(idata) 236 ALLOCATE( ztemp2d(jpi,jpj) ) 237 ztemp2d(:,:) = 0. 238 CALL iom_open ( 'ahmcoef.nc', inum ) 239 CALL iom_get ( inum, jpdom_data, 'icof', ztemp2d) 240 icof(:,:) = NINT(ztemp2d(:,:)) 241 CALL iom_close( inum ) 242 DEALLOCATE(ztemp2d) 275 243 276 244 ! Set ahm1 and ahm2 ( T- and F- points) (used for laplacian operator) … … 369 337 !!---------------------------------------------------------------------- 370 338 USE ldftra_oce, ONLY: aht0 339 USE iom 371 340 ! 372 341 LOGICAL, INTENT (in) :: ld_print ! If true, output arrays on numout … … 380 349 CHARACTER (len=15) :: clexp 381 350 INTEGER, POINTER, DIMENSION(:,:) :: icof 382 INTEGER, ALLOCATABLE, DIMENSION(:,:) :: idata351 REAL(wp), DIMENSION(:,:), ALLOCATABLE :: ztemp2d ! temporary array to read ahmcoef file 383 352 !!---------------------------------------------------------------------- 384 353 ! 385 354 CALL wrk_alloc( jpi , jpj , icof ) 386 355 ! 387 388 356 IF(lwp) WRITE(numout,*) 389 357 IF(lwp) WRITE(numout,*) 'inildf: 2d eddy viscosity coefficient' … … 464 432 ! Read 2d integer array to specify western boundary increase in the 465 433 ! ===================== equatorial strip (20N-20S) defined at t-points 466 467 ALLOCATE( idata(jpidta,jpjdta), STAT=ierror ) 468 IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'ldf_dyn_c2d_orca_R1: unable to allocate idata array' ) 469 ! 470 CALL ctl_opn( inum, 'ahmcoef', 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL', & 471 & 1, numout, lwp ) 472 REWIND inum 473 READ(inum,9101) clexp, iim, ijm 474 READ(inum,'(/)') 475 ifreq = 40 476 il1 = 1 477 DO jn = 1, jpidta/ifreq+1 478 READ(inum,'(/)') 479 il2 = MIN( jpidta, il1+ifreq-1 ) 480 READ(inum,9201) ( ii, ji = il1, il2, 5 ) 481 READ(inum,'(/)') 482 DO jj = jpjdta, 1, -1 483 READ(inum,9202) ij, ( idata(ji,jj), ji = il1, il2 ) 484 END DO 485 il1 = il1 + ifreq 486 END DO 487 488 DO jj = 1, nlcj 489 DO ji = 1, nlci 490 icof(ji,jj) = idata( mig(ji), mjg(jj) ) 491 END DO 492 END DO 493 DO jj = nlcj+1, jpj 494 DO ji = 1, nlci 495 icof(ji,jj) = icof(ji,nlcj) 496 END DO 497 END DO 498 DO jj = 1, jpj 499 DO ji = nlci+1, jpi 500 icof(ji,jj) = icof(nlci,jj) 501 END DO 502 END DO 503 504 9101 FORMAT(1x,a15,2i8) 505 9201 FORMAT(3x,13(i3,12x)) 506 9202 FORMAT(i3,41i3) 507 508 DEALLOCATE(idata) 434 ALLOCATE( ztemp2d(jpi,jpj) ) 435 ztemp2d(:,:) = 0. 436 CALL iom_open ( 'ahmcoef.nc', inum ) 437 CALL iom_get ( inum, jpdom_data, 'icof', ztemp2d) 438 icof(:,:) = NINT(ztemp2d(:,:)) 439 CALL iom_close( inum ) 440 DEALLOCATE(ztemp2d) 509 441 510 442 ! Set ahm1 and ahm2 ( T- and F- points) (used for laplacian operator) -
branches/UKMO/dev_r5107_mld_zint/NEMOGCM/NEMO/OPA_SRC/LDF/ldfdyn_c3d.h90
r5247 r5447 27 27 !!---------------------------------------------------------------------- 28 28 USE ldftra_oce, ONLY : aht0 29 USE iom 29 30 !! 30 31 LOGICAL, INTENT (in) :: ld_print ! If true, output arrays on numout … … 193 194 !!---------------------------------------------------------------------- 194 195 USE ldftra_oce, ONLY: aht0 196 USE iom 195 197 !! 196 198 LOGICAL, INTENT(in) :: ld_print ! If true, output arrays on numout … … 204 206 CHARACTER (len=15) :: clexp 205 207 INTEGER , POINTER, DIMENSION(:,:) :: icof 206 INTEGER , POINTER, DIMENSION(:,:) :: idata207 208 REAL(wp), POINTER, DIMENSION(: ) :: zcoef 208 209 REAL(wp), POINTER, DIMENSION(:,:) :: zahm0 210 ! 211 REAL(wp), DIMENSION(:,:), ALLOCATABLE :: ztemp2d ! temporary array to read ahmcoef file 209 212 !!---------------------------------------------------------------------- 210 213 ! 211 214 CALL wrk_alloc( jpi , jpj , icof ) 212 CALL wrk_alloc( jpidta, jpjdta, idata )213 215 CALL wrk_alloc( jpk , zcoef ) 214 216 CALL wrk_alloc( jpi , jpj , zahm0 ) … … 221 223 ! Read 2d integer array to specify western boundary increase in the 222 224 ! ===================== equatorial strip (20N-20S) defined at t-points 223 224 CALL ctl_opn( inum, 'ahmcoef', 'OLD', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp ) 225 READ(inum,9101) clexp, iim, ijm 226 READ(inum,'(/)') 227 ifreq = 40 228 il1 = 1 229 DO jn = 1, jpidta/ifreq+1 230 READ(inum,'(/)') 231 il2 = MIN( jpidta, il1+ifreq-1 ) 232 READ(inum,9201) ( ii, ji = il1, il2, 5 ) 233 READ(inum,'(/)') 234 DO jj = jpjdta, 1, -1 235 READ(inum,9202) ij, ( idata(ji,jj), ji = il1, il2 ) 236 END DO 237 il1 = il1 + ifreq 238 END DO 239 240 DO jj = 1, nlcj 241 DO ji = 1, nlci 242 icof(ji,jj) = idata( mig(ji), mjg(jj) ) 243 END DO 244 END DO 245 DO jj = nlcj+1, jpj 246 DO ji = 1, nlci 247 icof(ji,jj) = icof(ji,nlcj) 248 END DO 249 END DO 250 DO jj = 1, jpj 251 DO ji = nlci+1, jpi 252 icof(ji,jj) = icof(nlci,jj) 253 END DO 254 END DO 255 256 9101 FORMAT(1x,a15,2i8) 257 9201 FORMAT(3x,13(i3,12x)) 258 9202 FORMAT(i3,41i3) 259 225 ALLOCATE( ztemp2d(jpi,jpj) ) 226 ztemp2d(:,:) = 0. 227 CALL iom_open ( 'ahmcoef.nc', inum ) 228 CALL iom_get ( inum, jpdom_data, 'icof', ztemp2d) 229 icof(:,:) = NINT(ztemp2d(:,:)) 230 CALL iom_close( inum ) 231 DEALLOCATE(ztemp2d) 232 260 233 ! Set ahm1 and ahm2 261 234 ! ================= … … 455 428 ! 456 429 CALL wrk_dealloc( jpi , jpj , icof ) 457 CALL wrk_dealloc( jpidta, jpjdta, idata )458 430 CALL wrk_dealloc( jpk , zcoef ) 459 431 CALL wrk_dealloc( jpi , jpj , zahm0 ) -
branches/UKMO/dev_r5107_mld_zint/NEMOGCM/NEMO/OPA_SRC/LDF/ldfdyn_smag.F90
- Property svn:keywords set to Id
-
branches/UKMO/dev_r5107_mld_zint/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90
r5247 r5447 142 142 DO jj = 1, jpjm1 143 143 DO ji = 1, jpim1 144 ! IF should be useless check zpshde (PM) 145 IF ( mbku(ji,jj) > 1 ) zgru(ji,jj,mbku(ji,jj)) = gru(ji,jj) 146 IF ( mbkv(ji,jj) > 1 ) zgrv(ji,jj,mbkv(ji,jj)) = grv(ji,jj) 144 zgru(ji,jj,mbku(ji,jj)) = gru(ji,jj) 145 zgrv(ji,jj,mbkv(ji,jj)) = grv(ji,jj) 146 END DO 147 END DO 148 ENDIF 149 IF( ln_zps .AND. ln_isfcav ) THEN ! partial steps correction at the bottom ocean level 150 DO jj = 1, jpjm1 151 DO ji = 1, jpim1 147 152 IF ( miku(ji,jj) > 1 ) zgru(ji,jj,miku(ji,jj)) = grui(ji,jj) 148 153 IF ( mikv(ji,jj) > 1 ) zgrv(ji,jj,mikv(ji,jj)) = grvi(ji,jj) … … 151 156 ENDIF 152 157 ! 153 zdzr(:,:,1) = 0._wp !== Local vertical density gradient at T-point == ! (evaluated from N^2) 154 DO jk = 1, jpkm1 158 !== Local vertical density gradient at T-point == ! (evaluated from N^2) 159 ! interior value 160 DO jk = 2, jpkm1 155 161 ! ! zdzr = d/dz(prd)= - ( prd ) / grav * mk(pn2) -- at t point 156 162 ! ! trick: tmask(ik ) = 0 => all pn2 = 0 => zdzr = 0 … … 162 168 END DO 163 169 ! surface initialisation 164 DO jj = 1, jpjm1 165 DO ji = 1, jpim1 166 zdzr(ji,jj,1:mikt(ji,jj)) = 0._wp 167 END DO 168 END DO 170 zdzr(:,:,1) = 0._wp 171 IF ( ln_isfcav ) THEN 172 ! if isf need to overwrite the interior value at at the first ocean point 173 DO jj = 1, jpjm1 174 DO ji = 1, jpim1 175 zdzr(ji,jj,1:mikt(ji,jj)) = 0._wp 176 END DO 177 END DO 178 END IF 169 179 ! 170 180 ! !== Slopes just below the mixed layer ==! … … 175 185 ! =========================== | vslp = d/dj( prd ) / d/dz( prd ) 176 186 ! 177 DO jj = 2, jpjm1 178 DO ji = fs_2, fs_jpim1 ! vector opt. 179 IF (miku(ji,jj) .GT. miku(ji+1,jj)) zhmlpu(ji,jj) = hmlpt(ji ,jj) 180 IF (miku(ji,jj) .LT. miku(ji+1,jj)) zhmlpu(ji,jj) = hmlpt(ji+1,jj) 181 IF (miku(ji,jj) .EQ. miku(ji+1,jj)) zhmlpu(ji,jj) = MAX(hmlpt(ji ,jj), hmlpt(ji+1,jj)) 182 IF (mikv(ji,jj) .GT. miku(ji,jj+1)) zhmlpv(ji,jj) = hmlpt(ji ,jj) 183 IF (mikv(ji,jj) .LT. miku(ji,jj+1)) zhmlpv(ji,jj) = hmlpt(ji,jj+1) 184 IF (mikv(ji,jj) .EQ. miku(ji,jj+1)) zhmlpv(ji,jj) = MAX(hmlpt(ji,jj), hmlpt(ji,jj+1)) 187 IF ( ln_isfcav ) THEN 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) = MAX(hmlpt(ji ,jj ), 5._wp) 191 IF (miku(ji,jj) .LT. miku(ji+1,jj)) zhmlpu(ji,jj) = MAX(hmlpt(ji+1,jj ), 5._wp) 192 IF (miku(ji,jj) .EQ. miku(ji+1,jj)) zhmlpu(ji,jj) = MAX(hmlpt(ji ,jj ), hmlpt(ji+1,jj ), 5._wp) 193 IF (mikv(ji,jj) .GT. miku(ji,jj+1)) zhmlpv(ji,jj) = MAX(hmlpt(ji ,jj ), 5._wp) 194 IF (mikv(ji,jj) .LT. miku(ji,jj+1)) zhmlpv(ji,jj) = MAX(hmlpt(ji ,jj+1), 5._wp) 195 IF (mikv(ji,jj) .EQ. miku(ji,jj+1)) zhmlpv(ji,jj) = MAX(hmlpt(ji ,jj ), hmlpt(ji ,jj+1), 5._wp) 196 ENDDO 185 197 ENDDO 186 ENDDO 198 ELSE 199 DO jj = 2, jpjm1 200 DO ji = fs_2, fs_jpim1 ! vector opt. 201 zhmlpu(ji,jj) = MAX(hmlpt(ji,jj), hmlpt(ji+1,jj ), 5._wp) 202 zhmlpv(ji,jj) = MAX(hmlpt(ji,jj), hmlpt(ji ,jj+1), 5._wp) 203 ENDDO 204 ENDDO 205 END IF 187 206 DO jk = 2, jpkm1 !* Slopes at u and v points 188 207 DO jj = 2, jpjm1 … … 198 217 zbv = MIN( zbv, -100._wp* ABS( zav ) , -7.e+3_wp/fse3v(ji,jj,jk)* ABS( zav ) ) 199 218 ! ! uslp and vslp output in zwz and zww, resp. 200 zfi = MAX( omlmask(ji,jj,jk), omlmask(ji+1,jj ,jk) )201 zfj = MAX( omlmask(ji,jj,jk), omlmask(ji ,jj+1,jk) )219 zfi = MAX( omlmask(ji,jj,jk), omlmask(ji+1,jj ,jk) ) 220 zfj = MAX( omlmask(ji,jj,jk), omlmask(ji ,jj+1,jk) ) 202 221 ! thickness of water column between surface and level k at u/v point 203 zdepu = 0.5_wp * (( fsdept(ji,jj,jk) + fsdept(ji+1,jj ,jk) ) & 204 - 2 * MAX( risfdep(ji,jj), risfdep(ji+1,jj ) ) & 205 - fse3u(ji,jj,miku(ji,jj)) ) 206 zdepv = 0.5_wp * (( fsdept(ji,jj,jk) + fsdept(ji ,jj+1,jk) ) & 207 - 2 * MAX( risfdep(ji,jj), risfdep(ji,jj+1) ) & 208 - fse3v(ji,jj,mikv(ji,jj)) ) 209 zwz(ji,jj,jk) = ( 1. - zfi) * zau / ( zbu - zeps ) & 210 & + zfi * uslpml(ji,jj) & 211 & * zdepu / MAX( zhmlpu(ji,jj), 5._wp ) 212 zwz(ji,jj,jk) = zwz(ji,jj,jk) * umask(ji,jj,jk) * umask(ji,jj,jk-1) 213 zww(ji,jj,jk) = ( 1. - zfj) * zav / ( zbv - zeps ) & 214 & + zfj * vslpml(ji,jj) & 215 & * zdepv / MAX( zhmlpv(ji,jj), 5._wp ) 216 zww(ji,jj,jk) = zww(ji,jj,jk) * vmask(ji,jj,jk) * vmask(ji,jj,jk-1) 222 zdepu = 0.5_wp * ( ( fsdept(ji,jj,jk) + fsdept(ji+1,jj ,jk) ) & 223 - 2 * MAX( risfdep(ji,jj), risfdep(ji+1,jj ) ) - fse3u(ji,jj,miku(ji,jj)) ) 224 zdepv = 0.5_wp * ( ( fsdept(ji,jj,jk) + fsdept(ji ,jj+1,jk) ) & 225 - 2 * MAX( risfdep(ji,jj), risfdep(ji ,jj+1) ) - fse3v(ji,jj,mikv(ji,jj)) ) 226 ! 227 zwz(ji,jj,jk) = ( 1. - zfi) * zau / ( zbu - zeps ) & 228 & + zfi * uslpml(ji,jj) * zdepu / zhmlpu(ji,jj) 229 zwz(ji,jj,jk) = zwz(ji,jj,jk) * wumask(ji,jj,jk) 230 zww(ji,jj,jk) = ( 1. - zfj) * zav / ( zbv - zeps ) & 231 & + zfj * vslpml(ji,jj) * zdepv / zhmlpv(ji,jj) 232 zww(ji,jj,jk) = zww(ji,jj,jk) * wvmask(ji,jj,jk) 217 233 218 234 … … 266 282 uslp(ji,jj,jk) = uslp(ji,jj,jk) * ( umask(ji,jj+1,jk) + umask(ji,jj-1,jk ) ) * 0.5_wp & 267 283 & * ( umask(ji,jj ,jk) + umask(ji,jj ,jk+1) ) * 0.5_wp & 268 & * umask(ji,jj,jk-1) !* umask(ji,jj,jk) * umask(ji,jj,jk+1)284 & * umask(ji,jj,jk-1) 269 285 vslp(ji,jj,jk) = vslp(ji,jj,jk) * ( vmask(ji+1,jj,jk) + vmask(ji-1,jj,jk ) ) * 0.5_wp & 270 286 & * ( vmask(ji ,jj,jk) + vmask(ji ,jj,jk+1) ) * 0.5_wp & 271 & * vmask(ji,jj,jk-1) !* vmask(ji,jj,jk) * vmask(ji,jj,jk+1)287 & * vmask(ji,jj,jk-1) 272 288 END DO 273 289 END DO … … 282 298 DO ji = fs_2, fs_jpim1 ! vector opt. 283 299 ! !* Local vertical density gradient evaluated from N^2 284 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)300 zbw = zm1_2g * pn2 (ji,jj,jk) * ( prd (ji,jj,jk) + prd (ji,jj,jk-1) + 2. ) * wmask(ji,jj,jk) 285 301 ! !* Slopes at w point 286 302 ! ! i- & j-gradient of density at w-points … … 298 314 zbj = MIN( zbw , -100._wp* ABS( zaj ) , -7.e+3_wp/fse3w(ji,jj,jk)* ABS( zaj ) ) 299 315 ! ! wslpi and wslpj with ML flattening (output in zwz and zww, resp.) 300 zfk = MAX( omlmask(ji,jj,jk), omlmask(ji,jj,jk-1) ) 316 zfk = MAX( omlmask(ji,jj,jk), omlmask(ji,jj,jk-1) ) ! zfk=1 in the ML otherwise zfk=0 301 317 zck = ( fsdepw(ji,jj,jk) - fsdepw(ji,jj,mikt(ji,jj) ) ) / MAX( hmlp(ji,jj), 10._wp ) 302 318 zwz(ji,jj,jk) = ( zai / ( zbi - zeps ) * ( 1._wp - zfk ) & 303 & + zck * wslpiml(ji,jj) * zfk ) * tmask(ji,jj,jk) * tmask(ji,jj,jk-1)319 & + zck * wslpiml(ji,jj) * zfk ) * wmask(ji,jj,jk) 304 320 zww(ji,jj,jk) = ( zaj / ( zbj - zeps ) * ( 1._wp - zfk ) & 305 & + zck * wslpjml(ji,jj) * zfk ) * tmask(ji,jj,jk) * tmask(ji,jj,jk-1)321 & + zck * wslpjml(ji,jj) * zfk ) * wmask(ji,jj,jk) 306 322 307 323 !!gm modif to suppress omlmask.... (as in Griffies operator) … … 356 372 zck = ( umask(ji,jj,jk) + umask(ji-1,jj,jk) ) & 357 373 & * ( vmask(ji,jj,jk) + vmask(ji,jj-1,jk) ) * 0.25 358 wslpi(ji,jj,jk) = wslpi(ji,jj,jk) * zck * tmask(ji,jj,jk-1) * tmask(ji,jj,jk)359 wslpj(ji,jj,jk) = wslpj(ji,jj,jk) * zck * tmask(ji,jj,jk-1) * tmask(ji,jj,jk)374 wslpi(ji,jj,jk) = wslpi(ji,jj,jk) * zck * wmask(ji,jj,jk) 375 wslpj(ji,jj,jk) = wslpj(ji,jj,jk) * zck * wmask(ji,jj,jk) 360 376 END DO 361 377 END DO … … 423 439 vslp(ji,jj,jk) = -1./e2v(ji,jj) * ( fsdept_b(ji,jj+1,jk) - fsdept_b(ji ,jj ,jk) ) * vmask(ji,jj,jk) 424 440 wslpi(ji,jj,jk) = -1./e1t(ji,jj) * ( fsdepw_b(ji+1,jj,jk) - fsdepw_b(ji-1,jj,jk) ) & 425 & * tmask(ji,jj,jk) * tmask(ji,jj,jk-1) * 0.5441 & * wmask(ji,jj,jk) * 0.5 426 442 wslpj(ji,jj,jk) = -1./e2t(ji,jj) * ( fsdepw_b(ji,jj+1,jk) - fsdepw_b(ji,jj-1,jk) ) & 427 & * tmask(ji,jj,jk) * tmask(ji,jj,jk-1) * 0.5443 & * wmask(ji,jj,jk) * 0.5 428 444 END DO 429 445 END DO … … 736 752 DO ji = 1, jpi 737 753 ik = nmln(ji,jj) - 1 738 IF( jk <= ik .AND. jk >= mikt(ji,jj) ) THEN ; omlmask(ji,jj,jk) = 1._wp 739 ELSE ; omlmask(ji,jj,jk) = 0._wp 754 IF( jk <= ik .AND. jk >= mikt(ji,jj) ) THEN 755 omlmask(ji,jj,jk) = 1._wp 756 ELSE 757 omlmask(ji,jj,jk) = 0._wp 740 758 ENDIF 741 759 END DO … … 794 812 zbj = MIN( zbw , -100._wp* ABS( zaj ) , -7.e+3_wp/fse3w(ji,jj,ik)* ABS( zaj ) ) 795 813 ! !- i- & j-slope at w-points (wslpiml, wslpjml) 796 wslpiml(ji,jj) = zai / ( zbi - zeps ) * tmask (ji,jj,ik)797 wslpjml(ji,jj) = zaj / ( zbj - zeps ) * tmask (ji,jj,ik)814 wslpiml(ji,jj) = zai / ( zbi - zeps ) * wmask (ji,jj,ik) 815 wslpjml(ji,jj) = zaj / ( zbj - zeps ) * wmask (ji,jj,ik) 798 816 END DO 799 817 END DO -
branches/UKMO/dev_r5107_mld_zint/NEMOGCM/NEMO/OPA_SRC/LDF/ldftra_smag.F90
- Property svn:keywords set to Id
Note: See TracChangeset
for help on using the changeset viewer.