Fixed unspecified interface in prcaply1.

psblas3-type-indexed
Salvatore Filippone 19 years ago
parent ae9c37efae
commit a3151bb4de

@ -198,7 +198,8 @@ Subroutine psb_dgmresr(a,prec,b,x,eps,desc_a,info,&
call psb_errpush(info,name) call psb_errpush(info,name)
goto 9999 goto 9999
End If End If
if (debug) write(0,*) 'Size of V,W ',size(v),size(v,1),&
&size(w),size(w,1), size(v(:,1))
! ensure global coherence for convergence checks. ! ensure global coherence for convergence checks.
Call blacs_get(icontxt,16,isvch) Call blacs_get(icontxt,16,isvch)
ich = 1 ich = 1
@ -226,6 +227,11 @@ Subroutine psb_dgmresr(a,prec,b,x,eps,desc_a,info,&
If (debug) Write(0,*) 'restart: ',itx,it If (debug) Write(0,*) 'restart: ',itx,it
it = 0 it = 0
Call psb_axpby(one,b,zero,v(:,1),desc_a,info) Call psb_axpby(one,b,zero,v(:,1),desc_a,info)
if (info.ne.0) Then
info=4011
call psb_errpush(info,name)
goto 9999
End If
Call psb_spmm(-one,a,x,one,v(:,1),desc_a,info,work=aux) Call psb_spmm(-one,a,x,one,v(:,1),desc_a,info,work=aux)
call psb_prcaply(prec,v(:,1),desc_a,info) call psb_prcaply(prec,v(:,1),desc_a,info)

@ -305,10 +305,10 @@ contains
write (0,'("input argument n. ",i0," has an invalid value")')i_e_d(1) write (0,'("input argument n. ",i0," has an invalid value")')i_e_d(1)
write (0,'("current value is ",a)')a_e_d(2:2) write (0,'("current value is ",a)')a_e_d(2:2)
case(50) case(50)
write (0,'("input argument n. ",i0," must be equal or greater than input argument n. ",i0)') i_e_d(1), i_e_d(2) write (0,'("input argument n. ",i0," must be equal or greater than input argument n. ",i0)') i_e_d(1), i_e_d(3)
write (0,'("current values are ",i0," < ",i0)') i_e_d(3),i_e_d(4) write (0,'("current values are ",i0," < ",i0)') i_e_d(2),i_e_d(5)
case(60) case(60)
write (0,'("input argument n. ",i0," must be equal or greater than ",i0)')i_e_d(1),i_e_d(2) write (0,'("input argument n. ",i0," must be greater than or equal to ",i0)')i_e_d(1),i_e_d(2)
write (0,'("current value is ",i0," < ",i0)')i_e_d(3), i_e_d(2) write (0,'("current value is ",i0," < ",i0)')i_e_d(3), i_e_d(2)
case(70) case(70)
write (0,'("input argument n. ",i0," in entry # ",i0," has an invalid value")')i_e_d(1:2) write (0,'("input argument n. ",i0," in entry # ",i0," has an invalid value")')i_e_d(1:2)

@ -115,6 +115,7 @@ subroutine psb_dprecaply(prec,x,y,desc_data,info,trans, work)
write(0,*) 'Inconsistent preconditioner: neither SMTH nor BASE?' write(0,*) 'Inconsistent preconditioner: neither SMTH nor BASE?'
end if end if
if (size(prec%baseprecv) >1) then if (size(prec%baseprecv) >1) then
if (debug) write(0,*) 'Into mlprcaply',size(x),size(y)
call psb_dmlprcaply(prec%baseprecv,x,zero,y,desc_data,trans_,work_,info) call psb_dmlprcaply(prec%baseprecv,x,zero,y,desc_data,trans_,work_,info)
if(info /= 0) then if(info /= 0) then
call psb_errpush(4010,name,a_err='psb_dmlprcaply') call psb_errpush(4010,name,a_err='psb_dmlprcaply')
@ -335,7 +336,7 @@ subroutine psb_dbaseprcaply(prec,x,beta,y,desc_data,trans,work,info)
endif endif
if (debug) write(0,*)' vdiag: ',prec%d(:) if (debugprt) write(0,*)' vdiag: ',prec%d(:)
if (debug) write(0,*) 'Bi-CGSTAB with Additive Schwarz prec' if (debug) write(0,*) 'Bi-CGSTAB with Additive Schwarz prec'
tx(1:desc_data%matrix_data(psb_n_row_)) = x(1:desc_data%matrix_data(psb_n_row_)) tx(1:desc_data%matrix_data(psb_n_row_)) = x(1:desc_data%matrix_data(psb_n_row_))
@ -631,7 +632,7 @@ subroutine psb_dbjacaply(prec,x,beta,y,desc_data,trans,work,info)
case default case default
write(0,*) 'Unknown factorization type in bjac_prcaply',prec%iprcparm(f_type_) write(0,*) 'Unknown factorization type in bjac_prcaply',prec%iprcparm(f_type_)
end select end select
if (debug) write(0,*)' Y: ',y(:) if (debugprt) write(0,*)' Y: ',y(:)
else if (prec%iprcparm(jac_sweeps_) > 1) then else if (prec%iprcparm(jac_sweeps_) > 1) then
@ -933,7 +934,7 @@ subroutine psb_dmlprcaply(baseprecv,x,beta,y,desc_data,trans,work,info)
end if end if
if (debug) write(0,*)' Y2: ',Y(:) if (debugprt) write(0,*)' Y2: ',Y(:)
deallocate(t2l,w2l) deallocate(t2l,w2l)
@ -968,10 +969,13 @@ subroutine psb_dmlprcaply(baseprecv,x,beta,y,desc_data,trans,work,info)
! !
if (debug) write(0,*)' mult_ml_apply omega ',omega if (debug) write(0,*)' mult_ml_apply omega ',omega
if (debug) write(0,*)' mult_ml_apply X: ',X(:) if (debugprt) write(0,*)' mult_ml_apply X: ',X(:)
call psb_axpby(one,x,zero,tx,desc_data,info) call psb_axpby(one,x,zero,tx,desc_data,info)
if(info /=0) goto 9999 if(info /=0) then
if (debug) write(0,*)' From axpby1 ',size(x),size(tx),n_row,n_col,nr2l,nrg
call psb_errpush(4010,name,a_err='axpby post_smooth 1')
goto 9999
endif
if (ismth /= no_smth_) then if (ismth /= no_smth_) then
! !
! Smoothed aggregation ! Smoothed aggregation
@ -1315,6 +1319,21 @@ subroutine psb_dprecaply1(prec,x,desc_data,info,trans)
logical,parameter :: debug=.false., debugprt=.false. logical,parameter :: debug=.false., debugprt=.false.
real(kind(1.d0)), parameter :: one=1.d0, zero=0.d0 real(kind(1.d0)), parameter :: one=1.d0, zero=0.d0
interface
subroutine psb_dprecaply(prec,x,y,desc_data,info,trans, work)
use psb_descriptor_type
use psb_prec_type
implicit none
type(psb_desc_type),intent(in) :: desc_data
type(psb_dprec_type), intent(in) :: prec
real(kind(0.d0)),intent(inout) :: x(:), y(:)
integer, intent(out) :: info
character(len=1), optional :: trans
real(kind(0.d0)), optional, target :: work(:)
end subroutine psb_dprecaply
end interface
! Local variables ! Local variables
character :: trans_ character :: trans_
@ -1339,8 +1358,8 @@ subroutine psb_dprecaply1(prec,x,desc_data,info,trans)
call psb_errpush(4010,name,a_err='Allocate') call psb_errpush(4010,name,a_err='Allocate')
goto 9999 goto 9999
end if end if
if (debug) write(0,*) 'Prcaply1 Size(x) ',size(x), size(ww),size(w1)
call psb_dprecaply(prec,x,ww,desc_data,info,trans_,w1) call psb_dprecaply(prec,x,ww,desc_data,info,trans_,work=w1)
if(info /=0) goto 9999 if(info /=0) goto 9999
x(:) = ww(:) x(:) = ww(:)
deallocate(ww,W1) deallocate(ww,W1)

@ -925,7 +925,8 @@ subroutine psb_mlprec_bld(a,desc_a,p,info)
call psb_ipcoo2csr(p%av(ac_),info) call psb_ipcoo2csr(p%av(ac_),info)
if(info /= 0) then if(info /= 0) then
info=4011 info=4011
call psb_errpush(info,name) ch_err='psb_ipcoo2csr'
call psb_errpush(info,name,a_err=ch_err)
goto 9999 goto 9999
end if end if
@ -942,7 +943,8 @@ subroutine psb_mlprec_bld(a,desc_a,p,info)
call psb_splu(p%av(ac_),p%av(l_pr_),p%av(u_pr_),p%d,info) call psb_splu(p%av(ac_),p%av(l_pr_),p%av(u_pr_),p%d,info)
if(info /= 0) then if(info /= 0) then
info=4011 info=4011
call psb_errpush(info,name) ch_err='psb_splu'
call psb_errpush(info,name,a_err=ch_err)
goto 9999 goto 9999
end if end if
@ -952,7 +954,8 @@ subroutine psb_mlprec_bld(a,desc_a,p,info)
call psb_ipcsr2coo(p%av(ac_),info) call psb_ipcsr2coo(p%av(ac_),info)
if(info /= 0) then if(info /= 0) then
info=4011 info=4011
call psb_errpush(info,name) ch_err='psb_ipcsr2coo'
call psb_errpush(info,name,a_err=ch_err)
goto 9999 goto 9999
end if end if
k=0 k=0
@ -971,7 +974,8 @@ subroutine psb_mlprec_bld(a,desc_a,p,info)
& p%av(ac_)%aspk,p%av(ac_)%ia2,p%av(ac_)%ia1,p%iprcparm(slu_ptr_),info) & p%av(ac_)%aspk,p%av(ac_)%ia2,p%av(ac_)%ia1,p%iprcparm(slu_ptr_),info)
if(info /= 0) then if(info /= 0) then
info=4011 info=4011
call psb_errpush(info,name) ch_err='psb_fort_slu_factor'
call psb_errpush(info,name,a_err=ch_err)
goto 9999 goto 9999
end if end if
@ -981,7 +985,8 @@ subroutine psb_mlprec_bld(a,desc_a,p,info)
call psb_ipcsr2coo(p%av(ac_),info) call psb_ipcsr2coo(p%av(ac_),info)
if(info /= 0) then if(info /= 0) then
info=4011 info=4011
call psb_errpush(info,name) ch_err='psb_ipcsr2coo'
call psb_errpush(info,name,a_err=ch_err)
goto 9999 goto 9999
end if end if
k=0 k=0
@ -1001,7 +1006,8 @@ subroutine psb_mlprec_bld(a,desc_a,p,info)
& p%iprcparm(umf_symptr_),p%iprcparm(umf_numptr_),info) & p%iprcparm(umf_symptr_),p%iprcparm(umf_numptr_),info)
if(info /= 0) then if(info /= 0) then
info=4011 info=4011
call psb_errpush(info,name) ch_err='psb_fort_umf_factor'
call psb_errpush(info,name,a_err=ch_err)
goto 9999 goto 9999
end if end if

