Merge branch 'development' into documentation

documentation
Salvatore Filippone 4 years ago
commit 77fd14b89b

@ -120,7 +120,9 @@ module amg_base_prec_type
procedure, pass(pm) :: printout => d_ml_parms_printout
end type amg_dml_parms
type amg_saggr_data
type amg_iaggr_data
!
! Aggregation data and defaults:
!
@ -129,28 +131,21 @@ module amg_base_prec_type
! We are assuming that the coarse size fits in
! integer range of psb_ipk_, but this is
! not very restrictive
integer(psb_ipk_) :: min_coarse_size = izero
integer(psb_ipk_) :: min_coarse_size = -ione
integer(psb_ipk_) :: min_coarse_size_per_process = -ione
integer(psb_lpk_) :: target_coarse_size
! 2. maximum number of levels. Defaults to 20
integer(psb_ipk_) :: max_levs = 20_psb_ipk_
end type amg_iaggr_data
type, extends(amg_iaggr_data) :: amg_saggr_data
! 3. min_cr_ratio = 1.5
real(psb_spk_) :: min_cr_ratio = 1.5_psb_spk_
real(psb_spk_) :: op_complexity = szero
real(psb_spk_) :: avg_cr = szero
end type amg_saggr_data
type amg_daggr_data
!
! Aggregation data and defaults:
!
!
! 1. min_coarse_size = 0 Default target size will be computed as
! 40*(N_fine)**(1./3.)
! We are assuming that the coarse size fits in
! integer range of psb_ipk_, but this is
! not very restrictive
integer(psb_ipk_) :: min_coarse_size = izero
! 2. maximum number of levels. Defaults to 20
integer(psb_ipk_) :: max_levs = 20_psb_ipk_
type, extends(amg_iaggr_data) :: amg_daggr_data
! 3. min_cr_ratio = 1.5
real(psb_dpk_) :: min_cr_ratio = 1.5_psb_dpk_
real(psb_dpk_) :: op_complexity = dzero

@ -425,11 +425,22 @@ contains
end if
end function amg_c_get_nzeros
function amg_cprec_sizeof(prec) result(val)
function amg_cprec_sizeof(prec, global) result(val)
implicit none
class(amg_cprec_type), intent(in) :: prec
logical, intent(in), optional :: global
integer(psb_epk_) :: val
integer(psb_ipk_) :: i
type(psb_ctxt_type) :: ctxt
logical :: global_
if (present(global)) then
global_ = global
else
global_ = .false.
end if
val = 0
val = val + psb_sizeof_ip
if (allocated(prec%precv)) then
@ -437,6 +448,11 @@ contains
val = val + prec%precv(i)%sizeof()
end do
end if
if (global_) then
ctxt = prec%ctxt
call psb_sum(ctxt,val)
end if
end function amg_cprec_sizeof
!

@ -425,11 +425,22 @@ contains
end if
end function amg_d_get_nzeros
function amg_dprec_sizeof(prec) result(val)
function amg_dprec_sizeof(prec, global) result(val)
implicit none
class(amg_dprec_type), intent(in) :: prec
logical, intent(in), optional :: global
integer(psb_epk_) :: val
integer(psb_ipk_) :: i
type(psb_ctxt_type) :: ctxt
logical :: global_
if (present(global)) then
global_ = global
else
global_ = .false.
end if
val = 0
val = val + psb_sizeof_ip
if (allocated(prec%precv)) then
@ -437,6 +448,11 @@ contains
val = val + prec%precv(i)%sizeof()
end do
end if
if (global_) then
ctxt = prec%ctxt
call psb_sum(ctxt,val)
end if
end function amg_dprec_sizeof
!

