FrontISTR  5.7.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
130  subroutine hecmw_solve_isend_irecv &
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)))
157  allocate (wr(m*stack_import(neibpetot)))
158  allocate (req1(neibpetot))
159  allocate (req2(neibpetot))
160  !C
161  !C-- SEND
162  nreq1=0
163  do neib= 1, neibpetot
164  istart= stack_export(neib-1)
165  inum = stack_export(neib ) - istart
166  if (inum==0) cycle
167  nreq1=nreq1+1
168  do k= istart+1, istart+inum
169  ii = nod_export(k)
170  do kk = 1, m
171  ws(m*(k-1)+kk)= x(m*(ii-1)+kk)
172  enddo
173  enddo
174  call mpi_isend (ws(m*istart+1), m*inum,mpi_double_precision, &
175  & neibpe(neib), 0, solver_comm, req1(nreq1), ierr)
176  enddo
177  !C
178  !C-- RECEIVE
179  nreq2=0
180  do neib= 1, neibpetot
181  istart= stack_import(neib-1)
182  inum = stack_import(neib ) - istart
183  if (inum==0) cycle
184  nreq2=nreq2+1
185  call mpi_irecv (wr(m*istart+1), m*inum, mpi_double_precision, &
186  & neibpe(neib), 0, solver_comm, req2(nreq2), ierr)
187  enddo
188  !C
189  !C-- Find empty abuf
190  ireq = 0
191  do i = 1, max_nreq
192  if (abuf(i)%NEIBPETOT == 0) then
193  ireq = i
194  exit
195  endif
196  enddo
197  if (ireq == 0) then
198  stop 'Error: HECMW_SOLVE_ISEND_IRECV: exceeded maximum num of requests'
199  endif
200  !C
201  !C-- Save in abuf
202  abuf(ireq)%NEIBPETOT = neibpetot
203  abuf(ireq)%STACK_IMPORT=> stack_import
204  abuf(ireq)%NOD_IMPORT => nod_import
205  abuf(ireq)%WS => ws
206  abuf(ireq)%WR => wr
207  abuf(ireq)%X => x
208  abuf(ireq)%req1 => req1
209  abuf(ireq)%req2 => req2
210  abuf(ireq)%nreq1 = nreq1
211  abuf(ireq)%nreq2 = nreq2
212 #endif
213  end subroutine hecmw_solve_isend_irecv
214 
215  !C
216  !C*** SOLVER_ISEND_IRECV_WAIT
217  !C
218  subroutine hecmw_solve_isend_irecv_wait( m, ireq )
220  implicit none
221  integer(kind=kint ), intent(in) :: m, ireq
222 
223 #ifndef HECMW_SERIAL
224  ! local valiables
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
236  !C-- Check ireq
237  if (ireq < 0 .or. ireq > max_nreq) then
238  stop 'ERROR: HECMW_SOLVE_ISEND_IRECV_WAIT: invalid ireq'
239  endif
240  !C-- Restore from abuf
241  neibpetot = abuf(ireq)%NEIBPETOT
242  stack_import=> abuf(ireq)%STACK_IMPORT
243  nod_import => abuf(ireq)%NOD_IMPORT
244  ws => abuf(ireq)%WS
245  wr => abuf(ireq)%WR
246  x => abuf(ireq)%X
247  req1 => abuf(ireq)%req1
248  req2 => abuf(ireq)%req2
249  nreq1 = abuf(ireq)%nreq1
250  nreq2 = abuf(ireq)%nreq2
251  !C-- Free abuf
252  abuf(ireq)%NEIBPETOT = 0
253  !C-- Allocate
254  allocate (sta1(mpi_status_size,neibpetot))
255  allocate (sta2(mpi_status_size,neibpetot))
256  !C-- Wait irecv
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
262  ii = nod_import(k)
263  do j = 1, m
264  x(m*(ii-1)+j)= wr(m*(k-1)+j)
265  enddo
266  enddo
267  enddo
268  !C-- Wait isend
269  call mpi_waitall (nreq1, req1, sta1, ierr)
270  !C-- Deallocate
271  deallocate (sta1, sta2)
272  deallocate (req1, req2)
273  deallocate (ws, wr)
274 #endif
275  end subroutine hecmw_solve_isend_irecv_wait
276 
277  !C
278  !C*** SOLVER_REVERSE_SEND_RECV
279  !C
280  subroutine hecmw_solve_rev_send_recv &
281  & ( n, m, neibpetot, neibpe, stack_import, nod_import, &
282  & stack_export, nod_export, &
283  & ws, wr, x, solver_comm,my_rank)
284 
285  implicit none
286 
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
299 
300 #ifndef HECMW_SERIAL
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
305 
306  ! local valiables
307  integer(kind=kint ) :: neib,istart,inum,k,kk,ii,ierr,nreq1,nreq2
308  !C
309  !C-- INIT.
310  allocate (sta1(mpi_status_size,neibpetot))
311  allocate (sta2(mpi_status_size,neibpetot))
312  allocate (req1(neibpetot))
313  allocate (req2(neibpetot))
314 
315  !C
316  !C-- SEND
317  nreq1=0
318  do neib= 1, neibpetot
319  istart= stack_import(neib-1)
320  inum = stack_import(neib ) - istart
321  if (inum==0) cycle
322  nreq1=nreq1+1
323  do k= istart+1, istart+inum
324  ii = nod_import(k)
325  do kk = 1, m
326  ws(m*(k-1)+kk)= x(m*(ii-1)+kk)
327  enddo
328  enddo
329 
330  call mpi_isend (ws(m*istart+1), m*inum,mpi_double_precision, &
331  & neibpe(neib), 0, solver_comm, req1(nreq1), ierr)
332  enddo
333 
334  !C
335  !C-- RECEIVE
336  nreq2=0
337  do neib= 1, neibpetot
338  istart= stack_export(neib-1)
339  inum = stack_export(neib ) - istart
340  if (inum==0) cycle
341  nreq2=nreq2+1
342  call mpi_irecv (wr(m*istart+1), m*inum, mpi_double_precision, &
343  & neibpe(neib), 0, solver_comm, req2(nreq2), ierr)
344  enddo
345 
346  call mpi_waitall (nreq2, req2, sta2, ierr)
347 
348  do neib= 1, neibpetot
349  istart= stack_export(neib-1)
350  inum = stack_export(neib ) - istart
351  do k= istart+1, istart+inum
352  ii = nod_export(k)
353  do kk = 1, m
354  x(m*(ii-1)+kk)= x(m*(ii-1)+kk)+wr(m*(k-1)+kk)
355  enddo
356  enddo
357  enddo
358 
359  call mpi_waitall (nreq1, req1, sta1, ierr)
360  deallocate (sta1, sta2, req1, req2)
361 #endif
362  end subroutine hecmw_solve_rev_send_recv
363 
364 end module hecmw_solver_sr
hecmw_solver_sr::max_nreq
integer(kind=kint), parameter max_nreq
Definition: hecmw_solver_SR.F90:28
hecmw_solver_sr
Definition: hecmw_solver_SR.F90:11
hecmw_solver_sr::hecmw_solve_rev_send_recv
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)
Definition: hecmw_solver_SR.F90:284
hecmw_solver_sr::abuf
type(async_buf), dimension(max_nreq), save abuf
Definition: hecmw_solver_SR.F90:29
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
hecmw_solver_sr::hecmw_solve_send_recv
subroutine hecmw_solve_send_recv(N, m, NEIBPETOT, NEIBPE, STACK_IMPORT, NOD_IMPORT, STACK_EXPORT, NOD_EXPORT, WS, WR, X, SOLVER_COMM, my_rank)
Definition: hecmw_solver_SR.F90:41
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
hecmw_solver_sr::hecmw_solve_isend_irecv_wait
subroutine hecmw_solve_isend_irecv_wait(m, ireq)
Definition: hecmw_solver_SR.F90:219
hecmw_solver_sr::async_buf
Definition: hecmw_solver_SR.F90:15
hecmw_solver_sr::hecmw_solve_isend_irecv
subroutine hecmw_solve_isend_irecv(N, M, NEIBPETOT, NEIBPE, STACK_IMPORT, NOD_IMPORT, STACK_EXPORT, NOD_EXPORT, X, SOLVER_COMM, my_rank, ireq)
Definition: hecmw_solver_SR.F90:134