@ -99,28 +99,28 @@ subroutine psb_chkglobvect( m, n, lldx, ix, jx, desc_dec, info)
int_err(3) = 6 int_err(3) = 6
int_err(4) = psb_n_col_ int_err(4) = psb_n_col_
int_err(5) = desc_dec(psb_n_col_) int_err(5) = desc_dec(psb_n_col_)
else if (desc_dec(n_).lt.m) then else if (desc_dec(psb_n_).lt.m) then
info=60 info=60
int_err(1) = 1 int_err(1) = 1
int_err(2) = m int_err(2) = m
int_err(3) = 6 int_err(3) = 6
int_err(4) = n_ int_err(4) = psb_n_
int_err(5) = desc_dec(n_) int_err(5) = desc_dec(psb_n_)
else if (desc_dec(n_).lt.ix) then else if (desc_dec(psb_n_).lt.ix) then
info=60 info=60
int_err(1) = 4 int_err(1) = 4
int_err(2) = ix int_err(2) = ix
int_err(3) = 6 int_err(3) = 6
int_err(4) = n_ int_err(4) = psb_n_
int_err(5) = desc_dec(n_) int_err(5) = desc_dec(psb_n_)
else if (desc_dec(m_).lt.jx) then else if (desc_dec(psb_m_).lt.jx) then
info=60 info=60
int_err(1) = 5 int_err(1) = 5
int_err(2) = jx int_err(2) = jx
int_err(3) = 6 int_err(3) = 6
int_err(4) = m_ int_err(4) = psb_m_
int_err(5) = desc_dec(m_) int_err(5) = desc_dec(psb_m_)
else if (desc_dec(n_).lt.(ix+m-1)) then else if (desc_dec(psb_n_).lt.(ix+m-1)) then
info=80 info=80
int_err(1) = 1 int_err(1) = 1
int_err(2) = m int_err(2) = m

