Fixed build and apply to actually use degree

Poly-novrl
Cirdans-Home 11 months ago
parent 54d608d2dd
commit ea8974f88c

@ -49,7 +49,7 @@ subroutine amg_d_poly_smoother_apply_vect(alpha,sm,x,beta,y,desc_data,trans,&
type(psb_d_vect_type),intent(inout) :: y
real(psb_dpk_),intent(in) :: alpha,beta
character(len=1),intent(in) :: trans
integer(psb_ipk_), intent(in) :: sweeps
integer(psb_ipk_), intent(in) :: sweeps ! this is ignored here, the polynomial degree dictates the value
real(psb_dpk_),target, intent(inout) :: work(:)
type(psb_d_vect_type),intent(inout) :: wv(:)
integer(psb_ipk_), intent(out) :: info
@ -115,22 +115,22 @@ subroutine amg_d_poly_smoother_apply_vect(alpha,sm,x,beta,y,desc_data,trans,&
& a_err='invalid wv size in smoother_apply')
goto 9999
end if
sm%pdegree = sweeps
associate(tx => wv(1), ty => wv(2), tz => wv(3), r => wv(4))
call psb_geaxpby(done,x,dzero,r,desc_data,info)
call tx%zero()
call ty%zero()
call tz%zero()
call tz%zero()
select case(sm%variant)
case(amg_poly_lottes_)
block
block
real(psb_dpk_) :: cz, cr
! b == x
! b == x
! x == tx
!
do i=1, sweeps
!
do i=1, sm%pdegree
! B r_{k-1}
call sm%sv%apply(done,r,dzero,ty,desc_data,trans_,aux,wv(5:),info,init='Z')
cz = (2*i*done-3)/(2*i*done+done)
@ -153,20 +153,20 @@ subroutine amg_d_poly_smoother_apply_vect(alpha,sm,x,beta,y,desc_data,trans,&
case(amg_poly_lottes_beta_)
block
block
real(psb_dpk_) :: cz, cr
! b == x
! b == x
! x == tx
!
if (allocated(sm%poly_beta)) then
if (size(sm%poly_beta) /= sweeps) deallocate(sm%poly_beta)
if (size(sm%poly_beta) /= sm%pdegree) deallocate(sm%poly_beta)
end if
if (.not.allocated(sm%poly_beta)) then
call psb_realloc(sweeps,sm%poly_beta,info)
sm%poly_beta(1:sweeps) = amg_d_poly_beta_mat(1:sweeps,sweeps)
call psb_realloc(sm%pdegree,sm%poly_beta,info)
sm%poly_beta(1:sm%pdegree) = amg_d_poly_beta_mat(1:sm%pdegree,sm%pdegree)
end if
do i=1, sweeps
do i=1, sm%pdegree
! B r_{k-1}
call sm%sv%apply(done,r,dzero,ty,desc_data,trans_,aux,wv(5:),info,init='Z')
cz = (2*i*done-3)/(2*i*done+done)
@ -186,14 +186,14 @@ subroutine amg_d_poly_smoother_apply_vect(alpha,sm,x,beta,y,desc_data,trans,&
! x_k = x_{k-1} + z_k
end do
end block
case(amg_poly_new_)
block
block
real(psb_dpk_) :: sigma, theta, delta, rho_old, rho
! b == x
! b == x
! x == tx
!
sm%cf_a = amg_d_poly_a_vect(sweeps)
theta = (done+sm%cf_a)/2
delta = (done-sm%cf_a)/2
@ -203,10 +203,10 @@ subroutine amg_d_poly_smoother_apply_vect(alpha,sm,x,beta,y,desc_data,trans,&
call psb_geaxpby((done/sm%rho_ba),ty,dzero,r,desc_data,info)
call psb_geaxpby((done/theta),r,dzero,tz,desc_data,info)
! tz == d
do i=1, sweeps
do i=1, sm%pdegree
! x_{k+1} = x_k + d_k
call psb_geaxpby(done,tz,done,tx,desc_data,info)
!
!
! r_{k-1} = r_k - (1/rho(BA)) B A d_k
call psb_spmm(done,sm%pa,tz,dzero,ty,desc_data,info,work=aux,trans=trans_)
call sm%sv%apply(-done,ty,done,r,desc_data,trans_,aux,wv(5:),info,init='Z')

@ -56,7 +56,7 @@ subroutine amg_d_poly_smoother_bld(a,desc_a,sm,info,amold,vmold,imold)
type(psb_dspmat_type) :: tmpa
integer(psb_ipk_) :: n_row,n_col, nrow_a, nztota, nzeros
type(psb_ctxt_type) :: ctxt
real(psb_dpk_), allocatable :: da(:), dsv(:)
real(psb_dpk_), allocatable :: da(:), dsv(:)
integer(psb_ipk_) :: np, me, i, err_act, debug_unit, debug_level
character(len=20) :: name='d_poly_smoother_bld', ch_err
@ -74,8 +74,7 @@ subroutine amg_d_poly_smoother_bld(a,desc_a,sm,info,amold,vmold,imold)
n_col = desc_a%get_local_cols()
nrow_a = a%get_nrows()
nztota = a%get_nzeros()
if (.false.) then
select case(sm%variant)
select case(sm%variant)
case(amg_poly_lottes_)
! do nothing
case(amg_poly_lottes_beta_)
@ -89,7 +88,7 @@ subroutine amg_d_poly_smoother_bld(a,desc_a,sm,info,amold,vmold,imold)
goto 9999
end if
case(amg_poly_new_)
write(psb_out_unit,*) "Nella fase di build sm%pdegree = ",sm%pdegree
if ((1<=sm%pdegree).and.(sm%pdegree<=30)) then
!Ok
sm%cf_a = amg_d_poly_a_vect(sm%pdegree)
@ -100,15 +99,15 @@ subroutine amg_d_poly_smoother_bld(a,desc_a,sm,info,amold,vmold,imold)
goto 9999
end if
case default
case default
info = psb_err_internal_error_
call psb_errpush(info,name,&
& a_err='invalid sm%variant')
goto 9999
end select
end if
sm%pa => a
if (.not.allocated(sm%sv)) then
if (.not.allocated(sm%sv)) then
info = psb_err_internal_error_
call psb_errpush(info,name,&
& a_err='unallocated sm%sv')
@ -121,7 +120,7 @@ subroutine amg_d_poly_smoother_bld(a,desc_a,sm,info,amold,vmold,imold)
goto 9999
end if
!!$ if (.false.) then
!!$ if (.false.) then
!!$ select type(ssv => sm%sv)
!!$ class is(amg_d_l1_diag_solver_type)
!!$ da = a%arwsum(info)
@ -129,7 +128,7 @@ subroutine amg_d_poly_smoother_bld(a,desc_a,sm,info,amold,vmold,imold)
!!$ sm%rho_ba = maxval(da(1:n_row)*dsv(1:n_row))
!!$ class default
!!$ write(0,*) 'PolySmoother BUILD: only L1-Jacobi/L1-DIAG for now ',ssv%get_fmt()
!!$ sm%rho_ba = done
!!$ sm%rho_ba = done
!!$ end select
!!$ else
if (sm%rho_ba <= dzero) then
@ -153,9 +152,9 @@ subroutine amg_d_poly_smoother_bld(a,desc_a,sm,info,amold,vmold,imold)
call sm%sv%apply_v(done,tt,dzero,tz,desc_a,'NoTrans',work,wv,info) ! z_{k+1} = BA q_k
do i=1,sm%rho_estimate_iterations
znrm = psb_genrm2(tz,desc_a,info) ! znrm = |z_k|_2
call psb_geaxpby((done/znrm),tz,dzero,tq,desc_a,info) ! q_k = z_k/znrm
call psb_geaxpby((done/znrm),tz,dzero,tq,desc_a,info) ! q_k = z_k/znrm
call psb_spmm(done,a,tq,dzero,tt,desc_a,info) ! t_{k+1} = BA q_k
call sm%sv%apply_v(done,tt,dzero,tz,desc_a,'NoTrans',work,wv,info) ! z_{k+1} = B t_{k+1}
call sm%sv%apply_v(done,tt,dzero,tz,desc_a,'NoTrans',work,wv,info) ! z_{k+1} = B t_{k+1}
lambda = psb_gedot(tq,tz,desc_a,info) ! lambda = q_k^T z_{k+1} = q_k^T BA q_k
!write(0,*) 'BLD: lambda estimate ',i,lambda
end do

@ -1,15 +1,15 @@
!
!
!
!
! AMG4PSBLAS version 1.0
! Algebraic Multigrid Package
! based on PSBLAS (Parallel Sparse BLAS version 3.7)
!
! (C) Copyright 2021
!
! Salvatore Filippone
! Pasqua D'Ambra
! Fabio Durastante
!
!
! (C) Copyright 2021
!
! Salvatore Filippone
! Pasqua D'Ambra
! Fabio Durastante
!
! Redistribution and use in source and binary forms, with or without
! modification, are permitted provided that the following conditions
! are met:
@ -21,7 +21,7 @@
! 3. The name of the AMG4PSBLAS group or the names of its contributors may
! not be used to endorse or promote products derived from this
! software without specific written permission.
!
!
! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
! ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
! TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
@ -33,8 +33,8 @@
! CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
! ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
! POSSIBILITY OF SUCH DAMAGE.
!
!
!
!
subroutine amg_d_poly_smoother_descr(sm,info,iout,coarse,prefix)
use psb_base_mod
@ -59,13 +59,13 @@ subroutine amg_d_poly_smoother_descr(sm,info,iout,coarse,prefix)
call psb_erractionsave(err_act)
info = psb_success_
if (present(coarse)) then
if (present(coarse)) then
coarse_ = coarse
else
coarse_ = .false.
end if
if (present(iout)) then
iout_ = iout
if (present(iout)) then
iout_ = iout
else
iout_ = psb_out_unit
endif
@ -78,19 +78,19 @@ subroutine amg_d_poly_smoother_descr(sm,info,iout,coarse,prefix)
write(iout_,*) trim(prefix_), ' Polynomial smoother '
select case(sm%variant)
case(amg_poly_lottes_)
write(iout_,*) trim(prefix_), ' variant: ','POLY_LOTTES'
!write(iout_,*) trim(prefix_), ' Degree: ',sm%pdegree
write(iout_,*) trim(prefix_), ' variant: ','POLY_LOTTES'
write(iout_,*) trim(prefix_), ' Degree: ',sm%pdegree
write(iout_,*) trim(prefix_), ' rho_ba: ',sm%rho_ba
case(amg_poly_lottes_beta_)
write(iout_,*) trim(prefix_), ' variant: ','POLY_LOTTES_BETA'
!write(iout_,*) trim(prefix_), ' Degree: ',sm%pdegree
write(iout_,*) trim(prefix_), ' Degree: ',sm%pdegree
write(iout_,*) trim(prefix_), ' rho_ba: ',sm%rho_ba
!if (allocated(sm%poly_beta)) write(iout_,*) trim(prefix_), ' Coefficients: ',sm%poly_beta(1:sm%pdegree)
if (allocated(sm%poly_beta)) write(iout_,*) trim(prefix_), ' Coefficients: ',sm%poly_beta(1:sm%pdegree)
case(amg_poly_new_)
write(iout_,*) trim(prefix_), ' variant: ','POLY_NEW'
!write(iout_,*) trim(prefix_), ' Degree: ',sm%pdegree
write(iout_,*) trim(prefix_), ' Degree: ',sm%pdegree
write(iout_,*) trim(prefix_), ' rho_ba: ',sm%rho_ba
!write(iout_,*) trim(prefix_), ' Coefficient: ',sm%cf_a
write(iout_,*) trim(prefix_), ' Coefficient: ',sm%cf_a
case default
write(iout_,*) trim(prefix_), ' variant: ','UNKNOWN???'
end select

@ -49,7 +49,7 @@ subroutine amg_s_poly_smoother_apply_vect(alpha,sm,x,beta,y,desc_data,trans,&
type(psb_s_vect_type),intent(inout) :: y
real(psb_spk_),intent(in) :: alpha,beta
character(len=1),intent(in) :: trans
integer(psb_ipk_), intent(in) :: sweeps
integer(psb_ipk_), intent(in) :: sweeps ! this is ignored here, the polynomial degree dictates the value
real(psb_spk_),target, intent(inout) :: work(:)
type(psb_s_vect_type),intent(inout) :: wv(:)
integer(psb_ipk_), intent(out) :: info
@ -115,22 +115,22 @@ subroutine amg_s_poly_smoother_apply_vect(alpha,sm,x,beta,y,desc_data,trans,&
& a_err='invalid wv size in smoother_apply')
goto 9999
end if
sm%pdegree = sweeps
associate(tx => wv(1), ty => wv(2), tz => wv(3), r => wv(4))
call psb_geaxpby(sone,x,szero,r,desc_data,info)
call tx%zero()
call ty%zero()
call tz%zero()
call tz%zero()
select case(sm%variant)
case(amg_poly_lottes_)
block
block
real(psb_spk_) :: cz, cr
! b == x
! b == x
! x == tx
!
do i=1, sweeps
!
do i=1, sm%pdegree
! B r_{k-1}
call sm%sv%apply(sone,r,szero,ty,desc_data,trans_,aux,wv(5:),info,init='Z')
cz = (2*i*sone-3)/(2*i*sone+sone)
@ -153,20 +153,20 @@ subroutine amg_s_poly_smoother_apply_vect(alpha,sm,x,beta,y,desc_data,trans,&
case(amg_poly_lottes_beta_)
block
block
real(psb_spk_) :: cz, cr
! b == x
! b == x
! x == tx
!
if (allocated(sm%poly_beta)) then
if (size(sm%poly_beta) /= sweeps) deallocate(sm%poly_beta)
if (size(sm%poly_beta) /= sm%pdegree) deallocate(sm%poly_beta)
end if
if (.not.allocated(sm%poly_beta)) then
call psb_realloc(sweeps,sm%poly_beta,info)
sm%poly_beta(1:sweeps) = amg_d_poly_beta_mat(1:sweeps,sweeps)
call psb_realloc(sm%pdegree,sm%poly_beta,info)
sm%poly_beta(1:sm%pdegree) = amg_d_poly_beta_mat(1:sm%pdegree,sm%pdegree)
end if
do i=1, sweeps
do i=1, sm%pdegree
! B r_{k-1}
call sm%sv%apply(sone,r,szero,ty,desc_data,trans_,aux,wv(5:),info,init='Z')
cz = (2*i*sone-3)/(2*i*sone+sone)
@ -186,14 +186,14 @@ subroutine amg_s_poly_smoother_apply_vect(alpha,sm,x,beta,y,desc_data,trans,&
! x_k = x_{k-1} + z_k
end do
end block
case(amg_poly_new_)
block
block
real(psb_spk_) :: sigma, theta, delta, rho_old, rho
! b == x
! b == x
! x == tx
!
sm%cf_a = amg_d_poly_a_vect(sweeps)
theta = (sone+sm%cf_a)/2
delta = (sone-sm%cf_a)/2
@ -203,10 +203,10 @@ subroutine amg_s_poly_smoother_apply_vect(alpha,sm,x,beta,y,desc_data,trans,&
call psb_geaxpby((sone/sm%rho_ba),ty,szero,r,desc_data,info)
call psb_geaxpby((sone/theta),r,szero,tz,desc_data,info)
! tz == d
do i=1, sweeps
do i=1, sm%pdegree
! x_{k+1} = x_k + d_k
call psb_geaxpby(sone,tz,sone,tx,desc_data,info)
!
!
! r_{k-1} = r_k - (1/rho(BA)) B A d_k
call psb_spmm(sone,sm%pa,tz,szero,ty,desc_data,info,work=aux,trans=trans_)
call sm%sv%apply(-sone,ty,sone,r,desc_data,trans_,aux,wv(5:),info,init='Z')

@ -56,7 +56,7 @@ subroutine amg_s_poly_smoother_bld(a,desc_a,sm,info,amold,vmold,imold)
type(psb_sspmat_type) :: tmpa
integer(psb_ipk_) :: n_row,n_col, nrow_a, nztota, nzeros
type(psb_ctxt_type) :: ctxt
real(psb_spk_), allocatable :: da(:), dsv(:)
real(psb_spk_), allocatable :: da(:), dsv(:)
integer(psb_ipk_) :: np, me, i, err_act, debug_unit, debug_level
character(len=20) :: name='d_poly_smoother_bld', ch_err
@ -74,8 +74,7 @@ subroutine amg_s_poly_smoother_bld(a,desc_a,sm,info,amold,vmold,imold)
n_col = desc_a%get_local_cols()
nrow_a = a%get_nrows()
nztota = a%get_nzeros()
if (.false.) then
select case(sm%variant)
select case(sm%variant)
case(amg_poly_lottes_)
! do nothing
case(amg_poly_lottes_beta_)
@ -89,7 +88,7 @@ subroutine amg_s_poly_smoother_bld(a,desc_a,sm,info,amold,vmold,imold)
goto 9999
end if
case(amg_poly_new_)
write(psb_out_unit,*) "Nella fase di build sm%pdegree = ",sm%pdegree
if ((1<=sm%pdegree).and.(sm%pdegree<=30)) then
!Ok
sm%cf_a = amg_d_poly_a_vect(sm%pdegree)
@ -100,15 +99,15 @@ subroutine amg_s_poly_smoother_bld(a,desc_a,sm,info,amold,vmold,imold)
goto 9999
end if
case default
case default
info = psb_err_internal_error_
call psb_errpush(info,name,&
& a_err='invalid sm%variant')
goto 9999
end select
end if
sm%pa => a
if (.not.allocated(sm%sv)) then
if (.not.allocated(sm%sv)) then
info = psb_err_internal_error_
call psb_errpush(info,name,&
& a_err='unallocated sm%sv')
@ -121,7 +120,7 @@ subroutine amg_s_poly_smoother_bld(a,desc_a,sm,info,amold,vmold,imold)
goto 9999
end if
!!$ if (.false.) then
!!$ if (.false.) then
!!$ select type(ssv => sm%sv)
!!$ class is(amg_s_l1_diag_solver_type)
!!$ da = a%arwsum(info)
@ -129,7 +128,7 @@ subroutine amg_s_poly_smoother_bld(a,desc_a,sm,info,amold,vmold,imold)
!!$ sm%rho_ba = maxval(da(1:n_row)*dsv(1:n_row))
!!$ class default
!!$ write(0,*) 'PolySmoother BUILD: only L1-Jacobi/L1-DIAG for now ',ssv%get_fmt()
!!$ sm%rho_ba = sone
!!$ sm%rho_ba = sone
!!$ end select
!!$ else
if (sm%rho_ba <= szero) then
@ -153,9 +152,9 @@ subroutine amg_s_poly_smoother_bld(a,desc_a,sm,info,amold,vmold,imold)
call sm%sv%apply_v(sone,tt,szero,tz,desc_a,'NoTrans',work,wv,info) ! z_{k+1} = BA q_k
do i=1,sm%rho_estimate_iterations
znrm = psb_genrm2(tz,desc_a,info) ! znrm = |z_k|_2
call psb_geaxpby((sone/znrm),tz,szero,tq,desc_a,info) ! q_k = z_k/znrm
call psb_geaxpby((sone/znrm),tz,szero,tq,desc_a,info) ! q_k = z_k/znrm
call psb_spmm(sone,a,tq,szero,tt,desc_a,info) ! t_{k+1} = BA q_k
call sm%sv%apply_v(sone,tt,szero,tz,desc_a,'NoTrans',work,wv,info) ! z_{k+1} = B t_{k+1}
call sm%sv%apply_v(sone,tt,szero,tz,desc_a,'NoTrans',work,wv,info) ! z_{k+1} = B t_{k+1}
lambda = psb_gedot(tq,tz,desc_a,info) ! lambda = q_k^T z_{k+1} = q_k^T BA q_k
!write(0,*) 'BLD: lambda estimate ',i,lambda
end do

@ -1,15 +1,15 @@
!
!
!
!
! AMG4PSBLAS version 1.0
! Algebraic Multigrid Package
! based on PSBLAS (Parallel Sparse BLAS version 3.7)
!
! (C) Copyright 2021
!
! Salvatore Filippone
! Pasqua D'Ambra
! Fabio Durastante
!
!
! (C) Copyright 2021
!
! Salvatore Filippone
! Pasqua D'Ambra
! Fabio Durastante
!
! Redistribution and use in source and binary forms, with or without
! modification, are permitted provided that the following conditions
! are met:
@ -21,7 +21,7 @@
! 3. The name of the AMG4PSBLAS group or the names of its contributors may
! not be used to endorse or promote products derived from this
! software without specific written permission.
!
!
! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
! ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
! TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
@ -33,8 +33,8 @@
! CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
! ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
! POSSIBILITY OF SUCH DAMAGE.
!
!
!
!
subroutine amg_s_poly_smoother_descr(sm,info,iout,coarse,prefix)
use psb_base_mod
@ -59,13 +59,13 @@ subroutine amg_s_poly_smoother_descr(sm,info,iout,coarse,prefix)
call psb_erractionsave(err_act)
info = psb_success_
if (present(coarse)) then
if (present(coarse)) then
coarse_ = coarse
else
coarse_ = .false.
end if
if (present(iout)) then
iout_ = iout
if (present(iout)) then
iout_ = iout
else
iout_ = psb_out_unit
endif
@ -78,19 +78,19 @@ subroutine amg_s_poly_smoother_descr(sm,info,iout,coarse,prefix)
write(iout_,*) trim(prefix_), ' Polynomial smoother '
select case(sm%variant)
case(amg_poly_lottes_)
write(iout_,*) trim(prefix_), ' variant: ','POLY_LOTTES'
!write(iout_,*) trim(prefix_), ' Degree: ',sm%pdegree
write(iout_,*) trim(prefix_), ' variant: ','POLY_LOTTES'
write(iout_,*) trim(prefix_), ' Degree: ',sm%pdegree
write(iout_,*) trim(prefix_), ' rho_ba: ',sm%rho_ba
case(amg_poly_lottes_beta_)
write(iout_,*) trim(prefix_), ' variant: ','POLY_LOTTES_BETA'
!write(iout_,*) trim(prefix_), ' Degree: ',sm%pdegree
write(iout_,*) trim(prefix_), ' Degree: ',sm%pdegree
write(iout_,*) trim(prefix_), ' rho_ba: ',sm%rho_ba
!if (allocated(sm%poly_beta)) write(iout_,*) trim(prefix_), ' Coefficients: ',sm%poly_beta(1:sm%pdegree)
if (allocated(sm%poly_beta)) write(iout_,*) trim(prefix_), ' Coefficients: ',sm%poly_beta(1:sm%pdegree)
case(amg_poly_new_)
write(iout_,*) trim(prefix_), ' variant: ','POLY_NEW'
!write(iout_,*) trim(prefix_), ' Degree: ',sm%pdegree
write(iout_,*) trim(prefix_), ' Degree: ',sm%pdegree
write(iout_,*) trim(prefix_), ' rho_ba: ',sm%rho_ba
!write(iout_,*) trim(prefix_), ' Coefficient: ',sm%cf_a
write(iout_,*) trim(prefix_), ' Coefficient: ',sm%cf_a
case default
write(iout_,*) trim(prefix_), ' variant: ','UNKNOWN???'
end select

Loading…
Cancel
Save