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)), stat=ierr)
157 if( ierr /= 0 ) stop
"Allocation error: WS in hecmw_solver_SR"
158 allocate (wr(m*stack_import(neibpetot)), stat=ierr)
159 if( ierr /= 0 ) stop
"Allocation error: WR in hecmw_solver_SR"
160 allocate (req1(neibpetot))
161 allocate (req2(neibpetot))
165 do neib= 1, neibpetot
166 istart= stack_export(neib-1)
167 inum = stack_export(neib ) - istart
170 do k= istart+1, istart+inum
173 ws(m*(k-1)+kk)= x(m*(ii-1)+kk)
176 call mpi_isend (ws(m*istart+1), m*inum,mpi_double_precision, &
177 & neibpe(neib), 0, solver_comm, req1(nreq1), ierr)
182 do neib= 1, neibpetot
183 istart= stack_import(neib-1)
184 inum = stack_import(neib ) - istart
187 call mpi_irecv (wr(m*istart+1), m*inum, mpi_double_precision, &
188 & neibpe(neib), 0, solver_comm, req2(nreq2), ierr)
194 if (
abuf(i)%NEIBPETOT == 0)
then
200 stop
'Error: HECMW_SOLVE_ISEND_IRECV: exceeded maximum num of requests'
204 abuf(ireq)%NEIBPETOT = neibpetot
205 abuf(ireq)%STACK_IMPORT=> stack_import
206 abuf(ireq)%NOD_IMPORT => nod_import
210 abuf(ireq)%req1 => req1
211 abuf(ireq)%req2 => req2
212 abuf(ireq)%nreq1 = nreq1
213 abuf(ireq)%nreq2 = nreq2
223 integer(kind=kint ),
intent(in) :: m, ireq
227 integer(kind=kint ) :: NEIBPETOT
228 integer(kind=kint ),
pointer :: STACK_IMPORT(:)
229 integer(kind=kint ),
pointer :: NOD_IMPORT (:)
230 real (kind=
kreal),
pointer :: ws(:)
231 real (kind=
kreal),
pointer :: wr(:)
232 real (kind=
kreal),
pointer :: x(:)
233 integer(kind=kint ),
pointer :: req1(:)
234 integer(kind=kint ),
pointer :: req2(:)
235 integer(kind=kint ),
dimension(:,:),
allocatable :: sta1
236 integer(kind=kint ),
dimension(:,:),
allocatable :: sta2
237 integer(kind=kint ) :: neib,istart,inum,k,j,ii,ierr,nreq1,nreq2
239 if (ireq < 0 .or. ireq >
max_nreq)
then
240 stop
'ERROR: HECMW_SOLVE_ISEND_IRECV_WAIT: invalid ireq'
243 neibpetot =
abuf(ireq)%NEIBPETOT
244 stack_import=>
abuf(ireq)%STACK_IMPORT
245 nod_import =>
abuf(ireq)%NOD_IMPORT
249 req1 =>
abuf(ireq)%req1
250 req2 =>
abuf(ireq)%req2
251 nreq1 =
abuf(ireq)%nreq1
252 nreq2 =
abuf(ireq)%nreq2
254 abuf(ireq)%NEIBPETOT = 0
256 allocate (sta1(mpi_status_size,neibpetot))
257 allocate (sta2(mpi_status_size,neibpetot))
259 call mpi_waitall (nreq2, req2, sta2, ierr)
260 do neib= 1, neibpetot
261 istart= stack_import(neib-1)
262 inum = stack_import(neib ) - istart
263 do k= istart+1, istart+inum
266 x(m*(ii-1)+j)= wr(m*(k-1)+j)
271 call mpi_waitall (nreq1, req1, sta1, ierr)
273 deallocate (sta1, sta2)
274 deallocate (req1, req2)
283 & ( n, m, neibpetot, neibpe, stack_import, nod_import, &
284 & stack_export, nod_export, &
285 & ws, wr, x, solver_comm,my_rank)
289 integer(kind=kint ) ,
intent(in) :: N, M
290 integer(kind=kint ) ,
intent(in) :: NEIBPETOT
291 integer(kind=kint ),
pointer :: NEIBPE (:)
292 integer(kind=kint ),
pointer :: STACK_IMPORT(:)
293 integer(kind=kint ),
pointer :: NOD_IMPORT (:)
294 integer(kind=kint ),
pointer :: STACK_EXPORT(:)
295 integer(kind=kint ),
pointer :: NOD_EXPORT (:)
296 real (kind=
kreal),
dimension(: ),
intent(inout):: ws
297 real (kind=
kreal),
dimension(: ),
intent(inout):: wr
298 real (kind=
kreal),
dimension(: ),
intent(inout):: x
299 integer(kind=kint ) ,
intent(in) ::SOLVER_COMM
300 integer(kind=kint ) ,
intent(in) :: my_rank
303 integer(kind=kint ),
dimension(:,:),
allocatable :: sta1
304 integer(kind=kint ),
dimension(:,:),
allocatable :: sta2
305 integer(kind=kint ),
dimension(: ),
allocatable :: req1
306 integer(kind=kint ),
dimension(: ),
allocatable :: req2
309 integer(kind=kint ) :: neib,istart,inum,k,kk,ii,ierr,nreq1,nreq2
312 allocate (sta1(mpi_status_size,neibpetot))
313 allocate (sta2(mpi_status_size,neibpetot))
314 allocate (req1(neibpetot))
315 allocate (req2(neibpetot))
320 do neib= 1, neibpetot
321 istart= stack_import(neib-1)
322 inum = stack_import(neib ) - istart
325 do k= istart+1, istart+inum
328 ws(m*(k-1)+kk)= x(m*(ii-1)+kk)
332 call mpi_isend (ws(m*istart+1), m*inum,mpi_double_precision, &
333 & neibpe(neib), 0, solver_comm, req1(nreq1), ierr)
339 do neib= 1, neibpetot
340 istart= stack_export(neib-1)
341 inum = stack_export(neib ) - istart
344 call mpi_irecv (wr(m*istart+1), m*inum, mpi_double_precision, &
345 & neibpe(neib), 0, solver_comm, req2(nreq2), ierr)
348 call mpi_waitall (nreq2, req2, sta2, ierr)
350 do neib= 1, neibpetot
351 istart= stack_export(neib-1)
352 inum = stack_export(neib ) - istart
353 do k= istart+1, istart+inum
356 x(m*(ii-1)+kk)= x(m*(ii-1)+kk)+wr(m*(k-1)+kk)
361 call mpi_waitall (nreq1, req1, sta1, ierr)
362 deallocate (sta1, sta2, req1, req2)
subroutine hecmw_solve_isend_irecv_wait(m, ireq)
type(async_buf), dimension(max_nreq), save abuf
integer(kind=kint), parameter max_nreq
subroutine hecmw_solve_isend_irecv(N, M, NEIBPETOT, NEIBPE, STACK_IMPORT, NOD_IMPORT, STACK_EXPORT, NOD_EXPORT, X, SOLVER_COMM, my_rank, ireq)
subroutine hecmw_solve_send_recv(N, m, NEIBPETOT, NEIBPE, STACK_IMPORT, NOD_IMPORT, STACK_EXPORT, NOD_EXPORT, WS, WR, X, SOLVER_COMM, my_rank)
subroutine hecmw_solve_rev_send_recv(N, M, NEIBPETOT, NEIBPE, STACK_IMPORT, NOD_IMPORT, STACK_EXPORT, NOD_EXPORT, WS, WR, X, SOLVER_COMM, my_rank)
integer(kind=4), parameter kreal