@ -425,11 +425,22 @@ contains
end if
end function amg_s_get_nzeros
function amg_sprec_sizeof(prec) result(val)
function amg_sprec_sizeof(prec, global) result(val)
implicit none
class(amg_sprec_type), intent(in) :: prec
logical, intent(in), optional :: global
integer(psb_epk_) :: val
integer(psb_ipk_) :: i
type(psb_ctxt_type) :: ctxt
logical :: global_
if (present(global)) then
global_ = global
else
global_ = .false.
end if
val = 0
val = val + psb_sizeof_ip
if (allocated(prec%precv)) then
@ -437,6 +448,11 @@ contains
val = val + prec%precv(i)%sizeof()
end do
end if
if (global_) then
ctxt = prec%ctxt
call psb_sum(ctxt,val)
end if
end function amg_sprec_sizeof
!

@ -425,11 +425,22 @@ contains
end if
end function amg_z_get_nzeros
function amg_zprec_sizeof(prec) result(val)
function amg_zprec_sizeof(prec, global) result(val)
implicit none
class(amg_zprec_type), intent(in) :: prec
logical, intent(in), optional :: global
integer(psb_epk_) :: val
integer(psb_ipk_) :: i
type(psb_ctxt_type) :: ctxt
logical :: global_
if (present(global)) then
global_ = global
else
global_ = .false.
end if
val = 0
val = val + psb_sizeof_ip
if (allocated(prec%precv)) then
@ -437,6 +448,11 @@ contains
val = val + prec%precv(i)%sizeof()
end do
end if
if (global_) then
ctxt = prec%ctxt
call psb_sum(ctxt,val)
end if
end function amg_zprec_sizeof
!

@ -82,7 +82,7 @@ subroutine amg_c_hierarchy_bld(a,desc_a,prec,info)
integer(psb_ipk_) :: me,np
integer(psb_ipk_) :: err,i,k, err_act, iszv, newsz,&
& nplevs, mxplevs
integer(psb_lpk_) :: iaggsize, casize
integer(psb_lpk_) :: iaggsize, casize, mncsize, mncszpp
real(psb_spk_) :: mnaggratio, sizeratio, athresh, aomega
class(amg_c_base_smoother_type), allocatable :: coarse_sm, med_sm, &
& med_sm2, coarse_sm2
@ -132,17 +132,24 @@ subroutine amg_c_hierarchy_bld(a,desc_a,prec,info)
newsz = -1
mxplevs = prec%ag_data%max_levs
mnaggratio = prec%ag_data%min_cr_ratio
casize = prec%ag_data%min_coarse_size
mncsize = prec%ag_data%min_coarse_size
mncszpp = prec%ag_data%min_coarse_size_per_process
iszv = size(prec%precv)
call psb_bcast(ctxt,iszv)
call psb_bcast(ctxt,casize)
call psb_bcast(ctxt,mncsize)
call psb_bcast(ctxt,mncszpp)
call psb_bcast(ctxt,mxplevs)
call psb_bcast(ctxt,mnaggratio)
if (casize /= prec%ag_data%min_coarse_size) then
if (mncsize /= prec%ag_data%min_coarse_size) then
info=psb_err_internal_error_
call psb_errpush(info,name,a_err='Inconsistent min_coarse_size')
goto 9999
end if
if (mncszpp /= prec%ag_data%min_coarse_size_per_process) then
info=psb_err_internal_error_
call psb_errpush(info,name,a_err='Inconsistent min_coarse_size_per_process')
goto 9999
end if
if (mxplevs /= prec%ag_data%max_levs) then
info=psb_err_internal_error_
call psb_errpush(info,name,a_err='Inconsistent max_levs')
@ -192,26 +199,23 @@ subroutine amg_c_hierarchy_bld(a,desc_a,prec,info)
! coarse size is hit, or the gain falls below the min_cr_ratio
! threshold.
!
if (casize < 0) then
!
! Default to the cubic root of the size at base level.
!
casize = desc_a%get_global_rows()
casize = int((sone*casize)**(sone/(sone*3)),psb_lpk_)
casize = max(casize,lone)
casize = casize*40_psb_lpk_
call psb_bcast(ctxt,casize)
if (casize > huge(prec%ag_data%min_coarse_size)) then
!
! computed coarse size does not fit in IPK_.
! This is very unlikely, but make sure to put a positive number
!
prec%ag_data%min_coarse_size = huge(prec%ag_data%min_coarse_size)
if ((mncszpp < 0).and.(mncsize<0)) then
mncszpp = 200
prec%ag_data%min_coarse_size_per_process = mncszpp
end if
if (mncszpp > 0) then
casize = mncszpp*np
if (casize > huge(ione)) casize = huge(ione)
else
prec%ag_data%min_coarse_size = casize
if (mncsize < np) then
if (me == 0) write(0,*) &
& 'Warning: resetting coarse size to NP (1 variable per process)'
mncsize = np
end if
casize = mncsize
end if
prec%ag_data%target_coarse_size = casize
nplevs = max(itwo,mxplevs)
!

