FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
hecmw_solver_GMRESR.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_GMRESR
8 !C***
9 !
11 
12  public :: hecmw_solve_gmresr
13 
14 contains
15  !C
16  !C*** hecmw_solve_GMRESR
17  !C
18  subroutine hecmw_solve_gmresr( 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  integer(kind=kint),dimension(:), allocatable :: idxbfgs
43  real(kind=kreal), dimension(:) , allocatable :: vecr,workpc,tmpvecbfgs,rho,alpha
44  real(kind=kreal), dimension(:,:), allocatable :: u,c,uin,cin,sbfgs,ybfgs
45 
46  integer(kind=kint ) :: maxit, nrest,nbfgs
47 
48  real (kind=kreal) :: tol
49 
50  real (kind=kreal) zero, one
51  parameter( zero = 0.0d+0, one = 1.0d+0 )
52 
53  integer(kind=kint ) :: nrk,i,k,kk,jj,info,ik,iorth,idx,tmpidx,ibfgs
54  integer(kind=kint ) :: irow
55  real (kind=kreal) :: s_time,e_time,s1_time,e1_time
56  real (kind=kreal) :: ldh,ldw,bnrm2,dnrm2,rnorm
57  real (kind=kreal) :: commtime,comptime, coef,coef2,val,vcs,vsn,dtemp,aa,bb,r0,scale,rr
58  integer(kind=kint ) :: estcond
59  real (kind=kreal) :: t_max,t_min,t_avg,t_sd
60 
61 
62  call hecmw_barrier(hecmesh)
63  s_time= hecmw_wtime()
64  !C
65  !C-- INIT.
66  n = hecmat%N
67  np = hecmat%NP
68  ndof = hecmat%NDOF
69  nndof = n * ndof
70  my_rank = hecmesh%my_rank
71  x => hecmat%X
72  b => hecmat%B
73 
74  iterlog = hecmw_mat_get_iterlog( hecmat )
75  timelog = hecmw_mat_get_timelog( hecmat )
76  maxit = hecmw_mat_get_iter( hecmat )
77  tol = hecmw_mat_get_resid( hecmat )
78  nrest = hecmw_mat_get_nrest( hecmat )
79  nbfgs = hecmw_mat_get_nbfgs( hecmat )
80  estcond = hecmw_mat_get_estcond( hecmat )
81 
82  error= 0
83 
84  allocate (vecr(ndof*np))
85  allocate (workpc(ndof*np))
86  allocate (u(ndof*np,nrest))
87  allocate (c(ndof*np,nrest))
88  allocate (uin(ndof*np,nrest))
89  allocate (cin(ndof*np,nrest))
90 
91  if(nbfgs>0)then
92  allocate (tmpvecbfgs(ndof*np))
93  allocate (sbfgs(ndof*np,nbfgs))
94  allocate (ybfgs(ndof*np,nbfgs))
95  allocate (idxbfgs(nbfgs))
96  allocate (rho(nbfgs))
97  allocate (alpha(nbfgs))
98  endif
99 
100  do idx = 1, nbfgs
101  idxbfgs(idx) = idx
102  enddo
103 
104  commtime= 0.d0
105  comptime= 0.d0
106 
107  !C
108  !C-- SCALING
109  call hecmw_solver_scaling_fw(hecmesh, hecmat, tcomm)
110 
111  !C
112  !C-- matrix integration for OpenACC
113  !C
114  !C @note:
115  !C Combine hecMAT%AL, D, and AU into a single matrix for GPU execution.
116  !C This is a no-op for CPU builds.
117  call hecmw_mat_integrate(hecmat)
118 
119  !C===
120  !C +----------------------+
121  !C | SETUP PRECONDITIONER |
122  !C +----------------------+
123  !C===
124  call hecmw_precond_setup(hecmat, hecmesh, 0)
125 
126 
127  call hecmw_innerproduct_r(hecmesh, ndof, b, b, bnrm2, tcomm)
128  if (bnrm2.eq.0.d0) then
129  iter = 0
130  maxit = 0
131  resid = 0.d0
132  x = 0.d0
133  endif
134 
135  e_time= hecmw_wtime()
136  if (timelog.eq.2) then
137  call hecmw_time_statistics(hecmesh, e_time - s_time, &
138  t_max, t_min, t_avg, t_sd)
139  if (hecmesh%my_rank.eq.0) then
140  write(*,*) 'Time solver setup'
141  write(*,*) ' Max :',t_max
142  write(*,*) ' Min :',t_min
143  write(*,*) ' Avg :',t_avg
144  write(*,*) ' Std Dev :',t_sd
145  endif
146  tset = t_max
147  else
148  tset = e_time - s_time
149  endif
150  !C===
151 
152 
153  call hecmw_barrier(hecmesh)
154  s1_time= hecmw_wtime()
155  iter= 0
156  ibfgs=0
157  outer: do
158 
159  call hecmw_matresid(hecmesh, hecmat, x, b, vecr, tcomm)
160  do i = 1, nrest
161  iter= iter + 1
162 
163  !C Solve M*r = uin(:,1)
164  if (ibfgs == 0)then
165  call hecmw_precond_apply(hecmesh, hecmat, vecr, uin(:,1), workpc, tcomm)
166  else
167  call hecmw_copy_r(nndof, vecr, tmpvecbfgs)
168  do k = 1,ibfgs
169  idx = idxbfgs(k)
170  call hecmw_innerproduct_r(hecmesh, ndof, sbfgs(:,idx), ybfgs(:,idx), coef, tcomm)
171  rho(k) = 1.0d0 / coef
172  call hecmw_innerproduct_r(hecmesh, ndof, sbfgs(:,idx), tmpvecbfgs, coef2, tcomm)
173  alpha(k) = rho(k)*coef2
174  call hecmw_axpy_r(nndof, -alpha(k), ybfgs(:,idx), tmpvecbfgs)
175  enddo
176  call hecmw_precond_apply(hecmesh, hecmat, tmpvecbfgs, uin(:,1), workpc, tcomm)
177  do k = ibfgs,1,-1
178  idx = idxbfgs(k)
179  call hecmw_innerproduct_r(hecmesh, ndof, ybfgs(:,idx), uin(:,1), coef, tcomm)
180  coef2 = rho(k) * coef
181  call hecmw_axpy_r(nndof, alpha(k)-coef2, sbfgs(:,idx), uin(:,1))
182  enddo
183  endif
184  !C cin(:,1) = A*uin(:,1)
185  call hecmw_matvec(hecmesh, hecmat, uin(:,1), cin(:,1), tcomm)
186 
187  do iorth = 1, i-1
188  !C c_{i}^T cin_{i}
189  call hecmw_innerproduct_r(hecmesh, ndof, c(:,iorth), cin(:,iorth), coef, tcomm)
190 
191  call hecmw_axpyz_r(nndof, -coef, c(:,iorth), cin(:,iorth), cin(:,iorth+1))
192  call hecmw_axpyz_r(nndof, -coef, u(:,iorth), uin(:,iorth), uin(:,iorth+1))
193  enddo
194  call hecmw_innerproduct_r(hecmesh, ndof, cin(:,i), cin(:,i), coef, tcomm)
195  coef = 1.0d0 / dsqrt(coef)
196  call hecmw_axpby_r(nndof, coef, 0.0d0, cin(:,i), c(:,i))
197  call hecmw_axpby_r(nndof, coef, 0.0d0, uin(:,i), u(:,i))
198 
199  call hecmw_innerproduct_r(hecmesh, ndof, c(:,i), vecr, coef, tcomm)
200  call hecmw_axpy_r(nndof, coef, u(:,i), x)
201  call hecmw_axpy_r(nndof, -coef, c(:,i), vecr)
202 
203  if (nbfgs > 0)then
204  ibfgs = ibfgs + 1
205  if (ibfgs == nbfgs+1)then
206  tmpidx = idxbfgs(1)
207  do kk = 1, nbfgs-1
208  idxbfgs(kk) = idxbfgs(kk+1)
209  enddo
210  idxbfgs(nbfgs) = tmpidx
211  ibfgs = ibfgs - 1
212  endif
213  call hecmw_axpby_r(nndof, coef, 0.0d0, c(:,i), ybfgs(:,idxbfgs(ibfgs)))
214  call hecmw_axpby_r(nndof, coef, 0.0d0, u(:,i), sbfgs(:,idxbfgs(ibfgs)))
215  endif
216 
217 
218  call hecmw_innerproduct_r(hecmesh, ndof, vecr, vecr, dnrm2, tcomm)
219  resid= dsqrt(dnrm2/bnrm2)
220 
221  !C##### ITERATION HISTORY
222  if (my_rank.eq.0.and.iterlog.eq.1) write (*,'(i7, 1pe16.6)') iter, resid
223  !C#####
224 
225  if ( resid.le.tol ) exit outer
226  if ( iter.gt.maxit ) then
228  exit outer
229  end if
230  end do
231 
232  end do outer
233 
234  call hecmw_solver_scaling_bk(hecmat)
235 
236  !C
237  !C-- INTERFACE data EXCHANGE
238  s_time = hecmw_wtime()
239  !call hecmw_update_m_R (hecMESH, X, hecMAT%NP, hecMAT%NDOF)
240  call hecmw_update_r (hecmesh, x, hecmat%NP, hecmat%NDOF)
241 
242  e_time = hecmw_wtime()
243  tcomm = tcomm + e_time - s_time
244 
245  !deallocate (H, WW, SS)
246  deallocate (vecr)
247  deallocate (workpc)
248  deallocate (u )
249  deallocate (c )
250  deallocate (uin)
251  deallocate (cin)
252  call hecmw_precond_clear(hecmat)
253 
254  e1_time= hecmw_wtime()
255  if (timelog.eq.2) then
256  call hecmw_time_statistics(hecmesh, e1_time - s1_time, &
257  t_max, t_min, t_avg, t_sd)
258  if (hecmesh%my_rank.eq.0) then
259  write(*,*) 'Time solver iterations'
260  write(*,*) ' Max :',t_max
261  write(*,*) ' Min :',t_min
262  write(*,*) ' Avg :',t_avg
263  write(*,*) ' Std Dev :',t_sd
264  endif
265  tsol = t_max
266  else
267  tsol = e1_time - s1_time
268  endif
269 
270  end subroutine hecmw_solve_gmresr
271 
272 end module hecmw_solver_gmresr
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_nbfgs(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_clear(hecMAT)
subroutine, public hecmw_precond_setup(hecMAT, hecMESH, sym)
subroutine, public hecmw_precond_apply(hecMESH, hecMAT, R, Z, ZP, COMMtime)
subroutine, public hecmw_solve_gmresr(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_axpyz_r(n, alpha, X, Y, Z)
subroutine hecmw_innerproduct_r(hecMESH, ndof, X, Y, sum, COMMtime)
subroutine hecmw_axpby_r(n, alpha, beta, X, Y)
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=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