FrontISTR  5.9.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-- matrix integration for OpenACC
118  !C
119  !C @note:
120  !C Combine hecMAT%AL, D, and AU into a single matrix for GPU execution.
121  !C This is a no-op for CPU builds.
122  call hecmw_mat_integrate(hecmat)
123 
124  !C===
125  !C +----------------------+
126  !C | SETUP PRECONDITIONER |
127  !C +----------------------+
128  !C===
129  call hecmw_precond_setup(hecmat, hecmesh, 0)
130 
131  !C
132  !C
133  !C +--------------------+
134  !C | {r}= {b} - [A]{x0} |
135  !C +--------------------+
136  !C===
137  call hecmw_matresid(hecmesh, hecmat, x, b, ww(:,r), tcomm)
138  !C===
139 
140  call hecmw_innerproduct_r(hecmesh, ndof, b, b, bnrm2, tcomm)
141  if (bnrm2.eq.0.d0) then
142  iter = 0
143  maxit = 0
144  resid = 0.d0
145  x = 0.d0
146  endif
147 
148  e_time= hecmw_wtime()
149  if (timelog.eq.2) then
150  call hecmw_time_statistics(hecmesh, e_time - s_time, &
151  t_max, t_min, t_avg, t_sd)
152  if (hecmesh%my_rank.eq.0) then
153  write(*,*) 'Time solver setup'
154  write(*,*) ' Max :',t_max
155  write(*,*) ' Min :',t_min
156  write(*,*) ' Avg :',t_avg
157  write(*,*) ' Std Dev :',t_sd
158  endif
159  tset = t_max
160  else
161  tset = e_time - s_time
162  endif
163  !C===
164 
165  call hecmw_barrier(hecmesh)
166  s1_time= hecmw_wtime()
167  iter= 0
168 
169  outer: do
170 
171  !C
172  !C************************************************ GMRES Iteration
173  !C
174  i= 0
175  !C
176  !C +---------------+
177  !C | {v1}= {r}/|r| |
178  !C +---------------+
179  !C===
180  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,r), ww(:,r), dnrm2, tcomm)
181  if (dnrm2 == 0.d0) exit ! converged
182 
183  rnorm= dsqrt(dnrm2)
184  coef= one/rnorm
185  call hecmw_axpby_r(nndof, coef, 0.d0, ww(:,r), ww(:,v))
186  !C===
187 
188  !C
189  !C +--------------+
190  !C | {s}= |r|{e1} |
191  !C +--------------+
192  !C===
193  call hecmw_scale_r(nndof, zero, ww(:,s))
194  ww(1 ,s) = rnorm
195  !C===
196 
197  !C************************************************ GMRES(m) restart
198  do i = 1, nrest
199  iter= iter + 1
200 
201  !C
202  !C +-------------------+
203  !C | {w}= [A][Minv]{v} |
204  !C +-------------------+
205  !C===
206  call hecmw_precond_apply(hecmesh, hecmat, ww(:,v+i-1), ww(:,zq), ww(:,zp), tcomm)
207 
208  call hecmw_matvec(hecmesh, hecmat, ww(:,zq), ww(:,w), tcomm)
209  !C===
210 
211  !C
212  !C +------------------------------+
213  !C | ORTH. BASIS by GRAMM-SCHMIDT |
214  !C +------------------------------+
215  !C Construct the I-th column of the upper Hessenberg matrix H
216  !C using the Gram-Schmidt process on V and W.
217  !C===
218  do k= 1, i
219  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,w), ww(:,v+k-1), val, tcomm)
220 
221  call hecmw_axpy_r(nndof, -val, ww(:,v+k-1), ww(:,w))
222  h(k,i)= val
223  enddo
224 
225  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,w), ww(:,w), val, tcomm)
226  if (val == 0.d0) exit ! converged
227 
228  h(i+1,i)= dsqrt(val)
229  coef= one / h(i+1,i)
230  call hecmw_axpby_r(nndof, coef, 0.d0, ww(:,w), ww(:,v+i+1-1))
231  !C===
232 
233  !C
234  !C +-----------------+
235  !C | GIVENS ROTARION |
236  !C +-----------------+
237  !C===
238 
239  !C
240  !C-- Plane Rotation
241  do k = 1, i-1
242  vcs= h(k,cs)
243  vsn= h(k,sn)
244  dtemp = vcs*h(k ,i) + vsn*h(k+1,i)
245  h(k+1,i)= vcs*h(k+1,i) - vsn*h(k ,i)
246  h(k ,i)= dtemp
247  enddo
248 
249  !C
250  !C-- Construct Givens Plane Rotation
251  aa = h(i ,i)
252  bb = h(i+1,i)
253  r0= bb
254  if (dabs(aa).gt.dabs(bb)) r0= aa
255  scale= dabs(aa) + dabs(bb)
256 
257  if (scale.ne.0.d0) then
258  rr= scale * dsqrt((aa/scale)**2+(bb/scale)**2)
259  rr= dsign(1.d0,r0)*rr
260  h(i,cs)= aa/rr
261  h(i,sn)= bb/rr
262  else
263  h(i,cs)= 1.d0
264  h(i,sn)= 0.d0
265  rr = 0.d0
266  endif
267 
268  !C
269  !C-- Plane Rotation
270  vcs= h(i,cs)
271  vsn= h(i,sn)
272  dtemp = vcs*h(i ,i) + vsn*h(i+1,i)
273  h(i+1,i)= vcs*h(i+1,i) - vsn*h(i ,i)
274  h(i ,i)= dtemp
275 
276  dtemp = vcs*ww(i ,s) + vsn*ww(i+1,s)
277  ww(i+1,s)= vcs*ww(i+1,s) - vsn*ww(i ,s)
278  ww(i ,s)= dtemp
279 
280  resid = dabs( ww(i+1,s))/dsqrt(bnrm2)
281 
282  if (my_rank.eq.0 .and. iterlog.eq.1) &
283  & write (*, '(2i8, 1pe16.6)') iter,i+1, resid
284 
285  if (estcond /= 0 .and. hecmesh%my_rank == 0) then
286  if (mod(iter,estcond) == 0) call hecmw_estimate_condition_gmres(i, h)
287  endif
288 
289  if ( resid.le.tol ) then
290  !C-- [H]{y}= {s_tld}
291  do ik= 1, i
292  ss(ik)= ww(ik,s)
293  enddo
294  irow= i
295  ww(irow,y)= ss(irow) / h(irow,irow)
296 
297  do kk= irow-1, 1, -1
298  do jj= irow, kk+1, -1
299  ss(kk)= ss(kk) - h(kk,jj)*ww(jj,y)
300  enddo
301  ww(kk,y)= ss(kk) / h(kk,kk)
302  enddo
303 
304  !C-- {x}= {x} + {y}{V}
305  call hecmw_scale_r(nndof, 0.d0, ww(:,av))
306 
307  jj= irow
308  do jj= 1, irow
309  call hecmw_axpy_r(nndof, ww(jj,y), ww(:,v+jj-1), ww(:,av))
310  enddo
311 
312  call hecmw_precond_apply(hecmesh, hecmat, ww(:,av), ww(:,zq), ww(:,zp), tcomm)
313 
314  call hecmw_axpy_r(nndof, 1.d0, ww(:,zq), x)
315 
316  exit outer
317  endif
318 
319  if ( iter.gt.maxit ) then
321  exit outer
322  end if
323  end do
324  !C===
325 
326  !C
327  !C +------------------+
328  !C | CURRENT SOLUTION |
329  !C +------------------+
330  !C===
331 
332  !C-- [H]{y}= {s_tld}
333  do ik= 1, nrest
334  ss(ik)= ww(ik,s)
335  enddo
336  irow= nrest
337  ww(irow,y)= ss(irow) / h(irow,irow)
338 
339  do kk= irow-1, 1, -1
340  do jj= irow, kk+1, -1
341  ss(kk)= ss(kk) - h(kk,jj)*ww(jj,y)
342  enddo
343  ww(kk,y)= ss(kk) / h(kk,kk)
344  enddo
345 
346  !C-- {x}= {x} + {y}{V}
347  call hecmw_scale_r(nndof, 0.d0, ww(:,av))
348 
349  jj= irow
350  do jj= 1, irow
351  call hecmw_axpy_r(nndof, ww(jj,y), ww(:,v+jj-1), ww(:,av))
352  enddo
353 
354  call hecmw_precond_apply(hecmesh, hecmat, ww(:,av), ww(:,zq), ww(:,zp), tcomm)
355 
356  call hecmw_axpy_r(nndof, 1.d0, ww(:,zq), x)
357 
358  !C
359  !C-- Compute residual vector R, find norm, then check for tolerance.
360  call hecmw_matresid(hecmesh, hecmat, x, b, ww(:,r), tcomm)
361 
362  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,r), ww(:,r), dnrm2, tcomm)
363 
364  ww(i+1,s)= dsqrt(dnrm2/bnrm2)
365  resid = ww( i+1,s )
366 
367  if ( resid.le.tol ) exit outer
368  if ( iter .gt.maxit ) then
370  exit outer
371  end if
372  !C
373  !C-- RESTART
374  end do outer
375 
376  !C
377  !C-- iteration FAILED
378 
379  if (error == hecmw_solver_error_noconv_maxit) then
380  info = iter
381 
382  !C-- [H]{y}= {s_tld}
383  do ik= 1, i
384  ss(ik)= ww(ik,s)
385  enddo
386  irow= i
387  ww(irow,y)= ss(irow) / h(irow,irow)
388 
389  do kk= irow-1, 1, -1
390  do jj= irow, kk+1, -1
391  ss(kk)= ss(kk) - h(kk,jj)*ww(jj,y)
392  enddo
393  ww(kk,y)= ss(kk) / h(kk,kk)
394  enddo
395 
396  !C-- {x}= {x} + {y}{V}
397  call hecmw_scale_r(nndof, 0.d0, ww(:,av))
398 
399  jj= irow
400  do jj= 1, irow
401  call hecmw_axpy_r(nndof, ww(jj,y), ww(:,v+jj-1), ww(:,av))
402  enddo
403 
404  call hecmw_precond_apply(hecmesh, hecmat, ww(:,av), ww(:,zq), ww(:,zp), tcomm)
405 
406  call hecmw_axpy_r(nndof, 1.d0, ww(:,zq), x)
407  end if
408 
409  call hecmw_solver_scaling_bk(hecmat)
410 
411  if (estcond /= 0 .and. hecmesh%my_rank == 0) then
413  endif
414  !C
415  !C-- INTERFACE data EXCHANGE
416  s_time = hecmw_wtime()
417  call hecmw_update_r (hecmesh, x, hecmat%NP, hecmat%NDOF)
418  e_time = hecmw_wtime()
419  tcomm = tcomm + e_time - s_time
420 
421  deallocate (h, ww, ss)
422  !call hecmw_precond_clear(hecMAT)
423 
424  e1_time= hecmw_wtime()
425  if (timelog.eq.2) then
426  call hecmw_time_statistics(hecmesh, e1_time - s1_time, &
427  t_max, t_min, t_avg, t_sd)
428  if (hecmesh%my_rank.eq.0) then
429  write(*,*) 'Time solver iterations'
430  write(*,*) ' Max :',t_max
431  write(*,*) ' Min :',t_min
432  write(*,*) ' Avg :',t_avg
433  write(*,*) ' Std Dev :',t_sd
434  endif
435  tsol = t_max
436  else
437  tsol = e1_time - s1_time
438  endif
439 
440  end subroutine hecmw_solve_gmres
441 
442 end module hecmw_solver_gmres
subroutine, public hecmw_estimate_condition_gmres(I, H)
subroutine, public hecmw_mat_integrate(hecMAT)
Integrate matrix components into a single array for efficient access.
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