AbstractLinAlgPack: C++ Interfaces For Vectors, Matrices And Related Linear Algebra Objects  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Types | Public Member Functions | List of all members
AbstractLinAlgPack::MatrixConvertToSparse Class Referenceabstract

Mix-in interface for extracing explicit elements from a sparse matrix in one of several Fortran compatible formats. More...

#include <AbstractLinAlgPack_MatrixConvertToSparse.hpp>

Inheritance diagram for AbstractLinAlgPack::MatrixConvertToSparse:
Inheritance graph
[legend]

Public Types

enum  EExtractRegion { EXTRACT_FULL_MATRIX, EXTRACT_UPPER_TRIANGULAR, EXTRACT_LOWER_TRIANGULAR }
 
enum  EElementUniqueness { ELEMENTS_FORCE_UNIQUE, ELEMENTS_ALLOW_DUPLICATES_SUM }
 

Public Member Functions

virtual index_type num_nonzeros (EExtractRegion extract_region, EElementUniqueness element_uniqueness) const =0
 Returns the number of nonzeros in the matrix. More...
 
virtual 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=0, const index_type col_offset=0) const =0
 Extract sparse elements in a coordinate data structure. More...
 
- Public Member Functions inherited from AbstractLinAlgPack::MatrixBase
virtual ~MatrixBase ()
 Virtual destructor. More...
 
virtual const VectorSpacespace_cols () const =0
 Vector space for vectors that are compatible with the columns of the matrix. More...
 
virtual const VectorSpacespace_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...
 

Detailed Description

Mix-in interface for extracing explicit elements from a 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 simplest interface for clients to use to extract nonzero elements from matrices.

The following matrix is used in the documentation for the following extraction functions.

   [  1.1       1.3   1.4 ] 1
   [      2.2       2.4 ] 2
 A =  [     3.2         ] 3
   [  4.1       4.3   4.4 ] 4
   [  5.1       5.3   5.4 ] 5
       1      2     3     4

Note that above, A has:

Definition at line 88 of file AbstractLinAlgPack_MatrixConvertToSparse.hpp.

Member Enumeration Documentation

Enumerator
EXTRACT_FULL_MATRIX 

Extract the nonzero elements from the full matrix.

EXTRACT_UPPER_TRIANGULAR 

Extract the nonzero elements from the upper triangular region including the diagonal.

EXTRACT_LOWER_TRIANGULAR 

Extract the nonzero elements from the lower triangular region including the diagonal.

Definition at line 94 of file AbstractLinAlgPack_MatrixConvertToSparse.hpp.

Enumerator
ELEMENTS_FORCE_UNIQUE 

Entries must have unique row and column indexes.

ELEMENTS_ALLOW_DUPLICATES_SUM 

Entries with duplicate row and column indexes are allowed with the understanding that the values are summed.

Definition at line 100 of file AbstractLinAlgPack_MatrixConvertToSparse.hpp.

Member Function Documentation

virtual index_type AbstractLinAlgPack::MatrixConvertToSparse::num_nonzeros ( EExtractRegion  extract_region,
EElementUniqueness  element_uniqueness 
) const
pure virtual

Returns the number of nonzeros in the matrix.

Parameters
extract_region[in] Determines what matix region to extract:
  • EXTRACT_FULL_MATRIX: Extract the nonzero elements for the full matrix.
  • EXTRACT_UPPER_TRIANGULAR: Extract the nonzero elements for the upper triangular region (including the diagonal).
  • EXTRACT_LOWER_TRIANGULAR: Extract the nonzero elements for the lower triangular region (including the diagonal).
element_uniqueness[in] Determines if element row and column indexes must be unique.
  • ELEMENTS_FORCE_UNIQUE: The row and column indexes must be unique.
  • ELEMENTS_ALLOW_DUPLICATES_SUM: Entries with duplicate row and column indexes are allowed with the understanding that the values will be summed.

Implemented in AbstractLinAlgPack::MatrixExtractSparseElements, AbstractLinAlgPack::MatrixConvertToSparseEncap, and AbstractLinAlgPack::MatrixSymDiagSparse.

virtual void AbstractLinAlgPack::MatrixConvertToSparse::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 = 0,
const index_type  col_offset = 0 
) const
pure virtual

Extract sparse elements in a coordinate data structure.

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->num_nonzeros(extract_region,element_uniqueness) and then Arow[] and Acol[] are filled with the row and column indexes. If the client wants the nonzero values theit sets len_Aval = this->num_nonzeros(extract_region,element_uniqueness) and then Aval[] will be set on output.

The client can choose to extract the nonzeros for the full matrix (extract_region == EXTRACT_FULL_MATRIX) or the upper (extract_region == EXTRACT_UPPER_TRIANGULAR) or lower (extract_region == EXTRACT_LOWER_TRIANGULAR) triangular regions.

