Changeset 7508 for branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/DYN
- Timestamp:
- 2016-12-19T13:15:59+01:00 (7 years ago)
- Location:
- branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/DYN
- Files:
-
- 12 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/DYN/dynbfr.F90
r6748 r7508 47 47 INTEGER, INTENT(in) :: kt ! ocean time-step index 48 48 !! 49 INTEGER :: j i, jj ! dummy loop indexes49 INTEGER :: jk, ji, jj ! dummy loop indexes 50 50 INTEGER :: ikbu, ikbv ! local integers 51 51 REAL(wp) :: zm1_2dt ! local scalar … … 65 65 IF( l_trddyn ) THEN ! trends: store the input trends 66 66 CALL wrk_alloc( jpi,jpj,jpk, ztrdu, ztrdv ) 67 !$OMP PARALLEL WORKSHARE 68 ztrdu(:,:,:) = ua(:,:,:) 69 ztrdv(:,:,:) = va(:,:,:) 70 !$OMP END PARALLEL WORKSHARE 67 !$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 68 DO jk = 1, jpk 69 DO jj = 1, jpj 70 DO ji = 1, jpi 71 ztrdu(ji,jj,jk) = ua(ji,jj,jk) 72 ztrdv(ji,jj,jk) = va(ji,jj,jk) 73 END DO 74 END DO 75 END DO 71 76 ENDIF 72 77 … … 102 107 ! 103 108 IF( l_trddyn ) THEN ! trends: send trends to trddyn for further diagnostics 104 !$OMP PARALLEL WORKSHARE 105 ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:) 106 ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:) 107 !$OMP END PARALLEL WORKSHARE 109 !$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 110 DO jk = 1, jpk 111 DO jj = 1, jpj 112 DO ji = 1, jpi 113 ztrdu(ji,jj,jk) = ua(ji,jj,jk) - ztrdu(ji,jj,jk) 114 ztrdv(ji,jj,jk) = va(ji,jj,jk) - ztrdv(ji,jj,jk) 115 END DO 116 END DO 117 END DO 108 118 CALL trd_dyn( ztrdu(:,:,:), ztrdv(:,:,:), jpdyn_bfr, kt ) 109 119 CALL wrk_dealloc( jpi,jpj,jpk, ztrdu, ztrdv ) -
branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/DYN/dynhpg.F90
r7037 r7508 84 84 !!---------------------------------------------------------------------- 85 85 INTEGER, INTENT(in) :: kt ! ocean time-step index 86 INTEGER :: jk, jj, ji 86 87 REAL(wp), POINTER, DIMENSION(:,:,:) :: ztrdu, ztrdv 87 88 !!---------------------------------------------------------------------- … … 90 91 ! 91 92 IF( l_trddyn ) THEN ! Temporary saving of ua and va trends (l_trddyn) 92 CALL wrk_alloc( jpi,jpj,jpk, ztrdu, ztrdv ) 93 !$OMP PARALLEL WORKSHARE 94 ztrdu(:,:,:) = ua(:,:,:) 95 ztrdv(:,:,:) = va(:,:,:) 96 !$OMP END PARALLEL WORKSHARE 93 CALL wrk_alloc( jpi,jpj,jpk, ztrdu, ztrdv ) 94 !$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 95 DO jk = 1, jpk 96 DO jj = 1, jpj 97 DO ji = 1, jpi 98 ztrdu(ji,jj,jk) = ua(ji,jj,jk) 99 ztrdv(ji,jj,jk) = va(ji,jj,jk) 100 END DO 101 END DO 102 END DO 97 103 ENDIF 98 104 ! … … 107 113 ! 108 114 IF( l_trddyn ) THEN ! save the hydrostatic pressure gradient trends for momentum trend diagnostics 109 !$OMP PARALLEL WORKSHARE 110 ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:) 111 ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:) 112 !$OMP END PARALLEL WORKSHARE 115 !$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 116 DO jk = 1, jpk 117 DO jj = 1, jpj 118 DO ji = 1, jpi 119 ztrdu(ji,jj,jk) = ua(ji,jj,jk) - ztrdu(ji,jj,jk) 120 ztrdv(ji,jj,jk) = va(ji,jj,jk) - ztrdv(ji,jj,jk) 121 END DO 122 END DO 123 END DO 113 124 CALL trd_dyn( ztrdu, ztrdv, jpdyn_hpg, kt ) 114 125 CALL wrk_dealloc( jpi,jpj,jpk, ztrdu, ztrdv ) … … 202 213 ! 203 214 ! initialisation of ice shelf load 204 IF ( .NOT. ln_isfcav ) riceload(:,:)=0.0 215 IF ( .NOT. ln_isfcav ) THEN 216 !$OMP PARALLEL DO schedule(static) private(jj, ji) 217 DO jj = 1, jpj 218 DO ji = 1, jpi 219 riceload(ji,jj)=0.0 220 END DO 221 END DO 222 END IF 205 223 IF ( ln_isfcav ) THEN 206 224 CALL wrk_alloc( jpi,jpj, 2, ztstop) … … 216 234 217 235 ! assume water displaced by the ice shelf is at T=-1.9 and S=34.4 (rude) 218 ztstop(:,:,1)=-1.9_wp ; ztstop(:,:,2)=34.4_wp 236 !$OMP PARALLEL DO schedule(static) private(jj, ji) 237 DO jj = 1, jpj 238 DO ji = 1, jpi 239 ztstop(ji,jj,1)=-1.9_wp 240 ztstop(ji,jj,2)=34.4_wp 241 END DO 242 END DO 219 243 220 244 ! compute density of the water displaced by the ice shelf … … 244 268 END DO 245 269 END DO 246 !$OMP WORKSHARE 247 riceload(:,:)=ziceload(:,:) ! need to be saved for diaar5 248 !$OMP END WORKSHARE NOWAIT 270 !$OMP DO schedule(static) private(jj, ji) 271 DO jj = 1, jpj 272 DO ji = 1, jpi 273 riceload(ji,jj)=ziceload(ji,jj) ! need to be saved for diaar5 274 END DO 275 END DO 249 276 !$OMP END PARALLEL 250 277 -
branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/DYN/dynkeg.F90
r6748 r7508 92 92 IF( l_trddyn ) THEN ! Save ua and va trends 93 93 CALL wrk_alloc( jpi,jpj,jpk, ztrdu, ztrdv ) 94 !$OMP PARALLEL WORKSHARE 95 ztrdu(:,:,:) = ua(:,:,:) 96 ztrdv(:,:,:) = va(:,:,:) 97 !$OMP END PARALLEL WORKSHARE 94 !$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 95 DO jk = 1, jpk 96 DO jj = 1, jpj 97 DO ji = 1, jpi 98 ztrdu(ji,jj,jk) = ua(ji,jj,jk) 99 ztrdv(ji,jj,jk) = va(ji,jj,jk) 100 END DO 101 END DO 102 END DO 98 103 ENDIF 99 100 zhke(:,:,jpk) = 0._wp 104 !$OMP PARALLEL DO schedule(static) private(jj, ji) 105 DO jj = 1, jpj 106 DO ji = 1, jpi 107 zhke(ji,jj,jpk) = 0._wp 108 END DO 109 END DO 101 110 102 111 SELECT CASE ( kscheme ) !== Horizontal kinetic energy at T-point ==! … … 149 158 ! 150 159 IF( l_trddyn ) THEN ! save the Kinetic Energy trends for diagnostic 151 !$OMP PARALLEL WORKSHARE 152 ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:) 153 ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:) 154 !$OMP END PARALLEL WORKSHARE 160 !$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 161 DO jk = 1, jpk 162 DO jj = 1, jpj 163 DO ji = 1, jpi 164 ztrdu(ji,jj,jk) = ua(ji,jj,jk) - ztrdu(ji,jj,jk) 165 ztrdv(ji,jj,jk) = va(ji,jj,jk) - ztrdv(ji,jj,jk) 166 END DO 167 END DO 168 END DO 155 169 CALL trd_dyn( ztrdu, ztrdv, jpdyn_keg, kt ) 156 170 CALL wrk_dealloc( jpi,jpj,jpk, ztrdu, ztrdv ) -
branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/DYN/dynldf.F90
r6748 r7508 61 61 !!---------------------------------------------------------------------- 62 62 INTEGER, INTENT(in) :: kt ! ocean time-step index 63 INTEGER :: jk, jj, ji 63 64 ! 64 65 REAL(wp), POINTER, DIMENSION(:,:,:) :: ztrdu, ztrdv … … 69 70 IF( l_trddyn ) THEN ! temporary save of momentum trends 70 71 CALL wrk_alloc( jpi,jpj,jpk, ztrdu, ztrdv ) 71 !$OMP PARALLEL WORKSHARE 72 ztrdu(:,:,:) = ua(:,:,:) 73 ztrdv(:,:,:) = va(:,:,:) 74 !$OMP END PARALLEL WORKSHARE 72 !$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 73 DO jk = 1, jpk 74 DO jj = 1, jpj 75 DO ji = 1, jpi 76 ztrdu(ji,jj,jk) = ua(ji,jj,jk) 77 ztrdv(ji,jj,jk) = va(ji,jj,jk) 78 END DO 79 END DO 80 END DO 75 81 ENDIF 76 82 … … 84 90 85 91 IF( l_trddyn ) THEN ! save the horizontal diffusive trends for further diagnostics 86 !$OMP PARALLEL WORKSHARE 87 ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:) 88 ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:) 89 !$OMP END PARALLEL WORKSHARE 92 !$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 93 DO jk = 1, jpk 94 DO jj = 1, jpj 95 DO ji = 1, jpi 96 ztrdu(ji,jj,jk) = ua(ji,jj,jk) - ztrdu(ji,jj,jk) 97 ztrdv(ji,jj,jk) = va(ji,jj,jk) - ztrdv(ji,jj,jk) 98 END DO 99 END DO 100 END DO 90 101 CALL trd_dyn( ztrdu, ztrdv, jpdyn_ldf, kt ) 91 102 CALL wrk_dealloc( jpi,jpj,jpk, ztrdu, ztrdv ) -
branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/DYN/dynldf_lap_blp.F90
r7037 r7508 132 132 !!---------------------------------------------------------------------- 133 133 INTEGER , INTENT(in ) :: kt ! ocean time-step index 134 INTEGER :: jk, jj, ji 134 135 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in ) :: pub, pvb ! before velocity fields 135 136 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) :: pua, pva ! momentum trend … … 148 149 ENDIF 149 150 ! 150 !$OMP PARALLEL WORKSHARE 151 zulap(:,:,:) = 0._wp 152 zvlap(:,:,:) = 0._wp 153 !$OMP END PARALLEL WORKSHARE 151 !$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 152 DO jk = 1, jpk 153 DO jj = 1, jpj 154 DO ji = 1, jpi 155 zulap(ji,jj,jk) = 0._wp 156 zvlap(ji,jj,jk) = 0._wp 157 END DO 158 END DO 159 END DO 154 160 ! 155 161 CALL dyn_ldf_lap( kt, pub, pvb, zulap, zvlap, 1 ) ! rotated laplacian applied to ptb (output in zlap) -
branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/DYN/dynnxt.F90
r7037 r7508 116 116 ! Compute actual transport and replace it with ts estimate at "after" time step 117 117 !$OMP PARALLEL 118 !$OMP WORKSHARE 119 zue(:,:) = e3u_a(:,:,1) * ua(:,:,1) * umask(:,:,1) 120 zve(:,:) = e3v_a(:,:,1) * va(:,:,1) * vmask(:,:,1) 121 !$OMP END WORKSHARE 118 !$OMP DO schedule(static) private(jj, ji) 119 DO jj = 1, jpj 120 DO ji = 1, jpi 121 zue(ji,jj) = e3u_a(ji,jj,1) * ua(ji,jj,1) * umask(ji,jj,1) 122 zve(ji,jj) = e3v_a(ji,jj,1) * va(ji,jj,1) * vmask(ji,jj,1) 123 END DO 124 END DO 122 125 DO jk = 2, jpkm1 123 126 !$OMP DO schedule(static) private(jj,ji) … … 176 179 ! 177 180 IF( ln_dyn_trd ) THEN ! 3D output: total momentum trends 178 !$OMP PARALLEL WORKSHARE 179 zua(:,:,:) = ( ua(:,:,:) - ub(:,:,:) ) * z1_2dt 180 zva(:,:,:) = ( va(:,:,:) - vb(:,:,:) ) * z1_2dt 181 !$OMP END PARALLEL WORKSHARE 182 CALL iom_put( "utrd_tot", zua ) ! total momentum trends, except the asselin time filter 183 CALL iom_put( "vtrd_tot", zva ) 181 !$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 182 DO jk = 1, jpk 183 DO jj = 1, jpj 184 DO ji = 1, jpi 185 zua(ji,jj,jk) = ( ua(ji,jj,jk) - ub(ji,jj,jk) ) * z1_2dt 186 zva(ji,jj,jk) = ( va(ji,jj,jk) - vb(ji,jj,jk) ) * z1_2dt 187 END DO 188 END DO 189 END DO 190 CALL iom_put( "utrd_tot", zua ) ! total momentum trends, except the asselin time filter 191 CALL iom_put( "vtrd_tot", zva ) 184 192 ENDIF 185 193 ! 186 !$OMP PARALLEL WORKSHARE 187 zua(:,:,:) = un(:,:,:) ! save the now velocity before the asselin filter 188 zva(:,:,:) = vn(:,:,:) ! (caution: there will be a shift by 1 timestep in the 189 !$OMP END PARALLEL WORKSHARE 194 !$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 195 DO jk = 1, jpk 196 DO jj = 1, jpj 197 DO ji = 1, jpi 198 zua(ji,jj,jk) = un(ji,jj,jk) ! save the now velocity before the asselin filter 199 zva(ji,jj,jk) = vn(ji,jj,jk) ! (caution: there will be a shift by 1 timestep in the 200 END DO 201 END DO 202 END DO 190 203 ! ! computation of the asselin filter trends) 191 204 ENDIF … … 318 331 END DO 319 332 END DO 320 !$OMP WORKSHARE 321 e3u_b(:,:,1:jpkm1) = ze3u_f(:,:,1:jpkm1) ! e3u_b <-- filtered scale factor 322 e3v_b(:,:,1:jpkm1) = ze3v_f(:,:,1:jpkm1) 323 !$OMP END WORKSHARE NOWAIT 333 !$OMP DO schedule(static) private(jj, ji) 334 DO jj = 1, jpj 335 DO ji = 1, jpi 336 e3u_b(ji,jj,1:jpkm1) = ze3u_f(ji,jj,1:jpkm1) ! e3u_b <-- filtered scale factor 337 e3v_b(ji,jj,1:jpkm1) = ze3v_f(ji,jj,1:jpkm1) 338 END DO 339 END DO 324 340 !$OMP END PARALLEL 325 341 ! … … 333 349 ! Doing it here also means that asselin filter contribution is removed 334 350 !$OMP PARALLEL 335 !$OMP WORKSHARE 336 zue(:,:) = e3u_b(:,:,1) * ub(:,:,1) * umask(:,:,1) 337 zve(:,:) = e3v_b(:,:,1) * vb(:,:,1) * vmask(:,:,1) 338 !$OMP END WORKSHARE 351 !$OMP DO schedule(static) private(jj, ji) 352 DO jj = 1, jpj 353 DO ji = 1, jpi 354 zue(ji,jj) = e3u_b(ji,jj,1) * ub(ji,jj,1) * umask(ji,jj,1) 355 zve(ji,jj) = e3v_b(ji,jj,1) * vb(ji,jj,1) * vmask(ji,jj,1) 356 END DO 357 END DO 339 358 DO jk = 2, jpkm1 340 359 !$OMP DO schedule(static) private(jj, ji) … … 364 383 IF(.NOT.ln_linssh ) THEN 365 384 !$OMP PARALLEL 366 !$OMP WORKSHARE 367 hu_b(:,:) = e3u_b(:,:,1) * umask(:,:,1) 368 hv_b(:,:) = e3v_b(:,:,1) * vmask(:,:,1) 369 !$OMP END WORKSHARE 385 !$OMP DO schedule(static) private(jj, ji) 386 DO jj = 1, jpj 387 DO ji = 1, jpi 388 hu_b(ji,jj) = e3u_b(ji,jj,1) * umask(ji,jj,1) 389 hv_b(ji,jj) = e3v_b(ji,jj,1) * vmask(ji,jj,1) 390 END DO 391 END DO 370 392 DO jk = 2, jpkm1 371 393 !$OMP DO schedule(static) private(jj, ji) … … 378 400 !$OMP END DO 379 401 END DO 380 !$OMP WORKSHARE 381 r1_hu_b(:,:) = ssumask(:,:) / ( hu_b(:,:) + 1._wp - ssumask(:,:) ) 382 r1_hv_b(:,:) = ssvmask(:,:) / ( hv_b(:,:) + 1._wp - ssvmask(:,:) ) 383 !$OMP END WORKSHARE 402 !$OMP DO schedule(static) private(jj, ji) 403 DO jj = 1, jpj 404 DO ji = 1, jpi 405 r1_hu_b(ji,jj) = ssumask(ji,jj) / ( hu_b(ji,jj) + 1._wp - ssumask(ji,jj) ) 406 r1_hv_b(ji,jj) = ssvmask(ji,jj) / ( hv_b(ji,jj) + 1._wp - ssvmask(ji,jj) ) 407 END DO 408 END DO 384 409 !$OMP END PARALLEL 385 410 ENDIF 386 411 ! 387 412 !$OMP PARALLEL 388 !$OMP WORKSHARE 389 un_b(:,:) = e3u_a(:,:,1) * un(:,:,1) * umask(:,:,1) 390 ub_b(:,:) = e3u_b(:,:,1) * ub(:,:,1) * umask(:,:,1) 391 vn_b(:,:) = e3v_a(:,:,1) * vn(:,:,1) * vmask(:,:,1) 392 vb_b(:,:) = e3v_b(:,:,1) * vb(:,:,1) * vmask(:,:,1) 393 !$OMP END WORKSHARE 413 !$OMP DO schedule(static) private(jj, ji) 414 DO jj = 1, jpj 415 DO ji = 1, jpi 416 un_b(ji,jj) = e3u_a(ji,jj,1) * un(ji,jj,1) * umask(ji,jj,1) 417 ub_b(ji,jj) = e3u_b(ji,jj,1) * ub(ji,jj,1) * umask(ji,jj,1) 418 vn_b(ji,jj) = e3v_a(ji,jj,1) * vn(ji,jj,1) * vmask(ji,jj,1) 419 vb_b(ji,jj) = e3v_b(ji,jj,1) * vb(ji,jj,1) * vmask(ji,jj,1) 420 END DO 421 END DO 394 422 DO jk = 2, jpkm1 395 423 !$OMP DO schedule(static) private(jj, ji) … … 403 431 END DO 404 432 END DO 405 !$OMP WORKSHARE 406 un_b(:,:) = un_b(:,:) * r1_hu_a(:,:) 407 vn_b(:,:) = vn_b(:,:) * r1_hv_a(:,:) 408 ub_b(:,:) = ub_b(:,:) * r1_hu_b(:,:) 409 vb_b(:,:) = vb_b(:,:) * r1_hv_b(:,:) 410 !$OMP END WORKSHARE NOWAIT 433 !$OMP DO schedule(static) private(jj, ji) 434 DO jj = 1, jpj 435 DO ji = 1, jpi 436 un_b(ji,jj) = un_b(ji,jj) * r1_hu_a(ji,jj) 437 vn_b(ji,jj) = vn_b(ji,jj) * r1_hv_a(ji,jj) 438 ub_b(ji,jj) = ub_b(ji,jj) * r1_hu_b(ji,jj) 439 vb_b(ji,jj) = vb_b(ji,jj) * r1_hv_b(ji,jj) 440 END DO 441 END DO 411 442 !$OMP END PARALLEL 412 443 ! … … 416 447 ENDIF 417 448 IF( l_trddyn ) THEN ! 3D output: asselin filter trends on momentum 418 !$OMP PARALLEL WORKSHARE 419 zua(:,:,:) = ( ub(:,:,:) - zua(:,:,:) ) * z1_2dt 420 zva(:,:,:) = ( vb(:,:,:) - zva(:,:,:) ) * z1_2dt 421 !$OMP END PARALLEL WORKSHARE 422 CALL trd_dyn( zua, zva, jpdyn_atf, kt ) 449 !$OMP DO schedule(static) private(jk, jj, ji) 450 DO jk = 1, jpkm1 451 DO jj = 1, jpj 452 DO ji = 1, jpi 453 zua(ji,jj,jk) = ( ub(ji,jj,jk) - zua(ji,jj,jk) ) * z1_2dt 454 zva(ji,jj,jk) = ( vb(ji,jj,jk) - zva(ji,jj,jk) ) * z1_2dt 455 END DO 456 END DO 457 END DO 458 CALL trd_dyn( zua, zva, jpdyn_atf, kt ) 423 459 ENDIF 424 460 ! -
branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg.F90
r7037 r7508 83 83 IF( l_trddyn ) THEN ! temporary save of ta and sa trends 84 84 CALL wrk_alloc( jpi,jpj,jpk, ztrdu, ztrdv ) 85 !$OMP PARALLEL WORKSHARE 86 ztrdu(:,:,:) = ua(:,:,:) 87 ztrdv(:,:,:) = va(:,:,:) 88 !$OMP END PARALLEL WORKSHARE 85 !$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 86 DO jk = 1, jpk 87 DO jj = 1, jpj 88 DO ji = 1, jpi 89 ztrdu(ji,jj,jk) = ua(ji,jj,jk) 90 ztrdv(ji,jj,jk) = va(ji,jj,jk) 91 END DO 92 END DO 93 END DO 89 94 ENDIF 90 95 ! … … 134 139 zgrau0r = - grav * r1_rau0 135 140 !$OMP PARALLEL 136 !$OMP WORKSHARE 137 zpice(:,:) = ( zintp * snwice_mass(:,:) + ( 1.- zintp ) * snwice_mass_b(:,:) ) * zgrau0r 138 !$OMP END WORKSHARE 141 !$OMP DO schedule(static) private(jj, ji) 142 DO jj = 1, jpj 143 DO ji = 1, jpi 144 zpice(ji,jj) = ( zintp * snwice_mass(ji,jj) + ( 1.- zintp ) * snwice_mass_b(ji,jj) ) * zgrau0r 145 END DO 146 END DO 139 147 !$OMP DO schedule(static) private(jj, ji) 140 148 DO jj = 2, jpjm1 … … 170 178 ! 171 179 IF( l_trddyn ) THEN ! save the surface pressure gradient trends for further diagnostics 172 !$OMP PARALLEL WORKSHARE 173 ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:) 174 ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:) 175 !$OMP END PARALLEL WORKSHARE 180 !$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 181 DO jk = 1, jpk 182 DO jj = 1, jpj 183 DO ji = 1, jpi 184 ztrdu(ji,jj,jk) = ua(ji,jj,jk) - ztrdu(ji,jj,jk) 185 ztrdv(ji,jj,jk) = va(ji,jj,jk) - ztrdv(ji,jj,jk) 186 END DO 187 END DO 188 END DO 176 189 CALL trd_dyn( ztrdu, ztrdv, jpdyn_spg, kt ) 177 190 CALL wrk_dealloc( jpi,jpj,jpk, ztrdu, ztrdv ) -
branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg_ts.F90
r7037 r7508 244 244 ! 245 245 !$OMP PARALLEL 246 !$OMP WORKSHARE247 ftne(1,:) = 0._wp ; ftnw(1,:) = 0._wp ; ftse(1,:) = 0._wp ; ftsw(1,:) = 0._wp248 !$OMP END WORKSHARE 249 246 !$OMP DO schedule(static) private(jj) 247 DO jj = 1, jpj 248 ftne(1,jj) = 0._wp ; ftnw(1,jj) = 0._wp ; ftse(1,jj) = 0._wp ; ftsw(1,jj) = 0._wp 249 END DO 250 250 !$OMP DO schedule(static) private(jj, ji) 251 251 DO jj = 2, jpj … … 261 261 ! 262 262 ELSE !== all other schemes (ENE, ENS, MIX) 263 !$OMP PARALLEL WORKSHARE 264 zwz(:,:) = 0._wp 265 zhf(:,:) = 0._wp 266 !$OMP END PARALLEL WORKSHARE 263 !$OMP PARALLEL DO schedule(static) private(jj, ji) 264 DO jj = 1, jpj 265 DO ji = 1, jpi 266 zwz(ji,jj) = 0._wp 267 zhf(ji,jj) = 0._wp 268 END DO 269 END DO 267 270 IF ( .not. ln_sco ) THEN 268 271 … … 276 279 ! zhf(:,:) = gdepw_0(:,:,jk+1) 277 280 ELSE 278 !$OMP PARALLEL WORKSHARE 279 zhf(:,:) = hbatf(:,:) 280 !$OMP END PARALLEL WORKSHARE 281 !$OMP PARALLEL DO schedule(static) private(jj, ji) 282 DO jj = 1, jpj 283 DO ji = 1, jpi 284 zhf(ji,jj) = hbatf(ji,jj) 285 END DO 286 END DO 281 287 END IF 282 288 … … 308 314 END DO 309 315 END DO 310 !$OMP WORKSHARE 311 zwz(:,:) = ff(:,:) * zwz(:,:) 312 !$OMP END WORKSHARE NOWAIT 316 !$OMP DO schedule(static) private(jj, ji) 317 DO jj = 1, jpj 318 DO ji = 1, jpi 319 zwz(ji,jj) = ff(ji,jj) * zwz(ji,jj) 320 END DO 321 END DO 313 322 !$OMP END PARALLEL 314 323 ENDIF … … 330 339 ! ! -------------------------------------------------- 331 340 !$OMP PARALLEL 332 !$OMP WORKSHARE 333 zu_frc(:,:) = 0._wp 334 zv_frc(:,:) = 0._wp 335 !$OMP END WORKSHARE 341 !$OMP DO schedule(static) private(jj, ji) 342 DO jj = 1, jpj 343 DO ji = 1, jpi 344 zu_frc(ji,jj) = 0._wp 345 zv_frc(ji,jj) = 0._wp 346 END DO 347 END DO 336 348 ! 337 349 DO jk = 1, jpkm1 … … 345 357 END DO 346 358 ! 347 !$OMP WORKSHARE 348 zu_frc(:,:) = zu_frc(:,:) * r1_hu_n(:,:) 349 zv_frc(:,:) = zv_frc(:,:) * r1_hv_n(:,:) 350 !$OMP END WORKSHARE 351 ! 359 !$OMP DO schedule(static) private(jj, ji) 360 DO jj = 1, jpj 361 DO ji = 1, jpi 362 zu_frc(ji,jj) = zu_frc(ji,jj) * r1_hu_n(ji,jj) 363 zv_frc(ji,jj) = zv_frc(ji,jj) * r1_hv_n(ji,jj) 364 END DO 365 END DO 352 366 ! 353 367 ! !* baroclinic momentum trend (remove the vertical mean trend) … … 364 378 ! !* barotropic Coriolis trends (vorticity scheme dependent) 365 379 ! ! -------------------------------------------------------- 366 !$OMP WORKSHARE 367 zwx(:,:) = un_b(:,:) * hu_n(:,:) * e2u(:,:) ! now fluxes 368 zwy(:,:) = vn_b(:,:) * hv_n(:,:) * e1v(:,:) 369 !$OMP END WORKSHARE NOWAIT 380 !$OMP DO schedule(static) private(jj, ji) 381 DO jj = 1, jpj 382 DO ji = 1, jpi 383 zwx(ji,jj) = un_b(ji,jj) * hu_n(ji,jj) * e2u(ji,jj) ! now fluxes 384 zwy(ji,jj) = vn_b(ji,jj) * hv_n(ji,jj) * e1v(ji,jj) 385 END DO 386 END DO 370 387 !$OMP END PARALLEL 371 388 ! … … 419 436 IF( ln_wd ) THEN ! Calculating and applying W/D gravity filters 420 437 !$OMP PARALLEL 421 !$OMP WORKSHARE 422 wduflt1(:,:) = 1.0_wp 423 wdvflt1(:,:) = 1.0_wp 424 !$OMP END WORKSHARE 438 !$OMP DO schedule(static) private(jj, ji) 439 DO jj = 1, jpj 440 DO ji = 1, jpi 441 wduflt1(ji,jj) = 1.0_wp 442 wdvflt1(ji,jj) = 1.0_wp 443 END DO 444 END DO 425 445 !$OMP DO schedule(static) private(jj,ji,ll_tmp1,ll_tmp2) 426 446 DO jj = 2, jpjm1 … … 529 549 END DO 530 550 ELSE 531 !$OMP PARALLEL WORKSHARE 532 zu_frc(:,:) = zu_frc(:,:) + r1_hu_n(:,:) * bfrua(:,:) * zwx(:,:) 533 zv_frc(:,:) = zv_frc(:,:) + r1_hv_n(:,:) * bfrva(:,:) * zwy(:,:) 534 !$OMP END PARALLEL WORKSHARE 551 !$OMP PARALLEL DO schedule(static) private(jj,ji) 552 DO jj = 1, jpj 553 DO ji = 1, jpi 554 zu_frc(ji,jj) = zu_frc(ji,jj) + r1_hu_n(ji,jj) * bfrua(ji,jj) * zwx(ji,jj) 555 zv_frc(ji,jj) = zv_frc(ji,jj) + r1_hv_n(ji,jj) * bfrva(ji,jj) * zwy(ji,jj) 556 END DO 557 END DO 535 558 END IF 536 559 ! … … 559 582 ! 560 583 ! Note that the "unclipped" top friction parameter is used even with explicit drag 561 !$OMP PARALLEL WORKSHARE 562 zu_frc(:,:) = zu_frc(:,:) + r1_hu_n(:,:) * tfrua(:,:) * zwx(:,:) 563 zv_frc(:,:) = zv_frc(:,:) + r1_hv_n(:,:) * tfrva(:,:) * zwy(:,:) 564 !$OMP END PARALLEL WORKSHARE 584 !$OMP PARALLEL DO schedule(static) private(jj,ji) 585 DO jj = 1, jpj 586 DO ji = 1, jpi 587 zu_frc(ji,jj) = zu_frc(ji,jj) + r1_hu_n(ji,jj) * tfrua(ji,jj) * zwx(ji,jj) 588 zv_frc(ji,jj) = zv_frc(ji,jj) + r1_hv_n(ji,jj) * tfrva(ji,jj) * zwy(ji,jj) 589 END DO 590 END DO 565 591 ! 566 592 IF (ln_bt_fw) THEN ! Add wind forcing 567 !$OMP PARALLEL WORKSHARE 568 zu_frc(:,:) = zu_frc(:,:) + zraur * utau(:,:) * r1_hu_n(:,:) 569 zv_frc(:,:) = zv_frc(:,:) + zraur * vtau(:,:) * r1_hv_n(:,:) 570 !$OMP END PARALLEL WORKSHARE 593 !$OMP PARALLEL DO schedule(static) private(jj,ji) 594 DO jj = 1, jpj 595 DO ji = 1, jpi 596 zu_frc(ji,jj) = zu_frc(ji,jj) + zraur * utau(ji,jj) * r1_hu_n(ji,jj) 597 zv_frc(ji,jj) = zv_frc(ji,jj) + zraur * vtau(ji,jj) * r1_hv_n(ji,jj) 598 END DO 599 END DO 571 600 ELSE 572 !$OMP PARALLEL WORKSHARE 573 zu_frc(:,:) = zu_frc(:,:) + zraur * z1_2 * ( utau_b(:,:) + utau(:,:) ) * r1_hu_n(:,:) 574 zv_frc(:,:) = zv_frc(:,:) + zraur * z1_2 * ( vtau_b(:,:) + vtau(:,:) ) * r1_hv_n(:,:) 575 !$OMP END PARALLEL WORKSHARE 601 !$OMP PARALLEL DO schedule(static) private(jj,ji) 602 DO jj = 1, jpj 603 DO ji = 1, jpi 604 zu_frc(ji,jj) = zu_frc(ji,jj) + zraur * z1_2 * ( utau_b(ji,jj) + utau(ji,jj) ) * r1_hu_n(ji,jj) 605 zv_frc(ji,jj) = zv_frc(ji,jj) + zraur * z1_2 * ( vtau_b(ji,jj) + vtau(ji,jj) ) * r1_hv_n(ji,jj) 606 END DO 607 END DO 576 608 ENDIF 577 609 ! … … 605 637 ! ! Surface net water flux and rivers 606 638 IF (ln_bt_fw) THEN 607 !$OMP PARALLEL WORKSHARE 608 zssh_frc(:,:) = zraur * ( emp(:,:) - rnf(:,:) + fwfisf(:,:) ) 609 !$OMP END PARALLEL WORKSHARE 639 !$OMP PARALLEL DO schedule(static) private(jj,ji) 640 DO jj = 1, jpj 641 DO ji = 1, jpi 642 zssh_frc(ji,jj) = zraur * ( emp(ji,jj) - rnf(ji,jj) + fwfisf(ji,jj) ) 643 END DO 644 END DO 610 645 ELSE 611 !$OMP PARALLEL WORKSHARE 612 zssh_frc(:,:) = zraur * z1_2 * ( emp(:,:) + emp_b(:,:) - rnf(:,:) - rnf_b(:,:) & 613 & + fwfisf(:,:) + fwfisf_b(:,:) ) 614 !$OMP END PARALLEL WORKSHARE 646 !$OMP PARALLEL DO schedule(static) private(jj,ji) 647 DO jj = 1, jpj 648 DO ji = 1, jpi 649 zssh_frc(ji,jj) = zraur * z1_2 * ( emp(ji,jj) + emp_b(ji,jj) - rnf(ji,jj) - rnf_b(ji,jj) & 650 & + fwfisf(ji,jj) + fwfisf_b(ji,jj) ) 651 END DO 652 END DO 615 653 ENDIF 616 654 #if defined key_asminc 617 655 ! ! Include the IAU weighted SSH increment 618 656 IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN 619 !$OMP PARALLEL WORKSHARE 620 zssh_frc(:,:) = zssh_frc(:,:) - ssh_iau(:,:) 621 !$OMP END PARALLEL WORKSHARE 657 !$OMP PARALLEL DO schedule(static) private(jj,ji) 658 DO jj = 1, jpj 659 DO ji = 1, jpi 660 zssh_frc(ji,jj) = zssh_frc(ji,jj) - ssh_iau(ji,jj) 661 END DO 662 END DO 622 663 ENDIF 623 664 #endif … … 637 678 ! Initialize barotropic variables: 638 679 IF( ll_init )THEN 639 !$OMP PARALLEL WORKSHARE 640 sshbb_e(:,:) = 0._wp 641 ubb_e (:,:) = 0._wp 642 vbb_e (:,:) = 0._wp 643 sshb_e (:,:) = 0._wp 644 ub_e (:,:) = 0._wp 645 vb_e (:,:) = 0._wp 646 !$OMP END PARALLEL WORKSHARE 680 !$OMP PARALLEL DO schedule(static) private(jj,ji) 681 DO jj = 1, jpj 682 DO ji = 1, jpi 683 sshbb_e(ji,jj) = 0._wp 684 ubb_e (ji,jj) = 0._wp 685 vbb_e (ji,jj) = 0._wp 686 sshb_e (ji,jj) = 0._wp 687 ub_e (ji,jj) = 0._wp 688 vb_e (ji,jj) = 0._wp 689 END DO 690 END DO 647 691 ENDIF 648 692 … … 659 703 ! 660 704 IF (ln_bt_fw) THEN ! FORWARD integration: start from NOW fields 661 !$OMP PARALLEL WORKSHARE 662 sshn_e(:,:) = sshn(:,:) 663 un_e (:,:) = un_b(:,:) 664 vn_e (:,:) = vn_b(:,:) 665 ! 666 hu_e (:,:) = hu_n(:,:) 667 hv_e (:,:) = hv_n(:,:) 668 hur_e (:,:) = r1_hu_n(:,:) 669 hvr_e (:,:) = r1_hv_n(:,:) 670 !$OMP END PARALLEL WORKSHARE 705 !$OMP PARALLEL DO schedule(static) private(jj,ji) 706 DO jj = 1, jpj 707 DO ji = 1, jpi 708 sshn_e(ji,jj) = sshn(ji,jj) 709 un_e (ji,jj) = un_b(ji,jj) 710 vn_e (ji,jj) = vn_b(ji,jj) 711 ! 712 hu_e (ji,jj) = hu_n(ji,jj) 713 hv_e (ji,jj) = hv_n(ji,jj) 714 hur_e (ji,jj) = r1_hu_n(ji,jj) 715 hvr_e (ji,jj) = r1_hv_n(ji,jj) 716 END DO 717 END DO 671 718 ELSE ! CENTRED integration: start from BEFORE fields 672 !$OMP PARALLEL WORKSHARE 673 sshn_e(:,:) = sshb(:,:) 674 un_e (:,:) = ub_b(:,:) 675 vn_e (:,:) = vb_b(:,:) 676 ! 677 hu_e (:,:) = hu_b(:,:) 678 hv_e (:,:) = hv_b(:,:) 679 hur_e (:,:) = r1_hu_b(:,:) 680 hvr_e (:,:) = r1_hv_b(:,:) 681 !$OMP END PARALLEL WORKSHARE 719 !$OMP PARALLEL DO schedule(static) private(jj,ji) 720 DO jj = 1, jpj 721 DO ji = 1, jpi 722 sshn_e(ji,jj) = sshb(ji,jj) 723 un_e (ji,jj) = ub_b(ji,jj) 724 vn_e (ji,jj) = vb_b(ji,jj) 725 ! 726 hu_e (ji,jj) = hu_b(ji,jj) 727 hv_e (ji,jj) = hv_b(ji,jj) 728 hur_e (ji,jj) = r1_hu_b(ji,jj) 729 hvr_e (ji,jj) = r1_hv_b(ji,jj) 730 END DO 731 END DO 682 732 ENDIF 683 733 ! … … 685 735 ! 686 736 ! Initialize sums: 687 !$OMP PARALLEL WORKSHARE 688 ua_b (:,:) = 0._wp ! After barotropic velocities (or transport if flux form) 689 va_b (:,:) = 0._wp 690 ssha (:,:) = 0._wp ! Sum for after averaged sea level 691 un_adv(:,:) = 0._wp ! Sum for now transport issued from ts loop 692 vn_adv(:,:) = 0._wp 693 !$OMP END PARALLEL WORKSHARE 737 !$OMP PARALLEL DO schedule(static) private(jj,ji) 738 DO jj = 1, jpj 739 DO ji = 1, jpi 740 ua_b (ji,jj) = 0._wp ! After barotropic velocities (or transport if flux form) 741 va_b (ji,jj) = 0._wp 742 ssha (ji,jj) = 0._wp ! Sum for after averaged sea level 743 un_adv(ji,jj) = 0._wp ! Sum for now transport issued from ts loop 744 vn_adv(ji,jj) = 0._wp 745 END DO 746 END DO 694 747 ! ! ==================== ! 695 748 DO jn = 1, icycle ! sub-time-step loop ! … … 715 768 716 769 ! Extrapolate barotropic velocities at step jit+0.5: 717 !$OMP PARALLEL WORKSHARE 718 ua_e(:,:) = za1 * un_e(:,:) + za2 * ub_e(:,:) + za3 * ubb_e(:,:) 719 va_e(:,:) = za1 * vn_e(:,:) + za2 * vb_e(:,:) + za3 * vbb_e(:,:) 720 !$OMP END PARALLEL WORKSHARE 770 !$OMP PARALLEL DO schedule(static) private(jj,ji) 771 DO jj = 1, jpj 772 DO ji = 1, jpi 773 ua_e(ji,jj) = za1 * un_e(ji,jj) + za2 * ub_e(ji,jj) + za3 * ubb_e(ji,jj) 774 va_e(ji,jj) = za1 * vn_e(ji,jj) + za2 * vb_e(ji,jj) + za3 * vbb_e(ji,jj) 775 END DO 776 END DO 721 777 722 778 IF( .NOT.ln_linssh ) THEN !* Update ocean depth (variable volume case only) … … 724 780 ! Extrapolate Sea Level at step jit+0.5: 725 781 !$OMP PARALLEL 726 !$OMP WORKSHARE 727 zsshp2_e(:,:) = za1 * sshn_e(:,:) + za2 * sshb_e(:,:) + za3 * sshbb_e(:,:) 728 !$OMP END WORKSHARE 782 !$OMP DO schedule(static) private(jj,ji) 783 DO jj = 1, jpj 784 DO ji = 1, jpi 785 zsshp2_e(ji,jj) = za1 * sshn_e(ji,jj) + za2 * sshb_e(ji,jj) + za3 * sshbb_e(ji,jj) 786 END DO 787 END DO 729 788 ! 730 789 !$OMP DO schedule(static) private(jj,ji) … … 743 802 CALL lbc_lnk_multi( zwx, 'U', 1._wp, zwy, 'V', 1._wp ) 744 803 ! 745 !$OMP PARALLEL WORKSHARE 746 zhup2_e (:,:) = hu_0(:,:) + zwx(:,:) ! Ocean depth at U- and V-points 747 zhvp2_e (:,:) = hv_0(:,:) + zwy(:,:) 748 !$OMP END PARALLEL WORKSHARE 804 !$OMP PARALLEL DO schedule(static) private(jj,ji) 805 DO jj = 1, jpj 806 DO ji = 1, jpi 807 zhup2_e (ji,jj) = hu_0(ji,jj) + zwx(ji,jj) ! Ocean depth at U- and V-points 808 zhvp2_e (ji,jj) = hv_0(ji,jj) + zwy(ji,jj) 809 END DO 810 END DO 749 811 IF( ln_wd ) THEN 750 812 !$OMP PARALLEL DO schedule(static) private(jj,ji) 751 DO jj = 1, jpj752 DO ji = 1, jpi ! vector opt.753 zhup2_e(ji,jj) = MAX(zhup2_e (ji,jj), rn_wdmin1)754 zhvp2_e(ji,jj) = MAX(zhvp2_e (ji,jj), rn_wdmin1)755 END DO756 END DO813 DO jj = 1, jpj 814 DO ji = 1, jpi ! vector opt. 815 zhup2_e(ji,jj) = MAX(zhup2_e (ji,jj), rn_wdmin1) 816 zhvp2_e(ji,jj) = MAX(zhvp2_e (ji,jj), rn_wdmin1) 817 END DO 818 END DO 757 819 END IF 758 820 ELSE 759 !$OMP PARALLEL WORKSHARE 760 zhup2_e (:,:) = hu_n(:,:) 761 zhvp2_e (:,:) = hv_n(:,:) 762 !$OMP END PARALLEL WORKSHARE 821 !$OMP PARALLEL DO schedule(static) private(jj,ji) 822 DO jj = 1, jpj 823 DO ji = 1, jpi 824 zhup2_e (ji,jj) = hu_n(ji,jj) 825 zhvp2_e (ji,jj) = hv_n(ji,jj) 826 END DO 827 END DO 763 828 ENDIF 764 829 ! !* after ssh … … 767 832 ! considering fluxes below: 768 833 ! 769 !$OMP PARALLEL WORKSHARE 770 zwx(:,:) = e2u(:,:) * ua_e(:,:) * zhup2_e(:,:) ! fluxes at jn+0.5 771 zwy(:,:) = e1v(:,:) * va_e(:,:) * zhvp2_e(:,:) 772 !$OMP END PARALLEL WORKSHARE 834 !$OMP PARALLEL DO schedule(static) private(jj,ji) 835 DO jj = 1, jpj 836 DO ji = 1, jpi 837 zwx(ji,jj) = e2u(ji,jj) * ua_e(ji,jj) * zhup2_e(ji,jj) ! fluxes at jn+0.5 838 zwy(ji,jj) = e1v(ji,jj) * va_e(ji,jj) * zhvp2_e(ji,jj) 839 END DO 840 END DO 773 841 ! 774 842 #if defined key_agrif … … 802 870 za2 = wgtbtp2(jn) 803 871 !$OMP PARALLEL 804 !$OMP WORKSHARE 805 un_adv(:,:) = un_adv(:,:) + za2 * zwx(:,:) * r1_e2u(:,:) 806 vn_adv(:,:) = vn_adv(:,:) + za2 * zwy(:,:) * r1_e1v(:,:) 807 !$OMP END WORKSHARE NOWAIT 872 !$OMP DO schedule(static) private(jj,ji) 873 DO jj = 1, jpj 874 DO ji = 1, jpi 875 un_adv(ji,jj) = un_adv(ji,jj) + za2 * zwx(ji,jj) * r1_e2u(ji,jj) 876 vn_adv(ji,jj) = vn_adv(ji,jj) + za2 * zwy(ji,jj) * r1_e1v(ji,jj) 877 END DO 878 END DO 879 !$OMP END DO NOWAIT 808 880 ! 809 881 ! Set next sea level: … … 815 887 END DO 816 888 END DO 817 !$OMP WORKSHARE 818 ssha_e(:,:) = ( sshn_e(:,:) - rdtbt * ( zssh_frc(:,:) + zhdiv(:,:) ) ) * ssmask(:,:) 819 !$OMP END WORKSHARE NOWAIT 889 !$OMP DO schedule(static) private(jj,ji) 890 DO jj = 1, jpj 891 DO ji = 1, jpi 892 ssha_e(ji,jj) = ( sshn_e(ji,jj) - rdtbt * ( zssh_frc(ji,jj) + zhdiv(ji,jj) ) ) * ssmask(ji,jj) 893 END DO 894 END DO 895 !$OMP END DO NOWAIT 820 896 !$OMP END PARALLEL 821 897 IF( ln_wd ) THEN … … 872 948 ENDIF 873 949 ! 874 zsshp2_e(:,:) = za0 * ssha_e(:,:) + za1 * sshn_e (:,:) & 875 & + za2 * sshb_e(:,:) + za3 * sshbb_e(:,:) 950 !$OMP PARALLEL DO schedule(static) private(jj,ji) 951 DO jj = 1, jpj 952 DO ji = 1, jpi 953 zsshp2_e(ji,jj) = za0 * ssha_e(ji,jj) + za1 * sshn_e (ji,jj) & 954 & + za2 * sshb_e(ji,jj) + za3 * sshbb_e(ji,jj) 955 END DO 956 END DO 876 957 IF( ln_wd ) THEN ! Calculating and applying W/D gravity filters 877 958 !$OMP PARALLEL 878 !$OMP WORKSHARE 879 wduflt1(:,:) = 1._wp 880 wdvflt1(:,:) = 1._wp 881 !$OMP END WORKSHARE 959 !$OMP DO schedule(static) private(jj,ji) 960 DO jj = 1, jpj 961 DO ji = 1, jpi 962 wduflt1(ji,jj) = 1._wp 963 wdvflt1(ji,jj) = 1._wp 964 END DO 965 END DO 882 966 !$OMP DO schedule(static) private(jj,ji,ll_tmp1,ll_tmp2) 883 967 DO jj = 2, jpjm1 … … 1014 1098 ! 1015 1099 ! Add bottom stresses: 1016 !$OMP PARALLEL WORKSHARE 1017 zu_trd(:,:) = zu_trd(:,:) + bfrua(:,:) * un_e(:,:) * hur_e(:,:) 1018 zv_trd(:,:) = zv_trd(:,:) + bfrva(:,:) * vn_e(:,:) * hvr_e(:,:) 1019 ! 1020 ! Add top stresses: 1021 zu_trd(:,:) = zu_trd(:,:) + tfrua(:,:) * un_e(:,:) * hur_e(:,:) 1022 zv_trd(:,:) = zv_trd(:,:) + tfrva(:,:) * vn_e(:,:) * hvr_e(:,:) 1023 !$OMP END PARALLEL WORKSHARE 1100 !$OMP PARALLEL DO schedule(static) private(jj,ji) 1101 DO jj = 1, jpj 1102 DO ji = 1, jpi 1103 zu_trd(ji,jj) = zu_trd(ji,jj) + bfrua(ji,jj) * un_e(ji,jj) * hur_e(ji,jj) 1104 zv_trd(ji,jj) = zv_trd(ji,jj) + bfrva(ji,jj) * vn_e(ji,jj) * hvr_e(ji,jj) 1105 ! 1106 ! Add top stresses: 1107 zu_trd(ji,jj) = zu_trd(ji,jj) + tfrua(ji,jj) * un_e(ji,jj) * hur_e(ji,jj) 1108 zv_trd(ji,jj) = zv_trd(ji,jj) + tfrva(ji,jj) * vn_e(ji,jj) * hvr_e(ji,jj) 1109 END DO 1110 END DO 1024 1111 ! 1025 1112 ! Surface pressure trend: … … 1109 1196 END DO 1110 1197 ELSE 1111 !$OMP PARALLEL WORKSHARE 1112 hu_e (:,:) = hu_0(:,:) + zsshu_a(:,:) 1113 hv_e (:,:) = hv_0(:,:) + zsshv_a(:,:) 1114 !$OMP END PARALLEL WORKSHARE 1198 !$OMP PARALLEL DO schedule(static) private(jj,ji) 1199 DO jj = 1, jpj 1200 DO ji = 1, jpi 1201 hu_e (ji,jj) = hu_0(ji,jj) + zsshu_a(ji,jj) 1202 hv_e (ji,jj) = hv_0(ji,jj) + zsshv_a(ji,jj) 1203 END DO 1204 END DO 1115 1205 END IF 1116 !$OMP PARALLEL WORKSHARE 1117 hur_e(:,:) = ssumask(:,:) / ( hu_e(:,:) + 1._wp - ssumask(:,:) ) 1118 hvr_e(:,:) = ssvmask(:,:) / ( hv_e(:,:) + 1._wp - ssvmask(:,:) ) 1119 !$OMP END PARALLEL WORKSHARE 1206 !$OMP PARALLEL DO schedule(static) private(jj,ji) 1207 DO jj = 1, jpj 1208 DO ji = 1, jpi 1209 hur_e(ji,jj) = ssumask(ji,jj) / ( hu_e(ji,jj) + 1._wp - ssumask(ji,jj) ) 1210 hvr_e(ji,jj) = ssvmask(ji,jj) / ( hv_e(ji,jj) + 1._wp - ssvmask(ji,jj) ) 1211 END DO 1212 END DO 1120 1213 ! 1121 1214 ENDIF … … 1132 1225 ! !* Swap 1133 1226 ! ! ---- 1134 !$OMP PARALLEL WORKSHARE 1135 ubb_e (:,:) = ub_e (:,:) 1136 ub_e (:,:) = un_e (:,:) 1137 un_e (:,:) = ua_e (:,:) 1138 ! 1139 vbb_e (:,:) = vb_e (:,:) 1140 vb_e (:,:) = vn_e (:,:) 1141 vn_e (:,:) = va_e (:,:) 1142 ! 1143 sshbb_e(:,:) = sshb_e(:,:) 1144 sshb_e (:,:) = sshn_e(:,:) 1145 sshn_e (:,:) = ssha_e(:,:) 1146 !$OMP END PARALLEL WORKSHARE 1227 !$OMP PARALLEL DO schedule(static) private(jj,ji) 1228 DO jj = 1, jpj 1229 DO ji = 1, jpi 1230 ubb_e (ji,jj) = ub_e (ji,jj) 1231 ub_e (ji,jj) = un_e (ji,jj) 1232 un_e (ji,jj) = ua_e (ji,jj) 1233 ! 1234 vbb_e (ji,jj) = vb_e (ji,jj) 1235 vb_e (ji,jj) = vn_e (ji,jj) 1236 vn_e (ji,jj) = va_e (ji,jj) 1237 ! 1238 sshbb_e(ji,jj) = sshb_e(ji,jj) 1239 sshb_e (ji,jj) = sshn_e(ji,jj) 1240 sshn_e (ji,jj) = ssha_e(ji,jj) 1241 END DO 1242 END DO 1147 1243 1148 1244 ! !* Sum over whole bt loop … … 1150 1246 za1 = wgtbtp1(jn) 1151 1247 IF( ln_dynadv_vec .OR. ln_linssh ) THEN ! Sum velocities 1152 !$OMP PARALLEL WORKSHARE 1153 ua_b (:,:) = ua_b (:,:) + za1 * ua_e (:,:) 1154 va_b (:,:) = va_b (:,:) + za1 * va_e (:,:) 1155 !$OMP END PARALLEL WORKSHARE 1248 !$OMP PARALLEL DO schedule(static) private(jj,ji) 1249 DO jj = 1, jpj 1250 DO ji = 1, jpi 1251 ua_b (ji,jj) = ua_b (ji,jj) + za1 * ua_e (ji,jj) 1252 va_b (ji,jj) = va_b (ji,jj) + za1 * va_e (ji,jj) 1253 END DO 1254 END DO 1156 1255 ELSE ! Sum transports 1157 !$OMP PARALLEL WORKSHARE 1158 ua_b (:,:) = ua_b (:,:) + za1 * ua_e (:,:) * hu_e (:,:) 1159 va_b (:,:) = va_b (:,:) + za1 * va_e (:,:) * hv_e (:,:) 1160 !$OMP END PARALLEL WORKSHARE 1256 !$OMP PARALLEL DO schedule(static) private(jj,ji) 1257 DO jj = 1, jpj 1258 DO ji = 1, jpi 1259 ua_b (ji,jj) = ua_b (ji,jj) + za1 * ua_e (ji,jj) * hu_e (ji,jj) 1260 va_b (ji,jj) = va_b (ji,jj) + za1 * va_e (ji,jj) * hv_e (ji,jj) 1261 END DO 1262 END DO 1161 1263 ENDIF 1162 1264 ! ! Sum sea level 1163 !$OMP PARALLEL WORKSHARE 1164 ssha(:,:) = ssha(:,:) + za1 * ssha_e(:,:) 1165 !$OMP END PARALLEL WORKSHARE 1265 !$OMP PARALLEL DO schedule(static) private(jj,ji) 1266 DO jj = 1, jpj 1267 DO ji = 1, jpi 1268 ssha(ji,jj) = ssha(ji,jj) + za1 * ssha_e(ji,jj) 1269 END DO 1270 END DO 1166 1271 ! ! ==================== ! 1167 1272 END DO ! end loop ! … … 1172 1277 ! 1173 1278 ! Set advection velocity correction: 1174 !$OMP PARALLEL WORKSHARE 1175 zwx(:,:) = un_adv(:,:) 1176 zwy(:,:) = vn_adv(:,:) 1177 !$OMP END PARALLEL WORKSHARE 1279 !$OMP PARALLEL DO schedule(static) private(jj,ji) 1280 DO jj = 1, jpj 1281 DO ji = 1, jpi 1282 zwx(ji,jj) = un_adv(ji,jj) 1283 zwy(ji,jj) = vn_adv(ji,jj) 1284 END DO 1285 END DO 1178 1286 IF( ( kt == nit000 .AND. neuler==0 ) .OR. .NOT.ln_bt_fw ) THEN 1179 !$OMP PARALLEL WORKSHARE 1180 un_adv(:,:) = zwx(:,:) * r1_hu_n(:,:) 1181 vn_adv(:,:) = zwy(:,:) * r1_hv_n(:,:) 1182 !$OMP END PARALLEL WORKSHARE 1287 !$OMP PARALLEL DO schedule(static) private(jj,ji) 1288 DO jj = 1, jpj 1289 DO ji = 1, jpi 1290 un_adv(ji,jj) = zwx(ji,jj) * r1_hu_n(ji,jj) 1291 vn_adv(ji,jj) = zwy(ji,jj) * r1_hv_n(ji,jj) 1292 END DO 1293 END DO 1183 1294 ELSE 1184 !$OMP PARALLEL WORKSHARE 1185 un_adv(:,:) = z1_2 * ( ub2_b(:,:) + zwx(:,:) ) * r1_hu_n(:,:) 1186 vn_adv(:,:) = z1_2 * ( vb2_b(:,:) + zwy(:,:) ) * r1_hv_n(:,:) 1187 !$OMP END PARALLEL WORKSHARE 1295 !$OMP PARALLEL DO schedule(static) private(jj,ji) 1296 DO jj = 1, jpj 1297 DO ji = 1, jpi 1298 un_adv(ji,jj) = z1_2 * ( ub2_b(ji,jj) + zwx(ji,jj) ) * r1_hu_n(ji,jj) 1299 vn_adv(ji,jj) = z1_2 * ( vb2_b(ji,jj) + zwy(ji,jj) ) * r1_hv_n(ji,jj) 1300 END DO 1301 END DO 1188 1302 END IF 1189 1303 1190 1304 IF( ln_bt_fw ) THEN ! Save integrated transport for next computation 1191 !$OMP PARALLEL WORKSHARE 1192 ub2_b(:,:) = zwx(:,:) 1193 vb2_b(:,:) = zwy(:,:) 1194 !$OMP END PARALLEL WORKSHARE 1305 !$OMP PARALLEL DO schedule(static) private(jj,ji) 1306 DO jj = 1, jpj 1307 DO ji = 1, jpi 1308 ub2_b(ji,jj) = zwx(ji,jj) 1309 vb2_b(ji,jj) = zwy(ji,jj) 1310 END DO 1311 END DO 1195 1312 ENDIF 1196 1313 ! … … 1225 1342 !$OMP END DO NOWAIT 1226 1343 ! Save barotropic velocities not transport: 1227 !$OMP WORKSHARE 1228 ua_b(:,:) = ua_b(:,:) / ( hu_0(:,:) + zsshu_a(:,:) + 1._wp - ssumask(:,:) ) 1229 va_b(:,:) = va_b(:,:) / ( hv_0(:,:) + zsshv_a(:,:) + 1._wp - ssvmask(:,:) ) 1230 !$OMP END WORKSHARE NOWAIT 1344 !$OMP DO schedule(static) private(jj,ji) 1345 DO jj = 1, jpj 1346 DO ji = 1, jpi 1347 ua_b(ji,jj) = ua_b(ji,jj) / ( hu_0(ji,jj) + zsshu_a(ji,jj) + 1._wp - ssumask(ji,jj) ) 1348 va_b(ji,jj) = va_b(ji,jj) / ( hv_0(ji,jj) + zsshv_a(ji,jj) + 1._wp - ssvmask(ji,jj) ) 1349 END DO 1350 END DO 1231 1351 !$OMP END PARALLEL 1232 1352 ENDIF … … 1249 1369 IF( .NOT.Agrif_Root() .AND. ln_bt_fw ) THEN 1250 1370 IF( Agrif_NbStepint() == 0 ) THEN 1251 !$OMP PARALLEL WORKSHARE 1252 ub2_i_b(:,:) = 0._wp 1253 vb2_i_b(:,:) = 0._wp 1254 !$OMP END PARALLEL WORKSHARE 1371 !$OMP PARALLEL DO schedule(static) private(jj,ji) 1372 DO jj = 1, jpj 1373 DO ji = 1, jpi 1374 ub2_i_b(ji,jj) = 0._wp 1375 vb2_i_b(ji,jj) = 0._wp 1376 END DO 1377 END DO 1255 1378 END IF 1256 1379 ! 1257 1380 za1 = 1._wp / REAL(Agrif_rhot(), wp) 1258 !$OMP PARALLEL WORKSHARE 1259 ub2_i_b(:,:) = ub2_i_b(:,:) + za1 * ub2_b(:,:) 1260 vb2_i_b(:,:) = vb2_i_b(:,:) + za1 * vb2_b(:,:) 1261 !$OMP END PARALLEL WORKSHARE 1381 !$OMP PARALLEL DO schedule(static) private(jj,ji) 1382 DO jj = 1, jpj 1383 DO ji = 1, jpi 1384 ub2_i_b(ji,jj) = ub2_i_b(ji,jj) + za1 * ub2_b(ji,jj) 1385 vb2_i_b(ji,jj) = vb2_i_b(ji,jj) + za1 * vb2_b(ji,jj) 1386 END DO 1387 END DO 1262 1388 ENDIF 1263 1389 #endif -
branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/DYN/dynvor.F90
r7037 r7508 94 94 !!---------------------------------------------------------------------- 95 95 INTEGER, INTENT( in ) :: kt ! ocean time-step index 96 INTEGER :: jk, jj, ji 96 97 ! 97 98 REAL(wp), POINTER, DIMENSION(:,:,:) :: ztrdu, ztrdv … … 106 107 CASE ( np_ENE ) !* energy conserving scheme 107 108 IF( l_trddyn ) THEN ! trend diagnostics: split the trend in two 108 !$OMP PARALLEL WORKSHARE 109 ztrdu(:,:,:) = ua(:,:,:) 110 ztrdv(:,:,:) = va(:,:,:) 111 !$OMP END PARALLEL WORKSHARE 109 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 110 DO jk = 1, jpk 111 DO jj = 1, jpj 112 DO ji = 1, jpi 113 ztrdu(ji,jj,jk) = ua(ji,jj,jk) 114 ztrdv(ji,jj,jk) = va(ji,jj,jk) 115 END DO 116 END DO 117 END DO 112 118 CALL vor_ene( kt, nrvm, ua, va ) ! relative vorticity or metric trend 113 !$OMP PARALLEL WORKSHARE 114 ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:) 115 ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:) 116 !$OMP END PARALLEL WORKSHARE 119 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 120 DO jk = 1, jpk 121 DO jj = 1, jpj 122 DO ji = 1, jpi 123 ztrdu(ji,jj,jk) = ua(ji,jj,jk) - ztrdu(ji,jj,jk) 124 ztrdv(ji,jj,jk) = va(ji,jj,jk) - ztrdv(ji,jj,jk) 125 END DO 126 END DO 127 END DO 117 128 CALL trd_dyn( ztrdu, ztrdv, jpdyn_rvo, kt ) 118 !$OMP PARALLEL WORKSHARE 119 ztrdu(:,:,:) = ua(:,:,:) 120 ztrdv(:,:,:) = va(:,:,:) 121 !$OMP END PARALLEL WORKSHARE 129 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 130 DO jk = 1, jpk 131 DO jj = 1, jpj 132 DO ji = 1, jpi 133 ztrdu(ji,jj,jk) = ua(ji,jj,jk) 134 ztrdv(ji,jj,jk) = va(ji,jj,jk) 135 END DO 136 END DO 137 END DO 122 138 CALL vor_ene( kt, ncor, ua, va ) ! planetary vorticity trend 123 !$OMP PARALLEL WORKSHARE 124 ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:) 125 ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:) 126 !$OMP END PARALLEL WORKSHARE 139 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 140 DO jk = 1, jpk 141 DO jj = 1, jpj 142 DO ji = 1, jpi 143 ztrdu(ji,jj,jk) = ua(ji,jj,jk) - ztrdu(ji,jj,jk) 144 ztrdv(ji,jj,jk) = va(ji,jj,jk) - ztrdv(ji,jj,jk) 145 END DO 146 END DO 147 END DO 127 148 CALL trd_dyn( ztrdu, ztrdv, jpdyn_pvo, kt ) 128 149 ELSE … … 132 153 CASE ( np_ENS ) !* enstrophy conserving scheme 133 154 IF( l_trddyn ) THEN ! trend diagnostics: splitthe trend in two 134 !$OMP PARALLEL WORKSHARE 135 ztrdu(:,:,:) = ua(:,:,:) 136 ztrdv(:,:,:) = va(:,:,:) 137 !$OMP END PARALLEL WORKSHARE 155 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 156 DO jk = 1, jpk 157 DO jj = 1, jpj 158 DO ji = 1, jpi 159 ztrdu(ji,jj,jk) = ua(ji,jj,jk) 160 ztrdv(ji,jj,jk) = va(ji,jj,jk) 161 END DO 162 END DO 163 END DO 138 164 CALL vor_ens( kt, nrvm, ua, va ) ! relative vorticity or metric trend 139 !$OMP PARALLEL WORKSHARE 140 ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:) 141 ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:) 142 !$OMP END PARALLEL WORKSHARE 165 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 166 DO jk = 1, jpk 167 DO jj = 1, jpj 168 DO ji = 1, jpi 169 ztrdu(ji,jj,jk) = ua(ji,jj,jk) - ztrdu(ji,jj,jk) 170 ztrdv(ji,jj,jk) = va(ji,jj,jk) - ztrdv(ji,jj,jk) 171 END DO 172 END DO 173 END DO 143 174 CALL trd_dyn( ztrdu, ztrdv, jpdyn_rvo, kt ) 144 !$OMP PARALLEL WORKSHARE 145 ztrdu(:,:,:) = ua(:,:,:) 146 ztrdv(:,:,:) = va(:,:,:) 147 !$OMP END PARALLEL WORKSHARE 175 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 176 DO jk = 1, jpk 177 DO jj = 1, jpj 178 DO ji = 1, jpi 179 ztrdu(ji,jj,jk) = ua(ji,jj,jk) 180 ztrdv(ji,jj,jk) = va(ji,jj,jk) 181 END DO 182 END DO 183 END DO 148 184 CALL vor_ens( kt, ncor, ua, va ) ! planetary vorticity trend 149 !$OMP PARALLEL WORKSHARE 150 ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:) 151 ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:) 152 !$OMP END PARALLEL WORKSHARE 185 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 186 DO jk = 1, jpk 187 DO jj = 1, jpj 188 DO ji = 1, jpi 189 ztrdu(ji,jj,jk) = ua(ji,jj,jk) - ztrdu(ji,jj,jk) 190 ztrdv(ji,jj,jk) = va(ji,jj,jk) - ztrdv(ji,jj,jk) 191 END DO 192 END DO 193 END DO 153 194 CALL trd_dyn( ztrdu, ztrdv, jpdyn_pvo, kt ) 154 195 ELSE … … 158 199 CASE ( np_MIX ) !* mixed ene-ens scheme 159 200 IF( l_trddyn ) THEN ! trend diagnostics: split the trend in two 160 !$OMP PARALLEL WORKSHARE 161 ztrdu(:,:,:) = ua(:,:,:) 162 ztrdv(:,:,:) = va(:,:,:) 163 !$OMP END PARALLEL WORKSHARE 201 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 202 DO jk = 1, jpk 203 DO jj = 1, jpj 204 DO ji = 1, jpi 205 ztrdu(ji,jj,jk) = ua(ji,jj,jk) 206 ztrdv(ji,jj,jk) = va(ji,jj,jk) 207 END DO 208 END DO 209 END DO 164 210 CALL vor_ens( kt, nrvm, ua, va ) ! relative vorticity or metric trend (ens) 165 !$OMP PARALLEL WORKSHARE 166 ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:) 167 ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:) 168 !$OMP END PARALLEL WORKSHARE 211 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 212 DO jk = 1, jpk 213 DO jj = 1, jpj 214 DO ji = 1, jpi 215 ztrdu(ji,jj,jk) = ua(ji,jj,jk) - ztrdu(ji,jj,jk) 216 ztrdv(ji,jj,jk) = va(ji,jj,jk) - ztrdv(ji,jj,jk) 217 END DO 218 END DO 219 END DO 169 220 CALL trd_dyn( ztrdu, ztrdv, jpdyn_rvo, kt ) 170 !$OMP PARALLEL WORKSHARE 171 ztrdu(:,:,:) = ua(:,:,:) 172 ztrdv(:,:,:) = va(:,:,:) 173 !$OMP END PARALLEL WORKSHARE 221 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 222 DO jk = 1, jpk 223 DO jj = 1, jpj 224 DO ji = 1, jpi 225 ztrdu(ji,jj,jk) = ua(ji,jj,jk) 226 ztrdv(ji,jj,jk) = va(ji,jj,jk) 227 END DO 228 END DO 229 END DO 174 230 CALL vor_ene( kt, ncor, ua, va ) ! planetary vorticity trend (ene) 175 !$OMP PARALLEL WORKSHARE 176 ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:) 177 ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:) 178 !$OMP END PARALLEL WORKSHARE 231 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 232 DO jk = 1, jpk 233 DO jj = 1, jpj 234 DO ji = 1, jpi 235 ztrdu(ji,jj,jk) = ua(ji,jj,jk) - ztrdu(ji,jj,jk) 236 ztrdv(ji,jj,jk) = va(ji,jj,jk) - ztrdv(ji,jj,jk) 237 END DO 238 END DO 239 END DO 179 240 CALL trd_dyn( ztrdu, ztrdv, jpdyn_pvo, kt ) 180 241 ELSE … … 185 246 CASE ( np_EEN ) !* energy and enstrophy conserving scheme 186 247 IF( l_trddyn ) THEN ! trend diagnostics: split the trend in two 187 !$OMP PARALLEL WORKSHARE 188 ztrdu(:,:,:) = ua(:,:,:) 189 ztrdv(:,:,:) = va(:,:,:) 190 !$OMP END PARALLEL WORKSHARE 248 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 249 DO jk = 1, jpk 250 DO jj = 1, jpj 251 DO ji = 1, jpi 252 ztrdu(ji,jj,jk) = ua(ji,jj,jk) 253 ztrdv(ji,jj,jk) = va(ji,jj,jk) 254 END DO 255 END DO 256 END DO 191 257 CALL vor_een( kt, nrvm, ua, va ) ! relative vorticity or metric trend 192 !$OMP PARALLEL WORKSHARE 193 ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:) 194 ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:) 195 !$OMP END PARALLEL WORKSHARE 258 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 259 DO jk = 1, jpk 260 DO jj = 1, jpj 261 DO ji = 1, jpi 262 ztrdu(ji,jj,jk) = ua(ji,jj,jk) - ztrdu(ji,jj,jk) 263 ztrdv(ji,jj,jk) = va(ji,jj,jk) - ztrdv(ji,jj,jk) 264 END DO 265 END DO 266 END DO 196 267 CALL trd_dyn( ztrdu, ztrdv, jpdyn_rvo, kt ) 197 !$OMP PARALLEL WORKSHARE 198 ztrdu(:,:,:) = ua(:,:,:) 199 ztrdv(:,:,:) = va(:,:,:) 200 !$OMP END PARALLEL WORKSHARE 268 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 269 DO jk = 1, jpk 270 DO jj = 1, jpj 271 DO ji = 1, jpi 272 ztrdu(ji,jj,jk) = ua(ji,jj,jk) 273 ztrdv(ji,jj,jk) = va(ji,jj,jk) 274 END DO 275 END DO 276 END DO 201 277 CALL vor_een( kt, ncor, ua, va ) ! planetary vorticity trend 202 !$OMP PARALLEL WORKSHARE 203 ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:) 204 ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:) 205 !$OMP END PARALLEL WORKSHARE 278 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 279 DO jk = 1, jpk 280 DO jj = 1, jpj 281 DO ji = 1, jpi 282 ztrdu(ji,jj,jk) = ua(ji,jj,jk) - ztrdu(ji,jj,jk) 283 ztrdv(ji,jj,jk) = va(ji,jj,jk) - ztrdv(ji,jj,jk) 284 END DO 285 END DO 286 END DO 206 287 CALL trd_dyn( ztrdu, ztrdv, jpdyn_pvo, kt ) 207 288 ELSE … … 269 350 SELECT CASE( kvor ) !== vorticity considered ==! 270 351 CASE ( np_COR ) !* Coriolis (planetary vorticity) 271 !$OMP PARALLEL WORKSHARE 272 zwz(:,:) = ff(:,:) 273 !$OMP END PARALLEL WORKSHARE 352 !$OMP PARALLEL DO schedule(static) private(jj,ji) 353 DO jj = 1, jpj 354 DO ji = 1, jpi 355 zwz(ji,jj) = ff(ji,jj) 356 END DO 357 END DO 274 358 CASE ( np_RVO ) !* relative vorticity 275 !$OMP PARALLEL DO private(jj,ji)359 !$OMP PARALLEL DO schedule(static) private(jj,ji) 276 360 DO jj = 1, jpjm1 277 361 DO ji = 1, fs_jpim1 ! vector opt. … … 281 365 END DO 282 366 CASE ( np_MET ) !* metric term 283 !$OMP PARALLEL DO private(jj,ji)367 !$OMP PARALLEL DO schedule(static) private(jj,ji) 284 368 DO jj = 1, jpjm1 285 369 DO ji = 1, fs_jpim1 ! vector opt. … … 290 374 END DO 291 375 CASE ( np_CRV ) !* Coriolis + relative vorticity 292 !$OMP PARALLEL DO private(jj,ji)376 !$OMP PARALLEL DO schedule(static) private(jj,ji) 293 377 DO jj = 1, jpjm1 294 378 DO ji = 1, fs_jpim1 ! vector opt. … … 299 383 END DO 300 384 CASE ( np_CME ) !* Coriolis + metric 301 !$OMP PARALLEL DO private(jj,ji)385 !$OMP PARALLEL DO schedule(static) private(jj,ji) 302 386 DO jj = 1, jpjm1 303 387 DO ji = 1, fs_jpim1 ! vector opt. … … 313 397 ! 314 398 IF( ln_dynvor_msk ) THEN !== mask/unmask vorticity ==! 315 !$OMP PARALLEL DO private(jj,ji)399 !$OMP PARALLEL DO schedule(static) private(jj,ji) 316 400 DO jj = 1, jpjm1 317 401 DO ji = 1, fs_jpim1 ! vector opt. … … 322 406 323 407 IF( ln_sco ) THEN 324 !$OMP PARALLEL WORKSHARE 325 zwz(:,:) = zwz(:,:) / e3f_n(:,:,jk) 326 zwx(:,:) = e2u(:,:) * e3u_n(:,:,jk) * un(:,:,jk) 327 zwy(:,:) = e1v(:,:) * e3v_n(:,:,jk) * vn(:,:,jk) 328 !$OMP END PARALLEL WORKSHARE 408 !$OMP PARALLEL DO schedule(static) private(jj,ji) 409 DO jj = 1, jpj 410 DO ji = 1, jpi 411 zwz(ji,jj) = zwz(ji,jj) / e3f_n(ji,jj,jk) 412 zwx(ji,jj) = e2u(ji,jj) * e3u_n(ji,jj,jk) * un(ji,jj,jk) 413 zwy(ji,jj) = e1v(ji,jj) * e3v_n(ji,jj,jk) * vn(ji,jj,jk) 414 END DO 415 END DO 329 416 ELSE 330 !$OMP PARALLEL WORKSHARE 331 zwx(:,:) = e2u(:,:) * un(:,:,jk) 332 zwy(:,:) = e1v(:,:) * vn(:,:,jk) 333 !$OMP END PARALLEL WORKSHARE 417 !$OMP PARALLEL DO schedule(static) private(jj,ji) 418 DO jj = 1, jpj 419 DO ji = 1, jpi 420 zwx(ji,jj) = e2u(ji,jj) * un(ji,jj,jk) 421 zwy(ji,jj) = e1v(ji,jj) * vn(ji,jj,jk) 422 END DO 423 END DO 334 424 ENDIF 335 425 ! !== compute and add the vorticity term trend =! 336 !$OMP PARALLEL DO private(jj, ji, zy1, zy2, zx1, zx2)426 !$OMP PARALLEL DO schedule(static) private(jj, ji, zy1, zy2, zx1, zx2) 337 427 DO jj = 2, jpjm1 338 428 DO ji = fs_2, fs_jpim1 ! vector opt. … … 523 613 SELECT CASE( nn_een_e3f ) ! == reciprocal of e3 at F-point 524 614 CASE ( 0 ) ! original formulation (masked averaging of e3t divided by 4) 525 !$OMP PARALLEL DO private(jj,ji,ze3)615 !$OMP PARALLEL DO schedule(static) private(jj,ji,ze3) 526 616 DO jj = 1, jpjm1 527 617 DO ji = 1, fs_jpim1 ! vector opt. … … 534 624 END DO 535 625 CASE ( 1 ) ! new formulation (masked averaging of e3t divided by the sum of mask) 536 !$OMP PARALLEL DO private(jj,ji,ze3,zmsk)626 !$OMP PARALLEL DO schedule(static) private(jj,ji,ze3,zmsk) 537 627 DO jj = 1, jpjm1 538 628 DO ji = 1, fs_jpim1 ! vector opt. … … 550 640 SELECT CASE( kvor ) !== vorticity considered ==! 551 641 CASE ( np_COR ) !* Coriolis (planetary vorticity) 552 !$OMP PARALLEL DO private(jj,ji)642 !$OMP PARALLEL DO schedule(static) private(jj,ji) 553 643 DO jj = 1, jpjm1 554 644 DO ji = 1, fs_jpim1 ! vector opt. … … 557 647 END DO 558 648 CASE ( np_RVO ) !* relative vorticity 559 !$OMP PARALLEL DO private(jj,ji)649 !$OMP PARALLEL DO schedule(static) private(jj,ji) 560 650 DO jj = 1, jpjm1 561 651 DO ji = 1, fs_jpim1 ! vector opt. … … 566 656 END DO 567 657 CASE ( np_MET ) !* metric term 568 !$OMP PARALLEL DO private(jj,ji)658 !$OMP PARALLEL DO schedule(static) private(jj,ji) 569 659 DO jj = 1, jpjm1 570 660 DO ji = 1, fs_jpim1 ! vector opt. … … 575 665 END DO 576 666 CASE ( np_CRV ) !* Coriolis + relative vorticity 577 !$OMP PARALLEL DO private(jj,ji)667 !$OMP PARALLEL DO schedule(static) private(jj,ji) 578 668 DO jj = 1, jpjm1 579 669 DO ji = 1, fs_jpim1 ! vector opt. … … 584 674 END DO 585 675 CASE ( np_CME ) !* Coriolis + metric 586 !$OMP PARALLEL DO private(jj,ji)676 !$OMP PARALLEL DO schedule(static) private(jj,ji) 587 677 DO jj = 1, jpjm1 588 678 DO ji = 1, fs_jpim1 ! vector opt. … … 598 688 ! 599 689 IF( ln_dynvor_msk ) THEN !== mask/unmask vorticity ==! 600 !$OMP PARALLEL DO private(jj,ji)690 !$OMP PARALLEL DO schedule(static) private(jj,ji) 601 691 DO jj = 1, jpjm1 602 692 DO ji = 1, fs_jpim1 ! vector opt. … … 609 699 ! 610 700 ! !== horizontal fluxes ==! 611 !$OMP PARALLEL DO private(jj,ji) 701 !$OMP PARALLEL 702 !$OMP DO schedule(static) private(jj,ji) 612 703 DO jj = 1, jpj 613 704 DO ji = 1, jpi … … 618 709 619 710 ! !== compute and add the vorticity term trend =! 711 !$OMP DO schedule(static) private(jj) 712 DO jj = 1, jpj 713 ztne(1,jj) = 0 ; ztnw(1,jj) = 0 ; ztse(1,jj) = 0 ; ztsw(1,jj) = 0 714 END DO 620 715 jj = 2 621 ztne(1,:) = 0 ; ztnw(1,:) = 0 ; ztse(1,:) = 0 ; ztsw(1,:) = 0 716 !$OMP DO schedule(static) private(ji) 622 717 DO ji = 2, jpi ! split in 2 parts due to vector opt. 623 718 ztne(ji,jj) = zwz(ji-1,jj ) + zwz(ji ,jj ) + zwz(ji ,jj-1) … … 626 721 ztsw(ji,jj) = zwz(ji ,jj-1) + zwz(ji-1,jj-1) + zwz(ji-1,jj ) 627 722 END DO 628 !$OMP PARALLEL 629 !$OMP DO private(jj,ji) 723 !$OMP DO schedule(static) private(jj,ji) 630 724 DO jj = 3, jpj 631 725 DO ji = fs_2, jpi ! vector opt. ok because we start at jj = 3 … … 636 730 END DO 637 731 END DO 638 !$OMP DO private(jj,ji,zua,zva)732 !$OMP DO schedule(static) private(jj,ji,zua,zva) 639 733 DO jj = 2, jpjm1 640 734 DO ji = fs_2, fs_jpim1 ! vector opt. … … 703 797 IF(lwp) WRITE(numout,*) ' namlbc: change fmask value in the angles (T) ln_vorlat = ', ln_vorlat 704 798 IF( ln_vorlat .AND. ( ln_dynvor_ene .OR. ln_dynvor_ens .OR. ln_dynvor_mix ) ) THEN 705 !$OMP PARALLEL DO private(jk,jj,ji)799 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 706 800 DO jk = 1, jpk 707 801 DO jj = 2, jpjm1 -
branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/DYN/dynzad.F90
r7037 r7508 77 77 IF( l_trddyn ) THEN ! Save ua and va trends 78 78 CALL wrk_alloc( jpi, jpj, jpk, ztrdu, ztrdv ) 79 !$OMP PARALLEL WORKSHARE 80 ztrdu(:,:,:) = ua(:,:,:) 81 ztrdv(:,:,:) = va(:,:,:) 82 !$OMP END PARALLEL WORKSHARE 79 !$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 80 DO jk = 1, jpk 81 DO jj = 1, jpj 82 DO ji = 1, jpi 83 ztrdu(ji,jj,jk) = ua(ji,jj,jk) 84 ztrdv(ji,jj,jk) = va(ji,jj,jk) 85 END DO 86 END DO 87 END DO 83 88 ENDIF 84 89 … … 139 144 140 145 IF( l_trddyn ) THEN ! save the vertical advection trends for diagnostic 141 !$OMP PARALLEL WORKSHARE 142 ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:) 143 ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:) 144 !$OMP END PARALLEL WORKSHARE 146 !$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 147 DO jk = 1, jpk 148 DO jj = 1, jpj 149 DO ji = 1, jpi 150 ztrdu(ji,jj,jk) = ua(ji,jj,jk) - ztrdu(ji,jj,jk) 151 ztrdv(ji,jj,jk) = va(ji,jj,jk) - ztrdv(ji,jj,jk) 152 END DO 153 END DO 154 END DO 145 155 CALL trd_dyn( ztrdu, ztrdv, jpdyn_zad, kt ) 146 156 CALL wrk_dealloc( jpi, jpj, jpk, ztrdu, ztrdv ) -
branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/DYN/dynzdf.F90
r6748 r7508 53 53 !! 54 54 INTEGER, INTENT( in ) :: kt ! ocean time-step index 55 INTEGER :: ji, jj, jk ! dummy loop indices 55 56 ! 56 57 REAL(wp), POINTER, DIMENSION(:,:,:) :: ztrdu, ztrdv … … 66 67 IF( l_trddyn ) THEN ! temporary save of ta and sa trends 67 68 CALL wrk_alloc( jpi, jpj, jpk, ztrdu, ztrdv ) 68 !$OMP PARALLEL WORKSHARE 69 ztrdu(:,:,:) = ua(:,:,:) 70 ztrdv(:,:,:) = va(:,:,:) 71 !$OMP END PARALLEL WORKSHARE 69 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 70 DO jk = 1, jpk 71 DO jj = 1, jpj 72 DO ji = 1, jpi 73 ztrdu(ji,jj,jk) = ua(ji,jj,jk) 74 ztrdv(ji,jj,jk) = va(ji,jj,jk) 75 END DO 76 END DO 77 END DO 72 78 ENDIF 73 79 … … 80 86 81 87 IF( l_trddyn ) THEN ! save the vertical diffusive trends for further diagnostics 82 !$OMP PARALLEL WORKSHARE 83 ztrdu(:,:,:) = ( ua(:,:,:) - ub(:,:,:) ) / r2dt - ztrdu(:,:,:) 84 ztrdv(:,:,:) = ( va(:,:,:) - vb(:,:,:) ) / r2dt - ztrdv(:,:,:) 85 !$OMP END PARALLEL WORKSHARE 88 !$OMP PARALLEL DO schedule(static) private(jk,jj,ji) 89 DO jk = 1, jpk 90 DO jj = 1, jpj 91 DO ji = 1, jpi 92 ztrdu(ji,jj,jk) = ( ua(ji,jj,jk) - ub(ji,jj,jk) ) / r2dt - ztrdu(ji,jj,jk) 93 ztrdv(ji,jj,jk) = ( va(ji,jj,jk) - vb(ji,jj,jk) ) / r2dt - ztrdv(ji,jj,jk) 94 END DO 95 END DO 96 END DO 86 97 CALL trd_dyn( ztrdu, ztrdv, jpdyn_zdf, kt ) 87 98 CALL wrk_dealloc( jpi, jpj, jpk, ztrdu, ztrdv ) -
branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90
r7037 r7508 97 97 ! !------------------------------! 98 98 !$OMP PARALLEL 99 !$OMP WORKSHARE 100 zhdiv(:,:) = 0._wp 101 !$OMP END WORKSHARE 99 !$OMP DO schedule(static) private(jj, ji) 100 DO jj = 1, jpj 101 DO ji = 1, jpi 102 zhdiv(ji,jj) = 0._wp 103 END DO 104 END DO 102 105 DO jk = 1, jpkm1 ! Horizontal divergence of barotropic transports 103 106 !$OMP DO schedule(static) private(jj, ji) 104 107 DO jj = 1, jpj 105 108 DO ji = 1, jpi ! vector opt. 106 zhdiv(ji,jj) = zhdiv(ji,jj) + e3t_n(ji,jj,jk) * hdivn(ji,jj,jk)109 zhdiv(ji,jj) = zhdiv(ji,jj) + e3t_n(ji,jj,jk) * hdivn(ji,jj,jk) 107 110 END DO 108 111 END DO … … 116 119 117 120 IF(ln_wd) CALL wad_lmt(sshb, zcoef * (emp_b(:,:) + emp(:,:)), z2dt) 118 !$OMP PARALLEL WORKSHARE 119 ssha(:,:) = ( sshb(:,:) - z2dt * ( zcoef * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) ) ) * ssmask(:,:) 120 !$OMP END PARALLEL WORKSHARE 121 !$OMP PARALLEL DO schedule(static) private(jj, ji) 122 DO jj = 1, jpj 123 DO ji = 1, jpi 124 ssha(ji,jj) = ( sshb(ji,jj) - z2dt * ( zcoef * ( emp_b(ji,jj) + emp(ji,jj) ) + zhdiv(ji,jj) ) ) * ssmask(ji,jj) 125 END DO 126 END DO 121 127 IF ( .NOT.ln_dynspg_ts ) THEN 122 128 ! These lines are not necessary with time splitting since … … 136 142 IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN ! Include the IAU weighted SSH increment 137 143 CALL ssh_asm_inc( kt ) 138 !$OMP PARALLEL WORKSHARE 139 ssha(:,:) = ssha(:,:) + z2dt * ssh_iau(:,:) 140 !$OMP END PARALLEL WORKSHARE 144 !$OMP PARALLEL DO schedule(static) private(jj, ji) 145 DO jj = 1, jpj 146 DO ji = 1, jpi 147 ssha(ji,jj) = ssha(ji,jj) + z2dt * ssh_iau(ji,jj) 148 END DO 149 END DO 141 150 ENDIF 142 151 #endif … … 184 193 IF(lwp) WRITE(numout,*) '~~~~~ ' 185 194 ! 186 wn( :,:,jpk) = 0._wp ! bottom boundary condition: w=0 (set once for all)195 wn(ji,jj,jpk) = 0._wp ! bottom boundary condition: w=0 (set once for all) 187 196 ENDIF 188 197 ! !------------------------------! … … 267 276 INTEGER, INTENT(in) :: kt ! ocean time-step index 268 277 ! 278 INTEGER :: ji, jj, jk ! dummy loop indices 269 279 REAL(wp) :: zcoef ! local scalar 270 280 !!---------------------------------------------------------------------- … … 280 290 IF( ( neuler == 0 .AND. kt == nit000 ) .OR. & 281 291 & ( ln_bt_fw .AND. ln_dynspg_ts ) ) THEN 282 !$OMP PARALLEL WORKSHARE 283 sshb(:,:) = sshn(:,:) ! before <-- now 284 sshn(:,:) = ssha(:,:) ! now <-- after (before already = now) 285 !$OMP END PARALLEL WORKSHARE 292 !$OMP PARALLEL DO schedule(static) private(jj, ji) 293 DO jj = 1, jpj 294 DO ji = 1, jpi 295 sshb(ji,jj) = sshn(ji,jj) ! before <-- now 296 sshn(ji,jj) = ssha(ji,jj) ! now <-- after (before already = now) 297 END DO 298 END DO 286 299 ! 287 300 ELSE !== Leap-Frog time-stepping: Asselin filter + swap ==! 288 301 ! ! before <-- now filtered 289 !$OMP PARALLEL WORKSHARE 290 sshb(:,:) = sshn(:,:) + atfp * ( sshb(:,:) - 2 * sshn(:,:) + ssha(:,:) ) 291 !$OMP END PARALLEL WORKSHARE 302 !$OMP PARALLEL DO schedule(static) private(jj, ji) 303 DO jj = 1, jpj 304 DO ji = 1, jpi 305 sshb(ji,jj) = sshn(ji,jj) + atfp * ( sshb(ji,jj) - 2 * sshn(ji,jj) + ssha(ji,jj) ) 306 END DO 307 END DO 292 308 IF( .NOT.ln_linssh ) THEN ! before <-- with forcing removed 293 309 zcoef = atfp * rdt * r1_rau0 294 !$OMP PARALLEL WORKSHARE 295 sshb(:,:) = sshb(:,:) - zcoef * ( emp_b(:,:) - emp (:,:) & 296 & - rnf_b(:,:) + rnf (:,:) & 297 & + fwfisf_b(:,:) - fwfisf(:,:) ) * ssmask(:,:) 298 !$OMP END PARALLEL WORKSHARE 310 !$OMP PARALLEL DO schedule(static) private(jj, ji) 311 DO jj = 1, jpj 312 DO ji = 1, jpi 313 sshb(ji,jj) = sshb(ji,jj) - zcoef * ( emp_b(ji,jj) - emp (ji,jj) & 314 & - rnf_b(ji,jj) + rnf (ji,jj) & 315 & + fwfisf_b(ji,jj) - fwfisf(ji,jj) ) * ssmask(ji,jj) 316 END DO 317 END DO 299 318 ENDIF 300 !$OMP PARALLEL WORKSHARE 301 sshn(:,:) = ssha(:,:) ! now <-- after 302 !$OMP END PARALLEL WORKSHARE 319 !$OMP PARALLEL DO schedule(static) private(jj, ji) 320 DO jj = 1, jpj 321 DO ji = 1, jpi 322 sshn(ji,jj) = ssha(ji,jj) ! now <-- after 323 END DO 324 END DO 303 325 ENDIF 304 326 !
Note: See TracChangeset
for help on using the changeset viewer.