@ -96,41 +96,41 @@ subroutine psb_chkmat( m, n, ia, ja, desc_dec, info, iia, jja)
int_err(1) = 6 int_err(1) = 6
int_err(2) = psb_n_row_ int_err(2) = psb_n_row_
int_err(3) = desc_dec(psb_n_row_) int_err(3) = desc_dec(psb_n_row_)
else if (desc_dec(m_).lt.m) then else if (desc_dec(psb_m_).lt.m) then
info=60 info=60
int_err(1) = 1 int_err(1) = 1
int_err(2) = m int_err(2) = m
int_err(3) = 5 int_err(3) = 5
int_err(4) = m_ int_err(4) = psb_m_
int_err(5) = desc_dec(m_) int_err(5) = desc_dec(psb_m_)
else if (desc_dec(m_).lt.m) then else if (desc_dec(psb_n_).lt.m) then
info=60 info=60
int_err(1) = 2 int_err(1) = 2
int_err(2) = n int_err(2) = n
int_err(3) = 5 int_err(3) = 5
int_err(4) = n_ int_err(4) = psb_n_
int_err(5) = desc_dec(n_) int_err(5) = desc_dec(psb_n_)
else if (desc_dec(m_).lt.ia) then else if (desc_dec(psb_m_).lt.ia) then
info=60 info=60
int_err(1) = 3 int_err(1) = 3
int_err(2) = ia int_err(2) = ia
int_err(3) = 5 int_err(3) = 5
int_err(4) = m_ int_err(4) = psb_m_
int_err(5) = desc_dec(m_) int_err(5) = desc_dec(psb_m_)
else if (desc_dec(n_).lt.ja) then else if (desc_dec(psb_n_).lt.ja) then
info=60 info=60
int_err(1) = 4 int_err(1) = 4
int_err(2) = ja int_err(2) = ja
int_err(3) = 5 int_err(3) = 5
int_err(4) = n_ int_err(4) = psb_n_
int_err(5) = desc_dec(n_) int_err(5) = desc_dec(psb_n_)
else if (desc_dec(m_).lt.(ia+m-1)) then else if (desc_dec(psb_m_).lt.(ia+m-1)) then
info=80 info=80
int_err(1) = 1 int_err(1) = 1
int_err(2) = m int_err(2) = m
int_err(3) = 3 int_err(3) = 3
int_err(4) = ia int_err(4) = ia
else if (desc_dec(n_).lt.(ja+n-1)) then else if (desc_dec(psb_n_).lt.(ja+n-1)) then
info=80 info=80
int_err(1) = 2 int_err(1) = 2
int_err(2) = n int_err(2) = n