@ -152,6 +152,9 @@ subroutine amg_ccprecseti(p,what,val,info,ilev,ilmax,pos,idx)
case ('MIN_COARSE_SIZE')
p%ag_data%min_coarse_size = max(val,-1)
return
case ('MIN_COARSE_SIZE_PER_PROCESS')
p%ag_data%min_coarse_size_per_process = max(val,-1)
return
case('MAX_LEVS')
p%ag_data%max_levs = max(val,1)
return

@ -125,6 +125,7 @@ subroutine amg_cprecinit(ctxt,prec,ptype,info)
endif
prec%ctxt = ctxt
prec%ag_data%min_coarse_size = -1
prec%ag_data%min_coarse_size_per_process = -1
select case(psb_toupper(trim(ptype)))
case ('NOPREC','NONE')

@ -82,7 +82,7 @@ subroutine amg_d_hierarchy_bld(a,desc_a,prec,info)
integer(psb_ipk_) :: me,np
integer(psb_ipk_) :: err,i,k, err_act, iszv, newsz,&
& nplevs, mxplevs
integer(psb_lpk_) :: iaggsize, casize
integer(psb_lpk_) :: iaggsize, casize, mncsize, mncszpp
real(psb_dpk_) :: mnaggratio, sizeratio, athresh, aomega
class(amg_d_base_smoother_type), allocatable :: coarse_sm, med_sm, &
& med_sm2, coarse_sm2
@ -132,17 +132,24 @@ subroutine amg_d_hierarchy_bld(a,desc_a,prec,info)
newsz = -1
mxplevs = prec%ag_data%max_levs
mnaggratio = prec%ag_data%min_cr_ratio
casize = prec%ag_data%min_coarse_size
mncsize = prec%ag_data%min_coarse_size
mncszpp = prec%ag_data%min_coarse_size_per_process
iszv = size(prec%precv)
call psb_bcast(ctxt,iszv)
call psb_bcast(ctxt,casize)
call psb_bcast(ctxt,mncsize)
call psb_bcast(ctxt,mncszpp)
call psb_bcast(ctxt,mxplevs)
call psb_bcast(ctxt,mnaggratio)
if (casize /= prec%ag_data%min_coarse_size) then
if (mncsize /= prec%ag_data%min_coarse_size) then
info=psb_err_internal_error_
call psb_errpush(info,name,a_err='Inconsistent min_coarse_size')
goto 9999
end if
if (mncszpp /= prec%ag_data%min_coarse_size_per_process) then
info=psb_err_internal_error_
call psb_errpush(info,name,a_err='Inconsistent min_coarse_size_per_process')
goto 9999
end if
if (mxplevs /= prec%ag_data%max_levs) then
info=psb_err_internal_error_
call psb_errpush(info,name,a_err='Inconsistent max_levs')
@ -192,26 +199,23 @@ subroutine amg_d_hierarchy_bld(a,desc_a,prec,info)
! coarse size is hit, or the gain falls below the min_cr_ratio
! threshold.
!
if (casize < 0) then
!
! Default to the cubic root of the size at base level.
!
casize = desc_a%get_global_rows()
casize = int((done*casize)**(done/(done*3)),psb_lpk_)
casize = max(casize,lone)
casize = casize*40_psb_lpk_
call psb_bcast(ctxt,casize)
if (casize > huge(prec%ag_data%min_coarse_size)) then
!
! computed coarse size does not fit in IPK_.
! This is very unlikely, but make sure to put a positive number
!
prec%ag_data%min_coarse_size = huge(prec%ag_data%min_coarse_size)
if ((mncszpp < 0).and.(mncsize<0)) then
mncszpp = 200
prec%ag_data%min_coarse_size_per_process = mncszpp
end if
if (mncszpp > 0) then
casize = mncszpp*np
if (casize > huge(ione)) casize = huge(ione)
else
prec%ag_data%min_coarse_size = casize
if (mncsize < np) then
if (me == 0) write(0,*) &
& 'Warning: resetting coarse size to NP (1 variable per process)'
mncsize = np
end if
casize = mncsize
end if
prec%ag_data%target_coarse_size = casize
nplevs = max(itwo,mxplevs)
!

