Integrated INVT and INVK in BJAC prec

implement-ainv
Cirdans-Home 4 years ago
parent 34ffbc3845
commit 0c9098065a

@ -195,7 +195,8 @@ subroutine psb_c_bjac_apply_vect(alpha,prec,x,beta,y,desc_data,info,trans,work)
goto 9999 goto 9999
end if end if
case(psb_f_ainv_) case(psb_f_ainv_,psb_f_invt_,psb_f_invk_)
! Application of approximate inverse preconditioner, just some spmm
select case(trans_) select case(trans_)
case('N') case('N')
@ -373,11 +374,8 @@ subroutine psb_c_bjac_apply(alpha,prec,x,beta,y,desc_data,info,trans,work)
goto 9999 goto 9999
end if end if
case(psb_f_ainv_) case(psb_f_ainv_,psb_f_invt_,psb_f_invk_)
! Application of approximate inverse preconditioner, just some spmm ! Application of approximate inverse preconditioner, just some spmm
! call psb_spmm(alpha, a, x, beta, y,desc_a, info, &
! & trans, work)
select case(trans_) select case(trans_)
@ -866,7 +864,7 @@ subroutine psb_c_bjac_precbld(a,desc_a,prec,info,amold,vmold,imold)
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate') call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
goto 9999 goto 9999
endif endif
! Computin the factorization ! Computing the factorization
call psb_ainv_fact(a,iinvalg,inv_fill,inv_thresh,lf,dd,uf,desc_a,info,iscale=iscale) call psb_ainv_fact(a,iinvalg,inv_fill,inv_thresh,lf,dd,uf,desc_a,info,iscale=iscale)
if(info == psb_success_) then if(info == psb_success_) then
@ -885,6 +883,151 @@ subroutine psb_c_bjac_precbld(a,desc_a,prec,info,amold,vmold,imold)
goto 9999 goto 9999
end if end if
case(psb_f_invk_)
! Approximate Inverse Factorizations based on the sparse inversion of
! triangular factors of an ILU(n) factorization
if (allocated(prec%av)) then
if (size(prec%av) < psb_bp_ilu_avsz) then
do i=1,size(prec%av)
call prec%av(i)%free()
enddo
deallocate(prec%av,stat=info)
endif
end if
if (.not.allocated(prec%av)) then
allocate(prec%av(psb_max_avsz),stat=info)
if (info /= psb_success_) then
call psb_errpush(psb_err_alloc_dealloc_,name)
goto 9999
end if
endif
nrow_a = desc_a%get_local_rows()
nztota = a%get_nzeros()
n_col = desc_a%get_local_cols()
nhalo = n_col-nrow_a
n_row = nrow_a
allocate(lf,uf,stat=info)
if (info == psb_success_) call lf%allocate(n_row,n_row,nztota)
if (info == psb_success_) call uf%allocate(n_row,n_row,nztota)
if(info /= psb_success_) then
info=psb_err_from_subroutine_
ch_err='psb_sp_all'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
allocate(dd(n_row),stat=info)
if (info == psb_success_) then
allocate(prec%dv, stat=info)
if (info == 0) then
if (present(vmold)) then
allocate(prec%dv%v,mold=vmold,stat=info)
else
allocate(psb_c_base_vect_type :: prec%dv%v,stat=info)
end if
end if
end if
if (info /= psb_success_) then
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
goto 9999
endif
! Computing the factorization
call psb_invk_fact(a,fill_in, inv_fill,lf,dd,uf,desc_a,info)
if(info == psb_success_) then
call prec%av(psb_l_pr_)%mv_from(lf%a)
call prec%av(psb_u_pr_)%mv_from(uf%a)
call prec%av(psb_l_pr_)%set_asb()
call prec%av(psb_u_pr_)%set_asb()
call prec%av(psb_l_pr_)%trim()
call prec%av(psb_u_pr_)%trim()
call prec%dv%bld(dd)
! call move_alloc(dd,prec%d)
else
info=psb_err_from_subroutine_
ch_err='psb_ilut_fact'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
case(psb_f_invt_)
! Approximate Inverse Factorizations based on the sparse inversion of
! triangular factors of an ILU(eps) factorization
if (allocated(prec%av)) then
if (size(prec%av) < psb_bp_ilu_avsz) then
do i=1,size(prec%av)
call prec%av(i)%free()
enddo
deallocate(prec%av,stat=info)
endif
end if
if (.not.allocated(prec%av)) then
allocate(prec%av(psb_max_avsz),stat=info)
if (info /= psb_success_) then
call psb_errpush(psb_err_alloc_dealloc_,name)
goto 9999
end if
endif
nrow_a = desc_a%get_local_rows()
nztota = a%get_nzeros()
n_col = desc_a%get_local_cols()
nhalo = n_col-nrow_a
n_row = nrow_a
allocate(lf,uf,stat=info)
if (info == psb_success_) call lf%allocate(n_row,n_row,nztota)
if (info == psb_success_) call uf%allocate(n_row,n_row,nztota)
if(info /= psb_success_) then
info=psb_err_from_subroutine_
ch_err='psb_sp_all'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
allocate(dd(n_row),stat=info)
if (info == psb_success_) then
allocate(prec%dv, stat=info)
if (info == 0) then
if (present(vmold)) then
allocate(prec%dv%v,mold=vmold,stat=info)
else
allocate(psb_c_base_vect_type :: prec%dv%v,stat=info)
end if
end if
end if
if (info /= psb_success_) then
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
goto 9999
endif
! Computing the factorization
call psb_invt_fact(a,fill_in,inv_fill,fact_eps,inv_thresh,lf,dd,uf,&
& desc_a,info)
if(info == psb_success_) then
call prec%av(psb_l_pr_)%mv_from(lf%a)
call prec%av(psb_u_pr_)%mv_from(uf%a)
call prec%av(psb_l_pr_)%set_asb()
call prec%av(psb_u_pr_)%set_asb()
call prec%av(psb_l_pr_)%trim()
call prec%av(psb_u_pr_)%trim()
call prec%dv%bld(dd)
! call move_alloc(dd,prec%d)
else
info=psb_err_from_subroutine_
ch_err='psb_ilut_fact'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
case(psb_f_none_) case(psb_f_none_)
info=psb_err_from_subroutine_ info=psb_err_from_subroutine_
ch_err='Inconsistent prec psb_f_none_' ch_err='Inconsistent prec psb_f_none_'

