Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
Epetra_CrsSingletonFilter Class Reference

Epetra_CrsSingletonFilter: A class for explicitly eliminating matrix rows and columns. More...

#include <Epetra_CrsSingletonFilter.h>

Protected Member Functions

Epetra_CrsMatrixFullCrsMatrix () const
 
const Epetra_MapFullMatrixRowMap () const
 
const Epetra_MapFullMatrixColMap () const
 
const Epetra_MapFullMatrixDomainMap () const
 
const Epetra_MapFullMatrixRangeMap () 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)
 

Protected Attributes

Epetra_LinearProblemFullProblem_
 
Epetra_LinearProblemReducedProblem_
 
Epetra_RowMatrixFullMatrix_
 
Epetra_CrsMatrixFullCrsMatrix_
 
Epetra_CrsMatrixReducedMatrix_
 
Epetra_MultiVectorReducedRHS_
 
Epetra_MultiVectorReducedLHS_
 
Epetra_MapReducedMatrixRowMap_
 
Epetra_MapReducedMatrixColMap_
 
Epetra_MapReducedMatrixDomainMap_
 
Epetra_MapReducedMatrixRangeMap_
 
Epetra_MapOrigReducedMatrixDomainMap_
 
Epetra_ImportFull2ReducedRHSImporter_
 
Epetra_ImportFull2ReducedLHSImporter_
 
Epetra_ExportRedistributeDomainExporter_
 
int * ColSingletonRowLIDs_
 
int * ColSingletonColLIDs_
 
int * ColSingletonPivotLIDs_
 
double * ColSingletonPivots_
 
int AbsoluteThreshold_
 
double RelativeThreshold_
 
int NumMyRowSingletons_
 
int NumMyColSingletons_
 
int NumGlobalRowSingletons_
 
int NumGlobalColSingletons_
 
double RatioOfDimensions_
 
double RatioOfNonzeros_
 
bool HaveReducedProblem_
 
bool UserDefinedEliminateMaps_
 
bool AnalysisDone_
 
bool SymmetricElimination_
 
Epetra_MultiVectortempExportX_
 
Epetra_MultiVectortempX_
 
Epetra_MultiVectortempB_
 
Epetra_MultiVectorRedistributeReducedLHS_
 
int * Indices_
 
Epetra_SerialDenseVector Values_
 
Epetra_MapColoringRowMapColors_
 
Epetra_MapColoringColMapColors_
 
bool FullMatrixIsCrsMatrix_
 
int MaxNumMyEntries_
 

Private Member Functions

 Epetra_CrsSingletonFilter (const Epetra_CrsSingletonFilter &Problem)
 Copy constructor (defined as private so it is unavailable to user). More...
 
Epetra_CrsSingletonFilteroperator= (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_LinearProblemFullProblem () const
 Returns pointer to the original unreduced Epetra_LinearProblem. More...
 
Epetra_LinearProblemReducedProblem () const
 Returns pointer to the derived reduced Epetra_LinearProblem. More...
 
Epetra_RowMatrixFullMatrix () const
 Returns pointer to Epetra_CrsMatrix from full problem. More...
 
Epetra_CrsMatrixReducedMatrix () const
 Returns pointer to Epetra_CrsMatrix from full problem. More...
 
Epetra_MapColoringRowMapColors () const
 Returns pointer to Epetra_MapColoring object: color 0 rows are part of reduced system. More...
 
Epetra_MapColoringColMapColors () const
 Returns pointer to Epetra_MapColoring object: color 0 columns are part of reduced system. More...
 
Epetra_MapReducedMatrixRowMap () const
 Returns pointer to Epetra_Map describing the reduced system row distribution. More...
 
Epetra_MapReducedMatrixColMap () const
 Returns pointer to Epetra_Map describing the reduced system column distribution. More...
 
Epetra_MapReducedMatrixDomainMap () const
 Returns pointer to Epetra_Map describing the domain map for the reduced system. More...
 
Epetra_MapReducedMatrixRangeMap () const
 Returns pointer to Epetra_Map describing the range map for the reduced system. More...
 

Detailed Description

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:
  1. Construct full problem: Construct and Epetra_LinearProblem containing the "full" matrix, RHS and LHS. This is done outside of Epetra_CrsSingletonFilter class. Presumably, you have some reason to believe that this system may contain singletons.
  2. Construct an Epetra_CrsSingletonFilter instance: Constructor needs no arguments.
  3. Analyze matrix: Invoke the Analyze() method, passing in the Epetra_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 Epetra_LinearProblem), you should call the ComputeFullSolution() method to compute the LHS values assocaited with the singleton columns.
  6. Solve reduced problem: Obtain a pointer 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 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 structural 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 pointer generated by ReducedProblem() will not change when UpdateReducedProblem() is called.