@ -103,28 +103,28 @@ subroutine psb_chkvect( m, n, lldx, ix, jx, desc_dec, info, iix, jjx)
int_err(3) = 6 int_err(3) = 6
int_err(4) = psb_n_col_ int_err(4) = psb_n_col_
int_err(5) = desc_dec(psb_n_col_) int_err(5) = desc_dec(psb_n_col_)
else if (desc_dec(n_).lt.m) then else if (desc_dec(psb_n_).lt.m) then
info=60 info=60
int_err(1) = 1 int_err(1) = 1
int_err(2) = m int_err(2) = m
int_err(3) = 6 int_err(3) = 6
int_err(4) = n_ int_err(4) = psb_n_
int_err(5) = desc_dec(n_) int_err(5) = desc_dec(psb_n_)
else if (desc_dec(n_).lt.ix) then else if (desc_dec(psb_n_).lt.ix) then
info=60 info=60
int_err(1) = 4 int_err(1) = 4
int_err(2) = ix int_err(2) = ix
int_err(3) = 6 int_err(3) = 6
int_err(4) = n_ int_err(4) = psb_n_
int_err(5) = desc_dec(n_) int_err(5) = desc_dec(psb_n_)
else if (desc_dec(m_).lt.jx) then else if (desc_dec(psb_m_).lt.jx) then
info=60 info=60
int_err(1) = 5 int_err(1) = 5
int_err(2) = jx int_err(2) = jx
int_err(3) = 6 int_err(3) = 6
int_err(4) = m_ int_err(4) = psb_m_
int_err(5) = desc_dec(m_) int_err(5) = desc_dec(psb_m_)
else if (desc_dec(n_).lt.(ix+m-1)) then else if (desc_dec(psb_n_).lt.(ix+m-1)) then
info=80 info=80
int_err(1) = 1 int_err(1) = 1
int_err(2) = m int_err(2) = m
@ -140,7 +140,7 @@ subroutine psb_chkvect( m, n, lldx, ix, jx, desc_dec, info, iix, jjx)
! Compute local indices for submatrix starting ! Compute local indices for submatrix starting
! at global indices ix and jx ! at global indices ix and jx
if(present(iix)) iix=ix ! (for our applications iix=ix)) if(present(iix)) iix=ix ! (for our applications iix=ix))
if(present(jjx)) iix=ix ! (for our applications jjx=jx)) if(present(jjx)) jjx=jx ! (for our applications jjx=jx))
call psb_erractionrestore(err_act) call psb_erractionrestore(err_act)
return return

