FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_comm_f.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 
7 contains
8 
9  !C
10  !C***
11  !C*** hecmw_barrier
12  !C***
13  !C
14  subroutine hecmw_barrier (hecMESH)
16  implicit none
17  type (hecmwST_local_mesh) :: hecMESH
18 #ifndef HECMW_SERIAL
19  integer(kind=kint):: ierr
20  call mpi_barrier (hecmesh%MPI_COMM, ierr)
21 #endif
22  end subroutine hecmw_barrier
23 
24  subroutine hecmw_scatterv_dp(sbuf, sc, disp, rbuf, rc, root, comm)
26  implicit none
27  integer(kind=kint) :: sc !send counts
28  double precision :: sbuf(sc) !send buffer
29  integer(kind=kint) :: disp !displacement
30  integer(kind=kint) :: rc !receive counts
31  double precision :: rbuf(rc) !receive buffer
32  integer(kind=kint) :: root
33  integer(kind=kint) :: comm
34 #ifndef HECMW_SERIAL
35  integer(kind=kint) :: ierr
36  call mpi_scatterv( sbuf, sc, disp, mpi_double_precision, &
37  rbuf, rc, mpi_double_precision, &
38  root, comm, ierr )
39 #else
40  rbuf(1) = sbuf(1)
41 #endif
42  end subroutine hecmw_scatterv_dp
43 
44 #ifndef HECMW_SERIAL
45  function hecmw_operation_hec2mpi(operation)
47  implicit none
48  integer(kind=kint) :: hecmw_operation_hec2mpi
49  integer(kind=kint) :: operation
51  if (operation == hecmw_sum) then
52  hecmw_operation_hec2mpi = mpi_sum
53  elseif (operation == hecmw_prod) then
54  hecmw_operation_hec2mpi = mpi_prod
55  elseif (operation == hecmw_max) then
56  hecmw_operation_hec2mpi = mpi_max
57  elseif (operation == hecmw_min) then
58  hecmw_operation_hec2mpi = mpi_min
59  endif
60  end function hecmw_operation_hec2mpi
61 #endif
62 
63  subroutine hecmw_scatter_int_1(sbuf, rval, root, comm)
65  implicit none
66  integer(kind=kint) :: sbuf(*) !send buffer
67  integer(kind=kint) :: rval !receive value
68  integer(kind=kint) :: root
69  integer(kind=kint) :: comm
70 #ifndef HECMW_SERIAL
71  integer(kind=kint) :: ierr
72  call mpi_scatter( sbuf, 1, mpi_integer, &
73  & rval, 1, mpi_integer, root, comm, ierr )
74 #else
75  rval=sbuf(1)
76 #endif
77  end subroutine hecmw_scatter_int_1
78 
79  subroutine hecmw_gather_int_1(sval, rbuf, root, comm)
81  implicit none
82  integer(kind=kint) :: sval !send buffer
83  integer(kind=kint) :: rbuf(*) !receive buffer
84  integer(kind=kint) :: root
85  integer(kind=kint) :: comm
86 #ifndef HECMW_SERIAL
87  integer(kind=kint) :: ierr
88  call mpi_gather( sval, 1, mpi_integer, &
89  & rbuf, 1, mpi_integer, root, comm, ierr )
90 #else
91  rbuf(1)=sval
92 #endif
93  end subroutine hecmw_gather_int_1
94 
95  subroutine hecmw_allgather_int_1(sval, rbuf, comm)
97  implicit none
98  integer(kind=kint) :: sval !send buffer
99  integer(kind=kint) :: rbuf(*) !receive buffer
100  integer(kind=kint) :: comm
101 #ifndef HECMW_SERIAL
102  integer(kind=kint) :: ierr
103  call mpi_allgather( sval, 1, mpi_integer, &
104  & rbuf, 1, mpi_integer, comm, ierr )
105 #else
106  rbuf(1)=sval
107 #endif
108  end subroutine hecmw_allgather_int_1
109 
110  subroutine hecmw_scatterv_real(sbuf, scs, disp, &
111  & rbuf, rc, root, comm)
113  implicit none
114  real(kind=kreal) :: sbuf(*) !send buffer
115  integer(kind=kint) :: scs(*) !send counts
116  integer(kind=kint) :: disp(*) !displacement
117  real(kind=kreal) :: rbuf(*) !receive buffer
118  integer(kind=kint) :: rc !receive counts
119  integer(kind=kint) :: root
120  integer(kind=kint) :: comm
121 #ifndef HECMW_SERIAL
122  integer(kind=kint) :: ierr
123  call mpi_scatterv( sbuf, scs, disp, mpi_real8, &
124  & rbuf, rc, mpi_real8, root, comm, ierr )
125 #else
126  rbuf(1:rc)=sbuf(1:rc)
127 #endif
128  end subroutine hecmw_scatterv_real
129 
130  subroutine hecmw_gatherv_real(sbuf, sc, &
131  & rbuf, rcs, disp, root, comm)
133  implicit none
134  real(kind=kreal) :: sbuf(*) !send buffer
135  integer(kind=kint) :: sc !send counts
136  real(kind=kreal) :: rbuf(*) !receive buffer
137  integer(kind=kint) :: rcs(*) !receive counts
138  integer(kind=kint) :: disp(*) !displacement
139  integer(kind=kint) :: root
140  integer(kind=kint) :: comm
141 #ifndef HECMW_SERIAL
142  integer(kind=kint) :: ierr
143  call mpi_gatherv( sbuf, sc, mpi_real8, &
144  & rbuf, rcs, disp, mpi_real8, root, comm, ierr )
145 #else
146  rbuf(1:sc)=sbuf(1:sc)
147 #endif
148  end subroutine hecmw_gatherv_real
149 
150  subroutine hecmw_gatherv_int(sbuf, sc, &
151  & rbuf, rcs, disp, root, comm)
153  implicit none
154  integer(kind=kint) :: sbuf(*) !send buffer
155  integer(kind=kint) :: sc !send counts
156  integer(kind=kint) :: rbuf(*) !receive buffer
157  integer(kind=kint) :: rcs(*) !receive counts
158  integer(kind=kint) :: disp(*) !displacement
159  integer(kind=kint) :: root
160  integer(kind=kint) :: comm
161 #ifndef HECMW_SERIAL
162  integer(kind=kint) :: ierr
163  call mpi_gatherv( sbuf, sc, mpi_integer, &
164  & rbuf, rcs, disp, mpi_integer, root, comm, ierr )
165 #else
166  rbuf(1:sc)=sbuf(1:sc)
167 #endif
168  end subroutine hecmw_gatherv_int
169 
170  subroutine hecmw_allreduce_int_1(sval, rval, op, comm)
172  implicit none
173  integer(kind=kint) :: sval !send val
174  integer(kind=kint) :: rval !receive val
175  integer(kind=kint):: op, comm, ierr
176 #ifndef HECMW_SERIAL
177  call mpi_allreduce(sval, rval, 1, mpi_integer, &
178  & hecmw_operation_hec2mpi(op), comm, ierr)
179 #else
180  rval=sval
181 #endif
182  end subroutine hecmw_allreduce_int_1
183 
184  subroutine hecmw_alltoall_int(sbuf, sc, rbuf, rc, comm)
186  implicit none
187  integer(kind=kint) :: sbuf(*)
188  integer(kind=kint) :: sc
189  integer(kind=kint) :: rbuf(*)
190  integer(kind=kint) :: rc
191  integer(kind=kint) :: comm
192 #ifndef HECMW_SERIAL
193  integer(kind=kint) :: ierr
194  call mpi_alltoall( sbuf, sc, mpi_integer, &
195  & rbuf, rc, mpi_integer, comm, ierr )
196 #else
197  rbuf(1:sc)=sbuf(1:sc)
198 #endif
199  end subroutine hecmw_alltoall_int
200 
201  subroutine hecmw_isend_int(sbuf, sc, dest, &
202  & tag, comm, req)
204  implicit none
205  integer(kind=kint) :: sbuf(*)
206  integer(kind=kint) :: sc
207  integer(kind=kint) :: dest
208  integer(kind=kint) :: tag
209  integer(kind=kint) :: comm
210  integer(kind=kint) :: req
211 #ifndef HECMW_SERIAL
212  integer(kind=kint) :: ierr
213  call mpi_isend(sbuf, sc, mpi_integer, &
214  & dest, tag, comm, req, ierr)
215 #endif
216  end subroutine hecmw_isend_int
217 
218  subroutine hecmw_isend_r(sbuf, sc, dest, &
219  & tag, comm, req)
221  implicit none
222  integer(kind=kint) :: sc
223  double precision, dimension(sc) :: sbuf
224  integer(kind=kint) :: dest
225  integer(kind=kint) :: tag
226  integer(kind=kint) :: comm
227  integer(kind=kint) :: req
228 #ifndef HECMW_SERIAL
229  integer(kind=kint) :: ierr
230  call mpi_isend(sbuf, sc, mpi_double_precision, &
231  & dest, tag, comm, req, ierr)
232 #endif
233  end subroutine hecmw_isend_r
234 
235  subroutine hecmw_irecv_int(rbuf, rc, source, &
236  & tag, comm, req)
238  implicit none
239  integer(kind=kint) :: rbuf(*)
240  integer(kind=kint) :: rc
241  integer(kind=kint) :: source
242  integer(kind=kint) :: tag
243  integer(kind=kint) :: comm
244  integer(kind=kint) :: req
245 #ifndef HECMW_SERIAL
246  integer(kind=kint) :: ierr
247  call mpi_irecv(rbuf, rc, mpi_integer, &
248  & source, tag, comm, req, ierr)
249 #endif
250  end subroutine hecmw_irecv_int
251 
252  subroutine hecmw_irecv_r(rbuf, rc, source, &
253  & tag, comm, req)
255  implicit none
256  integer(kind=kint) :: rc
257  double precision, dimension(rc) :: rbuf
258  integer(kind=kint) :: source
259  integer(kind=kint) :: tag
260  integer(kind=kint) :: comm
261  integer(kind=kint) :: req
262 #ifndef HECMW_SERIAL
263  integer(kind=kint) :: ierr
264  call mpi_irecv(rbuf, rc, mpi_double_precision, &
265  & source, tag, comm, req, ierr)
266 #endif
267  end subroutine hecmw_irecv_r
268 
269  subroutine hecmw_waitall(cnt, reqs, stats)
271  implicit none
272  integer(kind=kint) :: cnt
273  integer(kind=kint) :: reqs(*)
274  integer(kind=kint) :: stats(HECMW_STATUS_SIZE,*)
275 #ifndef HECMW_SERIAL
276  integer(kind=kint) :: ierr
277  call mpi_waitall(cnt, reqs, stats, ierr)
278 #endif
279  end subroutine hecmw_waitall
280 
281  subroutine hecmw_recv_int(rbuf, rc, source, &
282  & tag, comm, stat)
284  implicit none
285  integer(kind=kint) :: rbuf(*)
286  integer(kind=kint) :: rc
287  integer(kind=kint) :: source
288  integer(kind=kint) :: tag
289  integer(kind=kint) :: comm
290  integer(kind=kint) :: stat(HECMW_STATUS_SIZE)
291 #ifndef HECMW_SERIAL
292  integer(kind=kint) :: ierr
293  call mpi_recv(rbuf, rc, mpi_integer, &
294  & source, tag, comm, stat, ierr)
295 #endif
296  end subroutine hecmw_recv_int
297 
298  subroutine hecmw_recv_r(rbuf, rc, source, &
299  & tag, comm, stat)
301  implicit none
302  integer(kind=kint) :: rc
303  double precision, dimension(rc) :: rbuf
304  integer(kind=kint) :: source
305  integer(kind=kint) :: tag
306  integer(kind=kint) :: comm
307  integer(kind=kint) :: stat(HECMW_STATUS_SIZE)
308 #ifndef HECMW_SERIAL
309  integer(kind=kint) :: ierr
310  call mpi_recv(rbuf, rc, mpi_double_precision, &
311  & source, tag, comm, stat, ierr)
312 #endif
313  end subroutine hecmw_recv_r
314  !C
315  !C***
316  !C*** hecmw_allREDUCE
317  !C***
318  !C
319  subroutine hecmw_allreduce_dp(val,VALM,n,hec_op,comm )
321  implicit none
322  integer(kind=kint) :: n, hec_op,op, comm, ierr
323  double precision, dimension(n) :: val
324  double precision, dimension(n) :: VALM
325 #ifndef HECMW_SERIAL
326  select case( hec_op )
327  case ( hecmw_sum )
328  op = mpi_sum
329  case ( hecmw_prod )
330  op = mpi_prod
331  case ( hecmw_max )
332  op = mpi_max
333  case ( hecmw_min )
334  op = mpi_min
335  end select
336  call mpi_allreduce(val,valm,n,mpi_double_precision,op,comm,ierr)
337 #else
338  integer(kind=kint) :: i
339  do i=1,n
340  valm(i) = val(i)
341  end do
342 #endif
343  end subroutine hecmw_allreduce_dp
344 
345  subroutine hecmw_allreduce_dp1(s1,s2,hec_op,comm )
347  implicit none
348  integer(kind=kint) :: hec_op, comm
349  double precision :: s1, s2
350  double precision, dimension(1) :: val
351  double precision, dimension(1) :: VALM
352 #ifndef HECMW_SERIAL
353  val(1) = s1
354  valm(1) = s2
355  call hecmw_allreduce_dp(val,valm,1,hec_op,comm )
356  s1 = val(1)
357  s2 = valm(1)
358 #else
359  s2 = s1
360 #endif
361  end subroutine hecmw_allreduce_dp1
362  !C
363  !C***
364  !C*** hecmw_allREDUCE_R
365  !C***
366  !C
367  subroutine hecmw_allreduce_r (hecMESH, val, n, ntag)
369  implicit none
370  integer(kind=kint):: n, ntag
371  real(kind=kreal), dimension(n) :: val
372  type (hecmwST_local_mesh) :: hecMESH
373 #ifndef HECMW_SERIAL
374  integer(kind=kint):: ierr
375  real(kind=kreal), dimension(:), allocatable :: valm
376 
377  allocate (valm(n))
378  valm= 0.d0
379  if (ntag .eq. hecmw_sum) then
380  call mpi_allreduce &
381  & (val, valm, n, mpi_double_precision, mpi_sum, &
382  & hecmesh%MPI_COMM, ierr)
383  endif
384 
385  if (ntag .eq. hecmw_max) then
386  call mpi_allreduce &
387  & (val, valm, n, mpi_double_precision, mpi_max, &
388  & hecmesh%MPI_COMM, ierr)
389  endif
390 
391  if (ntag .eq. hecmw_min) then
392  call mpi_allreduce &
393  & (val, valm, n, mpi_double_precision, mpi_min, &
394  & hecmesh%MPI_COMM, ierr)
395  endif
396 
397  val= valm
398  deallocate (valm)
399 #endif
400  end subroutine hecmw_allreduce_r
401 
402  subroutine hecmw_allreduce_r1 (hecMESH, s, ntag)
404  implicit none
405  integer(kind=kint):: ntag
406  real(kind=kreal) :: s
407  type (hecmwST_local_mesh) :: hecMESH
408 #ifndef HECMW_SERIAL
409  real(kind=kreal), dimension(1) :: val
410  val(1) = s
411  call hecmw_allreduce_r(hecmesh, val, 1, ntag )
412  s = val(1)
413 #endif
414  end subroutine hecmw_allreduce_r1
415 
416  !C
417  !C***
418  !C*** hecmw_allREDUCE_I
419  !C***
420  !C
421  subroutine hecmw_allreduce_i(hecMESH, val, n, ntag)
423  implicit none
424  integer(kind=kint):: n, ntag
425  integer(kind=kint), dimension(n) :: val
426  type (hecmwST_local_mesh) :: hecMESH
427 #ifndef HECMW_SERIAL
428  integer(kind=kint):: ierr
429  integer(kind=kint), dimension(:), allocatable :: VALM
430 
431  allocate (valm(n))
432  valm= 0
433  if (ntag .eq. hecmw_sum) then
434  call mpi_allreduce &
435  & (val, valm, n, mpi_integer, mpi_sum, &
436  & hecmesh%MPI_COMM, ierr)
437  endif
438 
439  if (ntag .eq. hecmw_max) then
440  call mpi_allreduce &
441  & (val, valm, n, mpi_integer, mpi_max, &
442  & hecmesh%MPI_COMM, ierr)
443  endif
444 
445  if (ntag .eq. hecmw_min) then
446  call mpi_allreduce &
447  & (val, valm, n, mpi_integer, mpi_min, &
448  & hecmesh%MPI_COMM, ierr)
449  endif
450 
451 
452  val= valm
453  deallocate (valm)
454 #endif
455  end subroutine hecmw_allreduce_i
456 
457  subroutine hecmw_allreduce_i1 (hecMESH, s, ntag)
459  implicit none
460  integer(kind=kint):: ntag, s
461  type (hecmwST_local_mesh) :: hecMESH
462 #ifndef HECMW_SERIAL
463  integer(kind=kint), dimension(1) :: val
464 
465  val(1) = s
466  call hecmw_allreduce_i(hecmesh, val, 1, ntag )
467  s = val(1)
468 #endif
469  end subroutine hecmw_allreduce_i1
470 
471  !C
472  !C***
473  !C*** hecmw_bcast_R
474  !C***
475  !C
476  subroutine hecmw_bcast_r (hecMESH, val, n, nbase)
478  implicit none
479  integer(kind=kint):: n, nbase
480  real(kind=kreal), dimension(n) :: val
481  type (hecmwST_local_mesh) :: hecMESH
482 #ifndef HECMW_SERIAL
483  integer(kind=kint):: ierr
484  call mpi_bcast (val, n, mpi_double_precision, nbase, hecmesh%MPI_COMM, ierr)
485 #endif
486  end subroutine hecmw_bcast_r
487 
488  subroutine hecmw_bcast_r_comm (val, n, nbase, comm)
490  implicit none
491  integer(kind=kint):: n, nbase
492  real(kind=kreal), dimension(n) :: val
493  integer(kind=kint):: comm
494 #ifndef HECMW_SERIAL
495  integer(kind=kint):: ierr
496  call mpi_bcast (val, n, mpi_double_precision, nbase, comm, ierr)
497 #endif
498  end subroutine hecmw_bcast_r_comm
499 
500  subroutine hecmw_bcast_r1 (hecMESH, s, nbase)
502  implicit none
503  integer(kind=kint):: nbase, ierr
504  real(kind=kreal) :: s
505  type (hecmwST_local_mesh) :: hecMESH
506 #ifndef HECMW_SERIAL
507  real(kind=kreal), dimension(1) :: val
508  val(1)=s
509  call mpi_bcast (val, 1, mpi_double_precision, nbase, hecmesh%MPI_COMM, ierr)
510  s = val(1)
511 #endif
512  end subroutine hecmw_bcast_r1
513 
514  subroutine hecmw_bcast_r1_comm (s, nbase, comm)
516  implicit none
517  integer(kind=kint):: nbase
518  real(kind=kreal) :: s
519  integer(kind=kint):: comm
520 #ifndef HECMW_SERIAL
521  integer(kind=kint):: ierr
522  real(kind=kreal), dimension(1) :: val
523  val(1)=s
524  call mpi_bcast (val, 1, mpi_double_precision, nbase, comm, ierr)
525  s = val(1)
526 #endif
527  end subroutine hecmw_bcast_r1_comm
528  !C
529  !C***
530  !C*** hecmw_bcast_I
531  !C***
532  !C
533  subroutine hecmw_bcast_i (hecMESH, val, n, nbase)
535  implicit none
536  integer(kind=kint):: n, nbase
537  integer(kind=kint), dimension(n) :: val
538  type (hecmwST_local_mesh) :: hecMESH
539 #ifndef HECMW_SERIAL
540  integer(kind=kint):: ierr
541  call mpi_bcast (val, n, mpi_integer, nbase, hecmesh%MPI_COMM, ierr)
542 #endif
543  end subroutine hecmw_bcast_i
544 
545  subroutine hecmw_bcast_i_comm (val, n, nbase, comm)
547  implicit none
548  integer(kind=kint):: n, nbase
549  integer(kind=kint), dimension(n) :: val
550  integer(kind=kint):: comm
551 #ifndef HECMW_SERIAL
552  integer(kind=kint):: ierr
553  call mpi_bcast (val, n, mpi_integer, nbase, comm, ierr)
554 #endif
555  end subroutine hecmw_bcast_i_comm
556 
557  subroutine hecmw_bcast_i1 (hecMESH, s, nbase)
559  implicit none
560  integer(kind=kint):: nbase, s
561  type (hecmwST_local_mesh) :: hecMESH
562 #ifndef HECMW_SERIAL
563  integer(kind=kint):: ierr
564  integer(kind=kint), dimension(1) :: val
565  val(1) = s
566  call mpi_bcast (val, 1, mpi_integer, nbase, hecmesh%MPI_COMM, ierr)
567  s = val(1)
568 #endif
569  end subroutine hecmw_bcast_i1
570 
571  subroutine hecmw_bcast_i1_comm (s, nbase, comm)
573  implicit none
574  integer(kind=kint):: nbase, s
575  integer(kind=kint):: comm
576 #ifndef HECMW_SERIAL
577  integer(kind=kint):: ierr
578  integer(kind=kint), dimension(1) :: val
579  val(1) = s
580  call mpi_bcast (val, 1, mpi_integer, nbase, comm, ierr)
581  s = val(1)
582 #endif
583  end subroutine hecmw_bcast_i1_comm
584  !C
585  !C***
586  !C*** hecmw_bcast_C
587  !C***
588  !C
589  subroutine hecmw_bcast_c (hecMESH, val, n, nn, nbase)
591  implicit none
592  integer(kind=kint):: n, nn, nbase
593  character(len=n) :: val(nn)
594  type (hecmwST_local_mesh) :: hecMESH
595 #ifndef HECMW_SERIAL
596  integer(kind=kint):: ierr
597  call mpi_bcast (val, n*nn, mpi_character, nbase, hecmesh%MPI_COMM,&
598  & ierr)
599 #endif
600  end subroutine hecmw_bcast_c
601 
602  subroutine hecmw_bcast_c_comm (val, n, nn, nbase, comm)
604  implicit none
605  integer(kind=kint):: n, nn, nbase
606  character(len=n) :: val(nn)
607  integer(kind=kint):: comm
608 #ifndef HECMW_SERIAL
609  integer(kind=kint):: ierr
610  call mpi_bcast (val, n*nn, mpi_character, nbase, comm,&
611  & ierr)
612 #endif
613  end subroutine hecmw_bcast_c_comm
614 
615  !C
616  !C***
617  !C*** hecmw_assemble_R
618  !C***
619  subroutine hecmw_assemble_r (hecMESH, val, n, m)
621  use hecmw_solver_sr
622 
623  implicit none
624  integer(kind=kint):: n, m
625  real(kind=kreal), dimension(m*n) :: val
626  type (hecmwST_local_mesh) :: hecMESH
627 #ifndef HECMW_SERIAL
628  integer(kind=kint):: ns, nr
629  real(kind=kreal), dimension(:), allocatable :: ws, wr
630 
631  if( hecmesh%n_neighbor_pe == 0 ) return
632 
633  ns = hecmesh%import_index(hecmesh%n_neighbor_pe)
634  nr = hecmesh%export_index(hecmesh%n_neighbor_pe)
635 
636  allocate (ws(m*ns), wr(m*nr))
638  & ( n, m, hecmesh%n_neighbor_pe, hecmesh%neighbor_pe, &
639  & hecmesh%import_index, hecmesh%import_item, &
640  & hecmesh%export_index, hecmesh%export_item, &
641  & ws, wr, val , hecmesh%MPI_COMM, hecmesh%my_rank)
642  deallocate (ws, wr)
643 #endif
644  end subroutine hecmw_assemble_r
645 
646  !C
647  !C***
648  !C*** hecmw_assemble_I
649  !C***
650  subroutine hecmw_assemble_i (hecMESH, val, n, m)
653 
654  implicit none
655  integer(kind=kint):: n, m
656  integer(kind=kint), dimension(m*n) :: val
657  type (hecmwST_local_mesh) :: hecMESH
658 #ifndef HECMW_SERIAL
659  integer(kind=kint):: ns, nr
660  integer(kind=kint), dimension(:), allocatable :: WS, WR
661 
662  if( hecmesh%n_neighbor_pe == 0 ) return
663 
664  ns = hecmesh%import_index(hecmesh%n_neighbor_pe)
665  nr = hecmesh%export_index(hecmesh%n_neighbor_pe)
666 
667  allocate (ws(m*ns), wr(m*nr))
669  & ( n, m, hecmesh%n_neighbor_pe, hecmesh%neighbor_pe, &
670  & hecmesh%import_index, hecmesh%import_item, &
671  & hecmesh%export_index, hecmesh%export_item, &
672  & ws, wr, val , hecmesh%MPI_COMM, hecmesh%my_rank)
673  deallocate (ws, wr)
674 #endif
675  end subroutine hecmw_assemble_i
676 
677  !C
678  !C***
679  !C*** hecmw_update_R
680  !C***
681  !C
682  subroutine hecmw_update_r (hecMESH, val, n, m)
684  use hecmw_solver_sr
685 
686  implicit none
687  integer(kind=kint):: n, m
688  real(kind=kreal), dimension(m*n) :: val
689  type (hecmwST_local_mesh) :: hecMESH
690 #ifndef HECMW_SERIAL
691  integer(kind=kint):: ns, nr
692  real(kind=kreal), dimension(:), allocatable :: ws, wr
693 
694  if( hecmesh%n_neighbor_pe == 0 ) return
695 
696  ns = hecmesh%export_index(hecmesh%n_neighbor_pe)
697  nr = hecmesh%import_index(hecmesh%n_neighbor_pe)
698 
699  allocate (ws(m*ns), wr(m*nr))
700  call hecmw_solve_send_recv &
701  & ( n, m, hecmesh%n_neighbor_pe, hecmesh%neighbor_pe, &
702  & hecmesh%import_index, hecmesh%import_item, &
703  & hecmesh%export_index, hecmesh%export_item, &
704  & ws, wr, val , hecmesh%MPI_COMM, hecmesh%my_rank)
705  deallocate (ws, wr)
706 #endif
707  end subroutine hecmw_update_r
708 
709  !C
710  !C***
711  !C*** hecmw_update_R_async
712  !C***
713  !C
714  !C REAL
715  !C
716  subroutine hecmw_update_r_async (hecMESH, val, n, m, ireq)
718  use hecmw_solver_sr
719  implicit none
720  integer(kind=kint):: n, m, ireq
721  real(kind=kreal), dimension(m*n) :: val
722  type (hecmwST_local_mesh) :: hecMESH
723 #ifndef HECMW_SERIAL
724  if( hecmesh%n_neighbor_pe == 0 ) return
725 
727  & ( n, m, hecmesh%n_neighbor_pe, hecmesh%neighbor_pe, &
728  & hecmesh%import_index, hecmesh%import_item, &
729  & hecmesh%export_index, hecmesh%export_item, &
730  & val , hecmesh%MPI_COMM, hecmesh%my_rank, ireq)
731 #endif
732  end subroutine hecmw_update_r_async
733 
734  !C
735  !C***
736  !C*** hecmw_update_R_wait
737  !C***
738  !C
739  !C REAL
740  !C
741  subroutine hecmw_update_r_wait (hecMESH, ireq)
743  use hecmw_solver_sr
744  implicit none
745  integer(kind=kint):: ireq
746  type (hecmwST_local_mesh) :: hecMESH
747 #ifndef HECMW_SERIAL
748  if( hecmesh%n_neighbor_pe == 0 ) return
749 
750  call hecmw_solve_isend_irecv_wait(hecmesh%n_dof, ireq )
751 #endif
752  end subroutine hecmw_update_r_wait
753 
754 
755  !C
756  !C***
757  !C*** hecmw_update_I
758  !C***
759  !C
760  subroutine hecmw_update_i (hecMESH, val, n, m)
763 
764  implicit none
765  integer(kind=kint):: n, m
766  integer(kind=kint), dimension(m*n) :: val
767  type (hecmwST_local_mesh) :: hecMESH
768 #ifndef HECMW_SERIAL
769  integer(kind=kint):: ns, nr
770  integer(kind=kint), dimension(:), allocatable :: WS, WR
771 
772  if( hecmesh%n_neighbor_pe == 0 ) return
773 
774  ns = hecmesh%export_index(hecmesh%n_neighbor_pe)
775  nr = hecmesh%import_index(hecmesh%n_neighbor_pe)
776 
777  allocate (ws(m*ns), wr(m*nr))
779  & ( n, m, hecmesh%n_neighbor_pe, hecmesh%neighbor_pe, &
780  & hecmesh%import_index, hecmesh%import_item, &
781  & hecmesh%export_index, hecmesh%export_item, &
782  & ws, wr, val , hecmesh%MPI_COMM, hecmesh%my_rank)
783  deallocate (ws, wr)
784 #endif
785  end subroutine hecmw_update_i
786 
787 end module m_hecmw_comm_f
hecmw_util::hecmw_max
integer(kind=kint), parameter hecmw_max
Definition: hecmw_util_f.F90:25
m_hecmw_comm_f::hecmw_allreduce_i
subroutine hecmw_allreduce_i(hecMESH, val, n, ntag)
Definition: hecmw_comm_f.F90:422
m_hecmw_comm_f::hecmw_operation_hec2mpi
integer(kind=kint) function hecmw_operation_hec2mpi(operation)
Definition: hecmw_comm_f.F90:46
m_hecmw_comm_f::hecmw_allgather_int_1
subroutine hecmw_allgather_int_1(sval, rbuf, comm)
Definition: hecmw_comm_f.F90:96
m_hecmw_comm_f::hecmw_update_r_async
subroutine hecmw_update_r_async(hecMESH, val, n, m, ireq)
Definition: hecmw_comm_f.F90:717
m_hecmw_comm_f::hecmw_bcast_r
subroutine hecmw_bcast_r(hecMESH, val, n, nbase)
Definition: hecmw_comm_f.F90:477
m_hecmw_comm_f::hecmw_update_r_wait
subroutine hecmw_update_r_wait(hecMESH, ireq)
Definition: hecmw_comm_f.F90:742
m_hecmw_comm_f::hecmw_assemble_i
subroutine hecmw_assemble_i(hecMESH, val, n, m)
Definition: hecmw_comm_f.F90:651
m_hecmw_comm_f::hecmw_irecv_int
subroutine hecmw_irecv_int(rbuf, rc, source, tag, comm, req)
Definition: hecmw_comm_f.F90:237
hecmw_util::hecmw_sum
integer(kind=kint), parameter hecmw_sum
Definition: hecmw_util_f.F90:23
m_hecmw_comm_f::hecmw_irecv_r
subroutine hecmw_irecv_r(rbuf, rc, source, tag, comm, req)
Definition: hecmw_comm_f.F90:254
m_hecmw_comm_f::hecmw_scatter_int_1
subroutine hecmw_scatter_int_1(sbuf, rval, root, comm)
Definition: hecmw_comm_f.F90:64
hecmw_solver_sr
Definition: hecmw_solver_SR.F90:11
hecmw_solver_sr_i::hecmw_solve_send_recv_i
subroutine hecmw_solve_send_recv_i(N, m, NEIBPETOT, NEIBPE, STACK_IMPORT, NOD_IMPORT, STACK_EXPORT, NOD_EXPORT, WS, WR, X, SOLVER_COMM, my_rank)
Definition: hecmw_solver_SR_i.F90:21
m_hecmw_comm_f::hecmw_allreduce_r
subroutine hecmw_allreduce_r(hecMESH, val, n, ntag)
Definition: hecmw_comm_f.F90:368
m_hecmw_comm_f::hecmw_barrier
subroutine hecmw_barrier(hecMESH)
Definition: hecmw_comm_f.F90:15
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
m_hecmw_comm_f::hecmw_gatherv_real
subroutine hecmw_gatherv_real(sbuf, sc, rbuf, rcs, disp, root, comm)
Definition: hecmw_comm_f.F90:132
m_hecmw_comm_f::hecmw_allreduce_dp
subroutine hecmw_allreduce_dp(val, VALM, n, hec_op, comm)
Definition: hecmw_comm_f.F90:320
m_hecmw_comm_f::hecmw_scatterv_dp
subroutine hecmw_scatterv_dp(sbuf, sc, disp, rbuf, rc, root, comm)
Definition: hecmw_comm_f.F90:25
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
m_hecmw_comm_f::hecmw_bcast_i
subroutine hecmw_bcast_i(hecMESH, val, n, nbase)
Definition: hecmw_comm_f.F90:534
m_hecmw_comm_f::hecmw_gatherv_int
subroutine hecmw_gatherv_int(sbuf, sc, rbuf, rcs, disp, root, comm)
Definition: hecmw_comm_f.F90:152
m_hecmw_comm_f::hecmw_allreduce_r1
subroutine hecmw_allreduce_r1(hecMESH, s, ntag)
Definition: hecmw_comm_f.F90:403
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
m_hecmw_comm_f
Definition: hecmw_comm_f.F90:6
hecmw_util::hecmw_min
integer(kind=kint), parameter hecmw_min
Definition: hecmw_util_f.F90:26
m_hecmw_comm_f::hecmw_update_i
subroutine hecmw_update_i(hecMESH, val, n, m)
Definition: hecmw_comm_f.F90:761
m_hecmw_comm_f::hecmw_bcast_c
subroutine hecmw_bcast_c(hecMESH, val, n, nn, nbase)
Definition: hecmw_comm_f.F90:590
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
m_hecmw_comm_f::hecmw_allreduce_i1
subroutine hecmw_allreduce_i1(hecMESH, s, ntag)
Definition: hecmw_comm_f.F90:458
m_hecmw_comm_f::hecmw_assemble_r
subroutine hecmw_assemble_r(hecMESH, val, n, m)
Definition: hecmw_comm_f.F90:620
hecmw_solver_sr_i::hecmw_solve_rev_send_recv_i
subroutine hecmw_solve_rev_send_recv_i(N, M, NEIBPETOT, NEIBPE, STACK_IMPORT, NOD_IMPORT, STACK_EXPORT, NOD_EXPORT, WS, WR, X, SOLVER_COMM, my_rank)
Definition: hecmw_solver_SR_i.F90:114
m_hecmw_comm_f::hecmw_bcast_i1_comm
subroutine hecmw_bcast_i1_comm(s, nbase, comm)
Definition: hecmw_comm_f.F90:572
hecmw_solver_sr::hecmw_solve_isend_irecv_wait
subroutine hecmw_solve_isend_irecv_wait(m, ireq)
Definition: hecmw_solver_SR.F90:219
m_hecmw_comm_f::hecmw_alltoall_int
subroutine hecmw_alltoall_int(sbuf, sc, rbuf, rc, comm)
Definition: hecmw_comm_f.F90:185
m_hecmw_comm_f::hecmw_recv_r
subroutine hecmw_recv_r(rbuf, rc, source, tag, comm, stat)
Definition: hecmw_comm_f.F90:300
m_hecmw_comm_f::hecmw_isend_r
subroutine hecmw_isend_r(sbuf, sc, dest, tag, comm, req)
Definition: hecmw_comm_f.F90:220
hecmw_util::hecmw_prod
integer(kind=kint), parameter hecmw_prod
Definition: hecmw_util_f.F90:24
m_hecmw_comm_f::hecmw_allreduce_int_1
subroutine hecmw_allreduce_int_1(sval, rval, op, comm)
Definition: hecmw_comm_f.F90:171
m_hecmw_comm_f::hecmw_bcast_r1
subroutine hecmw_bcast_r1(hecMESH, s, nbase)
Definition: hecmw_comm_f.F90:501
m_hecmw_comm_f::hecmw_bcast_r1_comm
subroutine hecmw_bcast_r1_comm(s, nbase, comm)
Definition: hecmw_comm_f.F90:515
hecmw_solver_sr_i
Definition: hecmw_solver_SR_i.F90:11
m_hecmw_comm_f::hecmw_bcast_i_comm
subroutine hecmw_bcast_i_comm(val, n, nbase, comm)
Definition: hecmw_comm_f.F90:546
m_hecmw_comm_f::hecmw_recv_int
subroutine hecmw_recv_int(rbuf, rc, source, tag, comm, stat)
Definition: hecmw_comm_f.F90:283
m_hecmw_comm_f::hecmw_allreduce_dp1
subroutine hecmw_allreduce_dp1(s1, s2, hec_op, comm)
Definition: hecmw_comm_f.F90:346
m_hecmw_comm_f::hecmw_bcast_c_comm
subroutine hecmw_bcast_c_comm(val, n, nn, nbase, comm)
Definition: hecmw_comm_f.F90:603
m_hecmw_comm_f::hecmw_bcast_r_comm
subroutine hecmw_bcast_r_comm(val, n, nbase, comm)
Definition: hecmw_comm_f.F90:489
m_hecmw_comm_f::hecmw_waitall
subroutine hecmw_waitall(cnt, reqs, stats)
Definition: hecmw_comm_f.F90:270
m_hecmw_comm_f::hecmw_update_r
subroutine hecmw_update_r(hecMESH, val, n, m)
Definition: hecmw_comm_f.F90:683
m_hecmw_comm_f::hecmw_bcast_i1
subroutine hecmw_bcast_i1(hecMESH, s, nbase)
Definition: hecmw_comm_f.F90:558
m_hecmw_comm_f::hecmw_isend_int
subroutine hecmw_isend_int(sbuf, sc, dest, tag, comm, req)
Definition: hecmw_comm_f.F90:203
m_hecmw_comm_f::hecmw_gather_int_1
subroutine hecmw_gather_int_1(sval, rbuf, root, comm)
Definition: hecmw_comm_f.F90:80
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
m_hecmw_comm_f::hecmw_scatterv_real
subroutine hecmw_scatterv_real(sbuf, scs, disp, rbuf, rc, root, comm)
Definition: hecmw_comm_f.F90:112