Definition at line 115 of file Epetra_CrsSingletonFilter.h.

Constructor & Destructor Documentation

Epetra_CrsSingletonFilter::Epetra_CrsSingletonFilter ( )

Epetra_CrsSingletonFilter default constructor.

Definition at line 66 of file Epetra_CrsSingletonFilter.cpp.

Epetra_CrsSingletonFilter::~Epetra_CrsSingletonFilter ( )
virtual

Epetra_CrsSingletonFilter Destructor.

Definition at line 70 of file Epetra_CrsSingletonFilter.cpp.

Epetra_CrsSingletonFilter::Epetra_CrsSingletonFilter ( const Epetra_CrsSingletonFilter Problem)
private

Copy constructor (defined as private so it is unavailable to user).

Member Function Documentation

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.

bool Epetra_CrsSingletonFilter::SingletonsDetected ( ) const
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().

Returns
Error code, set to 0 if no error.

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.

Returns
Error code, set to 0 if no error.

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.

int Epetra_CrsSingletonFilter::NumRowSingletons ( ) const
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.

int Epetra_CrsSingletonFilter::NumColSingletons ( ) const
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.

int Epetra_CrsSingletonFilter::NumSingletons ( ) const
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.

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

Definition at line 188 of file Epetra_CrsSingletonFilter.h.

double Epetra_CrsSingletonFilter::RatioOfDimensions ( ) const
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.

double Epetra_CrsSingletonFilter::RatioOfNonzeros ( ) const
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.

Epetra_LinearProblem* Epetra_CrsSingletonFilter::FullProblem ( ) const
inline

Returns pointer to the original unreduced Epetra_LinearProblem.

Definition at line 201 of file Epetra_CrsSingletonFilter.h.

Epetra_LinearProblem* Epetra_CrsSingletonFilter::ReducedProblem ( ) const
inline

Returns pointer to the derived reduced Epetra_LinearProblem.

Definition at line 204 of file Epetra_CrsSingletonFilter.h.

Epetra_RowMatrix* Epetra_CrsSingletonFilter::FullMatrix ( ) const
inline

Returns pointer to Epetra_CrsMatrix from full problem.

Definition at line 207 of file Epetra_CrsSingletonFilter.h.

Epetra_CrsMatrix* Epetra_CrsSingletonFilter::ReducedMatrix ( ) const
inline

Returns pointer to Epetra_CrsMatrix from full problem.

Definition at line 210 of file Epetra_CrsSingletonFilter.h.

Epetra_MapColoring* Epetra_CrsSingletonFilter::RowMapColors ( ) const
inline

Returns pointer to Epetra_MapColoring object: color 0 rows are part of reduced system.

Definition at line 213 of file Epetra_CrsSingletonFilter.h.

Epetra_MapColoring* Epetra_CrsSingletonFilter::ColMapColors ( ) const
inline

Returns pointer to Epetra_MapColoring object: color 0 columns are part of reduced system.

Definition at line 216 of file Epetra_CrsSingletonFilter.h.

Epetra_Map* Epetra_CrsSingletonFilter::ReducedMatrixRowMap ( ) const
inline

Returns pointer to Epetra_Map describing the reduced system row distribution.

Definition at line 219 of file Epetra_CrsSingletonFilter.h.

Epetra_Map* Epetra_CrsSingletonFilter::ReducedMatrixColMap ( ) const
inline

Returns pointer to Epetra_Map describing the reduced system column distribution.

Definition at line 222 of file Epetra_CrsSingletonFilter.h.

Epetra_Map* Epetra_CrsSingletonFilter::ReducedMatrixDomainMap ( ) const
inline

Returns pointer to Epetra_Map describing the domain map for the reduced system.

Definition at line 225 of file Epetra_CrsSingletonFilter.h.

Epetra_Map* Epetra_CrsSingletonFilter::ReducedMatrixRangeMap ( ) const
inline

Returns pointer to Epetra_Map describing the range map for the reduced system.

Definition at line 228 of file Epetra_CrsSingletonFilter.h.

Epetra_CrsMatrix* Epetra_CrsSingletonFilter::FullCrsMatrix ( ) const
inlineprotected

Definition at line 236 of file Epetra_CrsSingletonFilter.h.

const Epetra_Map& Epetra_CrsSingletonFilter::FullMatrixRowMap ( ) const
inlineprotected

Definition at line 238 of file Epetra_CrsSingletonFilter.h.

const Epetra_Map& Epetra_CrsSingletonFilter::FullMatrixColMap ( ) const
inlineprotected

