Make test program to dump preconditioner controlled from input file.

fix-sludist7
Salvatore Filippone 3 years ago
parent af178daa84
commit 92f7cde375

@ -174,6 +174,17 @@ program amg_d_pde2d
real(psb_dpk_) :: cthres ! threshold for ILUT factorization
integer(psb_ipk_) :: cjswp ! sweeps for GS or JAC coarsest-lev subsolver
! Dump data
logical :: dump = .false.
integer(psb_ipk_) :: dlmin ! Minimum level to dump
integer(psb_ipk_) :: dlmax ! Maximum level to dump
logical :: dump_ac = .false.
logical :: dump_rp = .false.
logical :: dump_tprol = .false.
logical :: dump_smoother = .false.
logical :: dump_solver = .false.
logical :: dump_global_num = .false.
end type precdata
type(precdata) :: p_choice
@ -392,6 +403,12 @@ program amg_d_pde2d
write(psb_out_unit,'(" ")')
end if
if (p_choice%dump) then
call prec%dump(info,istart=p_choice%dlmin,iend=p_choice%dlmax,&
& ac=p_choice%dump_ac,rp=p_choice%dump_rp,tprol=p_choice%dump_tprol,&
& smoother=p_choice%dump_smoother, solver=p_choice%dump_solver, &
& global_num=p_choice%dump_global_num)
end if
!
! iterative method parameters
!
@ -530,6 +547,7 @@ contains
! preconditioner type
call read_data(prec%descr,inp_unit) ! verbose description of the prec
call read_data(prec%ptype,inp_unit) ! preconditioner type
! First smoother / 1-lev preconditioner
call read_data(prec%smther,inp_unit) ! smoother type
call read_data(prec%jsweeps,inp_unit) ! (pre-)smoother / 1-lev prec sweeps
@ -580,6 +598,17 @@ contains
call read_data(prec%cfill,inp_unit) ! fill-in for incompl LU
call read_data(prec%cthres,inp_unit) ! Threshold for ILUT
call read_data(prec%cjswp,inp_unit) ! sweeps for GS/JAC subsolver
! dump
call read_data(prec%dump,inp_unit) ! Dump on file?
call read_data(prec%dlmin,inp_unit) ! Minimum level to dump
call read_data(prec%dlmax,inp_unit) ! Maximum level to dump
call read_data(prec%dump_ac,inp_unit)
call read_data(prec%dump_rp,inp_unit)
call read_data(prec%dump_tprol,inp_unit)
call read_data(prec%dump_smoother,inp_unit)
call read_data(prec%dump_solver,inp_unit)
call read_data(prec%dump_global_num,inp_unit)
if (inp_unit /= psb_inp_unit) then
close(inp_unit)
end if
@ -626,6 +655,7 @@ contains
call psb_bcast(ctxt,prec%mlcycle)
call psb_bcast(ctxt,prec%outer_sweeps)
call psb_bcast(ctxt,prec%maxlevs)
call psb_bcast(ctxt,prec%csizepp)
call psb_bcast(ctxt,prec%aggr_prol)
call psb_bcast(ctxt,prec%par_aggr_alg)
@ -641,13 +671,23 @@ contains
end if
call psb_bcast(ctxt,prec%athres)
call psb_bcast(ctxt,prec%csizepp)
call psb_bcast(ctxt,prec%cmat)
call psb_bcast(ctxt,prec%csolve)
call psb_bcast(ctxt,prec%csbsolve)
call psb_bcast(ctxt,prec%cfill)
call psb_bcast(ctxt,prec%cthres)
call psb_bcast(ctxt,prec%cjswp)
! dump
call psb_bcast(ctxt,prec%dump)
call psb_bcast(ctxt,prec%dlmin)
call psb_bcast(ctxt,prec%dlmax)
call psb_bcast(ctxt,prec%dump_ac)
call psb_bcast(ctxt,prec%dump_rp)
call psb_bcast(ctxt,prec%dump_tprol)
call psb_bcast(ctxt,prec%dump_smoother)
call psb_bcast(ctxt,prec%dump_solver)
call psb_bcast(ctxt,prec%dump_global_num)
end subroutine get_parms

