16 integer(kind=kint ) :: neibpetot = 0
17 integer(kind=kint ),
pointer :: stack_import(:)
18 integer(kind=kint ),
pointer :: nod_import (:)
19 real (kind=
kreal),
pointer :: ws(:)
20 real (kind=
kreal),
pointer :: wr(:)
21 real (kind=
kreal),
pointer :: x(:)
22 integer(kind=kint ),
pointer :: req1(:)
23 integer(kind=kint ),
pointer :: req2(:)
24 integer(kind=kint ) :: nreq1
25 integer(kind=kint ) :: nreq2
37 & ( n, m, neibpetot, neibpe, &
38 & stack_import, nod_import, &
39 & stack_export, nod_export, ws, wr, x, &
40 & solver_comm,my_rank)
47 integer(kind=kint ) ,
intent(in) :: N, m
48 integer(kind=kint ) ,
intent(in) :: NEIBPETOT
49 integer(kind=kint ),
pointer :: NEIBPE (:)
50 integer(kind=kint ),
pointer :: STACK_IMPORT(:)
51 integer(kind=kint ),
pointer :: NOD_IMPORT (:)
52 integer(kind=kint ),
pointer :: STACK_EXPORT(:)
53 integer(kind=kint ),
pointer :: NOD_EXPORT (:)
54 real (kind=
kreal),
dimension(: ),
intent(inout):: ws
55 real (kind=
kreal),
dimension(: ),
intent(inout):: wr
56 real (kind=
kreal),
dimension(: ),
intent(inout):: x
57 integer(kind=kint ) ,
intent(in) ::SOLVER_COMM
58 integer(kind=kint ) ,
intent(in) :: my_rank
61 integer(kind=kint ),
dimension(:,:),
allocatable :: sta1
62 integer(kind=kint ),
dimension(:,:),
allocatable :: sta2
63 integer(kind=kint ),
dimension(: ),
allocatable :: req1
64 integer(kind=kint ),
dimension(: ),
allocatable :: req2
66 integer(kind=kint ),
save :: NFLAG
69 integer(kind=kint ) :: neib,istart,inum,k,kk,ii,ierr,nreq1,nreq2
72 allocate (sta1(mpi_status_size,neibpetot))
73 allocate (sta2(mpi_status_size,neibpetot))
74 allocate (req1(neibpetot))
75 allocate (req2(neibpetot))
81 istart= stack_export(neib-1)
82 inum = stack_export(neib ) - istart
85 do k= istart+1, istart+inum
88 ws(m*(k-1)+kk)= x(m*(ii-1)+kk)
92 call mpi_isend (ws(m*istart+1), m*inum,mpi_double_precision, &
93 & neibpe(neib), 0, solver_comm, req1(nreq1), ierr)
100 istart= stack_import(neib-1)
101 inum = stack_import(neib ) - istart
104 call mpi_irecv (wr(m*istart+1), m*inum, mpi_double_precision, &
105 & neibpe(neib), 0, solver_comm, req2(nreq2), ierr)
108 call mpi_waitall (nreq2, req2, sta2, ierr)
110 do neib= 1, neibpetot
111 istart= stack_import(neib-1)
112 inum = stack_import(neib ) - istart
113 do k= istart+1, istart+inum
116 x(m*(ii-1)+kk)= wr(m*(k-1)+kk)
121 call mpi_waitall (nreq1, req1, sta1, ierr)
123 deallocate (sta1, sta2, req1, req2)
131 & ( n, m, neibpetot, neibpe, stack_import, nod_import, &
132 & stack_export, nod_export, &
133 & x, solver_comm,my_rank,ireq)
135 integer(kind=kint ) ,
intent(in) :: N, M
136 integer(kind=kint ) ,
intent(in) :: NEIBPETOT
137 integer(kind=kint ),
pointer :: NEIBPE (:)
138 integer(kind=kint ),
pointer :: STACK_IMPORT(:)
139 integer(kind=kint ),
pointer :: NOD_IMPORT (:)
140 integer(kind=kint ),
pointer :: STACK_EXPORT(:)
141 integer(kind=kint ),
pointer :: NOD_EXPORT (:)
142 real (kind=
kreal),
target,
intent(inout):: x(:)
143 integer(kind=kint ) ,
intent(in) ::SOLVER_COMM
144 integer(kind=kint ) ,
intent(in) :: my_rank
145 integer(kind=kint ) ,
intent(out) :: ireq
149 real (kind=
kreal),
pointer :: ws(:)
150 real (kind=
kreal),
pointer :: wr(:)
151 integer(kind=kint ),
pointer :: req1(:)
152 integer(kind=kint ),
pointer :: req2(:)
153 integer(kind=kint ) :: neib,istart,inum,k,kk,ii,ierr,i,nreq1,nreq2
156 allocate (ws(m*stack_export(neibpetot)))
157 allocate (wr(m*stack_import(neibpetot)))
158 allocate (req1(neibpetot))
159 allocate (req2(neibpetot))
163 do neib= 1, neibpetot
164 istart= stack_export(neib-1)
165 inum = stack_export(neib ) - istart
168 do k= istart+1, istart+inum
171 ws(m*(k-1)+kk)= x(m*(ii-1)+kk)
174 call mpi_isend (ws(m*istart+1), m*inum,mpi_double_precision, &
175 & neibpe(neib), 0, solver_comm, req1(nreq1), ierr)
180 do neib= 1, neibpetot
181 istart= stack_import(neib-1)
182 inum = stack_import(neib ) - istart
185 call mpi_irecv (wr(m*istart+1), m*inum, mpi_double_precision, &
186 & neibpe(neib), 0, solver_comm, req2(nreq2), ierr)
192 if (
abuf(i)%NEIBPETOT == 0)
then
198 stop
'Error: HECMW_SOLVE_ISEND_IRECV: exceeded maximum num of requests'
202 abuf(ireq)%NEIBPETOT = neibpetot
203 abuf(ireq)%STACK_IMPORT=> stack_import
204 abuf(ireq)%NOD_IMPORT => nod_import
208 abuf(ireq)%req1 => req1
209 abuf(ireq)%req2 => req2
210 abuf(ireq)%nreq1 = nreq1
211 abuf(ireq)%nreq2 = nreq2
221 integer(kind=kint ),
intent(in) :: m, ireq
225 integer(kind=kint ) :: NEIBPETOT
226 integer(kind=kint ),
pointer :: STACK_IMPORT(:)
227 integer(kind=kint ),
pointer :: NOD_IMPORT (:)
228 real (kind=
kreal),
pointer :: ws(:)
229 real (kind=
kreal),
pointer :: wr(:)
230 real (kind=
kreal),
pointer :: x(:)
231 integer(kind=kint ),
pointer :: req1(:)
232 integer(kind=kint ),
pointer :: req2(:)
233 integer(kind=kint ),
dimension(:,:),
allocatable :: sta1
234 integer(kind=kint ),
dimension(:,:),
allocatable :: sta2
235 integer(kind=kint ) :: neib,istart,inum,k,j,ii,ierr,nreq1,nreq2
237 if (ireq < 0 .or. ireq >
max_nreq)
then
238 stop
'ERROR: HECMW_SOLVE_ISEND_IRECV_WAIT: invalid ireq'
241 neibpetot =
abuf(ireq)%NEIBPETOT
242 stack_import=>
abuf(ireq)%STACK_IMPORT
243 nod_import =>
abuf(ireq)%NOD_IMPORT
247 req1 =>
abuf(ireq)%req1
248 req2 =>
abuf(ireq)%req2
249 nreq1 =
abuf(ireq)%nreq1
250 nreq2 =
abuf(ireq)%nreq2
252 abuf(ireq)%NEIBPETOT = 0
254 allocate (sta1(mpi_status_size,neibpetot))
255 allocate (sta2(mpi_status_size,neibpetot))
257 call mpi_waitall (nreq2, req2, sta2, ierr)
258 do neib= 1, neibpetot
259 istart= stack_import(neib-1)
260 inum = stack_import(neib ) - istart
261 do k= istart+1, istart+inum
264 x(m*(ii-1)+j)= wr(m*(k-1)+j)
269 call mpi_waitall (nreq1, req1, sta1, ierr)
271 deallocate (sta1, sta2)
272 deallocate (req1, req2)
281 & ( n, m, neibpetot, neibpe, stack_import, nod_import, &
282 & stack_export, nod_export, &
283 & ws, wr, x, solver_comm,my_rank)
287 integer(kind=kint ) ,
intent(in) :: N, M
288 integer(kind=kint ) ,
intent(in) :: NEIBPETOT
289 integer(kind=kint ),
pointer :: NEIBPE (:)
290 integer(kind=kint ),
pointer :: STACK_IMPORT(:)
291 integer(kind=kint ),
pointer :: NOD_IMPORT (:)
292 integer(kind=kint ),
pointer :: STACK_EXPORT(:)
293 integer(kind=kint ),
pointer :: NOD_EXPORT (:)
294 real (kind=
kreal),
dimension(: ),
intent(inout):: ws
295 real (kind=
kreal),
dimension(: ),
intent(inout):: wr
296 real (kind=
kreal),
dimension(: ),
intent(inout):: x
297 integer(kind=kint ) ,
intent(in) ::SOLVER_COMM
298 integer(kind=kint ) ,
intent(in) :: my_rank
301 integer(kind=kint ),
dimension(:,:),
allocatable :: sta1
302 integer(kind=kint ),
dimension(:,:),
allocatable :: sta2
303 integer(kind=kint ),
dimension(: ),
allocatable :: req1
304 integer(kind=kint ),
dimension(: ),
allocatable :: req2
307 integer(kind=kint ) :: neib,istart,inum,k,kk,ii,ierr,nreq1,nreq2
310 allocate (sta1(mpi_status_size,neibpetot))
311 allocate (sta2(mpi_status_size,neibpetot))
312 allocate (req1(neibpetot))
313 allocate (req2(neibpetot))
318 do neib= 1, neibpetot
319 istart= stack_import(neib-1)
320 inum = stack_import(neib ) - istart
323 do k= istart+1, istart+inum
326 ws(m*(k-1)+kk)= x(m*(ii-1)+kk)
330 call mpi_isend (ws(m*istart+1), m*inum,mpi_double_precision, &
331 & neibpe(neib), 0, solver_comm, req1(nreq1), ierr)
337 do neib= 1, neibpetot
338 istart= stack_export(neib-1)
339 inum = stack_export(neib ) - istart
342 call mpi_irecv (wr(m*istart+1), m*inum, mpi_double_precision, &
343 & neibpe(neib), 0, solver_comm, req2(nreq2), ierr)
346 call mpi_waitall (nreq2, req2, sta2, ierr)
348 do neib= 1, neibpetot
349 istart= stack_export(neib-1)
350 inum = stack_export(neib ) - istart
351 do k= istart+1, istart+inum
354 x(m*(ii-1)+kk)= x(m*(ii-1)+kk)+wr(m*(k-1)+kk)
359 call mpi_waitall (nreq1, req1, sta1, ierr)
360 deallocate (sta1, sta2, req1, req2)