stopcriterion
Salvatore Filippone 17 years ago
parent 2c558043f1
commit 3f6ff2a28e

@ -1,7 +1,13 @@
!!$
!!$
!!$ Parallel Sparse BLAS v2.0
!!$ (C) Copyright 2006 Salvatore Filippone University of Rome Tor Vergata
!!$ Alfredo Buttari University of Rome Tor Vergata
!!$ MLD2P4
!!$ MultiLevel Domain Decomposition Parallel Preconditioners Package
!!$ based on PSBLAS (Parallel Sparse BLAS v.2.0)
!!$
!!$ (C) Copyright 2006 Alfredo Buttari University of Rome Tor Vergata
!!$ Pasqua D'Ambra ICAR-CNR, Naples
!!$ Daniela di Serafino Second University of Naples
!!$ Salvatore Filippone University of Rome Tor Vergata
!!$
!!$ Redistribution and use in source and binary forms, with or without
!!$ modification, are permitted provided that the following conditions
@ -47,7 +53,7 @@ Module mld_krylov_mod
interface psb_cg
subroutine psb_dcg(a,prec,b,x,eps,&
& desc_a,info,itmax,iter,err,itrace,istop)
& desc_a,info,itmax,iter,err,itrace,istop)
use psb_base_mod
use psb_prec_mod
type(psb_dspmat_type), intent(in) :: a
@ -65,7 +71,7 @@ Module mld_krylov_mod
interface psb_bicg
subroutine psb_dbicg(a,prec,b,x,eps,&
& desc_a,info,itmax,iter,err,itrace,istop)
& desc_a,info,itmax,iter,err,itrace,istop)
use psb_base_mod
use psb_prec_mod
type(psb_dspmat_type), intent(in) :: a
@ -83,7 +89,7 @@ Module mld_krylov_mod
interface psb_bicgstab
subroutine psb_dcgstab(a,prec,b,x,eps,&
& desc_a,info,itmax,iter,err,itrace,istop)
& desc_a,info,itmax,iter,err,itrace,istop)
use psb_base_mod
use psb_prec_mod
type(psb_dspmat_type), intent(in) :: a
@ -98,7 +104,7 @@ Module mld_krylov_mod
real(kind(1.d0)), optional, intent(out) :: err
end subroutine psb_dcgstab
subroutine psb_zcgstab(a,prec,b,x,eps,&
& desc_a,info,itmax,iter,err,itrace,istop)
& desc_a,info,itmax,iter,err,itrace,istop)
use psb_base_mod
use psb_prec_mod
type(psb_zspmat_type), intent(in) :: a
@ -182,7 +188,7 @@ Module mld_krylov_mod
real(kind(1.d0)), optional, intent(out) :: err
end subroutine psb_dcgs
subroutine psb_zcgs(a,prec,b,x,eps,&
& desc_a,info,itmax,iter,err,itrace,istop)
& desc_a,info,itmax,iter,err,itrace,istop)
use psb_base_mod
use psb_prec_mod
type(psb_zspmat_type), intent(in) :: a

@ -161,6 +161,168 @@ subroutine mld_dprecseti(p,what,val,info,ilev)
endif
end subroutine mld_dprecseti
subroutine mld_dprecsetc(p,what,string,info,ilev)
use psb_base_mod
use mld_prec_mod, mld_protect_name => mld_dprecsetc
implicit none
type(mld_dprec_type), intent(inout) :: p
integer, intent(in) :: what
character(len=*), intent(in) :: string
integer, intent(out) :: info
integer, optional, intent(in) :: ilev
integer :: ilev_, nlev_,val
info = 0
if (.not.allocated(p%baseprecv)) then
write(0,*) 'Error: trying to call PRECSET on an uninitialized preconditioner'
info = -1
return
endif
nlev_ = size(p%baseprecv)
if (present(ilev)) then
ilev_ = ilev
else
ilev_ = 1
end if
if ((ilev_<1).or.(ilev_ > nlev_)) then
write(0,*) 'PRECSET ERRROR: ilev out of bounds'
info = -1
return
endif
if (.not.allocated(p%baseprecv(ilev_)%iprcparm)) then
write(0,*) 'Error: trying to call PRECSET on an uninitialized preconditioner'
info = -1
return
endif
if (present(ilev)) then
if (ilev_ == 1) then
! Rules for fine level are slightly different.
select case(what)
case(mld_prec_type_,mld_sub_solve_,mld_sub_restr_,mld_sub_prol_)
call get_stringval(string,val,info)
p%baseprecv(ilev_)%iprcparm(what) = val
case default
write(0,*) 'Error: trying to call PRECSET with an invalid WHAT'
info = -2
end select
else if (ilev_ > 1) then
select case(what)
case(mld_prec_type_,mld_sub_solve_,mld_sub_restr_,mld_sub_prol_,&
& mld_ml_type_,mld_aggr_alg_,mld_aggr_kind_,mld_coarse_mat_,&
& mld_smooth_pos_,mld_aggr_eig_)
call get_stringval(string,val,info)
p%baseprecv(ilev_)%iprcparm(what) = val
case(mld_coarse_solve_)
call get_stringval(string,val,info)
if (ilev_ /= nlev_) then
write(0,*) 'Inconsistent specification of WHAT vs. ILEV'
info = -2
return
end if
p%baseprecv(ilev_)%iprcparm(mld_sub_solve_) = val
case default
write(0,*) 'Error: trying to call PRECSET with an invalid WHAT'
info = -2
end select
endif
else if (.not.present(ilev)) then
select case(what)
case(mld_prec_type_,mld_sub_solve_,mld_sub_restr_,mld_sub_prol_,mld_sub_ren_,&
& mld_smooth_sweeps_,mld_ml_type_,mld_aggr_alg_,mld_aggr_kind_,mld_coarse_mat_,&
& mld_smooth_pos_)
call get_stringval(string,val,info)
do ilev_=1,nlev_-1
if (.not.allocated(p%baseprecv(ilev_)%iprcparm)) then
write(0,*) 'Error: trying to call PRECSET on an uninitialized preconditioner'
info = -1
return
endif
p%baseprecv(ilev_)%iprcparm(what) = val
end do
case(mld_coarse_solve_)
if (.not.allocated(p%baseprecv(nlev_)%iprcparm)) then
write(0,*) 'Error: trying to call PRECSET on an uninitialized preconditioner'
info = -1
return
endif
call get_stringval(string,val,info)
p%baseprecv(nlev_)%iprcparm(mld_sub_solve_) = val
case default
write(0,*) 'Error: trying to call PRECSET with an invalid WHAT'
info = -2
end select
endif
contains
subroutine get_stringval(string,val,info)
character(len=*), intent(in) :: string
integer, intent(out) :: val, info
info = 0
select case(toupper(trim(string)))
case('NONE')
val = 0
case('HALO')
val = psb_halo_
case('SUM')
val = psb_sum_
case('AVG')
val = psb_avg_
case('ILU')
val = mld_ilu_n_
case('MILU')
val = mld_milu_n_
case('ILUT')
val = mld_ilu_t_
case('SLU')
val = mld_slu_
case('UMFP')
val = mld_umf_
case('ADD')
val = mld_add_ml_
case('MULT')
val = mld_mult_ml_
case('DEC')
val = mld_dec_aggr_
case('REPL')
val = mld_repl_mat_
case('DIST')
val = mld_distr_mat_
case('SYMDEC')
val = mld_sym_dec_aggr_
case('GLB')
val = mld_glb_aggr_
case('SMOOTH')
val = mld_smooth_prol_
case('PRE')
val = mld_pre_smooth_
case('POST')
val = mld_post_smooth_
case('TWOSIDE','BOTH')
val = mld_twoside_smooth_
case default
val = -1
info = -1
end select
if (info /= 0) then
write(0,*) 'Error in get_Stringval: unknown: "',trim(string),'"'
end if
end subroutine get_stringval
end subroutine mld_dprecsetc
subroutine mld_dprecsetd(p,what,val,info,ilev)
use psb_base_mod

@ -101,6 +101,15 @@ module mld_prec_mod
integer, intent(out) :: info
integer, optional, intent(in) :: ilev
end subroutine mld_dprecsetd
subroutine mld_dprecsetc(p,what,string,info,ilev)
use psb_base_mod
use mld_prec_type
type(mld_dprec_type), intent(inout) :: p
integer, intent(in) :: what
character(len=*), intent(in) :: string
integer, intent(out) :: info
integer, optional, intent(in) :: ilev
end subroutine mld_dprecsetc
subroutine mld_zprecseti(p,what,val,info,ilev)
use psb_base_mod
use mld_prec_type

Loading…
Cancel
Save