Given the sparsity pattern of an m
by n
matrix a
,
this subroutine determines the degree sequence for
the intersection graph of the columns of a
.
In graph-theory terminology, the intersection graph of
the columns of a
is the loopless graph g
with vertices
a(j), j = 1,2,...,n
where a(j)
is the j
-th column of a
and with edge (a(i),a(j))
if and only if columns i
and j
have a non-zero in the same row position.
Note
The value of m
is not needed by degr
and is
therefore not present in the subroutine statement.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n |
a positive integer input variable set to the number
of columns of |
||
integer, | intent(in), | dimension(*) | :: | indrow |
an integer input array which contains the row
indices for the non-zeroes in the matrix |
|
integer, | intent(in), | dimension(n+1) | :: | jpntr |
an integer input array of length |
|
integer, | intent(in), | dimension(*) | :: | indcol |
an integer input array which contains the
column indices for the non-zeroes in the matrix |
|
integer, | intent(in), | dimension(*) | :: | ipntr |
an integer input array of length |
|
integer, | intent(out), | dimension(n) | :: | ndeg |
an integer output array of length |
|
integer, | dimension(n) | :: | iwa |
an integer work array of length |
subroutine degr(n,Indrow,Jpntr,Indcol,Ipntr,Ndeg,Iwa) implicit none integer,intent(in) :: n !! a positive integer input variable set to the number !! of columns of `a`. integer,dimension(*),intent(in) :: indrow !! an integer input array which contains the row !! indices for the non-zeroes in the matrix `a`. integer,dimension(n+1),intent(in) :: jpntr !! an integer input array of length `n + 1` which !! specifies the locations of the row indices in `indrow`. !! the row indices for column `j` are !! `indrow(k), k = jpntr(j),...,jpntr(j+1)-1`. !! **note** that `jpntr(n+1)-1` is then the number of non-zero !! elements of the matrix `a`. integer,dimension(*),intent(in) :: indcol !! an integer input array which contains the !! column indices for the non-zeroes in the matrix `a`. integer,dimension(*),intent(in) :: ipntr !! an integer input array of length `m + 1` which !! specifies the locations of the column indices in `indcol`. !! the column indices for row `i` are !! `indcol(k), k = ipntr(i),...,ipntr(i+1)-1`. !! **note** that `ipntr(m+1)-1` is then the number of non-zero !! elements of the matrix `a`. integer,dimension(n),intent(out) :: ndeg !! an integer output array of length `n` which !! specifies the degree sequence. the degree of the !! `j`-th column of `a` is `ndeg(j)`. integer,dimension(n) :: iwa !! an integer work array of length `n` integer :: ic , ip , ir , jcol , jp ! initialization block. do jp = 1 , n Ndeg(jp) = 0 Iwa(jp) = 0 enddo ! compute the degree sequence by determining the contributions ! to the degrees from the current(jcol) column and further ! columns which have not yet been considered. do jcol = 2 , n Iwa(jcol) = n ! determine all positions (ir,jcol) which correspond ! to non-zeroes in the matrix. do jp = Jpntr(jcol) , Jpntr(jcol+1) - 1 ir = Indrow(jp) ! for each row ir, determine all positions (ir,ic) ! which correspond to non-zeroes in the matrix. do ip = Ipntr(ir) , Ipntr(ir+1) - 1 ic = Indcol(ip) ! array iwa marks columns which have contributed to ! the degree count of column jcol. update the degree ! counts of these columns as well as column jcol. if ( Iwa(ic)<jcol ) then Iwa(ic) = jcol Ndeg(ic) = Ndeg(ic) + 1 Ndeg(jcol) = Ndeg(jcol) + 1 endif enddo enddo enddo end subroutine degr