34 integer(kind=kint ),
intent(inout):: iter, error
35 real (kind=
kreal),
intent(inout):: resid, tset, tsol, tcomm
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(:)
42 real(kind=
kreal),
dimension(:,:),
allocatable :: ww
44 integer(kind=kint ) :: maxit, nrest
46 real (kind=
kreal) :: tol
48 real (kind=
kreal),
dimension(:),
allocatable :: ss
49 real (kind=
kreal),
dimension(:,:),
allocatable :: h
51 integer(kind=kint ) :: cs, sn
53 real (kind=
kreal) zero, one
54 parameter( zero = 0.0d+0, one = 1.0d+0 )
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
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
81 my_rank = hecmesh%my_rank
92 if (nrest >= ndof*np-1) nrest = ndof*np-2
98 allocate (ww(ndof*np,nrk))
133 if (bnrm2.eq.0.d0)
then
141 if (timelog.eq.2)
then
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
153 tset = e_time - s_time
173 if (dnrm2 == 0.d0)
exit
200 call hecmw_matvec(hecmesh, hecmat, ww(:,zq), ww(:,w), tcomm)
218 if (val == 0.d0)
exit
222 call hecmw_axpby_r(nndof, coef, 0.d0, ww(:,w), ww(:,v+i+1-1))
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)
246 if (dabs(aa).gt.dabs(bb)) r0= aa
247 scale= dabs(aa) + dabs(bb)
249 if (scale.ne.0.d0)
then
250 rr= scale * dsqrt((aa/scale)**2+(bb/scale)**2)
251 rr= dsign(1.d0,r0)*rr
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)
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)
272 resid = dabs( ww(i+1,s))/dsqrt(bnrm2)
274 if (my_rank.eq.0 .and. iterlog.eq.1) &
275 &
write (*,
'(2i8, 1pe16.6)') iter,i+1, resid
277 if (estcond /= 0 .and. hecmesh%my_rank == 0)
then
281 if ( resid.le.tol )
then
287 ww(irow,y)= ss(irow) / h(irow,irow)
290 do jj= irow, kk+1, -1
291 ss(kk)= ss(kk) - h(kk,jj)*ww(jj,y)
293 ww(kk,y)= ss(kk) / h(kk,kk)
301 call hecmw_axpy_r(nndof, ww(jj,y), ww(:,v+jj-1), ww(:,av))
311 if ( iter.gt.maxit )
then
329 ww(irow,y)= ss(irow) / h(irow,irow)
332 do jj= irow, kk+1, -1
333 ss(kk)= ss(kk) - h(kk,jj)*ww(jj,y)
335 ww(kk,y)= ss(kk) / h(kk,kk)
343 call hecmw_axpy_r(nndof, ww(jj,y), ww(:,v+jj-1), ww(:,av))
356 ww(i+1,s)= dsqrt(dnrm2/bnrm2)
359 if ( resid.le.tol )
exit outer
360 if ( iter .gt.maxit )
then
379 ww(irow,y)= ss(irow) / h(irow,irow)
382 do jj= irow, kk+1, -1
383 ss(kk)= ss(kk) - h(kk,jj)*ww(jj,y)
385 ww(kk,y)= ss(kk) / h(kk,kk)
393 call hecmw_axpy_r(nndof, ww(jj,y), ww(:,v+jj-1), ww(:,av))
403 if (estcond /= 0 .and. hecmesh%my_rank == 0)
then
411 tcomm = tcomm + e_time - s_time
413 deallocate (h, ww, ss)
417 if (timelog.eq.2)
then
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
429 tsol = e1_time - s1_time
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)
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