@ -175,6 +175,17 @@ program amg_d_pde3d
real(psb_dpk_) :: cthres ! threshold for ILUT factorization
integer(psb_ipk_) :: cjswp ! sweeps for GS or JAC coarsest-lev subsolver
! Dump data
logical :: dump = .false.
integer(psb_ipk_) :: dlmin ! Minimum level to dump
integer(psb_ipk_) :: dlmax ! Maximum level to dump
logical :: dump_ac = .false.
logical :: dump_rp = .false.
logical :: dump_tprol = .false.
logical :: dump_smoother = .false.
logical :: dump_solver = .false.
logical :: dump_global_num = .false.
end type precdata
type(precdata) :: p_choice
@ -396,6 +407,12 @@ program amg_d_pde3d
write(psb_out_unit,'(" ")')
end if
if (p_choice%dump) then
call prec%dump(info,istart=p_choice%dlmin,iend=p_choice%dlmax,&
& ac=p_choice%dump_ac,rp=p_choice%dump_rp,tprol=p_choice%dump_tprol,&
& smoother=p_choice%dump_smoother, solver=p_choice%dump_solver, &
& global_num=p_choice%dump_global_num)
end if
!
! iterative method parameters
!
@ -534,6 +551,7 @@ contains
! preconditioner type
call read_data(prec%descr,inp_unit) ! verbose description of the prec
call read_data(prec%ptype,inp_unit) ! preconditioner type
! First smoother / 1-lev preconditioner
call read_data(prec%smther,inp_unit) ! smoother type
call read_data(prec%jsweeps,inp_unit) ! (pre-)smoother / 1-lev prec sweeps
@ -584,6 +602,17 @@ contains
call read_data(prec%cfill,inp_unit) ! fill-in for incompl LU
call read_data(prec%cthres,inp_unit) ! Threshold for ILUT
call read_data(prec%cjswp,inp_unit) ! sweeps for GS/JAC subsolver
! dump
call read_data(prec%dump,inp_unit) ! Dump on file?
call read_data(prec%dlmin,inp_unit) ! Minimum level to dump
call read_data(prec%dlmax,inp_unit) ! Maximum level to dump
call read_data(prec%dump_ac,inp_unit)
call read_data(prec%dump_rp,inp_unit)
call read_data(prec%dump_tprol,inp_unit)
call read_data(prec%dump_smoother,inp_unit)
call read_data(prec%dump_solver,inp_unit)
call read_data(prec%dump_global_num,inp_unit)
if (inp_unit /= psb_inp_unit) then
close(inp_unit)
end if
@ -652,6 +681,18 @@ contains
call psb_bcast(ctxt,prec%cfill)
call psb_bcast(ctxt,prec%cthres)
call psb_bcast(ctxt,prec%cjswp)
! dump
call psb_bcast(ctxt,prec%dump)
call psb_bcast(ctxt,prec%dlmin)
call psb_bcast(ctxt,prec%dlmax)
call psb_bcast(ctxt,prec%dump_ac)
call psb_bcast(ctxt,prec%dump_rp)
call psb_bcast(ctxt,prec%dump_tprol)
call psb_bcast(ctxt,prec%dump_smoother)
call psb_bcast(ctxt,prec%dump_solver)
call psb_bcast(ctxt,prec%dump_global_num)
end subroutine get_parms

