module huti_cg
use huti_aux
implicit none
#include "huti_fdefs.h"
#define X xvec
#define B rhsvec
#define Z work(:,1)
#define Z_ind 1
#define P work(:,2)
#define P_ind 2
#define Q work(:,3)
#define Q_ind 3
#define R work(:,4)
#define R_ind 4
contains
recursive subroutine huti_scgsolv ( ndim, wrkdim, xvec, rhsvec, &
ipar, dpar, work, matvecsubr, pcondlsubr, &
pcondrsubr, dotprodfun, normfun, stopcfun )
use huti_interfaces
implicit none
procedure( mv_iface_s ), pointer :: matvecsubr
procedure( pc_iface_s ), pointer :: pcondlsubr
procedure( pc_iface_s ), pointer :: pcondrsubr
procedure( dotp_iface_s ), pointer :: dotprodfun
procedure( norm_iface_s ), pointer :: normfun
procedure( stopc_iface_s ), pointer :: stopcfun
integer :: ndim, wrkdim
real, dimension(ndim) :: xvec, rhsvec
integer, dimension(HUTI_IPAR_DFLTSIZE) :: ipar
double precision, dimension(HUTI_DPAR_DFLTSIZE) :: dpar
real, dimension(ndim,wrkdim) :: work
real :: alpha, beta, rho, oldrho
integer iter_count
real :: residual, rhsnorm, precrhsnorm
iter_count = 1
if ( HUTI_STOPC .eq. HUTI_TRESID_SCALED_BYB .or. &
HUTI_STOPC .eq. HUTI_PRESID_SCALED_BYB ) then
rhsnorm = normfun( HUTI_NDIM, B, 1 )
end if
if ( HUTI_STOPC .eq. HUTI_PRESID_SCALED_BYPRECB ) then
call pcondlsubr( P, B, ipar )
precrhsnorm = normfun( HUTI_NDIM, P, 1 )
end if
HUTI_EXTOP_MATTYPE = HUTI_MAT_NOTTRPSED
if ( HUTI_INITIALX .eq. HUTI_RANDOMX ) then
call huti_srandvec ( X, ipar )
else if ( HUTI_INITIALX .ne. HUTI_USERSUPPLIEDX ) then
X = 1
end if
call matvecsubr( X, R, ipar )
R = B - R
300 continue
call pcondlsubr( Q, R, ipar )
call pcondrsubr( Z, Q, ipar )
rho = dotprodfun( HUTI_NDIM, R, 1, Z, 1 )
if ( rho .eq. 0 ) then
HUTI_INFO = HUTI_CG_RHO
go to 1000
end if
if ( iter_count .eq. 1 ) then
P = Z
else
beta = rho / oldrho
P = Z + beta * P
end if
call matvecsubr( P, Q, ipar )
alpha = rho / dotprodfun( HUTI_NDIM, P, 1, Q, 1 )
X = X + alpha * P
R = R - alpha * Q
select case (HUTI_STOPC)
case (HUTI_TRUERESIDUAL)
call matvecsubr( X, Z, ipar )
Z = Z - B
residual = normfun( HUTI_NDIM, Z, 1 )
case (HUTI_TRESID_SCALED_BYB)
call matvecsubr( X, Z, ipar )
Z = Z - B
residual = normfun( HUTI_NDIM, Z, 1 ) / rhsnorm
case (HUTI_PSEUDORESIDUAL)
residual = normfun( HUTI_NDIM, R, 1 )
case (HUTI_PRESID_SCALED_BYB)
residual = normfun( HUTI_NDIM, R, 1 ) / rhsnorm
case (HUTI_PRESID_SCALED_BYPRECB)
residual = normfun( HUTI_NDIM, R, 1 ) / precrhsnorm
case (HUTI_XDIFF_NORM)
Z = alpha * P
residual = normfun( HUTI_NDIM, Z, 1 )
case (HUTI_USUPPLIED_STOPC)
residual = stopcfun( X, B, R, ipar, dpar )
case default
call matvecsubr( X, Z, ipar )
Z = Z - B
residual = normfun( HUTI_NDIM, Z, 1 )
end select
if ( HUTI_DBUGLVL .ne. HUTI_NO_DEBUG ) then
if ( mod(iter_count, HUTI_DBUGLVL) .eq. 0 ) then
WRITE (*, '(I8, E11.4)') iter_count, residual
CALL FLUSH(6)
end if
end if
IF ( residual .LT. HUTI_TOLERANCE ) THEN
HUTI_INFO = HUTI_CONVERGENCE
go to 1000
end if
IF( residual /= residual .OR. residual > HUTI_MAXTOLERANCE ) THEN
HUTI_INFO = HUTI_DIVERGENCE
GOTO 1000
END IF
oldrho = rho
iter_count = iter_count + 1
if ( iter_count .gt. HUTI_MAXIT ) then
HUTI_INFO = HUTI_MAXITER
go to 1000
end if
go to 300
1000 continue
if ( HUTI_DBUGLVL .ne. HUTI_NO_DEBUG ) then
WRITE (*, '(I8, E11.4)') iter_count, residual
CALL FLUSH(6)
end if
HUTI_ITERS = iter_count
return
end subroutine huti_scgsolv
recursive subroutine huti_dcgsolv ( ndim, wrkdim, xvec, rhsvec, &
ipar, dpar, work, matvecsubr, pcondlsubr, &
pcondrsubr, dotprodfun, normfun, stopcfun )
use huti_interfaces
implicit none
procedure( mv_iface_d ), pointer :: matvecsubr
procedure( pc_iface_d ), pointer :: pcondlsubr
procedure( pc_iface_d ), pointer :: pcondrsubr
procedure( dotp_iface_d ), pointer :: dotprodfun
procedure( norm_iface_d ), pointer :: normfun
procedure( stopc_iface_d ), pointer :: stopcfun
integer :: ndim, wrkdim
double precision, dimension(ndim) :: xvec, rhsvec
integer, dimension(HUTI_IPAR_DFLTSIZE) :: ipar
double precision, dimension(HUTI_DPAR_DFLTSIZE) :: dpar
double precision, dimension(ndim,wrkdim) :: work
double precision :: alpha, beta, rho, oldrho
integer iter_count, i
double precision :: residual, rhsnorm, precrhsnorm
iter_count = 1
if ( HUTI_STOPC .eq. HUTI_TRESID_SCALED_BYB .or. &
HUTI_STOPC .eq. HUTI_PRESID_SCALED_BYB ) then
rhsnorm = normfun( HUTI_NDIM, B, 1 )
end if
if ( HUTI_STOPC .eq. HUTI_PRESID_SCALED_BYPRECB ) then
call pcondlsubr( P, B, ipar )
precrhsnorm = normfun( HUTI_NDIM, P, 1 )
end if
HUTI_EXTOP_MATTYPE = HUTI_MAT_NOTTRPSED
if ( HUTI_INITIALX .eq. HUTI_RANDOMX ) then
call huti_drandvec ( X, ipar )
else if ( HUTI_INITIALX .ne. HUTI_USERSUPPLIEDX ) then
X = 1
end if
call matvecsubr( X, R, ipar )
#ifdef _OPENMP
do i=1,ndim
work(i,R_ind)=rhsvec(i)-work(i,R_ind)
end do
#else
R = B - R
#endif
300 continue
call pcondlsubr( Q, R, ipar )
call pcondrsubr( Z, Q, ipar )
rho = dotprodfun( HUTI_NDIM, R, 1, Z, 1 )
if ( rho .eq. 0 ) then
HUTI_INFO = HUTI_CG_RHO
go to 1000
end if
if ( iter_count .eq. 1 ) then
#ifdef _OPENMP
DO i=1,ndim
work(i,P_ind)=work(i,Z_ind)
END DO
#else
P = Z
#endif
else
beta = rho / oldrho
#ifdef _OPENMP
DO i=1,ndim
work(i,P_ind)=work(i,Z_ind) + beta * work(i,P_ind)
END DO
#else
P = Z + beta * P
#endif
end if
call matvecsubr( P, Q, ipar )
alpha = rho / dotprodfun( HUTI_NDIM, P, 1, Q, 1 )
#ifdef _OPENMP
DO i=1,ndim
xvec(i) = xvec(i) + alpha * work(i,P_ind)
END DO
DO i=1,ndim
work(i,R_ind) = work(i,R_ind) - alpha * work(i,Q_ind)
END DO
#else
X = X + alpha * P
R = R - alpha * Q
#endif
select case (HUTI_STOPC)
case (HUTI_TRUERESIDUAL)
call matvecsubr( X, Z, ipar )
#ifdef _OPENMP
DO i=1,ndim
work(i,Z_ind) = work(i,Z_ind) - rhsvec(i)
END DO
#else
Z = Z - B
#endif
residual = normfun( HUTI_NDIM, Z, 1 )
case (HUTI_TRESID_SCALED_BYB)
call matvecsubr( X, Z, ipar )
#ifdef _OPENMP
DO i=1,ndim
work(i,Z_ind) = work(i,Z_ind) - rhsvec(i)
END DO
#else
Z = Z - B
#endif
residual = normfun( HUTI_NDIM, Z, 1 ) / rhsnorm
case (HUTI_PSEUDORESIDUAL)
residual = normfun( HUTI_NDIM, R, 1 )
case (HUTI_PRESID_SCALED_BYB)
residual = normfun( HUTI_NDIM, R, 1 ) / rhsnorm
case (HUTI_PRESID_SCALED_BYPRECB)
residual = normfun( HUTI_NDIM, R, 1 ) / precrhsnorm
case (HUTI_XDIFF_NORM)
#ifdef _OPENMP
DO i=1,ndim
work(i,Z_ind) = alpha * work(i,P_ind)
END DO
#else
Z = alpha * P
#endif
residual = normfun( HUTI_NDIM, Z, 1 )
case (HUTI_USUPPLIED_STOPC)
residual = stopcfun( X, B, R, ipar, dpar )
case default
call matvecsubr( X, Z, ipar )
#ifdef _OPENMP
DO i=1,ndim
work(i,Z_ind) = work(i,Z_ind) - rhsvec(i)
END DO
#else
Z = Z - B
#endif
residual = normfun( HUTI_NDIM, Z, 1 )
end select
if ( HUTI_DBUGLVL .ne. HUTI_NO_DEBUG ) then
if ( mod(iter_count, HUTI_DBUGLVL) .eq. 0 ) then
WRITE (*, '(I8, E11.4)') iter_count, residual
CALL FLUSH(6)
end if
end if
if ( residual .lt. HUTI_TOLERANCE ) then
HUTI_INFO = HUTI_CONVERGENCE
go to 1000
end if
IF( residual /= residual .OR. residual > HUTI_MAXTOLERANCE ) THEN
HUTI_INFO = HUTI_DIVERGENCE
GOTO 1000
END IF
oldrho = rho
iter_count = iter_count + 1
if ( iter_count .gt. HUTI_MAXIT ) then
HUTI_INFO = HUTI_MAXITER
go to 1000
end if
go to 300
1000 continue
IF ( HUTI_DBUGLVL .NE. HUTI_NO_DEBUG ) THEN
WRITE (*, '(I8, E11.4)') iter_count, residual
CALL FLUSH(6)
END IF
HUTI_ITERS = iter_count
return
end subroutine huti_dcgsolv
recursive subroutine huti_ccgsolv ( ndim, wrkdim, xvec, rhsvec, &
ipar, dpar, work, matvecsubr, pcondlsubr, &
pcondrsubr, dotprodfun, normfun, stopcfun )
use huti_interfaces
implicit none
procedure( mv_iface_c ), pointer :: matvecsubr
procedure( pc_iface_c ), pointer :: pcondlsubr
procedure( pc_iface_c ), pointer :: pcondrsubr
procedure( dotp_iface_c ), pointer :: dotprodfun
procedure( norm_iface_c ), pointer :: normfun
procedure( stopc_iface_c ), pointer :: stopcfun
integer :: ndim, wrkdim
complex, dimension(ndim) :: xvec, rhsvec
integer, dimension(HUTI_IPAR_DFLTSIZE) :: ipar
double precision, dimension(HUTI_DPAR_DFLTSIZE) :: dpar
complex, dimension(ndim,wrkdim) :: work
complex :: alpha, beta, rho, oldrho
integer iter_count
real :: residual, rhsnorm, precrhsnorm
iter_count = 1
if ( HUTI_STOPC .eq. HUTI_TRESID_SCALED_BYB .or. &
HUTI_STOPC .eq. HUTI_PRESID_SCALED_BYB ) then
rhsnorm = normfun( HUTI_NDIM, B, 1 )
end if
if ( HUTI_STOPC .eq. HUTI_PRESID_SCALED_BYPRECB ) then
call pcondlsubr( P, B, ipar )
precrhsnorm = normfun( HUTI_NDIM, P, 1 )
end if
HUTI_EXTOP_MATTYPE = HUTI_MAT_NOTTRPSED
if ( HUTI_INITIALX .eq. HUTI_RANDOMX ) then
call huti_crandvec ( X, ipar )
else if ( HUTI_INITIALX .ne. HUTI_USERSUPPLIEDX ) then
X = 1
end if
call matvecsubr( X, R, ipar )
R = B - R
300 continue
call pcondlsubr( Q, R, ipar )
call pcondrsubr( Z, Q, ipar )
rho = dotprodfun( HUTI_NDIM, R, 1, Z, 1 )
if ( rho .eq. 0 ) then
HUTI_INFO = HUTI_CG_RHO
go to 1000
end if
if ( iter_count .eq. 1 ) then
P = Z
else
beta = rho / oldrho
P = Z + beta * P
end if
call matvecsubr( P, Q, ipar )
alpha = rho / dotprodfun( HUTI_NDIM, P, 1, Q, 1 )
X = X + alpha * P
R = R - alpha * Q
select case (HUTI_STOPC)
case (HUTI_TRUERESIDUAL)
call matvecsubr( X, Z, ipar )
Z = Z - B
residual = normfun( HUTI_NDIM, Z, 1 )
case (HUTI_TRESID_SCALED_BYB)
call matvecsubr( X, Z, ipar )
Z = Z - B
residual = normfun( HUTI_NDIM, Z, 1 ) / rhsnorm
case (HUTI_PSEUDORESIDUAL)
residual = normfun( HUTI_NDIM, R, 1 )
case (HUTI_PRESID_SCALED_BYB)
residual = normfun( HUTI_NDIM, R, 1 ) / rhsnorm
case (HUTI_PRESID_SCALED_BYPRECB)
residual = normfun( HUTI_NDIM, R, 1 ) / precrhsnorm
case (HUTI_XDIFF_NORM)
Z = alpha * P
residual = normfun( HUTI_NDIM, Z, 1 )
case (HUTI_USUPPLIED_STOPC)
residual = stopcfun( X, B, R, ipar, dpar )
case default
call matvecsubr( X, Z, ipar )
Z = Z - B
residual = normfun( HUTI_NDIM, Z, 1 )
end select
if ( HUTI_DBUGLVL .ne. HUTI_NO_DEBUG ) then
if ( mod(iter_count, HUTI_DBUGLVL) .eq. 0 ) then
WRITE (*, '(I8, E11.4)') iter_count, residual
CALL FLUSH(6)
end if
end if
if ( residual .lt. HUTI_TOLERANCE ) then
HUTI_INFO = HUTI_CONVERGENCE
go to 1000
end if
IF( residual /= residual .OR. residual > HUTI_MAXTOLERANCE ) THEN
HUTI_INFO = HUTI_DIVERGENCE
GOTO 1000
END IF
oldrho = rho
iter_count = iter_count + 1
if ( iter_count .gt. HUTI_MAXIT ) then
HUTI_INFO = HUTI_MAXITER
go to 1000
end if
go to 300
1000 continue
if ( HUTI_DBUGLVL .ne. HUTI_NO_DEBUG ) then
WRITE (*, '(I8, E11.4)') iter_count, residual
call flush(6)
end if
HUTI_ITERS = iter_count
return
end subroutine huti_ccgsolv
recursive subroutine huti_zcgsolv ( ndim, wrkdim, xvec, rhsvec, &
ipar, dpar, work, matvecsubr, pcondlsubr, &
pcondrsubr, dotprodfun, normfun, stopcfun )
use huti_interfaces
implicit none
procedure( mv_iface_z ), pointer :: matvecsubr
procedure( pc_iface_z ), pointer :: pcondlsubr
procedure( pc_iface_z ), pointer :: pcondrsubr
procedure( dotp_iface_z ), pointer :: dotprodfun
procedure( norm_iface_z ), pointer :: normfun
procedure( stopc_iface_z ), pointer :: stopcfun
integer :: ndim, wrkdim
double complex, dimension(ndim) :: xvec, rhsvec
integer, dimension(HUTI_IPAR_DFLTSIZE) :: ipar
double precision, dimension(HUTI_DPAR_DFLTSIZE) :: dpar
double complex, dimension(ndim,wrkdim) :: work
double complex :: alpha, beta, rho, oldrho
integer iter_count
double precision :: residual, rhsnorm, precrhsnorm
iter_count = 1
if ( HUTI_STOPC .eq. HUTI_TRESID_SCALED_BYB .or. &
HUTI_STOPC .eq. HUTI_PRESID_SCALED_BYB ) then
rhsnorm = normfun( HUTI_NDIM, B, 1 )
end if
if ( HUTI_STOPC .eq. HUTI_PRESID_SCALED_BYPRECB ) then
call pcondlsubr( P, B, ipar )
precrhsnorm = normfun( HUTI_NDIM, P, 1 )
end if
HUTI_EXTOP_MATTYPE = HUTI_MAT_NOTTRPSED
if ( HUTI_INITIALX .eq. HUTI_RANDOMX ) then
call huti_zrandvec ( X, ipar )
else if ( HUTI_INITIALX .ne. HUTI_USERSUPPLIEDX ) then
X = 1
end if
call matvecsubr( X, R, ipar )
R = B - R
300 continue
call pcondlsubr( Q, R, ipar )
call pcondrsubr( Z, Q, ipar )
rho = dotprodfun( HUTI_NDIM, R, 1, Z, 1 )
if ( rho .eq. 0 ) then
HUTI_INFO = HUTI_CG_RHO
go to 1000
end if
if ( iter_count .eq. 1 ) then
P = Z
else
beta = rho / oldrho
P = Z + beta * P
end if
call matvecsubr( P, Q, ipar )
alpha = rho / dotprodfun( HUTI_NDIM, P, 1, Q, 1 )
X = X + alpha * P
R = R - alpha * Q
select case (HUTI_STOPC)
case (HUTI_TRUERESIDUAL)
call matvecsubr( X, Z, ipar )
Z = Z - B
residual = normfun( HUTI_NDIM, Z, 1 )
case (HUTI_TRESID_SCALED_BYB)
call matvecsubr( X, Z, ipar )
Z = Z - B
residual = normfun( HUTI_NDIM, Z, 1 ) / rhsnorm
case (HUTI_PSEUDORESIDUAL)
residual = normfun( HUTI_NDIM, R, 1 )
case (HUTI_PRESID_SCALED_BYB)
residual = normfun( HUTI_NDIM, R, 1 ) / rhsnorm
case (HUTI_PRESID_SCALED_BYPRECB)
residual = normfun( HUTI_NDIM, R, 1 ) / precrhsnorm
case (HUTI_XDIFF_NORM)
Z = alpha * P
residual = normfun( HUTI_NDIM, Z, 1 )
case (HUTI_USUPPLIED_STOPC)
residual = stopcfun( X, B, R, ipar, dpar )
case default
call matvecsubr( X, Z, ipar )
Z = Z - B
residual = normfun( HUTI_NDIM, Z, 1 )
end select
if ( HUTI_DBUGLVL .ne. HUTI_NO_DEBUG ) then
if ( mod(iter_count, HUTI_DBUGLVL) .eq. 0 ) then
WRITE (*, '(I8, E11.4)') iter_count, residual
call flush(6)
end if
end if
if ( residual .lt. HUTI_TOLERANCE ) then
HUTI_INFO = HUTI_CONVERGENCE
go to 1000
end if
IF( residual /= residual .OR. residual > HUTI_MAXTOLERANCE ) THEN
HUTI_INFO = HUTI_DIVERGENCE
GOTO 1000
END IF
oldrho = rho
iter_count = iter_count + 1
if ( iter_count .gt. HUTI_MAXIT ) then
HUTI_INFO = HUTI_MAXITER
go to 1000
end if
go to 300
1000 continue
if ( HUTI_DBUGLVL .ne. HUTI_NO_DEBUG ) then
WRITE (*, '(I8, E11.4)') iter_count, residual
CALL FLUSH(6)
end if
HUTI_ITERS = iter_count
return
end subroutine huti_zcgsolv
end module huti_cg