Updated precset to have a non-optional INFO argument.

Defined a new PSB_WTIME function, and use it in user-level code.
psblas3-type-indexed
Salvatore Filippone 19 years ago
parent acd573ca17
commit 817221494f

@ -51,6 +51,10 @@ module psb_penv_mod
module procedure psb_barrier module procedure psb_barrier
end interface end interface
interface psb_wtime
module procedure psb_wtime
end interface
interface psb_bcast interface psb_bcast
module procedure psb_ibcasts, psb_ibcastv, psb_ibcastm,& module procedure psb_ibcasts, psb_ibcastv, psb_ibcastm,&
& psb_dbcasts, psb_dbcastv, psb_dbcastm,& & psb_dbcasts, psb_dbcastv, psb_dbcastm,&
@ -192,6 +196,13 @@ contains
end subroutine psb_barrier end subroutine psb_barrier
function psb_wtime()
real(kind(1.d0)) :: psb_wtime
real(kind(1.d0)), external :: mpi_wtime
psb_wtime = mpi_wtime()
end function psb_wtime
subroutine psb_abort(ictxt) subroutine psb_abort(ictxt)
integer, intent(in) :: ictxt integer, intent(in) :: ictxt
@ -202,6 +213,7 @@ contains
subroutine psb_info(ictxt,iam,np) subroutine psb_info(ictxt,iam,np)
integer, intent(in) :: ictxt integer, intent(in) :: ictxt
integer, intent(out) :: iam, np integer, intent(out) :: iam, np
integer :: nprow, npcol, myprow, mypcol integer :: nprow, npcol, myprow, mypcol

@ -57,29 +57,29 @@ module psb_prec_mod
end interface end interface
interface psb_precset interface psb_precset
subroutine psb_dprecset(prec,ptype,iv,rs,rv,ierr) subroutine psb_dprecset(prec,ptype,info,iv,rs,rv)
use psb_serial_mod use psb_serial_mod
use psb_descriptor_type use psb_descriptor_type
use psb_prec_type use psb_prec_type
implicit none implicit none
type(psb_dprec_type), intent(inout) :: prec type(psb_dprec_type), intent(inout) :: prec
character(len=*), intent(in) :: ptype character(len=*), intent(in) :: ptype
integer, intent(out) :: info
integer, optional, intent(in) :: iv(:) integer, optional, intent(in) :: iv(:)
real(kind(1.d0)), optional, intent(in) :: rs real(kind(1.d0)), optional, intent(in) :: rs
real(kind(1.d0)), optional, intent(in) :: rv(:) real(kind(1.d0)), optional, intent(in) :: rv(:)
integer, optional, intent(out) :: ierr
end subroutine psb_dprecset end subroutine psb_dprecset
subroutine psb_zprecset(prec,ptype,iv,rs,rv,ierr) subroutine psb_zprecset(prec,ptype,info,iv,rs,rv)
use psb_serial_mod use psb_serial_mod
use psb_descriptor_type use psb_descriptor_type
use psb_prec_type use psb_prec_type
implicit none implicit none
type(psb_zprec_type), intent(inout) :: prec type(psb_zprec_type), intent(inout) :: prec
character(len=*), intent(in) :: ptype character(len=*), intent(in) :: ptype
integer, intent(out) :: info
integer, optional, intent(in) :: iv(:) integer, optional, intent(in) :: iv(:)
real(kind(1.d0)), optional, intent(in) :: rs real(kind(1.d0)), optional, intent(in) :: rs
real(kind(1.d0)), optional, intent(in) :: rv(:) real(kind(1.d0)), optional, intent(in) :: rv(:)
integer, optional, intent(out) :: ierr
end subroutine psb_zprecset end subroutine psb_zprecset
end interface end interface