@ -158,6 +158,9 @@ subroutine amg_dcprecseti(p,what,val,info,ilev,ilmax,pos,idx)
case ('MIN_COARSE_SIZE')
p%ag_data%min_coarse_size = max(val,-1)
return
case ('MIN_COARSE_SIZE_PER_PROCESS')
p%ag_data%min_coarse_size_per_process = max(val,-1)
return
case('MAX_LEVS')
p%ag_data%max_levs = max(val,1)
return

@ -128,6 +128,7 @@ subroutine amg_dprecinit(ctxt,prec,ptype,info)
endif
prec%ctxt = ctxt
prec%ag_data%min_coarse_size = -1
prec%ag_data%min_coarse_size_per_process = -1
select case(psb_toupper(trim(ptype)))
case ('NOPREC','NONE')

@ -82,7 +82,7 @@ subroutine amg_s_hierarchy_bld(a,desc_a,prec,info)
integer(psb_ipk_) :: me,np
integer(psb_ipk_) :: err,i,k, err_act, iszv, newsz,&
& nplevs, mxplevs
integer(psb_lpk_) :: iaggsize, casize
integer(psb_lpk_) :: iaggsize, casize, mncsize, mncszpp
real(psb_spk_) :: mnaggratio, sizeratio, athresh, aomega
class(amg_s_base_smoother_type), allocatable :: coarse_sm, med_sm, &
& med_sm2, coarse_sm2
@ -132,17 +132,24 @@ subroutine amg_s_hierarchy_bld(a,desc_a,prec,info)
newsz = -1
mxplevs = prec%ag_data%max_levs
mnaggratio = prec%ag_data%min_cr_ratio
casize = prec%ag_data%min_coarse_size
mncsize = prec%ag_data%min_coarse_size
mncszpp = prec%ag_data%min_coarse_size_per_process
iszv = size(prec%precv)
call psb_bcast(ctxt,iszv)
call psb_bcast(ctxt,casize)
call psb_bcast(ctxt,mncsize)
call psb_bcast(ctxt,mncszpp)
call psb_bcast(ctxt,mxplevs)
call psb_bcast(ctxt,mnaggratio)
if (casize /= prec%ag_data%min_coarse_size) then
if (mncsize /= prec%ag_data%min_coarse_size) then
info=psb_err_internal_error_
call psb_errpush(info,name,a_err='Inconsistent min_coarse_size')
goto 9999
end if
if (mncszpp /= prec%ag_data%min_coarse_size_per_process) then
info=psb_err_internal_error_
call psb_errpush(info,name,a_err='Inconsistent min_coarse_size_per_process')
goto 9999
end if
if (mxplevs /= prec%ag_data%max_levs) then
info=psb_err_internal_error_
call psb_errpush(info,name,a_err='Inconsistent max_levs')
@ -192,26 +199,23 @@ subroutine amg_s_hierarchy_bld(a,desc_a,prec,info)
! coarse size is hit, or the gain falls below the min_cr_ratio
! threshold.
!
if (casize < 0) then
!
! Default to the cubic root of the size at base level.
!
casize = desc_a%get_global_rows()
casize = int((sone*casize)**(sone/(sone*3)),psb_lpk_)
casize = max(casize,lone)
casize = casize*40_psb_lpk_
call psb_bcast(ctxt,casize)
if (casize > huge(prec%ag_data%min_coarse_size)) then
!
! computed coarse size does not fit in IPK_.
! This is very unlikely, but make sure to put a positive number
!
prec%ag_data%min_coarse_size = huge(prec%ag_data%min_coarse_size)
if ((mncszpp < 0).and.(mncsize<0)) then
mncszpp = 200
prec%ag_data%min_coarse_size_per_process = mncszpp
end if
if (mncszpp > 0) then
casize = mncszpp*np
if (casize > huge(ione)) casize = huge(ione)
else
prec%ag_data%min_coarse_size = casize
if (mncsize < np) then
if (me == 0) write(0,*) &
& 'Warning: resetting coarse size to NP (1 variable per process)'
mncsize = np
end if
casize = mncsize
end if
prec%ag_data%target_coarse_size = casize
nplevs = max(itwo,mxplevs)
!

