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