@ -34,7 +34,7 @@
!!$ POSSIBILITY OF SUCH DAMAGE. !!$ POSSIBILITY OF SUCH DAMAGE.
!!$ !!$
!!$ !!$
subroutine psb_dprecset(p,ptype,iv,rs,rv,info) subroutine psb_dprecset(p,ptype,info,iv,rs,rv)
use psb_serial_mod use psb_serial_mod
use psb_descriptor_type use psb_descriptor_type
@ -43,16 +43,16 @@ subroutine psb_dprecset(p,ptype,iv,rs,rv,info)
type(psb_dprec_type), intent(inout) :: p type(psb_dprec_type), intent(inout) :: p
character(len=*), intent(in) :: ptype character(len=*), intent(in) :: ptype
integer, intent(out) :: info
integer, optional, intent(in) :: iv(:) integer, optional, intent(in) :: iv(:)
real(kind(1.d0)), optional, intent(in) :: rs real(kind(1.d0)), optional, intent(in) :: rs
real(kind(1.d0)), optional, intent(in) :: rv(:) real(kind(1.d0)), optional, intent(in) :: rv(:)
integer, optional, intent(out) :: info
type(psb_dbaseprc_type), pointer :: bpv(:)=>null() type(psb_dbaseprc_type), pointer :: bpv(:)=>null()
character(len=len(ptype)) :: typeup character(len=len(ptype)) :: typeup
integer :: isz, err integer :: isz, err
if (present(info)) info = 0 info = 0
if (.not.associated(p%baseprecv)) then if (.not.associated(p%baseprecv)) then
allocate(p%baseprecv(1),stat=err) allocate(p%baseprecv(1),stat=err)
@ -195,6 +195,6 @@ subroutine psb_dprecset(p,ptype,iv,rs,rv,info)
end select end select
if (present(info)) info = err info = err
end subroutine psb_dprecset end subroutine psb_dprecset

@ -70,7 +70,7 @@ subroutine psb_zilu_bld(a,desc_a,p,upd,info)
character, intent(in) :: upd character, intent(in) :: upd
! .. Local Scalars .. ! .. Local Scalars ..
integer :: i, j, jj, k, kk, m integer :: i, j, jj, k, kk, m, i1, i2, ia
integer :: int_err(5) integer :: int_err(5)
character :: trans, unitd character :: trans, unitd
type(psb_zspmat_type) :: blck, atmp type(psb_zspmat_type) :: blck, atmp

@ -34,7 +34,7 @@
!!$ POSSIBILITY OF SUCH DAMAGE. !!$ POSSIBILITY OF SUCH DAMAGE.
!!$ !!$
!!$ !!$
subroutine psb_zprecset(p,ptype,iv,rs,rv,info) subroutine psb_zprecset(p,ptype,info,iv,rs,rv)
use psb_serial_mod use psb_serial_mod
use psb_descriptor_type use psb_descriptor_type
@ -43,16 +43,16 @@ subroutine psb_zprecset(p,ptype,iv,rs,rv,info)
type(psb_zprec_type), intent(inout) :: p type(psb_zprec_type), intent(inout) :: p
character(len=*), intent(in) :: ptype character(len=*), intent(in) :: ptype
integer, intent(out) :: info
integer, optional, intent(in) :: iv(:) integer, optional, intent(in) :: iv(:)
real(kind(1.d0)), optional, intent(in) :: rs real(kind(1.d0)), optional, intent(in) :: rs
real(kind(1.d0)), optional, intent(in) :: rv(:) real(kind(1.d0)), optional, intent(in) :: rv(:)
integer, optional, intent(out) :: info
type(psb_zbaseprc_type), pointer :: bpv(:)=>null() type(psb_zbaseprc_type), pointer :: bpv(:)=>null()
character(len=len(ptype)) :: typeup character(len=len(ptype)) :: typeup
integer :: isz, err integer :: isz, err
if (present(info)) info = 0 info = 0
if (.not.associated(p%baseprecv)) then if (.not.associated(p%baseprecv)) then
allocate(p%baseprecv(1),stat=err) allocate(p%baseprecv(1),stat=err)
@ -196,6 +196,6 @@ subroutine psb_zprecset(p,ptype,iv,rs,rv,info)
end select end select
if (present(info)) info = err info = err
end subroutine psb_zprecset end subroutine psb_zprecset

