Epetra Package Browser (Single Doxygen Collection)
Development
|
Epetra_CrsSingletonFilter: A class for explicitly eliminating matrix rows and columns. More...
#include <Epetra_CrsSingletonFilter.h>
Protected Member Functions | |
Epetra_CrsMatrix * | FullCrsMatrix () const |
const Epetra_Map & | FullMatrixRowMap () const |
const Epetra_Map & | FullMatrixColMap () const |
const Epetra_Map & | FullMatrixDomainMap () const |
const Epetra_Map & | FullMatrixRangeMap () const |
void | InitializeDefaults () |
int | ComputeEliminateMaps () |
int | Setup (Epetra_LinearProblem *Problem) |
int | InitFullMatrixAccess () |
int | GetRow (int Row, int &NumIndices, int *&Indices) |
int | GetRowGCIDs (int Row, int &NumIndices, double *&Values, int *&GlobalIndices) |
int | GetRow (int Row, int &NumIndices, double *&Values, int *&Indices) |
int | CreatePostSolveArrays (const Epetra_IntVector &RowIDs, const Epetra_MapColoring &RowMapColors, const Epetra_IntVector &ColProfiles, const Epetra_IntVector &NewColProfiles, const Epetra_IntVector &ColHasRowWithSingleton) |
int | ConstructRedistributeExporter (Epetra_Map *SourceMap, Epetra_Map *TargetMap, Epetra_Export *&RedistributeExporter, Epetra_Map *&RedistributeMap) |
Private Member Functions | |
Epetra_CrsSingletonFilter (const Epetra_CrsSingletonFilter &Problem) | |
Copy constructor (defined as private so it is unavailable to user). More... | |
Epetra_CrsSingletonFilter & | operator= (const Epetra_CrsSingletonFilter &Problem) |
Constructors/Destructor | |
Epetra_CrsSingletonFilter () | |
Epetra_CrsSingletonFilter default constructor. More... | |
virtual | ~Epetra_CrsSingletonFilter () |
Epetra_CrsSingletonFilter Destructor. More... | |
Analyze methods | |
int | Analyze (Epetra_RowMatrix *FullMatrix) |
Analyze the input matrix, removing row/column pairs that have singletons. More... | |
bool | SingletonsDetected () const |
Returns true if singletons were detected in this matrix (must be called after Analyze() to be effective). More... | |
Reduce methods | |
int | ConstructReducedProblem (Epetra_LinearProblem *Problem) |
Return a reduced linear problem based on results of Analyze(). More... | |
int | UpdateReducedProblem (Epetra_LinearProblem *Problem) |
Update a reduced linear problem using new values. More... | |
Methods to construct Full System Solution | |
int | ComputeFullSolution () |
Compute a solution for the full problem using the solution of the reduced problem, put in LHS of FullProblem(). More... | |
Filter Statistics | |
int | NumRowSingletons () const |
Return number of rows that contain a single entry, returns -1 if Analysis not performed yet. More... | |
int | NumColSingletons () const |
Return number of columns that contain a single entry that are not associated with singleton row, returns -1 if Analysis not performed yet. More... | |
int | NumSingletons () const |
Return total number of singletons detected, returns -1 if Analysis not performed yet. More... | |
double | RatioOfDimensions () const |
Returns ratio of reduced system to full system dimensions, returns -1.0 if reduced problem not constructed. More... | |
double | RatioOfNonzeros () const |
Returns ratio of reduced system to full system nonzero count, returns -1.0 if reduced problem not constructed. More... | |
Attribute Access Methods | |
Epetra_LinearProblem * | FullProblem () const |
Returns pointer to the original unreduced Epetra_LinearProblem. More... | |
Epetra_LinearProblem * | ReducedProblem () const |
Returns pointer to the derived reduced Epetra_LinearProblem. More... | |
Epetra_RowMatrix * | FullMatrix () const |
Returns pointer to Epetra_CrsMatrix from full problem. More... | |
Epetra_CrsMatrix * | ReducedMatrix () const |
Returns pointer to Epetra_CrsMatrix from full problem. More... | |
Epetra_MapColoring * | RowMapColors () const |
Returns pointer to Epetra_MapColoring object: color 0 rows are part of reduced system. More... | |
Epetra_MapColoring * | ColMapColors () const |
Returns pointer to Epetra_MapColoring object: color 0 columns are part of reduced system. More... | |
Epetra_Map * | ReducedMatrixRowMap () const |
Returns pointer to Epetra_Map describing the reduced system row distribution. More... | |
Epetra_Map * | ReducedMatrixColMap () const |
Returns pointer to Epetra_Map describing the reduced system column distribution. More... | |
Epetra_Map * | ReducedMatrixDomainMap () const |
Returns pointer to Epetra_Map describing the domain map for the reduced system. More... | |
Epetra_Map * | ReducedMatrixRangeMap () const |
Returns pointer to Epetra_Map describing the range map for the reduced system. More... | |
Epetra_CrsSingletonFilter: A class for explicitly eliminating matrix rows and columns.
The Epetra_CrsSingletonFilter class takes an existing Epetra_LinearProblem object, analyzes it 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 Epetra_LinearProblem object. The reduced system can then be solved using any solver that is understands an Epetra_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 115 of file Epetra_CrsSingletonFilter.h.
Epetra_CrsSingletonFilter::Epetra_CrsSingletonFilter | ( | ) |
Epetra_CrsSingletonFilter default constructor.
Definition at line 66 of file Epetra_CrsSingletonFilter.cpp.
|
virtual |
Epetra_CrsSingletonFilter Destructor.
Definition at line 70 of file Epetra_CrsSingletonFilter.cpp.
|
private |
Copy constructor (defined as private so it is unavailable to user).
int Epetra_CrsSingletonFilter::Analyze | ( | Epetra_RowMatrix * | 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 153 of file Epetra_CrsSingletonFilter.cpp.
|
inline |
Returns true if singletons were detected in this matrix (must be called after Analyze() to be effective).
Definition at line 142 of file Epetra_CrsSingletonFilter.h.
int Epetra_CrsSingletonFilter::ConstructReducedProblem | ( | Epetra_LinearProblem * | Problem | ) |
Return a reduced linear problem based on results of Analyze().
Creates a new Epetra_LinearProblem object based on the results of the Analyze phase. A pointer to the reduced problem is obtained via a call to ReducedProblem().
Definition at line 288 of file Epetra_CrsSingletonFilter.cpp.
int Epetra_CrsSingletonFilter::UpdateReducedProblem | ( | Epetra_LinearProblem * | Problem | ) |
Update a reduced linear problem using new values.
Updates an existing Epetra_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 458 of file Epetra_CrsSingletonFilter.cpp.
int Epetra_CrsSingletonFilter::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 Epetra_LinearProblem.
Definition at line 643 of file Epetra_CrsSingletonFilter.cpp.
|
inline |
Return number of rows that contain a single entry, returns -1 if Analysis not performed yet.
Definition at line 177 of file Epetra_CrsSingletonFilter.h.
|
inline |
Return number of columns that contain a single entry that are not associated with singleton row, returns -1 if Analysis not performed yet.
Definition at line 180 of file Epetra_CrsSingletonFilter.h.
|
inline |
Return total number of singletons detected, returns -1 if Analysis not performed yet.
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 188 of file Epetra_CrsSingletonFilter.h.
|
inline |
Returns ratio of reduced system to full system dimensions, returns -1.0 if reduced problem not constructed.
Definition at line 191 of file Epetra_CrsSingletonFilter.h.
|
inline |
Returns ratio of reduced system to full system nonzero count, returns -1.0 if reduced problem not constructed.
Definition at line 194 of file Epetra_CrsSingletonFilter.h.
|
inline |
Returns pointer to the original unreduced Epetra_LinearProblem.
Definition at line 201 of file Epetra_CrsSingletonFilter.h.
|
inline |
Returns pointer to the derived reduced Epetra_LinearProblem.
Definition at line 204 of file Epetra_CrsSingletonFilter.h.
|
inline |
Returns pointer to Epetra_CrsMatrix from full problem.
Definition at line 207 of file Epetra_CrsSingletonFilter.h.
|
inline |
Returns pointer to Epetra_CrsMatrix from full problem.
Definition at line 210 of file Epetra_CrsSingletonFilter.h.
|
inline |
Returns pointer to Epetra_MapColoring object: color 0 rows are part of reduced system.
Definition at line 213 of file Epetra_CrsSingletonFilter.h.
|
inline |
Returns pointer to Epetra_MapColoring object: color 0 columns are part of reduced system.
Definition at line 216 of file Epetra_CrsSingletonFilter.h.
|
inline |
Returns pointer to Epetra_Map describing the reduced system row distribution.
Definition at line 219 of file Epetra_CrsSingletonFilter.h.
|
inline |
Returns pointer to Epetra_Map describing the reduced system column distribution.
Definition at line 222 of file Epetra_CrsSingletonFilter.h.
|
inline |
Returns pointer to Epetra_Map describing the domain map for the reduced system.
Definition at line 225 of file Epetra_CrsSingletonFilter.h.
|
inline |
Returns pointer to Epetra_Map describing the range map for the reduced system.
Definition at line 228 of file Epetra_CrsSingletonFilter.h.
|
inlineprotected |
Definition at line 236 of file Epetra_CrsSingletonFilter.h.
|
inlineprotected |
Definition at line 238 of file Epetra_CrsSingletonFilter.h.
|
inlineprotected |
Definition at line 239 of file Epetra_CrsSingletonFilter.h.
|
inlineprotected |
Definition at line 240 of file Epetra_CrsSingletonFilter.h.
|
inlineprotected |
Definition at line 241 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 100 of file Epetra_CrsSingletonFilter.cpp.
|
protected |
|
protected |
|
protected |
Definition at line 690 of file Epetra_CrsSingletonFilter.cpp.
|
protected |
Definition at line 703 of file Epetra_CrsSingletonFilter.cpp.
|
protected |
Definition at line 731 of file Epetra_CrsSingletonFilter.cpp.
|
protected |
Definition at line 716 of file Epetra_CrsSingletonFilter.cpp.
|
protected |
Definition at line 742 of file Epetra_CrsSingletonFilter.cpp.
|
protected |
Definition at line 552 of file Epetra_CrsSingletonFilter.cpp.
|
private |
|
protected |
Definition at line 259 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 260 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 261 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 262 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 263 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 264 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 265 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 267 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 268 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 269 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 270 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 271 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 272 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 273 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 274 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 276 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 277 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 278 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 279 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 282 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 283 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 285 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 286 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 287 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 288 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 289 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 290 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 292 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 293 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 294 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 295 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 297 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 298 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 299 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 300 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 301 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 302 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 304 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 305 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 306 of file Epetra_CrsSingletonFilter.h.
|
protected |
Definition at line 307 of file Epetra_CrsSingletonFilter.h.