New fnd_owner implementation, taking into account CD%REINIT

merge-paraggr-newops^2
Salvatore Filippone 5 years ago
parent 055e342253
commit 58b7489db9

@ -162,7 +162,7 @@ subroutine psi_graph_fnd_owner(idx,iprc,idxmap,info)
!
call psb_safe_ab_cpy(idxmap%p_adjcncy,ladj,info)
nadj = psb_size(ladj)
! This makes ladj allocated with size 0 just in case
! This makes ladj allocated with size 0 if needed, as opposed to unallocated
call psb_realloc(nadj,ladj,info)
!
! Throughout the subroutine, nreqst is the number of local inquiries

@ -74,7 +74,8 @@ subroutine psi_indx_map_fnd_owner(idx,iprc,idxmap,info)
integer(psb_ipk_), allocatable :: hhidx(:)
integer(psb_mpk_) :: icomm, minfo, iictxt
integer(psb_ipk_) :: i, err_act, hsize, nv
integer(psb_ipk_) :: i, err_act, hsize
integer(psb_lpk_) :: nv
integer(psb_lpk_) :: mglob
integer(psb_ipk_) :: ictxt,np,me, nresp
logical, parameter :: gettime=.false.
@ -140,7 +141,66 @@ subroutine psi_indx_map_fnd_owner(idx,iprc,idxmap,info)
else
call psi_graph_fnd_owner(idx,iprc,idxmap,info)
if (allocated(idxmap%halo_owner)) then
!
! Maybe we are coming here after a REINIT event.
! In this case, reuse the existing information as much as possible.
!
block
integer(psb_ipk_), allocatable :: tprc(:), lidx(:)
integer(psb_lpk_), allocatable :: tidx(:)
integer(psb_lpk_) :: k1, k2, nh
allocate(lidx(nv),stat=info)
if (info /= psb_success_) then
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
goto 9999
end if
!
! Get local answers, if any
!
call idxmap%g2l(idx,lidx,info,owned=.false.)
call idxmap%fnd_halo_owner(lidx,iprc,info)
nh = count(iprc<0)
!write(0,*) me,'Going through new impl from ',nv,' to ',nh
allocate(tidx(nh),tprc(nh),stat=info)
if (info /= psb_success_) then
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
goto 9999
end if
!
! Prepare remote queries
!
k2 = 0
do k1 = 1, nv
if (iprc(k1) < 0) then
k2 = k2 + 1
if (k2 > nh) then
info = psb_err_internal_error_
call psb_errpush(info,name,a_err='Wrong auxiliary count')
goto 9999
end if
tidx(k2) = idx(k1)
end if
end do
call psi_graph_fnd_owner(tidx,tprc,idxmap,info)
k2 = 0
do k1 = 1, nv
if (iprc(k1) < 0) then
k2 = k2 + 1
if (k2 > nh) then
info = psb_err_internal_error_
call psb_errpush(info,name,a_err='Wrong auxiliary count')
goto 9999
end if
iprc(k1) = tprc(k2)
end if
end do
end block
else
call psi_graph_fnd_owner(idx,iprc,idxmap,info)
end if
end if

@ -1503,6 +1503,7 @@ contains
subroutine base_fnd_halo_owner_s(idxmap,xin,xout,info)
use psb_penv_mod
use psb_error_mod
use psb_realloc_mod
implicit none
class(psb_indx_map), intent(inout) :: idxmap
integer(psb_ipk_), intent(in) :: xin
@ -1512,16 +1513,17 @@ contains
integer(psb_ipk_) :: i, j, nr, nc, nh
nr = idxmap%local_rows
nc = idxmap%local_cols
nc = min(idxmap%local_cols, (nr+psb_size(idxmap%halo_owner)))
xout = -1
if (.not.allocated(idxmap%halo_owner)) then
!write(0,*) 'Halo_owner not allocated!', nr, nc, xin
return
end if
if ((nr<xin).and.(xin <= nc)) then
if (size(idxmap%halo_owner)<(xin-nr)) then
!write(0,*) 'Halo_owner bad size',xin,nr,xin-nr,size(idxmap%halo_owner)
return
end if
!!$ if (size(idxmap%halo_owner)<(xin-nr)) then
!!$ !write(0,*) 'Halo_owner bad size',xin,nr,xin-nr,size(idxmap%halo_owner)
!!$ return
!!$ end if
xout = idxmap%halo_owner(xin-nr)
end if
@ -1530,6 +1532,7 @@ contains
subroutine base_fnd_halo_owner_v(idxmap,xin,xout,info)
use psb_penv_mod
use psb_error_mod
use psb_realloc_mod
implicit none
class(psb_indx_map), intent(inout) :: idxmap
integer(psb_ipk_), intent(in) :: xin(:)
@ -1538,13 +1541,15 @@ contains
integer(psb_ipk_) :: i, j, nr, nc, nh, sz
nr = idxmap%local_rows
nc = idxmap%local_cols
nc = min(idxmap%local_cols, (nr+psb_size(idxmap%halo_owner)))
sz = min(size(xin),size(xout))
xout = -1
do i = 1, sz
if ((nr<xin(i)).and.(xin(i) <= nc)) xout = idxmap%halo_owner(xin(i)-nr)
xout(i) = -1
if ((nr<xin(i)).and.(xin(i) <= nc)) xout(i) = idxmap%halo_owner(xin(i)-nr)
end do
do i=sz+1,size(xout)
xout(i) = -1
end do
end subroutine base_fnd_halo_owner_v
end module psb_indx_map_mod

Loading…
Cancel
Save