@ -195,7 +195,8 @@ subroutine psb_d_bjac_apply_vect(alpha,prec,x,beta,y,desc_data,info,trans,work)
goto 9999 goto 9999
end if end if
case(psb_f_ainv_) case(psb_f_ainv_,psb_f_invt_,psb_f_invk_)
! Application of approximate inverse preconditioner, just some spmm
select case(trans_) select case(trans_)
case('N') case('N')
@ -373,11 +374,8 @@ subroutine psb_d_bjac_apply(alpha,prec,x,beta,y,desc_data,info,trans,work)
goto 9999 goto 9999
end if end if
case(psb_f_ainv_) case(psb_f_ainv_,psb_f_invt_,psb_f_invk_)
! Application of approximate inverse preconditioner, just some spmm ! Application of approximate inverse preconditioner, just some spmm
! call psb_spmm(alpha, a, x, beta, y,desc_a, info, &
! & trans, work)
select case(trans_) select case(trans_)
@ -866,7 +864,7 @@ subroutine psb_d_bjac_precbld(a,desc_a,prec,info,amold,vmold,imold)
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate') call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
goto 9999 goto 9999
endif endif
! Computin the factorization ! Computing the factorization
call psb_ainv_fact(a,iinvalg,inv_fill,inv_thresh,lf,dd,uf,desc_a,info,iscale=iscale) call psb_ainv_fact(a,iinvalg,inv_fill,inv_thresh,lf,dd,uf,desc_a,info,iscale=iscale)
if(info == psb_success_) then if(info == psb_success_) then
@ -885,6 +883,151 @@ subroutine psb_d_bjac_precbld(a,desc_a,prec,info,amold,vmold,imold)
goto 9999 goto 9999
end if end if
case(psb_f_invk_)
! Approximate Inverse Factorizations based on the sparse inversion of
! triangular factors of an ILU(n) factorization
if (allocated(prec%av)) then
if (size(prec%av) < psb_bp_ilu_avsz) then
do i=1,size(prec%av)
call prec%av(i)%free()
enddo
deallocate(prec%av,stat=info)
endif
end if
if (.not.allocated(prec%av)) then
allocate(prec%av(psb_max_avsz),stat=info)
if (info /= psb_success_) then
call psb_errpush(psb_err_alloc_dealloc_,name)
goto 9999
end if
endif
nrow_a = desc_a%get_local_rows()
nztota = a%get_nzeros()
n_col = desc_a%get_local_cols()
nhalo = n_col-nrow_a
n_row = nrow_a
allocate(lf,uf,stat=info)
if (info == psb_success_) call lf%allocate(n_row,n_row,nztota)
if (info == psb_success_) call uf%allocate(n_row,n_row,nztota)
if(info /= psb_success_) then
info=psb_err_from_subroutine_
ch_err='psb_sp_all'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
allocate(dd(n_row),stat=info)
if (info == psb_success_) then
allocate(prec%dv, stat=info)
if (info == 0) then
if (present(vmold)) then
allocate(prec%dv%v,mold=vmold,stat=info)
else
allocate(psb_d_base_vect_type :: prec%dv%v,stat=info)
end if
end if
end if
if (info /= psb_success_) then
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
goto 9999
endif
! Computing the factorization
call psb_invk_fact(a,fill_in, inv_fill,lf,dd,uf,desc_a,info)
if(info == psb_success_) then
call prec%av(psb_l_pr_)%mv_from(lf%a)
call prec%av(psb_u_pr_)%mv_from(uf%a)
call prec%av(psb_l_pr_)%set_asb()
call prec%av(psb_u_pr_)%set_asb()
call prec%av(psb_l_pr_)%trim()
call prec%av(psb_u_pr_)%trim()
call prec%dv%bld(dd)
! call move_alloc(dd,prec%d)
else
info=psb_err_from_subroutine_
ch_err='psb_ilut_fact'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
case(psb_f_invt_)
! Approximate Inverse Factorizations based on the sparse inversion of
! triangular factors of an ILU(eps) factorization
if (allocated(prec%av)) then
if (size(prec%av) < psb_bp_ilu_avsz) then
do i=1,size(prec%av)
call prec%av(i)%free()
enddo
deallocate(prec%av,stat=info)
endif
end if
if (.not.allocated(prec%av)) then
allocate(prec%av(psb_max_avsz),stat=info)
if (info /= psb_success_) then
call psb_errpush(psb_err_alloc_dealloc_,name)
goto 9999
end if
endif
nrow_a = desc_a%get_local_rows()
nztota = a%get_nzeros()
n_col = desc_a%get_local_cols()
nhalo = n_col-nrow_a
n_row = nrow_a
allocate(lf,uf,stat=info)
if (info == psb_success_) call lf%allocate(n_row,n_row,nztota)
if (info == psb_success_) call uf%allocate(n_row,n_row,nztota)
if(info /= psb_success_) then
info=psb_err_from_subroutine_
ch_err='psb_sp_all'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
allocate(dd(n_row),stat=info)
if (info == psb_success_) then
allocate(prec%dv, stat=info)
if (info == 0) then
if (present(vmold)) then
allocate(prec%dv%v,mold=vmold,stat=info)
else
allocate(psb_d_base_vect_type :: prec%dv%v,stat=info)
end if
end if
end if
if (info /= psb_success_) then
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
goto 9999
endif
! Computing the factorization
call psb_invt_fact(a,fill_in,inv_fill,fact_eps,inv_thresh,lf,dd,uf,&
& desc_a,info)
if(info == psb_success_) then
call prec%av(psb_l_pr_)%mv_from(lf%a)
call prec%av(psb_u_pr_)%mv_from(uf%a)
call prec%av(psb_l_pr_)%set_asb()
call prec%av(psb_u_pr_)%set_asb()
call prec%av(psb_l_pr_)%trim()
call prec%av(psb_u_pr_)%trim()
call prec%dv%bld(dd)
! call move_alloc(dd,prec%d)
else
info=psb_err_from_subroutine_
ch_err='psb_ilut_fact'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
case(psb_f_none_) case(psb_f_none_)
info=psb_err_from_subroutine_ info=psb_err_from_subroutine_
ch_err='Inconsistent prec psb_f_none_' ch_err='Inconsistent prec psb_f_none_'

