Changeset 11533
- Timestamp:
- 2019-09-11T15:37:04+02:00 (5 years ago)
- Location:
- NEMO/releases/release-4.0
- Files:
-
- 10 edited
Legend:
- Unmodified
- Added
- Removed
-
NEMO/releases/release-4.0/cfgs/SHARED/namelist_ref
r11162 r11533 189 189 &namsbc ! Surface Boundary Condition manager (default: NO selection) 190 190 !----------------------------------------------------------------------- 191 nn_fsbc = 4! frequency of SBC module call191 nn_fsbc = 2 ! frequency of SBC module call 192 192 ! ! (control sea-ice & iceberg model call) 193 193 ! Type of air-sea fluxes -
NEMO/releases/release-4.0/src/OCE/DIA/diawri.F90
r10425 r11533 210 210 ENDIF 211 211 212 IF( ln_zad_Aimp ) wn = wn + wi ! Recombine explicit and implicit parts of vertical velocity for diagnostic output 213 ! 212 214 CALL iom_put( "woce", wn ) ! vertical velocity 213 215 IF( iom_use('w_masstr') .OR. iom_use('w_masstr2') ) THEN ! vertical mass transport & its square value … … 220 222 IF( iom_use('w_masstr2') ) CALL iom_put( "w_masstr2", z3d(:,:,:) * z3d(:,:,:) ) 221 223 ENDIF 224 ! 225 IF( ln_zad_Aimp ) wn = wn - wi ! Remove implicit part of vertical velocity that was added for diagnostic output 222 226 223 227 CALL iom_put( "avt" , avt ) ! T vert. eddy diff. coef. … … 842 846 CALL histwrite( nid_V, "sometauy", it, vtau , ndim_hV, ndex_hV ) ! j-wind stress 843 847 844 CALL histwrite( nid_W, "vovecrtz", it, wn , ndim_T, ndex_T ) ! vert. current 848 IF( ln_zad_Aimp ) THEN 849 CALL histwrite( nid_W, "vovecrtz", it, wn + wi , ndim_T, ndex_T ) ! vert. current 850 ELSE 851 CALL histwrite( nid_W, "vovecrtz", it, wn , ndim_T, ndex_T ) ! vert. current 852 ENDIF 845 853 CALL histwrite( nid_W, "votkeavt", it, avt , ndim_T, ndex_T ) ! T vert. eddy diff. coef. 846 854 CALL histwrite( nid_W, "votkeavm", it, avm , ndim_T, ndex_T ) ! T vert. eddy visc. coef. … … 903 911 CALL iom_rstput( 0, 0, inum, 'vozocrtx', un ) ! now i-velocity 904 912 CALL iom_rstput( 0, 0, inum, 'vomecrty', vn ) ! now j-velocity 905 CALL iom_rstput( 0, 0, inum, 'vovecrtz', wn ) ! now k-velocity 913 IF( ln_zad_Aimp ) THEN 914 CALL iom_rstput( 0, 0, inum, 'vovecrtz', wn + wi ) ! now k-velocity 915 ELSE 916 CALL iom_rstput( 0, 0, inum, 'vovecrtz', wn ) ! now k-velocity 917 ENDIF 906 918 IF( ALLOCATED(ahtu) ) THEN 907 919 CALL iom_rstput( 0, 0, inum, 'ahtu', ahtu ) ! aht at u-point -
NEMO/releases/release-4.0/src/OCE/DOM/domwri.F90
r10425 r11533 162 162 CALL iom_rstput( 0, 0, inum, 'isfdraft', zprt, ktype = jp_r8 ) ! ! nb of ocean T-points 163 163 ! ! vertical mesh 164 CALL iom_rstput( 0, 0, inum, 'e3t_0', e3t_0, ktype = jp_r8 ) ! ! scale factors 165 CALL iom_rstput( 0, 0, inum, 'e3u_0', e3u_0, ktype = jp_r8 ) 166 CALL iom_rstput( 0, 0, inum, 'e3v_0', e3v_0, ktype = jp_r8 ) 167 CALL iom_rstput( 0, 0, inum, 'e3w_0', e3w_0, ktype = jp_r8 ) 164 CALL iom_rstput( 0, 0, inum, 'e3t_1d', e3t_1d, ktype = jp_r8 ) ! ! scale factors 165 CALL iom_rstput( 0, 0, inum, 'e3w_1d', e3w_1d, ktype = jp_r8 ) 166 167 CALL iom_rstput( 0, 0, inum, 'e3t_0' , e3t_0 , ktype = jp_r8 ) 168 CALL iom_rstput( 0, 0, inum, 'e3u_0' , e3u_0 , ktype = jp_r8 ) 169 CALL iom_rstput( 0, 0, inum, 'e3v_0' , e3v_0 , ktype = jp_r8 ) 170 CALL iom_rstput( 0, 0, inum, 'e3f_0' , e3f_0 , ktype = jp_r8 ) 171 CALL iom_rstput( 0, 0, inum, 'e3w_0' , e3w_0 , ktype = jp_r8 ) 172 CALL iom_rstput( 0, 0, inum, 'e3uw_0', e3uw_0, ktype = jp_r8 ) 173 CALL iom_rstput( 0, 0, inum, 'e3vw_0', e3vw_0, ktype = jp_r8 ) 168 174 ! 169 175 CALL iom_rstput( 0, 0, inum, 'gdept_1d' , gdept_1d , ktype = jp_r8 ) ! stretched system -
NEMO/releases/release-4.0/src/OCE/DYN/dynhpg.F90
r10491 r11533 37 37 USE trd_oce ! trends: ocean variables 38 38 USE trddyn ! trend manager: dynamics 39 !jcUSE zpshde ! partial step: hor. derivative (zps_hde routine)39 USE zpshde ! partial step: hor. derivative (zps_hde routine) 40 40 ! 41 41 USE in_out_manager ! I/O manager … … 338 338 REAL(wp) :: zcoef0, zcoef1, zcoef2, zcoef3 ! temporary scalars 339 339 REAL(wp), DIMENSION(jpi,jpj,jpk) :: zhpi, zhpj 340 REAL(wp), DIMENSION(jpi,jpj) :: zgtsu, zgtsv, zgru, zgrv 340 341 !!---------------------------------------------------------------------- 341 342 ! … … 346 347 ENDIF 347 348 348 ! Partial steps: bottom beforehorizontal gradient of t, s, rd at the last ocean level349 !jc CALL zps_hde ( kt, jpts, tsn, gtsu, gtsv, rhd, gru ,grv )349 ! Partial steps: Compute NOW horizontal gradient of t, s, rd at the last ocean level 350 CALL zps_hde( kt, jpts, tsn, zgtsu, zgtsv, rhd, zgru , zgrv ) 350 351 351 352 ! Local constant initialization … … 385 386 END DO 386 387 387 ! partial steps correction at the last level (use gru &grv computed in zpshde.F90)388 ! partial steps correction at the last level (use zgru & zgrv computed in zpshde.F90) 388 389 DO jj = 2, jpjm1 389 390 DO ji = 2, jpim1 … … 395 396 ua (ji,jj,iku) = ua(ji,jj,iku) - zhpi(ji,jj,iku) ! subtract old value 396 397 zhpi(ji,jj,iku) = zhpi(ji,jj,iku-1) & ! compute the new one 397 & + zcoef2 * ( rhd(ji+1,jj,iku-1) - rhd(ji,jj,iku-1) + gru(ji,jj) ) * r1_e1u(ji,jj)398 & + zcoef2 * ( rhd(ji+1,jj,iku-1) - rhd(ji,jj,iku-1) + zgru(ji,jj) ) * r1_e1u(ji,jj) 398 399 ua (ji,jj,iku) = ua(ji,jj,iku) + zhpi(ji,jj,iku) ! add the new one to the general momentum trend 399 400 ENDIF … … 401 402 va (ji,jj,ikv) = va(ji,jj,ikv) - zhpj(ji,jj,ikv) ! subtract old value 402 403 zhpj(ji,jj,ikv) = zhpj(ji,jj,ikv-1) & ! compute the new one 403 & + zcoef3 * ( rhd(ji,jj+1,ikv-1) - rhd(ji,jj,ikv-1) + grv(ji,jj) ) * r1_e2v(ji,jj)404 & + zcoef3 * ( rhd(ji,jj+1,ikv-1) - rhd(ji,jj,ikv-1) + zgrv(ji,jj) ) * r1_e2v(ji,jj) 404 405 va (ji,jj,ikv) = va(ji,jj,ikv) + zhpj(ji,jj,ikv) ! add the new one to the general momentum trend 405 406 ENDIF -
NEMO/releases/release-4.0/src/OCE/DYN/sshwzv.F90
r11294 r11533 9 9 !! - ! 2010-09 (D.Storkey and E.O'Dea) bug fixes for BDY module 10 10 !! 3.3 ! 2011-10 (M. Leclair) split former ssh_wzv routine and remove all vvl related work 11 !! 4.0 ! 2018-12 (A. Coward) add mixed implicit/explicit advection 11 12 !!---------------------------------------------------------------------- 12 13 -
NEMO/releases/release-4.0/src/OCE/LBC/mppini.F90
r10616 r11533 669 669 ! 670 670 CALL mpp_init_ioipsl ! Prepare NetCDF output file (if necessary) 671 ! 672 IF ( ln_nnogather) THEN671 ! 672 IF (( jperio >= 3 .AND. jperio <= 6 .AND. jpni > 1 ).AND.( ln_nnogather )) THEN 673 673 CALL mpp_init_nfdcom ! northfold neighbour lists 674 674 IF (llwrtlay) THEN -
NEMO/releases/release-4.0/src/OCE/TRA/traadv_fct.F90
r10425 r11533 21 21 USE diaar5 ! AR5 diagnostics 22 22 USE phycst , ONLY : rau0_rcp 23 USE zdf_oce , ONLY : ln_zad_Aimp 23 24 ! 24 25 USE in_out_manager ! I/O manager … … 86 87 REAL(wp), DIMENSION(jpi,jpj,jpk) :: zwi, zwx, zwy, zwz, ztu, ztv, zltu, zltv, ztw 87 88 REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: ztrdx, ztrdy, ztrdz, zptry 89 REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: zwinf, zwdia, zwsup 90 LOGICAL :: ll_zAimp ! flag to apply adaptive implicit vertical advection 88 91 !!---------------------------------------------------------------------- 89 92 ! … … 97 100 l_hst = .FALSE. 98 101 l_ptr = .FALSE. 102 ll_zAimp = .FALSE. 99 103 IF( ( cdtype =='TRA' .AND. l_trdtra ) .OR. ( cdtype =='TRC' .AND. l_trdtrc ) ) l_trd = .TRUE. 100 104 IF( cdtype =='TRA' .AND. ln_diaptr ) l_ptr = .TRUE. … … 116 120 ! 117 121 zwi(:,:,:) = 0._wp 122 ! 123 ! If adaptive vertical advection, check if it is needed on this PE at this time 124 IF( ln_zad_Aimp ) THEN 125 IF( MAXVAL( ABS( wi(:,:,:) ) ) > 0._wp ) ll_zAimp = .TRUE. 126 END IF 127 ! If active adaptive vertical advection, build tridiagonal matrix 128 IF( ll_zAimp ) THEN 129 ALLOCATE(zwdia(jpi,jpj,jpk), zwinf(jpi,jpj,jpk),zwsup(jpi,jpj,jpk)) 130 DO jk = 1, jpkm1 131 DO jj = 2, jpjm1 132 DO ji = fs_2, fs_jpim1 ! vector opt. (ensure same order of calculation as below if wi=0.) 133 zwdia(ji,jj,jk) = 1._wp + p2dt * ( MAX( wi(ji,jj,jk ) , 0._wp ) - MIN( wi(ji,jj,jk+1) , 0._wp ) ) / e3t_a(ji,jj,jk) 134 zwinf(ji,jj,jk) = p2dt * MIN( wi(ji,jj,jk ) , 0._wp ) / e3t_a(ji,jj,jk) 135 zwsup(ji,jj,jk) = -p2dt * MAX( wi(ji,jj,jk+1) , 0._wp ) / e3t_a(ji,jj,jk) 136 END DO 137 END DO 138 END DO 139 END IF 118 140 ! 119 141 DO jn = 1, kjpt !== loop over the tracers ==! … … 169 191 END DO 170 192 END DO 193 194 IF ( ll_zAimp ) THEN 195 CALL tridia_solver( zwdia, zwsup, zwinf, zwi, zwi , 0 ) 196 ! 197 ztw(:,:,1) = 0._wp ; ztw(:,:,jpk) = 0._wp ; 198 DO jk = 2, jpkm1 ! Interior value ( multiplied by wmask) 199 DO jj = 2, jpjm1 200 DO ji = fs_2, fs_jpim1 ! vector opt. 201 zfp_wk = wi(ji,jj,jk) + ABS( wi(ji,jj,jk) ) 202 zfm_wk = wi(ji,jj,jk) - ABS( wi(ji,jj,jk) ) 203 ztw(ji,jj,jk) = 0.5 * e1e2t(ji,jj) * ( zfp_wk * zwi(ji,jj,jk) + zfm_wk * zwi(ji,jj,jk-1) ) * wmask(ji,jj,jk) 204 zwz(ji,jj,jk) = zwz(ji,jj,jk) + ztw(ji,jj,jk) ! update vertical fluxes 205 END DO 206 END DO 207 END DO 208 DO jk = 1, jpkm1 209 DO jj = 2, jpjm1 210 DO ji = fs_2, fs_jpim1 ! vector opt. 211 pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) - ( ztw(ji,jj,jk) - ztw(ji ,jj ,jk+1) ) & 212 & * r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk) 213 END DO 214 END DO 215 END DO 216 ! 217 END IF 171 218 ! 172 219 IF( l_trd .OR. l_hst ) THEN ! trend diagnostics (contribution of upstream fluxes) … … 277 324 zwz(:,:,1) = 0._wp ! only ocean surface as interior zwz values have been w-masked 278 325 ENDIF 326 ! 327 IF ( ll_zAimp ) THEN 328 DO jk = 1, jpkm1 !* trend and after field with monotonic scheme 329 DO jj = 2, jpjm1 330 DO ji = fs_2, fs_jpim1 ! vector opt. 331 ! ! total intermediate advective trends 332 ztra = - ( zwx(ji,jj,jk) - zwx(ji-1,jj ,jk ) & 333 & + zwy(ji,jj,jk) - zwy(ji ,jj-1,jk ) & 334 & + zwz(ji,jj,jk) - zwz(ji ,jj ,jk+1) ) * r1_e1e2t(ji,jj) 335 ztw(ji,jj,jk) = zwi(ji,jj,jk) + p2dt * ztra / e3t_a(ji,jj,jk) * tmask(ji,jj,jk) 336 END DO 337 END DO 338 END DO 339 ! 340 CALL tridia_solver( zwdia, zwsup, zwinf, ztw, ztw , 0 ) 341 ! 342 DO jk = 2, jpkm1 ! Interior value ( multiplied by wmask) 343 DO jj = 2, jpjm1 344 DO ji = fs_2, fs_jpim1 ! vector opt. 345 zfp_wk = wi(ji,jj,jk) + ABS( wi(ji,jj,jk) ) 346 zfm_wk = wi(ji,jj,jk) - ABS( wi(ji,jj,jk) ) 347 zwz(ji,jj,jk) = zwz(ji,jj,jk) + 0.5 * e1e2t(ji,jj) * ( zfp_wk * ztw(ji,jj,jk) + zfm_wk * ztw(ji,jj,jk-1) ) * wmask(ji,jj,jk) 348 END DO 349 END DO 350 END DO 351 END IF 279 352 ! 280 353 CALL lbc_lnk_multi( 'traadv_fct', zwi, 'T', 1., zwx, 'U', -1. , zwy, 'V', -1., zwz, 'W', 1. ) … … 289 362 DO jj = 2, jpjm1 290 363 DO ji = fs_2, fs_jpim1 ! vector opt. 291 pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) - ( zwx(ji,jj,jk) - zwx(ji-1,jj ,jk ) & 292 & + zwy(ji,jj,jk) - zwy(ji ,jj-1,jk ) & 293 & + zwz(ji,jj,jk) - zwz(ji ,jj ,jk+1) ) & 294 & * r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk) 295 END DO 296 END DO 297 END DO 364 ztra = - ( zwx(ji,jj,jk) - zwx(ji-1,jj ,jk ) & 365 & + zwy(ji,jj,jk) - zwy(ji ,jj-1,jk ) & 366 & + zwz(ji,jj,jk) - zwz(ji ,jj ,jk+1) ) * r1_e1e2t(ji,jj) 367 pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) + ztra / e3t_n(ji,jj,jk) 368 zwi(ji,jj,jk) = zwi(ji,jj,jk) + p2dt * ztra / e3t_a(ji,jj,jk) * tmask(ji,jj,jk) 369 END DO 370 END DO 371 END DO 372 ! 373 IF ( ll_zAimp ) THEN 374 ! 375 ztw(:,:,1) = 0._wp ; ztw(:,:,jpk) = 0._wp 376 DO jk = 2, jpkm1 ! Interior value ( multiplied by wmask) 377 DO jj = 2, jpjm1 378 DO ji = fs_2, fs_jpim1 ! vector opt. 379 zfp_wk = wi(ji,jj,jk) + ABS( wi(ji,jj,jk) ) 380 zfm_wk = wi(ji,jj,jk) - ABS( wi(ji,jj,jk) ) 381 ztw(ji,jj,jk) = - 0.5 * e1e2t(ji,jj) * ( zfp_wk * zwi(ji,jj,jk) + zfm_wk * zwi(ji,jj,jk-1) ) * wmask(ji,jj,jk) 382 zwz(ji,jj,jk) = zwz(ji,jj,jk) + ztw(ji,jj,jk) ! Update vertical fluxes for trend diagnostic 383 END DO 384 END DO 385 END DO 386 DO jk = 1, jpkm1 387 DO jj = 2, jpjm1 388 DO ji = fs_2, fs_jpim1 ! vector opt. 389 pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) - ( ztw(ji,jj,jk) - ztw(ji ,jj ,jk+1) ) & 390 & * r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk) 391 END DO 392 END DO 393 END DO 394 END IF 298 395 ! 299 396 IF( l_trd .OR. l_hst ) THEN ! trend diagnostics // heat/salt transport … … 318 415 END DO ! end of tracer loop 319 416 ! 417 IF ( ll_zAimp ) THEN 418 DEALLOCATE( zwdia, zwinf, zwsup ) 419 ENDIF 320 420 IF( l_trd .OR. l_hst ) THEN 321 421 DEALLOCATE( ztrdx, ztrdy, ztrdz ) -
NEMO/releases/release-4.0/src/OCE/TRA/traqsr.F90
r10425 r11533 168 168 DO jj = 2, jpjm1 ! Separation in R-G-B depending of the surface Chl 169 169 DO ji = fs_2, fs_jpim1 170 zchl = sf_chl(1)%fnow(ji,jj,1)170 zchl = MIN( 10. , MAX( 0.03, sf_chl(1)%fnow(ji,jj,1) ) ) 171 171 zCtot = 40.6 * zchl**0.459 172 172 zze = 568.2 * zCtot**(-0.746) -
NEMO/releases/release-4.0/src/OCE/step.F90
r11288 r11533 165 165 CALL eos ( tsn, rhd, rhop, gdept_n(:,:,:) ) ! now in situ density for hpg computation 166 166 167 !!jc: fs simplification168 !!jc: lines below are useless if ln_linssh=F. Keep them here (which maintains a bug if ln_linssh=T and ln_zps=T, cf ticket #1636)169 !! but ensures reproductible results170 !! with previous versions using split-explicit free surface171 IF( ln_zps .AND. .NOT. ln_isfcav ) &172 & CALL zps_hde ( kstp, jpts, tsn, gtsu, gtsv, & ! Partial steps: before horizontal gradient173 & rhd, gru , grv ) ! of t, s, rd at the last ocean level174 IF( ln_zps .AND. ln_isfcav ) &175 & CALL zps_hde_isf( kstp, jpts, tsn, gtsu, gtsv, gtui, gtvi, & ! Partial steps for top cell (ISF)176 & rhd, gru , grv , grui, grvi ) ! of t, s, rd at the first ocean level177 !!jc: fs simplification178 167 179 168 ua(:,:,:) = 0._wp ! set dynamics trends to zero -
NEMO/releases/release-4.0/src/OCE/stpctl.F90
r10570 r11533 96 96 IF( ln_zad_Aimp ) THEN 97 97 istatus = NF90_DEF_VAR( idrun, 'abs_wi_max', NF90_DOUBLE, (/ idtime /), idw1 ) 98 istatus = NF90_DEF_VAR( idrun, 'C u_max', NF90_DOUBLE, (/ idtime /), idc1 )98 istatus = NF90_DEF_VAR( idrun, 'Cf_max', NF90_DOUBLE, (/ idtime /), idc1 ) 99 99 ENDIF 100 100 istatus = NF90_ENDDEF(idrun) … … 123 123 IF( ln_zad_Aimp ) THEN 124 124 zmax(8) = MAXVAL( ABS( wi(:,:,:) ) , mask = wmask(:,:,:) == 1._wp ) ! implicit vertical vel. max 125 zmax(9) = MAXVAL( Cu_adv(:,:,:) , mask = tmask(:,:,:) == 1._wp ) ! cell Courant no. max125 zmax(9) = MAXVAL( Cu_adv(:,:,:) , mask = tmask(:,:,:) == 1._wp ) ! partitioning coeff. max 126 126 ENDIF 127 127 !
Note: See TracChangeset
for help on using the changeset viewer.