FrontISTR  5.8.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  sum = 0.0d0
61  do i = 1, hecmesh%nn_internal * ndof
62  sum = sum + x(i)*y(i)
63  end do
64 
65  start_time= hecmw_wtime()
66  call hecmw_allreduce_r1 (hecmesh, sum, hecmw_sum)
67  end_time= hecmw_wtime()
68  if (present(commtime)) commtime = commtime + end_time - start_time
69 
70  end subroutine hecmw_innerproduct_r
71 
72  !C
73  !C***
74  !C*** hecmw_absMax_R
75  !C***
76  !C
77  subroutine hecmw_absmax_r (hecMESH, ndof, X, absMax, COMMtime)
78  use hecmw_util
79  use m_hecmw_comm_f
80 
81  implicit none
82  type (hecmwST_local_mesh) :: hecMESH
83  integer(kind=kint) :: ndof
84  real(kind=kreal) :: x(:)
85  real(kind=kreal) :: absmax
86  real(kind=kreal), optional :: commtime
87 
88  integer(kind=kint) :: i
89  real(kind=kreal) :: start_time, end_time
90 
91  absmax = 0.0d0
92  do i = 1, hecmesh%nn_internal * ndof
93  absmax = max(absmax, abs(x(i)))
94  end do
95 
96  start_time= hecmw_wtime()
97  call hecmw_allreduce_r1 (hecmesh, absmax, hecmw_max)
98  end_time= hecmw_wtime()
99  if (present(commtime)) commtime = commtime + end_time - start_time
100  end subroutine hecmw_absmax_r
101 
102 
103  !C
104  !C***
105  !C*** hecmw_innerProduct_R_nocomm
106  !C***
107  !C
108  subroutine hecmw_innerproduct_r_nocomm (hecMESH, ndof, X, Y, sum)
109  use hecmw_util
110 
111  implicit none
112  type (hecmwST_local_mesh) :: hecMESH
113  integer(kind=kint) :: ndof
114  real(kind=kreal) :: x(:), y(:)
115  real(kind=kreal) :: sum
116 
117  integer(kind=kint) :: i
118 
119  sum = 0.0d0
120  do i = 1, hecmesh%nn_internal * ndof
121  sum = sum + x(i)*y(i)
122  end do
123 
124  end subroutine hecmw_innerproduct_r_nocomm
125 
126  !C
127  !C***
128  !C*** hecmw_axpy_R
129  !C***
130  !C
131  subroutine hecmw_axpy_r (n, alpha, X, Y)
132 
133  use hecmw_util
134 
135  implicit none
136  integer(kind=kint) :: n
137  real(kind=kreal) :: x(:), y(:)
138  real(kind=kreal) :: alpha
139 
140  integer(kind=kint) :: i
141 
142 
143  !$OMP PARALLEL
144  !$OMP DO
145  do i = 1, n
146  y(i) = alpha * x(i) + y(i)
147  end do
148  !$OMP END DO
149  !$OMP END PARALLEL
150 
151  end subroutine hecmw_axpy_r
152 
153  !C
154  !C***
155  !C*** hecmw_xpay_R
156  !C***
157  !C
158  subroutine hecmw_xpay_r (n, alpha, X, Y)
159 
160  use hecmw_util
161 
162  implicit none
163  integer(kind=kint) :: n
164  real(kind=kreal) :: x(:), y(:)
165  real(kind=kreal) :: alpha
166 
167  integer(kind=kint) :: i
168 
169  !$OMP PARALLEL
170  !$OMP DO
171  do i = 1, n
172  y(i) = x(i) + alpha * y(i)
173  end do
174  !$OMP END DO
175  !$OMP END PARALLEL
176 
177  end subroutine hecmw_xpay_r
178 
179  !C
180  !C***
181  !C*** hecmw_axpby_R
182  !C***
183  !C
184  subroutine hecmw_axpby_r (n, alpha, beta, X, Y)
185 
186  use hecmw_util
187 
188  implicit none
189  integer(kind=kint) :: n
190  real(kind=kreal) :: x(:), y(:)
191  real(kind=kreal) :: alpha, beta
192 
193  integer(kind=kint) :: i
194 
195 
196  if (beta == 0.d0) then
197  !$OMP PARALLEL
198  !$OMP DO
199  do i = 1, n
200  y(i) = alpha * x(i)
201  end do
202  !$OMP END DO
203  !$OMP END PARALLEL
204  else
205  !$OMP PARALLEL
206  !$OMP DO
207  do i = 1, n
208  y(i) = alpha * x(i) + beta * y(i)
209  end do
210  !$OMP END DO
211  !$OMP END PARALLEL
212  end if
213 
214  end subroutine hecmw_axpby_r
215 
216  !C
217  !C***
218  !C*** hecmw_axpyz_R
219  !C***
220  !C
221  subroutine hecmw_axpyz_r (n, alpha, X, Y, Z)
222 
223  use hecmw_util
224 
225  implicit none
226  integer(kind=kint) :: n
227  real(kind=kreal) :: x(:), y(:), z(:)
228  real(kind=kreal) :: alpha
229 
230  integer(kind=kint) :: i
231 
232 
233  !$OMP PARALLEL
234  !$OMP DO
235  do i = 1, n
236  z(i) = alpha * x(i) + y(i)
237  end do
238  !$OMP END DO
239  !$OMP END PARALLEL
240  end subroutine hecmw_axpyz_r
241 
242  !C
243  !C***
244  !C*** hecmw_scale_R
245  !C***
246  !C
247  subroutine hecmw_scale_r (n, alpha, X)
248 
249  use hecmw_util
250 
251  implicit none
252  integer(kind=kint) :: n
253  real(kind=kreal) :: x(:)
254  real(kind=kreal) :: alpha
255 
256  integer(kind=kint) :: i
257 
258  if (alpha == 0.d0) then
259  !$OMP PARALLEL
260  !$OMP DO
261  do i = 1, n
262  x(i) = 0.d0
263  end do
264  !$OMP END DO
265  !$OMP END PARALLEL
266  else
267  !$OMP PARALLEL
268  !$OMP DO
269  do i = 1, n
270  x(i) = alpha * x(i)
271  end do
272  !$OMP END DO
273  !$OMP END PARALLEL
274  end if
275 
276  end subroutine hecmw_scale_r
277 
278  !C
279  !C***
280  !C*** hecmw_copy_R
281  !C***
282  !C
283  subroutine hecmw_copy_r (n, X, Y)
284 
285  use hecmw_util
286 
287  implicit none
288  integer(kind=kint) :: n
289  real(kind=kreal) :: x(:), y(:)
290 
291  integer(kind=kint) :: i
292 
293 
294  !$OMP PARALLEL
295  !$OMP DO
296  do i = 1, n
297  y(i) = x(i)
298  end do
299  !$OMP END DO
300  !$OMP END PARALLEL
301 
302  end subroutine hecmw_copy_r
303 
304  !C
305  !C***
306  !C*** hecmw_time_statistics
307  !C***
308  !C
309  subroutine hecmw_time_statistics (hecMESH, time, &
310  t_max, t_min, t_avg, t_sd)
311  use hecmw_util
312  use m_hecmw_comm_f
313 
314  implicit none
315  type (hecmwST_local_mesh), intent(in) :: hecMESH
316  real(kind=kreal), intent(in) :: time
317  real(kind=kreal), intent(out) :: t_max
318  real(kind=kreal), intent(out), optional :: t_min, t_avg, t_sd
319  real(kind=kreal) :: t2_avg
320  integer(kind=kint) :: nprocs
321 
322  nprocs = hecmw_comm_get_size()
323 
324  t_max = time
325  call hecmw_allreduce_r1(hecmesh, t_max, hecmw_max)
326 
327  if (.not. present(t_min)) return
328  t_min = time
329  call hecmw_allreduce_r1(hecmesh, t_min, hecmw_min)
330 
331  if (.not. present(t_avg)) return
332  t_avg = time
333  call hecmw_allreduce_r1(hecmesh, t_avg, hecmw_sum)
334  t_avg = t_avg / nprocs
335 
336  if (.not. present(t_sd)) return
337  t2_avg = time*time
338  call hecmw_allreduce_r1(hecmesh, t2_avg, hecmw_sum)
339  t2_avg = t2_avg / nprocs
340 
341  t_sd = dsqrt(t2_avg - t_avg*t_avg)
342  end subroutine hecmw_time_statistics
343 
344 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)