FrontISTR  5.2.0
Large-scale structural analysis program with finit element method
hecmw_solver_SR_11i.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*** module hecmw_solver_SR_11i
8 !C***
9 !
11 contains
12  !C
13  !C*** SOLVER_SEND_RECV
14  !C
16  & ( n, neibpetot, neibpe, stack_import, nod_import, &
17  & stack_export, nod_export, &
18  & ws, wr, x, solver_comm,my_rank)
19 
20  use hecmw_util
21  implicit real*8 (a-h,o-z)
22  ! include 'mpif.h'
23  ! include 'hecmw_config_f.h'
24 
25  integer(kind=kint ) , intent(in) :: N
26  integer(kind=kint ) , intent(in) :: NEIBPETOT
27  integer(kind=kint ), pointer :: NEIBPE (:)
28  integer(kind=kint ), pointer :: STACK_IMPORT(:)
29  integer(kind=kint ), pointer :: NOD_IMPORT (:)
30  integer(kind=kint ), pointer :: STACK_EXPORT(:)
31  integer(kind=kint ), pointer :: NOD_EXPORT (:)
32  integer(kind=kint ), dimension(:) , intent(inout):: WS
33  integer(kind=kint ), dimension(:) , intent(inout):: WR
34  integer(kind=kint ), dimension(:) , intent(inout):: X
35  integer(kind=kint ) , intent(in) ::SOLVER_COMM
36  integer(kind=kint ) , intent(in) :: my_rank
37 
38 #ifndef HECMW_SERIAL
39  integer(kind=kint ), dimension(:,:), allocatable :: sta1
40  integer(kind=kint ), dimension(:,:), allocatable :: sta2
41  integer(kind=kint ), dimension(: ), allocatable :: req1
42  integer(kind=kint ), dimension(: ), allocatable :: req2
43 
44  integer(kind=kint ), save :: NFLAG
45  data nflag/0/
46  ! local valiables
47  integer(kind=kint ) :: neib, istart,inum,k,ierr,nreq1,nreq2
48 
49  !C
50  !C-- INIT.
51  allocate (sta1(mpi_status_size,neibpetot))
52  allocate (sta2(mpi_status_size,neibpetot))
53  allocate (req1(neibpetot))
54  allocate (req2(neibpetot))
55 
56  !C
57  !C-- SEND
58 
59  nreq1=0
60  do neib= 1, neibpetot
61  istart= stack_export(neib-1)
62  inum = stack_export(neib ) - istart
63  if (inum==0) cycle
64  nreq1=nreq1+1
65 
66  do k= istart+1, istart+inum
67  ws(k)= x(nod_export(k))
68  enddo
69  call mpi_isend (ws(istart+1), inum, mpi_integer, &
70  & neibpe(neib), 0, solver_comm, &
71  & req1(nreq1), ierr)
72  enddo
73 
74  !C
75  !C-- RECEIVE
76 
77  nreq2=0
78  do neib= 1, neibpetot
79  istart= stack_import(neib-1)
80  inum = stack_import(neib ) - istart
81  if (inum==0) cycle
82  nreq2=nreq2+1
83  call mpi_irecv (wr(istart+1), inum, mpi_integer, &
84  & neibpe(neib), 0, solver_comm, &
85  & req2(nreq2), ierr)
86  enddo
87 
88  call mpi_waitall (nreq2, req2, sta2, ierr)
89 
90  do neib= 1, neibpetot
91  istart= stack_import(neib-1)
92  inum = stack_import(neib ) - istart
93  do k= istart+1, istart+inum
94  x(nod_import(k))= wr(k)
95  enddo
96  enddo
97 
98  call mpi_waitall (nreq1, req1, sta1, ierr)
99 
100  deallocate (sta1, sta2, req1, req2)
101 #endif
102  end subroutine hecmw_solve_send_recv_11i
103 end module hecmw_solver_sr_11i
subroutine hecmw_solve_send_recv_11i(N, 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