@ -195,7 +195,8 @@ subroutine psb_s_bjac_apply_vect(alpha,prec,x,beta,y,desc_data,info,trans,work)
goto 9999 goto 9999
end if end if
case(psb_f_ainv_) case(psb_f_ainv_,psb_f_invt_,psb_f_invk_)
! Application of approximate inverse preconditioner, just some spmm
select case(trans_) select case(trans_)
case('N') case('N')
@ -373,11 +374,8 @@ subroutine psb_s_bjac_apply(alpha,prec,x,beta,y,desc_data,info,trans,work)
goto 9999 goto 9999
end if end if
case(psb_f_ainv_) case(psb_f_ainv_,psb_f_invt_,psb_f_invk_)
! Application of approximate inverse preconditioner, just some spmm ! Application of approximate inverse preconditioner, just some spmm
! call psb_spmm(alpha, a, x, beta, y,desc_a, info, &
! & trans, work)
select case(trans_) select case(trans_)
@ -866,7 +864,7 @@ subroutine psb_s_bjac_precbld(a,desc_a,prec,info,amold,vmold,imold)
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate') call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
goto 9999 goto 9999
endif endif
! Computin the factorization ! Computing the factorization
call psb_ainv_fact(a,iinvalg,inv_fill,inv_thresh,lf,dd,uf,desc_a,info,iscale=iscale) call psb_ainv_fact(a,iinvalg,inv_fill,inv_thresh,lf,dd,uf,desc_a,info,iscale=iscale)
if(info == psb_success_) then if(info == psb_success_) then
@ -885,6 +883,151 @@ subroutine psb_s_bjac_precbld(a,desc_a,prec,info,amold,vmold,imold)
goto 9999 goto 9999
end if end if
case(psb_f_invk_)
! Approximate Inverse Factorizations based on the sparse inversion of
! triangular factors of an ILU(n) factorization
if (allocated(prec%av)) then
if (size(prec%av) < psb_bp_ilu_avsz) then
do i=1,size(prec%av)
call prec%av(i)%free()
enddo
deallocate(prec%av,stat=info)
endif
end if
if (.not.allocated(prec%av)) then
allocate(prec%av(psb_max_avsz),stat=info)
if (info /= psb_success_) then
call psb_errpush(psb_err_alloc_dealloc_,name)
goto 9999
end if
endif
nrow_a = desc_a%get_local_rows()
nztota = a%get_nzeros()
n_col = desc_a%get_local_cols()
nhalo = n_col-nrow_a
n_row = nrow_a
allocate(lf,uf,stat=info)
if (info == psb_success_) call lf%allocate(n_row,n_row,nztota)
if (info == psb_success_) call uf%allocate(n_row,n_row,nztota)
if(info /= psb_success_) then
info=psb_err_from_subroutine_
ch_err='psb_sp_all'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
allocate(dd(n_row),stat=info)
if (info == psb_success_) then
allocate(prec%dv, stat=info)
if (info == 0) then
if (present(vmold)) then
allocate(prec%dv%v,mold=vmold,stat=info)
else
allocate(psb_s_base_vect_type :: prec%dv%v,stat=info)
end if
end if
end if
if (info /= psb_success_) then
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
goto 9999
endif
! Computing the factorization
call psb_invk_fact(a,fill_in, inv_fill,lf,dd,uf,desc_a,info)
if(info == psb_success_) then
call prec%av(psb_l_pr_)%mv_from(lf%a)
call prec%av(psb_u_pr_)%mv_from(uf%a)
call prec%av(psb_l_pr_)%set_asb()
call prec%av(psb_u_pr_)%set_asb()
call prec%av(psb_l_pr_)%trim()
call prec%av(psb_u_pr_)%trim()
call prec%dv%bld(dd)
! call move_alloc(dd,prec%d)
else
info=psb_err_from_subroutine_
ch_err='psb_ilut_fact'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
case(psb_f_invt_)
! Approximate Inverse Factorizations based on the sparse inversion of
! triangular factors of an ILU(eps) factorization
if (allocated(prec%av)) then
if (size(prec%av) < psb_bp_ilu_avsz) then
do i=1,size(prec%av)
call prec%av(i)%free()
enddo
deallocate(prec%av,stat=info)
endif
end if
if (.not.allocated(prec%av)) then
allocate(prec%av(psb_max_avsz),stat=info)
if (info /= psb_success_) then
call psb_errpush(psb_err_alloc_dealloc_,name)
goto 9999
end if
endif
nrow_a = desc_a%get_local_rows()
nztota = a%get_nzeros()
n_col = desc_a%get_local_cols()
nhalo = n_col-nrow_a
n_row = nrow_a
allocate(lf,uf,stat=info)
if (info == psb_success_) call lf%allocate(n_row,n_row,nztota)
if (info == psb_success_) call uf%allocate(n_row,n_row,nztota)
if(info /= psb_success_) then
info=psb_err_from_subroutine_
ch_err='psb_sp_all'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
allocate(dd(n_row),stat=info)
if (info == psb_success_) then
allocate(prec%dv, stat=info)
if (info == 0) then
if (present(vmold)) then
allocate(prec%dv%v,mold=vmold,stat=info)
else
allocate(psb_s_base_vect_type :: prec%dv%v,stat=info)
end if
end if
end if
if (info /= psb_success_) then
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
goto 9999
endif
! Computing the factorization
call psb_invt_fact(a,fill_in,inv_fill,fact_eps,inv_thresh,lf,dd,uf,&
& desc_a,info)
if(info == psb_success_) then
call prec%av(psb_l_pr_)%mv_from(lf%a)
call prec%av(psb_u_pr_)%mv_from(uf%a)
call prec%av(psb_l_pr_)%set_asb()
call prec%av(psb_u_pr_)%set_asb()
call prec%av(psb_l_pr_)%trim()
call prec%av(psb_u_pr_)%trim()
call prec%dv%bld(dd)
! call move_alloc(dd,prec%d)
else
info=psb_err_from_subroutine_
ch_err='psb_ilut_fact'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
case(psb_f_none_) case(psb_f_none_)
info=psb_err_from_subroutine_ info=psb_err_from_subroutine_
ch_err='Inconsistent prec psb_f_none_' ch_err='Inconsistent prec psb_f_none_'

