FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
hecmw_solver_BiCGSTAB.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 
6 !C
7 !C***
8 !C*** module hecmw_solver_BiCGSTAB
9 !C***
10 !C
12 contains
13  !C
14  !C*** BiCGSTAB_3
15  !C
16  subroutine hecmw_solve_bicgstab( hecMESH, hecMAT, ITER, RESID, error, &
17  & Tset, Tsol, Tcomm )
18 
19  use hecmw_util
21  use m_hecmw_comm_f
26  use hecmw_precond
27 
28  implicit none
29 
30  type(hecmwst_local_mesh) :: hecMESH
31  type(hecmwst_matrix) :: hecMAT
32  integer(kind=kint ), intent(inout):: ITER, error
33  real (kind=kreal), intent(inout):: resid, tset, tsol, tcomm
34 
35  integer(kind=kint ) :: N, NP, NDOF, NNDOF
36  integer(kind=kint ) :: my_rank
37  integer(kind=kint ) :: ITERlog, TIMElog
38  real(kind=kreal), pointer :: b(:), x(:)
39 
40  real(kind=kreal), dimension(:,:), allocatable :: ww
41  real(kind=kreal), dimension(2) :: cg
42 
43  integer(kind=kint ) :: MAXIT
44 
45  ! local variables
46  real (kind=kreal):: tol
47  integer(kind=kint )::i
48  real (kind=kreal)::s_time,s1_time,e_time,e1_time, start_time, end_time
49  real (kind=kreal)::bnrm2,c2
50  real (kind=kreal)::rho,rho1,beta,alpha,dnrm2
51  real (kind=kreal)::omega
52  real (kind=kreal)::t_max,t_min,t_avg,t_sd
53 
54  integer(kind=kint), parameter :: R = 1
55  integer(kind=kint), parameter :: RT= 2
56  integer(kind=kint), parameter :: P = 3
57  integer(kind=kint), parameter :: PT= 4
58  integer(kind=kint), parameter :: S = 5
59  integer(kind=kint), parameter :: ST= 1
60  integer(kind=kint), parameter :: T = 6
61  integer(kind=kint), parameter :: V = 7
62  integer(kind=kint), parameter :: WK= 8
63 
64  integer(kind=kint), parameter :: N_ITER_RECOMPUTE_R= 100
65 
66  call hecmw_barrier(hecmesh)
67  s_time= hecmw_wtime()
68 
69  !C===
70  !C +-------+
71  !C | INIT. |
72  !C +-------+
73  !C===
74  n = hecmat%N
75  np = hecmat%NP
76  ndof = hecmat%NDOF
77  nndof = n * ndof
78  my_rank = hecmesh%my_rank
79  x => hecmat%X
80  b => hecmat%B
81 
82  iterlog = hecmw_mat_get_iterlog( hecmat )
83  timelog = hecmw_mat_get_timelog( hecmat )
84  maxit = hecmw_mat_get_iter( hecmat )
85  tol = hecmw_mat_get_resid( hecmat )
86 
87  error = 0
88  rho1 = 0.0d0
89  alpha = 0.0d0
90  beta = 0.0d0
91  omega = 0.0d0
92 
93  allocate (ww(ndof*np, 8))
94  ww = 0.d0
95 
96  !C
97  !C-- SCALING
98  call hecmw_solver_scaling_fw(hecmesh, hecmat, tcomm)
99 
100  !C
101  !C-- matrix integration for OpenACC
102  !C
103  !C @note:
104  !C Combine hecMAT%AL, D, and AU into a single matrix for GPU execution.
105  !C This is a no-op for CPU builds.
106  call hecmw_mat_integrate(hecmat)
107 
108  !C===
109  !C +----------------------+
110  !C | SETUP PRECONDITIONER |
111  !C +----------------------+
112  !C===
113  call hecmw_precond_setup(hecmat, hecmesh, 0)
114 
115  !C===
116  !C +---------------------+
117  !C | {r0}= {b} - [A]{x0} |
118  !C +---------------------+
119  !C===
120  call hecmw_matresid(hecmesh, hecmat, x, b, ww(:,r), tcomm)
121 
122  !C-- set arbitrary {r_tld}
123  call hecmw_copy_r(nndof, ww(:,r), ww(:,rt))
124 
125  !C-- compute ||{b}||
126  call hecmw_innerproduct_r(hecmesh, ndof, b, b, bnrm2, tcomm)
127  if (bnrm2.eq.0.d0) then
128  iter = 0
129  maxit = 0
130  resid = 0.d0
131  x = 0.d0
132  endif
133 
134  e_time = hecmw_wtime()
135  if (timelog.eq.2) then
136  call hecmw_time_statistics(hecmesh, e_time - s_time, &
137  t_max, t_min, t_avg, t_sd)
138  if (hecmesh%my_rank.eq.0) then
139  write(*,*) 'Time solver setup'
140  write(*,*) ' Max :',t_max
141  write(*,*) ' Min :',t_min
142  write(*,*) ' Avg :',t_avg
143  write(*,*) ' Std Dev :',t_sd
144  endif
145  tset = t_max
146  else
147  tset = e_time - s_time
148  endif
149 
150  tcomm = 0.d0
151  call hecmw_barrier(hecmesh)
152  s1_time = hecmw_wtime()
153  !C
154  !C*************************************************************** iterative procedures start
155  !C
156  do iter = 1, maxit
157 
158  !C===
159  !C +-----------------+
160  !C | RHO= {r}{r_tld} |
161  !C +-----------------+
162  !C===
163  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,r), ww(:,rt), rho, tcomm)
164 
165  !C===
166  !C +----------------------------------------+
167  !C | BETA= (RHO/RHO1) * (ALPHA/OMEGA) |
168  !C | {p} = {r} + BETA * ( {p} - OMEGA*{v} ) |
169  !C +----------------------------------------+
170  !C===
171  if ( iter.gt.1 ) then
172  beta = (rho/rho1) * (alpha/omega)
173  call hecmw_axpy_r(nndof, -omega, ww(:,v), ww(:,p))
174  call hecmw_xpay_r(nndof, beta, ww(:,r), ww(:,p))
175  else
176  call hecmw_copy_r(nndof, ww(:,r), ww(:,p))
177  endif
178 
179  !C===
180  !C +--------------------+
181  !C | {p_tld}= [Minv]{p} |
182  !C +--------------------+
183  !C===
184  call hecmw_precond_apply(hecmesh, hecmat, ww(:, p), ww(:, pt), ww(:, wk), tcomm)
185 
186  !C===
187  !C +-------------------------+
188  !C | {v}= [A] {p_tld} |
189  !C +-------------------------+
190  !C===
191  call hecmw_matvec(hecmesh, hecmat, ww(:,pt), ww(:,v), tcomm)
192 
193  !C
194  !C-- calc. ALPHA
195  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,rt), ww(:,v), c2, tcomm)
196 
197  alpha = rho / c2
198 
199  !C
200  !C-- {s}= {r} - ALPHA*{V}
201  call hecmw_axpyz_r(nndof, -alpha, ww(:,v), ww(:,r), ww(:,s))
202 
203  !C===
204  !C +--------------------+
205  !C | {s_tld}= [Minv]{s} |
206  !C +--------------------+
207  !C===
208  call hecmw_precond_apply(hecmesh, hecmat, ww(:, s), ww(:, st), ww(:, wk), tcomm)
209 
210  !C===
211  !C +-------------------------+
212  !C | {t}= [A] {s_tld} |
213  !C +-------------------------+
214  !C===
215  call hecmw_matvec(hecmesh, hecmat, ww(:,st), ww(:,t), tcomm)
216 
217  !C===
218  !C +----------------------------+
219  !C | OMEGA= ({t}{s}) / ({t}{t}) |
220  !C +----------------------------+
221  !C===
222  call hecmw_innerproduct_r_nocomm(hecmesh, ndof, ww(:,t), ww(:,s), cg(1))
223  call hecmw_innerproduct_r_nocomm(hecmesh, ndof, ww(:,t), ww(:,t), cg(2))
224  s_time= hecmw_wtime()
225  call hecmw_allreduce_r(hecmesh, cg, 2, hecmw_sum)
226  e_time= hecmw_wtime()
227  tcomm = tcomm + e_time - s_time
228 
229  omega = cg(1) / cg(2)
230 
231  !C===
232  !C +----------------+
233  !C | update {x},{r} |
234  !C +----------------+
235  !C===
236  call hecmw_axpy_r(nndof, alpha, ww(:,pt), x)
237  call hecmw_axpy_r(nndof, omega, ww(:,st), x)
238  !C
239  !C--- recompute R sometimes
240  if ( mod(iter,n_iter_recompute_r)==0 ) then
241  call hecmw_matresid(hecmesh, hecmat, x, b, ww(:,r), tcomm)
242  else
243  call hecmw_axpyz_r(nndof, -omega, ww(:,t), ww(:,s), ww(:,r))
244  endif
245 
246  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,r), ww(:,r), dnrm2, tcomm)
247 
248  resid= dsqrt(dnrm2/bnrm2)
249 
250  !C##### ITERATION HISTORY
251  if (my_rank.eq.0.and.iterlog.eq.1) write (*,'(i7, 1pe16.6)') iter, resid
252  !C#####
253 
254  if ( resid.le.tol ) then
255  if ( mod(iter,n_iter_recompute_r)==0 ) exit
256  !C----- recompute R to make sure it is really converged
257  call hecmw_matresid(hecmesh, hecmat, x, b, ww(:,r), tcomm)
258  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,r), ww(:,r), dnrm2, tcomm)
259  resid= dsqrt(dnrm2/bnrm2)
260  if ( resid.le.tol ) exit
261  endif
262  if ( iter .eq.maxit ) error = hecmw_solver_error_noconv_maxit
263 
264  rho1 = rho
265 
266  enddo
267  !C
268  !C*************************************************************** iterative procedures end
269  !C
270 
271  call hecmw_solver_scaling_bk(hecmat)
272  !C
273  !C-- INTERFACE data EXCHANGE
274  !C
275  start_time = hecmw_wtime()
276  call hecmw_update_r (hecmesh, x, hecmat%NP, hecmat%NDOF)
277  end_time = hecmw_wtime()
278  tcomm = tcomm + end_time - start_time
279 
280  deallocate (ww)
281  !call hecmw_precond_clear(hecMAT)
282 
283  e1_time = hecmw_wtime()
284  if (timelog.eq.2) then
285  call hecmw_time_statistics(hecmesh, e1_time - s1_time, &
286  t_max, t_min, t_avg, t_sd)
287  if (hecmesh%my_rank.eq.0) then
288  write(*,*) 'Time solver iterations'
289  write(*,*) ' Max :',t_max
290  write(*,*) ' Min :',t_min
291  write(*,*) ' Avg :',t_avg
292  write(*,*) ' Std Dev :',t_sd
293  endif
294  tsol = t_max
295  else
296  tsol = e1_time - s1_time
297  endif
298 
299  end subroutine hecmw_solve_bicgstab
300 end module hecmw_solver_bicgstab
subroutine, public hecmw_mat_integrate(hecMAT)
Integrate matrix components into a single array for efficient access.
real(kind=kreal) function, public hecmw_mat_get_resid(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_iterlog(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_timelog(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_iter(hecMAT)
subroutine, public hecmw_precond_setup(hecMAT, hecMESH, sym)
subroutine, public hecmw_precond_apply(hecMESH, hecMAT, R, Z, ZP, COMMtime)
subroutine hecmw_solve_bicgstab(hecMESH, hecMAT, ITER, RESID, error, Tset, Tsol, Tcomm)
subroutine, public hecmw_matresid(hecMESH, hecMAT, X, B, R, COMMtime)
subroutine, public hecmw_matvec(hecMESH, hecMAT, X, Y, COMMtime)
subroutine hecmw_xpay_r(n, alpha, X, Y)
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_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, public hecmw_solver_scaling_fw(hecMESH, hecMAT, COMMtime)
subroutine, public hecmw_solver_scaling_bk(hecMAT)
I/O and Utility.
Definition: hecmw_util_f.F90:7
integer(kind=kint), parameter hecmw_sum
integer(kind=4), parameter kreal
real(kind=kreal) function hecmw_wtime()
subroutine hecmw_update_r(hecMESH, val, n, m)
subroutine hecmw_allreduce_r(hecMESH, val, n, ntag)
subroutine hecmw_barrier(hecMESH)
integer(kind=kint), parameter hecmw_solver_error_noconv_maxit