!!$ !!$ Parallel Sparse BLAS version 2.2 !!$ (C) Copyright 2006/2007/2008 !!$ Salvatore Filippone University of Rome Tor Vergata !!$ Alfredo Buttari University of Rome Tor Vergata !!$ !!$ 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 PSBLAS 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 PSBLAS 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: psb_csymbmm.f90 ! Subroutine: ! Arguments: ! ! ! Note: This subroutine performs the symbolic product of two sparse matrices. ! It is modeled after the SMMP package by R. Bank and C. Douglas, but is ! rewritten in Fortran 95 making use of our sparse matrix facilities. ! subroutine psb_csymbmm(a,b,c,info) use psb_sparse_mod, psb_protect_name => psb_csymbmm implicit none type(psb_c_sparse_mat), intent(in) :: a,b type(psb_c_sparse_mat), intent(out) :: c integer, intent(out) :: info type(psb_c_csr_sparse_mat), allocatable :: ccsr integer :: err_act character(len=*), parameter :: name='psb_symbmm' call psb_erractionsave(err_act) info = 0 if ((a%is_null()) .or.(b%is_null())) then info = 1121 call psb_errpush(info,name) goto 9999 endif allocate(ccsr, stat=info) if (info /= 0) then info = 4000 call psb_errpush(info,name) goto 9999 end if call psb_symbmm(a%a,b%a,ccsr,info) if (info /= 0) then call psb_errpush(info,name) goto 9999 end if call move_alloc(ccsr,c%a) call psb_erractionrestore(err_act) return 9999 continue call psb_erractionrestore(err_act) if (err_act == psb_act_abort_) then call psb_error() return end if return end subroutine psb_csymbmm subroutine psb_cbase_symbmm(a,b,c,info) use psb_mat_mod use psb_serial_mod, psb_protect_name => psb_cbase_symbmm implicit none class(psb_c_base_sparse_mat), intent(in) :: a,b type(psb_c_csr_sparse_mat), intent(out) :: c integer, intent(out) :: info integer, allocatable :: itemp(:) integer :: nze, ma,na,mb,nb character(len=20) :: name integer :: err_act name='psb_symbmm' call psb_erractionsave(err_act) info = 0 ma = a%get_nrows() na = a%get_ncols() mb = b%get_nrows() nb = b%get_ncols() if ( mb /= na ) then write(0,*) 'Mismatch in SYMBMM: ',ma,na,mb,nb endif allocate(itemp(max(ma,na,mb,nb)),stat=info) if (info /= 0) then info = 4000 call psb_Errpush(info,name) goto 9999 endif ! ! Note: we need to test whether there is a performance impact ! in not using the original Douglas & Bank code. ! select type(a) type is (psb_c_csr_sparse_mat) select type(b) type is (psb_c_csr_sparse_mat) call csr_symbmm(a,b,c,itemp,info) class default call gen_symbmm(a,b,c,itemp,info) end select class default call gen_symbmm(a,b,c,itemp,info) end select if (info /= 0) then call psb_errpush(info,name) goto 9999 end if call psb_realloc(size(c%ja),c%val,info) deallocate(itemp) call psb_erractionrestore(err_act) return 9999 continue call psb_erractionrestore(err_act) if (err_act == psb_act_abort_) then call psb_error() return end if return contains subroutine csr_symbmm(a,b,c,itemp,info) type(psb_c_csr_sparse_mat), intent(in) :: a,b type(psb_c_csr_sparse_mat), intent(out) :: c integer :: itemp(:) integer, intent(out) :: info interface subroutine symbmm (n, m, l, ia, ja, diaga, & & ib, jb, diagb, ic, jc, diagc, index) integer n,m,l, ia(*), ja(*), diaga, ib(*), jb(*), diagb,& & diagc, index(*) integer, allocatable :: ic(:),jc(:) end subroutine symbmm end interface integer :: nze, ma,na,mb,nb info = 0 ma = a%get_nrows() na = a%get_ncols() mb = b%get_nrows() nb = b%get_ncols() nze = max(ma+1,2*ma) call c%allocate(ma,nb,nze) call symbmm(ma,na,nb,a%irp,a%ja,0,& & b%irp,b%ja,0,& & c%irp,c%ja,0,itemp) end subroutine csr_symbmm subroutine gen_symbmm(a,b,c,index,info) class(psb_c_base_sparse_mat), intent(in) :: a,b type(psb_c_csr_sparse_mat), intent(out) :: c integer :: index(:),info integer, allocatable :: iarw(:), iacl(:),ibrw(:),ibcl(:) integer :: maxlmn,i,j,m,n,k,l,istart,length,nazr,nbzr,jj,minlm,minmn integer :: nze, ma,na,mb,nb ma = a%get_nrows() na = a%get_ncols() mb = b%get_nrows() nb = b%get_ncols() nze = max(ma+1,2*ma) call c%allocate(ma,nb,nze) n = ma m = na l = nb maxlmn = max(l,m,n) allocate(iarw(maxlmn),iacl(maxlmn),ibrw(maxlmn),ibcl(maxlmn),& & stat=info) if (info /= 0) then info = 4000 return endif do i=1,maxlmn index(i)=0 end do c%irp(1)=1 minlm = min(l,m) minmn = min(m,n) main: do i=1,n istart=-1 length=0 call a%csget(i,i,nazr,iarw,iacl,info) do jj=1, nazr j=iacl(jj) if ((j<1).or.(j>m)) then write(0,*) ' SymbMM: Problem with A ',i,jj,j,m info = 1 return endif call b%csget(j,j,nbzr,ibrw,ibcl,info) do k=1,nbzr if ((ibcl(k)<1).or.(ibcl(k)>maxlmn)) then write(0,*) 'Problem in SYMBMM 1:',j,k,ibcl(k),maxlmn info=2 return else if(index(ibcl(k)) == 0) then index(ibcl(k))=istart istart=ibcl(k) length=length+1 endif endif end do end do c%irp(i+1)=c%irp(i)+length if (c%irp(i+1) > size(c%ja)) then if (n > (2*i)) then nze = max(c%irp(i+1), c%irp(i)*((n+i-1)/i)) else nze = max(c%irp(i+1), nint((dble(c%irp(i))*(dble(n)/i))) ) endif call psb_realloc(nze,c%ja,info) end if do j= c%irp(i),c%irp(i+1)-1 c%ja(j)=istart istart=index(istart) index(c%ja(j))=0 end do call psb_msort(c%ja(c%irp(i):c%irp(i)+length-1)) index(i) = 0 end do main end subroutine gen_symbmm end subroutine psb_cbase_symbmm