@ -91,13 +91,11 @@ program df_sample
! other variables ! other variables
integer :: i,info,j,m_problem integer :: i,info,j,m_problem
integer :: internal, m,ii,nnzero integer :: internal, m,ii,nnzero
real(kind(1.d0)) :: mpi_wtime, t1, t2, tprec, r_amax, b_amax,& real(kind(1.d0)) :: t1, t2, tprec, r_amax, b_amax,&
&scale,resmx,resmxp &scale,resmx,resmxp
integer :: nrhs, nrow, n_row, dim, nv, ne integer :: nrhs, nrow, n_row, dim, nv, ne
integer, pointer :: ivg(:), ipv(:), neigh(:) integer, pointer :: ivg(:), ipv(:), neigh(:)
external mpi_wtime
call psb_init(ictxt) call psb_init(ictxt)
call psb_info(ictxt,iam,np) call psb_info(ictxt,iam,np)
@ -122,7 +120,7 @@ program df_sample
& ipart,afmt,istopc,itmax,itrace,novr,iprec,eps) & ipart,afmt,istopc,itmax,itrace,novr,iprec,eps)
call psb_barrier(ictxt) call psb_barrier(ictxt)
t1 = mpi_wtime() t1 = psb_wtime()
! read the input matrix to be processed and (possibly) the rhs ! read the input matrix to be processed and (possibly) the rhs
nrhs = 1 nrhs = 1
@ -213,7 +211,7 @@ program df_sample
call psb_geall(r_col,desc_a,info) call psb_geall(r_col,desc_a,info)
r_col(:) =0.0 r_col(:) =0.0
call psb_geasb(r_col,desc_a,info) call psb_geasb(r_col,desc_a,info)
t2 = mpi_wtime() - t1 t2 = psb_wtime() - t1
call psb_amx(ictxt, t2) call psb_amx(ictxt, t2)
@ -234,28 +232,28 @@ program df_sample
igsmth=-1 igsmth=-1
select case(iprec) select case(iprec)
case(noprec_) case(noprec_)
call psb_precset(pre,'noprec') call psb_precset(pre,'noprec',info)
case(diagsc_) case(diagsc_)
call psb_precset(pre,'diagsc') call psb_precset(pre,'diagsc',info)
case(bja_) case(bja_)
call psb_precset(pre,'ilu') call psb_precset(pre,'ilu',info)
case(asm_) case(asm_)
call psb_precset(pre,'asm',iv=(/novr,halo_,sum_/)) call psb_precset(pre,'asm',info,iv=(/novr,halo_,sum_/))
case(ash_) case(ash_)
call psb_precset(pre,'asm',iv=(/novr,nohalo_,sum_/)) call psb_precset(pre,'asm',info,iv=(/novr,nohalo_,sum_/))
case(ras_) case(ras_)
call psb_precset(pre,'asm',iv=(/novr,halo_,none_/)) call psb_precset(pre,'asm',info,iv=(/novr,halo_,none_/))
case(rash_) case(rash_)
call psb_precset(pre,'asm',iv=(/novr,nohalo_,none_/)) call psb_precset(pre,'asm',info,iv=(/novr,nohalo_,none_/))
case default case default
call psb_precset(pre,'ilu') call psb_precset(pre,'ilu',info)
end select end select
! building the preconditioner ! building the preconditioner
t1 = mpi_wtime() t1 = psb_wtime()
call psb_precbld(a,desc_a,pre,info) call psb_precbld(a,desc_a,pre,info)
tprec = mpi_wtime()-t1 tprec = psb_wtime()-t1
if (info /= 0) then if (info /= 0) then
call psb_errpush(4010,name,a_err='psb_precbld') call psb_errpush(4010,name,a_err='psb_precbld')
goto 9999 goto 9999
@ -271,11 +269,11 @@ program df_sample
iparm = 0 iparm = 0
call psb_barrier(ictxt) call psb_barrier(ictxt)
t1 = mpi_wtime() t1 = psb_wtime()
call psb_krylov(cmethd,a,pre,b_col,x_col,eps,desc_a,info,& call psb_krylov(cmethd,a,pre,b_col,x_col,eps,desc_a,info,&
& itmax=itmax,iter=iter,err=err,itrace=itrace,istop=istopc,irst=ml) & itmax=itmax,iter=iter,err=err,itrace=itrace,istop=istopc,irst=ml)
call psb_barrier(ictxt) call psb_barrier(ictxt)
t2 = mpi_wtime() - t1 t2 = psb_wtime() - t1
call psb_amx(ictxt,t2) call psb_amx(ictxt,t2)
call psb_geaxpby(1.d0,b_col,0.d0,r_col,desc_a,info) call psb_geaxpby(1.d0,b_col,0.d0,r_col,desc_a,info)
call psb_spmm(-1.d0,a,x_col,1.d0,r_col,desc_a,info) call psb_spmm(-1.d0,a,x_col,1.d0,r_col,desc_a,info)

