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.
388 lines
12 KiB
Fortran
388 lines
12 KiB
Fortran
17 years ago
|
!!$
|
||
18 years ago
|
!!$
|
||
15 years ago
|
!!$ MLD2P4 version 2.0
|
||
17 years ago
|
!!$ MultiLevel Domain Decomposition Parallel Preconditioners Package
|
||
15 years ago
|
!!$ based on PSBLAS (Parallel Sparse BLAS version 3.0)
|
||
17 years ago
|
!!$
|
||
15 years ago
|
!!$ (C) Copyright 2008,2009,2010
|
||
17 years ago
|
!!$
|
||
16 years ago
|
!!$ Salvatore Filippone University of Rome Tor Vergata
|
||
15 years ago
|
!!$ Alfredo Buttari CNRS-IRIT, Toulouse
|
||
17 years ago
|
!!$ Pasqua D'Ambra ICAR-CNR, Naples
|
||
|
!!$ Daniela di Serafino Second University of Naples
|
||
18 years ago
|
!!$
|
||
|
!!$ 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.
|
||
17 years ago
|
!!$ 3. The name of the MLD2P4 group or the names of its contributors may
|
||
18 years ago
|
!!$ 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
|
||
17 years ago
|
!!$ PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE MLD2P4 GROUP OR ITS CONTRIBUTORS
|
||
18 years ago
|
!!$ 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.
|
||
|
!!$
|
||
17 years ago
|
!!$
|
||
17 years ago
|
! File: mld_zsp_renum.f90
|
||
17 years ago
|
!
|
||
17 years ago
|
! Subroutine: mld_zsp_renum
|
||
17 years ago
|
! Version: complex
|
||
17 years ago
|
! Contains: gps_reduction
|
||
17 years ago
|
!
|
||
|
! This routine reorders the rows and the columns of the local part of a sparse
|
||
|
! distributed matrix, according to one of the following criteria:
|
||
|
! 1. the numbering of the global column indices,
|
||
|
! 2. the Gibbs-Poole-Stockmeyer (GPS) band reduction algorithm.
|
||
|
! NOTE: the GPS algorithm is disabled for the time being (see mld_prec_type.f90).
|
||
|
!
|
||
|
! The matrix to be reordered is stored into a and blck, as specified in the
|
||
|
! description of the arguments below.
|
||
|
!
|
||
|
! If required by the user (p%iprcparm(mld_sub_ren_) /= 0), the routine is
|
||
17 years ago
|
! used by mld_fact_bld in building the block-Jacobi and Additive Schwarz
|
||
17 years ago
|
! 'base preconditioners' corresponding to any level of a multilevel
|
||
|
! preconditioner.
|
||
|
!
|
||
|
!
|
||
|
! Arguments:
|
||
17 years ago
|
! a - type(psb_zspmat_type), input.
|
||
17 years ago
|
! The sparse matrix structure containing the 'original' local
|
||
|
! part of the matrix to be reordered, i.e. the rows of the matrix
|
||
|
! held by the calling process according to the initial data
|
||
|
! distribution.
|
||
17 years ago
|
! blck - type(psb_zspmat_type), input.
|
||
17 years ago
|
! The sparse matrix structure containing the remote rows of the
|
||
17 years ago
|
! matrix to be reordered, that have been retrieved by mld_as_bld
|
||
17 years ago
|
! to build an Additive Schwarz base preconditioner with overlap
|
||
|
! greater than 0.If the overlap is 0, then blck does not contain
|
||
|
! any row.
|
||
16 years ago
|
! p - type(mld_zbaseprec_type), input/output.
|
||
17 years ago
|
! The base preconditioner data structure containing the local
|
||
|
! part of the base preconditioner to be built. In input it
|
||
|
! contains information on the type of reordering to be applied
|
||
|
! and on the matrix to be reordered. In output it contains
|
||
|
! information on the reordering applied.
|
||
17 years ago
|
! atmp - type(psb_zspmat_type), output.
|
||
17 years ago
|
! The sparse matrix structure containing the whole local reordered
|
||
|
! matrix.
|
||
13 years ago
|
! info - integer, output.
|
||
17 years ago
|
! Error code.
|
||
|
!
|
||
17 years ago
|
subroutine mld_zsp_renum(a,blck,p,atmp,info)
|
||
17 years ago
|
|
||
14 years ago
|
use psb_base_mod
|
||
14 years ago
|
use mld_z_inner_mod, mld_protect_name => mld_zsp_renum
|
||
18 years ago
|
|
||
18 years ago
|
implicit none
|
||
|
|
||
17 years ago
|
! Arguments
|
||
18 years ago
|
type(psb_zspmat_type), intent(in) :: a,blck
|
||
17 years ago
|
type(psb_zspmat_type), intent(out) :: atmp
|
||
16 years ago
|
type(mld_zbaseprec_type), intent(inout) :: p
|
||
18 years ago
|
integer, intent(out) :: info
|
||
|
|
||
17 years ago
|
! Local variables
|
||
17 years ago
|
character(len=20) :: name, ch_err
|
||
13 years ago
|
integer :: nztota, nztotb, nztmp, nzt2, nnr, i,k, ma, mb
|
||
17 years ago
|
integer, allocatable :: itmp(:), itmp2(:)
|
||
|
integer :: ictxt,np,me, err_act
|
||
13 years ago
|
type(psb_z_coo_sparse_mat) :: cootmp, cootmp2
|
||
|
type(psb_z_csr_sparse_mat) :: csrtmp
|
||
17 years ago
|
real(psb_dpk_) :: t3,t4
|
||
18 years ago
|
|
||
|
if (psb_get_errstatus().ne.0) return
|
||
15 years ago
|
info=psb_success_
|
||
18 years ago
|
name='mld_zsp_renum'
|
||
18 years ago
|
call psb_erractionsave(err_act)
|
||
|
|
||
14 years ago
|
ictxt=p%desc_data%get_context()
|
||
18 years ago
|
call psb_info(ictxt, me, np)
|
||
|
|
||
18 years ago
|
!
|
||
17 years ago
|
! NOTE: the matrix to be reordered is converted into the COO format.
|
||
|
! If necessary it is converted from the COO to the CSR format.
|
||
|
! The output matrix is in COO format.
|
||
|
!
|
||
|
|
||
|
!
|
||
|
! Convert a into the COO format and extend it up to a%m+blck%m rows
|
||
|
! by adding null rows. The converted extended matrix is stored in atmp.
|
||
18 years ago
|
!
|
||
13 years ago
|
nztota=a%get_nzeros()
|
||
|
nztotb=blck%get_nzeros()
|
||
|
ma = a%get_nrows()
|
||
|
mb = blck%get_nrows()
|
||
|
|
||
|
|
||
15 years ago
|
if (p%iprcparm(mld_sub_ren_) == mld_renum_glb_) then
|
||
18 years ago
|
|
||
17 years ago
|
!
|
||
|
! Remember: we have switched IA1=COLS and IA2=ROWS.
|
||
|
! Now identify the set of distinct local column indices.
|
||
|
!
|
||
14 years ago
|
nnr = p%desc_data%get_local_rows()
|
||
18 years ago
|
allocate(p%perm(nnr),p%invperm(nnr),itmp2(nnr),stat=info)
|
||
15 years ago
|
if (info /= psb_success_) then
|
||
|
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
|
||
18 years ago
|
goto 9999
|
||
|
end if
|
||
17 years ago
|
do i=1, nnr
|
||
|
itmp2(i) = i
|
||
|
end do
|
||
17 years ago
|
call psb_loc_to_glob(itmp2(1:nnr),p%desc_data,info,iact='I')
|
||
18 years ago
|
!
|
||
17 years ago
|
! Compute reordering. We want new(i) = old(perm(i)).
|
||
18 years ago
|
!
|
||
18 years ago
|
call psb_msort(itmp2(1:nnr),ix=p%perm)
|
||
17 years ago
|
!
|
||
|
! Compute the inverse of the permutation stored in perm
|
||
|
!
|
||
18 years ago
|
do k=1, nnr
|
||
|
p%invperm(p%perm(k)) = k
|
||
|
enddo
|
||
18 years ago
|
t3 = psb_wtime()
|
||
18 years ago
|
|
||
15 years ago
|
else if (p%iprcparm(mld_sub_ren_) == mld_renum_gps_) then
|
||
17 years ago
|
|
||
|
!
|
||
18 years ago
|
! This is a renumbering with Gibbs-Poole-Stockmeyer
|
||
|
! band reduction. Switched off for now. To be fixed,
|
||
17 years ago
|
! gps_reduction should get p%perm.
|
||
|
!
|
||
|
|
||
|
!
|
||
|
! Convert atmp into the CSR format
|
||
|
!
|
||
13 years ago
|
call a%cscnv(atmp,info,type='coo',dupl=psb_dupl_add_)
|
||
|
call psb_rwextd(ma+mb,atmp,info,blck)
|
||
|
call atmp%mv_to(csrtmp)
|
||
|
|
||
|
nztmp = csrtmp%get_nzeros()
|
||
18 years ago
|
|
||
17 years ago
|
!
|
||
|
! Realloc the permutation arrays
|
||
|
!
|
||
13 years ago
|
call psb_realloc(csrtmp%get_nrows(),p%perm,info)
|
||
15 years ago
|
if(info /= psb_success_) then
|
||
|
info=psb_err_from_subroutine_
|
||
18 years ago
|
ch_err='psb_realloc'
|
||
|
call psb_errpush(info,name,a_err=ch_err)
|
||
|
goto 9999
|
||
|
end if
|
||
|
|
||
13 years ago
|
call psb_realloc(csrtmp%get_nrows(),p%invperm,info)
|
||
15 years ago
|
if(info /= psb_success_) then
|
||
|
info=psb_err_from_subroutine_
|
||
18 years ago
|
ch_err='psb_realloc'
|
||
|
call psb_errpush(info,name,a_err=ch_err)
|
||
|
goto 9999
|
||
|
end if
|
||
|
|
||
13 years ago
|
allocate(itmp(max(8,csrtmp%get_nrows()+2,nztmp+2)),stat=info)
|
||
15 years ago
|
if (info /= psb_success_) then
|
||
|
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
|
||
18 years ago
|
goto 9999
|
||
|
end if
|
||
|
|
||
|
itmp(1:8) = 0
|
||
|
|
||
17 years ago
|
!
|
||
|
! Renumber rows and columns according to the GPS algorithm
|
||
|
!
|
||
13 years ago
|
call gps_reduction(csrtmp%get_nrows(),csrtmp%irp,csrtmp%ja,p%perm,p%invperm,info)
|
||
15 years ago
|
if(info /= psb_success_) then
|
||
|
info=psb_err_from_subroutine_
|
||
18 years ago
|
ch_err='gps_reduction'
|
||
|
call psb_errpush(info,name,a_err=ch_err)
|
||
|
goto 9999
|
||
|
end if
|
||
|
|
||
17 years ago
|
!
|
||
|
! Compute the inverse permutation
|
||
|
!
|
||
13 years ago
|
do k=1, csrtmp%get_nrows()
|
||
18 years ago
|
p%invperm(p%perm(k)) = k
|
||
|
enddo
|
||
18 years ago
|
t3 = psb_wtime()
|
||
17 years ago
|
|
||
18 years ago
|
end if
|
||
|
|
||
17 years ago
|
!
|
||
|
! Rebuild atmp with the new numbering (COO format)
|
||
|
!
|
||
13 years ago
|
call a%cp_to(cootmp)
|
||
|
nztmp=cootmp%get_nzeros()
|
||
18 years ago
|
do i=1,nztmp
|
||
13 years ago
|
cootmp%ia(i) = p%perm(cootmp%ia(i))
|
||
|
cootmp%ja(i) = p%invperm(cootmp%ja(i))
|
||
|
end do
|
||
|
call blck%cp_to(cootmp2)
|
||
|
nzt2 = cootmp2%get_nzeros()
|
||
|
call psb_ensure_size(nztmp+nzt2,cootmp%ia,info)
|
||
|
call psb_ensure_size(nztmp+nzt2,cootmp%ja,info)
|
||
|
call psb_ensure_size(nztmp+nzt2,cootmp%val,info)
|
||
|
do i=1,nzt2
|
||
|
cootmp%ia(nztmp+i) = p%perm(cootmp2%ia(i))
|
||
|
cootmp%ja(nztmp+i) = p%invperm(cootmp2%ja(i))
|
||
|
cootmp%val(nztmp+i) = (cootmp2%val(i))
|
||
18 years ago
|
end do
|
||
13 years ago
|
call cootmp2%free()
|
||
|
call cootmp%set_nzeros(nztmp+nzt2)
|
||
|
call cootmp%set_dupl(psb_dupl_add_)
|
||
|
call cootmp%fix(info)
|
||
|
call atmp%mv_from(cootmp)
|
||
15 years ago
|
if (info /= psb_success_) then
|
||
|
call psb_errpush(psb_err_from_subroutine_,name,a_err='psb_fixcoo')
|
||
18 years ago
|
goto 9999
|
||
|
end if
|
||
17 years ago
|
|
||
18 years ago
|
t4 = psb_wtime()
|
||
17 years ago
|
|
||
18 years ago
|
call psb_erractionrestore(err_act)
|
||
|
return
|
||
|
|
||
|
9999 continue
|
||
|
call psb_erractionrestore(err_act)
|
||
18 years ago
|
if (err_act.eq.psb_act_abort_) then
|
||
18 years ago
|
call psb_error()
|
||
|
return
|
||
|
end if
|
||
|
return
|
||
|
|
||
|
contains
|
||
|
|
||
17 years ago
|
!
|
||
17 years ago
|
! Subroutine: gps_reduction
|
||
|
! Note: internal subroutine of mld_zsp_renum
|
||
17 years ago
|
!
|
||
|
! Compute a renumbering of the row and column indices of a sparse matrix
|
||
13 years ago
|
! according to the Gibbs-Poole-Stockmeyer band reduction algorithm. The
|
||
17 years ago
|
! matrix is stored in CSR format.
|
||
|
!
|
||
|
! This routine has been obtained by adapting ACM TOMS Algorithm 582.
|
||
|
!
|
||
|
!
|
||
|
! Arguments:
|
||
|
! m - integer, ...
|
||
17 years ago
|
! The number of rows of the matrix to which the renumbering
|
||
17 years ago
|
! is applied.
|
||
|
! ia - integer, dimension(:), ...
|
||
|
! The indices identifying the first nonzero entry of each row
|
||
|
! of the matrix, according to the CSR storage format.
|
||
|
! ja - integer, dimension(:), ...
|
||
|
! The column indices of the nonzero entries of the matrix,
|
||
|
! according to the CSR storage format.
|
||
|
! perm - integer, dimension(:), ...
|
||
|
! The row/column index permutation corresponding to the
|
||
|
! renumbering.
|
||
|
! iperm - integer, dimension(:),...
|
||
|
! The inverse of the row/column permutation stored in perm.
|
||
|
! info - integer, output.
|
||
|
! Error code
|
||
|
!
|
||
18 years ago
|
subroutine gps_reduction(m,ia,ja,perm,iperm,info)
|
||
17 years ago
|
|
||
17 years ago
|
! Arguments
|
||
|
integer :: m
|
||
|
integer,dimension(:) :: ia,ja,perm,iperm
|
||
18 years ago
|
integer, intent(out) :: info
|
||
|
|
||
17 years ago
|
! Local variables
|
||
|
integer :: i,j,dgConn,Npnt
|
||
|
integer :: n,idpth,ideg,ibw2,ipf2
|
||
18 years ago
|
integer,dimension(:,:),allocatable::NDstk
|
||
|
integer,dimension(:),allocatable::iOld,renum,ndeg,lvl,lvls1,lvls2,ccstor
|
||
17 years ago
|
character(len=20) :: name
|
||
18 years ago
|
|
||
|
if(psb_get_errstatus().ne.0) return
|
||
15 years ago
|
info=psb_success_
|
||
18 years ago
|
name='gps_reduction'
|
||
|
call psb_erractionsave(err_act)
|
||
|
|
||
17 years ago
|
! Compute the maximum connectivity degree
|
||
18 years ago
|
npnt = m
|
||
|
dgConn=0
|
||
|
do i=1,m
|
||
|
dgconn = max(dgconn,(ia(i+1)-ia(i)))
|
||
|
enddo
|
||
17 years ago
|
! The maximum connectivity value is dgConn
|
||
18 years ago
|
|
||
17 years ago
|
n=Npnt ! Max number of rows
|
||
|
iDeg=dgConn ! Max connectivity
|
||
|
! iDpth= ! Number of level (initialization not needed)
|
||
18 years ago
|
|
||
|
allocate(NDstk(Npnt,dgConn),stat=info)
|
||
15 years ago
|
if (info /= psb_success_) then
|
||
|
info=psb_err_alloc_dealloc_
|
||
18 years ago
|
call psb_errpush(info,name)
|
||
|
goto 9999
|
||
|
endif
|
||
|
allocate(iOld(Npnt),renum(Npnt+1),ndeg(Npnt),lvl(Npnt),lvls1(Npnt),&
|
||
|
&lvls2(Npnt),ccstor(Npnt),stat=info)
|
||
15 years ago
|
if (info /= psb_success_) then
|
||
|
info=psb_err_alloc_dealloc_
|
||
18 years ago
|
call psb_errpush(info,name)
|
||
|
goto 9999
|
||
|
endif
|
||
|
|
||
17 years ago
|
! Prepare the matrix graph
|
||
18 years ago
|
Ndstk(:,:)=0
|
||
|
do i=1,Npnt
|
||
|
k=0
|
||
|
do j = ia(i),ia(i+1) - 1
|
||
|
if ((1<=ja(j)).and.( ja( j ) /= i ).and.(ja(j)<=npnt)) then
|
||
|
k = k+1
|
||
|
Ndstk(i,k)=ja(j)
|
||
|
endif
|
||
|
enddo
|
||
|
ndeg(i)=k
|
||
|
enddo
|
||
|
|
||
17 years ago
|
! Numbering
|
||
18 years ago
|
do i=1,Npnt
|
||
|
iOld(i)=i
|
||
|
enddo
|
||
17 years ago
|
|
||
|
! Call gps_reduce
|
||
18 years ago
|
call psb_gps_reduce(Ndstk,Npnt,iOld,renum,ndeg,lvl,lvls1, lvls2,ccstor,&
|
||
|
& ibw2,ipf2,n,idpth,ideg)
|
||
17 years ago
|
|
||
|
! Build permutation vector
|
||
18 years ago
|
perm(1:Npnt)=renum(1:Npnt)
|
||
17 years ago
|
|
||
|
!Build inverse permutation vector
|
||
18 years ago
|
do i=1,Npnt
|
||
|
iperm(perm(i))=i
|
||
|
enddo
|
||
17 years ago
|
|
||
|
! Deallocate memory
|
||
18 years ago
|
deallocate(NDstk,iOld,renum,ndeg,lvl,lvls1,lvls2,ccstor)
|
||
|
|
||
|
call psb_erractionrestore(err_act)
|
||
|
return
|
||
|
|
||
|
9999 continue
|
||
|
call psb_erractionrestore(err_act)
|
||
18 years ago
|
if (err_act.eq.psb_act_abort_) then
|
||
18 years ago
|
call psb_error()
|
||
|
return
|
||
|
end if
|
||
|
return
|
||
|
|
||
|
end subroutine gps_reduction
|
||
|
|
||
18 years ago
|
end subroutine mld_zsp_renum
|