Fix sorting of dep_list: store in CSR-like, handle self-loops

pizdaint-runs
Salvatore Filippone 4 years ago
parent 781f0ef083
commit 62a67a0d0e

@ -191,11 +191,16 @@ subroutine psi_i_crea_index(desc_a,index_in,index_out,nxch,nsnd,nrcv,info)
if (do_timings) call psb_tic(idx_phase2)
call psi_bld_glb_dep_list(ictxt,&
& loc_dl,length_dl,c_dep_list,dl_ptr,info)
if (info /= 0) then
write(0,*) me,trim(name),' From bld_glb_list ',info
end if
!!$ call psi_dl_check(dep_list,dl_lda,np,length_dl)
!!$
!!$ ! ....now i can sort dependency lists.
call psi_sort_dl(dl_ptr,c_dep_list,length_dl,np,info)
if (info /= 0) then
write(0,*) me,trim(name),' From sort_dl ',info
end if
ldl = length_dl(me)
loc_dl = c_dep_list(dl_ptr(me):dl_ptr(me)+ldl-1)
@ -253,7 +258,7 @@ contains
logical :: val
val = .not.(((dlmax>(26*4)).or.((dlavg>=(26*2)).and.(np>=128))))
val = .true.
!val = .true.
end function choose_sorting
end subroutine psi_i_crea_index

@ -136,7 +136,7 @@ end subroutine psi_i_sort_dl
! node in the dependency list for the current one *
! *
!**********************************************************************
subroutine psi_i_csr_sort_dl(dl_ptr,c_dep_list,l_dep_list,np,info)
subroutine psi_i_csr_sort_dl(dl_ptr,c_dep_list,l_dep_list,ictxt,info)
use psi_mod, psb_protect_name => psi_i_csr_sort_dl
use psb_const_mod
use psb_error_mod
@ -144,14 +144,18 @@ subroutine psi_i_csr_sort_dl(dl_ptr,c_dep_list,l_dep_list,np,info)
implicit none
integer(psb_ipk_), intent(inout) :: c_dep_list(:), dl_ptr(0:), l_dep_list(0:)
integer(psb_ipk_), intent(in) :: np
integer(psb_ipk_), intent(in) :: ictxt
integer(psb_ipk_), intent(out) :: info
! Local variables
integer(psb_ipk_), allocatable :: dg(:), dgp(:),&
& idx(:), upd(:), edges(:,:), ich(:)
integer(psb_ipk_) :: i, j, nedges, ip1, ip2, nch, ip, iedge,&
& i1, ix, ist, iswap(2)
logical :: internal_error
integer(psb_ipk_) :: me, np
info = 0
call psb_info(ictxt,me,np)
nedges = size(c_dep_list)
allocate(dg(0:np-1),dgp(nedges),edges(2,nedges),upd(0:np-1),&
@ -173,7 +177,7 @@ subroutine psi_i_csr_sort_dl(dl_ptr,c_dep_list,l_dep_list,np,info)
do i = 0, np-1
do j = dl_ptr(i),dl_ptr(i+1) - 1
ip = c_dep_list(j)
if (i<ip) then
if (i<=ip) then
nedges = nedges + 1
edges(1,nedges) = i
edges(2,nedges) = ip
@ -243,15 +247,22 @@ subroutine psi_i_csr_sort_dl(dl_ptr,c_dep_list,l_dep_list,np,info)
dg(i) = dg(i) + upd(i)
end do
end do
internal_error = .false.
do i=0, np-1
if (dg(i) /= 0) then
write(psb_err_unit,*)&
& 'SRTLIST Error on exit:',i,dg(i)
internal_error = .true.
if (me == 0) write(psb_err_unit,*)&
& 'csr_SRTLIST Error on exit:',i,dg(i)
end if
dg(i) = 0
end do
if (internal_error .and. (me==0)) then
write(0,*) 'Error on srt_list. Input:'
do i = 0, np-1
write(0,*) 'Proc: ',i,' list: '
write(0,*) c_dep_list(dl_ptr(i):dl_ptr(i+1) - 1)
end do
end if
!
! 10. Scan the edge sequence;
! for each edge, take each one of its
@ -268,11 +279,16 @@ subroutine psi_i_csr_sort_dl(dl_ptr,c_dep_list,l_dep_list,np,info)
ix = dl_ptr(i)
c_dep_list(ix+dg(i)) = edges(1,j)
dg(i) = dg(i)+1
!
! If there are any self loops, adjust for error condition
! check
!
if (edges(1,j) == edges(2,j)) dg(i) = dg(i) -1
end do
do i=0, np-1
if (dg(i) /= l_dep_list(i)) then
write(psb_err_unit,*) &
if (me == 0) write(psb_err_unit,*) &
& 'SRTLIST Mismatch on output',i,dg(i),l_dep_list(i)
end if
end do

@ -114,7 +114,7 @@ subroutine srtlist(dep_list,dl_lda,ldl,np,dg,dgp,upd, edges,idx,ich,info)
do i=1, np
do j=1, dg(i)
ip = dep_list(j,i) + 1
if (ip.gt.i) then
if (ip >= i) then
iedge = iedge + 1
edges(1,iedge) = i
edges(2,iedge) = ip
@ -187,6 +187,7 @@ subroutine srtlist(dep_list,dl_lda,ldl,np,dg,dgp,upd, edges,idx,ich,info)
i = edges(2,j)
dg(i) = dg(i)+1
dep_list(dg(i),i) = edges(1,j)-1
if (edges(1,j) == edges(2,j)) dg(i) = dg(i) -1
enddo
do i=1, np
if (dg(i).ne.ldl(i)) then

@ -111,13 +111,13 @@ module psi_i_mod
integer(psb_ipk_) :: np
integer(psb_ipk_) :: info
end subroutine psi_i_sort_dl
subroutine psi_i_csr_sort_dl(dl_ptr,c_dep_list,l_dep_list,np,info)
subroutine psi_i_csr_sort_dl(dl_ptr,c_dep_list,l_dep_list,ictxt,info)
import
implicit none
integer(psb_ipk_), intent(in) :: c_dep_list(:), dl_ptr(0:)
integer(psb_ipk_), intent(inout) :: l_dep_list(0:)
integer(psb_ipk_), intent(in) :: np
integer(psb_ipk_), intent(out) :: info
integer(psb_ipk_) :: ictxt
integer(psb_ipk_) :: info
end subroutine psi_i_csr_sort_dl
end interface
@ -151,7 +151,7 @@ module psi_i_mod
integer(psb_ipk_), allocatable, intent(out) :: c_dep_list(:), dl_ptr(:)
integer(psb_ipk_), intent(out) :: info
end subroutine psi_i_bld_glb_csr_dep_list
end interface psi_bld_glb_dep_list
end interface
interface psi_extract_loc_dl
subroutine psi_i_xtr_loc_dl(ictxt,is_bld,is_upd,desc_str,loc_dl,length_dl,info)

Loading…
Cancel
Save