The client can force the extracted nonzero elements to have unique row and column indexes (element_uniqueness == ELEMENTS_FORCE_UNIQUE) or can allow elements with duplicate row and column indexes with the understanding that the values will be summed (element_uniqueness == ELEMENTS_ALLOW_DUPLICATES_SUM).

Parameters
extract_region[in] Determines what matix region to extract:
  • EXTRACT_FULL_MATRIX: Extract the nonzero elements for the full matrix.
  • EXTRACT_UPPER_TRIANGULAR: Extract the nonzero elements for the upper triangular region (including the diagonal).
  • EXTRACT_LOWER_TRIANGULAR: Extract the nonzero elements for the lower triangular region (including the diagonal).
element_uniqueness[in] Determines if element row and column indexes must be unique.
  • ELEMENTS_FORCE_UNIQUE: The row and column indexes must be unique.
  • ELEMENTS_ALLOW_DUPLICATES_SUM: Entries with duplicate row and column indexes are allowed with the understanding that the values will be summed.
len_Aval[in] If the client wants to extract the nonzero values of the matrix in Avar[] then len_Aval must be set to the return value from this->num_nonzeros(extract_region) 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 must 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 must be set to len_Aij = this->num_nonzeros(extract_region). 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 index 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 index 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.

Preconditions

  • [len_Aval != 0] len_Aval == this->num_nonzeros(etract_region,element_uniqueness) (throw std::invalid_argument).
  • [len_Aij != 0] len_Aij == this->num_nonzeros(etract_region,element_uniqueness) (throw std::invalid_argument).
  • [len_Aval == 0] Aval == NULL (throw std::invalid_argument).
  • [len_Aij == 0] Arow == NULL (throw std::invalid_argument).
  • [len_Aij == 0] Acol == NULL (throw std::invalid_argument).

To illustrate the behavior of this function consider the example matix A (shown above) in the following examples:

(1) To extract all of the nonzero elements and structure for the full matrix with the row and column indexes offset by 2 and 4 respectively the following input would yield the following output.

Input:

   extract_region     = EXTRACT_FULL_MATRIX
         element_uniqueness = ELEMENTS_FORCE_UNIQUE
   len_Aval           = A_nz = 12
   len_Aij            = A_nz = 12
   row_offset         = 2
   col_offset         = 4

Output:

   k   Aval[k-1]  Arow[k-1]  Acol[k-1]
   -  ---------  ---------  ---------
   1         1.1          3          5
   2         4.1          5          5
   3         5.1          7          5
   4         2.2          4          6
   5         3.2          5          6
   6         1.3          3          7
   7         4.3          6          7
   8         5.3          7          7
   9         1.4          3          8
   10        2.4          4          8
   11        4.4          6          8
   12        5.4          7          8

Note that in the example above that the elements are sorted by column but they do not have to be on output.

(2) To extract all of the nonzero elements and structure for the upper triangular region of the matrix with the row and column indexes offset by 0 and 0 respectively the following input would yield the following output.

Input:

   extract_region     = EXTRACT_UPPER_TRIANGULAR
         element_uniqueness = ELEMENTS_FORCE_UNIQUE
   len_Aval           = A_up_nz = 6
   len_Aij            = A_up_nz = 6
   row_offset         = 0
   col_offset         = 0

Output:

   k  Aval[k-1]  Arow[k-1]  Acol[k-1]
   -  ---------  ---------  ---------
   1        1.1          1          1
   2        2.2          2          2
   3        1.3          1          3
   4        1.4          1          4
   5        2.4          2          4
   6        4.4          4          4

Note that in the example above that the elements are sorted by column but they do not have to be on output.

(3) To extract all of the nonzero elements and structure for the lower triangular region of the matrix with the row and column indexes offset by 0 and 0 respectively the following input would yield the following output.

Input:

   extract_region     = EXTRACT_LOWER_TRIANGULAR
         element_uniqueness = ELEMENTS_FORCE_UNIQUE
   len_Aval           = A_lo_nz = 9
   len_Aij            = A_lo_nz = 9
   row_offset         = 0
   col_offset         = 0

Output:

   k  Aval[k-1]  Arow[k-1]  Acol[k-1]
   -  ---------  ---------  ---------
   1        1.1          1          1
   2        4.1          4          1
   3        5.1          5          1
   4        2.2          2          2
   5        3.2          3          2
   6        4.3          4          3
   7        5.3          5          3
   8        4.4          4          4
   9        5.4          5          4

Note that in the example above that the elements are sorted by column but they do not have to be on output.

Implemented in AbstractLinAlgPack::MatrixExtractSparseElements, AbstractLinAlgPack::MatrixConvertToSparseEncap, and AbstractLinAlgPack::MatrixSymDiagSparse.


The documentation for this class was generated from the following file: