mld2p4-2:

mlprec/impl/mld_c_dec_map_bld.f90
 mlprec/impl/mld_d_dec_map_bld.f90
 mlprec/impl/mld_s_dec_map_bld.f90
 mlprec/impl/mld_z_dec_map_bld.f90
 mlprec/impl/solver/mld_d_mumps_solver_apply.F90
 mlprec/impl/solver/mld_d_mumps_solver_bld.F90
 mlprec/mld_base_prec_type.F90

Fixed bug in decoupled aggregation  (how did it survive so long???)
Mumps integration
stopcriterion
Salvatore Filippone 9 years ago
parent 260a6aa3ce
commit 3d6e5fab87

@ -54,7 +54,7 @@ subroutine mld_c_dec_map_bld(theta,a,desc_a,nlaggr,ilaggr,info)
! Local variables ! Local variables
integer(psb_ipk_), allocatable :: ils(:), neigh(:), irow(:), icol(:) integer(psb_ipk_), allocatable :: ils(:), neigh(:), irow(:), icol(:)
complex(psb_spk_), allocatable :: val(:), diag(:) complex(psb_spk_), allocatable :: val(:), diag(:)
integer(psb_ipk_) :: icnt,nlp,k,n,ia,isz,nr, naggr,i,j,m, nz integer(psb_ipk_) :: icnt,nlp,k,n,ia,isz,nr, naggr,i,j,m, nz, ilg
real(psb_spk_) :: cpling, tcl real(psb_spk_) :: cpling, tcl
logical :: recovery logical :: recovery
integer(psb_ipk_) :: debug_level, debug_unit,err_act integer(psb_ipk_) :: debug_level, debug_unit,err_act
@ -125,7 +125,8 @@ subroutine mld_c_dec_map_bld(theta,a,desc_a,nlaggr,ilaggr,info)
do k=1, nz do k=1, nz
j = icol(k) j = icol(k)
if ((1<=j).and.(j<=nr).and.(i /= j)) then ilg = ilaggr(j)
if ((ilg<0).and.(1<=j).and.(j<=nr).and.(i /= j)) then
if (abs(val(k)) > theta*sqrt(abs(diag(i)*diag(j)))) then if (abs(val(k)) > theta*sqrt(abs(diag(i)*diag(j)))) then
ilaggr(j) = naggr ilaggr(j) = naggr
else else

@ -54,7 +54,7 @@ subroutine mld_d_dec_map_bld(theta,a,desc_a,nlaggr,ilaggr,info)
! Local variables ! Local variables
integer(psb_ipk_), allocatable :: ils(:), neigh(:), irow(:), icol(:) integer(psb_ipk_), allocatable :: ils(:), neigh(:), irow(:), icol(:)
real(psb_dpk_), allocatable :: val(:), diag(:) real(psb_dpk_), allocatable :: val(:), diag(:)
integer(psb_ipk_) :: icnt,nlp,k,n,ia,isz,nr, naggr,i,j,m, nz integer(psb_ipk_) :: icnt,nlp,k,n,ia,isz,nr, naggr,i,j,m, nz, ilg
real(psb_dpk_) :: cpling, tcl real(psb_dpk_) :: cpling, tcl
logical :: recovery logical :: recovery
integer(psb_ipk_) :: debug_level, debug_unit,err_act integer(psb_ipk_) :: debug_level, debug_unit,err_act
@ -125,7 +125,8 @@ subroutine mld_d_dec_map_bld(theta,a,desc_a,nlaggr,ilaggr,info)
do k=1, nz do k=1, nz
j = icol(k) j = icol(k)
if ((1<=j).and.(j<=nr).and.(i /= j)) then ilg = ilaggr(j)
if ((ilg<0).and.(1<=j).and.(j<=nr).and.(i /= j)) then
if (abs(val(k)) > theta*sqrt(abs(diag(i)*diag(j)))) then if (abs(val(k)) > theta*sqrt(abs(diag(i)*diag(j)))) then
ilaggr(j) = naggr ilaggr(j) = naggr
else else

@ -54,7 +54,7 @@ subroutine mld_s_dec_map_bld(theta,a,desc_a,nlaggr,ilaggr,info)
! Local variables ! Local variables
integer(psb_ipk_), allocatable :: ils(:), neigh(:), irow(:), icol(:) integer(psb_ipk_), allocatable :: ils(:), neigh(:), irow(:), icol(:)
real(psb_spk_), allocatable :: val(:), diag(:) real(psb_spk_), allocatable :: val(:), diag(:)
integer(psb_ipk_) :: icnt,nlp,k,n,ia,isz,nr, naggr,i,j,m, nz integer(psb_ipk_) :: icnt,nlp,k,n,ia,isz,nr, naggr,i,j,m, nz, ilg
real(psb_spk_) :: cpling, tcl real(psb_spk_) :: cpling, tcl
logical :: recovery logical :: recovery
integer(psb_ipk_) :: debug_level, debug_unit,err_act integer(psb_ipk_) :: debug_level, debug_unit,err_act
@ -125,7 +125,8 @@ subroutine mld_s_dec_map_bld(theta,a,desc_a,nlaggr,ilaggr,info)
do k=1, nz do k=1, nz
j = icol(k) j = icol(k)
if ((1<=j).and.(j<=nr).and.(i /= j)) then ilg = ilaggr(j)
if ((ilg<0).and.(1<=j).and.(j<=nr).and.(i /= j)) then
if (abs(val(k)) > theta*sqrt(abs(diag(i)*diag(j)))) then if (abs(val(k)) > theta*sqrt(abs(diag(i)*diag(j)))) then
ilaggr(j) = naggr ilaggr(j) = naggr
else else

