AbstractLinAlgPack: C++ Interfaces For Vectors, Matrices And Related Linear Algebra Objects
Version of the Day
|
Interface for extracting nonzero elements from a banded subregion of a permuted sparse matrix in one of several Fortran compatible formats. More...
#include <AbstractLinAlgPack_MatrixExtractSparseElements.hpp>
Public Member Functions | |
virtual index_type | count_nonzeros (EElementUniqueness element_uniqueness, const index_type inv_row_perm[], const index_type inv_col_perm[], const Range1D &row_rng, const Range1D &col_rng, index_type dl, index_type du) const =0 |
Returns the number of nonzeros in the banded submatrix of the permuted matrix. More... | |
virtual void | coor_extract_nonzeros (EElementUniqueness element_uniqueness, const index_type inv_row_perm[], const index_type inv_col_perm[], const Range1D &row_rng, const Range1D &col_rng, index_type dl, index_type du, value_type alpha, const index_type len_Aval, value_type Aval[], const index_type len_Aij, index_type Arow[], index_type Acol[], const index_type row_offset=0, const index_type col_offset=0) const =0 |
Extract elements in a coordinate data structure. More... | |
Public Member Functions inherited from AbstractLinAlgPack::MatrixBase | |
virtual | ~MatrixBase () |
Virtual destructor. More... | |
virtual const VectorSpace & | space_cols () const =0 |
Vector space for vectors that are compatible with the columns of the matrix. More... | |
virtual const VectorSpace & | space_rows () const =0 |
Vector space for vectors that are compatible with the rows of the matrix. More... | |
virtual size_type | rows () const |
Return the number of rows in the matrix. More... | |
virtual size_type | cols () const |
Return the number of columns in the matrix. More... | |
virtual size_type | nz () const |
Return the number of nonzero elements in the matrix. More... | |
Overridden from MatrixConvertToSparse | |
index_type | num_nonzeros (EExtractRegion extract_region, EElementUniqueness element_uniqueness) const |
void | coor_extract_nonzeros (EExtractRegion extract_region, EElementUniqueness element_uniqueness, const index_type len_Aval, value_type Aval[], const index_type len_Aij, index_type Arow[], index_type Acol[], const index_type row_offset, const index_type col_offset) const |
Additional Inherited Members | |
Public Types inherited from AbstractLinAlgPack::MatrixConvertToSparse | |
enum | EExtractRegion { EXTRACT_FULL_MATRIX, EXTRACT_UPPER_TRIANGULAR, EXTRACT_LOWER_TRIANGULAR } |
enum | EElementUniqueness { ELEMENTS_FORCE_UNIQUE, ELEMENTS_ALLOW_DUPLICATES_SUM } |
Interface for extracting nonzero elements from a banded subregion of a permuted sparse matrix in one of several Fortran compatible formats.
The formats supported are:
Coordiante:
Aval[k], Arow[k], Acol[k], k = 0..num_nonzeros(...)-1
Compressed Row (Column):
Aval[k], Acol[k], k = 0..num_nonzeros(...)-1 Arow_start[j], j = 0..rows()-1
This is meant to be the do-all interface for clients to use to extract nonzero elements for sparse matrices.
The idea is that given a matrix A, row and column permutations P
and Q
, a range of rows and columns (rl,ru)
and (cl,cu)
defining a square submatrix M and a range of lower and upper bands dl
and du
within the submatrix M, this interface is used to extract those nonzero elements.
In Matlab like terms we have:
let M = (P'*A*Q)(rl:ru,cl:cu)
This interface extracts nonzero elements from M within the banded region [dl,du]
where d = 0
is the diagonal of M, d < 0
is below the diagonal and d > 0
is above the diagonal.
The following matrix is used in the documentation for the following extraction functions.
[ 1 6 9 ] 1 [ 4 10 ] 2 A = [ 5 ] 3 [ 2 7 11 ] 4 [ 3 8 12 ] 5 1 2 3 4
Note that above, A has:
A_nz == this->num_nonzeros(EXTRACT_FULL_MATRIX,ELEMENTS_FORCE_UNIQUE) == 12
A_up_nz = this->num_nonzeros(EXTRACT_UPPER_TRIANGULAR,ELEMENTS_FORCE_UNIQUE) == 6
A_lo_nz = this->num_nonzeros(EXTRACT_LOWER_TRIANGULAR,ELEMENTS_FORCE_UNIQUE) == 9
Definition at line 101 of file AbstractLinAlgPack_MatrixExtractSparseElements.hpp.
|
pure virtual |
Returns the number of nonzeros in the banded submatrix of the permuted matrix.
Let B be the banded submatrix that is being specifed and let A be this full matrix object. Then the (dense) elemements of B are defined by:
/ 0 : for dl <= (j-i) <= du B(i,j) = | \ A(row_perm(i+rl-1),col_perm(j+cl-1)) : for (j-i) < dl || du < (j-i) for i = 1..ru-rl+1, j = 1..cu-cl+1
Above rl = row_rng.lbound()
, ru = row_rng.ubound(), cl = col_rng.lbound()
and cu = col_rng.ubound()
.
Preconditions:
1 <= row_perm(i) <= this->rows(), i = 1..rows()
(throw std::out_of_bounds
) 1 <= col_perm(i) <= this->cols(), i = 1..cols()
(throw std::out_of_bounds
) row_rng.ubound() <= this->rows()
(throw std::out_of_bounds
) col_rng.ubound() <= this->cols()
(throw std::out_of_bounds
) du >= dl
(throw std::range_error
) -(ru-rl) <= dl && du <= (cu-cl)
(throw std::out_of_bounds
) To illustrate the behavior of this function consider the example matix A (see intro) in the following example:
Define the permutations as:
row_perm = { 2, 4, 1, 3, 5 } col_perm = { 3, 2, 1, 4 }
The permuted matrix (P'*A*Q)
would then be:
2 | 1 [ 4 10 ] 4 | 2 [ 7 2 11 ] 1 | 3 [ 6 1 9 ] (P'*A*Q) = 3 | 4 [ 5 ] 5 | 5 [ 8 3 12 ] 1 2 3 4 - - - - 3 2 1 4
Now define the square submatrix in the range:
row_rng = [ rl, ru ] = [ 2, 5 ] col_rng = [ cl, cu ] = [ 2, 4 ]
The square submatrix is then:
4 | 1 [ 2 11 ] (P'*A*Q)(r1:r2,cl:cu) = 1 | 2 [ 1 9 ] 3 | 3 [ 5 ] 1 2 3 - - - 2 1 4
Now define the range of diagonals as:
dl = -1, du = 1
This finally gives us our matrix that we wish to extract nonzeros from as (the x's show elemements that where excluded out of the diagonal range:
4 | 1 [ 4 x ] B = 1 | 2 [ 1 9 ] 3 | 3 [ x ] 1 2 3 - - - 2 1 4
Now we can see that for this example that this->count_nonzeros() would return 3.
In summary, for the example A shown above we have:
Input:
element_uniqueness = ELEMENTS_FORCE_UNIQUE row_perm[] = { 2, 4, 1, 3, 5 } col_perm[] = { 3, 2, 1, 4 } row_rng = [ 2, 5 ] col_rng = [ 2, 4 ] dl = -1 du = +1
Output:
3 <- count_nonzeros(element_uniqueness,row_perm,col_perm,row_rng,col_rng,dl,du)
element_uniqueness | [in] Determines if element row and column indexes must be unique.
|
inv_row_perm | [in] Arary (length this->rows() ) Defines the row permutation P . Row i of A is row inv_row_perm [i-1] of . If inv_row_perm==NULL on input then the identity permutation P = I is used. |
inv_col_perm | [in] Arary (length this->cols() ) Defines the column permutation Q . Column j of A is column inv_col_perm [j-1] of (A*Q) . If inv_col_perm==NULL on input then the identity permutation Q = I is used. |
row_rng | [in] Defines the range of rows [rl,ru] that the submatrix of (P'*A*Q) is taken from (see preconditions). |
col_rng | [in] Defines the range of columns [cl,cu] that the submatrix of (P'*A*Q) is taken from (see preconditions). |
dl | [in] Lower diagonal to extract elements above (see preconditions). |
du | [in] Upper diagonal to extract elements below (see preconditions). |
Implemented in AbstractLinAlgPack::MatrixSparseCOORSerial.
|
pure virtual |
Extract elements in a coordinate data structure.
Let B be the scaled banded submatrix that is being specifed and let A be this matrix object. Then the elemements of B are defined by:
/ 0 : for dl <= (j-i) <= du B(i,j) = | \ alpha*A(row_perm(i+rl-1),col_perm(j+cl-1)) : for (j-i) < dl || du < (j-i) for i = 1..ru-rl+1, j = 1..cu-cl+1
were rl = row_rng.lbound()
, ru = row_rng.ubound()
, cl = col_rng.lbound()
and cu = col_rng.ubound()
.
The client can extract the structure in Arow
[] and Acol
[] and/or the nonzero elements in Aval
[]. If the client wants to extract the structure it sets len_Aij = this->count_nonzeros()
and then Arow
[] and Acol
[] are filled with the row and column indexes. If the client wants the nonzero values it sets len_Aval = this->count_nonzeros()
and then Aval
[] will be set on output.
The input arguments passed to this->count_nonzeros()
must correspond to the same arguments passed to this function.
To illustrate the behavior of this function consider the same example as outlined in the documentation for count_nonzeros()
. Let's assume that we want to scale the example banded matrix by alpha = 2.0
, make the row and column indexes start at (4,7) and extract the structure (len_Aij > 0
) and the nonzero values (len_Aval
). The input and output from this function for this example would then be:
Input:
elements_uniqueness = ELEMENTS_FORCE_UNIQUE row_perm = { 2, 4, 1, 3, 5 } col_perm = { 3, 2, 1, 4 } row_rng = [ 2, 5 ] col_rng = [ 2, 4 ] dl = -1 du = +1 alpha = 2.0 len_Aval = count_nonzeros(...) = 3 len_Aij = count_nonzeros(...) = 3 row_offset = 3 col_offset = 6
Output:
k A(i,j) B(i,j) Avar Arow Acol - ------ ------ ---- ---- ---- 1 4(4,1) 4(1,2) 8 4 8 2 1(1,1) 1(2,2) 2 5 8 3 9(1,5) 9(2,3) 18 5 9
Some of the most common uses of this interface are:
(bl = -(ru-rl+1) && bu = (cu-cl+1))
(bl = 0 && bu = (cu-cl+1))
(bl = -(ru-rl+1) && bu = 0)
element_uniqueness | [in] Same as passed to count_nonzeros (...). |
row_perm | [in] Same as passed to count_nonzeros (...). |
col_perm | [in] Same as passed to count_nonzeros (...). |
row_rng | [in] Same as passed to count_nonzeros (...). |
col_rng | [in] Same as passed to count_nonzeros (...). |
dl | [in] Same as passed to count_nonzeros (...). |
du | [in] Same as passed to count_nonzeros (...). |
alpha | [in] Scaling parameter (see above) |
len_Aval | [in] If the client wants to extract the nonzero values of the matrix in Aval [] then this should be set to len_Aval = this->count_nonzeros(...) . Otherwise len_Avar should be set to zero. |
Aval | [out] If len_Aval > 0 then Aval must point to the begining of an array of length len_Aval . If len_Aval == 0 then Aval may be NULL. On output, Aval [k] is the nonzero value of the kth nonzero element, for k = 0...len_Aval-1 . |
len_Aij | [in] If the client wants to extract the structure of the matrix in Arow [] and Acol [] then this should be set to len_Aij = this->num_nonzeros(...) . Otherwise len_Aij should be set to zero. |
Arow | [out] If len_Aij > 0 then Arow must point to the begining of an array of length len_Aij . If len_Aij == 0 then Arow must be NULL . On output, Arow [k] is the row index of the kth nonzero element, for k = 0...len_Aij-1 . |
Acol | [out] If len_Aij > 0 then Acol must point to the begining of an array of length len_Aij . If len_Aij == 0 then Acol must be NULL . On output, Acol [k] is the column index of the kth nonzero element, for k = 0...len_Aij-1 . |
row_offset | [in] The row indexes in Arow [] are offset by +row_offset on output. To leave the first row index as 1 set row_offset = 0 . This is to allow the client to place this matrix as a submatrix into a larger matrix stored in the coordinate format. Default value = 0. |
col_offset | [in] The column indexes in Acol [] are offset by +col_offset on output. To leave the first column index as 1 set col_offset = 0 . This is to allow the client to place this matrix as a submatrix into a larger matrix stored in the coordinate format. Default value = 0. |
Implemented in AbstractLinAlgPack::MatrixSparseCOORSerial.
|
virtual |
Implements AbstractLinAlgPack::MatrixConvertToSparse.
Definition at line 51 of file AbstractLinAlgPack_MatrixExtractSparseElements.cpp.
|
virtual |
Implements AbstractLinAlgPack::MatrixConvertToSparse.
Definition at line 62 of file AbstractLinAlgPack_MatrixExtractSparseElements.cpp.