Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
List of all members
Tpetra::CrsSingletonFilter_LinearProblem< Scalar, LocalOrdinal, GlobalOrdinal, Node > Class Template Reference

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_typeFullProblem () const
 Returns RCP to the original full (unreduced) Tpetra::LinearProblem.
 
Teuchos::RCP< linear_problem_typeReducedProblem () const
 Returns RCP to the derived reduced Tpetra::LinearProblem.
 
Teuchos::RCP< row_matrix_typeFullMatrix () const
 ! Returns RCP to Tpetra::RowMatrix from full problem.
 
Teuchos::RCP< crs_matrix_typeReducedMatrix () const
 Returns RCP to Tpetra::CrsMatrix from reduced problem.
 
Teuchos::RCP< const map_typeReducedMatrixRowMap () const
 Returns RCP to Tpetra::Map describing the reduced system row distribution.
 
Teuchos::RCP< const map_typeReducedMatrixColMap () const
 Returns RCP to Tpetra::Map describing the reduced system column distribution.
 
Teuchos::RCP< const map_typeReducedMatrixDomainMap () const
 Returns RCP to Tpetra::Map describing the domain map for the reduced system.
 
Teuchos::RCP< const map_typeReducedMatrixRangeMap () 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.
 

Detailed Description

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
class Tpetra::CrsSingletonFilter_LinearProblem< Scalar, LocalOrdinal, GlobalOrdinal, Node >

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:

  1. Construct full problem: Construct and Tpetra::LinearProblem containing the "full" matrix, RHS and LHS. This is done outside of Tpetra:CrsSingletonFilter class. Presumably, you have some reason to believe that this system may contain singletons.
  2. Construct an Tpetra::CrsSingletonFilter instance: Constructor needs no arguments.
  3. Analyze matrix: Invoke the Analyze() method, passing in the Tpetra::RowMatrix object from your full linear problem mentioned in the first step above.
  4. Go/No Go decision to construct reduced problem: Query the results of the Analyze method using the SingletonsDetected() method. This method returns "true" if there were singletons found in the matrix. You can also query any of the other methods in the Filter Statistics section to determine if you want to proceed with the construction of the reduced system.
  5. Construct reduced problem: If, in the previous step, you determine that you want to proceed with the construction of the reduced problem, you should next call the ConstructReducedProblem() method, passing in the full linear problem object from the first step. This method will use the information from the Analyze() method to construct a reduce problem that has explicitly eliminated the singleton rows, solved for the corresponding LHS values and updated the RHS. This step will also remove singleton columns from the reduced system. Once the solution of the reduced problem is is computed (via any solver that understands an Tpetra::LinearProblem), you should call the ComputeFullSolution() method to compute the LHS values assocaited with the singleton columns.
  6. Solve reduced problem: Obtain an RCP to the reduced problem using the ReducedProblem() method. Using the solver of your choice, solve the reduced system.
  7. Compute solution to full problem: Once the solution of the reduced problem is determined, the ComputeFullSolution() method will place the reduced solution values into the appropriate locations of the full solution LHS and then compute the values associated with column singletons. At this point, you have a complete solution to the original full problem.
  8. Solve a subsequent full problem that differs from the original problem only in values: It is often the case that the structure of a problem will be the same for a sequence of linear problems. In this case, the UpdateReducedProblem() method can be useful. After going through the above process one time, if you have a linear problem that is structurally identical to the previous problem, you can minimize memory and time costs by using the UpdateReducedProblem() method, passing in the subsequent problem. Once you have called the UpdateReducedProblem() method, you can then solve the reduce problem problem as you wish, and then compute the full solution as before. The RCP generated by ReducedProblem() will not change when UpdateReducedProblem() is called.

Definition at line 120 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.

Constructor & Destructor Documentation

◆ CrsSingletonFilter_LinearProblem()

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.

◆ ~CrsSingletonFilter_LinearProblem()

Destructor.

Member Function Documentation

◆ operator()()

Analysis of transform operation on original object and construction of new object.

Returns
Returns an RCP to the newly created object of type NewType.

Implements Tpetra::Transform< T, T >.

Definition at line 456 of file Tpetra_CrsSingletonFilter_LinearProblem_def.hpp.

◆ analyze()

void Tpetra::CrsSingletonFilter_LinearProblem< Scalar, LocalOrdinal, GlobalOrdinal, Node >::analyze ( const OriginalType &  orig)
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.

◆ construct()

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.

◆ fwd()

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.

◆ rvs()

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.

◆ Analyze()

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.

◆ SingletonsDetected()

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.

◆ ConstructReducedProblem()

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.

◆ UpdateReducedProblem()

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.

◆ 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.

◆ NumSingletonRows()

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.

◆ NumSingletonCols()

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.

◆ NumSingletons()

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.

Warning
This method returns -1 if Analyze() method has not been called.

Definition at line 269 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.

◆ RatioOfDimensions()

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.

◆ RatioOfNonzeros()

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.

◆ FullProblem()

Returns RCP to the original full (unreduced) Tpetra::LinearProblem.

Definition at line 284 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.

◆ ReducedProblem()

Teuchos::RCP< linear_problem_type > Tpetra::CrsSingletonFilter_LinearProblem< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ReducedProblem ( ) const
inline

Returns RCP to the derived reduced Tpetra::LinearProblem.

Definition at line 287 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.

◆ FullMatrix()

! Returns RCP to Tpetra::RowMatrix from full problem.

Definition at line 290 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.

◆ ReducedMatrix()

Teuchos::RCP< crs_matrix_type > Tpetra::CrsSingletonFilter_LinearProblem< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ReducedMatrix ( ) const
inline

Returns RCP to Tpetra::CrsMatrix from reduced problem.

Definition at line 293 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.

◆ ReducedMatrixRowMap()

Teuchos::RCP< const map_type > Tpetra::CrsSingletonFilter_LinearProblem< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ReducedMatrixRowMap ( ) const
inline

Returns RCP to Tpetra::Map describing the reduced system row distribution.

Definition at line 296 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.

◆ ReducedMatrixColMap()

Teuchos::RCP< const map_type > Tpetra::CrsSingletonFilter_LinearProblem< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ReducedMatrixColMap ( ) const
inline

Returns RCP to Tpetra::Map describing the reduced system column distribution.

Definition at line 299 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.

◆ ReducedMatrixDomainMap()

Teuchos::RCP< const map_type > Tpetra::CrsSingletonFilter_LinearProblem< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ReducedMatrixDomainMap ( ) const
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.

◆ ReducedMatrixRangeMap()

Teuchos::RCP< const map_type > Tpetra::CrsSingletonFilter_LinearProblem< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ReducedMatrixRangeMap ( ) const
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.

◆ isConstructed()

bool Tpetra::Transform< T, T >::isConstructed ( )
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.

Member Data Documentation

◆ localNumSingletonRows_

local_ordinal_type Tpetra::CrsSingletonFilter_LinearProblem< Scalar, LocalOrdinal, GlobalOrdinal, Node >::localNumSingletonRows_
protected

Number of singleton rows.

Definition at line 364 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.

◆ localNumSingletonCols_

local_ordinal_type Tpetra::CrsSingletonFilter_LinearProblem< Scalar, LocalOrdinal, GlobalOrdinal, Node >::localNumSingletonCols_
protected

Number of singleton columns not eliminated by singleton rows.

Definition at line 365 of file Tpetra_CrsSingletonFilter_LinearProblem_decl.hpp.


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