Merged internal comments for FND_OWNER methods.

psblas-3.6-maint
Salvatore Filippone 5 years ago
parent 5d6a380664
commit 3368b44174

@ -30,21 +30,28 @@
!
!
!
! File: psi_fnd_owner.f90
! File: psb_fnd_owner.f90
!
! Subroutine: psi_fnd_owner
! Subroutine: psb_fnd_owner
! Figure out who owns global indices.
!
! Arguments:
! nv - integer Number of indices required on the calling
! process
! idx(:) - integer Required indices on the calling process.
! Note: the indices should be unique!
! iprc(:) - integer(psb_ipk_), allocatable Output: process identifiers for the corresponding
! indices
! desc_a - type(psb_desc_type). The communication descriptor.
! idxmap - class(psb_indx_map). The index map
! info - integer. return code.
!
! This is the default implementation of the FND_OWNER method.
! If a particular index map class has additional information, it can override it
! (see e.g. the GEN_BLOCK_MAP class).
!
! 1. Check if IDXM%PARTS is available, and use it; or
! 2. Check if TEMPVG(:) is allocated, and use it; or
! 3. Call the general method PSI_GRAPH_FND_OWNER.
!
!
subroutine psb_indx_map_fnd_owner(idx,iprc,idxmap,info)
use psb_serial_mod
use psb_const_mod

@ -29,22 +29,24 @@
! POSSIBILITY OF SUCH DAMAGE.
!
!
! File: psi_a2a_fnd_owner.f90
!
! File: psi_fnd_owner.f90
!
! Subroutine: psi_fnd_owner
! Subroutine: psi_a2a_fnd_owner
! Figure out who owns global indices.
!
! Arguments:
! nv - integer Number of indices required on the calling
! process
! idx(:) - integer Required indices on the calling process.
! Note: the indices should be unique!
! iprc(:) - integer(psb_ipk_), allocatable Output: process identifiers for the corresponding
! indices
! desc_a - type(psb_desc_type). The communication descriptor.
! idxmap - class(psb_indx_map). The index map
! info - integer. return code.
!
! This version does not assume any prior knowledge about the process topology,
! so it goes for an all-to-all by building an auxiliary neighbours list and
! reusing the neighbour version.
!
!
subroutine psi_a2a_fnd_owner(idx,iprc,idxmap,info)
use psb_serial_mod
use psb_const_mod

@ -30,12 +30,28 @@
!
!
!
! File: psi_fnd_owner.f90
!
! Subroutine: psi_fnd_owner
! File: psi_adjcncy_fnd_owner.f90
!
! Subroutine: psi_adjcncy_fnd_owner
! Figure out who owns global indices.
!
! Arguments:
! idx(:) - integer Required indices on the calling process.
! Note: the indices should be unique!
! iprc(:) - integer(psb_ipk_), allocatable Output: process identifiers for the corresponding
! indices
! adj(:) - integer(psb_ipk_) Input: list of topological neighbours for current process.
!
! idxmap - class(psb_indx_map). The index map
! info - integer. return code.
!
! This version takes on input a list of processes that are assumed to
! be topological neighbours of the current one. Each process will send to all
! of its neighbours the list of indices for which it is trying to find the
! owner, prepare its own answers, and collect answers from others.
! There are three possibile implementations: using mpi_alltoallv, using mpi_isend/irecv,
! using psb_snd/psb_rcv. The default is mpi_alltoallv.
!
subroutine psi_adjcncy_fnd_owner(idx,iprc,adj,idxmap,info)
use psb_serial_mod
@ -69,8 +85,10 @@ subroutine psi_adjcncy_fnd_owner(idx,iprc,adj,idxmap,info)
& last_ih, last_j, nidx, nrecv, nadj
integer(psb_ipk_) :: mglob, ih
integer(psb_ipk_) :: ictxt,np,me
logical, parameter :: gettime=.false., new_impl=.true.
logical, parameter :: a2av_impl=.true., debug=.false.
logical, parameter :: gettime=.false., debug=.false.
logical, parameter :: a2av_impl=.true.
logical, parameter :: mpi_irecv_impl=.false.
logical, parameter :: psb_rcv_impl=.false.
real(psb_dpk_) :: t0, t1, t2, t3, t4, tamx, tidx
character(len=20) :: name
@ -175,6 +193,10 @@ subroutine psi_adjcncy_fnd_owner(idx,iprc,adj,idxmap,info)
call mpi_alltoallv(tproc,rvsz,rvidx,psb_mpi_ipk_integer,&
& lclidx,sdsz,sdidx,psb_mpi_ipk_integer,icomm,iret)
!
! Because IPRC has been initialized to -1, the MAX operation selects
! the answers.
!
do i=0, np-1
if (sdsz(i)>0) then
! Must be nidx == sdsz(i)
@ -183,8 +205,7 @@ subroutine psi_adjcncy_fnd_owner(idx,iprc,adj,idxmap,info)
end do
if (debug) write(0,*) me,' End of adjcncy_fnd ',iprc(1:nidx)
else
if (new_impl) then
else if (mpi_irecv_impl) then
!
! First simple minded version with auxiliary arrays
! dimensioned on NP.
@ -303,7 +324,7 @@ subroutine psi_adjcncy_fnd_owner(idx,iprc,adj,idxmap,info)
end do
if (debug) write(0,*) me,' End of adjcncy_fnd ',iprc(1:nidx)
else
else if (psb_rcv_impl) then
Allocate(hidx(0:np),hsz(np),&
& sdsz(0:np-1),rvsz(0:np-1),stat=info)
@ -371,7 +392,10 @@ subroutine psi_adjcncy_fnd_owner(idx,iprc,adj,idxmap,info)
if (nidx > 0) call psb_rcv(ictxt,tproc(1:nidx),adj(j))
iprc(1:nidx) = max(iprc(1:nidx), tproc(1:nidx))
end do
end if
else
info = psb_err_internal_error_
call psb_errpush(info,name,a_err='invalid exchange alg choice')
goto 9999
end if
call psb_erractionrestore(err_act)

Loading…
Cancel
Save