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
40 real (kind=
kreal),
pointer :: b(:)
41 real (kind=
kreal),
pointer :: x(:)
43 integer(kind=kint ) :: iterlog, timelog
45 real(kind=
kreal),
dimension(:,:),
allocatable :: ww
47 real(kind=
kreal),
dimension(2) :: rr
49 integer(kind=kint ) :: maxit
50 real (kind=
kreal) :: tol
51 integer(kind=kint ) :: i,j
52 real (kind=
kreal) :: s_time,s1_time,e_time,e1_time
53 real (kind=
kreal) :: bnrm2
54 real (kind=
kreal) :: rho,rho1,beta,alpha,dnrm2
55 real (kind=
kreal) :: qsi,eta,coef1
56 real (kind=
kreal) :: t_max,t_min,t_avg,t_sd
58 integer(kind=kint),
parameter :: r= 1
59 integer(kind=kint),
parameter ::rt= 2
60 integer(kind=kint),
parameter :: t= 3
61 integer(kind=kint),
parameter ::tt= 4
62 integer(kind=kint),
parameter ::t0= 5
63 integer(kind=kint),
parameter :: p= 6
64 integer(kind=kint),
parameter ::pt= 7
65 integer(kind=kint),
parameter :: u= 8
66 integer(kind=kint),
parameter ::w1= 9
67 integer(kind=kint),
parameter :: y=10
68 integer(kind=kint),
parameter :: z=11
69 integer(kind=kint),
parameter ::wk=12
70 integer(kind=kint),
parameter ::w2=13
71 integer(kind=kint),
parameter ::zq=14
73 integer(kind=kint),
parameter :: n_iter_recompute_r= 20
83 my_rank = hecmesh%my_rank
95 allocate (ww(ndof*np,14))
119 if (bnrm2.eq.0.d0)
then
129 if (timelog.eq.2)
then
131 t_max, t_min, t_avg, t_sd)
132 if (hecmesh%my_rank.eq.0)
then
133 write(*,*)
'Time solver setup'
134 write(*,*)
' Max :',t_max
135 write(*,*)
' Min :',t_min
136 write(*,*)
' Avg :',t_avg
137 write(*,*)
' Std Dev :',t_sd
141 tset = e_time - s_time
182 call hecmw_matvec(hecmesh, hecmat, ww(:,p), ww(:,pt), tcomm)
197 call hecmw_axpyz_r(nndof, -1.0d0, ww(:,w1), ww(:,pt), ww(:,y))
200 call hecmw_axpyz_r(nndof, -alpha, ww(:,pt), ww(:,wk), ww(:,t))
222 call hecmw_matvec(hecmesh, hecmat, ww(:,tt), ww(:,wk), tcomm)
234 call pol_coef_vanilla2(iter, ww, t, tt, y, qsi, eta)
254 call hecmw_axpyz_r(nndof, -1.0d0, ww(:,r), ww(:,t0), ww(:,u))
271 if ( mod(iter,n_iter_recompute_r)==0 )
then
283 tcomm = tcomm + e_time - s_time
287 beta = alpha*coef1 / (qsi*rho)
288 call hecmw_axpyz_r(nndof, beta, ww(:,pt), ww(:,tt), ww(:,w1))
290 resid= dsqrt(dnrm2/bnrm2)
294 if (my_rank.eq.0 .and. iterlog.eq.1) &
295 &
write (*, 1000) iter, resid
296 1000
format (i5, 1pe16.6)
299 if (resid.le.tol )
then
300 if ( mod(iter,n_iter_recompute_r)==0 )
exit
304 resid= dsqrt(dnrm2/bnrm2)
305 if (resid.le.tol)
exit
318 tcomm = tcomm + e_time - s_time
324 if (timelog.eq.2)
then
326 t_max, t_min, t_avg, t_sd)
327 if (hecmesh%my_rank.eq.0)
then
328 write(*,*)
'Time solver iterations'
329 write(*,*)
' Max :',t_max
330 write(*,*)
' Min :',t_min
331 write(*,*)
' Avg :',t_avg
332 write(*,*)
' Std Dev :',t_sd
336 tsol = e1_time - s1_time
346 integer(kind=kint),
intent(in) :: iter
347 real(kind=
kreal),
intent(in) :: ww(:,:)
348 integer(kind=kint),
intent(in) :: T, TT, Y
349 real(kind=
kreal),
intent(out) :: qsi, eta
351 real(kind=
kreal),
dimension(5) :: cg
352 real(kind=
kreal),
dimension(2) :: eq
353 real(kind=
kreal) :: delta
363 tcomm = tcomm + e_time - s_time
369 delta= (cg(5)*cg(1)-cg(4)*cg(4))
370 eq(1)= (cg(1)*cg(2)-cg(3)*cg(4)) / delta
371 eq(2)= (cg(5)*cg(3)-cg(4)*cg(2)) / delta
382 subroutine pol_coef_vanilla(iter, WW, T, TT, Y, QSI, ETA)
384 integer(kind=kint),
intent(in) :: iter
385 real(kind=
kreal),
intent(inout) :: ww(:,:)
386 integer(kind=kint),
intent(in) :: T, TT, Y
387 real(kind=
kreal),
intent(out) :: qsi, eta
389 real(kind=
kreal),
dimension(3) :: cg
390 real(kind=
kreal) :: gamma1, gamma2
391 real(kind=
kreal) :: c, c_abs
393 real(kind=
kreal),
parameter :: omega = 0.707106781d0
405 tcomm = tcomm + e_time - s_time
421 tcomm = tcomm + e_time - s_time
423 c = cg(3) / dsqrt(cg(1)*cg(2))
425 if (c_abs > omega)
then
426 qsi = c * dsqrt(cg(1)/cg(2))
430 qsi = omega * dsqrt(cg(1)/cg(2))
432 qsi = -omega * dsqrt(cg(1)/cg(2))
435 eta = gamma1 - qsi*gamma2
436 end subroutine pol_coef_vanilla
441 subroutine pol_coef_vanilla2(iter, WW, T, TT, Y, QSI, ETA)
443 integer(kind=kint),
intent(in) :: iter
444 real(kind=
kreal),
intent(inout) :: ww(:,:)
445 integer(kind=kint),
intent(in) :: T, TT, Y
446 real(kind=
kreal),
intent(out) :: qsi, eta
448 real(kind=
kreal),
dimension(6) :: cg
449 real(kind=
kreal) :: gamma1, gamma2
450 real(kind=
kreal) :: c, c_abs
452 real(kind=
kreal),
parameter :: omega = 0.707106781d0
465 tcomm = tcomm + e_time - s_time
472 tcomm = tcomm + e_time - s_time
477 c = cg(3) / dsqrt(cg(1)*cg(2))
479 if (c_abs > omega)
then
480 qsi = c * dsqrt(cg(1)/cg(2))
483 qsi = omega * dsqrt(cg(1)/cg(2))
485 qsi = -omega * dsqrt(cg(1)/cg(2))
488 eta = gamma1 - qsi*gamma2
489 end subroutine pol_coef_vanilla2
subroutine pol_coef(iter, WW, T, TT, Y, QSI, ETA)
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_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_gpbicg(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_xpay_r(n, alpha, X, Y)
subroutine hecmw_axpyz_r(n, alpha, X, Y, Z)
subroutine hecmw_innerproduct_r_nocomm(hecMESH, ndof, X, Y, sum)
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)
integer(kind=kint), parameter hecmw_sum
integer(kind=4), parameter kreal
real(kind=kreal) function hecmw_wtime()
subroutine hecmw_update_r(hecMESH, val, n, m)
subroutine hecmw_allreduce_r(hecMESH, val, n, ntag)
subroutine hecmw_barrier(hecMESH)
integer(kind=kint), parameter hecmw_solver_error_noconv_maxit