@ -152,6 +152,9 @@ subroutine amg_scprecseti(p,what,val,info,ilev,ilmax,pos,idx)
case ('MIN_COARSE_SIZE')
p%ag_data%min_coarse_size = max(val,-1)
return
case ('MIN_COARSE_SIZE_PER_PROCESS')
p%ag_data%min_coarse_size_per_process = max(val,-1)
return
case('MAX_LEVS')
p%ag_data%max_levs = max(val,1)
return

@ -125,6 +125,7 @@ subroutine amg_sprecinit(ctxt,prec,ptype,info)
endif
prec%ctxt = ctxt
prec%ag_data%min_coarse_size = -1
prec%ag_data%min_coarse_size_per_process = -1
select case(psb_toupper(trim(ptype)))
case ('NOPREC','NONE')

@ -82,7 +82,7 @@ subroutine amg_z_hierarchy_bld(a,desc_a,prec,info)
integer(psb_ipk_) :: me,np
integer(psb_ipk_) :: err,i,k, err_act, iszv, newsz,&
& nplevs, mxplevs
integer(psb_lpk_) :: iaggsize, casize
integer(psb_lpk_) :: iaggsize, casize, mncsize, mncszpp
real(psb_dpk_) :: mnaggratio, sizeratio, athresh, aomega
class(amg_z_base_smoother_type), allocatable :: coarse_sm, med_sm, &
& med_sm2, coarse_sm2
@ -132,17 +132,24 @@ subroutine amg_z_hierarchy_bld(a,desc_a,prec,info)
newsz = -1
mxplevs = prec%ag_data%max_levs
mnaggratio = prec%ag_data%min_cr_ratio
casize = prec%ag_data%min_coarse_size
mncsize = prec%ag_data%min_coarse_size
mncszpp = prec%ag_data%min_coarse_size_per_process
iszv = size(prec%precv)
call psb_bcast(ctxt,iszv)
call psb_bcast(ctxt,casize)
call psb_bcast(ctxt,mncsize)
call psb_bcast(ctxt,mncszpp)
call psb_bcast(ctxt,mxplevs)
call psb_bcast(ctxt,mnaggratio)
if (casize /= prec%ag_data%min_coarse_size) then
if (mncsize /= prec%ag_data%min_coarse_size) then
info=psb_err_internal_error_
call psb_errpush(info,name,a_err='Inconsistent min_coarse_size')
goto 9999
end if
if (mncszpp /= prec%ag_data%min_coarse_size_per_process) then
info=psb_err_internal_error_
call psb_errpush(info,name,a_err='Inconsistent min_coarse_size_per_process')
goto 9999
end if
if (mxplevs /= prec%ag_data%max_levs) then
info=psb_err_internal_error_
call psb_errpush(info,name,a_err='Inconsistent max_levs')
@ -192,26 +199,23 @@ subroutine amg_z_hierarchy_bld(a,desc_a,prec,info)
! coarse size is hit, or the gain falls below the min_cr_ratio
! threshold.
!
if (casize < 0) then
!
! Default to the cubic root of the size at base level.
!
casize = desc_a%get_global_rows()
casize = int((done*casize)**(done/(done*3)),psb_lpk_)
casize = max(casize,lone)
casize = casize*40_psb_lpk_
call psb_bcast(ctxt,casize)
if (casize > huge(prec%ag_data%min_coarse_size)) then
!
! computed coarse size does not fit in IPK_.
! This is very unlikely, but make sure to put a positive number
!
prec%ag_data%min_coarse_size = huge(prec%ag_data%min_coarse_size)
if ((mncszpp < 0).and.(mncsize<0)) then
mncszpp = 200
prec%ag_data%min_coarse_size_per_process = mncszpp
end if
if (mncszpp > 0) then
casize = mncszpp*np
if (casize > huge(ione)) casize = huge(ione)
else
prec%ag_data%min_coarse_size = casize
if (mncsize < np) then
if (me == 0) write(0,*) &
& 'Warning: resetting coarse size to NP (1 variable per process)'
mncsize = np
end if
casize = mncsize
end if
prec%ag_data%target_coarse_size = casize
nplevs = max(itwo,mxplevs)
!

