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