@ -174,6 +174,17 @@ program amg_s_pde2d
real(psb_spk_) :: cthres ! threshold for ILUT factorization
integer(psb_ipk_) :: cjswp ! sweeps for GS or JAC coarsest-lev subsolver
! Dump data
logical :: dump = .false.
integer(psb_ipk_) :: dlmin ! Minimum level to dump
integer(psb_ipk_) :: dlmax ! Maximum level to dump
logical :: dump_ac = .false.
logical :: dump_rp = .false.
logical :: dump_tprol = .false.
logical :: dump_smoother = .false.
logical :: dump_solver = .false.
logical :: dump_global_num = .false.
end type precdata
type(precdata) :: p_choice
@ -392,6 +403,12 @@ program amg_s_pde2d
write(psb_out_unit,'(" ")')
end if
if (p_choice%dump) then
call prec%dump(info,istart=p_choice%dlmin,iend=p_choice%dlmax,&
& ac=p_choice%dump_ac,rp=p_choice%dump_rp,tprol=p_choice%dump_tprol,&
& smoother=p_choice%dump_smoother, solver=p_choice%dump_solver, &
& global_num=p_choice%dump_global_num)
end if
!
! iterative method parameters
!
@ -530,6 +547,7 @@ contains
! preconditioner type
call read_data(prec%descr,inp_unit) ! verbose description of the prec
call read_data(prec%ptype,inp_unit) ! preconditioner type
! First smoother / 1-lev preconditioner
call read_data(prec%smther,inp_unit) ! smoother type
call read_data(prec%jsweeps,inp_unit) ! (pre-)smoother / 1-lev prec sweeps
@ -580,6 +598,17 @@ contains
call read_data(prec%cfill,inp_unit) ! fill-in for incompl LU
call read_data(prec%cthres,inp_unit) ! Threshold for ILUT
call read_data(prec%cjswp,inp_unit) ! sweeps for GS/JAC subsolver
! dump
call read_data(prec%dump,inp_unit) ! Dump on file?
call read_data(prec%dlmin,inp_unit) ! Minimum level to dump
call read_data(prec%dlmax,inp_unit) ! Maximum level to dump
call read_data(prec%dump_ac,inp_unit)
call read_data(prec%dump_rp,inp_unit)
call read_data(prec%dump_tprol,inp_unit)
call read_data(prec%dump_smoother,inp_unit)
call read_data(prec%dump_solver,inp_unit)
call read_data(prec%dump_global_num,inp_unit)
if (inp_unit /= psb_inp_unit) then
close(inp_unit)
end if
@ -626,6 +655,7 @@ contains
call psb_bcast(ctxt,prec%mlcycle)
call psb_bcast(ctxt,prec%outer_sweeps)
call psb_bcast(ctxt,prec%maxlevs)
call psb_bcast(ctxt,prec%csizepp)
call psb_bcast(ctxt,prec%aggr_prol)
call psb_bcast(ctxt,prec%par_aggr_alg)
@ -641,13 +671,23 @@ contains
end if
call psb_bcast(ctxt,prec%athres)
call psb_bcast(ctxt,prec%csizepp)
call psb_bcast(ctxt,prec%cmat)
call psb_bcast(ctxt,prec%csolve)
call psb_bcast(ctxt,prec%csbsolve)
call psb_bcast(ctxt,prec%cfill)
call psb_bcast(ctxt,prec%cthres)
call psb_bcast(ctxt,prec%cjswp)
! dump
call psb_bcast(ctxt,prec%dump)
call psb_bcast(ctxt,prec%dlmin)
call psb_bcast(ctxt,prec%dlmax)
call psb_bcast(ctxt,prec%dump_ac)
call psb_bcast(ctxt,prec%dump_rp)
call psb_bcast(ctxt,prec%dump_tprol)
call psb_bcast(ctxt,prec%dump_smoother)
call psb_bcast(ctxt,prec%dump_solver)
call psb_bcast(ctxt,prec%dump_global_num)
end subroutine get_parms