@ -158,6 +158,9 @@ subroutine amg_zcprecseti(p,what,val,info,ilev,ilmax,pos,idx)
case ('MIN_COARSE_SIZE')
p%ag_data%min_coarse_size = max(val,-1)
return
case ('MIN_COARSE_SIZE_PER_PROCESS')
p%ag_data%min_coarse_size_per_process = max(val,-1)
return
case('MAX_LEVS')
p%ag_data%max_levs = max(val,1)
return

@ -128,6 +128,7 @@ subroutine amg_zprecinit(ctxt,prec,ptype,info)
endif
prec%ctxt = ctxt
prec%ag_data%min_coarse_size = -1
prec%ag_data%min_coarse_size_per_process = -1
select case(psb_toupper(trim(ptype)))
case ('NOPREC','NONE')

@ -133,7 +133,7 @@ program amg_d_pde2d
real(psb_dpk_), allocatable :: athresv(:) ! smoothed aggregation threshold vector
integer(psb_ipk_) :: thrvsz ! size of threshold vector
real(psb_dpk_) :: athres ! smoothed aggregation threshold
integer(psb_ipk_) :: csize ! minimum size of coarsest matrix
integer(psb_ipk_) :: csizepp ! minimum size of coarsest matrix per process
! AMG smoother or pre-smoother; also 1-lev preconditioner
character(len=16) :: smther ! (pre-)smoother type: BJAC, AS
@ -274,8 +274,8 @@ program amg_d_pde2d
call prec%set('ml_cycle', p_choice%mlcycle, info)
call prec%set('outer_sweeps', p_choice%outer_sweeps,info)
if (p_choice%csize>0)&
& call prec%set('min_coarse_size', p_choice%csize, info)
if (p_choice%csizepp>0)&
& call prec%set('min_coarse_size_per_process', p_choice%csizepp, info)
if (p_choice%mncrratio>1)&
& call prec%set('min_cr_ratio', p_choice%mncrratio, info)
if (p_choice%maxlevs>0)&
@ -551,7 +551,7 @@ contains
call read_data(prec%mlcycle,inp_unit) ! AMG cycle type
call read_data(prec%outer_sweeps,inp_unit) ! number of 1lev/outer sweeps
call read_data(prec%maxlevs,inp_unit) ! max number of levels in AMG prec
call read_data(prec%csize,inp_unit) ! min size coarsest mat
call read_data(prec%csizepp,inp_unit) ! min size coarsest mat
! aggregation
call read_data(prec%aggr_prol,inp_unit) ! aggregation type
call read_data(prec%par_aggr_alg,inp_unit) ! parallel aggregation alg
@ -632,7 +632,7 @@ contains
end if
call psb_bcast(ctxt,prec%athres)
call psb_bcast(ctxt,prec%csize)
call psb_bcast(ctxt,prec%csizepp)
call psb_bcast(ctxt,prec%cmat)
call psb_bcast(ctxt,prec%csolve)
call psb_bcast(ctxt,prec%csbsolve)