@ -91,12 +91,11 @@ program zf_sample
! other variables ! other variables
integer :: i,info,j,m_problem integer :: i,info,j,m_problem
integer :: internal, m,ii,nnzero integer :: internal, m,ii,nnzero
real(kind(1.d0)) :: mpi_wtime, t1, t2, tprec, r_amax, b_amax,& real(kind(1.d0)) :: t1, t2, tprec, r_amax, b_amax,&
&scale,resmx,resmxp &scale,resmx,resmxp
integer :: nrhs, nrow, n_row, dim, nv, ne integer :: nrhs, nrow, n_row, dim, nv, ne
integer, pointer :: ivg(:), ipv(:), neigh(:) integer, pointer :: ivg(:), ipv(:), neigh(:)
external mpi_wtime
call psb_init(ictxt) call psb_init(ictxt)
@ -122,7 +121,7 @@ program zf_sample
& ipart,afmt,istopc,itmax,itrace,novr,iprec,eps) & ipart,afmt,istopc,itmax,itrace,novr,iprec,eps)
call psb_barrier(ictxt) call psb_barrier(ictxt)
t1 = mpi_wtime() t1 = psb_wtime()
! read the input matrix to be processed and (possibly) the rhs ! read the input matrix to be processed and (possibly) the rhs
nrhs = 1 nrhs = 1
@ -214,7 +213,7 @@ program zf_sample
call psb_geall(r_col,desc_a,info) call psb_geall(r_col,desc_a,info)
r_col(:) =0.0 r_col(:) =0.0
call psb_geasb(r_col,desc_a,info) call psb_geasb(r_col,desc_a,info)
t2 = mpi_wtime() - t1 t2 = psb_wtime() - t1
call psb_amx(ictxt, t2) call psb_amx(ictxt, t2)
@ -235,27 +234,27 @@ program zf_sample
igsmth=-1 igsmth=-1
select case(iprec) select case(iprec)
case(noprec_) case(noprec_)
call psb_precset(pre,'noprec') call psb_precset(pre,'noprec',info)
case(diagsc_) case(diagsc_)
call psb_precset(pre,'diagsc') call psb_precset(pre,'diagsc',info)
case(bja_) case(bja_)
call psb_precset(pre,'ilu') call psb_precset(pre,'ilu',info)
case(asm_) case(asm_)
call psb_precset(pre,'asm',iv=(/novr,halo_,sum_/)) call psb_precset(pre,'asm',info,iv=(/novr,halo_,sum_/))
case(ash_) case(ash_)
call psb_precset(pre,'asm',iv=(/novr,nohalo_,sum_/)) call psb_precset(pre,'asm',info,iv=(/novr,nohalo_,sum_/))
case(ras_) case(ras_)
call psb_precset(pre,'asm',iv=(/novr,halo_,none_/)) call psb_precset(pre,'asm',info,iv=(/novr,halo_,none_/))
case(rash_) case(rash_)
call psb_precset(pre,'asm',iv=(/novr,nohalo_,none_/)) call psb_precset(pre,'asm',info,iv=(/novr,nohalo_,none_/))
case default case default
call psb_precset(pre,'ilu') call psb_precset(pre,'ilu',info)
end select end select
! building the preconditioner ! building the preconditioner
t1 = mpi_wtime() t1 = psb_wtime()
call psb_precbld(a,desc_a,pre,info) call psb_precbld(a,desc_a,pre,info)
tprec = mpi_wtime()-t1 tprec = psb_wtime()-t1
if (info /= 0) then if (info /= 0) then
call psb_errpush(4010,name,a_err='psb_precbld') call psb_errpush(4010,name,a_err='psb_precbld')
goto 9999 goto 9999
@ -271,11 +270,11 @@ program zf_sample
iparm = 0 iparm = 0
call psb_barrier(ictxt) call psb_barrier(ictxt)
t1 = mpi_wtime() t1 = psb_wtime()
call psb_krylov(cmethd,a,pre,b_col,x_col,eps,desc_a,info,& call psb_krylov(cmethd,a,pre,b_col,x_col,eps,desc_a,info,&
& itmax=itmax,iter=iter,err=err,itrace=itrace,istop=istopc,irst=ml) & itmax=itmax,iter=iter,err=err,itrace=itrace,istop=istopc,irst=ml)
call psb_barrier(ictxt) call psb_barrier(ictxt)
t2 = mpi_wtime() - t1 t2 = psb_wtime() - t1
call psb_amx(ictxt,t2) call psb_amx(ictxt,t2)
call psb_geaxpby(zone,b_col,zzero,r_col,desc_a,info) call psb_geaxpby(zone,b_col,zzero,r_col,desc_a,info)
call psb_spmm(-zone,a,x_col,zone,r_col,desc_a,info) call psb_spmm(-zone,a,x_col,zone,r_col,desc_a,info)