@ -195,7 +195,8 @@ subroutine psb_z_bjac_apply_vect(alpha,prec,x,beta,y,desc_data,info,trans,work)
goto 9999 goto 9999
end if end if
case(psb_f_ainv_) case(psb_f_ainv_,psb_f_invt_,psb_f_invk_)
! Application of approximate inverse preconditioner, just some spmm
select case(trans_) select case(trans_)
case('N') case('N')
@ -373,11 +374,8 @@ subroutine psb_z_bjac_apply(alpha,prec,x,beta,y,desc_data,info,trans,work)
goto 9999 goto 9999
end if end if
case(psb_f_ainv_) case(psb_f_ainv_,psb_f_invt_,psb_f_invk_)
! Application of approximate inverse preconditioner, just some spmm ! Application of approximate inverse preconditioner, just some spmm
! call psb_spmm(alpha, a, x, beta, y,desc_a, info, &
! & trans, work)
select case(trans_) select case(trans_)
@ -866,7 +864,7 @@ subroutine psb_z_bjac_precbld(a,desc_a,prec,info,amold,vmold,imold)
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate') call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
goto 9999 goto 9999
endif endif
! Computin the factorization ! Computing the factorization
call psb_ainv_fact(a,iinvalg,inv_fill,inv_thresh,lf,dd,uf,desc_a,info,iscale=iscale) call psb_ainv_fact(a,iinvalg,inv_fill,inv_thresh,lf,dd,uf,desc_a,info,iscale=iscale)
if(info == psb_success_) then if(info == psb_success_) then
@ -885,6 +883,151 @@ subroutine psb_z_bjac_precbld(a,desc_a,prec,info,amold,vmold,imold)
goto 9999 goto 9999
end if end if
case(psb_f_invk_)
! Approximate Inverse Factorizations based on the sparse inversion of
! triangular factors of an ILU(n) factorization
if (allocated(prec%av)) then
if (size(prec%av) < psb_bp_ilu_avsz) then
do i=1,size(prec%av)
call prec%av(i)%free()
enddo
deallocate(prec%av,stat=info)
endif
end if
if (.not.allocated(prec%av)) then
allocate(prec%av(psb_max_avsz),stat=info)
if (info /= psb_success_) then
call psb_errpush(psb_err_alloc_dealloc_,name)
goto 9999
end if
endif
nrow_a = desc_a%get_local_rows()
nztota = a%get_nzeros()
n_col = desc_a%get_local_cols()
nhalo = n_col-nrow_a
n_row = nrow_a
allocate(lf,uf,stat=info)
if (info == psb_success_) call lf%allocate(n_row,n_row,nztota)
if (info == psb_success_) call uf%allocate(n_row,n_row,nztota)
if(info /= psb_success_) then
info=psb_err_from_subroutine_
ch_err='psb_sp_all'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
allocate(dd(n_row),stat=info)
if (info == psb_success_) then
allocate(prec%dv, stat=info)
if (info == 0) then
if (present(vmold)) then
allocate(prec%dv%v,mold=vmold,stat=info)
else
allocate(psb_z_base_vect_type :: prec%dv%v,stat=info)
end if
end if
end if
if (info /= psb_success_) then
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
goto 9999
endif
! Computing the factorization
call psb_invk_fact(a,fill_in, inv_fill,lf,dd,uf,desc_a,info)
if(info == psb_success_) then
call prec%av(psb_l_pr_)%mv_from(lf%a)
call prec%av(psb_u_pr_)%mv_from(uf%a)
call prec%av(psb_l_pr_)%set_asb()
call prec%av(psb_u_pr_)%set_asb()
call prec%av(psb_l_pr_)%trim()
call prec%av(psb_u_pr_)%trim()
call prec%dv%bld(dd)
! call move_alloc(dd,prec%d)
else
info=psb_err_from_subroutine_
ch_err='psb_ilut_fact'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
case(psb_f_invt_)
! Approximate Inverse Factorizations based on the sparse inversion of
! triangular factors of an ILU(eps) factorization
if (allocated(prec%av)) then
if (size(prec%av) < psb_bp_ilu_avsz) then
do i=1,size(prec%av)
call prec%av(i)%free()
enddo
deallocate(prec%av,stat=info)
endif
end if
if (.not.allocated(prec%av)) then
allocate(prec%av(psb_max_avsz),stat=info)
if (info /= psb_success_) then
call psb_errpush(psb_err_alloc_dealloc_,name)
goto 9999
end if
endif
nrow_a = desc_a%get_local_rows()
nztota = a%get_nzeros()
n_col = desc_a%get_local_cols()
nhalo = n_col-nrow_a
n_row = nrow_a
allocate(lf,uf,stat=info)
if (info == psb_success_) call lf%allocate(n_row,n_row,nztota)
if (info == psb_success_) call uf%allocate(n_row,n_row,nztota)
if(info /= psb_success_) then
info=psb_err_from_subroutine_
ch_err='psb_sp_all'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
allocate(dd(n_row),stat=info)
if (info == psb_success_) then
allocate(prec%dv, stat=info)
if (info == 0) then
if (present(vmold)) then
allocate(prec%dv%v,mold=vmold,stat=info)
else
allocate(psb_z_base_vect_type :: prec%dv%v,stat=info)
end if
end if
end if
if (info /= psb_success_) then
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
goto 9999
endif
! Computing the factorization
call psb_invt_fact(a,fill_in,inv_fill,fact_eps,inv_thresh,lf,dd,uf,&
& desc_a,info)
if(info == psb_success_) then
call prec%av(psb_l_pr_)%mv_from(lf%a)
call prec%av(psb_u_pr_)%mv_from(uf%a)
call prec%av(psb_l_pr_)%set_asb()
call prec%av(psb_u_pr_)%set_asb()
call prec%av(psb_l_pr_)%trim()
call prec%av(psb_u_pr_)%trim()
call prec%dv%bld(dd)
! call move_alloc(dd,prec%d)
else
info=psb_err_from_subroutine_
ch_err='psb_ilut_fact'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
case(psb_f_none_) case(psb_f_none_)
info=psb_err_from_subroutine_ info=psb_err_from_subroutine_
ch_err='Inconsistent prec psb_f_none_' ch_err='Inconsistent prec psb_f_none_'