@ -224,6 +224,7 @@ subroutine psb_daxpbyv(alpha, x, beta,y,desc_a,info)
integer :: int_err(5), icontxt, nprow, npcol, myrow, mycol,& integer :: int_err(5), icontxt, nprow, npcol, myrow, mycol,&
& err_act, n, iix, jjx, temp(2), ix, iy, ijx, m, iiy, in, jjy & err_act, n, iix, jjx, temp(2), ix, iy, ijx, m, iiy, in, jjy
character(len=20) :: name, ch_err character(len=20) :: name, ch_err
logical, parameter :: debug=.true.
name='psb_daxpby' name='psb_daxpby'
if(psb_get_errstatus().ne.0) return if(psb_get_errstatus().ne.0) return
@ -252,11 +253,18 @@ subroutine psb_daxpbyv(alpha, x, beta,y,desc_a,info)
! check vector correctness ! check vector correctness
call psb_chkvect(m,ione,size(x),ix,ione,desc_a%matrix_data,info,iix,jjx) call psb_chkvect(m,ione,size(x),ix,ione,desc_a%matrix_data,info,iix,jjx)
if(info.ne.0) then
info=4010
ch_err='psb_chkvect 1'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
call psb_chkvect(m,ione,size(y),iy,ione,desc_a%matrix_data,info,iiy,jjy) call psb_chkvect(m,ione,size(y),iy,ione,desc_a%matrix_data,info,iiy,jjy)
if(info.ne.0) then if(info.ne.0) then
info=4010 info=4010
ch_err='psb_chkvect' ch_err='psb_chkvect 2'
call psb_errpush(info,name,a_err=ch_err) call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if end if
if ((iix.ne.ione).or.(iiy.ne.ione)) then if ((iix.ne.ione).or.(iiy.ne.ione)) then

Loading…
Cancel
Save