FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
hecmw_solver_misc.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 
8 contains
9 
10  !C
11  !C***
12  !C*** hecmw_innerProduct_I
13  !C***
14  !C
15  subroutine hecmw_innerproduct_i (hecMESH, ndof, X, Y, sum, COMMtime )
16  use hecmw_util
17  use m_hecmw_comm_f
18 
19  implicit none
20  type (hecmwST_local_mesh) :: hecMESH
21  integer(kind=kint) :: ndof
22  integer(kind=kint) :: X(:), Y(:)
23  integer(kind=kint) :: sum
24  real(kind=kreal), optional :: commtime
25 
26  integer(kind=kint) :: i
27  real(kind=kreal) :: start_time, end_time
28 
29  sum = 0
30  do i = 1, hecmesh%nn_internal * ndof
31  sum = sum + x(i)*y(i)
32  end do
33 
34  start_time= hecmw_wtime()
35  call hecmw_allreduce_i1 (hecmesh, sum, hecmw_sum)
36  end_time= hecmw_wtime()
37  if (present(commtime)) commtime = commtime + end_time - start_time
38 
39  end subroutine hecmw_innerproduct_i
40 
41  !C
42  !C***
43  !C*** hecmw_innerProduct_R
44  !C***
45  !C
46  subroutine hecmw_innerproduct_r (hecMESH, ndof, X, Y, sum, COMMtime )
47  use hecmw_util
48  use m_hecmw_comm_f
49 
50  implicit none
51  type (hecmwST_local_mesh) :: hecMESH
52  integer(kind=kint) :: ndof
53  real(kind=kreal) :: x(:), y(:)
54  real(kind=kreal) :: sum
55  real(kind=kreal), optional :: commtime
56 
57  integer(kind=kint) :: i
58  real(kind=kreal) :: start_time, end_time
59 
60  !$acc serial
61  sum = 0.0d0
62  !$acc end serial
63 
64  !$acc kernels
65  !$acc loop reduction(+:sum)
66  do i = 1, hecmesh%nn_internal * ndof
67  sum = sum + x(i)*y(i)
68  end do
69  !$acc end kernels
70 
71  start_time= hecmw_wtime()
72  call hecmw_allreduce_r1 (hecmesh, sum, hecmw_sum)
73  end_time= hecmw_wtime()
74  if (present(commtime)) commtime = commtime + end_time - start_time
75 
76  end subroutine hecmw_innerproduct_r
77 
78  !C
79  !C***
80  !C*** hecmw_absMax_R
81  !C***
82  !C
83  subroutine hecmw_absmax_r (hecMESH, ndof, X, absMax, COMMtime)
84  use hecmw_util
85  use m_hecmw_comm_f
86 
87  implicit none
88  type (hecmwST_local_mesh) :: hecMESH
89  integer(kind=kint) :: ndof
90  real(kind=kreal) :: x(:)
91  real(kind=kreal) :: absmax
92  real(kind=kreal), optional :: commtime
93 
94  integer(kind=kint) :: i
95  real(kind=kreal) :: start_time, end_time
96 
97  absmax = 0.0d0
98  do i = 1, hecmesh%nn_internal * ndof
99  absmax = max(absmax, abs(x(i)))
100  end do
101 
102  start_time= hecmw_wtime()
103  call hecmw_allreduce_r1 (hecmesh, absmax, hecmw_max)
104  end_time= hecmw_wtime()
105  if (present(commtime)) commtime = commtime + end_time - start_time
106  end subroutine hecmw_absmax_r
107 
108 
109  !C
110  !C***
111  !C*** hecmw_innerProduct_R_nocomm
112  !C***
113  !C
114  subroutine hecmw_innerproduct_r_nocomm (hecMESH, ndof, X, Y, sum)
115  use hecmw_util
116 
117  implicit none
118  type (hecmwST_local_mesh) :: hecMESH
119  integer(kind=kint) :: ndof
120  real(kind=kreal) :: x(:), y(:)
121  real(kind=kreal) :: sum
122 
123  integer(kind=kint) :: i
124 
125  sum = 0.0d0
126  do i = 1, hecmesh%nn_internal * ndof
127  sum = sum + x(i)*y(i)
128  end do
129 
130  end subroutine hecmw_innerproduct_r_nocomm
131 
132  !C
133  !C***
134  !C*** hecmw_axpy_R
135  !C***
136  !C
137  subroutine hecmw_axpy_r (n, alpha, X, Y)
138 
139  use hecmw_util
140 
141  implicit none
142  integer(kind=kint) :: n
143  real(kind=kreal) :: x(:), y(:)
144  real(kind=kreal) :: alpha
145 
146  integer(kind=kint) :: i
147 
148 
149 #ifdef _OPENACC
150  !$acc kernels
151  !$acc loop independent
152 #else
153  !$OMP PARALLEL
154  !$OMP DO
155 #endif
156  do i = 1, n
157  y(i) = alpha * x(i) + y(i)
158  end do
159 #ifdef _OPENACC
160  !$acc end kernels
161 #else
162  !$OMP END DO
163  !$OMP END PARALLEL
164 #endif
165 
166  end subroutine hecmw_axpy_r
167 
168  !C
169  !C***
170  !C*** hecmw_xpay_R
171  !C***
172  !C
173  subroutine hecmw_xpay_r (n, alpha, X, Y)
174 
175  use hecmw_util
176 
177  implicit none
178  integer(kind=kint) :: n
179  real(kind=kreal) :: x(:), y(:)
180  real(kind=kreal) :: alpha
181 
182  integer(kind=kint) :: i
183 
184 #ifdef _OPENACC
185  !$acc kernels
186  !$acc loop independent
187 #else
188  !$OMP PARALLEL
189  !$OMP DO
190 #endif
191  do i = 1, n
192  y(i) = x(i) + alpha * y(i)
193  end do
194 #ifdef _OPENACC
195  !$acc end kernels
196 #else
197  !$OMP END DO
198  !$OMP END PARALLEL
199 #endif
200 
201  end subroutine hecmw_xpay_r
202 
203  !C
204  !C***
205  !C*** hecmw_axpby_R
206  !C***
207  !C
208  subroutine hecmw_axpby_r (n, alpha, beta, X, Y)
209 
210  use hecmw_util
211 
212  implicit none
213  integer(kind=kint) :: n
214  real(kind=kreal) :: x(:), y(:)
215  real(kind=kreal) :: alpha, beta
216 
217  integer(kind=kint) :: i
218 
219 
220  if (beta == 0.d0) then
221 #ifdef _OPENACC
222  !$acc kernels
223  !$acc loop independent
224 #else
225  !$OMP PARALLEL
226  !$OMP DO
227 #endif
228  do i = 1, n
229  y(i) = alpha * x(i)
230  end do
231 #ifdef _OPENACC
232  !$acc end kernels
233 #else
234  !$OMP END DO
235  !$OMP END PARALLEL
236 #endif
237  else
238 #ifdef _OPENACC
239  !$acc kernels
240  !$acc loop independent
241 #else
242  !$OMP PARALLEL
243  !$OMP DO
244 #endif
245  do i = 1, n
246  y(i) = alpha * x(i) + beta * y(i)
247  end do
248 #ifdef _OPENACC
249  !$acc end kernels
250 #else
251  !$OMP END DO
252  !$OMP END PARALLEL
253 #endif
254  end if
255 
256  end subroutine hecmw_axpby_r
257 
258  !C
259  !C***
260  !C*** hecmw_axpyz_R
261  !C***
262  !C
263  subroutine hecmw_axpyz_r (n, alpha, X, Y, Z)
264 
265  use hecmw_util
266 
267  implicit none
268  integer(kind=kint) :: n
269  real(kind=kreal) :: x(:), y(:), z(:)
270  real(kind=kreal) :: alpha
271 
272  integer(kind=kint) :: i
273 
274 
275 #ifdef _OPENACC
276  !$acc kernels
277  !$acc loop independent
278 #else
279  !$OMP PARALLEL
280  !$OMP DO
281 #endif
282  do i = 1, n
283  z(i) = alpha * x(i) + y(i)
284  end do
285 #ifdef _OPENACC
286  !$acc end kernels
287 #else
288  !$OMP END DO
289  !$OMP END PARALLEL
290 #endif
291  end subroutine hecmw_axpyz_r
292 
293  !C
294  !C***
295  !C*** hecmw_scale_R
296  !C***
297  !C
298  subroutine hecmw_scale_r (n, alpha, X)
299 
300  use hecmw_util
301 
302  implicit none
303  integer(kind=kint) :: n
304  real(kind=kreal) :: x(:)
305  real(kind=kreal) :: alpha
306 
307  integer(kind=kint) :: i
308 
309  if (alpha == 0.d0) then
310 #ifdef _OPENACC
311  !$acc kernels
312  !$acc loop independent
313 #else
314  !$OMP PARALLEL
315  !$OMP DO
316 #endif
317  do i = 1, n
318  x(i) = 0.d0
319  end do
320 #ifdef _OPENACC
321  !$acc end kernels
322 #else
323  !$OMP END DO
324  !$OMP END PARALLEL
325 #endif
326  else
327 #ifdef _OPENACC
328  !$acc kernels
329  !$acc loop independent
330 #else
331  !$OMP PARALLEL
332  !$OMP DO
333 #endif
334  do i = 1, n
335  x(i) = alpha * x(i)
336  end do
337 #ifdef _OPENACC
338  !$acc end kernels
339 #else
340  !$OMP END DO
341  !$OMP END PARALLEL
342 #endif
343  end if
344 
345  end subroutine hecmw_scale_r
346 
347  !C
348  !C***
349  !C*** hecmw_copy_R
350  !C***
351  !C
352  subroutine hecmw_copy_r (n, X, Y)
353 
354  use hecmw_util
355 
356  implicit none
357  integer(kind=kint) :: n
358  real(kind=kreal) :: x(:), y(:)
359 
360  integer(kind=kint) :: i
361 
362 
363 #ifdef _OPENACC
364  !$acc kernels
365  !$acc loop independent
366 #else
367  !$OMP PARALLEL
368  !$OMP DO
369 #endif
370  do i = 1, n
371  y(i) = x(i)
372  end do
373 #ifdef _OPENACC
374  !$acc end kernels
375 #else
376  !$OMP END DO
377  !$OMP END PARALLEL
378 #endif
379 
380  end subroutine hecmw_copy_r
381 
382  !C
383  !C***
384  !C*** hecmw_time_statistics
385  !C***
386  !C
387  subroutine hecmw_time_statistics (hecMESH, time, &
388  t_max, t_min, t_avg, t_sd)
389  use hecmw_util
390  use m_hecmw_comm_f
391 
392  implicit none
393  type (hecmwST_local_mesh), intent(in) :: hecMESH
394  real(kind=kreal), intent(in) :: time
395  real(kind=kreal), intent(out) :: t_max
396  real(kind=kreal), intent(out), optional :: t_min, t_avg, t_sd
397  real(kind=kreal) :: t2_avg
398  integer(kind=kint) :: nprocs
399 
400  nprocs = hecmw_comm_get_size()
401 
402  t_max = time
403  call hecmw_allreduce_r1(hecmesh, t_max, hecmw_max)
404 
405  if (.not. present(t_min)) return
406  t_min = time
407  call hecmw_allreduce_r1(hecmesh, t_min, hecmw_min)
408 
409  if (.not. present(t_avg)) return
410  t_avg = time
411  call hecmw_allreduce_r1(hecmesh, t_avg, hecmw_sum)
412  t_avg = t_avg / nprocs
413 
414  if (.not. present(t_sd)) return
415  t2_avg = time*time
416  call hecmw_allreduce_r1(hecmesh, t2_avg, hecmw_sum)
417  t2_avg = t2_avg / nprocs
418 
419  t_sd = dsqrt(t2_avg - t_avg*t_avg)
420  end subroutine hecmw_time_statistics
421 
422 end module hecmw_solver_misc
subroutine hecmw_xpay_r(n, alpha, X, Y)
subroutine hecmw_absmax_r(hecMESH, ndof, X, absMax, COMMtime)
subroutine hecmw_axpyz_r(n, alpha, X, Y, Z)
subroutine hecmw_innerproduct_r_nocomm(hecMESH, ndof, X, Y, sum)
subroutine hecmw_innerproduct_r(hecMESH, ndof, X, Y, sum, COMMtime)
subroutine hecmw_scale_r(n, alpha, X)
subroutine hecmw_axpby_r(n, alpha, beta, X, Y)
subroutine hecmw_axpy_r(n, alpha, X, Y)
subroutine hecmw_copy_r(n, X, Y)
subroutine hecmw_time_statistics(hecMESH, time, t_max, t_min, t_avg, t_sd)
subroutine hecmw_innerproduct_i(hecMESH, ndof, X, Y, sum, COMMtime)
I/O and Utility.
Definition: hecmw_util_f.F90:7
integer(kind=kint), parameter hecmw_sum
integer(kind=kint) function hecmw_comm_get_size()
integer(kind=kint), parameter hecmw_max
integer(kind=4), parameter kreal
integer(kind=kint), parameter hecmw_min
real(kind=kreal) function hecmw_wtime()
subroutine hecmw_allreduce_i1(hecMESH, s, ntag)
subroutine hecmw_allreduce_r1(hecMESH, s, ntag)