@ -34,6 +34,8 @@ module psb_c_bjacprec
use psb_c_base_prec_mod use psb_c_base_prec_mod
use psb_c_ilu_fact_mod use psb_c_ilu_fact_mod
use psb_c_ainv_fact_mod use psb_c_ainv_fact_mod
use psb_c_invk_fact_mod
use psb_c_invt_fact_mod
type, extends(psb_c_base_prec_type) :: psb_c_bjac_prec_type type, extends(psb_c_base_prec_type) :: psb_c_bjac_prec_type
integer(psb_ipk_), allocatable :: iprcparm(:) integer(psb_ipk_), allocatable :: iprcparm(:)
@ -59,9 +61,9 @@ module psb_c_bjacprec
end type psb_c_bjac_prec_type end type psb_c_bjac_prec_type
character(len=15), parameter, private :: & character(len=15), parameter, private :: &
& fact_names(0:5)=(/'None ','ILU(0) ',& & fact_names(0:6)=(/'None ','ILU(0) ',&
& 'ILU(n) ','ILU(eps) ','AINV(eps) ',& & 'ILU(n) ','ILU(eps) ','AINV(eps) ',&
& 'INVT '/) & 'INVT ','INVK '/)
private :: psb_c_bjac_sizeof, psb_c_bjac_precdescr, psb_c_bjac_get_nzeros private :: psb_c_bjac_sizeof, psb_c_bjac_precdescr, psb_c_bjac_get_nzeros

