FrontISTR  5.8.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 +----------------------+
113  !C | SETUP PRECONDITIONER |
114  !C +----------------------+
115  !C===
116  call hecmw_precond_setup(hecmat, hecmesh, 0)
117 
118 
119  call hecmw_innerproduct_r(hecmesh, ndof, b, b, bnrm2, tcomm)
120  if (bnrm2.eq.0.d0) then
121  iter = 0
122  maxit = 0
123  resid = 0.d0
124  x = 0.d0
125  endif
126 
127  e_time= hecmw_wtime()
128  if (timelog.eq.2) then
129  call hecmw_time_statistics(hecmesh, e_time - s_time, &
130  t_max, t_min, t_avg, t_sd)
131  if (hecmesh%my_rank.eq.0) then
132  write(*,*) 'Time solver setup'
133  write(*,*) ' Max :',t_max
134  write(*,*) ' Min :',t_min
135  write(*,*) ' Avg :',t_avg
136  write(*,*) ' Std Dev :',t_sd
137  endif
138  tset = t_max
139  else
140  tset = e_time - s_time
141  endif
142  !C===
143 
144 
145  call hecmw_barrier(hecmesh)
146  s1_time= hecmw_wtime()
147  iter= 0
148  ibfgs=0
149  outer: do
150 
151  call hecmw_matresid(hecmesh, hecmat, x, b, vecr, tcomm)
152  do i = 1, nrest
153  iter= iter + 1
154 
155  !C Solve M*r = uin(:,1)
156  if (ibfgs == 0)then
157  call hecmw_precond_apply(hecmesh, hecmat, vecr, uin(:,1), workpc, tcomm)
158  else
159  call hecmw_copy_r(nndof, vecr, tmpvecbfgs)
160  do k = 1,ibfgs
161  idx = idxbfgs(k)
162  call hecmw_innerproduct_r(hecmesh, ndof, sbfgs(:,idx), ybfgs(:,idx), coef, tcomm)
163  rho(k) = 1.0d0 / coef
164  call hecmw_innerproduct_r(hecmesh, ndof, sbfgs(:,idx), tmpvecbfgs, coef2, tcomm)
165  alpha(k) = rho(k)*coef2
166  call hecmw_axpy_r(nndof, -alpha(k), ybfgs(:,idx), tmpvecbfgs)
167  enddo
168  call hecmw_precond_apply(hecmesh, hecmat, tmpvecbfgs, uin(:,1), workpc, tcomm)
169  do k = ibfgs,1,-1
170  idx = idxbfgs(k)
171  call hecmw_innerproduct_r(hecmesh, ndof, ybfgs(:,idx), uin(:,1), coef, tcomm)
172  coef2 = rho(k) * coef
173  call hecmw_axpy_r(nndof, alpha(k)-coef2, sbfgs(:,idx), uin(:,1))
174  enddo
175  endif
176  !C cin(:,1) = A*uin(:,1)
177  call hecmw_matvec(hecmesh, hecmat, uin(:,1), cin(:,1), tcomm)
178 
179  do iorth = 1, i-1
180  !C c_{i}^T cin_{i}
181  call hecmw_innerproduct_r(hecmesh, ndof, c(:,iorth), cin(:,iorth), coef, tcomm)
182 
183  call hecmw_axpyz_r(nndof, -coef, c(:,iorth), cin(:,iorth), cin(:,iorth+1))
184  call hecmw_axpyz_r(nndof, -coef, u(:,iorth), uin(:,iorth), uin(:,iorth+1))
185  enddo
186  call hecmw_innerproduct_r(hecmesh, ndof, cin(:,i), cin(:,i), coef, tcomm)
187  coef = 1.0d0 / dsqrt(coef)
188  call hecmw_axpby_r(nndof, coef, 0.0d0, cin(:,i), c(:,i))
189  call hecmw_axpby_r(nndof, coef, 0.0d0, uin(:,i), u(:,i))
190 
191  call hecmw_innerproduct_r(hecmesh, ndof, c(:,i), vecr, coef, tcomm)
192  call hecmw_axpy_r(nndof, coef, u(:,i), x)
193  call hecmw_axpy_r(nndof, -coef, c(:,i), vecr)
194 
195  if (nbfgs > 0)then
196  ibfgs = ibfgs + 1
197  if (ibfgs == nbfgs+1)then
198  tmpidx = idxbfgs(1)
199  do kk = 1, nbfgs-1
200  idxbfgs(kk) = idxbfgs(kk+1)
201  enddo
202  idxbfgs(nbfgs) = tmpidx
203  ibfgs = ibfgs - 1
204  endif
205  call hecmw_axpby_r(nndof, coef, 0.0d0, c(:,i), ybfgs(:,idxbfgs(ibfgs)))
206  call hecmw_axpby_r(nndof, coef, 0.0d0, u(:,i), sbfgs(:,idxbfgs(ibfgs)))
207  endif
208 
209 
210  call hecmw_innerproduct_r(hecmesh, ndof, vecr, vecr, dnrm2, tcomm)
211  resid= dsqrt(dnrm2/bnrm2)
212 
213  !C##### ITERATION HISTORY
214  if (my_rank.eq.0.and.iterlog.eq.1) write (*,'(i7, 1pe16.6)') iter, resid
215  !C#####
216 
217  if ( resid.le.tol ) exit outer
218  if ( iter.gt.maxit ) then
220  exit outer
221  end if
222  end do
223 
224  end do outer
225 
226  call hecmw_solver_scaling_bk(hecmat)
227 
228  !C
229  !C-- INTERFACE data EXCHANGE
230  s_time = hecmw_wtime()
231  !call hecmw_update_m_R (hecMESH, X, hecMAT%NP, hecMAT%NDOF)
232  call hecmw_update_r (hecmesh, x, hecmat%NP, hecmat%NDOF)
233 
234  e_time = hecmw_wtime()
235  tcomm = tcomm + e_time - s_time
236 
237  !deallocate (H, WW, SS)
238  deallocate (vecr)
239  deallocate (workpc)
240  deallocate (u )
241  deallocate (c )
242  deallocate (uin)
243  deallocate (cin)
244  call hecmw_precond_clear(hecmat)
245 
246  e1_time= hecmw_wtime()
247  if (timelog.eq.2) then
248  call hecmw_time_statistics(hecmesh, e1_time - s1_time, &
249  t_max, t_min, t_avg, t_sd)
250  if (hecmesh%my_rank.eq.0) then
251  write(*,*) 'Time solver iterations'
252  write(*,*) ' Max :',t_max
253  write(*,*) ' Min :',t_min
254  write(*,*) ' Avg :',t_avg
255  write(*,*) ' Std Dev :',t_sd
256  endif
257  tsol = t_max
258  else
259  tsol = e1_time - s1_time
260  endif
261 
262  end subroutine hecmw_solve_gmresr
263 
264 end module hecmw_solver_gmresr
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