@ -134,7 +134,7 @@ program amg_d_pde3d
real(psb_dpk_), allocatable :: athresv(:) ! smoothed aggregation threshold vector
integer(psb_ipk_) :: thrvsz ! size of threshold vector
real(psb_dpk_) :: athres ! smoothed aggregation threshold
integer(psb_ipk_) :: csize ! minimum size of coarsest matrix
integer(psb_ipk_) :: csizepp ! minimum size of coarsest matrix per process
! AMG smoother or pre-smoother; also 1-lev preconditioner
character(len=16) :: smther ! (pre-)smoother type: BJAC, AS
@ -278,8 +278,8 @@ program amg_d_pde3d
call prec%set('ml_cycle', p_choice%mlcycle, info)
call prec%set('outer_sweeps', p_choice%outer_sweeps,info)
if (p_choice%csize>0)&
& call prec%set('min_coarse_size', p_choice%csize, info)
if (p_choice%csizepp>0)&
& call prec%set('min_coarse_size_per_process', p_choice%csizepp, info)
if (p_choice%mncrratio>1)&
& call prec%set('min_cr_ratio', p_choice%mncrratio, info)
if (p_choice%maxlevs>0)&
@ -555,7 +555,7 @@ contains
call read_data(prec%mlcycle,inp_unit) ! AMG cycle type
call read_data(prec%outer_sweeps,inp_unit) ! number of 1lev/outer sweeps
call read_data(prec%maxlevs,inp_unit) ! max number of levels in AMG prec
call read_data(prec%csize,inp_unit) ! min size coarsest mat
call read_data(prec%csizepp,inp_unit) ! min size coarsest mat
! aggregation
call read_data(prec%aggr_prol,inp_unit) ! aggregation type
call read_data(prec%par_aggr_alg,inp_unit) ! parallel aggregation alg
@ -636,7 +636,7 @@ contains
end if
call psb_bcast(ctxt,prec%athres)
call psb_bcast(ctxt,prec%csize)
call psb_bcast(ctxt,prec%csizepp)
call psb_bcast(ctxt,prec%cmat)
call psb_bcast(ctxt,prec%csolve)
call psb_bcast(ctxt,prec%csbsolve)

@ -133,7 +133,7 @@ program amg_s_pde2d
real(psb_spk_), allocatable :: athresv(:) ! smoothed aggregation threshold vector
integer(psb_ipk_) :: thrvsz ! size of threshold vector
real(psb_spk_) :: athres ! smoothed aggregation threshold
integer(psb_ipk_) :: csize ! minimum size of coarsest matrix
integer(psb_ipk_) :: csizepp ! minimum size of coarsest matrix per process
! AMG smoother or pre-smoother; also 1-lev preconditioner
character(len=16) :: smther ! (pre-)smoother type: BJAC, AS
@ -274,8 +274,8 @@ program amg_s_pde2d
call prec%set('ml_cycle', p_choice%mlcycle, info)
call prec%set('outer_sweeps', p_choice%outer_sweeps,info)
if (p_choice%csize>0)&
& call prec%set('min_coarse_size', p_choice%csize, info)
if (p_choice%csizepp>0)&
& call prec%set('min_coarse_size_per_process', p_choice%csizepp, info)
if (p_choice%mncrratio>1)&
& call prec%set('min_cr_ratio', p_choice%mncrratio, info)
if (p_choice%maxlevs>0)&
@ -551,7 +551,7 @@ contains
call read_data(prec%mlcycle,inp_unit) ! AMG cycle type
call read_data(prec%outer_sweeps,inp_unit) ! number of 1lev/outer sweeps
call read_data(prec%maxlevs,inp_unit) ! max number of levels in AMG prec
call read_data(prec%csize,inp_unit) ! min size coarsest mat
call read_data(prec%csizepp,inp_unit) ! min size coarsest mat
! aggregation
call read_data(prec%aggr_prol,inp_unit) ! aggregation type
call read_data(prec%par_aggr_alg,inp_unit) ! parallel aggregation alg
@ -632,7 +632,7 @@ contains
end if
call psb_bcast(ctxt,prec%athres)
call psb_bcast(ctxt,prec%csize)
call psb_bcast(ctxt,prec%csizepp)
call psb_bcast(ctxt,prec%cmat)
call psb_bcast(ctxt,prec%csolve)
call psb_bcast(ctxt,prec%csbsolve)

