psblas3-integer8:

base/internals/psi_crea_index.f90
 base/internals/psi_desc_index.F90
 base/internals/psi_sort_dl.f90
 base/internals/srtlist.f
 base/modules/psb_error_impl.F90
 base/modules/psi_penv_mod.F90
 test/kernel/pdgenspmv.f90
 test/kernel/runs/fspmv.inp
 test/kernel/runs/spmv.inp
 test/pargen/ppde.f90

Progress: now test/kernel/pdgenspmv begins to work
 (but cdbldext needs fixing still)
psblas3-type-indexed
Salvatore Filippone 13 years ago
parent d069551ca7
commit 9c8ada8c2b

@ -126,11 +126,14 @@ subroutine psi_desc_index(desc,index_in,dep_list,&
integer(psb_ipk_) :: ictxt
integer(psb_ipk_), parameter :: no_comm=-1
! ...local arrays..
integer(psb_ipk_),allocatable :: brvindx(:),rvsz(:),&
& bsdindx(:),sdsz(:), sndbuf(:), rcvbuf(:)
integer(psb_ipk_),allocatable :: sndbuf(:), rcvbuf(:)
integer(psb_mpik_),allocatable :: brvindx(:),rvsz(:),&
& bsdindx(:),sdsz(:)
integer(psb_ipk_) :: ihinsz,ntot,k,err_act,nidx,&
& idxr, idxs, iszs, iszr, nesd, nerv, icomm
& idxr, idxs, iszs, iszr, nesd, nerv
integer(psb_mpik_) :: icomm, minfo
logical,parameter :: usempi=.true.
integer(psb_ipk_) :: debug_level, debug_unit
@ -169,8 +172,8 @@ subroutine psi_desc_index(desc,index_in,dep_list,&
goto 9999
end if
sdsz(:) = 0
rvsz(:) = 0
sdsz(:) = 0
rvsz(:) = 0
bsdindx(:) = 0
brvindx(:) = 0
i = 1
@ -183,8 +186,8 @@ subroutine psi_desc_index(desc,index_in,dep_list,&
i = i + nerv + 1
end do
ihinsz=i
call mpi_alltoall(sdsz,1,psb_mpi_integer,rvsz,1,psb_mpi_integer,icomm,info)
if(info /= psb_success_) then
call mpi_alltoall(sdsz,1,mpi_integer,rvsz,1,mpi_integer,icomm,minfo)
if (minfo /= psb_success_) then
call psb_errpush(psb_err_from_subroutine_,name,a_err='mpi_alltoall')
goto 9999
end if
@ -290,8 +293,8 @@ subroutine psi_desc_index(desc,index_in,dep_list,&
end do
call mpi_alltoallv(sndbuf,sdsz,bsdindx,psb_mpi_integer,&
& rcvbuf,rvsz,brvindx,psb_mpi_integer,icomm,info)
if(info /= psb_success_) then
& rcvbuf,rvsz,brvindx,psb_mpi_integer,icomm,minfo)
if (minfo /= psb_success_) then
call psb_errpush(psb_err_from_subroutine_,name,a_err='mpi_alltoallv')
goto 9999
end if

@ -72,9 +72,8 @@ subroutine psi_sort_dl(dep_list,l_dep_list,np,info)
allocate(work(isz))
! call srtlist(dep_list, dl_lda, l_dep_list, np, info)
call srtlist(dep_list,size(dep_list,1),l_dep_list,np,work(idg),&
call srtlist(dep_list,size(dep_list,1,kind=psb_ipk_),l_dep_list,np,work(idg),&
& work(idgp),work(iupd),work(iedges),work(iidx),work(iich),info)
if (info /= psb_success_) then
call psb_errpush(psb_err_from_subroutine_,name,a_err='srtlist')
goto 9999

@ -102,6 +102,7 @@ C
DG(I)=LDL(I)
ENDDO
NEDGES = 0
DO I=1, NP
DO J=1, DG(I)

@ -4,7 +4,6 @@ subroutine psb_errcomm(ictxt, err)
use psb_penv_mod
integer(psb_mpik_), intent(in) :: ictxt
integer(psb_ipk_), intent(inout):: err
integer(psb_ipk_) :: temp(2)
call psb_amx(ictxt, err)

@ -137,23 +137,30 @@ contains
#if defined(LONG_INTEGERS)
subroutine psb_init_ipk(ictxt,np,basectxt,ids)
integer(psb_ipk_), intent(out) :: ictxt
integer(psb_ipk_), intent(in), optional :: np, basectxt
integer(psb_ipk_), intent(in), optional :: np, basectxt, ids(:)
integer(psb_mpik_) :: iictxt
integer(psb_mpik_) :: inp, ibasectxt
integer(psb_mpik_), allocatable :: ids_(:)
if (present(ids)) then
allocate(ids_(size(ids)))
ids_ = ids
else
allocate(ids_(0))
end if
if (present(np).and.present(basectxt)) then
inp = np
ibasectxt = basectxt
call psb_init(iictxt,np=inp,basectxt=ibasectxt)
call psb_init(iictxt,np=inp,basectxt=ibasectxt,ids=ids_)
else if (present(np)) then
inp = np
call psb_init(iictxt,np=inp)
call psb_init(iictxt,np=inp,ids=ids_)
else if (present(basectxt)) then
ibasectxt = basectxt
call psb_init(iictxt,basectxt=ibasectxt)
call psb_init(iictxt,basectxt=ibasectxt,ids=ids_)
else
call psb_init(iictxt)
call psb_init(iictxt,ids=ids_)
end if
ictxt = iictxt
end subroutine psb_init_ipk
@ -491,8 +498,8 @@ contains
! !!!!!!!!!!!!!!!!!!!!!!
subroutine psi_iamx_op(inv, outv,len,type)
integer(psb_ipk_) :: inv(*),outv(*)
integer(psb_ipk_) :: len,type
integer(psb_ipk_) :: i
integer(psb_mpik_) :: len,type
integer(psb_mpik_) :: i
do i=1, len
if (abs(inv(i)) > abs(outv(i))) outv(i) = inv(i)
@ -501,8 +508,8 @@ contains
subroutine psi_iamn_op(inv, outv,len,type)
integer(psb_ipk_) :: inv(*),outv(*)
integer(psb_ipk_) :: len,type
integer(psb_ipk_) :: i
integer(psb_mpik_) :: len,type
integer(psb_mpik_) :: i
do i=1, len
if (abs(inv(i)) < abs(outv(i))) outv(i) = inv(i)
end do
@ -510,8 +517,8 @@ contains
subroutine psi_i8amx_op(inv, outv,len,type)
integer(psb_long_int_k_) :: inv(*),outv(*)
integer(psb_ipk_) :: len,type
integer(psb_ipk_) :: i
integer(psb_mpik_) :: len,type
integer(psb_mpik_) :: i
do i=1, len
if (abs(inv(i)) > abs(outv(i))) outv(i) = inv(i)
@ -527,8 +534,8 @@ contains
include 'mpif.h'
#endif
integer(psb_long_int_k_) :: inv(*),outv(*)
integer(psb_ipk_) :: len,type
integer(psb_ipk_) :: i
integer(psb_mpik_) :: len,type
integer(psb_mpik_) :: i
if (type /= mpi_integer8) then
write(psb_err_unit,*) 'Invalid type !!!'
end if
@ -538,88 +545,88 @@ contains
end subroutine psi_i8amn_op
subroutine psi_samx_op(vin,vinout,len,itype)
integer(psb_ipk_), intent(in) :: len, itype
integer(psb_mpik_), intent(in) :: len, itype
real(psb_spk_), intent(in) :: vin(len)
real(psb_spk_), intent(inout) :: vinout(len)
integer(psb_ipk_) :: i
integer(psb_mpik_) :: i
do i=1, len
if (abs(vinout(i)) < abs(vin(i))) vinout(i) = vin(i)
end do
end subroutine psi_samx_op
subroutine psi_samn_op(vin,vinout,len,itype)
integer(psb_ipk_), intent(in) :: len, itype
integer(psb_mpik_), intent(in) :: len, itype
real(psb_spk_), intent(in) :: vin(len)
real(psb_spk_), intent(inout) :: vinout(len)
integer(psb_ipk_) :: i
integer(psb_mpik_) :: i
do i=1, len
if (abs(vinout(i)) > abs(vin(i))) vinout(i) = vin(i)
end do
end subroutine psi_samn_op
subroutine psi_damx_op(vin,vinout,len,itype)
integer(psb_ipk_), intent(in) :: len, itype
integer(psb_mpik_), intent(in) :: len, itype
real(psb_dpk_), intent(in) :: vin(len)
real(psb_dpk_), intent(inout) :: vinout(len)
integer(psb_ipk_) :: i
integer(psb_mpik_) :: i
do i=1, len
if (abs(vinout(i)) < abs(vin(i))) vinout(i) = vin(i)
end do
end subroutine psi_damx_op
subroutine psi_damn_op(vin,vinout,len,itype)
integer(psb_ipk_), intent(in) :: len, itype
integer(psb_mpik_), intent(in) :: len, itype
real(psb_dpk_), intent(in) :: vin(len)
real(psb_dpk_), intent(inout) :: vinout(len)
integer(psb_ipk_) :: i
integer(psb_mpik_) :: i
do i=1, len
if (abs(vinout(i)) > abs(vin(i))) vinout(i) = vin(i)
end do
end subroutine psi_damn_op
subroutine psi_camx_op(vin,vinout,len,itype)
integer(psb_ipk_), intent(in) :: len, itype
integer(psb_mpik_), intent(in) :: len, itype
complex(psb_spk_), intent(in) :: vin(len)
complex(psb_spk_), intent(inout) :: vinout(len)
integer(psb_ipk_) :: i
integer(psb_mpik_) :: i
do i=1, len
if (abs(vinout(i)) < abs(vin(i))) vinout(i) = vin(i)
end do
end subroutine psi_camx_op
subroutine psi_camn_op(vin,vinout,len,itype)
integer(psb_ipk_), intent(in) :: len, itype
integer(psb_mpik_), intent(in) :: len, itype
complex(psb_spk_), intent(in) :: vin(len)
complex(psb_spk_), intent(inout) :: vinout(len)
integer(psb_ipk_) :: i
integer(psb_mpik_) :: i
do i=1, len
if (abs(vinout(i)) > abs(vin(i))) vinout(i) = vin(i)
end do
end subroutine psi_camn_op
subroutine psi_zamx_op(vin,vinout,len,itype)
integer(psb_ipk_), intent(in) :: len, itype
integer(psb_mpik_), intent(in) :: len, itype
complex(psb_dpk_), intent(in) :: vin(len)
complex(psb_dpk_), intent(inout) :: vinout(len)
integer(psb_ipk_) :: i
integer(psb_mpik_) :: i
do i=1, len
if (abs(vinout(i)) < abs(vin(i))) vinout(i) = vin(i)
end do
end subroutine psi_zamx_op
subroutine psi_zamn_op(vin,vinout,len,itype)
integer(psb_ipk_), intent(in) :: len, itype
integer(psb_mpik_), intent(in) :: len, itype
complex(psb_dpk_), intent(in) :: vin(len)
complex(psb_dpk_), intent(inout) :: vinout(len)
integer(psb_ipk_) :: i
integer(psb_mpik_) :: i
do i=1, len
if (abs(vinout(i)) > abs(vin(i))) vinout(i) = vin(i)
end do
@ -627,11 +634,11 @@ contains
subroutine psi_snrm2_op(vin,vinout,len,itype)
implicit none
integer(psb_ipk_), intent(in) :: len, itype
integer(psb_mpik_), intent(in) :: len, itype
real(psb_spk_), intent(in) :: vin(len)
real(psb_spk_), intent(inout) :: vinout(len)
integer(psb_ipk_) :: i
integer(psb_mpik_) :: i
real(psb_spk_) :: w, z
do i=1, len
w = max( vin(i), vinout(i) )
@ -646,11 +653,11 @@ contains
subroutine psi_dnrm2_op(vin,vinout,len,itype)
implicit none
integer(psb_ipk_), intent(in) :: len, itype
integer(psb_mpik_), intent(in) :: len, itype
real(psb_dpk_), intent(in) :: vin(len)
real(psb_dpk_), intent(inout) :: vinout(len)
integer(psb_ipk_) :: i
integer(psb_mpik_) :: i
real(psb_dpk_) :: w, z
do i=1, len
w = max( vin(i), vinout(i) )

@ -152,12 +152,12 @@ program pdgen
call psb_abort(ictxt)
end if
call xxv%set(done)
call psb_barrier(ictxt)
t1 = psb_wtime()
do i=1,times
write(0,*) 'Iteration ',i,' of ',times
call psb_spmm(done,a,xxv,dzero,bv,desc_a,info,'n')
end do
call psb_barrier(ictxt)
@ -380,6 +380,7 @@ contains
nt = nr
call psb_sum(ictxt,nt)
write(0,*) iam,' Going to generate ',nt
if (nt /= m) write(psb_err_unit,*) iam, 'Initialization error ',nr,nt,m
call psb_barrier(ictxt)
t0 = psb_wtime()
@ -392,6 +393,7 @@ contains
call psb_barrier(ictxt)
talc = psb_wtime()-t0
write(0,*) iam,' Done allocation'
if (info /= psb_success_) then
info=psb_err_from_subroutine_
ch_err='allocation rout.'
@ -417,7 +419,7 @@ contains
call psb_loc_to_glob(myidx,desc_a,info)
write(0,*) iam,' Done loc_to_glob'
! loop over rows belonging to current process in a block
! distribution.
@ -426,6 +428,7 @@ contains
do ii=1, nlr,nb
ib = min(nb,nlr-ii+1)
element = 1
write(0,*) iam,' iteration ',ii
do k=1,ib
i=ii+k-1
! local matrix pointer
@ -528,7 +531,7 @@ contains
call psb_geins(ib,myidx(ii:ii+ib-1),zt(1:ib),xxv,desc_a,info)
if(info /= psb_success_) exit
end do
write(0,*) iam,' end loop'
tgen = psb_wtime()-t1
if(info /= psb_success_) then
info=psb_err_from_subroutine_
@ -538,10 +541,11 @@ contains
end if
deallocate(val,irow,icol)
write(0,*) iam,' Going for cdasb'
call psb_barrier(ictxt)
t1 = psb_wtime()
call psb_cdasb(desc_a,info)
write(0,*) iam,' Done cdasb'
if (info == psb_success_) &
& call psb_spasb(a,desc_a,info,dupl=psb_dupl_err_,afmt=afmt)
call psb_barrier(ictxt)
@ -551,6 +555,7 @@ contains
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
write(0,*) iam,' Done cdasb/spasb'
if (info == psb_success_) call psb_geasb(xxv,desc_a,info)
if (info == psb_success_) call psb_geasb(bv,desc_a,info)
if(info /= psb_success_) then

@ -0,0 +1,5 @@
ASIC_100ks.mtx
MM
0

@ -1,5 +1,3 @@
ASIC_100ks.mtx
MM
0
CSR
10

@ -92,7 +92,6 @@ program ppde
integer(psb_ipk_) :: iter, itmax,itrace, istopc, irst
integer(psb_long_int_k_) :: amatsize, precsize, descsize, d2size
real(psb_dpk_) :: err, eps
integer(psb_mpik_) :: iict
! other variables
integer(psb_ipk_) :: info, i
character(len=20) :: name,ch_err
@ -101,8 +100,7 @@ program ppde
info=psb_success_
call psb_init(iict)
ictxt = iict
call psb_init(ictxt)
call psb_info(ictxt,iam,np)
write(0,*) 'Fromt init/info',iam,np
if (iam < 0) then

Loading…
Cancel
Save