|
Tpetra parallel linear algebra Version of the Day
|
A class for explicitly eliminating matrix rows and columns from a LinearProblem. More...
#include <Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp>
Inherits Tpetra::SameTypeTransform< T >.
Typedefs | |
| local_ordinal_type | localNumSingletonRows_ |
| Number of singleton rows. | |
| local_ordinal_type | localNumSingletonCols_ |
| Number of singleton columns not eliminated by singleton rows. | |
| CrsSingletonFilter_LinearProblem (bool run_on_host=false, bool verbose=false) | |
| Constructor. | |
| virtual | ~CrsSingletonFilter_LinearProblem ()=default |
| Destructor. | |
| NewType | operator() (const OriginalType &originalLinearProblem) |
| Analysis of transform operation on original object and construction of new object. | |
| void | analyze (const OriginalType &originalLinearProblem) |
| Initial analysis phase of transform. | |
| NewType | construct () |
| Construction of new object as a result of the transform. | |
| void | fwd () |
| Forward transfer of data. | |
| void | rvs () |
| Reverse transfer of data. | |
| void | Analyze (const Teuchos::RCP< row_matrix_type > &FullMatrix) |
| Analyze the input matrix, removing row/column pairs that have singletons. | |
| bool | SingletonsDetected () const |
| Returns true if singletons were detected in this matrix (must be called after Analyze() to be effective). | |
| void | ConstructReducedProblem (const Teuchos::RCP< linear_problem_type > &Problem) |
| Return a reduced linear problem based on results of Analyze(). | |
| void | UpdateReducedProblem (const Teuchos::RCP< linear_problem_type > &Problem) |
| Update a reduced linear problem using new values. | |
| void | ComputeFullSolution () |
| Compute a solution for the full problem using the solution of the reduced problem, put in LHS of FullProblem(). | |
| int | NumSingletonRows () const |
| Return number of rows that contain a single entry, returns -1 if Analysis has not been performed yet. | |
| int | NumSingletonCols () const |
| Return number of columns that contain a single entry that are not associated with singleton row, returns -1 if Analysis has not been performed yet. | |
| int | NumSingletons () const |
| Return total number of singletons detected. | |
| double | RatioOfDimensions () const |
| Returns ratio of reduced system to full system dimensions, returns -1.0 if reduced problem not constructed. | |
| double | RatioOfNonzeros () const |
| Returns ratio of reduced system to full system nonzero count, returns -1.0 if reduced problem not constructed. | |
| Teuchos::RCP< linear_problem_type > | FullProblem () const |
| Returns RCP to the original full (unreduced) Tpetra::LinearProblem. | |
| Teuchos::RCP< linear_problem_type > | ReducedProblem () const |
| Returns RCP to the derived reduced Tpetra::LinearProblem. | |
| Teuchos::RCP< row_matrix_type > | FullMatrix () const |
| ! Returns RCP to Tpetra::RowMatrix from full problem. | |
| Teuchos::RCP< crs_matrix_type > | ReducedMatrix () const |
| Returns RCP to Tpetra::CrsMatrix from reduced problem. | |
| Teuchos::RCP< const map_type > | ReducedMatrixRowMap () const |
| Returns RCP to Tpetra::Map describing the reduced system row distribution. | |
| Teuchos::RCP< const map_type > | ReducedMatrixColMap () const |
| Returns RCP to Tpetra::Map describing the reduced system column distribution. | |
| Teuchos::RCP< const map_type > | ReducedMatrixDomainMap () const |
| Returns RCP to Tpetra::Map describing the domain map for the reduced system. | |
| Teuchos::RCP< const map_type > | ReducedMatrixRangeMap () const |
| Returns RCP to Tpetra::Map describing the range map for the reduced system. | |
Virtual functions with default implements allowing for optional | |
implementation by the Transform developer | |
| virtual bool | isConstructed () |
| Check if transformed object has been constructed. | |
A class for explicitly eliminating matrix rows and columns from a LinearProblem.
The Tpetra::CrsSingletonFilter class takes an existing Tpetra::LinearProblem object, analyzes its structure and explicitly eliminates singleton rows and columns from the matrix and appropriately modifies the RHS and LHS of the linear problem. The result of this process is a reduced system of equations that is itself an Tpetra::LinearProblem object. The reduced system can then be solved using any solver that understands a Tpetra::LinearProblem. The solution for the full system is obtained by calling ComputeFullSolution().
Singleton rows are defined to be rows that have a single nonzero entry in the matrix. The equation associated with this row can be explicitly eliminated because it involved only one variable. For example if row i has a single nonzero value in column j, call it A(i,j), we can explicitly solve for x(j) = b(i)/A(i,j), where b(i) is the ith entry of the RHS and x(j) is the jth entry of the LHS.
Singleton columns are defined to be columns that have a single nonzero entry in the matrix. The variable associated with this column is fully dependent, meaning that the solution for all other variables does not depend on it. If this entry is A(i,j) then the ith row and jth column can be removed from the system and x(j) can be solved after the solution for all other variables is determined.
By removing singleton rows and columns, we can often produce a reduced system that is smaller and far less dense, and in general having better numerical properties.
The basic procedure for using this class is as follows:
Definition at line 120 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.
| Tpetra::CrsSingletonFilter_LinearProblem< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsSingletonFilter_LinearProblem | ( | bool | run_on_host = false, |
| bool | verbose = false |
||
| ) |
Constructor.
Definition at line 437 of file Tpetra_CrsSingletonFilter_LinearProblem_def.hpp.
|
virtualdefault |
Destructor.
|
virtual |
Analysis of transform operation on original object and construction of new object.
Implements Tpetra::Transform< T, T >.
Definition at line 456 of file Tpetra_CrsSingletonFilter_LinearProblem_def.hpp.
|
virtual |
Initial analysis phase of transform.
Initial analysis phase of transform to confirm the transform is possible allowing methods construct(), fwd() and rvs() to be successfully utilized.
The default implementation calls method operator() and stores the resulting object in an internal attribute newObj_.
Reimplemented from Tpetra::Transform< T, T >.
Definition at line 465 of file Tpetra_CrsSingletonFilter_LinearProblem_def.hpp.
|
virtual |
Construction of new object as a result of the transform.
The default implementation returns internal attribute newObj_.
Reimplemented from Tpetra::Transform< T, T >.
Definition at line 497 of file Tpetra_CrsSingletonFilter_LinearProblem_def.hpp.
|
virtual |
Forward transfer of data.
Forward transfer of data from orig object input in the operator() method call to the new object created in this same call.
Implements Tpetra::Transform< T, T >.
Definition at line 522 of file Tpetra_CrsSingletonFilter_LinearProblem_def.hpp.
|
virtual |
Reverse transfer of data.
Reverse transfer of data from new object created in the operator() method call to the orig object input to this same method.
Implements Tpetra::Transform< T, T >.
Definition at line 530 of file Tpetra_CrsSingletonFilter_LinearProblem_def.hpp.
| void Tpetra::CrsSingletonFilter_LinearProblem< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Analyze | ( | const Teuchos::RCP< row_matrix_type > & | FullMatrix | ) |
Analyze the input matrix, removing row/column pairs that have singletons.
Analyzes the user's input matrix to determine rows and columns that should be explicitly eliminated to create the reduced system. Look for rows and columns that have single entries. These rows/columns can easily be removed from the problem. The results of calling this method are two MapColoring objects accessible via RowMapColors() and ColMapColors() accessor methods. All rows/columns that would be eliminated in the reduced system have a color of 1 in the corresponding RowMapColors/ColMapColors object. All kept rows/cols have a color of 0.
Definition at line 538 of file Tpetra_CrsSingletonFilter_LinearProblem_def.hpp.
|
inline |
Returns true if singletons were detected in this matrix (must be called after Analyze() to be effective).
Definition at line 208 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.
| void Tpetra::CrsSingletonFilter_LinearProblem< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ConstructReducedProblem | ( | const Teuchos::RCP< linear_problem_type > & | Problem | ) |
Return a reduced linear problem based on results of Analyze().
Creates a new Tpetra::LinearProblem object based on the results of the Analyze phase. An RCP to the reduced problem is obtained via a call to ReducedProblem().
Definition at line 949 of file Tpetra_CrsSingletonFilter_LinearProblem_def.hpp.
| void Tpetra::CrsSingletonFilter_LinearProblem< Scalar, LocalOrdinal, GlobalOrdinal, Node >::UpdateReducedProblem | ( | const Teuchos::RCP< linear_problem_type > & | Problem | ) |
Update a reduced linear problem using new values.
Updates an existing Tpetra::LinearProblem object using new matrix, LHS and RHS values. The matrix structure must be identical to the matrix that was used to construct the original reduced problem.
Definition at line 1248 of file Tpetra_CrsSingletonFilter_LinearProblem_def.hpp.
| void Tpetra::CrsSingletonFilter_LinearProblem< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ComputeFullSolution | ( | ) |
Compute a solution for the full problem using the solution of the reduced problem, put in LHS of FullProblem().
After solving the reduced linear system, this method can be called to compute the solution to the original problem, assuming the solution for the reduced system is valid. The solution of the unreduced, original problem will be in the LHS of the original Tpetra::LinearProblem.
Definition at line 1502 of file Tpetra_CrsSingletonFilter_LinearProblem_def.hpp.
|
inline |
Return number of rows that contain a single entry, returns -1 if Analysis has not been performed yet.
Definition at line 254 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.
|
inline |
Return number of columns that contain a single entry that are not associated with singleton row, returns -1 if Analysis has not been performed yet.
Definition at line 259 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.
|
inline |
Return total number of singletons detected.
Return total number of singletons detected across all processors. This method will not return a valid result until after the Analyze() method is called. The dimension of the reduced system can be computed by subtracting this number from dimension of full system.
Definition at line 269 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.
|
inline |
Returns ratio of reduced system to full system dimensions, returns -1.0 if reduced problem not constructed.
Definition at line 273 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.
|
inline |
Returns ratio of reduced system to full system nonzero count, returns -1.0 if reduced problem not constructed.
Definition at line 277 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.
|
inline |
Returns RCP to the original full (unreduced) Tpetra::LinearProblem.
Definition at line 284 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.
|
inline |
Returns RCP to the derived reduced Tpetra::LinearProblem.
Definition at line 287 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.
|
inline |
! Returns RCP to Tpetra::RowMatrix from full problem.
Definition at line 290 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.
|
inline |
Returns RCP to Tpetra::CrsMatrix from reduced problem.
Definition at line 293 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.
|
inline |
Returns RCP to Tpetra::Map describing the reduced system row distribution.
Definition at line 296 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.
|
inline |
Returns RCP to Tpetra::Map describing the reduced system column distribution.
Definition at line 299 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.
|
inline |
Returns RCP to Tpetra::Map describing the domain map for the reduced system.
Definition at line 302 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.
|
inline |
Returns RCP to Tpetra::Map describing the range map for the reduced system.
Definition at line 305 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.
|
virtualinherited |
Check if transformed object has been constructed.
The default implementation returns true if newObj_ != 0.
Definition at line 94 of file Tpetra_Transform.hpp.
|
protected |
Number of singleton rows.
Definition at line 364 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.
|
protected |
Number of singleton columns not eliminated by singleton rows.
Definition at line 365 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.