@ -34,6 +34,8 @@ module psb_d_bjacprec
use psb_d_base_prec_mod use psb_d_base_prec_mod
use psb_d_ilu_fact_mod use psb_d_ilu_fact_mod
use psb_d_ainv_fact_mod use psb_d_ainv_fact_mod
use psb_d_invk_fact_mod
use psb_d_invt_fact_mod
type, extends(psb_d_base_prec_type) :: psb_d_bjac_prec_type type, extends(psb_d_base_prec_type) :: psb_d_bjac_prec_type
integer(psb_ipk_), allocatable :: iprcparm(:) integer(psb_ipk_), allocatable :: iprcparm(:)
@ -59,9 +61,9 @@ module psb_d_bjacprec
end type psb_d_bjac_prec_type end type psb_d_bjac_prec_type
character(len=15), parameter, private :: & character(len=15), parameter, private :: &
& fact_names(0:5)=(/'None ','ILU(0) ',& & fact_names(0:6)=(/'None ','ILU(0) ',&
& 'ILU(n) ','ILU(eps) ','AINV(eps) ',& & 'ILU(n) ','ILU(eps) ','AINV(eps) ',&
& 'INVT '/) & 'INVT ','INVK '/)
private :: psb_d_bjac_sizeof, psb_d_bjac_precdescr, psb_d_bjac_get_nzeros private :: psb_d_bjac_sizeof, psb_d_bjac_precdescr, psb_d_bjac_get_nzeros