Definition at line 239 of file Epetra_CrsSingletonFilter.h.

const Epetra_Map& Epetra_CrsSingletonFilter::FullMatrixDomainMap ( ) const
inlineprotected

Definition at line 240 of file Epetra_CrsSingletonFilter.h.

const Epetra_Map& Epetra_CrsSingletonFilter::FullMatrixRangeMap ( ) const
inlineprotected

Definition at line 241 of file Epetra_CrsSingletonFilter.h.

void Epetra_CrsSingletonFilter::InitializeDefaults ( )
protected

Definition at line 100 of file Epetra_CrsSingletonFilter.cpp.

int Epetra_CrsSingletonFilter::ComputeEliminateMaps ( )
protected
int Epetra_CrsSingletonFilter::Setup ( Epetra_LinearProblem Problem)
protected
int Epetra_CrsSingletonFilter::InitFullMatrixAccess ( )
protected

Definition at line 690 of file Epetra_CrsSingletonFilter.cpp.

int Epetra_CrsSingletonFilter::GetRow ( int  Row,
int &  NumIndices,
int *&  Indices 
)
protected

Definition at line 703 of file Epetra_CrsSingletonFilter.cpp.

int Epetra_CrsSingletonFilter::GetRowGCIDs ( int  Row,
int &  NumIndices,
double *&  Values,
int *&  GlobalIndices 
)
protected

Definition at line 731 of file Epetra_CrsSingletonFilter.cpp.

int Epetra_CrsSingletonFilter::GetRow ( int  Row,
int &  NumIndices,
double *&  Values,
int *&  Indices 
)
protected

Definition at line 716 of file Epetra_CrsSingletonFilter.cpp.

int Epetra_CrsSingletonFilter::CreatePostSolveArrays ( const Epetra_IntVector RowIDs,
const Epetra_MapColoring RowMapColors,
const Epetra_IntVector ColProfiles,
const Epetra_IntVector NewColProfiles,
const Epetra_IntVector ColHasRowWithSingleton 
)
protected

Definition at line 742 of file Epetra_CrsSingletonFilter.cpp.

int Epetra_CrsSingletonFilter::ConstructRedistributeExporter ( Epetra_Map SourceMap,
Epetra_Map TargetMap,
Epetra_Export *&  RedistributeExporter,
Epetra_Map *&  RedistributeMap 
)
protected

Definition at line 552 of file Epetra_CrsSingletonFilter.cpp.

Epetra_CrsSingletonFilter& Epetra_CrsSingletonFilter::operator= ( const Epetra_CrsSingletonFilter Problem)
private

Member Data Documentation

Epetra_LinearProblem* Epetra_CrsSingletonFilter::FullProblem_
protected

Definition at line 259 of file Epetra_CrsSingletonFilter.h.

Epetra_LinearProblem* Epetra_CrsSingletonFilter::ReducedProblem_
protected

Definition at line 260 of file Epetra_CrsSingletonFilter.h.

Epetra_RowMatrix* Epetra_CrsSingletonFilter::FullMatrix_
protected

Definition at line 261 of file Epetra_CrsSingletonFilter.h.

Epetra_CrsMatrix* Epetra_CrsSingletonFilter::FullCrsMatrix_
protected

Definition at line 262 of file Epetra_CrsSingletonFilter.h.

Epetra_CrsMatrix* Epetra_CrsSingletonFilter::ReducedMatrix_
protected

Definition at line 263 of file Epetra_CrsSingletonFilter.h.

Epetra_MultiVector* Epetra_CrsSingletonFilter::ReducedRHS_
protected

Definition at line 264 of file Epetra_CrsSingletonFilter.h.

Epetra_MultiVector* Epetra_CrsSingletonFilter::ReducedLHS_
protected

Definition at line 265 of file Epetra_CrsSingletonFilter.h.

Epetra_Map* Epetra_CrsSingletonFilter::ReducedMatrixRowMap_
protected

Definition at line 267 of file Epetra_CrsSingletonFilter.h.

Epetra_Map* Epetra_CrsSingletonFilter::ReducedMatrixColMap_
protected

Definition at line 268 of file Epetra_CrsSingletonFilter.h.

Epetra_Map* Epetra_CrsSingletonFilter::ReducedMatrixDomainMap_
protected

Definition at line 269 of file Epetra_CrsSingletonFilter.h.

Epetra_Map* Epetra_CrsSingletonFilter::ReducedMatrixRangeMap_
protected

Definition at line 270 of file Epetra_CrsSingletonFilter.h.

Epetra_Map* Epetra_CrsSingletonFilter::OrigReducedMatrixDomainMap_
protected

Definition at line 271 of file Epetra_CrsSingletonFilter.h.

