FrontISTR  5.8.0
Large-scale structural analysis program with finit element method
hecmw_solver_GMRES.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*** module hecmw_solver_GMRES
8 !C***
9 !
11 
12  public :: hecmw_solve_gmres
13 
14 contains
15  !C
16  !C*** hecmw_solve_GMRES
17  !C
18  subroutine hecmw_solve_gmres( hecMESH, hecMAT, ITER, RESID, error, &
19  & Tset, Tsol, Tcomm )
20  use hecmw_util
22  use m_hecmw_comm_f
27  use hecmw_precond
29 
30  implicit none
31 
32  type(hecmwst_local_mesh) :: hecmesh
33  type(hecmwst_matrix) :: hecmat
34  integer(kind=kint ), intent(inout):: iter, error
35  real (kind=kreal), intent(inout):: resid, tset, tsol, tcomm
36 
37  integer(kind=kint ) :: n, np, ndof, nndof
38  integer(kind=kint ) :: my_rank
39  integer(kind=kint ) :: iterlog, timelog
40  real(kind=kreal), pointer :: b(:), x(:)
41 
42  real(kind=kreal), dimension(:,:), allocatable :: ww
43 
44  integer(kind=kint ) :: maxit, nrest
45 
46  real (kind=kreal) :: tol
47 
48  real (kind=kreal), dimension(:), allocatable :: ss
49  real (kind=kreal), dimension(:,:), allocatable :: h
50 
51  integer(kind=kint ) :: cs, sn
52 
53  real (kind=kreal) zero, one
54  parameter( zero = 0.0d+0, one = 1.0d+0 )
55 
56  integer(kind=kint ) :: nrk,i,k,kk,jj,info,ik
57  integer(kind=kint ) :: irow
58  real (kind=kreal) :: s_time,e_time,s1_time,e1_time
59  real (kind=kreal) :: ldh,ldw,bnrm2,dnrm2,rnorm
60  real (kind=kreal) :: commtime,comptime, coef,val,vcs,vsn,dtemp,aa,bb,r0,scale,rr
61  integer(kind=kint ) :: estcond
62  real (kind=kreal) :: t_max,t_min,t_avg,t_sd
63 
64  integer(kind=kint), parameter :: r = 1
65  integer(kind=kint), parameter :: zp = r + 1
66  integer(kind=kint), parameter :: zq = r + 2
67  integer(kind=kint), parameter :: s = r + 3
68  integer(kind=kint), parameter :: w = s + 1
69  integer(kind=kint), parameter :: y = w
70  integer(kind=kint), parameter :: av = y + 1
71  integer(kind=kint), parameter :: v = av + 1
72 
73  call hecmw_barrier(hecmesh)
74  s_time= hecmw_wtime()
75  !C
76  !C-- INIT.
77  n = hecmat%N
78  np = hecmat%NP
79  ndof = hecmat%NDOF
80  nndof = n * ndof
81  my_rank = hecmesh%my_rank
82  x => hecmat%X
83  b => hecmat%B
84 
85  iterlog = hecmw_mat_get_iterlog( hecmat )
86  timelog = hecmw_mat_get_timelog( hecmat )
87  maxit = hecmw_mat_get_iter( hecmat )
88  tol = hecmw_mat_get_resid( hecmat )
89  nrest = hecmw_mat_get_nrest( hecmat )
90  estcond = hecmw_mat_get_estcond( hecmat )
91 
92  if (nrest >= ndof*np-1) nrest = ndof*np-2
93 
94  error= 0
95  nrk= nrest + 7
96 
97  allocate (h(nrk,nrk))
98  allocate (ww(ndof*np,nrk))
99  allocate (ss(nrk))
100 
101  commtime= 0.d0
102  comptime= 0.d0
103 
104  ldh= nrest + 2
105  ldw= n
106 
107  !C
108  !C-- Store the Givens parameters in matrix H.
109  cs= nrest + 1
110  sn= cs + 1
111 
112  !C
113  !C-- SCALING
114  call hecmw_solver_scaling_fw(hecmesh, hecmat, tcomm)
115 
116  !C===
117  !C +----------------------+
118  !C | SETUP PRECONDITIONER |
119  !C +----------------------+
120  !C===
121  call hecmw_precond_setup(hecmat, hecmesh, 0)
122 
123  !C
124  !C
125  !C +--------------------+
126  !C | {r}= {b} - [A]{x0} |
127  !C +--------------------+
128  !C===
129  call hecmw_matresid(hecmesh, hecmat, x, b, ww(:,r), tcomm)
130  !C===
131 
132  call hecmw_innerproduct_r(hecmesh, ndof, b, b, bnrm2, tcomm)
133  if (bnrm2.eq.0.d0) then
134  iter = 0
135  maxit = 0
136  resid = 0.d0
137  x = 0.d0
138  endif
139 
140  e_time= hecmw_wtime()
141  if (timelog.eq.2) then
142  call hecmw_time_statistics(hecmesh, e_time - s_time, &
143  t_max, t_min, t_avg, t_sd)
144  if (hecmesh%my_rank.eq.0) then
145  write(*,*) 'Time solver setup'
146  write(*,*) ' Max :',t_max
147  write(*,*) ' Min :',t_min
148  write(*,*) ' Avg :',t_avg
149  write(*,*) ' Std Dev :',t_sd
150  endif
151  tset = t_max
152  else
153  tset = e_time - s_time
154  endif
155  !C===
156 
157  call hecmw_barrier(hecmesh)
158  s1_time= hecmw_wtime()
159  iter= 0
160 
161  outer: do
162 
163  !C
164  !C************************************************ GMRES Iteration
165  !C
166  i= 0
167  !C
168  !C +---------------+
169  !C | {v1}= {r}/|r| |
170  !C +---------------+
171  !C===
172  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,r), ww(:,r), dnrm2, tcomm)
173  if (dnrm2 == 0.d0) exit ! converged
174 
175  rnorm= dsqrt(dnrm2)
176  coef= one/rnorm
177  call hecmw_axpby_r(nndof, coef, 0.d0, ww(:,r), ww(:,v))
178  !C===
179 
180  !C
181  !C +--------------+
182  !C | {s}= |r|{e1} |
183  !C +--------------+
184  !C===
185  call hecmw_scale_r(nndof, zero, ww(:,s))
186  ww(1 ,s) = rnorm
187  !C===
188 
189  !C************************************************ GMRES(m) restart
190  do i = 1, nrest
191  iter= iter + 1
192 
193  !C
194  !C +-------------------+
195  !C | {w}= [A][Minv]{v} |
196  !C +-------------------+
197  !C===
198  call hecmw_precond_apply(hecmesh, hecmat, ww(:,v+i-1), ww(:,zq), ww(:,zp), tcomm)
199 
200  call hecmw_matvec(hecmesh, hecmat, ww(:,zq), ww(:,w), tcomm)
201  !C===
202 
203  !C
204  !C +------------------------------+
205  !C | ORTH. BASIS by GRAMM-SCHMIDT |
206  !C +------------------------------+
207  !C Construct the I-th column of the upper Hessenberg matrix H
208  !C using the Gram-Schmidt process on V and W.
209  !C===
210  do k= 1, i
211  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,w), ww(:,v+k-1), val, tcomm)
212 
213  call hecmw_axpy_r(nndof, -val, ww(:,v+k-1), ww(:,w))
214  h(k,i)= val
215  enddo
216 
217  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,w), ww(:,w), val, tcomm)
218  if (val == 0.d0) exit ! converged
219 
220  h(i+1,i)= dsqrt(val)
221  coef= one / h(i+1,i)
222  call hecmw_axpby_r(nndof, coef, 0.d0, ww(:,w), ww(:,v+i+1-1))
223  !C===
224 
225  !C
226  !C +-----------------+
227  !C | GIVENS ROTARION |
228  !C +-----------------+
229  !C===
230 
231  !C
232  !C-- Plane Rotation
233  do k = 1, i-1
234  vcs= h(k,cs)
235  vsn= h(k,sn)
236  dtemp = vcs*h(k ,i) + vsn*h(k+1,i)
237  h(k+1,i)= vcs*h(k+1,i) - vsn*h(k ,i)
238  h(k ,i)= dtemp
239  enddo
240 
241  !C
242  !C-- Construct Givens Plane Rotation
243  aa = h(i ,i)
244  bb = h(i+1,i)
245  r0= bb
246  if (dabs(aa).gt.dabs(bb)) r0= aa
247  scale= dabs(aa) + dabs(bb)
248 
249  if (scale.ne.0.d0) then
250  rr= scale * dsqrt((aa/scale)**2+(bb/scale)**2)
251  rr= dsign(1.d0,r0)*rr
252  h(i,cs)= aa/rr
253  h(i,sn)= bb/rr
254  else
255  h(i,cs)= 1.d0
256  h(i,sn)= 0.d0
257  rr = 0.d0
258  endif
259 
260  !C
261  !C-- Plane Rotation
262  vcs= h(i,cs)
263  vsn= h(i,sn)
264  dtemp = vcs*h(i ,i) + vsn*h(i+1,i)
265  h(i+1,i)= vcs*h(i+1,i) - vsn*h(i ,i)
266  h(i ,i)= dtemp
267 
268  dtemp = vcs*ww(i ,s) + vsn*ww(i+1,s)
269  ww(i+1,s)= vcs*ww(i+1,s) - vsn*ww(i ,s)
270  ww(i ,s)= dtemp
271 
272  resid = dabs( ww(i+1,s))/dsqrt(bnrm2)
273 
274  if (my_rank.eq.0 .and. iterlog.eq.1) &
275  & write (*, '(2i8, 1pe16.6)') iter,i+1, resid
276 
277  if (estcond /= 0 .and. hecmesh%my_rank == 0) then
278  if (mod(iter,estcond) == 0) call hecmw_estimate_condition_gmres(i, h)
279  endif
280 
281  if ( resid.le.tol ) then
282  !C-- [H]{y}= {s_tld}
283  do ik= 1, i
284  ss(ik)= ww(ik,s)
285  enddo
286  irow= i
287  ww(irow,y)= ss(irow) / h(irow,irow)
288 
289  do kk= irow-1, 1, -1
290  do jj= irow, kk+1, -1
291  ss(kk)= ss(kk) - h(kk,jj)*ww(jj,y)
292  enddo
293  ww(kk,y)= ss(kk) / h(kk,kk)
294  enddo
295 
296  !C-- {x}= {x} + {y}{V}
297  call hecmw_scale_r(nndof, 0.d0, ww(:,av))
298 
299  jj= irow
300  do jj= 1, irow
301  call hecmw_axpy_r(nndof, ww(jj,y), ww(:,v+jj-1), ww(:,av))
302  enddo
303 
304  call hecmw_precond_apply(hecmesh, hecmat, ww(:,av), ww(:,zq), ww(:,zp), tcomm)
305 
306  call hecmw_axpy_r(nndof, 1.d0, ww(:,zq), x)
307 
308  exit outer
309  endif
310 
311  if ( iter.gt.maxit ) then
313  exit outer
314  end if
315  end do
316  !C===
317 
318  !C
319  !C +------------------+
320  !C | CURRENT SOLUTION |
321  !C +------------------+
322  !C===
323 
324  !C-- [H]{y}= {s_tld}
325  do ik= 1, nrest
326  ss(ik)= ww(ik,s)
327  enddo
328  irow= nrest
329  ww(irow,y)= ss(irow) / h(irow,irow)
330 
331  do kk= irow-1, 1, -1
332  do jj= irow, kk+1, -1
333  ss(kk)= ss(kk) - h(kk,jj)*ww(jj,y)
334  enddo
335  ww(kk,y)= ss(kk) / h(kk,kk)
336  enddo
337 
338  !C-- {x}= {x} + {y}{V}
339  call hecmw_scale_r(nndof, 0.d0, ww(:,av))
340 
341  jj= irow
342  do jj= 1, irow
343  call hecmw_axpy_r(nndof, ww(jj,y), ww(:,v+jj-1), ww(:,av))
344  enddo
345 
346  call hecmw_precond_apply(hecmesh, hecmat, ww(:,av), ww(:,zq), ww(:,zp), tcomm)
347 
348  call hecmw_axpy_r(nndof, 1.d0, ww(:,zq), x)
349 
350  !C
351  !C-- Compute residual vector R, find norm, then check for tolerance.
352  call hecmw_matresid(hecmesh, hecmat, x, b, ww(:,r), tcomm)
353 
354  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,r), ww(:,r), dnrm2, tcomm)
355 
356  ww(i+1,s)= dsqrt(dnrm2/bnrm2)
357  resid = ww( i+1,s )
358 
359  if ( resid.le.tol ) exit outer
360  if ( iter .gt.maxit ) then
362  exit outer
363  end if
364  !C
365  !C-- RESTART
366  end do outer
367 
368  !C
369  !C-- iteration FAILED
370 
371  if (error == hecmw_solver_error_noconv_maxit) then
372  info = iter
373 
374  !C-- [H]{y}= {s_tld}
375  do ik= 1, i
376  ss(ik)= ww(ik,s)
377  enddo
378  irow= i
379  ww(irow,y)= ss(irow) / h(irow,irow)
380 
381  do kk= irow-1, 1, -1
382  do jj= irow, kk+1, -1
383  ss(kk)= ss(kk) - h(kk,jj)*ww(jj,y)
384  enddo
385  ww(kk,y)= ss(kk) / h(kk,kk)
386  enddo
387 
388  !C-- {x}= {x} + {y}{V}
389  call hecmw_scale_r(nndof, 0.d0, ww(:,av))
390 
391  jj= irow
392  do jj= 1, irow
393  call hecmw_axpy_r(nndof, ww(jj,y), ww(:,v+jj-1), ww(:,av))
394  enddo
395 
396  call hecmw_precond_apply(hecmesh, hecmat, ww(:,av), ww(:,zq), ww(:,zp), tcomm)
397 
398  call hecmw_axpy_r(nndof, 1.d0, ww(:,zq), x)
399  end if
400 
401  call hecmw_solver_scaling_bk(hecmat)
402 
403  if (estcond /= 0 .and. hecmesh%my_rank == 0) then
405  endif
406  !C
407  !C-- INTERFACE data EXCHANGE
408  s_time = hecmw_wtime()
409  call hecmw_update_r (hecmesh, x, hecmat%NP, hecmat%NDOF)
410  e_time = hecmw_wtime()
411  tcomm = tcomm + e_time - s_time
412 
413  deallocate (h, ww, ss)
414  !call hecmw_precond_clear(hecMAT)
415 
416  e1_time= hecmw_wtime()
417  if (timelog.eq.2) then
418  call hecmw_time_statistics(hecmesh, e1_time - s1_time, &
419  t_max, t_min, t_avg, t_sd)
420  if (hecmesh%my_rank.eq.0) then
421  write(*,*) 'Time solver iterations'
422  write(*,*) ' Max :',t_max
423  write(*,*) ' Min :',t_min
424  write(*,*) ' Avg :',t_avg
425  write(*,*) ' Std Dev :',t_sd
426  endif
427  tsol = t_max
428  else
429  tsol = e1_time - s1_time
430  endif
431 
432  end subroutine hecmw_solve_gmres
433 
434 end module hecmw_solver_gmres
subroutine, public hecmw_estimate_condition_gmres(I, H)
integer(kind=kint) function, public hecmw_mat_get_nrest(hecMAT)
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_estcond(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, public hecmw_solve_gmres(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_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_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=4), parameter kreal
real(kind=kreal) function hecmw_wtime()
subroutine hecmw_update_r(hecMESH, val, n, m)
subroutine hecmw_barrier(hecMESH)
integer(kind=kint), parameter hecmw_solver_error_noconv_maxit