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  implicit none
44  ! include 'mpif.h'
45  ! include 'hecmw_config_f.h'
46 
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
59 
60 #ifndef HECMW_SERIAL
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
65 
66  integer(kind=kint ), save :: NFLAG
67  data nflag/0/
68  ! local valiables
69  integer(kind=kint ) :: neib,istart,inum,k,kk,ii,ierr,nreq1,nreq2
70  !C
71  !C-- INIT.
72  allocate (sta1(mpi_status_size,neibpetot))
73  allocate (sta2(mpi_status_size,neibpetot))
74  allocate (req1(neibpetot))
75  allocate (req2(neibpetot))
76 
77  !C
78  !C-- SEND
79  nreq1=0
80  do neib= 1, neibpetot
81  istart= stack_export(neib-1)
82  inum = stack_export(neib ) - istart
83  if (inum==0) cycle
84  nreq1=nreq1+1
85  do k= istart+1, istart+inum
86  ii = nod_export(k)
87  do kk = 1, m
88  ws(m*(k-1)+kk)= x(m*(ii-1)+kk)
89  enddo
90  enddo
91 
92  call mpi_isend (ws(m*istart+1), m*inum,mpi_double_precision, &
93  & neibpe(neib), 0, solver_comm, req1(nreq1), ierr)
94  enddo
95 
96  !C
97  !C-- RECEIVE
98  nreq2=0
99  do neib= 1, neibpetot
100  istart= stack_import(neib-1)
101  inum = stack_import(neib ) - istart
102  if (inum==0) cycle
103  nreq2=nreq2+1
104  call mpi_irecv (wr(m*istart+1), m*inum, mpi_double_precision, &
105  & neibpe(neib), 0, solver_comm, req2(nreq2), ierr)
106  enddo
107 
108  call mpi_waitall (nreq2, req2, sta2, ierr)
109 
110  do neib= 1, neibpetot
111  istart= stack_import(neib-1)
112  inum = stack_import(neib ) - istart
113  do k= istart+1, istart+inum
114  ii = nod_import(k)
115  do kk= 1, m
116  x(m*(ii-1)+kk)= wr(m*(k-1)+kk)
117  enddo
118  enddo
119  enddo
120 
121  call mpi_waitall (nreq1, req1, sta1, ierr)
122 
123  deallocate (sta1, sta2, req1, req2)
124 #endif
125  end subroutine hecmw_solve_send_recv
126 
127  !C
128  !C*** SOLVER_ISEND_IRECV
129  !C
131  & ( n, m, neibpetot, neibpe, stack_import, nod_import, &
132  & stack_export, nod_export, &
133  & x, solver_comm,my_rank,ireq)
134  implicit none
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
146 
147 #ifndef HECMW_SERIAL
148  ! local valiables
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
154  !C
155  !C-- INIT.
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))
162  !C
163  !C-- SEND
164  nreq1=0
165  do neib= 1, neibpetot
166  istart= stack_export(neib-1)
167  inum = stack_export(neib ) - istart
168  if (inum==0) cycle
169  nreq1=nreq1+1
170  do k= istart+1, istart+inum
171  ii = nod_export(k)
172  do kk = 1, m
173  ws(m*(k-1)+kk)= x(m*(ii-1)+kk)
174  enddo
175  enddo
176  call mpi_isend (ws(m*istart+1), m*inum,mpi_double_precision, &
177  & neibpe(neib), 0, solver_comm, req1(nreq1), ierr)
178  enddo
179  !C
180  !C-- RECEIVE
181  nreq2=0
182  do neib= 1, neibpetot
183  istart= stack_import(neib-1)
184  inum = stack_import(neib ) - istart
185  if (inum==0) cycle
186  nreq2=nreq2+1
187  call mpi_irecv (wr(m*istart+1), m*inum, mpi_double_precision, &
188  & neibpe(neib), 0, solver_comm, req2(nreq2), ierr)
189  enddo
190  !C
191  !C-- Find empty abuf
192  ireq = 0
193  do i = 1, max_nreq
194  if (abuf(i)%NEIBPETOT == 0) then
195  ireq = i
196  exit
197  endif
198  enddo
199  if (ireq == 0) then
200  stop 'Error: HECMW_SOLVE_ISEND_IRECV: exceeded maximum num of requests'
201  endif
202  !C
203  !C-- Save in abuf
204  abuf(ireq)%NEIBPETOT = neibpetot
205  abuf(ireq)%STACK_IMPORT=> stack_import
206  abuf(ireq)%NOD_IMPORT => nod_import
207  abuf(ireq)%WS => ws
208  abuf(ireq)%WR => wr
209  abuf(ireq)%X => x
210  abuf(ireq)%req1 => req1
211  abuf(ireq)%req2 => req2
212  abuf(ireq)%nreq1 = nreq1
213  abuf(ireq)%nreq2 = nreq2
214 #endif
215  end subroutine hecmw_solve_isend_irecv
216 
217  !C
218  !C*** SOLVER_ISEND_IRECV_WAIT
219  !C
220  subroutine hecmw_solve_isend_irecv_wait( m, ireq )
221  use hecmw_util
222  implicit none
223  integer(kind=kint ), intent(in) :: m, ireq
224 
225 #ifndef HECMW_SERIAL
226  ! local valiables
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
238  !C-- Check ireq
239  if (ireq < 0 .or. ireq > max_nreq) then
240  stop 'ERROR: HECMW_SOLVE_ISEND_IRECV_WAIT: invalid ireq'
241  endif
242  !C-- Restore from abuf
243  neibpetot = abuf(ireq)%NEIBPETOT
244  stack_import=> abuf(ireq)%STACK_IMPORT
245  nod_import => abuf(ireq)%NOD_IMPORT
246  ws => abuf(ireq)%WS
247  wr => abuf(ireq)%WR
248  x => abuf(ireq)%X
249  req1 => abuf(ireq)%req1
250  req2 => abuf(ireq)%req2
251  nreq1 = abuf(ireq)%nreq1
252  nreq2 = abuf(ireq)%nreq2
253  !C-- Free abuf
254  abuf(ireq)%NEIBPETOT = 0
255  !C-- Allocate
256  allocate (sta1(mpi_status_size,neibpetot))
257  allocate (sta2(mpi_status_size,neibpetot))
258  !C-- Wait irecv
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
264  ii = nod_import(k)
265  do j = 1, m
266  x(m*(ii-1)+j)= wr(m*(k-1)+j)
267  enddo
268  enddo
269  enddo
270  !C-- Wait isend
271  call mpi_waitall (nreq1, req1, sta1, ierr)
272  !C-- Deallocate
273  deallocate (sta1, sta2)
274  deallocate (req1, req2)
275  deallocate (ws, wr)
276 #endif
277  end subroutine hecmw_solve_isend_irecv_wait
278 
279  !C
280  !C*** SOLVER_REVERSE_SEND_RECV
281  !C
283  & ( n, m, neibpetot, neibpe, stack_import, nod_import, &
284  & stack_export, nod_export, &
285  & ws, wr, x, solver_comm,my_rank)
286 
287  implicit none
288 
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
301 
302 #ifndef HECMW_SERIAL
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
307 
308  ! local valiables
309  integer(kind=kint ) :: neib,istart,inum,k,kk,ii,ierr,nreq1,nreq2
310  !C
311  !C-- INIT.
312  allocate (sta1(mpi_status_size,neibpetot))
313  allocate (sta2(mpi_status_size,neibpetot))
314  allocate (req1(neibpetot))
315  allocate (req2(neibpetot))
316 
317  !C
318  !C-- SEND
319  nreq1=0
320  do neib= 1, neibpetot
321  istart= stack_import(neib-1)
322  inum = stack_import(neib ) - istart
323  if (inum==0) cycle
324  nreq1=nreq1+1
325  do k= istart+1, istart+inum
326  ii = nod_import(k)
327  do kk = 1, m
328  ws(m*(k-1)+kk)= x(m*(ii-1)+kk)
329  enddo
330  enddo
331 
332  call mpi_isend (ws(m*istart+1), m*inum,mpi_double_precision, &
333  & neibpe(neib), 0, solver_comm, req1(nreq1), ierr)
334  enddo
335 
336  !C
337  !C-- RECEIVE
338  nreq2=0
339  do neib= 1, neibpetot
340  istart= stack_export(neib-1)
341  inum = stack_export(neib ) - istart
342  if (inum==0) cycle
343  nreq2=nreq2+1
344  call mpi_irecv (wr(m*istart+1), m*inum, mpi_double_precision, &
345  & neibpe(neib), 0, solver_comm, req2(nreq2), ierr)
346  enddo
347 
348  call mpi_waitall (nreq2, req2, sta2, ierr)
349 
350  do neib= 1, neibpetot
351  istart= stack_export(neib-1)
352  inum = stack_export(neib ) - istart
353  do k= istart+1, istart+inum
354  ii = nod_export(k)
355  do kk = 1, m
356  x(m*(ii-1)+kk)= x(m*(ii-1)+kk)+wr(m*(k-1)+kk)
357  enddo
358  enddo
359  enddo
360 
361  call mpi_waitall (nreq1, req1, sta1, ierr)
362  deallocate (sta1, sta2, req1, req2)
363 #endif
364  end subroutine hecmw_solve_rev_send_recv
365 
366 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