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 +----------------------+
102  !C | SETUP PRECONDITIONER |
103  !C +----------------------+
104  !C===
105  call hecmw_precond_setup(hecmat, hecmesh, 0)
106 
107  !C===
108  !C +---------------------+
109  !C | {r0}= {b} - [A]{x0} |
110  !C +---------------------+
111  !C===
112  call hecmw_matresid(hecmesh, hecmat, x, b, ww(:,r), tcomm)
113 
114  !C-- set arbitrary {r_tld}
115  call hecmw_copy_r(nndof, ww(:,r), ww(:,rt))
116 
117  !C-- compute ||{b}||
118  call hecmw_innerproduct_r(hecmesh, ndof, b, b, bnrm2, tcomm)
119  if (bnrm2.eq.0.d0) then
120  iter = 0
121  maxit = 0
122  resid = 0.d0
123  x = 0.d0
124  endif
125 
126  e_time = hecmw_wtime()
127  if (timelog.eq.2) then
128  call hecmw_time_statistics(hecmesh, e_time - s_time, &
129  t_max, t_min, t_avg, t_sd)
130  if (hecmesh%my_rank.eq.0) then
131  write(*,*) 'Time solver setup'
132  write(*,*) ' Max :',t_max
133  write(*,*) ' Min :',t_min
134  write(*,*) ' Avg :',t_avg
135  write(*,*) ' Std Dev :',t_sd
136  endif
137  tset = t_max
138  else
139  tset = e_time - s_time
140  endif
141 
142  tcomm = 0.d0
143  call hecmw_barrier(hecmesh)
144  s1_time = hecmw_wtime()
145  !C
146  !C*************************************************************** iterative procedures start
147  !C
148  do iter = 1, maxit
149 
150  !C===
151  !C +-----------------+
152  !C | RHO= {r}{r_tld} |
153  !C +-----------------+
154  !C===
155  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,r), ww(:,rt), rho, tcomm)
156 
157  !C===
158  !C +----------------------------------------+
159  !C | BETA= (RHO/RHO1) * (ALPHA/OMEGA) |
160  !C | {p} = {r} + BETA * ( {p} - OMEGA*{v} ) |
161  !C +----------------------------------------+
162  !C===
163  if ( iter.gt.1 ) then
164  beta = (rho/rho1) * (alpha/omega)
165  call hecmw_axpy_r(nndof, -omega, ww(:,v), ww(:,p))
166  call hecmw_xpay_r(nndof, beta, ww(:,r), ww(:,p))
167  else
168  call hecmw_copy_r(nndof, ww(:,r), ww(:,p))
169  endif
170 
171  !C===
172  !C +--------------------+
173  !C | {p_tld}= [Minv]{p} |
174  !C +--------------------+
175  !C===
176  call hecmw_precond_apply(hecmesh, hecmat, ww(:, p), ww(:, pt), ww(:, wk), tcomm)
177 
178  !C===
179  !C +-------------------------+
180  !C | {v}= [A] {p_tld} |
181  !C +-------------------------+
182  !C===
183  call hecmw_matvec(hecmesh, hecmat, ww(:,pt), ww(:,v), tcomm)
184 
185  !C
186  !C-- calc. ALPHA
187  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,rt), ww(:,v), c2, tcomm)
188 
189  alpha = rho / c2
190 
191  !C
192  !C-- {s}= {r} - ALPHA*{V}
193  call hecmw_axpyz_r(nndof, -alpha, ww(:,v), ww(:,r), ww(:,s))
194 
195  !C===
196  !C +--------------------+
197  !C | {s_tld}= [Minv]{s} |
198  !C +--------------------+
199  !C===
200  call hecmw_precond_apply(hecmesh, hecmat, ww(:, s), ww(:, st), ww(:, wk), tcomm)
201 
202  !C===
203  !C +-------------------------+
204  !C | {t}= [A] {s_tld} |
205  !C +-------------------------+
206  !C===
207  call hecmw_matvec(hecmesh, hecmat, ww(:,st), ww(:,t), tcomm)
208 
209  !C===
210  !C +----------------------------+
211  !C | OMEGA= ({t}{s}) / ({t}{t}) |
212  !C +----------------------------+
213  !C===
214  call hecmw_innerproduct_r_nocomm(hecmesh, ndof, ww(:,t), ww(:,s), cg(1))
215  call hecmw_innerproduct_r_nocomm(hecmesh, ndof, ww(:,t), ww(:,t), cg(2))
216  s_time= hecmw_wtime()
217  call hecmw_allreduce_r(hecmesh, cg, 2, hecmw_sum)
218  e_time= hecmw_wtime()
219  tcomm = tcomm + e_time - s_time
220 
221  omega = cg(1) / cg(2)
222 
223  !C===
224  !C +----------------+
225  !C | update {x},{r} |
226  !C +----------------+
227  !C===
228  call hecmw_axpy_r(nndof, alpha, ww(:,pt), x)
229  call hecmw_axpy_r(nndof, omega, ww(:,st), x)
230  !C
231  !C--- recompute R sometimes
232  if ( mod(iter,n_iter_recompute_r)==0 ) then
233  call hecmw_matresid(hecmesh, hecmat, x, b, ww(:,r), tcomm)
234  else
235  call hecmw_axpyz_r(nndof, -omega, ww(:,t), ww(:,s), ww(:,r))
236  endif
237 
238  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,r), ww(:,r), dnrm2, tcomm)
239 
240  resid= dsqrt(dnrm2/bnrm2)
241 
242  !C##### ITERATION HISTORY
243  if (my_rank.eq.0.and.iterlog.eq.1) write (*,'(i7, 1pe16.6)') iter, resid
244  !C#####
245 
246  if ( resid.le.tol ) then
247  if ( mod(iter,n_iter_recompute_r)==0 ) exit
248  !C----- recompute R to make sure it is really converged
249  call hecmw_matresid(hecmesh, hecmat, x, b, ww(:,r), tcomm)
250  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,r), ww(:,r), dnrm2, tcomm)
251  resid= dsqrt(dnrm2/bnrm2)
252  if ( resid.le.tol ) exit
253  endif
254  if ( iter .eq.maxit ) error = hecmw_solver_error_noconv_maxit
255 
256  rho1 = rho
257 
258  enddo
259  !C
260  !C*************************************************************** iterative procedures end
261  !C
262 
263  call hecmw_solver_scaling_bk(hecmat)
264  !C
265  !C-- INTERFACE data EXCHANGE
266  !C
267  start_time = hecmw_wtime()
268  call hecmw_update_r (hecmesh, x, hecmat%NP, hecmat%NDOF)
269  end_time = hecmw_wtime()
270  tcomm = tcomm + end_time - start_time
271 
272  deallocate (ww)
273  !call hecmw_precond_clear(hecMAT)
274 
275  e1_time = hecmw_wtime()
276  if (timelog.eq.2) then
277  call hecmw_time_statistics(hecmesh, e1_time - s1_time, &
278  t_max, t_min, t_avg, t_sd)
279  if (hecmesh%my_rank.eq.0) then
280  write(*,*) 'Time solver iterations'
281  write(*,*) ' Max :',t_max
282  write(*,*) ' Min :',t_min
283  write(*,*) ' Avg :',t_avg
284  write(*,*) ' Std Dev :',t_sd
285  endif
286  tsol = t_max
287  else
288  tsol = e1_time - s1_time
289  endif
290 
291  end subroutine hecmw_solve_bicgstab
292 end module hecmw_solver_bicgstab
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