Changeset 11071
- Timestamp:
- 2019-06-04T14:58:06+02:00 (6 years ago)
- Location:
- NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src
- Files:
-
- 9 edited
Legend:
- Unmodified
- Added
- Removed
-
NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/BDY/bdy_oce.F90
r11067 r11071 139 139 TYPE(OBC_DATA) , DIMENSION(jp_bdy), TARGET :: dta_bdy !: bdy external data (local process) 140 140 !$AGRIF_END_DO_NOT_TREAT 141 LOGICAL, ALLOCATABLE, DIMENSION(:,:,:) :: lsend_bdy !: mark needed communication for given boundary, grid and direction 142 LOGICAL, ALLOCATABLE, DIMENSION(:,:,:) :: lrecv_bdy !: mark needed communication for given boundary, grid and direction 141 LOGICAL, ALLOCATABLE, DIMENSION(:,:,:) :: lsend_bdy !: mark needed communication for given boundary, grid and neighbour 142 LOGICAL, ALLOCATABLE, DIMENSION(:,:,:) :: lrecv_bdy !: when searching in any direction 143 LOGICAL, ALLOCATABLE, DIMENSION(:,:,:) :: lsend_bdyint !: mark needed communication for given boundary, grid and neighbour 144 LOGICAL, ALLOCATABLE, DIMENSION(:,:,:) :: lrecv_bdyint !: when searching towards the interior of the computational domain 145 LOGICAL, ALLOCATABLE, DIMENSION(:,:,:) :: lsend_bdyext !: mark needed communication for given boundary, grid and neighbour 146 LOGICAL, ALLOCATABLE, DIMENSION(:,:,:) :: lrecv_bdyext !: when searching towards the exterior of the computational domain 143 147 !!---------------------------------------------------------------------- 144 148 !! NEMO/OCE 4.0 , NEMO Consortium (2018) -
NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/BDY/bdydyn2d.F90
r11067 r11071 51 51 !! 52 52 INTEGER :: ib_bdy ! Loop counter 53 LOGICAL, DIMENSION(4) :: l send2, lrecv2, lsend3,lrecv3 ! indicate how communications are to be carried out53 LOGICAL, DIMENSION(4) :: llsend2, llrecv2, llsend3, llrecv3 ! indicate how communications are to be carried out 54 54 55 55 DO ib_bdy=1, nb_bdy … … 72 72 ENDDO 73 73 ! 74 l send2(:) = .false.75 l recv2(:) = .false.76 l send3(:) = .false.77 l recv3(:) = .false.74 llsend2(:) = .false. 75 llrecv2(:) = .false. 76 llsend3(:) = .false. 77 llrecv3(:) = .false. 78 78 DO ib_bdy=1, nb_bdy 79 79 SELECT CASE( cn_dyn2d(ib_bdy) ) 80 80 CASE('flather') 81 lsend2(:) = lsend2(:) .OR. lsend_bdy(ib_bdy,2,:) ! to every bdy neighbour, U points 82 lrecv2(:) = lrecv2(:) .OR. lrecv_bdy(ib_bdy,2,:) ! from every bdy neighbour, U points 83 lsend3(:) = lsend3(:) .OR. lsend_bdy(ib_bdy,3,:) ! to every bdy neighbour, V points 84 lrecv3(:) = lrecv3(:) .OR. lrecv_bdy(ib_bdy,3,:) ! from every bdy neighbour, V points 85 CASE('orlanski') 86 lsend2(:) = lsend2(:) .OR. lsend_bdy(ib_bdy,2,:) ! to every bdy neighbour, U points 87 lrecv2(:) = lrecv2(:) .OR. lrecv_bdy(ib_bdy,2,:) ! from every bdy neighbour, U points 88 lsend3(:) = lsend3(:) .OR. lsend_bdy(ib_bdy,3,:) ! to every bdy neighbour, V points 89 lrecv3(:) = lrecv3(:) .OR. lrecv_bdy(ib_bdy,3,:) ! from every bdy neighbour, V points 90 CASE('orlanski_npo') 91 lsend2(:) = lsend2(:) .OR. lsend_bdy(ib_bdy,2,:) ! to every bdy neighbour, U points 92 lrecv2(:) = lrecv2(:) .OR. lrecv_bdy(ib_bdy,2,:) ! from every bdy neighbour, U points 93 lsend3(:) = lsend3(:) .OR. lsend_bdy(ib_bdy,3,:) ! to every bdy neighbour, V points 94 lrecv3(:) = lrecv3(:) .OR. lrecv_bdy(ib_bdy,3,:) ! from every bdy neighbour, V points 81 llsend2(1:2) = llsend2(1:2) .OR. lsend_bdyint(ib_bdy,2,1:2) ! west/east, U points 82 llsend2(1) = llsend2(1) .OR. lsend_bdyext(ib_bdy,2,1) ! neighbour might search point towards bdy on its east 83 llrecv2(1:2) = llrecv2(1:2) .OR. lrecv_bdyint(ib_bdy,2,1:2) ! west/east, U points 84 llrecv2(2) = llrecv2(2) .OR. lrecv_bdyext(ib_bdy,2,2) ! might search point towards bdy on the east 85 llsend3(3:4) = llsend3(3:4) .OR. lsend_bdyint(ib_bdy,3,3:4) ! north/south, V points 86 llsend3(3) = llsend3(3) .OR. lsend_bdyext(ib_bdy,3,3) ! neighbour might search point towards bdy on its north 87 llrecv3(3:4) = llrecv3(3:4) .OR. lrecv_bdyint(ib_bdy,3,3:4) ! north/south, V points 88 llrecv3(4) = llrecv3(4) .OR. lrecv_bdyext(ib_bdy,3,4) ! might search point towards bdy on the north 89 CASE('orlanski', 'orlanski_npo') 90 llsend2(:) = llsend2(:) .OR. lsend_bdy(ib_bdy,2,:) ! possibly every direction, U points 91 llrecv2(:) = llrecv2(:) .OR. lrecv_bdy(ib_bdy,2,:) ! possibly every direction, U points 92 llsend3(:) = llsend3(:) .OR. lsend_bdy(ib_bdy,3,:) ! possibly every direction, V points 93 llrecv3(:) = llrecv3(:) .OR. lrecv_bdy(ib_bdy,3,:) ! possibly every direction, V points 95 94 END SELECT 96 95 END DO 97 IF( ANY(l send2) .OR. ANY(lrecv2) ) THEN ! if need to send/recv in at least one direction98 CALL lbc_bdy_lnk( 'bdydyn2d', l send2,lrecv2, pua2d, 'U', -1. )96 IF( ANY(llsend2) .OR. ANY(llrecv2) ) THEN ! if need to send/recv in at least one direction 97 CALL lbc_bdy_lnk( 'bdydyn2d', llsend2, llrecv2, pua2d, 'U', -1. ) 99 98 END IF 100 IF( ANY(l send3) .OR. ANY(lrecv3) ) THEN ! if need to send/recv in at least one direction101 CALL lbc_bdy_lnk( 'bdydyn2d', l send3,lrecv3, pva2d, 'V', -1. )99 IF( ANY(llsend3) .OR. ANY(llrecv3) ) THEN ! if need to send/recv in at least one direction 100 CALL lbc_bdy_lnk( 'bdydyn2d', llsend3, llrecv3, pva2d, 'V', -1. ) 102 101 END IF 103 102 ! … … 197 196 ENDIF 198 197 END DO 199 200 198 ! 201 199 igrd = 2 ! Flather bc on u-velocity; … … 284 282 !! 285 283 INTEGER :: ib_bdy ! bdy index 286 LOGICAL, DIMENSION(4) :: l send1,lrecv1 ! indicate how communications are to be carried out287 !!---------------------------------------------------------------------- 288 l send1(:) = .false.289 l recv1(:) = .false.284 LOGICAL, DIMENSION(4) :: llsend1, llrecv1 ! indicate how communications are to be carried out 285 !!---------------------------------------------------------------------- 286 llsend1(:) = .false. 287 llrecv1(:) = .false. 290 288 DO ib_bdy = 1, nb_bdy 291 289 CALL bdy_nmn( idx_bdy(ib_bdy), 1, zssh ) ! zssh is masked 292 l send1(:) = lsend1(:) .OR. lsend_bdy(ib_bdy,1,:) ! to every bdy neighbour, T points293 l recv1(:) = lrecv1(:) .OR. lrecv_bdy(ib_bdy,1,:) ! from every bdy neighbour, T points294 END DO 295 IF( ANY(l send1) .OR. ANY(lrecv1) ) THEN ! if need to send/recv in at least one direction296 CALL lbc_bdy_lnk( 'bdydyn2d', l send1,lrecv1, zssh(:,:,1), 'T', 1. )290 llsend1(:) = llsend1(:) .OR. lsend_bdyint(ib_bdy,1,:) ! possibly every direction, T points 291 llrecv1(:) = llrecv1(:) .OR. lrecv_bdyint(ib_bdy,1,:) ! possibly every direction, T points 292 END DO 293 IF( ANY(llsend1) .OR. ANY(llrecv1) ) THEN ! if need to send/recv in at least one direction 294 CALL lbc_bdy_lnk( 'bdydyn2d', llsend1, llrecv1, zssh(:,:,1), 'T', 1. ) 297 295 END IF 298 296 ! -
NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/BDY/bdydyn3d.F90
r11067 r11071 43 43 ! 44 44 INTEGER :: ib_bdy ! loop index 45 LOGICAL, DIMENSION(4) :: l send2, lrecv2, lsend3,lrecv3 ! indicate how communications are to be carried out45 LOGICAL, DIMENSION(4) :: llsend2, llrecv2, llsend3, llrecv3 ! indicate how communications are to be carried out 46 46 47 47 !!---------------------------------------------------------------------- … … 61 61 END DO 62 62 ! 63 l send2(:) = .false.64 l recv2(:) = .false.65 l send3(:) = .false.66 l recv3(:) = .false.63 llsend2(:) = .false. 64 llrecv2(:) = .false. 65 llsend3(:) = .false. 66 llrecv3(:) = .false. 67 67 DO ib_bdy=1, nb_bdy 68 68 SELECT CASE( cn_dyn3d(ib_bdy) ) 69 CASE('orlanski') 70 lsend2(:) = lsend2(:) .OR. lsend_bdy(ib_bdy,2,:) ! to every bdy neighbour, U points 71 lrecv2(:) = lrecv2(:) .OR. lrecv_bdy(ib_bdy,2,:) ! from every bdy neighbour, U points 72 lsend3(:) = lsend3(:) .OR. lsend_bdy(ib_bdy,3,:) ! to every bdy neighbour, V points 73 lrecv3(:) = lrecv3(:) .OR. lrecv_bdy(ib_bdy,3,:) ! from every bdy neighbour, V points 74 CASE('orlanski_npo') 75 lsend2(:) = lsend2(:) .OR. lsend_bdy(ib_bdy,2,:) ! to every bdy neighbour, U points 76 lrecv2(:) = lrecv2(:) .OR. lrecv_bdy(ib_bdy,2,:) ! from every bdy neighbour, U points 77 lsend3(:) = lsend3(:) .OR. lsend_bdy(ib_bdy,3,:) ! to every bdy neighbour, V points 78 lrecv3(:) = lrecv3(:) .OR. lrecv_bdy(ib_bdy,3,:) ! from every bdy neighbour, V points 69 CASE('orlanski', 'orlanski_npo') 70 llsend2(:) = llsend2(:) .OR. lsend_bdy(ib_bdy,2,:) ! possibly every direction, U points 71 llrecv2(:) = llrecv2(:) .OR. lrecv_bdy(ib_bdy,2,:) ! possibly every direction, U points 72 llsend3(:) = llsend3(:) .OR. lsend_bdy(ib_bdy,3,:) ! possibly every direction, V points 73 llrecv3(:) = llrecv3(:) .OR. lrecv_bdy(ib_bdy,3,:) ! possibly every direction, V points 79 74 CASE('zerograd') 80 l send2(:) = lsend2(:) .OR. lsend_bdy(ib_bdy,2,:) ! to every bdy neighbour, U points81 l recv2(:) = lrecv2(:) .OR. lrecv_bdy(ib_bdy,2,:) ! from every bdy neighbour, U points82 l send3(:) = lsend3(:) .OR. lsend_bdy(ib_bdy,3,:) ! to every bdy neighbour, V points83 l recv3(:) = lrecv3(:) .OR. lrecv_bdy(ib_bdy,3,:) ! from every bdy neighbour, V points75 llsend2(3:4) = llsend2(3:4) .OR. lsend_bdyint(ib_bdy,2,3:4) ! north/south, U points 76 llrecv2(3:4) = llrecv2(3:4) .OR. lrecv_bdyint(ib_bdy,2,3:4) ! north/south, U points 77 llsend3(1:2) = llsend3(1:2) .OR. lsend_bdyint(ib_bdy,3,1:2) ! west/east, V points 78 llrecv3(1:2) = llrecv3(1:2) .OR. lrecv_bdyint(ib_bdy,3,1:2) ! west/east, V points 84 79 CASE('neumann') 85 l send2(:) = lsend2(:) .OR. lsend_bdy(ib_bdy,2,:) ! to every bdy neighbour, U points86 l recv2(:) = lrecv2(:) .OR. lrecv_bdy(ib_bdy,2,:) ! from every bdy neighbour, U points87 l send3(:) = lsend3(:) .OR. lsend_bdy(ib_bdy,3,:) ! to every bdy neighbour, V points88 l recv3(:) = lrecv3(:) .OR. lrecv_bdy(ib_bdy,3,:) ! from every bdy neighbour, V points80 llsend2(:) = llsend2(:) .OR. lsend_bdyint(ib_bdy,2,:) ! possibly every direction, U points 81 llrecv2(:) = llrecv2(:) .OR. lrecv_bdyint(ib_bdy,2,:) ! possibly every direction, U points 82 llsend3(:) = llsend3(:) .OR. lsend_bdyint(ib_bdy,3,:) ! possibly every direction, V points 83 llrecv3(:) = llrecv3(:) .OR. lrecv_bdyint(ib_bdy,3,:) ! possibly every direction, V points 89 84 END SELECT 90 85 END DO 91 86 ! 92 IF( ANY(l send2) .OR. ANY(lrecv2) ) THEN ! if need to send/recv in at least one direction93 CALL lbc_bdy_lnk( 'bdydyn2d', l send2,lrecv2, ua, 'U', -1. )87 IF( ANY(llsend2) .OR. ANY(llrecv2) ) THEN ! if need to send/recv in at least one direction 88 CALL lbc_bdy_lnk( 'bdydyn2d', llsend2, llrecv2, ua, 'U', -1. ) 94 89 END IF 95 IF( ANY(l send3) .OR. ANY(lrecv3) ) THEN ! if need to send/recv in at least one direction96 CALL lbc_bdy_lnk( 'bdydyn2d', l send3,lrecv3, va, 'V', -1. )90 IF( ANY(llsend3) .OR. ANY(llrecv3) ) THEN ! if need to send/recv in at least one direction 91 CALL lbc_bdy_lnk( 'bdydyn2d', llsend3, llrecv3, va, 'V', -1. ) 97 92 END IF 98 93 ! … … 363 358 TYPE(OBC_INDEX), INTENT(in) :: idx ! OBC indices 364 359 INTEGER, INTENT(in) :: ib_bdy ! BDY set index 365 366 INTEGER :: jb, igrd ! dummy loop indices 367 LOGICAL, DIMENSION(4) :: lsend2, lrecv2, lsend3, lrecv3 ! indicate how communications are to be carried out 360 INTEGER :: igrd ! dummy indice 368 361 !!---------------------------------------------------------------------- 369 362 ! -
NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/BDY/bdyice.F90
r11067 r11071 56 56 ! 57 57 INTEGER :: jbdy ! BDY set index 58 LOGICAL, DIMENSION(4) :: l send1,lrecv1 ! indicate how communications are to be carried out58 LOGICAL, DIMENSION(4) :: llsend1, llrecv1 ! indicate how communications are to be carried out 59 59 !!---------------------------------------------------------------------- 60 60 ! controls … … 76 76 ! 77 77 ! Update bdy points 78 l send1(:) = .false.79 l recv1(:) = .false.78 llsend1(:) = .false. 79 llrecv1(:) = .false. 80 80 DO jbdy = 1, nb_bdy 81 81 IF( cn_ice(jbdy) == 'frs' ) THEN 82 l send1(:) = lsend1(:) .OR. lsend_bdy(jbdy,1,:) ! to every neighbour, T points83 l recv1(:) = lrecv1(:) .OR. lrecv_bdy(jbdy,1,:) ! from every neighbour, T points82 llsend1(:) = llsend1(:) .OR. lsend_bdyint(jbdy,1,:) ! possibly every direction, T points 83 llrecv1(:) = llrecv1(:) .OR. lrecv_bdyint(jbdy,1,:) ! possibly every direction, T points 84 84 END IF 85 85 END DO 86 IF( ANY(l send1) .OR. ANY(lrecv1) ) THEN ! if need to send/recv in at least one direction86 IF( ANY(llsend1) .OR. ANY(llrecv1) ) THEN ! if need to send/recv in at least one direction 87 87 ! exchange 3d arrays 88 CALL lbc_bdy_lnk_multi( 'bdyice', l send1,lrecv1, a_i , 'T', 1., h_i , 'T', 1., h_s , 'T', 1. &89 & , oa_i, 'T', 1., a_ip, 'T', 1., v_ip, 'T', 1. &90 & , s_i , 'T', 1., t_su, 'T', 1., v_i , 'T', 1. &91 & , v_s , 'T', 1., sv_i, 'T', 1. )88 CALL lbc_bdy_lnk_multi( 'bdyice', llsend1, llrecv1, a_i , 'T', 1., h_i , 'T', 1., h_s , 'T', 1. & 89 & , oa_i, 'T', 1., a_ip, 'T', 1., v_ip, 'T', 1. & 90 & , s_i , 'T', 1., t_su, 'T', 1., v_i , 'T', 1. & 91 & , v_s , 'T', 1., sv_i, 'T', 1. ) 92 92 ! exchange 4d arrays 93 CALL lbc_bdy_lnk_multi( 'bdyice', l send1,lrecv1, t_s , 'T', 1., e_s , 'T', 1. ) ! third dimension = 194 CALL lbc_bdy_lnk_multi( 'bdyice', l send1, lrecv1, t_i , 'T', 1., e_i , 'T', 1. ) ! third dimension = 293 CALL lbc_bdy_lnk_multi( 'bdyice', llsend1, llrecv1, t_s , 'T', 1., e_s , 'T', 1. ) ! third dimension = 1 94 CALL lbc_bdy_lnk_multi( 'bdyice', llsend1, llrecv1, t_i , 'T', 1., e_i , 'T', 1. ) ! third dimension = jpk 95 95 END IF 96 96 ! … … 172 172 ! 173 173 IF( flagu == 1. ) THEN 174 IF( ji+1 > jpi 174 IF( ji+1 > jpi ) CYCLE 175 175 IF( u_ice(ji ,jj ) < 0. ) jpbound = 1 ; ib = ji+1 176 176 END IF 177 177 IF( flagu == -1. ) THEN 178 IF( ji-1 < 1 178 IF( ji-1 < 1 ) CYCLE 179 179 IF( u_ice(ji-1,jj ) < 0. ) jpbound = 1 ; ib = ji-1 180 180 END IF 181 181 IF( flagv == 1. ) THEN 182 IF( j i+1 > jpj ) CYCLE182 IF( jj+1 > jpj ) CYCLE 183 183 IF( v_ice(ji ,jj ) < 0. ) jpbound = 1 ; jb = jj+1 184 184 END IF … … 299 299 INTEGER :: jbdy ! BDY set index 300 300 REAL(wp) :: zmsk1, zmsk2, zflag 301 LOGICAL, DIMENSION(4) :: l send2, lrecv2, lsend3,lrecv3 ! indicate how communications are to be carried out301 LOGICAL, DIMENSION(4) :: llsend2, llrecv2, llsend3, llrecv3 ! indicate how communications are to be carried out 302 302 !!------------------------------------------------------------------------------ 303 303 IF( ln_timing ) CALL timing_start('bdy_ice_dyn') … … 384 384 SELECT CASE ( cd_type ) 385 385 CASE ( 'U' ) 386 l send2(:) = .false. ;lrecv2(:) = .false.386 llsend2(:) = .false. ; llrecv2(:) = .false. 387 387 DO jbdy = 1, nb_bdy 388 388 IF( cn_ice(jbdy) == 'frs' .AND. nn_ice_dta(jbdy) /= 0 ) THEN 389 lsend2(:) = lsend2(:) .OR. lsend_bdy(jbdy,2,:) ! to every bdy neighbour, U points 390 lrecv2(:) = lrecv2(:) .OR. lrecv_bdy(jbdy,2,:) ! from every bdy neighbour, U points 389 llsend2(:) = llsend2(:) .OR. lsend_bdyint(jbdy,2,:) ! possibly every direction, U points 390 llsend2(1) = llsend2(1) .OR. lsend_bdyext(jbdy,2,1) ! neighbour might search point towards bdy on its east 391 llrecv2(:) = llrecv2(:) .OR. lrecv_bdyint(jbdy,2,:) ! possibly every direction, U points 392 llrecv2(2) = llrecv2(2) .OR. lrecv_bdyext(jbdy,2,2) ! might search point towards bdy on the east 391 393 END IF 392 394 END DO 393 IF( ANY(l send2) .OR. ANY(lrecv2) ) THEN ! if need to send/recv in at least one direction394 CALL lbc_bdy_lnk( 'bdyice', l send2,lrecv2, u_ice, 'U', -1. )395 IF( ANY(llsend2) .OR. ANY(llrecv2) ) THEN ! if need to send/recv in at least one direction 396 CALL lbc_bdy_lnk( 'bdyice', llsend2, llrecv2, u_ice, 'U', -1. ) 395 397 END IF 396 398 CASE ( 'V' ) 397 l send3(:) = .false. ;lrecv3(:) = .false.399 llsend3(:) = .false. ; llrecv3(:) = .false. 398 400 DO jbdy = 1, nb_bdy 399 401 IF( cn_ice(jbdy) == 'frs' .AND. nn_ice_dta(jbdy) /= 0 ) THEN 400 lsend3(:) = lsend3(:) .OR. lsend_bdy(jbdy,3,:) ! to every bdy neighbour, V points 401 lrecv3(:) = lrecv3(:) .OR. lrecv_bdy(jbdy,3,:) ! from every bdy neighbour, V points 402 llsend3(:) = llsend3(:) .OR. lsend_bdyint(jbdy,3,:) ! possibly every direction, V points 403 llsend3(3) = llsend3(3) .OR. lsend_bdyext(jbdy,3,3) ! neighbour might search point towards bdy on its north 404 llrecv3(:) = llrecv3(:) .OR. lrecv_bdyint(jbdy,3,:) ! possibly every direction, V points 405 llrecv3(4) = llrecv3(4) .OR. lrecv_bdyext(jbdy,3,4) ! might search point towards bdy on the north 402 406 END IF 403 407 END DO 404 IF( ANY(l send3) .OR. ANY(lrecv3) ) THEN ! if need to send/recv in at least one direction405 CALL lbc_bdy_lnk( 'bdyice', l send3,lrecv3, v_ice, 'V', -1. )408 IF( ANY(llsend3) .OR. ANY(llrecv3) ) THEN ! if need to send/recv in at least one direction 409 CALL lbc_bdy_lnk( 'bdyice', llsend3, llrecv3, v_ice, 'V', -1. ) 406 410 END IF 407 411 END SELECT -
NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/BDY/bdyini.F90
r11067 r11071 132 132 INTEGER :: ib_bdy1, ib_bdy2, ib1, ib2 ! - - 133 133 INTEGER :: i_offset, j_offset, inbdy, itreat ! - - 134 INTEGER :: ii1, ii2, ii3, ij1, ij2, ij3, iibe, ijbe ! - - 135 INTEGER :: flagu, flagv ! short cuts 134 136 INTEGER , POINTER :: nbi, nbj, nbr ! short cuts 135 137 REAL(wp), POINTER, DIMENSION(:,:) :: pmask ! pointer to 2D mask fields … … 145 147 REAL(wp), TARGET, DIMENSION(jpi,jpj) :: zfmask ! temporary fmask array excluding coastal boundary condition (shlat) 146 148 REAL(wp) , DIMENSION(jpi,jpj) :: ztmp 147 REAL(wp), POINTER :: flagu, flagv ! short cuts148 149 LOGICAL :: llnobdy, llsobdy, lleabdy, llwebdy ! local logicals 149 150 !! … … 867 868 ! Initialize array indicating communications in bdy 868 869 ! ------------------------------------------------- 869 870 ! Allocate array indicating if a send instruction is needed in bdy treatment 871 ALLOCATE( nbondi_bdy(nb_bdy) ) 872 ALLOCATE( nbondj_bdy(nb_bdy) ) 873 nbondi_bdy(:)=2 874 nbondj_bdy(:)=2 875 ! Allocate array indicating if a receive instruction is needed in bdy treatment 876 ALLOCATE( nbondi_bdy_b(nb_bdy)) 877 ALLOCATE( nbondj_bdy_b(nb_bdy)) 878 nbondi_bdy_b(:)=2 879 nbondj_bdy_b(:)=2 870 ALLOCATE( lsend_bdy(nb_bdy,jpbgrd,4), lrecv_bdy(nb_bdy,jpbgrd,4) ) 871 lsend_bdy(:,:,:) = .false. 872 lrecv_bdy(:,:,:) = .false. 880 873 881 874 DO ib_bdy = 1, nb_bdy 882 ! default : no send883 llsend_ea = .false.884 llsend_we = .false.885 llsend_so = .false.886 llsend_no = .false.887 ! default : no receive888 llrecv_ea = .false.889 llrecv_we = .false.890 llrecv_so = .false.891 llrecv_no = .false.892 875 DO igrd = 1, jpbgrd 893 876 DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd) ! only the rim triggers communications, see bdy routines … … 896 879 ! 897 880 ! check if point has to be sent to a neighbour 881 ! W neighbour and on the inner left side 882 IF( ii == 2 .and. (nbondi == 0 .or. nbondi == 1) ) lsend_bdy(ib_bdy,igrd,1) = .true. 898 883 ! E neighbour and on the inner right side 899 IF( ii == nlci-1 .and. (nbondi == 0 .or. nbondi == -1) ) llsend_ea= .true.900 ! W neighbour and on the inner leftside901 IF( i i == 2 .and. (nbondi == 0 .or. nbondi == 1) ) llsend_we= .true.884 IF( ii == jpi-1 .and. (nbondi == 0 .or. nbondi == -1) ) lsend_bdy(ib_bdy,igrd,2) = .true. 885 ! S neighbour and on the inner down side 886 IF( ij == 2 .and. (nbondj == 0 .or. nbondj == 1) ) lsend_bdy(ib_bdy,igrd,3) = .true. 902 887 ! N neighbour and on the inner up side 903 IF( ij == nlcj-1 .and. (nbondj == 0 .or. nbondj == -1) ) llsend_no = .true. 904 ! S neighbour and on the inner down side 905 IF( ij == 2 .and. (nbondj == 0 .or. nbondj == 1) ) llsend_so = .true. 888 IF( ij == jpj-1 .and. (nbondj == 0 .or. nbondj == -1) ) lsend_bdy(ib_bdy,igrd,4) = .true. 906 889 ! 907 890 ! check if point has to be received from a neighbour 891 ! W neighbour and on the outter left side 892 IF( ii == 1 .and. (nbondi == 0 .or. nbondi == 1) ) lrecv_bdy(ib_bdy,igrd,1) = .true. 908 893 ! E neighbour and on the outter right side 909 IF( ii == nlci .and. (nbondi == 0 .or. nbondi == -1) ) llrecv_ea= .true.910 ! W neighbour and on the outter leftside911 IF( i i == 1 .and. (nbondi == 0 .or. nbondi == 1) ) llrecv_we= .true.894 IF( ii == jpi .and. (nbondi == 0 .or. nbondi == -1) ) lrecv_bdy(ib_bdy,igrd,2) = .true. 895 ! S neighbour and on the outter down side 896 IF( ij == 1 .and. (nbondj == 0 .or. nbondj == 1) ) lrecv_bdy(ib_bdy,igrd,3) = .true. 912 897 ! N neighbour and on the outter up side 913 IF( ij == nlcj .and. (nbondj == 0 .or. nbondj == -1) ) llrecv_no = .true. 914 ! S neighbour and on the outter down side 915 IF( ij == 1 .and. (nbondj == 0 .or. nbondj == 1) ) llrecv_so = .true. 898 IF( ij == jpj .and. (nbondj == 0 .or. nbondj == -1) ) lrecv_bdy(ib_bdy,igrd,4) = .true. 916 899 ! 917 900 END DO 918 901 END DO ! igrd 919 920 ! definition of the i- and j- direction local boundaries arrays used for sending the boundaries921 IF( llsend_ea .and. llsend_we ) THEN ; nbondi_bdy(ib_bdy) = 0922 ELSEIF( llsend_ea .and. .not. llsend_we ) THEN ; nbondi_bdy(ib_bdy) = -1923 ELSEIF( .not. llsend_ea .and. llsend_we ) THEN ; nbondi_bdy(ib_bdy) = 1924 ENDIF925 IF( llsend_no .and. llsend_so ) THEN ; nbondj_bdy(ib_bdy) = 0926 ELSEIF( llsend_no .and. .not. llsend_so ) THEN ; nbondj_bdy(ib_bdy) = -1927 ELSEIF( .not. llsend_no .and. llsend_so ) THEN ; nbondj_bdy(ib_bdy) = 1928 ENDIF929 930 ! definition of the i- and j- direction local boundaries arrays used for receiving the boundaries931 IF( llrecv_ea .and. llrecv_we ) THEN ; nbondi_bdy_b(ib_bdy) = 0932 ELSEIF( llrecv_ea .and. .not. llrecv_we ) THEN ; nbondi_bdy_b(ib_bdy) = -1933 ELSEIF( .not. llrecv_ea .and. llrecv_we ) THEN ; nbondi_bdy_b(ib_bdy) = 1934 ENDIF935 IF( llrecv_no .and. llrecv_so ) THEN ; nbondj_bdy_b(ib_bdy) = 0936 ELSEIF( llrecv_no .and. .not. llrecv_so ) THEN ; nbondj_bdy_b(ib_bdy) = -1937 ELSEIF( .not. llrecv_no .and. llrecv_so ) THEN ; nbondj_bdy_b(ib_bdy) = 1938 ENDIF939 902 940 903 ! Compute rim weights for FRS scheme … … 1134 1097 ! 1 | o ! 2 o | ! 3 | x ! 4 x | 1135 1098 ! |_x_ _ ! _ _x_| ! | o ! o | 1136 IF( pmask(ii+1,ij+1) == 1. ) ztmp(ii,ij) = 1 1137 IF( pmask(ii-1,ij+1) == 1. ) ztmp(ii,ij) = 2 1138 IF( pmask(ii+1,ij-1) == 1. ) ztmp(ii,ij) = 3 1139 IF( pmask(ii-1,ij-1) == 1. ) ztmp(ii,ij) = 4 1099 IF( pmask(ii+1,ij+1) == 1. ) ztmp(ii,ij) = 1. 1100 IF( pmask(ii-1,ij+1) == 1. ) ztmp(ii,ij) = 2. 1101 IF( pmask(ii+1,ij-1) == 1. ) ztmp(ii,ij) = 3. 1102 IF( pmask(ii-1,ij-1) == 1. ) ztmp(ii,ij) = 4. 1140 1103 END IF 1141 1104 IF( inbdy == 1 ) THEN ! middle of linear bdy 1142 ztmp(ii,ij) = 0 ! regular treatment with flags1105 ztmp(ii,ij) = 0. ! regular treatment with flags 1143 1106 END IF 1144 1107 IF( inbdy == 2 ) THEN ! exterior of a corner … … 1146 1109 ! 5 ____x o ! 6 o x___ ! 7 x o ! 8 o x 1147 1110 ! | ! | ! o ! o 1148 IF( llnobdy .AND. lleabdy ) ztmp(ii,ij) = 5 1149 IF( llnobdy .AND. llwebdy ) ztmp(ii,ij) = 6 1150 IF( llsobdy .AND. lleabdy ) ztmp(ii,ij) = 7 1151 IF( llsobdy .AND. llwebdy ) ztmp(ii,ij) = 8 1111 IF( llnobdy .AND. lleabdy ) ztmp(ii,ij) = 5. 1112 IF( llnobdy .AND. llwebdy ) ztmp(ii,ij) = 6. 1113 IF( llsobdy .AND. lleabdy ) ztmp(ii,ij) = 7. 1114 IF( llsobdy .AND. llwebdy ) ztmp(ii,ij) = 8. 1152 1115 END IF 1153 1116 IF( inbdy == 3 ) THEN ! 3 neighbours __ __ … … 1155 1118 ! 9 _| x o ! 10 o x |_ ! 11 o x o ! 12 o x o 1156 1119 ! | o ! o | ! o ! __|¨|__ 1157 IF( llnobdy .AND. lleabdy .AND. llsobdy ) ztmp(ii,ij) = 9 1158 IF( llnobdy .AND. llwebdy .AND. llsobdy ) ztmp(ii,ij) = 10 1159 IF( llwebdy .AND. llsobdy .AND. lleabdy ) ztmp(ii,ij) = 11 1160 IF( llwebdy .AND. llnobdy .AND. lleabdy ) ztmp(ii,ij) = 12 1120 IF( llnobdy .AND. lleabdy .AND. llsobdy ) ztmp(ii,ij) = 9. 1121 IF( llnobdy .AND. llwebdy .AND. llsobdy ) ztmp(ii,ij) = 10. 1122 IF( llwebdy .AND. llsobdy .AND. lleabdy ) ztmp(ii,ij) = 11. 1123 IF( llwebdy .AND. llnobdy .AND. lleabdy ) ztmp(ii,ij) = 12. 1161 1124 END IF 1162 1125 IF( inbdy == 4 ) THEN … … 1171 1134 ii = idx_bdy(ib_bdy)%nbi(ib,igrd) 1172 1135 ij = idx_bdy(ib_bdy)%nbj(ib,igrd) 1173 idx_bdy(ib_bdy)%ntreat(ib,igrd) = ztmp(ii,ij)1136 idx_bdy(ib_bdy)%ntreat(ib,igrd) = NINT(ztmp(ii,ij)) 1174 1137 END DO 1175 1138 END DO … … 1177 1140 1178 1141 1179 ALLOCATE( lsend_bdy(nb_bdy,jpbgrd,4), lrecv_bdy(nb_bdy,jpbgrd,4) ) 1180 lsend_bdy(:,:,:) = .false. 1181 lrecv_bdy(:,:,:) = .false. 1142 ALLOCATE( lsend_bdyint(nb_bdy,jpbgrd,4), lrecv_bdyint(nb_bdy,jpbgrd,4) ) 1143 lsend_bdyint(:,:,:) = .false. 1144 lrecv_bdyint(:,:,:) = .false. 1145 ALLOCATE( lsend_bdyext(nb_bdy,jpbgrd,4), lrecv_bdyext(nb_bdy,jpbgrd,4) ) 1146 lsend_bdyext(:,:,:) = .false. 1147 lrecv_bdyext(:,:,:) = .false. 1182 1148 ! 1183 1149 ! Check which boundaries might need communication 1184 1150 DO igrd = 1, jpbgrd 1185 1151 DO ib_bdy = 1, nb_bdy 1186 IF ( nbondi_bdy (ib_bdy) == 0 ) THEN 1187 lsend_bdy(ib_bdy,igrd,1) = .true. 1188 lsend_bdy(ib_bdy,igrd,2) = .true. 1189 ELSE IF( nbondi_bdy (ib_bdy) == 1 ) THEN 1190 lsend_bdy(ib_bdy,igrd,1) = .true. 1191 ELSE IF( nbondi_bdy (ib_bdy) == -1 ) THEN 1192 lsend_bdy(ib_bdy,igrd,2) = .true. 1193 END IF 1194 IF ( nbondi_bdy_b(ib_bdy) == 0 ) THEN 1195 lrecv_bdy(ib_bdy,igrd,1) = .true. 1196 lrecv_bdy(ib_bdy,igrd,2) = .true. 1197 ELSE IF( nbondi_bdy_b(ib_bdy) == 1 ) THEN 1198 lrecv_bdy(ib_bdy,igrd,1) = .true. 1199 ELSE IF( nbondi_bdy_b(ib_bdy) == -1 ) THEN 1200 lrecv_bdy(ib_bdy,igrd,2) = .true. 1201 END IF 1202 IF( nbondj_bdy (ib_bdy) == 0 ) THEN 1203 lsend_bdy(ib_bdy,igrd,3) = .true. 1204 lsend_bdy(ib_bdy,igrd,4) = .true. 1205 ELSE IF( nbondj_bdy (ib_bdy) == 1 ) THEN 1206 lsend_bdy(ib_bdy,igrd,3) = .true. 1207 ELSE IF( nbondj_bdy (ib_bdy) == -1 ) THEN 1208 lsend_bdy(ib_bdy,igrd,4) = .true. 1209 END IF 1210 IF( nbondj_bdy_b(ib_bdy) == 0 ) THEN 1211 lrecv_bdy(ib_bdy,igrd,3) = .true. 1212 lrecv_bdy(ib_bdy,igrd,4) = .true. 1213 ELSE IF( nbondj_bdy_b(ib_bdy) == 1 ) THEN 1214 lrecv_bdy(ib_bdy,igrd,3) = .true. 1215 ELSE IF( nbondj_bdy_b(ib_bdy) == -1 ) THEN 1216 lrecv_bdy(ib_bdy,igrd,4) = .true. 1217 END IF 1152 DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd) 1153 ii = idx_bdy(ib_bdy)%nbi(ib,igrd) 1154 ij = idx_bdy(ib_bdy)%nbj(ib,igrd) 1155 flagu = NINT(idx_bdy(ib_bdy)%flagu(ib,igrd)) 1156 flagv = NINT(idx_bdy(ib_bdy)%flagv(ib,igrd)) 1157 iibe = ii - flagu ! neighbouring point towards the exterior of the computational domain 1158 ijbe = ij - flagv ! neighbouring point towards the exterior of the computational domain 1159 SELECT CASE( idx_bdy(ib_bdy)%ntreat(ib,igrd) ) ! points that will be used by bdy routines, -1 will be discarded 1160 CASE( 0 ) ; ii1 = ii + flagu ; ii2 = -1 ; ii3 = -1 1161 ij1 = ij + flagv ; ij2 = -1 ; ij3 = -1 1162 CASE( 1 ) ; ii1 = ii+1 ; ij1 = ij+1 ; ii2 = -1 ; ij2 = -1 ; ii3 = -1 ; ij3 = -1 1163 CASE( 2 ) ; ii1 = ii-1 ; ij1 = ij+1 ; ii2 = -1 ; ij2 = -1 ; ii3 = -1 ; ij3 = -1 1164 CASE( 3 ) ; ii1 = ii+1 ; ij1 = ij-1 ; ii2 = -1 ; ij2 = -1 ; ii3 = -1 ; ij3 = -1 1165 CASE( 4 ) ; ii1 = ii-1 ; ij1 = ij-1 ; ii2 = -1 ; ij2 = -1 ; ii3 = -1 ; ij3 = -1 1166 CASE( 5 ) ; ii1 = ii ; ij1 = ij+1 ; ii2 = ii+1 ; ij2 = ij ; ii3 = -1 ; ij3 = -1 1167 CASE( 6 ) ; ii1 = ii ; ij1 = ij+1 ; ii2 = ii-1 ; ij2 = ij ; ii3 = -1 ; ij3 = -1 1168 CASE( 7 ) ; ii1 = ii ; ij1 = ij-1 ; ii2 = ii+1 ; ij2 = ij ; ii3 = -1 ; ij3 = -1 1169 CASE( 8 ) ; ii1 = ii ; ij1 = ij-1 ; ii2 = ii-1 ; ij2 = ij ; ii3 = -1 ; ij3 = -1 1170 CASE( 9 ) ; ii1 = ii ; ij1 = ij+1 ; ii2 = ii+1 ; ij2 = ij ; ii3 = ii ; ij3 = ij-1 1171 CASE( 10 ) ; ii1 = ii ; ij1 = ij+1 ; ii2 = ii-1 ; ij2 = ij ; ii3 = ii ; ij3 = ij-1 1172 CASE( 11 ) ; ii1 = ii-1 ; ij1 = ij ; ii2 = ii ; ij2 = ij-1 ; ii3 = ii+1 ; ij3 = ij 1173 CASE( 12 ) ; ii1 = ii-1 ; ij1 = ij ; ii2 = ii ; ij2 = ij+1 ; ii3 = ii+1 ; ij3 = ij 1174 END SELECT 1175 ! 1176 ! search neighbour in the west/east direction 1177 ! Rim is on the halo and computed ocean is towards exterior of mpi domain 1178 ! <-- (o exterior) --> 1179 ! (1) o|x OR (2) x|o 1180 ! |___ ___| 1181 IF( ii1 == 0 .OR. ii2 == 0 .OR. ii3 == 0 ) lrecv_bdyint(ib_bdy,igrd,1) = .true. 1182 IF( ii1 == jpi+1 .OR. ii2 == jpi+1 .OR. ii3 == jpi+1 ) lrecv_bdyint(ib_bdy,igrd,2) = .true. 1183 IF( iibe == 0 ) lrecv_bdyext(ib_bdy,igrd,1) = .true. 1184 IF( iibe == jpi+1 ) lrecv_bdyext(ib_bdy,igrd,2) = .true. 1185 ! Check if neighbour has its rim parallel to its mpi subdomain border and located next to its halo 1186 ! :¨¨¨¨¨|¨¨--> | | <--¨¨|¨¨¨¨¨: 1187 ! : | x:o | neighbour limited by ... would need o | o:x | : 1188 ! :.....|_._:_____| (1) W neighbour E neighbour (2) |_____:_._|.....: 1189 IF( ii == 2 .AND. ( nbondi == 1 .OR. nbondi == 0 ) & 1190 & .AND. ( ii1 == 3 .OR. ii2 == 3 .OR. ii3 == 3 ) ) lsend_bdyint(ib_bdy,igrd,1) = .true. 1191 IF( ii == jpi-1 .AND. ( nbondi == -1 .OR. nbondi == 0 ) & 1192 & .AND. ( ii1 == jpi-2 .OR. ii2 == jpi-2 .OR. ii3 == jpi-2 ) ) lsend_bdyint(ib_bdy,igrd,2) = .true. 1193 IF( ii == 2 .AND. ( nbondi == 1 .OR. nbondi == 0 ) .AND. iibe == 3 ) lsend_bdyext(ib_bdy,igrd,1) = .true. 1194 IF( ii == jpi-1 .AND. ( nbondi == -1 .OR. nbondi == 0 ) .AND. iibe == jpi-2 ) lsend_bdyext(ib_bdy,igrd,2) = .true. 1195 ! 1196 ! search neighbour in the north/south direction 1197 ! Rim is on the halo and computed ocean is towards exterior of mpi domain 1198 !(3) | | ^ ___o___ 1199 ! | |___x___| OR | | x | 1200 ! v o (4) | | 1201 IF( ij1 == 0 .OR. ij2 == 0 .OR. ij3 == 0 ) lrecv_bdyint(ib_bdy,igrd,3) = .true. 1202 IF( ij1 == jpj+1 .OR. ij2 == jpj+1 .OR. ij3 == jpj+1 ) lrecv_bdyint(ib_bdy,igrd,4) = .true. 1203 IF( ijbe == 0 ) lrecv_bdyext(ib_bdy,igrd,3) = .true. 1204 IF( ijbe == jpj+1 ) lrecv_bdyext(ib_bdy,igrd,4) = .true. 1205 ! Check if neighbour has its rim parallel to its mpi subdomain _________ border and next to its halo 1206 ! ^ | o | : : 1207 ! | |¨¨¨¨x¨¨¨¨| neighbour limited by ... would need o | |....x....| 1208 ! :_________: (3) S neighbour N neighbour (4) v | o | 1209 IF( ij == 2 .AND. ( nbondj == 1 .OR. nbondj == 0 ) & 1210 & .AND. ( ij1 == 3 .OR. ij2 == 3 .OR. ij3 == 3 ) ) lsend_bdyint(ib_bdy,igrd,3) = .true. 1211 IF( ij == jpj-1 .AND. ( nbondj == -1 .OR. nbondj == 0 ) & 1212 & .AND. ( ij1 == jpj-2 .OR. ij2 == jpj-2 .OR. ij3 == jpj-2 ) ) lsend_bdyint(ib_bdy,igrd,4) = .true. 1213 IF( ij == 2 .AND. ( nbondj == 1 .OR. nbondj == 0 ) .AND. ijbe == 3 ) lsend_bdyext(ib_bdy,igrd,3) = .true. 1214 IF( ij == jpj-1 .AND. ( nbondj == -1 .OR. nbondj == 0 ) .AND. ijbe == jpj-2 ) lsend_bdyext(ib_bdy,igrd,4) = .true. 1215 END DO 1218 1216 END DO 1219 1217 END DO 1218 1219 DO ib_bdy = 1,nb_bdy 1220 IF( cn_dyn2d(ib_bdy) == 'orlanski' .OR. cn_dyn2d(ib_bdy) == 'orlanski_npo' .OR. & 1221 & cn_dyn3d(ib_bdy) == 'orlanski' .OR. cn_dyn3d(ib_bdy) == 'orlanski_npo' .OR. & 1222 & cn_tra(ib_bdy) == 'orlanski' .OR. cn_tra(ib_bdy) == 'orlanski_npo' ) THEN 1223 DO igrd = 1, jpbgrd 1224 DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd) 1225 ii = idx_bdy(ib_bdy)%nbi(ib,igrd) 1226 ij = idx_bdy(ib_bdy)%nbj(ib,igrd) 1227 IF( mig(ii) > 2 .AND. mig(ii) < jpiglo-2 .AND. mjg(ij) > 2 .AND. mjg(ij) < jpjglo-2 ) THEN 1228 WRITE(ctmp1,*) ' Orlanski can not be used when the open boundaries are on the interior of the computational domain' 1229 WRITE(ctmp2,*) ' ========== ' 1230 CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' ) 1231 END IF 1232 END DO 1233 END DO 1234 END IF 1235 END DO 1236 1220 1237 ! 1221 1238 ! Tidy up -
NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/BDY/bdylib.F90
r11049 r11071 187 187 ii = idx%nbi(jb,igrd) 188 188 ij = idx%nbj(jb,igrd) 189 IF( ii == 1 .OR. ii == jpi .OR. ij == 1 .OR. ij == jpj ) CYCLE 189 190 flagu = int( idx%flagu(jb,igrd) ) 190 191 flagv = int( idx%flagv(jb,igrd) ) … … 344 345 ii = idx%nbi(jb,igrd) 345 346 ij = idx%nbj(jb,igrd) 347 IF( ii == 1 .OR. ii == jpi .OR. ij == 1 .OR. ij == jpj ) CYCLE 346 348 flagu = int( idx%flagu(jb,igrd) ) 347 349 flagv = int( idx%flagv(jb,igrd) ) -
NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/BDY/bdytra.F90
r11067 r11071 51 51 INTEGER :: ib_bdy, jn, igrd ! Loop indeces 52 52 TYPE(ztrabdy), DIMENSION(jpts) :: zdta ! Temporary data structure 53 LOGICAL, DIMENSION(4) :: l send1,lrecv1 ! indicate how communications are to be carried out53 LOGICAL, DIMENSION(4) :: llsend1, llrecv1 ! indicate how communications are to be carried out 54 54 !!---------------------------------------------------------------------- 55 55 igrd = 1 … … 76 76 END DO 77 77 ! 78 l send1(:) = .false.79 l recv1(:) = .false.78 llsend1(:) = .false. 79 llrecv1(:) = .false. 80 80 DO ib_bdy=1, nb_bdy 81 81 SELECT CASE( TRIM(cn_tra(ib_bdy)) ) 82 CASE('neumann') 83 lsend1(:) = lsend1(:) .OR. lsend_bdy(ib_bdy,1,:) ! to every bdy neighbour, T points 84 lrecv1(:) = lrecv1(:) .OR. lrecv_bdy(ib_bdy,1,:) ! from every bdy neighbour, T points 85 CASE('orlanski') 86 lsend1(:) = lsend1(:) .OR. lsend_bdy(ib_bdy,1,:) ! to every bdy neighbour, T points 87 lrecv1(:) = lrecv1(:) .OR. lrecv_bdy(ib_bdy,1,:) ! from every bdy neighbour, T points 88 CASE('orlanski_npo') 89 lsend1(:) = lsend1(:) .OR. lsend_bdy(ib_bdy,1,:) ! to every bdy neighbour, T points 90 lrecv1(:) = lrecv1(:) .OR. lrecv_bdy(ib_bdy,1,:) ! from every bdy neighbour, T points 91 CASE('runoff') 92 lsend1(:) = lsend1(:) .OR. lsend_bdy(ib_bdy,1,:) ! to every bdy neighbour, T points 93 lrecv1(:) = lrecv1(:) .OR. lrecv_bdy(ib_bdy,1,:) ! from every bdy neighbour, T points 82 CASE('neumann','runoff') 83 llsend1(:) = llsend1(:) .OR. lsend_bdyint(ib_bdy,1,:) ! possibly every direction, T points 84 llrecv1(:) = llrecv1(:) .OR. lrecv_bdyint(ib_bdy,1,:) ! possibly every direction, T points 85 CASE('orlanski', 'orlanski_npo') 86 llsend1(:) = llsend1(:) .OR. lsend_bdy(ib_bdy,1,:) ! possibly every direction, T points 87 llrecv1(:) = llrecv1(:) .OR. lrecv_bdy(ib_bdy,1,:) ! possibly every direction, T points 94 88 END SELECT 95 89 END DO 96 IF( ANY(l send1) .OR. ANY(lrecv1) ) THEN ! if need to send/recv in at least one direction97 CALL lbc_bdy_lnk( 'bdytra', l send1,lrecv1, tsa, 'T', 1. )90 IF( ANY(llsend1) .OR. ANY(llrecv1) ) THEN ! if need to send/recv in at least one direction 91 CALL lbc_bdy_lnk( 'bdytra', llsend1, llrecv1, tsa, 'T', 1. ) 98 92 END IF 99 93 ! -
NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/DYN/dynkeg.F90
r11067 r11071 80 80 REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: ztrdu, ztrdv 81 81 REAL(wp) :: zweightu, zweightv 82 LOGICAL, DIMENSION(4) :: l send1,lrecv1 ! indicate how bdy communications are to be carried out82 LOGICAL, DIMENSION(4) :: llsend1, llrecv1 ! indicate how bdy communications are to be carried out 83 83 !!---------------------------------------------------------------------- 84 84 ! … … 136 136 END IF 137 137 END DO 138 ! send 2 and recv jpi, jpjused in the computation of the speed tendencies139 l send1(:) = .false.140 l recv1(:) = .false.138 ! send jpi-1, jpj-1 and receive 1 used in the computation of the speed tendencies 139 llsend1(:) = .false. 140 llrecv1(:) = .false. 141 141 DO ib_bdy = 1, nb_bdy 142 lsend1(:) = lsend1(:) .OR. lsend_bdy(ib_bdy,1,:) ! to every bdy neighbour, T points 143 lrecv1(:) = lrecv1(:) .OR. lrecv_bdy(ib_bdy,1,:) ! from every bdy neighbour, T points 144 END DO 145 IF( COUNT(lsend1) > 0 .OR. COUNT(lrecv1) > 0 ) THEN ! if need to send/recv in at least one direction 146 CALL lbc_bdy_lnk( 'bdydyn2d', lsend1, lrecv1, zhke, 'T', 1. ) 142 llsend1(2) = llsend1(2) .OR. lsend_bdy(ib_bdy,igrd,2) ! send east 143 llsend1(4) = llsend1(4) .OR. lsend_bdy(ib_bdy,igrd,4) ! send north 144 llrecv1(1) = llrecv1(1) .OR. lrecv_bdy(ib_bdy,igrd,1) ! receive west 145 llrecv1(3) = llrecv1(3) .OR. lrecv_bdy(ib_bdy,igrd,3) ! receive south 146 END DO 147 148 IF( ANY(llsend1) .OR. ANY(llrecv1) ) THEN ! if need to send/recv in at least one direction 149 CALL lbc_bdy_lnk( 'bdydyn2d', llsend1, llrecv1, zhke, 'T', 1. ) 147 150 END IF 148 151 END IF -
NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/TOP/trcbdy.F90
r11067 r11071 46 46 INTEGER, INTENT( in ) :: kt ! Main time step counter 47 47 !! 48 INTEGER :: ib_bdy ,jn ,igrd ! Loop ind eces48 INTEGER :: ib_bdy ,jn ,igrd ! Loop indices 49 49 REAL(wp), POINTER, DIMENSION(:,:) :: ztrc 50 50 REAL(wp), POINTER :: zfac 51 LOGICAL, DIMENSION(4) :: l send1,lrecv1 ! indicate how communications are to be carried out51 LOGICAL, DIMENSION(4) :: llsend1, llrecv1 ! indicate how communications are to be carried out 52 52 !!---------------------------------------------------------------------- 53 53 ! … … 75 75 END DO 76 76 ! 77 l send1(:) = .false.78 l recv1(:) = .false.77 llsend1(:) = .false. 78 llrecv1(:) = .false. 79 79 DO ib_bdy=1, nb_bdy 80 80 SELECT CASE( TRIM(cn_tra(ib_bdy)) ) 81 81 CASE('neumann') 82 lsend1(:) = lsend1(:) .OR. lsend_bdy(ib_bdy,1,:) ! to every bdy neighbour, T points 83 lrecv1(:) = lrecv1(:) .OR. lrecv_bdy(ib_bdy,1,:) ! from every bdy neighbour, T points 84 CASE('orlanski') 85 lsend1(:) = lsend1(:) .OR. lsend_bdy(ib_bdy,1,:) ! to every bdy neighbour, T points 86 lrecv1(:) = lrecv1(:) .OR. lrecv_bdy(ib_bdy,1,:) ! from every bdy neighbour, T points 87 CASE('orlanski_npo') 88 lsend1(:) = lsend1(:) .OR. lsend_bdy(ib_bdy,1,:) ! to every bdy neighbour, T points 89 lrecv1(:) = lrecv1(:) .OR. lrecv_bdy(ib_bdy,1,:) ! from every bdy neighbour, T points 82 llsend1(:) = llsend1(:) .OR. lsend_bdyint(ib_bdy,1,:) ! possibly every direction, T points 83 llrecv1(:) = llrecv1(:) .OR. lrecv_bdyint(ib_bdy,1,:) ! possibly every direction, T points 84 CASE('orlanski','orlanski_npo') 85 llsend1(:) = llsend1(:) .OR. lsend_bdy(ib_bdy,1,:) ! possibly every direction, T points 86 llrecv1(:) = llrecv1(:) .OR. lrecv_bdy(ib_bdy,1,:) ! possibly every direction, T points 90 87 END SELECT 91 88 END DO 92 IF( ANY(l send1) .OR. ANY(lrecv1) ) THEN ! if need to send/recv in at least one direction93 CALL lbc_bdy_lnk( 'bdytra', l send1,lrecv1, tsa, 'T', 1. )89 IF( ANY(llsend1) .OR. ANY(llrecv1) ) THEN ! if need to send/recv in at least one direction 90 CALL lbc_bdy_lnk( 'bdytra', llsend1, llrecv1, tsa, 'T', 1. ) 94 91 END IF 95 92 !
Note: See TracChangeset
for help on using the changeset viewer.