@ -54,7 +54,7 @@ subroutine mld_z_dec_map_bld(theta,a,desc_a,nlaggr,ilaggr,info)
! Local variables ! Local variables
integer(psb_ipk_), allocatable :: ils(:), neigh(:), irow(:), icol(:) integer(psb_ipk_), allocatable :: ils(:), neigh(:), irow(:), icol(:)
complex(psb_dpk_), allocatable :: val(:), diag(:) complex(psb_dpk_), allocatable :: val(:), diag(:)
integer(psb_ipk_) :: icnt,nlp,k,n,ia,isz,nr, naggr,i,j,m, nz integer(psb_ipk_) :: icnt,nlp,k,n,ia,isz,nr, naggr,i,j,m, nz, ilg
real(psb_dpk_) :: cpling, tcl real(psb_dpk_) :: cpling, tcl
logical :: recovery logical :: recovery
integer(psb_ipk_) :: debug_level, debug_unit,err_act integer(psb_ipk_) :: debug_level, debug_unit,err_act
@ -125,7 +125,8 @@ subroutine mld_z_dec_map_bld(theta,a,desc_a,nlaggr,ilaggr,info)
do k=1, nz do k=1, nz
j = icol(k) j = icol(k)
if ((1<=j).and.(j<=nr).and.(i /= j)) then ilg = ilaggr(j)
if ((ilg<0).and.(1<=j).and.(j<=nr).and.(i /= j)) then
if (abs(val(k)) > theta*sqrt(abs(diag(i)*diag(j)))) then if (abs(val(k)) > theta*sqrt(abs(diag(i)*diag(j)))) then
ilaggr(j) = naggr ilaggr(j) = naggr
else else

@ -108,10 +108,10 @@ subroutine d_mumps_solver_apply(alpha,sv,x,beta,y,desc_data,trans,work,info)
sv%id%rhs => gx sv%id%rhs => gx
sv%id%nrhs = 1 sv%id%nrhs = 1
sv%id%icntl(1)=-1 sv%id%icntl(1) = -1
sv%id%icntl(2)=-1 sv%id%icntl(2) = -1
sv%id%icntl(3)=-1 sv%id%icntl(3) = -1
sv%id%icntl(4)=-1 sv%id%icntl(4) = -1
sv%id%job = 3 sv%id%job = 3
call dmumps(sv%id) call dmumps(sv%id)
call psb_scatter(gx, ww, desc_data, info, root=0) call psb_scatter(gx, ww, desc_data, info, root=0)

@ -107,7 +107,7 @@
sv%id%comm = icomm sv%id%comm = icomm
sv%id%job = -1 sv%id%job = -1
sv%id%par=1 sv%id%par = 1
call dmumps(sv%id) call dmumps(sv%id)
!WARNING: CALLING dMUMPS WITH JOB=-1 DESTROY THE SETTING OF DEFAULT:TO FIX !WARNING: CALLING dMUMPS WITH JOB=-1 DESTROY THE SETTING OF DEFAULT:TO FIX
sv%id%icntl(3)=sv%ipar(2) sv%id%icntl(3)=sv%ipar(2)
@ -131,9 +131,9 @@
call psb_loc_to_glob(acoo%ia(1:nztota), desc_a, info, iact='I') call psb_loc_to_glob(acoo%ia(1:nztota), desc_a, info, iact='I')
end if end if
sv%id%irn_loc=> acoo%ia sv%id%irn_loc => acoo%ia
sv%id%jcn_loc=> acoo%ja sv%id%jcn_loc => acoo%ja
sv%id%a_loc=> acoo%val sv%id%a_loc => acoo%val
sv%id%icntl(18)=3 sv%id%icntl(18)=3
if(acoo%is_upper() .or. acoo%is_lower()) then if(acoo%is_upper() .or. acoo%is_lower()) then
sv%id%sym = 2 sv%id%sym = 2

@ -331,7 +331,7 @@ module mld_base_prec_type
& mld_fact_names(0:mld_max_sub_solve_)=(/& & mld_fact_names(0:mld_max_sub_solve_)=(/&
& 'none ','none ',& & 'none ','none ',&
& 'none ','none ',& & 'none ','none ',&
& 'none ', 'Point Jacobi ',& & 'none ','Point Jacobi ',&
& 'Gauss-Seidel ','ILU(n) ',& & 'Gauss-Seidel ','ILU(n) ',&
& 'MILU(n) ','ILU(t,n) ',& & 'MILU(n) ','ILU(t,n) ',&
& 'SuperLU ','UMFPACK LU ',& & 'SuperLU ','UMFPACK LU ',&

Loading…
Cancel
Save