FrontISTR  5.9.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)
15  use hecmw_util
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)
25  use hecmw_util
26  implicit none
27  double precision :: sbuf(*) !send buffer
28  integer(kind=kint) :: sc(*) !send counts
29  integer(kind=kint) :: disp(*) !displacement
30  double precision :: rbuf(*) !receive buffer
31  integer(kind=kint) :: rc !receive counts
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)
46  use hecmw_util
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)
64  use hecmw_util
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)
80  use hecmw_util
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)
96  use hecmw_util
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)
112  use hecmw_util
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)
132  use hecmw_util
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)
152  use hecmw_util
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)
171  use hecmw_util
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)
185  use hecmw_util
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)
203  use hecmw_util
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)
220  use hecmw_util
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)
237  use hecmw_util
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)
254  use hecmw_util
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)
270  use hecmw_util
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)
283  use hecmw_util
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)
300  use hecmw_util
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 )
320  use hecmw_util
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 )
346  use hecmw_util
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)
368  use hecmw_util
369 #ifdef _OPENACC
370  use openacc
371 #endif
372  implicit none
373  integer(kind=kint):: n, ntag
374  real(kind=kreal), dimension(n) :: val
375  type (hecmwST_local_mesh) :: hecMESH
376 #ifndef HECMW_SERIAL
377  integer(kind=kint):: ierr
378  real(kind=kreal), dimension(:), allocatable :: valm
379 #ifdef _OPENACC
380  integer(kind=kint) :: i
381  real(kind=kreal), dimension(n) :: tmp
382  type(c_devptr) :: tmp_dev
383 
384  tmp_dev = acc_malloc(kreal * n)
385  call acc_map_data(tmp, tmp_dev, kreal * n)
386  !$acc serial
387  do i = 1, n
388  tmp(i) = val(i)
389  enddo
390  !$acc end serial
391 #endif
392 
393  allocate (valm(n))
394  valm= 0.d0
395  if (ntag .eq. hecmw_sum) then
396 #ifdef _OPENACC
397  !$acc host_data use_device(tmp)
398  call mpi_allreduce &
399  & (tmp, valm, n, mpi_double_precision, mpi_sum, &
400  & hecmesh%MPI_COMM, ierr)
401  !$acc end host_data
402 #else
403  call mpi_allreduce &
404  & (val, valm, n, mpi_double_precision, mpi_sum, &
405  & hecmesh%MPI_COMM, ierr)
406 #endif
407  endif
408 
409  if (ntag .eq. hecmw_max) then
410  call mpi_allreduce &
411  & (val, valm, n, mpi_double_precision, mpi_max, &
412  & hecmesh%MPI_COMM, ierr)
413  endif
414 
415  if (ntag .eq. hecmw_min) then
416  call mpi_allreduce &
417  & (val, valm, n, mpi_double_precision, mpi_min, &
418  & hecmesh%MPI_COMM, ierr)
419  endif
420 
421 #ifdef _OPENACC
422  call acc_unmap_data(tmp)
423  call acc_free(tmp_dev)
424 #endif
425 
426  val= valm
427  deallocate (valm)
428 #endif
429  end subroutine hecmw_allreduce_r
430 
431  subroutine hecmw_allreduce_r1 (hecMESH, s, ntag)
432  use hecmw_util
433  implicit none
434  integer(kind=kint):: ntag
435  real(kind=kreal) :: s
436  type (hecmwst_local_mesh) :: hecmesh
437 #ifndef HECMW_SERIAL
438  real(kind=kreal), dimension(1) :: val
439  !$acc serial
440  val(1) = s
441  !$acc end serial
442  call hecmw_allreduce_r(hecmesh, val, 1, ntag )
443  s = val(1)
444 #endif
445  end subroutine hecmw_allreduce_r1
446 
447  !C
448  !C***
449  !C*** hecmw_allREDUCE_I
450  !C***
451  !C
452  subroutine hecmw_allreduce_i(hecMESH, val, n, ntag)
453  use hecmw_util
454  implicit none
455  integer(kind=kint):: n, ntag
456  integer(kind=kint), dimension(n) :: val
457  type (hecmwst_local_mesh) :: hecMESH
458 #ifndef HECMW_SERIAL
459  integer(kind=kint):: ierr
460  integer(kind=kint), dimension(:), allocatable :: VALM
461 
462  allocate (valm(n))
463  valm= 0
464  if (ntag .eq. hecmw_sum) then
465  call mpi_allreduce &
466  & (val, valm, n, mpi_integer, mpi_sum, &
467  & hecmesh%MPI_COMM, ierr)
468  endif
469 
470  if (ntag .eq. hecmw_max) then
471  call mpi_allreduce &
472  & (val, valm, n, mpi_integer, mpi_max, &
473  & hecmesh%MPI_COMM, ierr)
474  endif
475 
476  if (ntag .eq. hecmw_min) then
477  call mpi_allreduce &
478  & (val, valm, n, mpi_integer, mpi_min, &
479  & hecmesh%MPI_COMM, ierr)
480  endif
481 
482 
483  val= valm
484  deallocate (valm)
485 #endif
486  end subroutine hecmw_allreduce_i
487 
488  subroutine hecmw_allreduce_i1 (hecMESH, s, ntag)
489  use hecmw_util
490  implicit none
491  integer(kind=kint):: ntag, s
492  type (hecmwST_local_mesh) :: hecMESH
493 #ifndef HECMW_SERIAL
494  integer(kind=kint), dimension(1) :: val
495 
496  val(1) = s
497  call hecmw_allreduce_i(hecmesh, val, 1, ntag )
498  s = val(1)
499 #endif
500  end subroutine hecmw_allreduce_i1
501 
502  subroutine hecmw_allreduce_l1 (hecMESH, flag, ntag)
503  use hecmw_util
504  implicit none
505  integer(kind=kint):: ntag
506  logical :: flag
507  type (hecmwST_local_mesh) :: hecMESH
508 #ifndef HECMW_SERIAL
509  integer(kind=kint):: ierr
510  logical :: flagM
511 
512  flagm = .false.
513  if (ntag .eq. hecmw_lor) then
514  call mpi_allreduce &
515  & (flag, flagm, 1, mpi_logical, mpi_lor, &
516  & hecmesh%MPI_COMM, ierr)
517  endif
518 
519  if (ntag .eq. hecmw_land) then
520  call mpi_allreduce &
521  & (flag, flagm, 1, mpi_logical, mpi_land, &
522  & hecmesh%MPI_COMM, ierr)
523  endif
524 
525  flag = flagm
526 #endif
527  end subroutine hecmw_allreduce_l1
528 
529  !C
530  !C***
531  !C*** hecmw_bcast_R
532  !C***
533  !C
534  subroutine hecmw_bcast_r (hecMESH, val, n, nbase)
535  use hecmw_util
536  implicit none
537  integer(kind=kint):: n, nbase
538  real(kind=kreal), dimension(n) :: val
539  type (hecmwST_local_mesh) :: hecMESH
540 #ifndef HECMW_SERIAL
541  integer(kind=kint):: ierr
542  call mpi_bcast (val, n, mpi_double_precision, nbase, hecmesh%MPI_COMM, ierr)
543 #endif
544  end subroutine hecmw_bcast_r
545 
546  subroutine hecmw_bcast_r_comm (val, n, nbase, comm)
547  use hecmw_util
548  implicit none
549  integer(kind=kint):: n, nbase
550  real(kind=kreal), dimension(n) :: val
551  integer(kind=kint):: comm
552 #ifndef HECMW_SERIAL
553  integer(kind=kint):: ierr
554  call mpi_bcast (val, n, mpi_double_precision, nbase, comm, ierr)
555 #endif
556  end subroutine hecmw_bcast_r_comm
557 
558  subroutine hecmw_bcast_r1 (hecMESH, s, nbase)
559  use hecmw_util
560  implicit none
561  integer(kind=kint):: nbase, ierr
562  real(kind=kreal) :: s
563  type (hecmwST_local_mesh) :: hecMESH
564 #ifndef HECMW_SERIAL
565  real(kind=kreal), dimension(1) :: val
566  val(1)=s
567  call mpi_bcast (val, 1, mpi_double_precision, nbase, hecmesh%MPI_COMM, ierr)
568  s = val(1)
569 #endif
570  end subroutine hecmw_bcast_r1
571 
572  subroutine hecmw_bcast_r1_comm (s, nbase, comm)
573  use hecmw_util
574  implicit none
575  integer(kind=kint):: nbase
576  real(kind=kreal) :: s
577  integer(kind=kint):: comm
578 #ifndef HECMW_SERIAL
579  integer(kind=kint):: ierr
580  real(kind=kreal), dimension(1) :: val
581  val(1)=s
582  call mpi_bcast (val, 1, mpi_double_precision, nbase, comm, ierr)
583  s = val(1)
584 #endif
585  end subroutine hecmw_bcast_r1_comm
586  !C
587  !C***
588  !C*** hecmw_bcast_I
589  !C***
590  !C
591  subroutine hecmw_bcast_i (hecMESH, val, n, nbase)
592  use hecmw_util
593  implicit none
594  integer(kind=kint):: n, nbase
595  integer(kind=kint), dimension(n) :: val
596  type (hecmwST_local_mesh) :: hecMESH
597 #ifndef HECMW_SERIAL
598  integer(kind=kint):: ierr
599  call mpi_bcast (val, n, mpi_integer, nbase, hecmesh%MPI_COMM, ierr)
600 #endif
601  end subroutine hecmw_bcast_i
602 
603  subroutine hecmw_bcast_i_comm (val, n, nbase, comm)
604  use hecmw_util
605  implicit none
606  integer(kind=kint):: n, nbase
607  integer(kind=kint), dimension(n) :: val
608  integer(kind=kint):: comm
609 #ifndef HECMW_SERIAL
610  integer(kind=kint):: ierr
611  call mpi_bcast (val, n, mpi_integer, nbase, comm, ierr)
612 #endif
613  end subroutine hecmw_bcast_i_comm
614 
615  subroutine hecmw_bcast_i1 (hecMESH, s, nbase)
616  use hecmw_util
617  implicit none
618  integer(kind=kint):: nbase, s
619  type (hecmwST_local_mesh) :: hecMESH
620 #ifndef HECMW_SERIAL
621  integer(kind=kint):: ierr
622  integer(kind=kint), dimension(1) :: val
623  val(1) = s
624  call mpi_bcast (val, 1, mpi_integer, nbase, hecmesh%MPI_COMM, ierr)
625  s = val(1)
626 #endif
627  end subroutine hecmw_bcast_i1
628 
629  subroutine hecmw_bcast_i1_comm (s, nbase, comm)
630  use hecmw_util
631  implicit none
632  integer(kind=kint):: nbase, s
633  integer(kind=kint):: comm
634 #ifndef HECMW_SERIAL
635  integer(kind=kint):: ierr
636  integer(kind=kint), dimension(1) :: val
637  val(1) = s
638  call mpi_bcast (val, 1, mpi_integer, nbase, comm, ierr)
639  s = val(1)
640 #endif
641  end subroutine hecmw_bcast_i1_comm
642  !C
643  !C***
644  !C*** hecmw_bcast_C
645  !C***
646  !C
647  subroutine hecmw_bcast_c (hecMESH, val, n, nn, nbase)
648  use hecmw_util
649  implicit none
650  integer(kind=kint):: n, nn, nbase
651  character(len=n) :: val(nn)
652  type (hecmwST_local_mesh) :: hecMESH
653 #ifndef HECMW_SERIAL
654  integer(kind=kint):: ierr
655  call mpi_bcast (val, n*nn, mpi_character, nbase, hecmesh%MPI_COMM,&
656  & ierr)
657 #endif
658  end subroutine hecmw_bcast_c
659 
660  subroutine hecmw_bcast_c_comm (val, n, nn, nbase, comm)
661  use hecmw_util
662  implicit none
663  integer(kind=kint):: n, nn, nbase
664  character(len=n) :: val(nn)
665  integer(kind=kint):: comm
666 #ifndef HECMW_SERIAL
667  integer(kind=kint):: ierr
668  call mpi_bcast (val, n*nn, mpi_character, nbase, comm,&
669  & ierr)
670 #endif
671  end subroutine hecmw_bcast_c_comm
672 
673  !C
674  !C***
675  !C*** hecmw_assemble_R
676  !C***
677  subroutine hecmw_assemble_r (hecMESH, val, n, m)
678  use hecmw_util
679  use hecmw_solver_sr
680 
681  implicit none
682  integer(kind=kint):: n, m
683  real(kind=kreal), dimension(m*n) :: val
684  type (hecmwST_local_mesh) :: hecMESH
685 #ifndef HECMW_SERIAL
686  integer(kind=kint):: ns, nr
687  real(kind=kreal), dimension(:), allocatable :: ws, wr
688 
689  if( hecmesh%n_neighbor_pe == 0 ) return
690 
691  ns = hecmesh%import_index(hecmesh%n_neighbor_pe)
692  nr = hecmesh%export_index(hecmesh%n_neighbor_pe)
693 
694  allocate (ws(m*ns), wr(m*nr))
696  & ( n, m, hecmesh%n_neighbor_pe, hecmesh%neighbor_pe, &
697  & hecmesh%import_index, hecmesh%import_item, &
698  & hecmesh%export_index, hecmesh%export_item, &
699  & ws, wr, val , hecmesh%MPI_COMM, hecmesh%my_rank)
700  deallocate (ws, wr)
701 #endif
702  end subroutine hecmw_assemble_r
703 
704  !C
705  !C***
706  !C*** hecmw_assemble_I
707  !C***
708  subroutine hecmw_assemble_i (hecMESH, val, n, m)
709  use hecmw_util
711 
712  implicit none
713  integer(kind=kint):: n, m
714  integer(kind=kint), dimension(m*n) :: val
715  type (hecmwST_local_mesh) :: hecMESH
716 #ifndef HECMW_SERIAL
717  integer(kind=kint):: ns, nr
718  integer(kind=kint), dimension(:), allocatable :: WS, WR
719 
720  if( hecmesh%n_neighbor_pe == 0 ) return
721 
722  ns = hecmesh%import_index(hecmesh%n_neighbor_pe)
723  nr = hecmesh%export_index(hecmesh%n_neighbor_pe)
724 
725  allocate (ws(m*ns), wr(m*nr))
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  & ws, wr, val , hecmesh%MPI_COMM, hecmesh%my_rank)
731  deallocate (ws, wr)
732 #endif
733  end subroutine hecmw_assemble_i
734 
735  !C
736  !C***
737  !C*** hecmw_update_R
738  !C***
739  !C
740  subroutine hecmw_update_r (hecMESH, val, n, m)
741  use hecmw_util
742  use hecmw_solver_sr
743 
744  implicit none
745  integer(kind=kint):: n, m
746  real(kind=kreal), dimension(m*n) :: val
747  type (hecmwST_local_mesh) :: hecMESH
748 #ifndef HECMW_SERIAL
749  integer(kind=kint):: ns, nr
750  real(kind=kreal), dimension(:), allocatable :: ws, wr
751 
752  if( hecmesh%n_neighbor_pe == 0 ) return
753 
754  ns = hecmesh%export_index(hecmesh%n_neighbor_pe)
755  nr = hecmesh%import_index(hecmesh%n_neighbor_pe)
756 
757  allocate (ws(m*ns), wr(m*nr))
758  call hecmw_solve_send_recv &
759  & ( n, m, hecmesh%n_neighbor_pe, hecmesh%neighbor_pe, &
760  & hecmesh%import_index, hecmesh%import_item, &
761  & hecmesh%export_index, hecmesh%export_item, &
762  & ws, wr, val , hecmesh%MPI_COMM, hecmesh%my_rank)
763  deallocate (ws, wr)
764 #endif
765  end subroutine hecmw_update_r
766 
767  !C
768  !C***
769  !C*** hecmw_update_R_async
770  !C***
771  !C
772  !C REAL
773  !C
774  subroutine hecmw_update_r_async (hecMESH, val, n, m, ireq)
775  use hecmw_util
776  use hecmw_solver_sr
777  implicit none
778  integer(kind=kint):: n, m, ireq
779  real(kind=kreal), dimension(m*n) :: val
780  type (hecmwST_local_mesh) :: hecMESH
781 #ifndef HECMW_SERIAL
782  if( hecmesh%n_neighbor_pe == 0 ) return
783 
785  & ( n, m, hecmesh%n_neighbor_pe, hecmesh%neighbor_pe, &
786  & hecmesh%import_index, hecmesh%import_item, &
787  & hecmesh%export_index, hecmesh%export_item, &
788  & val , hecmesh%MPI_COMM, hecmesh%my_rank, ireq)
789 #endif
790  end subroutine hecmw_update_r_async
791 
792  !C
793  !C***
794  !C*** hecmw_update_R_wait
795  !C***
796  !C
797  !C REAL
798  !C
799  subroutine hecmw_update_r_wait (hecMESH, ireq)
800  use hecmw_util
801  use hecmw_solver_sr
802  implicit none
803  integer(kind=kint):: ireq
804  type (hecmwST_local_mesh) :: hecMESH
805 #ifndef HECMW_SERIAL
806  if( hecmesh%n_neighbor_pe == 0 ) return
807 
808  call hecmw_solve_isend_irecv_wait(hecmesh%n_dof, ireq )
809 #endif
810  end subroutine hecmw_update_r_wait
811 
812 
813  !C
814  !C***
815  !C*** hecmw_update_I
816  !C***
817  !C
818  subroutine hecmw_update_i (hecMESH, val, n, m)
819  use hecmw_util
821 
822  implicit none
823  integer(kind=kint):: n, m
824  integer(kind=kint), dimension(m*n) :: val
825  type (hecmwST_local_mesh) :: hecMESH
826 #ifndef HECMW_SERIAL
827  integer(kind=kint):: ns, nr
828  integer(kind=kint), dimension(:), allocatable :: WS, WR
829 
830  if( hecmesh%n_neighbor_pe == 0 ) return
831 
832  ns = hecmesh%export_index(hecmesh%n_neighbor_pe)
833  nr = hecmesh%import_index(hecmesh%n_neighbor_pe)
834 
835  allocate (ws(m*ns), wr(m*nr))
837  & ( n, m, hecmesh%n_neighbor_pe, hecmesh%neighbor_pe, &
838  & hecmesh%import_index, hecmesh%import_item, &
839  & hecmesh%export_index, hecmesh%export_item, &
840  & ws, wr, val , hecmesh%MPI_COMM, hecmesh%my_rank)
841  deallocate (ws, wr)
842 #endif
843  end subroutine hecmw_update_i
844 
845 end module m_hecmw_comm_f
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)
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)
subroutine hecmw_solve_isend_irecv_wait(m, ireq)
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=kint), parameter hecmw_land
integer(kind=kint), parameter hecmw_sum
integer(kind=kint), parameter hecmw_prod
integer(kind=kint), parameter hecmw_max
integer(kind=4), parameter kreal
integer(kind=kint), parameter hecmw_lor
integer(kind=kint), parameter hecmw_min
subroutine hecmw_assemble_r(hecMESH, val, n, m)
subroutine hecmw_bcast_c_comm(val, n, nn, nbase, comm)
subroutine hecmw_gatherv_real(sbuf, sc, rbuf, rcs, disp, root, comm)
subroutine hecmw_update_r(hecMESH, val, n, m)
subroutine hecmw_isend_int(sbuf, sc, dest, tag, comm, req)
subroutine hecmw_scatterv_real(sbuf, scs, disp, rbuf, rc, root, comm)
subroutine hecmw_recv_int(rbuf, rc, source, tag, comm, stat)
subroutine hecmw_bcast_r_comm(val, n, nbase, comm)
subroutine hecmw_bcast_c(hecMESH, val, n, nn, nbase)
subroutine hecmw_isend_r(sbuf, sc, dest, tag, comm, req)
subroutine hecmw_allreduce_i1(hecMESH, s, ntag)
subroutine hecmw_gather_int_1(sval, rbuf, root, comm)
subroutine hecmw_bcast_i_comm(val, n, nbase, comm)
subroutine hecmw_allreduce_i(hecMESH, val, n, ntag)
subroutine hecmw_allreduce_r(hecMESH, val, n, ntag)
subroutine hecmw_bcast_i1(hecMESH, s, nbase)
subroutine hecmw_allgather_int_1(sval, rbuf, comm)
subroutine hecmw_allreduce_dp1(s1, s2, hec_op, comm)
subroutine hecmw_waitall(cnt, reqs, stats)
subroutine hecmw_assemble_i(hecMESH, val, n, m)
subroutine hecmw_allreduce_dp(val, VALM, n, hec_op, comm)
subroutine hecmw_bcast_r(hecMESH, val, n, nbase)
subroutine hecmw_irecv_int(rbuf, rc, source, tag, comm, req)
subroutine hecmw_bcast_r1_comm(s, nbase, comm)
subroutine hecmw_scatter_int_1(sbuf, rval, root, comm)
subroutine hecmw_irecv_r(rbuf, rc, source, tag, comm, req)
subroutine hecmw_allreduce_l1(hecMESH, flag, ntag)
subroutine hecmw_bcast_i(hecMESH, val, n, nbase)
subroutine hecmw_update_r_async(hecMESH, val, n, m, ireq)
subroutine hecmw_update_i(hecMESH, val, n, m)
subroutine hecmw_recv_r(rbuf, rc, source, tag, comm, stat)
subroutine hecmw_allreduce_r1(hecMESH, s, ntag)
subroutine hecmw_alltoall_int(sbuf, sc, rbuf, rc, comm)
subroutine hecmw_update_r_wait(hecMESH, ireq)
subroutine hecmw_bcast_r1(hecMESH, s, nbase)
subroutine hecmw_scatterv_dp(sbuf, sc, disp, rbuf, rc, root, comm)
subroutine hecmw_allreduce_int_1(sval, rval, op, comm)
subroutine hecmw_gatherv_int(sbuf, sc, rbuf, rcs, disp, root, comm)
subroutine hecmw_barrier(hecMESH)
integer(kind=kint) function hecmw_operation_hec2mpi(operation)
subroutine hecmw_bcast_i1_comm(s, nbase, comm)