@ -57,8 +57,8 @@ module psb_prec_const_mod
integer(psb_ipk_), parameter :: psb_rfpsz=8 integer(psb_ipk_), parameter :: psb_rfpsz=8
! Factorization types: none, ILU(0), ILU(N), ILU(N,E) ! Factorization types: none, ILU(0), ILU(N), ILU(N,E)
integer(psb_ipk_), parameter :: psb_f_none_=0,psb_f_ilu_n_=1,psb_f_ilu_k_=2,psb_f_ilu_t_=3 integer(psb_ipk_), parameter :: psb_f_none_=0,psb_f_ilu_n_=1,psb_f_ilu_k_=2,psb_f_ilu_t_=3
! Approximate Inverse factorization type: AINV ! Approximate Inverse factorization type: AINV, INVT, INVK
integer(psb_ipk_), parameter :: psb_f_ainv_=4, psb_f_invt_=5 integer(psb_ipk_), parameter :: psb_f_ainv_=4, psb_f_invt_=5, psb_f_invk_=6
! Fields for sparse matrices ensembles: ! Fields for sparse matrices ensembles:
integer(psb_ipk_), parameter :: psb_l_pr_=1, psb_u_pr_=2, psb_bp_ilu_avsz=2 integer(psb_ipk_), parameter :: psb_l_pr_=1, psb_u_pr_=2, psb_bp_ilu_avsz=2
integer(psb_ipk_), parameter :: psb_max_avsz=psb_bp_ilu_avsz integer(psb_ipk_), parameter :: psb_max_avsz=psb_bp_ilu_avsz

