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