FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
hecmw_solver_SR.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2019 FrontISTR Commons
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
5 
6 !C
7 !C***
8 !C*** module hecmw_solver_SR_mm
9 !C***
10 !C
12  use hecmw_util
13 
14 #ifndef HECMW_SERIAL
15  type async_buf
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
26  end type async_buf
27 
28  integer(kind=kint), parameter :: max_nreq = 8
29  type(async_buf), save :: abuf(max_nreq)
30 #endif
31 
32 contains
33  !C
34  !C*** SOLVER_SEND_RECV
35  !C
36  subroutine hecmw_solve_send_recv &
37  & ( n, m, neibpetot, neibpe, &
38  & stack_import, nod_import, &
39  & stack_export, nod_export, ws, wr, x, &
40  & solver_comm,my_rank)
41 
42  use hecmw_util
43 #ifdef _OPENACC
44  use openacc
45 #endif
46  implicit none
47  ! include 'mpif.h'
48  ! include 'hecmw_config_f.h'
49 
50  integer(kind=kint ) , intent(in) :: N, m
51  integer(kind=kint ) , intent(in) :: NEIBPETOT
52  integer(kind=kint ), pointer :: NEIBPE (:)
53  integer(kind=kint ), pointer :: STACK_IMPORT(:)
54  integer(kind=kint ), pointer :: NOD_IMPORT (:)
55  integer(kind=kint ), pointer :: STACK_EXPORT(:)
56  integer(kind=kint ), pointer :: NOD_EXPORT (:)
57  real (kind=kreal), dimension(: ), intent(inout):: ws
58  real (kind=kreal), dimension(: ), intent(inout):: wr
59  real (kind=kreal), dimension(: ), intent(inout):: x
60  integer(kind=kint ) , intent(in) ::SOLVER_COMM
61  integer(kind=kint ) , intent(in) :: my_rank
62 
63 #ifndef HECMW_SERIAL
64  integer(kind=kint ), dimension(:,:), allocatable :: sta1
65  integer(kind=kint ), dimension(:,:), allocatable :: sta2
66  integer(kind=kint ), dimension(: ), allocatable :: req1
67  integer(kind=kint ), dimension(: ), allocatable :: req2
68 
69  integer(kind=kint ), save :: NFLAG
70  data nflag/0/
71  ! local valiables
72  integer(kind=kint ) :: neib,istart,inum,k,kk,ii,ierr,nreq1,nreq2
73 #ifdef _OPENACC
74  integer(kind=kint):: ns, nr
75  type(c_devptr) :: WS_dev, WR_dev
76 #endif
77  !C
78  !C-- INIT.
79  allocate (sta1(mpi_status_size,neibpetot))
80  allocate (sta2(mpi_status_size,neibpetot))
81  allocate (req1(neibpetot))
82  allocate (req2(neibpetot))
83 
84 #ifdef _OPENACC
85  ns = stack_export(neibpetot)
86  nr = stack_import(neibpetot)
87 
88  ws_dev = acc_malloc(kreal * m * ns)
89  wr_dev = acc_malloc(kreal * m * nr)
90 
91  call acc_map_data(ws, ws_dev, kreal * m * ns)
92  call acc_map_data(wr, wr_dev, kreal * m * nr)
93 #endif
94 
95  !C
96  !C-- SEND
97  nreq1=0
98  do neib= 1, neibpetot
99  istart= stack_export(neib-1)
100  inum = stack_export(neib ) - istart
101  if (inum==0) cycle
102  nreq1=nreq1+1
103 #ifdef _OPENACC
104  !$acc kernels
105  !$acc loop independent collapse(2)
106  do k= istart+1, istart+inum
107  do kk = 1, m
108  ii = nod_export(k)
109  ws(m*(k-1)+kk)= x(m*(ii-1)+kk)
110  enddo
111  enddo
112  !$acc end kernels
113 #else
114  do k= istart+1, istart+inum
115  ii = nod_export(k)
116  do kk = 1, m
117  ws(m*(k-1)+kk)= x(m*(ii-1)+kk)
118  enddo
119  enddo
120 #endif
121 
122  !$acc host_data use_device(WS)
123  call mpi_isend (ws(m*istart+1), m*inum,mpi_double_precision, &
124  & neibpe(neib), 0, solver_comm, req1(nreq1), ierr)
125  !$acc end host_data
126  enddo
127 
128  !C
129  !C-- RECEIVE
130  nreq2=0
131  do neib= 1, neibpetot
132  istart= stack_import(neib-1)
133  inum = stack_import(neib ) - istart
134  if (inum==0) cycle
135  nreq2=nreq2+1
136  !$acc host_data use_device(WR)
137  call mpi_irecv (wr(m*istart+1), m*inum, mpi_double_precision, &
138  & neibpe(neib), 0, solver_comm, req2(nreq2), ierr)
139  !$acc end host_data
140  enddo
141 
142  call mpi_waitall (nreq2, req2, sta2, ierr)
143 
144  do neib= 1, neibpetot
145  istart= stack_import(neib-1)
146  inum = stack_import(neib ) - istart
147 #ifdef _OPENACC
148  !$acc kernels
149  !$acc loop independent collapse(2)
150  do k= istart+1, istart+inum
151  do kk= 1, m
152  ii = nod_import(k)
153  x(m*(ii-1)+kk)= wr(m*(k-1)+kk)
154  enddo
155  enddo
156  !$acc end kernels
157 #else
158  do k= istart+1, istart+inum
159  ii = nod_import(k)
160  do kk= 1, m
161  x(m*(ii-1)+kk)= wr(m*(k-1)+kk)
162  enddo
163  enddo
164 #endif
165  enddo
166 
167  call mpi_waitall (nreq1, req1, sta1, ierr)
168 
169  deallocate (sta1, sta2, req1, req2)
170 
171 #ifdef _OPENACC
172  call acc_unmap_data(ws)
173  call acc_unmap_data(wr)
174 
175  call acc_free(ws_dev)
176  call acc_free(wr_dev)
177 #endif
178 #endif
179  end subroutine hecmw_solve_send_recv
180 
181  !C
182  !C*** SOLVER_ISEND_IRECV
183  !C
185  & ( n, m, neibpetot, neibpe, stack_import, nod_import, &
186  & stack_export, nod_export, &
187  & x, solver_comm,my_rank,ireq)
188  implicit none
189  integer(kind=kint ) , intent(in) :: N, M
190  integer(kind=kint ) , intent(in) :: NEIBPETOT
191  integer(kind=kint ), pointer :: NEIBPE (:)
192  integer(kind=kint ), pointer :: STACK_IMPORT(:)
193  integer(kind=kint ), pointer :: NOD_IMPORT (:)
194  integer(kind=kint ), pointer :: STACK_EXPORT(:)
195  integer(kind=kint ), pointer :: NOD_EXPORT (:)
196  real (kind=kreal), target, intent(inout):: x(:)
197  integer(kind=kint ) , intent(in) ::SOLVER_COMM
198  integer(kind=kint ) , intent(in) :: my_rank
199  integer(kind=kint ) , intent(out) :: ireq
200 
201 #ifndef HECMW_SERIAL
202  ! local valiables
203  real (kind=kreal), pointer :: ws(:)
204  real (kind=kreal), pointer :: wr(:)
205  integer(kind=kint ), pointer :: req1(:)
206  integer(kind=kint ), pointer :: req2(:)
207  integer(kind=kint ) :: neib,istart,inum,k,kk,ii,ierr,i,nreq1,nreq2
208  !C
209  !C-- INIT.
210  allocate (ws(m*stack_export(neibpetot)), stat=ierr)
211  if( ierr /= 0 ) stop "Allocation error: WS in hecmw_solver_SR"
212  allocate (wr(m*stack_import(neibpetot)), stat=ierr)
213  if( ierr /= 0 ) stop "Allocation error: WR in hecmw_solver_SR"
214  allocate (req1(neibpetot))
215  allocate (req2(neibpetot))
216  !C
217  !C-- SEND
218  nreq1=0
219  do neib= 1, neibpetot
220  istart= stack_export(neib-1)
221  inum = stack_export(neib ) - istart
222  if (inum==0) cycle
223  nreq1=nreq1+1
224  do k= istart+1, istart+inum
225  ii = nod_export(k)
226  do kk = 1, m
227  ws(m*(k-1)+kk)= x(m*(ii-1)+kk)
228  enddo
229  enddo
230  call mpi_isend (ws(m*istart+1), m*inum,mpi_double_precision, &
231  & neibpe(neib), 0, solver_comm, req1(nreq1), ierr)
232  enddo
233  !C
234  !C-- RECEIVE
235  nreq2=0
236  do neib= 1, neibpetot
237  istart= stack_import(neib-1)
238  inum = stack_import(neib ) - istart
239  if (inum==0) cycle
240  nreq2=nreq2+1
241  call mpi_irecv (wr(m*istart+1), m*inum, mpi_double_precision, &
242  & neibpe(neib), 0, solver_comm, req2(nreq2), ierr)
243  enddo
244  !C
245  !C-- Find empty abuf
246  ireq = 0
247  do i = 1, max_nreq
248  if (abuf(i)%NEIBPETOT == 0) then
249  ireq = i
250  exit
251  endif
252  enddo
253  if (ireq == 0) then
254  stop 'Error: HECMW_SOLVE_ISEND_IRECV: exceeded maximum num of requests'
255  endif
256  !C
257  !C-- Save in abuf
258  abuf(ireq)%NEIBPETOT = neibpetot
259  abuf(ireq)%STACK_IMPORT=> stack_import
260  abuf(ireq)%NOD_IMPORT => nod_import
261  abuf(ireq)%WS => ws
262  abuf(ireq)%WR => wr
263  abuf(ireq)%X => x
264  abuf(ireq)%req1 => req1
265  abuf(ireq)%req2 => req2
266  abuf(ireq)%nreq1 = nreq1
267  abuf(ireq)%nreq2 = nreq2
268 #endif
269  end subroutine hecmw_solve_isend_irecv
270 
271  !C
272  !C*** SOLVER_ISEND_IRECV_WAIT
273  !C
274  subroutine hecmw_solve_isend_irecv_wait( m, ireq )
275  use hecmw_util
276  implicit none
277  integer(kind=kint ), intent(in) :: m, ireq
278 
279 #ifndef HECMW_SERIAL
280  ! local valiables
281  integer(kind=kint ) :: NEIBPETOT
282  integer(kind=kint ), pointer :: STACK_IMPORT(:)
283  integer(kind=kint ), pointer :: NOD_IMPORT (:)
284  real (kind=kreal), pointer :: ws(:)
285  real (kind=kreal), pointer :: wr(:)
286  real (kind=kreal), pointer :: x(:)
287  integer(kind=kint ), pointer :: req1(:)
288  integer(kind=kint ), pointer :: req2(:)
289  integer(kind=kint ), dimension(:,:), allocatable :: sta1
290  integer(kind=kint ), dimension(:,:), allocatable :: sta2
291  integer(kind=kint ) :: neib,istart,inum,k,j,ii,ierr,nreq1,nreq2
292  !C-- Check ireq
293  if (ireq < 0 .or. ireq > max_nreq) then
294  stop 'ERROR: HECMW_SOLVE_ISEND_IRECV_WAIT: invalid ireq'
295  endif
296  !C-- Restore from abuf
297  neibpetot = abuf(ireq)%NEIBPETOT
298  stack_import=> abuf(ireq)%STACK_IMPORT
299  nod_import => abuf(ireq)%NOD_IMPORT
300  ws => abuf(ireq)%WS
301  wr => abuf(ireq)%WR
302  x => abuf(ireq)%X
303  req1 => abuf(ireq)%req1
304  req2 => abuf(ireq)%req2
305  nreq1 = abuf(ireq)%nreq1
306  nreq2 = abuf(ireq)%nreq2
307  !C-- Free abuf
308  abuf(ireq)%NEIBPETOT = 0
309  !C-- Allocate
310  allocate (sta1(mpi_status_size,neibpetot))
311  allocate (sta2(mpi_status_size,neibpetot))
312  !C-- Wait irecv
313  call mpi_waitall (nreq2, req2, sta2, ierr)
314  do neib= 1, neibpetot
315  istart= stack_import(neib-1)
316  inum = stack_import(neib ) - istart
317  do k= istart+1, istart+inum
318  ii = nod_import(k)
319  do j = 1, m
320  x(m*(ii-1)+j)= wr(m*(k-1)+j)
321  enddo
322  enddo
323  enddo
324  !C-- Wait isend
325  call mpi_waitall (nreq1, req1, sta1, ierr)
326  !C-- Deallocate
327  deallocate (sta1, sta2)
328  deallocate (req1, req2)
329  deallocate (ws, wr)
330 #endif
331  end subroutine hecmw_solve_isend_irecv_wait
332 
333  !C
334  !C*** SOLVER_REVERSE_SEND_RECV
335  !C
337  & ( n, m, neibpetot, neibpe, stack_import, nod_import, &
338  & stack_export, nod_export, &
339  & ws, wr, x, solver_comm,my_rank)
340 
341  implicit none
342 
343  integer(kind=kint ) , intent(in) :: N, M
344  integer(kind=kint ) , intent(in) :: NEIBPETOT
345  integer(kind=kint ), pointer :: NEIBPE (:)
346  integer(kind=kint ), pointer :: STACK_IMPORT(:)
347  integer(kind=kint ), pointer :: NOD_IMPORT (:)
348  integer(kind=kint ), pointer :: STACK_EXPORT(:)
349  integer(kind=kint ), pointer :: NOD_EXPORT (:)
350  real (kind=kreal), dimension(: ), intent(inout):: ws
351  real (kind=kreal), dimension(: ), intent(inout):: wr
352  real (kind=kreal), dimension(: ), intent(inout):: x
353  integer(kind=kint ) , intent(in) ::SOLVER_COMM
354  integer(kind=kint ) , intent(in) :: my_rank
355 
356 #ifndef HECMW_SERIAL
357  integer(kind=kint ), dimension(:,:), allocatable :: sta1
358  integer(kind=kint ), dimension(:,:), allocatable :: sta2
359  integer(kind=kint ), dimension(: ), allocatable :: req1
360  integer(kind=kint ), dimension(: ), allocatable :: req2
361 
362  ! local valiables
363  integer(kind=kint ) :: neib,istart,inum,k,kk,ii,ierr,nreq1,nreq2
364  !C
365  !C-- INIT.
366  allocate (sta1(mpi_status_size,neibpetot))
367  allocate (sta2(mpi_status_size,neibpetot))
368  allocate (req1(neibpetot))
369  allocate (req2(neibpetot))
370 
371  !C
372  !C-- SEND
373  nreq1=0
374  do neib= 1, neibpetot
375  istart= stack_import(neib-1)
376  inum = stack_import(neib ) - istart
377  if (inum==0) cycle
378  nreq1=nreq1+1
379  do k= istart+1, istart+inum
380  ii = nod_import(k)
381  do kk = 1, m
382  ws(m*(k-1)+kk)= x(m*(ii-1)+kk)
383  enddo
384  enddo
385 
386  call mpi_isend (ws(m*istart+1), m*inum,mpi_double_precision, &
387  & neibpe(neib), 0, solver_comm, req1(nreq1), ierr)
388  enddo
389 
390  !C
391  !C-- RECEIVE
392  nreq2=0
393  do neib= 1, neibpetot
394  istart= stack_export(neib-1)
395  inum = stack_export(neib ) - istart
396  if (inum==0) cycle
397  nreq2=nreq2+1
398  call mpi_irecv (wr(m*istart+1), m*inum, mpi_double_precision, &
399  & neibpe(neib), 0, solver_comm, req2(nreq2), ierr)
400  enddo
401 
402  call mpi_waitall (nreq2, req2, sta2, ierr)
403 
404  do neib= 1, neibpetot
405  istart= stack_export(neib-1)
406  inum = stack_export(neib ) - istart
407  do k= istart+1, istart+inum
408  ii = nod_export(k)
409  do kk = 1, m
410  x(m*(ii-1)+kk)= x(m*(ii-1)+kk)+wr(m*(k-1)+kk)
411  enddo
412  enddo
413  enddo
414 
415  call mpi_waitall (nreq1, req1, sta1, ierr)
416  deallocate (sta1, sta2, req1, req2)
417 #endif
418  end subroutine hecmw_solve_rev_send_recv
419 
420 end module hecmw_solver_sr
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)
I/O and Utility.
Definition: hecmw_util_f.F90:7
integer(kind=4), parameter kreal