@ -34,6 +34,8 @@ module psb_s_bjacprec
use psb_s_base_prec_mod use psb_s_base_prec_mod
use psb_s_ilu_fact_mod use psb_s_ilu_fact_mod
use psb_s_ainv_fact_mod use psb_s_ainv_fact_mod
use psb_s_invk_fact_mod
use psb_s_invt_fact_mod
type, extends(psb_s_base_prec_type) :: psb_s_bjac_prec_type type, extends(psb_s_base_prec_type) :: psb_s_bjac_prec_type
integer(psb_ipk_), allocatable :: iprcparm(:) integer(psb_ipk_), allocatable :: iprcparm(:)
@ -59,9 +61,9 @@ module psb_s_bjacprec
end type psb_s_bjac_prec_type end type psb_s_bjac_prec_type
character(len=15), parameter, private :: & character(len=15), parameter, private :: &
& fact_names(0:5)=(/'None ','ILU(0) ',& & fact_names(0:6)=(/'None ','ILU(0) ',&
& 'ILU(n) ','ILU(eps) ','AINV(eps) ',& & 'ILU(n) ','ILU(eps) ','AINV(eps) ',&
& 'INVT '/) & 'INVT ','INVK '/)
private :: psb_s_bjac_sizeof, psb_s_bjac_precdescr, psb_s_bjac_get_nzeros private :: psb_s_bjac_sizeof, psb_s_bjac_precdescr, psb_s_bjac_get_nzeros

@ -34,6 +34,8 @@ module psb_z_bjacprec
use psb_z_base_prec_mod use psb_z_base_prec_mod
use psb_z_ilu_fact_mod use psb_z_ilu_fact_mod
use psb_z_ainv_fact_mod use psb_z_ainv_fact_mod
use psb_z_invk_fact_mod
use psb_z_invt_fact_mod
type, extends(psb_z_base_prec_type) :: psb_z_bjac_prec_type type, extends(psb_z_base_prec_type) :: psb_z_bjac_prec_type
integer(psb_ipk_), allocatable :: iprcparm(:) integer(psb_ipk_), allocatable :: iprcparm(:)
@ -59,9 +61,9 @@ module psb_z_bjacprec
end type psb_z_bjac_prec_type end type psb_z_bjac_prec_type
character(len=15), parameter, private :: & character(len=15), parameter, private :: &
& fact_names(0:5)=(/'None ','ILU(0) ',& & fact_names(0:6)=(/'None ','ILU(0) ',&
& 'ILU(n) ','ILU(eps) ','AINV(eps) ',& & 'ILU(n) ','ILU(eps) ','AINV(eps) ',&
& 'INVT '/) & 'INVT ','INVK '/)
private :: psb_z_bjac_sizeof, psb_z_bjac_precdescr, psb_z_bjac_get_nzeros private :: psb_z_bjac_sizeof, psb_z_bjac_precdescr, psb_z_bjac_get_nzeros

Loading…
Cancel
Save