@ -89,8 +89,7 @@ program pde90
! miscellaneous ! miscellaneous
integer :: iargc,convert_descr,dim, check_descr integer :: iargc,convert_descr,dim, check_descr
real(kind(1.d0)), parameter :: one = 1.d0 real(kind(1.d0)), parameter :: one = 1.d0
real(kind(1.d0)) :: mpi_wtime, t1, t2, tprec, tsolve, t3, t4 real(kind(1.d0)) :: t1, t2, tprec, tsolve, t3, t4
external mpi_wtime
! sparse matrix and preconditioner ! sparse matrix and preconditioner
type(psb_dspmat_type) :: a, l, u, h type(psb_dspmat_type) :: a, l, u, h
@ -140,9 +139,9 @@ program pde90
! !
call psb_barrier(ictxt) call psb_barrier(ictxt)
t1 = mpi_wtime() t1 = psb_wtime()
call create_matrix(idim,a,b,x,desc_a,part_block,ictxt,afmt,info) call create_matrix(idim,a,b,x,desc_a,part_block,ictxt,afmt,info)
t2 = mpi_wtime() - t1 t2 = psb_wtime() - t1
if(info.ne.0) then if(info.ne.0) then
info=4010 info=4010
ch_err='create_matrix' ch_err='create_matrix'
@ -160,25 +159,25 @@ program pde90
if(iam == psb_root_) write(0,'("Setting preconditioner to : ",a)')pr_to_str(iprec) if(iam == psb_root_) write(0,'("Setting preconditioner to : ",a)')pr_to_str(iprec)
select case(iprec) select case(iprec)
case(noprec_) case(noprec_)
call psb_precset(pre,'noprec') call psb_precset(pre,'noprec',info)
case(diagsc_) case(diagsc_)
call psb_precset(pre,'diagsc') call psb_precset(pre,'diagsc',info)
case(bja_) case(bja_)
call psb_precset(pre,'ilu') call psb_precset(pre,'ilu',info)
case(asm_) case(asm_)
call psb_precset(pre,'asm',iv=(/novr,halo_,sum_/)) call psb_precset(pre,'asm',info,iv=(/novr,halo_,sum_/))
case(ash_) case(ash_)
call psb_precset(pre,'asm',iv=(/novr,nohalo_,sum_/)) call psb_precset(pre,'asm',info,iv=(/novr,nohalo_,sum_/))
case(ras_) case(ras_)
call psb_precset(pre,'asm',iv=(/novr,halo_,none_/)) call psb_precset(pre,'asm',info,iv=(/novr,halo_,none_/))
case(rash_) case(rash_)
call psb_precset(pre,'asm',iv=(/novr,nohalo_,none_/)) call psb_precset(pre,'asm',info,iv=(/novr,nohalo_,none_/))
case default case default
call psb_precset(pre,'ilu') call psb_precset(pre,'ilu',info)
end select end select
call psb_barrier(ictxt) call psb_barrier(ictxt)
t1 = mpi_wtime() t1 = psb_wtime()
call psb_precbld(a,desc_a,pre,info) call psb_precbld(a,desc_a,pre,info)
if(info.ne.0) then if(info.ne.0) then
info=4010 info=4010
@ -187,7 +186,7 @@ program pde90
goto 9999 goto 9999
end if end if
tprec = mpi_wtime()-t1 tprec = psb_wtime()-t1
call psb_amx(ictxt,tprec) call psb_amx(ictxt,tprec)
@ -199,7 +198,7 @@ program pde90
! !
if(iam == psb_root_) write(*,'("Calling iterative method ",a)')cmethd if(iam == psb_root_) write(*,'("Calling iterative method ",a)')cmethd
call psb_barrier(ictxt) call psb_barrier(ictxt)
t1 = mpi_wtime() t1 = psb_wtime()
eps = 1.d-9 eps = 1.d-9
call psb_krylov(cmethd,a,pre,b,x,eps,desc_a,info,& call psb_krylov(cmethd,a,pre,b,x,eps,desc_a,info,&
& itmax=itmax,iter=iter,err=err,itrace=itrace,istop=istopc,irst=ml) & itmax=itmax,iter=iter,err=err,itrace=itrace,istop=istopc,irst=ml)
@ -212,7 +211,7 @@ program pde90
end if end if
call psb_barrier(ictxt) call psb_barrier(ictxt)
t2 = mpi_wtime() - t1 t2 = psb_wtime() - t1
call psb_amx(ictxt,t2) call psb_amx(ictxt,t2)
if (iam == 0) then if (iam == 0) then
@ -427,9 +426,9 @@ contains
! deltat discretization time ! deltat discretization time
real(kind(1.d0)) :: deltah real(kind(1.d0)) :: deltah
real(kind(1.d0)),parameter :: rhs=0.d0,one=1.d0,zero=0.d0 real(kind(1.d0)),parameter :: rhs=0.d0,one=1.d0,zero=0.d0
real(kind(1.d0)) :: mpi_wtime, t1, t2, t3, tins, tasb real(kind(1.d0)) :: t1, t2, t3, tins, tasb
real(kind(1.d0)) :: a1, a2, a3, a4, b1, b2, b3 real(kind(1.d0)) :: a1, a2, a3, a4, b1, b2, b3
external mpi_wtime,a1, a2, a3, a4, b1, b2, b3 external :: a1, a2, a3, a4, b1, b2, b3
integer :: nb, ir1, ir2, ipr, err_act integer :: nb, ir1, ir2, ipr, err_act
logical :: own logical :: own
! common area ! common area
@ -478,7 +477,7 @@ contains
tins = 0.d0 tins = 0.d0
call psb_barrier(ictxt) call psb_barrier(ictxt)
t1 = mpi_wtime() t1 = psb_wtime()
! loop over rows belonging to current process in a block ! loop over rows belonging to current process in a block
! distribution. ! distribution.
@ -606,10 +605,10 @@ contains
irow(1:element-1)=glob_row irow(1:element-1)=glob_row
ia=glob_row ia=glob_row
t3 = mpi_wtime() t3 = psb_wtime()
call psb_spins(element-1,irow,icol,val,a,desc_a,info) call psb_spins(element-1,irow,icol,val,a,desc_a,info)
if(info.ne.0) exit if(info.ne.0) exit
tins = tins + (mpi_wtime()-t3) tins = tins + (psb_wtime()-t3)
call psb_geins(1,(/ia/),zt(1:1),b,desc_a,info) call psb_geins(1,(/ia/),zt(1:1),b,desc_a,info)
if(info.ne.0) exit if(info.ne.0) exit
zt(1)=0.d0 zt(1)=0.d0
@ -620,7 +619,7 @@ contains
end do end do
call psb_barrier(ictxt) call psb_barrier(ictxt)
t2 = mpi_wtime()-t1 t2 = psb_wtime()-t1
if(info.ne.0) then if(info.ne.0) then
info=4010 info=4010
@ -631,11 +630,11 @@ contains
deallocate(val,irow,icol) deallocate(val,irow,icol)
t1 = mpi_wtime() t1 = psb_wtime()
call psb_cdasb(desc_a,info) call psb_cdasb(desc_a,info)
call psb_spasb(a,desc_a,info,dupl=psb_dupl_err_,afmt=afmt) call psb_spasb(a,desc_a,info,dupl=psb_dupl_err_,afmt=afmt)
call psb_barrier(ictxt) call psb_barrier(ictxt)
tasb = mpi_wtime()-t1 tasb = psb_wtime()-t1
if(info.ne.0) then if(info.ne.0) then
info=4010 info=4010
ch_err='asb rout.' ch_err='asb rout.'

Loading…
Cancel
Save