Epetra_Import* Epetra_CrsSingletonFilter::Full2ReducedRHSImporter_
protected

Definition at line 272 of file Epetra_CrsSingletonFilter.h.

Epetra_Import* Epetra_CrsSingletonFilter::Full2ReducedLHSImporter_
protected

Definition at line 273 of file Epetra_CrsSingletonFilter.h.

Epetra_Export* Epetra_CrsSingletonFilter::RedistributeDomainExporter_
protected

Definition at line 274 of file Epetra_CrsSingletonFilter.h.

int* Epetra_CrsSingletonFilter::ColSingletonRowLIDs_
protected

Definition at line 276 of file Epetra_CrsSingletonFilter.h.

int* Epetra_CrsSingletonFilter::ColSingletonColLIDs_
protected

Definition at line 277 of file Epetra_CrsSingletonFilter.h.

int* Epetra_CrsSingletonFilter::ColSingletonPivotLIDs_
protected

Definition at line 278 of file Epetra_CrsSingletonFilter.h.

double* Epetra_CrsSingletonFilter::ColSingletonPivots_
protected

Definition at line 279 of file Epetra_CrsSingletonFilter.h.

int Epetra_CrsSingletonFilter::AbsoluteThreshold_
protected

Definition at line 282 of file Epetra_CrsSingletonFilter.h.

double Epetra_CrsSingletonFilter::RelativeThreshold_
protected

Definition at line 283 of file Epetra_CrsSingletonFilter.h.

int Epetra_CrsSingletonFilter::NumMyRowSingletons_
protected

Definition at line 285 of file Epetra_CrsSingletonFilter.h.

int Epetra_CrsSingletonFilter::NumMyColSingletons_
protected

Definition at line 286 of file Epetra_CrsSingletonFilter.h.

int Epetra_CrsSingletonFilter::NumGlobalRowSingletons_
protected

Definition at line 287 of file Epetra_CrsSingletonFilter.h.

int Epetra_CrsSingletonFilter::NumGlobalColSingletons_
protected

Definition at line 288 of file Epetra_CrsSingletonFilter.h.

double Epetra_CrsSingletonFilter::RatioOfDimensions_
protected

Definition at line 289 of file Epetra_CrsSingletonFilter.h.

double Epetra_CrsSingletonFilter::RatioOfNonzeros_
protected

Definition at line 290 of file Epetra_CrsSingletonFilter.h.

bool Epetra_CrsSingletonFilter::HaveReducedProblem_
protected

Definition at line 292 of file Epetra_CrsSingletonFilter.h.

bool Epetra_CrsSingletonFilter::UserDefinedEliminateMaps_
protected

Definition at line 293 of file Epetra_CrsSingletonFilter.h.

bool Epetra_CrsSingletonFilter::AnalysisDone_
protected

Definition at line 294 of file Epetra_CrsSingletonFilter.h.

bool Epetra_CrsSingletonFilter::SymmetricElimination_
protected

Definition at line 295 of file Epetra_CrsSingletonFilter.h.

Epetra_MultiVector* Epetra_CrsSingletonFilter::tempExportX_
protected

Definition at line 297 of file Epetra_CrsSingletonFilter.h.

Epetra_MultiVector* Epetra_CrsSingletonFilter::tempX_
protected

Definition at line 298 of file Epetra_CrsSingletonFilter.h.

Epetra_MultiVector* Epetra_CrsSingletonFilter::tempB_
protected

Definition at line 299 of file Epetra_CrsSingletonFilter.h.

Epetra_MultiVector* Epetra_CrsSingletonFilter::RedistributeReducedLHS_
protected

Definition at line 300 of file Epetra_CrsSingletonFilter.h.

int* Epetra_CrsSingletonFilter::Indices_
protected

Definition at line 301 of file Epetra_CrsSingletonFilter.h.

Epetra_SerialDenseVector Epetra_CrsSingletonFilter::Values_
protected

Definition at line 302 of file Epetra_CrsSingletonFilter.h.

Epetra_MapColoring* Epetra_CrsSingletonFilter::RowMapColors_
protected

Definition at line 304 of file Epetra_CrsSingletonFilter.h.

Epetra_MapColoring* Epetra_CrsSingletonFilter::ColMapColors_
protected

Definition at line 305 of file Epetra_CrsSingletonFilter.h.

bool Epetra_CrsSingletonFilter::FullMatrixIsCrsMatrix_
protected

Definition at line 306 of file Epetra_CrsSingletonFilter.h.

int Epetra_CrsSingletonFilter::MaxNumMyEntries_
protected

Definition at line 307 of file Epetra_CrsSingletonFilter.h.


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