You cannot select more than 25 topics
Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.
334 lines
9.6 KiB
Fortran
334 lines
9.6 KiB
Fortran
!
|
|
!
|
|
! MLD2P4 version 2.2
|
|
! MultiLevel Domain Decomposition Parallel Preconditioners Package
|
|
! based on PSBLAS (Parallel Sparse BLAS version 3.5)
|
|
!
|
|
! (C) Copyright 2008-2018
|
|
!
|
|
! Salvatore Filippone
|
|
! Pasqua D'Ambra
|
|
! Daniela di Serafino
|
|
!
|
|
! Redistribution and use in source and binary forms, with or without
|
|
! modification, are permitted provided that the following conditions
|
|
! are met:
|
|
! 1. Redistributions of source code must retain the above copyright
|
|
! notice, this list of conditions and the following disclaimer.
|
|
! 2. Redistributions in binary form must reproduce the above copyright
|
|
! notice, this list of conditions, and the following disclaimer in the
|
|
! documentation and/or other materials provided with the distribution.
|
|
! 3. The name of the MLD2P4 group or the names of its contributors may
|
|
! not be used to endorse or promote products derived from this
|
|
! software without specific written permission.
|
|
!
|
|
! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
|
|
! ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
|
|
! TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
|
|
! PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE MLD2P4 GROUP OR ITS CONTRIBUTORS
|
|
! BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
|
|
! CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
|
|
! SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
|
|
! INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
|
|
! CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
|
|
! ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
|
|
! POSSIBILITY OF SUCH DAMAGE.
|
|
!
|
|
!
|
|
!
|
|
! File: mld_s_soc2_map__bld.f90
|
|
!
|
|
! Subroutine: mld_s_soc2_map_bld
|
|
! Version: real
|
|
!
|
|
! The aggregator object hosts the aggregation method for building
|
|
! the multilevel hierarchy. This variant is based on the method
|
|
! presented in
|
|
!
|
|
! S. Gratton, P. Henon, P. Jiranek and X. Vasseur:
|
|
! Reducing complexity of algebraic multigrid by aggregation
|
|
! Numerical Lin. Algebra with Applications, 2016, 23:501-518
|
|
!
|
|
! Note: upon exit
|
|
!
|
|
! Arguments:
|
|
! a - type(psb_sspmat_type).
|
|
! The sparse matrix structure containing the local part of the
|
|
! matrix to be preconditioned.
|
|
! desc_a - type(psb_desc_type), input.
|
|
! The communication descriptor of a.
|
|
! p - type(mld_sprec_type), input/output.
|
|
! The preconditioner data structure; upon exit it contains
|
|
! the multilevel hierarchy of prolongators, restrictors
|
|
! and coarse matrices.
|
|
! info - integer, output.
|
|
! Error code.
|
|
!
|
|
!
|
|
!
|
|
subroutine mld_s_soc2_map_bld(iorder,theta,a,desc_a,nlaggr,ilaggr,info)
|
|
|
|
use psb_base_mod
|
|
use mld_base_prec_type
|
|
use mld_s_inner_mod!, mld_protect_name => mld_s_soc2_map_bld
|
|
|
|
implicit none
|
|
|
|
! Arguments
|
|
integer(psb_ipk_), intent(in) :: iorder
|
|
type(psb_sspmat_type), intent(in) :: a
|
|
type(psb_desc_type), intent(in) :: desc_a
|
|
real(psb_spk_), intent(in) :: theta
|
|
integer(psb_lpk_), allocatable, intent(out) :: ilaggr(:),nlaggr(:)
|
|
integer(psb_ipk_), intent(out) :: info
|
|
|
|
! Local variables
|
|
integer(psb_ipk_), allocatable :: ils(:), neigh(:), irow(:), icol(:),&
|
|
& ideg(:), idxs(:)
|
|
integer(psb_lpk_), allocatable :: tmpaggr(:)
|
|
real(psb_spk_), allocatable :: val(:), diag(:)
|
|
integer(psb_ipk_) :: icnt,nlp,k,n,ia,isz,nr,nc,naggr,i,j,m, nz, ilg, ii, ip, ip1,nzcnt
|
|
type(psb_s_csr_sparse_mat) :: acsr, muij, s_neigh
|
|
type(psb_s_coo_sparse_mat) :: s_neigh_coo
|
|
real(psb_spk_) :: cpling, tcl
|
|
logical :: disjoint
|
|
integer(psb_ipk_) :: debug_level, debug_unit,err_act
|
|
integer(psb_ipk_) :: ictxt,np,me
|
|
integer(psb_ipk_) :: nrow, ncol, n_ne
|
|
character(len=20) :: name, ch_err
|
|
|
|
info=psb_success_
|
|
name = 'mld_soc2_map_bld'
|
|
call psb_erractionsave(err_act)
|
|
if (psb_errstatus_fatal()) then
|
|
info = psb_err_internal_error_; goto 9999
|
|
end if
|
|
debug_unit = psb_get_debug_unit()
|
|
debug_level = psb_get_debug_level()
|
|
!
|
|
ictxt=desc_a%get_context()
|
|
call psb_info(ictxt,me,np)
|
|
nrow = desc_a%get_local_rows()
|
|
ncol = desc_a%get_local_cols()
|
|
|
|
nr = a%get_nrows()
|
|
nc = a%get_ncols()
|
|
allocate(ilaggr(nr),neigh(nr),ideg(nr),idxs(nr),icol(nr),stat=info)
|
|
if(info /= psb_success_) then
|
|
info=psb_err_alloc_request_
|
|
call psb_errpush(info,name,i_err=(/2*nr,izero,izero,izero,izero/),&
|
|
& a_err='integer')
|
|
goto 9999
|
|
end if
|
|
|
|
diag = a%get_diag(info)
|
|
if(info /= psb_success_) then
|
|
info=psb_err_from_subroutine_
|
|
call psb_errpush(info,name,a_err='psb_sp_getdiag')
|
|
goto 9999
|
|
end if
|
|
|
|
!
|
|
! Phase zero: compute muij
|
|
!
|
|
call a%cp_to(muij)
|
|
do i=1, nr
|
|
do k=muij%irp(i),muij%irp(i+1)-1
|
|
j = muij%ja(k)
|
|
if (j<= nr) muij%val(k) = abs(muij%val(k))/sqrt(abs(diag(i)*diag(j)))
|
|
end do
|
|
end do
|
|
!write(*,*) 'murows/cols ',maxval(mu_rows),maxval(mu_cols)
|
|
!
|
|
! Compute the 1-neigbour; mark strong links with +1, weak links with -1
|
|
!
|
|
call s_neigh_coo%allocate(nr,nr,muij%get_nzeros())
|
|
ip = 0
|
|
do i=1, nr
|
|
do k=muij%irp(i),muij%irp(i+1)-1
|
|
j = muij%ja(k)
|
|
if (j<=nr) then
|
|
ip = ip + 1
|
|
s_neigh_coo%ia(ip) = i
|
|
s_neigh_coo%ja(ip) = j
|
|
if (real(muij%val(k)) >= theta) then
|
|
s_neigh_coo%val(ip) = sone
|
|
else
|
|
s_neigh_coo%val(ip) = -sone
|
|
end if
|
|
end if
|
|
end do
|
|
end do
|
|
!write(*,*) 'S_NEIGH: ',nr,ip
|
|
call s_neigh_coo%set_nzeros(ip)
|
|
call s_neigh%mv_from_coo(s_neigh_coo,info)
|
|
|
|
if (iorder == mld_aggr_ord_nat_) then
|
|
do i=1, nr
|
|
ilaggr(i) = -(nr+1)
|
|
idxs(i) = i
|
|
end do
|
|
else
|
|
do i=1, nr
|
|
ilaggr(i) = -(nr+1)
|
|
ideg(i) = muij%irp(i+1) - muij%irp(i)
|
|
end do
|
|
call psb_msort(ideg,ix=idxs,dir=psb_sort_down_)
|
|
end if
|
|
|
|
|
|
!
|
|
! Phase one: Start with disjoint groups.
|
|
!
|
|
naggr = 0
|
|
icnt = 0
|
|
step1: do ii=1, nr
|
|
i = idxs(ii)
|
|
|
|
if (ilaggr(i) == -(nr+1)) then
|
|
!
|
|
! Get the 1-neighbourhood of I
|
|
!
|
|
ip1 = s_neigh%irp(i)
|
|
nz = s_neigh%irp(i+1)-ip1
|
|
!
|
|
! If the neighbourhood only contains I, skip it
|
|
!
|
|
if (nz ==0) then
|
|
ilaggr(i) = 0
|
|
cycle step1
|
|
end if
|
|
if ((nz==1).and.(s_neigh%ja(ip1)==i)) then
|
|
ilaggr(i) = 0
|
|
cycle step1
|
|
end if
|
|
!
|
|
! If the whole strongly coupled neighborhood of I is
|
|
! as yet unconnected, turn it into the next aggregate.
|
|
!
|
|
nzcnt = count(real(s_neigh%val(ip1:ip1+nz-1)) > 0)
|
|
icol(1:nzcnt) = pack(s_neigh%ja(ip1:ip1+nz-1),(real(s_neigh%val(ip1:ip1+nz-1)) > 0))
|
|
disjoint = all(ilaggr(icol(1:nzcnt)) == -(nr+1))
|
|
if (disjoint) then
|
|
icnt = icnt + 1
|
|
naggr = naggr + 1
|
|
do k=1, nzcnt
|
|
ilaggr(icol(k)) = naggr
|
|
end do
|
|
ilaggr(i) = naggr
|
|
end if
|
|
endif
|
|
enddo step1
|
|
|
|
if (debug_level >= psb_debug_outer_) then
|
|
write(debug_unit,*) me,' ',trim(name),&
|
|
& ' Check 1:',count(ilaggr == -(nr+1))
|
|
end if
|
|
|
|
!
|
|
! Phase two: join the neighbours
|
|
!
|
|
tmpaggr = ilaggr
|
|
step2: do ii=1,nr
|
|
i = idxs(ii)
|
|
|
|
if (ilaggr(i) == -(nr+1)) then
|
|
!
|
|
! Find the most strongly connected neighbour that is
|
|
! already aggregated, if any, and join its aggregate
|
|
!
|
|
cpling = szero
|
|
ip = 0
|
|
do k=s_neigh%irp(i), s_neigh%irp(i+1)-1
|
|
j = s_neigh%ja(k)
|
|
if ((1<=j).and.(j<=nr)) then
|
|
if ( (tmpaggr(j) > 0).and. (real(muij%val(k)) > cpling)&
|
|
& .and.(real(s_neigh%val(k))>0)) then
|
|
ip = k
|
|
cpling = muij%val(k)
|
|
end if
|
|
end if
|
|
enddo
|
|
if (ip > 0) then
|
|
ilaggr(i) = ilaggr(s_neigh%ja(ip))
|
|
end if
|
|
end if
|
|
end do step2
|
|
|
|
|
|
!
|
|
! Phase three: sweep over leftovers, if any
|
|
!
|
|
step3: do ii=1,nr
|
|
i = idxs(ii)
|
|
|
|
if (ilaggr(i) < 0) then
|
|
!
|
|
! Find its strongly connected neighbourhood not
|
|
! already aggregated, and make it into a new aggregate.
|
|
!
|
|
ip = 0
|
|
do k=s_neigh%irp(i), s_neigh%irp(i+1)-1
|
|
j = s_neigh%ja(k)
|
|
if ((1<=j).and.(j<=nr)) then
|
|
if (ilaggr(j) < 0) then
|
|
ip = ip + 1
|
|
icol(ip) = j
|
|
end if
|
|
end if
|
|
enddo
|
|
if (ip > 0) then
|
|
icnt = icnt + 1
|
|
naggr = naggr + 1
|
|
ilaggr(i) = naggr
|
|
do k=1, ip
|
|
ilaggr(icol(k)) = naggr
|
|
end do
|
|
end if
|
|
end if
|
|
end do step3
|
|
|
|
|
|
if (count(ilaggr<0) >0) then
|
|
info=psb_err_internal_error_
|
|
call psb_errpush(info,name,a_err='Fatal error: some leftovers')
|
|
goto 9999
|
|
endif
|
|
|
|
if (naggr > ncol) then
|
|
write(0,*) name,'Error : naggr > ncol',naggr,ncol
|
|
info=psb_err_internal_error_
|
|
call psb_errpush(info,name,a_err='Fatal error: naggr>ncol')
|
|
goto 9999
|
|
end if
|
|
|
|
call psb_realloc(ncol,ilaggr,info)
|
|
if (info /= psb_success_) then
|
|
info=psb_err_from_subroutine_
|
|
ch_err='psb_realloc'
|
|
call psb_errpush(info,name,a_err=ch_err)
|
|
goto 9999
|
|
end if
|
|
|
|
allocate(nlaggr(np),stat=info)
|
|
if (info /= psb_success_) then
|
|
info=psb_err_alloc_request_
|
|
call psb_errpush(info,name,i_err=(/np,izero,izero,izero,izero/),&
|
|
& a_err='integer')
|
|
goto 9999
|
|
end if
|
|
|
|
nlaggr(:) = 0
|
|
nlaggr(me+1) = naggr
|
|
call psb_sum(ictxt,nlaggr(1:np))
|
|
|
|
call psb_erractionrestore(err_act)
|
|
return
|
|
|
|
9999 call psb_error_handler(err_act)
|
|
|
|
return
|
|
|
|
end subroutine mld_s_soc2_map_bld
|
|
|