@ -134,7 +134,7 @@ program amg_s_pde3d
real(psb_spk_), allocatable :: athresv(:) ! smoothed aggregation threshold vector
integer(psb_ipk_) :: thrvsz ! size of threshold vector
real(psb_spk_) :: athres ! smoothed aggregation threshold
integer(psb_ipk_) :: csize ! minimum size of coarsest matrix
integer(psb_ipk_) :: csizepp ! minimum size of coarsest matrix per process
! AMG smoother or pre-smoother; also 1-lev preconditioner
character(len=16) :: smther ! (pre-)smoother type: BJAC, AS
@ -278,8 +278,8 @@ program amg_s_pde3d
call prec%set('ml_cycle', p_choice%mlcycle, info)
call prec%set('outer_sweeps', p_choice%outer_sweeps,info)
if (p_choice%csize>0)&
& call prec%set('min_coarse_size', p_choice%csize, info)
if (p_choice%csizepp>0)&
& call prec%set('min_coarse_size_per_process', p_choice%csizepp, info)
if (p_choice%mncrratio>1)&
& call prec%set('min_cr_ratio', p_choice%mncrratio, info)
if (p_choice%maxlevs>0)&
@ -555,7 +555,7 @@ contains
call read_data(prec%mlcycle,inp_unit) ! AMG cycle type
call read_data(prec%outer_sweeps,inp_unit) ! number of 1lev/outer sweeps
call read_data(prec%maxlevs,inp_unit) ! max number of levels in AMG prec
call read_data(prec%csize,inp_unit) ! min size coarsest mat
call read_data(prec%csizepp,inp_unit) ! min size coarsest mat
! aggregation
call read_data(prec%aggr_prol,inp_unit) ! aggregation type
call read_data(prec%par_aggr_alg,inp_unit) ! parallel aggregation alg
@ -636,7 +636,7 @@ contains
end if
call psb_bcast(ctxt,prec%athres)
call psb_bcast(ctxt,prec%csize)
call psb_bcast(ctxt,prec%csizepp)
call psb_bcast(ctxt,prec%cmat)
call psb_bcast(ctxt,prec%csolve)
call psb_bcast(ctxt,prec%csbsolve)

@ -37,7 +37,7 @@ LLK ! AINV variant, ignored otherwise
VCYCLE ! Type of multilevel CYCLE: VCYCLE WCYCLE KCYCLE MULT ADD
4 ! Number of outer sweeps for ML
-3 ! Max Number of levels in a multilevel preconditioner; if <0, lib default
-3 ! Target coarse matrix size; if <0, lib default
-3 ! Target coarse matrix size per process; if <0, lib default
SMOOTHED ! Type of aggregation: SMOOTHED UNSMOOTHED
DEC ! Parallel aggregation: DEC, SYMDEC
NATURAL ! Ordering of aggregation NATURAL DEGREE

@ -37,7 +37,7 @@ LLK ! AINV variant
VCYCLE ! Type of multilevel CYCLE: VCYCLE WCYCLE KCYCLE MULT ADD
1 ! Number of outer sweeps for ML
-3 ! Max Number of levels in a multilevel preconditioner; if <0, lib default
-3 ! Target coarse matrix size; if <0, lib default
-3 ! Target coarse matrix size per process; if <0, lib default
SMOOTHED ! Type of aggregation: SMOOTHED UNSMOOTHED
DEC ! Parallel aggregation: DEC, SYMDEC
NATURAL ! Ordering of aggregation NATURAL DEGREE

Loading…
Cancel
Save