@ -175,6 +175,17 @@ program amg_s_pde3d
real(psb_spk_) :: cthres ! threshold for ILUT factorization
integer(psb_ipk_) :: cjswp ! sweeps for GS or JAC coarsest-lev subsolver
! Dump data
logical :: dump = .false.
integer(psb_ipk_) :: dlmin ! Minimum level to dump
integer(psb_ipk_) :: dlmax ! Maximum level to dump
logical :: dump_ac = .false.
logical :: dump_rp = .false.
logical :: dump_tprol = .false.
logical :: dump_smoother = .false.
logical :: dump_solver = .false.
logical :: dump_global_num = .false.
end type precdata
type(precdata) :: p_choice
@ -396,6 +407,12 @@ program amg_s_pde3d
write(psb_out_unit,'(" ")')
end if
if (p_choice%dump) then
call prec%dump(info,istart=p_choice%dlmin,iend=p_choice%dlmax,&
& ac=p_choice%dump_ac,rp=p_choice%dump_rp,tprol=p_choice%dump_tprol,&
& smoother=p_choice%dump_smoother, solver=p_choice%dump_solver, &
& global_num=p_choice%dump_global_num)
end if
!
! iterative method parameters
!
@ -534,6 +551,7 @@ contains
! preconditioner type
call read_data(prec%descr,inp_unit) ! verbose description of the prec
call read_data(prec%ptype,inp_unit) ! preconditioner type
! First smoother / 1-lev preconditioner
call read_data(prec%smther,inp_unit) ! smoother type
call read_data(prec%jsweeps,inp_unit) ! (pre-)smoother / 1-lev prec sweeps
@ -584,6 +602,17 @@ contains
call read_data(prec%cfill,inp_unit) ! fill-in for incompl LU
call read_data(prec%cthres,inp_unit) ! Threshold for ILUT
call read_data(prec%cjswp,inp_unit) ! sweeps for GS/JAC subsolver
! dump
call read_data(prec%dump,inp_unit) ! Dump on file?
call read_data(prec%dlmin,inp_unit) ! Minimum level to dump
call read_data(prec%dlmax,inp_unit) ! Maximum level to dump
call read_data(prec%dump_ac,inp_unit)
call read_data(prec%dump_rp,inp_unit)
call read_data(prec%dump_tprol,inp_unit)
call read_data(prec%dump_smoother,inp_unit)
call read_data(prec%dump_solver,inp_unit)
call read_data(prec%dump_global_num,inp_unit)
if (inp_unit /= psb_inp_unit) then
close(inp_unit)
end if
@ -652,6 +681,18 @@ contains
call psb_bcast(ctxt,prec%cfill)
call psb_bcast(ctxt,prec%cthres)
call psb_bcast(ctxt,prec%cjswp)
! dump
call psb_bcast(ctxt,prec%dump)
call psb_bcast(ctxt,prec%dlmin)
call psb_bcast(ctxt,prec%dlmax)
call psb_bcast(ctxt,prec%dump_ac)
call psb_bcast(ctxt,prec%dump_rp)
call psb_bcast(ctxt,prec%dump_tprol)
call psb_bcast(ctxt,prec%dump_smoother)
call psb_bcast(ctxt,prec%dump_solver)
call psb_bcast(ctxt,prec%dump_global_num)
end subroutine get_parms

@ -55,3 +55,13 @@ DIST ! Coarsest-level matrix distribution: DIST REPL
1 ! Coarsest-level fillin P for ILU(P) and ILU(T,P)
1.d-4 ! Coarsest-level threshold T for ILU(T,P)
1 ! Number of sweeps for JACOBI/GS/BJAC coarsest-level solver
%%%%%%%%%%% Dump parms %%%%%%%%%%%%%%%%%%%%%%%%%%
F ! Dump preconditioner on file
1 ! Min level
20 ! Max level
T ! Dump AC
T ! Dump RP
F ! Dump TPROL
F ! Dump SMOOTHER
F ! Dump SOLVER
F ! Global numering ?

@ -41,7 +41,7 @@ VCYCLE ! Type of multilevel CYCLE: VCYCLE WCYCLE KCYCLE MUL
SMOOTHED ! Type of aggregation: SMOOTHED UNSMOOTHED
COUPLED ! Parallel aggregation: DEC, SYMDEC, COUPLED
MATCHBOXP ! aggregation measure SOC1, MATCHBOXP
4 ! Requested size of the aggregates for MATCHBOXP
8 ! Requested size of the aggregates for MATCHBOXP
NATURAL ! Ordering of aggregation NATURAL DEGREE
NOFILTER ! Filtering of matrix: FILTER NOFILTER
-1.5 ! Coarsening ratio, if < 0 use library default
@ -55,3 +55,13 @@ DIST ! Coarsest-level matrix distribution: DIST REPL
1 ! Coarsest-level fillin P for ILU(P) and ILU(T,P)
1.d-4 ! Coarsest-level threshold T for ILU(T,P)
1 ! Number of sweeps for JACOBI/GS/BJAC coarsest-level solver
%%%%%%%%%%% Dump parms %%%%%%%%%%%%%%%%%%%%%%%%%%
F ! Dump preconditioner on file
1 ! Min level
20 ! Max level
T ! Dump AC
T ! Dump RP
F ! Dump TPROL
F ! Dump SMOOTHER
F ! Dump SOLVER
F ! Global numering ?

Loading…
Cancel
Save