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  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)
403  use hecmw_util
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)
422  use hecmw_util
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)
458  use hecmw_util
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  subroutine hecmw_allreduce_l1 (hecMESH, flag, ntag)
472  use hecmw_util
473  implicit none
474  integer(kind=kint):: ntag
475  logical :: flag
476  type (hecmwST_local_mesh) :: hecMESH
477 #ifndef HECMW_SERIAL
478  integer(kind=kint):: ierr
479  logical :: flagM
480 
481  flagm = .false.
482  if (ntag .eq. hecmw_lor) then
483  call mpi_allreduce &
484  & (flag, flagm, 1, mpi_logical, mpi_lor, &
485  & hecmesh%MPI_COMM, ierr)
486  endif
487 
488  if (ntag .eq. hecmw_land) then
489  call mpi_allreduce &
490  & (flag, flagm, 1, mpi_logical, mpi_land, &
491  & hecmesh%MPI_COMM, ierr)
492  endif
493 
494  flag = flagm
495 #endif
496  end subroutine hecmw_allreduce_l1
497 
498  !C
499  !C***
500  !C*** hecmw_bcast_R
501  !C***
502  !C
503  subroutine hecmw_bcast_r (hecMESH, val, n, nbase)
504  use hecmw_util
505  implicit none
506  integer(kind=kint):: n, nbase
507  real(kind=kreal), dimension(n) :: val
508  type (hecmwST_local_mesh) :: hecMESH
509 #ifndef HECMW_SERIAL
510  integer(kind=kint):: ierr
511  call mpi_bcast (val, n, mpi_double_precision, nbase, hecmesh%MPI_COMM, ierr)
512 #endif
513  end subroutine hecmw_bcast_r
514 
515  subroutine hecmw_bcast_r_comm (val, n, nbase, comm)
516  use hecmw_util
517  implicit none
518  integer(kind=kint):: n, nbase
519  real(kind=kreal), dimension(n) :: val
520  integer(kind=kint):: comm
521 #ifndef HECMW_SERIAL
522  integer(kind=kint):: ierr
523  call mpi_bcast (val, n, mpi_double_precision, nbase, comm, ierr)
524 #endif
525  end subroutine hecmw_bcast_r_comm
526 
527  subroutine hecmw_bcast_r1 (hecMESH, s, nbase)
528  use hecmw_util
529  implicit none
530  integer(kind=kint):: nbase, ierr
531  real(kind=kreal) :: s
532  type (hecmwST_local_mesh) :: hecMESH
533 #ifndef HECMW_SERIAL
534  real(kind=kreal), dimension(1) :: val
535  val(1)=s
536  call mpi_bcast (val, 1, mpi_double_precision, nbase, hecmesh%MPI_COMM, ierr)
537  s = val(1)
538 #endif
539  end subroutine hecmw_bcast_r1
540 
541  subroutine hecmw_bcast_r1_comm (s, nbase, comm)
542  use hecmw_util
543  implicit none
544  integer(kind=kint):: nbase
545  real(kind=kreal) :: s
546  integer(kind=kint):: comm
547 #ifndef HECMW_SERIAL
548  integer(kind=kint):: ierr
549  real(kind=kreal), dimension(1) :: val
550  val(1)=s
551  call mpi_bcast (val, 1, mpi_double_precision, nbase, comm, ierr)
552  s = val(1)
553 #endif
554  end subroutine hecmw_bcast_r1_comm
555  !C
556  !C***
557  !C*** hecmw_bcast_I
558  !C***
559  !C
560  subroutine hecmw_bcast_i (hecMESH, val, n, nbase)
561  use hecmw_util
562  implicit none
563  integer(kind=kint):: n, nbase
564  integer(kind=kint), dimension(n) :: val
565  type (hecmwST_local_mesh) :: hecMESH
566 #ifndef HECMW_SERIAL
567  integer(kind=kint):: ierr
568  call mpi_bcast (val, n, mpi_integer, nbase, hecmesh%MPI_COMM, ierr)
569 #endif
570  end subroutine hecmw_bcast_i
571 
572  subroutine hecmw_bcast_i_comm (val, n, nbase, comm)
573  use hecmw_util
574  implicit none
575  integer(kind=kint):: n, nbase
576  integer(kind=kint), dimension(n) :: val
577  integer(kind=kint):: comm
578 #ifndef HECMW_SERIAL
579  integer(kind=kint):: ierr
580  call mpi_bcast (val, n, mpi_integer, nbase, comm, ierr)
581 #endif
582  end subroutine hecmw_bcast_i_comm
583 
584  subroutine hecmw_bcast_i1 (hecMESH, s, nbase)
585  use hecmw_util
586  implicit none
587  integer(kind=kint):: nbase, s
588  type (hecmwST_local_mesh) :: hecMESH
589 #ifndef HECMW_SERIAL
590  integer(kind=kint):: ierr
591  integer(kind=kint), dimension(1) :: val
592  val(1) = s
593  call mpi_bcast (val, 1, mpi_integer, nbase, hecmesh%MPI_COMM, ierr)
594  s = val(1)
595 #endif
596  end subroutine hecmw_bcast_i1
597 
598  subroutine hecmw_bcast_i1_comm (s, nbase, comm)
599  use hecmw_util
600  implicit none
601  integer(kind=kint):: nbase, s
602  integer(kind=kint):: comm
603 #ifndef HECMW_SERIAL
604  integer(kind=kint):: ierr
605  integer(kind=kint), dimension(1) :: val
606  val(1) = s
607  call mpi_bcast (val, 1, mpi_integer, nbase, comm, ierr)
608  s = val(1)
609 #endif
610  end subroutine hecmw_bcast_i1_comm
611  !C
612  !C***
613  !C*** hecmw_bcast_C
614  !C***
615  !C
616  subroutine hecmw_bcast_c (hecMESH, val, n, nn, nbase)
617  use hecmw_util
618  implicit none
619  integer(kind=kint):: n, nn, nbase
620  character(len=n) :: val(nn)
621  type (hecmwST_local_mesh) :: hecMESH
622 #ifndef HECMW_SERIAL
623  integer(kind=kint):: ierr
624  call mpi_bcast (val, n*nn, mpi_character, nbase, hecmesh%MPI_COMM,&
625  & ierr)
626 #endif
627  end subroutine hecmw_bcast_c
628 
629  subroutine hecmw_bcast_c_comm (val, n, nn, nbase, comm)
630  use hecmw_util
631  implicit none
632  integer(kind=kint):: n, nn, nbase
633  character(len=n) :: val(nn)
634  integer(kind=kint):: comm
635 #ifndef HECMW_SERIAL
636  integer(kind=kint):: ierr
637  call mpi_bcast (val, n*nn, mpi_character, nbase, comm,&
638  & ierr)
639 #endif
640  end subroutine hecmw_bcast_c_comm
641 
642  !C
643  !C***
644  !C*** hecmw_assemble_R
645  !C***
646  subroutine hecmw_assemble_r (hecMESH, val, n, m)
647  use hecmw_util
648  use hecmw_solver_sr
649 
650  implicit none
651  integer(kind=kint):: n, m
652  real(kind=kreal), dimension(m*n) :: val
653  type (hecmwST_local_mesh) :: hecMESH
654 #ifndef HECMW_SERIAL
655  integer(kind=kint):: ns, nr
656  real(kind=kreal), dimension(:), allocatable :: ws, wr
657 
658  if( hecmesh%n_neighbor_pe == 0 ) return
659 
660  ns = hecmesh%import_index(hecmesh%n_neighbor_pe)
661  nr = hecmesh%export_index(hecmesh%n_neighbor_pe)
662 
663  allocate (ws(m*ns), wr(m*nr))
665  & ( n, m, hecmesh%n_neighbor_pe, hecmesh%neighbor_pe, &
666  & hecmesh%import_index, hecmesh%import_item, &
667  & hecmesh%export_index, hecmesh%export_item, &
668  & ws, wr, val , hecmesh%MPI_COMM, hecmesh%my_rank)
669  deallocate (ws, wr)
670 #endif
671  end subroutine hecmw_assemble_r
672 
673  !C
674  !C***
675  !C*** hecmw_assemble_I
676  !C***
677  subroutine hecmw_assemble_i (hecMESH, val, n, m)
678  use hecmw_util
680 
681  implicit none
682  integer(kind=kint):: n, m
683  integer(kind=kint), dimension(m*n) :: val
684  type (hecmwST_local_mesh) :: hecMESH
685 #ifndef HECMW_SERIAL
686  integer(kind=kint):: ns, nr
687  integer(kind=kint), 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_i
703 
704  !C
705  !C***
706  !C*** hecmw_update_R
707  !C***
708  !C
709  subroutine hecmw_update_r (hecMESH, val, n, m)
710  use hecmw_util
711  use hecmw_solver_sr
712 
713  implicit none
714  integer(kind=kint):: n, m
715  real(kind=kreal), dimension(m*n) :: val
716  type (hecmwST_local_mesh) :: hecMESH
717 #ifndef HECMW_SERIAL
718  integer(kind=kint):: ns, nr
719  real(kind=kreal), dimension(:), allocatable :: ws, wr
720 
721  if( hecmesh%n_neighbor_pe == 0 ) return
722 
723  ns = hecmesh%export_index(hecmesh%n_neighbor_pe)
724  nr = hecmesh%import_index(hecmesh%n_neighbor_pe)
725 
726  allocate (ws(m*ns), wr(m*nr))
727  call hecmw_solve_send_recv &
728  & ( n, m, hecmesh%n_neighbor_pe, hecmesh%neighbor_pe, &
729  & hecmesh%import_index, hecmesh%import_item, &
730  & hecmesh%export_index, hecmesh%export_item, &
731  & ws, wr, val , hecmesh%MPI_COMM, hecmesh%my_rank)
732  deallocate (ws, wr)
733 #endif
734  end subroutine hecmw_update_r
735 
736  !C
737  !C***
738  !C*** hecmw_update_R_async
739  !C***
740  !C
741  !C REAL
742  !C
743  subroutine hecmw_update_r_async (hecMESH, val, n, m, ireq)
744  use hecmw_util
745  use hecmw_solver_sr
746  implicit none
747  integer(kind=kint):: n, m, ireq
748  real(kind=kreal), dimension(m*n) :: val
749  type (hecmwST_local_mesh) :: hecMESH
750 #ifndef HECMW_SERIAL
751  if( hecmesh%n_neighbor_pe == 0 ) return
752 
754  & ( n, m, hecmesh%n_neighbor_pe, hecmesh%neighbor_pe, &
755  & hecmesh%import_index, hecmesh%import_item, &
756  & hecmesh%export_index, hecmesh%export_item, &
757  & val , hecmesh%MPI_COMM, hecmesh%my_rank, ireq)
758 #endif
759  end subroutine hecmw_update_r_async
760 
761  !C
762  !C***
763  !C*** hecmw_update_R_wait
764  !C***
765  !C
766  !C REAL
767  !C
768  subroutine hecmw_update_r_wait (hecMESH, ireq)
769  use hecmw_util
770  use hecmw_solver_sr
771  implicit none
772  integer(kind=kint):: ireq
773  type (hecmwST_local_mesh) :: hecMESH
774 #ifndef HECMW_SERIAL
775  if( hecmesh%n_neighbor_pe == 0 ) return
776 
777  call hecmw_solve_isend_irecv_wait(hecmesh%n_dof, ireq )
778 #endif
779  end subroutine hecmw_update_r_wait
780 
781 
782  !C
783  !C***
784  !C*** hecmw_update_I
785  !C***
786  !C
787  subroutine hecmw_update_i (hecMESH, val, n, m)
788  use hecmw_util
790 
791  implicit none
792  integer(kind=kint):: n, m
793  integer(kind=kint), dimension(m*n) :: val
794  type (hecmwST_local_mesh) :: hecMESH
795 #ifndef HECMW_SERIAL
796  integer(kind=kint):: ns, nr
797  integer(kind=kint), dimension(:), allocatable :: WS, WR
798 
799  if( hecmesh%n_neighbor_pe == 0 ) return
800 
801  ns = hecmesh%export_index(hecmesh%n_neighbor_pe)
802  nr = hecmesh%import_index(hecmesh%n_neighbor_pe)
803 
804  allocate (ws(m*ns), wr(m*nr))
806  & ( n, m, hecmesh%n_neighbor_pe, hecmesh%neighbor_pe, &
807  & hecmesh%import_index, hecmesh%import_item, &
808  & hecmesh%export_index, hecmesh%export_item, &
809  & ws, wr, val , hecmesh%MPI_COMM, hecmesh%my_rank)
810  deallocate (ws, wr)
811 #endif
812  end subroutine hecmw_update_i
813 
814 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)