FrontISTR  5.7.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 )
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 )
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)
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)
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 (hecMESH, ndof, alpha, X, Y)
132 
133  use hecmw_util
134 
135  implicit none
136  type (hecmwST_local_mesh) :: hecMESH
137  integer(kind=kint) :: ndof
138  real(kind=kreal) :: x(:), y(:)
139  real(kind=kreal) :: alpha
140 
141  integer(kind=kint) :: i
142 
143 
144  !$OMP PARALLEL
145  !$OMP DO
146  do i = 1, hecmesh%nn_internal * ndof
147  y(i) = alpha * x(i) + y(i)
148  end do
149  !$OMP END DO
150  !$OMP END PARALLEL
151 
152  end subroutine hecmw_axpy_r
153 
154  !C
155  !C***
156  !C*** hecmw_xpay_R
157  !C***
158  !C
159  subroutine hecmw_xpay_r (hecMESH, ndof, alpha, X, Y)
160 
161  use hecmw_util
162 
163  implicit none
164  type (hecmwST_local_mesh) :: hecMESH
165  integer(kind=kint) :: ndof
166  real(kind=kreal) :: x(:), y(:)
167  real(kind=kreal) :: alpha
168 
169  integer(kind=kint) :: i
170 
171  !$OMP PARALLEL
172  !$OMP DO
173  do i = 1, hecmesh%nn_internal * ndof
174  y(i) = x(i) + alpha * y(i)
175  end do
176  !$OMP END DO
177  !$OMP END PARALLEL
178 
179  end subroutine hecmw_xpay_r
180 
181  !C
182  !C***
183  !C*** hecmw_axpyz_R
184  !C***
185  !C
186  subroutine hecmw_axpyz_r (hecMESH, ndof, alpha, X, Y, Z)
187 
188  use hecmw_util
189 
190  implicit none
191  type (hecmwST_local_mesh) :: hecMESH
192  integer(kind=kint) :: ndof
193  real(kind=kreal) :: x(:), y(:), z(:)
194  real(kind=kreal) :: alpha
195 
196  integer(kind=kint) :: i
197 
198 
199  !$OMP PARALLEL
200  !$OMP DO
201  do i = 1, hecmesh%nn_internal * ndof
202  z(i) = alpha * x(i) + y(i)
203  end do
204  !$OMP END DO
205  !$OMP END PARALLEL
206  end subroutine hecmw_axpyz_r
207 
208  !C
209  !C***
210  !C*** hecmw_scale_R
211  !C***
212  !C
213  subroutine hecmw_scale_r (hecMESH, ndof, alpha, X)
214 
215  use hecmw_util
216 
217  implicit none
218  type (hecmwST_local_mesh) :: hecMESH
219  integer(kind=kint) :: ndof
220  real(kind=kreal) :: x(:)
221  real(kind=kreal) :: alpha
222 
223  integer(kind=kint) :: i
224 
225  !$OMP PARALLEL
226  !$OMP DO
227  do i = 1, hecmesh%nn_internal * ndof
228  x(i) = alpha * x(i)
229  end do
230  !$OMP END DO
231  !$OMP END PARALLEL
232 
233  end subroutine hecmw_scale_r
234 
235  !C
236  !C***
237  !C*** hecmw_copy_R
238  !C***
239  !C
240  subroutine hecmw_copy_r (hecMESH, ndof, X, Y)
241 
242  use hecmw_util
243 
244  implicit none
245  type (hecmwST_local_mesh) :: hecMESH
246  integer(kind=kint) :: ndof
247  real(kind=kreal) :: x(:), y(:)
248 
249  integer(kind=kint) :: i
250 
251 
252  !$OMP PARALLEL
253  !$OMP DO
254  do i = 1, hecmesh%nn_internal * ndof
255  y(i) = x(i)
256  end do
257  !$OMP END DO
258  !$OMP END PARALLEL
259 
260  end subroutine hecmw_copy_r
261 
262  !C
263  !C***
264  !C*** hecmw_time_statistics
265  !C***
266  !C
267  subroutine hecmw_time_statistics (hecMESH, time, &
268  t_max, t_min, t_avg, t_sd)
270  use m_hecmw_comm_f
271 
272  implicit none
273  type (hecmwST_local_mesh), intent(in) :: hecMESH
274  real(kind=kreal), intent(in) :: time
275  real(kind=kreal), intent(out) :: t_max
276  real(kind=kreal), intent(out), optional :: t_min, t_avg, t_sd
277  real(kind=kreal) :: t2_avg
278  integer(kind=kint) :: nprocs
279 
280  nprocs = hecmw_comm_get_size()
281 
282  t_max = time
283  call hecmw_allreduce_r1(hecmesh, t_max, hecmw_max)
284 
285  if (.not. present(t_min)) return
286  t_min = time
287  call hecmw_allreduce_r1(hecmesh, t_min, hecmw_min)
288 
289  if (.not. present(t_avg)) return
290  t_avg = time
291  call hecmw_allreduce_r1(hecmesh, t_avg, hecmw_sum)
292  t_avg = t_avg / nprocs
293 
294  if (.not. present(t_sd)) return
295  t2_avg = time*time
296  call hecmw_allreduce_r1(hecmesh, t2_avg, hecmw_sum)
297  t2_avg = t2_avg / nprocs
298 
299  t_sd = dsqrt(t2_avg - t_avg*t_avg)
300  end subroutine hecmw_time_statistics
301 
302 end module hecmw_solver_misc
hecmw_util::hecmw_max
integer(kind=kint), parameter hecmw_max
Definition: hecmw_util_f.F90:25
hecmw_util::hecmw_wtime
real(kind=kreal) function hecmw_wtime()
Definition: hecmw_util_f.F90:549
hecmw_solver_misc::hecmw_absmax_r
subroutine hecmw_absmax_r(hecMESH, ndof, X, absMax, COMMtime)
Definition: hecmw_solver_misc.f90:78
hecmw_util::hecmw_sum
integer(kind=kint), parameter hecmw_sum
Definition: hecmw_util_f.F90:23
hecmw_solver_misc::hecmw_innerproduct_i
subroutine hecmw_innerproduct_i(hecMESH, ndof, X, Y, sum, COMMtime)
Definition: hecmw_solver_misc.f90:16
hecmw_solver_misc::hecmw_copy_r
subroutine hecmw_copy_r(hecMESH, ndof, X, Y)
Definition: hecmw_solver_misc.f90:241
hecmw_solver_misc::hecmw_innerproduct_r
subroutine hecmw_innerproduct_r(hecMESH, ndof, X, Y, sum, COMMtime)
Definition: hecmw_solver_misc.f90:47
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
m_hecmw_comm_f::hecmw_allreduce_r1
subroutine hecmw_allreduce_r1(hecMESH, s, ntag)
Definition: hecmw_comm_f.F90:403
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
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
hecmw_solver_misc::hecmw_xpay_r
subroutine hecmw_xpay_r(hecMESH, ndof, alpha, X, Y)
Definition: hecmw_solver_misc.f90:160
m_hecmw_comm_f::hecmw_allreduce_i1
subroutine hecmw_allreduce_i1(hecMESH, s, ntag)
Definition: hecmw_comm_f.F90:458
hecmw_solver_misc::hecmw_innerproduct_r_nocomm
subroutine hecmw_innerproduct_r_nocomm(hecMESH, ndof, X, Y, sum)
Definition: hecmw_solver_misc.f90:109
hecmw_solver_misc::hecmw_time_statistics
subroutine hecmw_time_statistics(hecMESH, time, t_max, t_min, t_avg, t_sd)
Definition: hecmw_solver_misc.f90:269
hecmw_solver_misc::hecmw_axpy_r
subroutine hecmw_axpy_r(hecMESH, ndof, alpha, X, Y)
Definition: hecmw_solver_misc.f90:132
hecmw_solver_misc::hecmw_scale_r
subroutine hecmw_scale_r(hecMESH, ndof, alpha, X)
Definition: hecmw_solver_misc.f90:214
hecmw_solver_misc::hecmw_axpyz_r
subroutine hecmw_axpyz_r(hecMESH, ndof, alpha, X, Y, Z)
Definition: hecmw_solver_misc.f90:187
hecmw_util::hecmw_comm_get_size
integer(kind=kint) function hecmw_comm_get_size()
Definition: hecmw_util_f.F90:593
hecmw_solver_misc
Definition: hecmw_solver_misc.f90:6