UPP (develop)
Loading...
Searching...
No Matches
EXCH.f
Go to the documentation of this file.
1
20!----------------------------------------------------------------------
25!----------------------------------------------------------------------
26 SUBROUTINE exch(A)
27
28 use ctlblk_mod, only: num_procs, jend, iup, jsta, idn, mpi_comm_comp, im,&
29 icoords,ibcoords,bufs,ibufs,me,numx, &
30 jsta_2l, jend_2u,ileft,iright,ista_2l,iend_2u,ista,iend,jm,modelname
31!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
32 implicit none
33!
34 include 'mpif.h'
35!
36 real,intent(inout) :: a ( ista_2l:iend_2u,jsta_2l:jend_2u )
37 real, allocatable :: coll(:), colr(:)
38 integer, allocatable :: icoll(:), icolr(:)
39 integer status(MPI_STATUS_SIZE)
40 integer ierr, jstam1, jendp1,j
41 integer size,ubound,lbound
42 integer msglenl, msglenr
43 integer i,ii,jj, ibl,ibu,jbl,jbu,icc,jcc
44 integer iwest,ieast
45 integer ifirst
46
47 logical, parameter :: checkcoords = .false.
48
49 data ifirst/0/
50 allocate(coll(jm))
51 allocate(colr(jm))
52 allocate(icolr(jm))
53 allocate(icoll(jm))
54 ibl=max(ista-1,1)
55 ibu=min(im,iend+1)
56 jbu=min(jm,jend+1)
57 jbl=max(jsta-1,1)
58!
59
60! write(*,*) 'mype=',me,'num_procs=',num_procs,'im=',im,'jsta_2l=', &
61! jsta_2l,'jend_2u=',jend_2u,'jend=',jend,'iup=',iup,'jsta=', &
62! jsta,'idn=',idn
63 if ( num_procs <= 1 ) return
64!
65! for global model apply cyclic boundary condition
66
67 IF(modelname == 'GFS') then
68 if(ifirst .le. 0 .and. me .eq. 0) print *,' CYCLIC BC APPLIED'
69 if(ileft .eq. mpi_proc_null) iwest=1 ! get eastern bc from western boundary of full domain
70 if(iright .eq. mpi_proc_null) ieast=1 ! get western bc from eastern boundary of full domain
71 if(ileft .eq. mpi_proc_null) ileft=me+(numx-1)
72 if(iright .eq. mpi_proc_null) iright=(me-numx) +1
73 endif
74
75 jstam1 = max(jsta_2l,jsta-1) ! Moorthi
76
77! send last row to iup's first row+ and receive first row- from idn's last row
78
79 call mpi_sendrecv(a(ista,jend),iend-ista+1,mpi_real,iup,1, &
80 & a(ista,jstam1),iend-ista+1,mpi_real,idn,1, &
81 & mpi_comm_comp,status,ierr)
82
83 if ( ierr /= 0 ) then
84 print *, ' problem with first sendrecv in exch, ierr = ',ierr
85 stop 6661
86 endif
87
88 if (checkcoords) then
89 if(ifirst .le. 0) then !IFIRST ONLY
90 call mpi_sendrecv(ibcoords(ista,jend),iend-ista+1,mpi_integer,iup,1, &
91 & ibcoords(ista,jstam1),iend-ista+1,mpi_integer,idn,1, &
92 & mpi_comm_comp,status,ierr)
93 if ( ierr /= 0 ) then
94 print *, ' problem with second sendrecv in exch, ierr = ',ierr
95 stop 7661
96 endif
97 do i=ista,iend
98 ii=ibcoords(i,jstam1)/10000
99 jj=ibcoords(i,jstam1)-(ii*10000)
100 if(ii .ne. i .or. jj .ne. jstam1 ) print *,' GWVX JEXCH CHECK FAIL ',ii,jj,ibcoords(i,jstam1),i
101 end do
102 endif !IFIRST
103 endif !checkcoords
104
105! build the I columns to send and receive
106
107 msglenl=jend-jsta+1
108 msglenr=jend-jsta+1
109 if(iright .lt. 0) msglenr=1
110 if(ileft .lt. 0) msglenl=1
111
112 do j=jsta,jend
113 coll(j)=a(ista,j)
114 end do
115
116 call mpi_barrier(mpi_comm_comp,ierr)
117
118! send first col to ileft last col+ and receive last col+ from ileft first col
119
120 call mpi_sendrecv(coll(jsta),msglenl ,mpi_real,ileft,1, &
121 & colr(jsta),msglenr ,mpi_real,iright,1, &
122 & mpi_comm_comp,status,ierr)
123
124 if ( ierr /= 0 ) then
125 print *, ' problem with third sendrecv in exch, ierr = ',ierr
126 stop 6662
127 endif
128
129 if(ifirst .le. 0) then ! IFIRST ONLY
130 call mpi_sendrecv(icoll(jsta),msglenl ,mpi_integer,ileft,1, &
131 & icolr(jsta),msglenr ,mpi_integer,iright,1, &
132 & mpi_comm_comp,status,ierr)
133 if ( ierr /= 0 ) then
134 print *, ' problem with fourth sendrecv in exch, ierr = ',ierr
135 stop 7662
136 endif
137 endif !IFIRST
138
139 if(iright .ge. 0) then
140 do j=jsta,jend
141 a(iend+1,j)=colr(j)
142 if(checkcoords) then
143 if(ifirst .le. 0) then !IFIRST ONLY
144 ibcoords(iend+1,j)=icolr(j)
145 ii=ibcoords(iend+1,j)/10000
146 jj=ibcoords( iend+1,j)-(ii*10000)
147 if( j .ne. jj .or. ii .ne. iend+1 .and. ii .ne. im .and. ii .ne. 1) &
148 write(*,921) j,iend+1,ii,jj,ibcoords(iend+1,j),'IEXCH COORD FAIL j,iend+1,ii,jj,ibcoord '
149 endif !IFIRST
150 endif !checkcoords
151 end do
152 endif ! for iright
153
154 921 format(5i10,a50)
155
156! print *,'mype=',me,'in EXCH, after first mpi_sendrecv'
157
158 if ( ierr /= 0 ) then
159 print *, ' problem with fifth sendrecv in exch, ierr = ',ierr
160 stop 6663
161 end if
162 jendp1 = min(jend+1,jend_2u) ! Moorthi
163
164!GWV. change from full im row exchange to iend-ista+1 subrow exchange,
165
166 do j=jsta,jend
167 colr(j)=a(iend,j)
168 end do
169
170! send first row to idown's last row+ and receive last row+ from iup's first row
171
172 call mpi_sendrecv(a(ista,jsta),iend-ista+1,mpi_real,idn,1, &
173 & a(ista,jendp1),iend-ista+1,mpi_real,iup,1, &
174 & mpi_comm_comp,status,ierr)
175 if ( ierr /= 0 ) then
176 print *, ' problem with sixth sendrecv in exch, ierr = ',ierr
177 stop 6664
178 endif
179
180 if (checkcoords) then
181 if (ifirst .le. 0) then
182 call mpi_sendrecv(ibcoords(ista,jsta),iend-ista+1,mpi_integer,idn,1, &
183 & ibcoords(ista,jendp1),iend-ista+1,mpi_integer,iup,1, &
184 & mpi_comm_comp,status,ierr)
185 if ( ierr /= 0 ) then
186 print *, ' problem with seventh sendrecv in exch, ierr = ',ierr
187 stop 7664
188 endif
189 endif ! IFIRST
190 endif ! checkcoords
191
192! send last col to iright first col- and receive first col- from ileft last col
193
194 call mpi_sendrecv(colr(jsta),msglenr ,mpi_real,iright,1 , &
195 & coll(jsta),msglenl ,mpi_real,ileft ,1, &
196 & mpi_comm_comp,status,ierr)
197
198 if ( ierr /= 0 ) then
199 print *, ' problem with eighth sendrecv in exch, ierr = ',ierr
200 stop 6665
201 endif
202
203 if (ifirst .le. 0) then
204 call mpi_sendrecv(icolr(jsta),msglenr ,mpi_integer,iright,1 , &
205 & icoll(jsta),msglenl ,mpi_integer,ileft ,1, &
206 & mpi_comm_comp,status,ierr)
207 if ( ierr /= 0 ) then
208 print *, ' problem with ninth sendrecv in exch, ierr = ',ierr
209 stop 7665
210 endif
211 endif !IFIRST
212
213 if(ileft .ge. 0) then
214 do j=jsta,jend
215 a(ista-1,j)=coll(j)
216 if(checkcoords) then
217 if(ifirst .le. 0) then
218 ibcoords(ista-1,j)=icoll(j)
219 ii=ibcoords(ista-1,j)/10000
220 jj=ibcoords( ista-1,j)-(ii*10000)
221 if( j .ne. jj .or. ii .ne. ista-1 .and. ii .ne. im .and. ii .ne. 1) &
222 write(*,921) j,ista-1,ii,jj,ibcoords(ista-1,j),'EXCH COORD FAIL j,ista-1,ii,jj,ibcoord '
223 endif !IFIRST
224 endif !checkcoords
225 end do
226 endif
227
228! interior check
229
230 if(checkcoords) then
231 if(ifirst .le. 0) then
232 do j=jsta,jend
233 do i=ista,iend
234 ii=ibcoords(i,j)/10000
235 jj=ibcoords( i,j)-(ii*10000)
236 if(ii .ne. i .or. jj .ne. j) write(*,151) 'INFAILED IJ ',i,j,ibcoords(i,j),ibl,jbl,ibu,jbu
237 end do
238 end do
239 endif !IFIRST
240 endif !checkcoords
241
242!! corner points. After the exchanges above, corner points are replicated in
243! neighbour halos so we can get them from the neighbors rather than
244! calculating more corner neighbor numbers
245! A(ista-1,jsta-1) is in the ileft a(iend,jsta-1) location
246! A(ista-1,jend+1) is in the ileft a(iend,jend+1) location
247! A(iend+1,jsta-1) is in the iright a(ista,jsta-1) location
248! A(iend+1,jend+1) is in the iright a(ista,jend+1) location
249!GWVx ibl=max(ista-1,1)
250!GWVx ibu=min(im,iend+1)
251
252 ibl=max(ista-1,1)
253 ibu=min(im,iend+1)
254 if(modelname == 'GFS') then
255 ibl=max(ista-1,0)
256 ibu=min(im+1,iend+1)
257 endif
258
259 jbu=min(jm,jend+1)
260 jbl=max(jsta-1,1)
261
262 call mpi_sendrecv(a(iend,jbl ),1, mpi_real,iright,1 , &
263 & a(ibl ,jbl ),1, mpi_real,ileft ,1, &
264 & mpi_comm_comp,status,ierr)
265 if ( ierr /= 0 ) then
266 print *, ' problem with tenth sendrecv in exch, ierr = ',ierr
267 stop 6771
268 endif
269
270 call mpi_sendrecv(a(iend,jbu ),1, mpi_real,iright,1 , &
271 & a(ibl ,jbu ),1, mpi_real,ileft ,1, &
272 & mpi_comm_comp,status,ierr)
273
274 if ( ierr /= 0 ) then
275 print *, ' problem with eleventh sendrecv in exch, ierr = ',ierr
276 stop 6772
277 endif
278
279 call mpi_sendrecv(a(ista,jbl ),1, mpi_real,ileft ,1, &
280 & a(ibu ,jbl ),1, mpi_real,iright,1, &
281 & mpi_comm_comp,status,ierr)
282
283 if ( ierr /= 0 ) then
284 print *, ' problem with twelft sendrecv in exch, ierr = ',ierr
285 stop 6773
286 endif
287
288 call mpi_sendrecv(a(ista,jbu ),1, mpi_real,ileft ,1 , &
289 & a(ibu ,jbu ),1, mpi_real,iright,1, &
290 & mpi_comm_comp,status,ierr)
291
292 if ( ierr /= 0 ) then
293 print *, ' problem with thirteenth sendrecv in exch, ierr = ',ierr
294 stop 6774
295 endif
296
297 139 format(a20,5(i10,i6,i6,'<>'))
298
299 if(checkcoords) then
300 if(ifirst .le. 0) then
301 call mpi_sendrecv(ibcoords(iend,jbl ),1 ,mpi_integer,iright,1 , &
302 & ibcoords(ibl ,jbl ),1 ,mpi_integer,ileft ,1, &
303 & mpi_comm_comp,status,ierr)
304
305 call mpi_sendrecv(ibcoords(iend,jbu ),1 ,mpi_integer,iright,1, &
306 & ibcoords(ibl ,jbu ),1 ,mpi_integer,ileft ,1, &
307 & mpi_comm_comp,status,ierr)
308 call mpi_sendrecv(ibcoords(ista,jbl ),1 ,mpi_integer,ileft ,1, &
309 & ibcoords(ibu ,jbl ),1 ,mpi_integer,iright,1, &
310 & mpi_comm_comp,status,ierr)
311 call mpi_sendrecv(ibcoords(ista,jbu ),1 ,mpi_integer,ileft ,1 , &
312 & ibcoords(ibu ,jbu ),1 ,mpi_integer,iright,1, &
313 mpi_comm_comp,status,ierr)
314
315! corner check for coordnates
316
317 icc=ibl
318 jcc=jbl
319 ii=ibcoords(icc,jcc)/10000
320 jj=ibcoords(icc,jcc)-(ii*10000)
321
322 if(ii .ne. icc .and. icc .ne. 0) write(*,151) ' CORNER FAILI ilb ll ',icc,jcc,ibcoords(icc,jcc),ii,jj
323 if( jj .ne. jcc) write(*,151) ' CORNER FAILJ ilb ll ',icc,jcc,ibcoords(icc,jcc),ii,jj
324
325 icc=ibu
326 jcc=jbl
327 ii=ibcoords(icc,jcc)/10000
328 jj=ibcoords(icc,jcc)-(ii*10000)
329 if(ii .ne. icc .and. icc .ne. im+1 ) write(*,151) ' CORNER FAILI ilb ul ',icc,jcc,ibcoords(icc,jcc),ii,jj
330 if( jj .ne. jcc ) write(*,151) ' CORNER FAILJ ilb ul ',icc,jcc,ibcoords(icc,jcc),ii,jj
331
332 icc=ibu
333 jcc=jbu
334 ii=ibcoords(icc,jcc)/10000
335 jj=ibcoords(icc,jcc)-(ii*10000)
336 if(ii .ne. icc .and. icc .ne. im+1) write(*,151) ' CORNER FAILI ilb uu ',icc,jcc,ibcoords(icc,jcc),ii,jj
337 if( jj .ne. jcc ) write(*,151) ' CORNER FAILJ ilb ul ',icc,jcc,ibcoords(icc,jcc),ii,jj
338
339 icc=ibl
340 jcc=jbu
341 ii=ibcoords(icc,jcc)/10000.
342 jj=ibcoords(icc,jcc)-(ii*10000)
343 if(ii .ne. icc .and. icc .ne. 0 ) write(*,151) ' CORNER FAILI ilb lu ',icc,jcc,ibcoords(icc,jcc),ii,jj
344 if( jj .ne. jcc ) write(*,151) ' CORNER FAILJ ilb ul ',icc,jcc,ibcoords(icc,jcc),ii,jj
345
346! if(ileft .ge. 0) then
347!119 format(' GWX LEFT EXCHANGE ileft,me,ibcoords(ista-1,jend+1),ibcoords(ista-1,jend-1),ista-1,jend-1,jend+1', &
348! 10i10)
349! endif
350
351! if(iright .ge. 0) then
352!! write(*,129) iright,me,ibcoords(ista+1,jend+1),ibcoords(ista+1,jend-1),ista-1,jend-1,jend+1 !GWVX
353!129 format(' GWX RIGHT EXCHANGE iright,me,ibcoords(ista+1,jend+1),ibcoords(ista-1,jend+1),ista-1,jend-1,jend+1', &
354! 10i10)
355! endif
356
357! interior check
358
359 do j=jsta,jend
360 do i=ista,iend
361 ii=ibcoords(i,j)/10000
362 jj=ibcoords( i,j)-(ii*10000)
363 if(ii .ne. i .or. jj .ne. j) write(*,151) 'GWVX FAILED IJ ',i,j,ibcoords(i,j),ibl,jbl,ibu,jbu
364 end do
365 end do
366
367 151 format(a70,10i10)
368
369! bounds check
370! first check top and bottom halo rows
371
372 j=jbu
373 do i=ista,iend
374 ii=ibcoords(i,j)/10000
375 jj=ibcoords( i,j)-(ii*10000)
376 if(ii .ne. i .or. jj .ne. j) write(*,151) 'GWVX FAILEDI JBU IJ ',i,j,ibcoords(i,j),ibl,jbl,ibu,jbu
377 end do
378
379 j=jbl
380 do i=ista,iend
381 ii=ibcoords(i,j)/10000
382 jj=ibcoords( i,j)-(ii*10000)
383 if(ii .ne. i .or. jj .ne. j) write(*,151) 'GWVX FAILEDI JBL IJ ',i,j,ibcoords(i,j),ibl,jbl,ibu,jbu
384 end do
385
386! second and last, check left and right halo columns
387
388 i=ibl
389 do j=jsta,jend
390 ii=ibcoords(i,j)/10000
391 jj=ibcoords( i,j)-(ii*10000)
392 if(ii .ne. i .and. ii .ne. im .or. jj .ne. j) write(*,151) 'GWVX FAILED IBL IJ ',ii,i,j,ibcoords(i,j),ibl,jbl,ibu,jbu
393 end do
394
395 i=ibu
396 do j=jsta,jend
397 ii=ibcoords(i,j)/10000
398 jj=ibcoords( i,j)-(ii*10000)
399 if(ii .ne. i .and. ii .ne. 1 .or. jj .ne. j) write(*,151) 'GWVX FAILED IBU ii i j ibcoords ibl,jbl,ibu,jbu',ii,i,j,ibcoords(i,j),ibl,jbl,ibu,jbu
400 end do
401
402 if(me .eq. 0) write(*,*) ' IFIRST CHECK'
403
404 endif ! IFIRST
405 endif !checkcoords
406
407! end halo checks
408 if ( ierr /= 0 ) then
409 print *, ' problem with second sendrecv in exch, ierr = ',ierr
410 stop
411 end if
412 call mpi_barrier(mpi_comm_comp,ierr)
413 ifirst=ifirst+1
414 end
415
416!!@PROCESS NOCHECK
423!----------------------------------------------------------------------
428!----------------------------------------------------------------------
429 subroutine exch_f(a)
430
431 use ctlblk_mod, only: num_procs, jend, iup, jsta, idn, &
432 & mpi_comm_comp, im, jsta_2l, jend_2u
433!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
434 implicit none
435!
436 include 'mpif.h'
437!
438 real,intent(inout) :: a ( im,jsta_2l:jend_2u )
439 integer status(MPI_STATUS_SIZE)
440 integer ierr, jstam1, jendp1
441!
442 if ( num_procs == 1 ) return
443!
444 jstam1 = max(jsta_2l,jsta-1) ! Moorthi
445 call mpi_sendrecv(a(1,jend),im,mpi_real,iup,1, &
446 & a(1,jstam1),im,mpi_real,idn,1, &
447 & mpi_comm_comp,status,ierr)
448 if ( ierr /= 0 ) then
449 print *, ' problem with first sendrecv in exch, ierr = ',ierr
450 stop
451 end if
452 jendp1=min(jend+1,jend_2u) ! Moorthi
453 call mpi_sendrecv(a(1,jsta),im,mpi_real,idn,1, &
454 & a(1,jendp1),im,mpi_real,iup,1, &
455 & mpi_comm_comp,status,ierr)
456 if ( ierr /= 0 ) then
457 print *, ' problem with second sendrecv in exch, ierr = ',ierr
458 stop
459 end if
460!
461 end
462
subroutine exch(a)
exch() Subroutine that exchanges one halo row.
Definition EXCH.f:27
subroutine exch_f(a)
exch_f() Subroutine that exchanges one halo row.
Definition EXCH.f:430
integer, dimension(:,:), allocatable ibcoords
Arrays that store the coordinates of their elements; used to validate communications; when scattered ...
Definition CTLBLK.f:177
integer, dimension(:,:), allocatable icoords
Arrays that store the coordinates of their elements; used to validate communications; when scattered ...
Definition CTLBLK.f:176