EpetraExt::LinearProblem_CrsSingletonFilter Class Reference

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

#include <EpetraExt_CrsSingletonFilter_LinearProblem.h>

Inheritance diagram for EpetraExt::LinearProblem_CrsSingletonFilter:

[legend]
List of all members.

Constructors/Destructor.

 LinearProblem_CrsSingletonFilter (bool verbose=false)
 Epetra_CrsSingletonFilter default constructor.
virtual ~LinearProblem_CrsSingletonFilter ()
 Epetra_CrsSingletonFilter Destructor.

Analyze methods.

int Analyze (Epetra_RowMatrix *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).

Reduce methods.

int ConstructReducedProblem (Epetra_LinearProblem *Problem)
 Return a reduced linear problem based on results of Analyze().
int UpdateReducedProblem (Epetra_LinearProblem *Problem)
 Update a reduced linear problem using new values.

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().

Filter Statistics.

int NumRowSingletons () const
 Return number of rows that contain a single entry, returns -1 if Analysis not performed yet.
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.
int NumSingletons () const
 Return total number of singletons detected, returns -1 if Analysis not performed yet.
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.

Attribute Access Methods.

Epetra_LinearProblemFullProblem () const
 Returns pointer to the original unreduced Epetra_LinearProblem.
Epetra_LinearProblemReducedProblem () const
 Returns pointer to the derived reduced Epetra_LinearProblem.
Epetra_RowMatrixFullMatrix () const
 Returns pointer to Epetra_CrsMatrix from full problem.
Epetra_CrsMatrixReducedMatrix () const
 Returns pointer to Epetra_CrsMatrix from full problem.
Epetra_MapColoringRowMapColors () const
 Returns pointer to Epetra_MapColoring object: color 0 rows are part of reduced system.
Epetra_MapColoringColMapColors () const
 Returns pointer to Epetra_MapColoring object: color 0 columns are part of reduced system.
Epetra_MapReducedMatrixRowMap () const
 Returns pointer to Epetra_Map describing the reduced system row distribution.
Epetra_MapReducedMatrixColMap () const
 Returns pointer to Epetra_Map describing the reduced system column distribution.
Epetra_MapReducedMatrixDomainMap () const
 Returns pointer to Epetra_Map describing the domain map for the reduced system.
Epetra_MapReducedMatrixRangeMap () const
 Returns pointer to Epetra_Map describing the range map for the reduced system.

Public Member Functions

NewTypeRef operator() (OriginalTypeRef orig)
bool analyze (OriginalTypeRef orig)
NewTypeRef construct ()
 Construction of new object as a result of the transform.
bool fwd ()
 Forward transfer of data from orig object input in the operator() method call to the new object created in this same call.
bool rvs ()
 Reverse transfer of data from new object created in the operator() method call to the orig object input to this same method.

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_
bool verbose_

Private Member Functions

 LinearProblem_CrsSingletonFilter (const LinearProblem_CrsSingletonFilter &Problem)
 Copy constructor (defined as private so it is unavailable to user).

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 105 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.


Constructor & Destructor Documentation

EpetraExt::LinearProblem_CrsSingletonFilter::LinearProblem_CrsSingletonFilter bool  verbose = false  ) 
 

Epetra_CrsSingletonFilter default constructor.

Definition at line 46 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

EpetraExt::LinearProblem_CrsSingletonFilter::~LinearProblem_CrsSingletonFilter  )  [virtual]
 

Epetra_CrsSingletonFilter Destructor.

Definition at line 53 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

EpetraExt::LinearProblem_CrsSingletonFilter::LinearProblem_CrsSingletonFilter const LinearProblem_CrsSingletonFilter Problem  )  [inline, private]
 

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

Definition at line 314 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.


Member Function Documentation

LinearProblem_CrsSingletonFilter::NewTypeRef EpetraExt::LinearProblem_CrsSingletonFilter::operator() OriginalTypeRef  orig  ) 
 

Definition at line 85 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

bool EpetraExt::LinearProblem_CrsSingletonFilter::analyze OriginalTypeRef  orig  ) 
 

Definition at line 94 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

LinearProblem_CrsSingletonFilter::NewTypeRef EpetraExt::LinearProblem_CrsSingletonFilter::construct  )  [virtual]
 

Construction of new object as a result of the transform.

Preconditions:

Invariants:

Postconditions:

The default implementation returns internal attribute newObj_.

Reimplemented from EpetraExt::Transform< T, U >.

Definition at line 130 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

bool EpetraExt::LinearProblem_CrsSingletonFilter::fwd  )  [virtual]
 

Forward transfer of data from orig object input in the operator() method call to the new object created in this same call.

Returns true is operation is successful.

Preconditions:

Invariants:

Postconditions:

Implements EpetraExt::Transform< T, U >.

Definition at line 152 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

bool EpetraExt::LinearProblem_CrsSingletonFilter::rvs  )  [virtual]
 

Reverse transfer of data from new object created in the operator() method call to the orig object input to this same method.

Returns true if operation is successful.

Preconditions:

Invariants:

Postconditions:

Implements EpetraExt::Transform< T, U >.

Definition at line 161 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

int EpetraExt::LinearProblem_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 223 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

bool EpetraExt::LinearProblem_CrsSingletonFilter::SingletonsDetected  )  const [inline]
 

Returns true if singletons were detected in this matrix (must be called after Analyze() to be effective).

Definition at line 146 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

int EpetraExt::LinearProblem_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 362 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

int EpetraExt::LinearProblem_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 538 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

int EpetraExt::LinearProblem_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 670 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

int EpetraExt::LinearProblem_CrsSingletonFilter::NumRowSingletons  )  const [inline]
 

Return number of rows that contain a single entry, returns -1 if Analysis not performed yet.

Definition at line 178 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

int EpetraExt::LinearProblem_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 181 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

int EpetraExt::LinearProblem_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 189 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

double EpetraExt::LinearProblem_CrsSingletonFilter::RatioOfDimensions  )  const [inline]
 

Returns ratio of reduced system to full system dimensions, returns -1.0 if reduced problem not constructed.

Definition at line 192 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

double EpetraExt::LinearProblem_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 195 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_LinearProblem* EpetraExt::LinearProblem_CrsSingletonFilter::FullProblem  )  const [inline]
 

Returns pointer to the original unreduced Epetra_LinearProblem.

Definition at line 201 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_LinearProblem* EpetraExt::LinearProblem_CrsSingletonFilter::ReducedProblem  )  const [inline]
 

Returns pointer to the derived reduced Epetra_LinearProblem.

Definition at line 204 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_RowMatrix* EpetraExt::LinearProblem_CrsSingletonFilter::FullMatrix  )  const [inline]
 

Returns pointer to Epetra_CrsMatrix from full problem.

Definition at line 207 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_CrsMatrix* EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrix  )  const [inline]
 

Returns pointer to Epetra_CrsMatrix from full problem.

Definition at line 210 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_MapColoring* EpetraExt::LinearProblem_CrsSingletonFilter::RowMapColors  )  const [inline]
 

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

Definition at line 213 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_MapColoring* EpetraExt::LinearProblem_CrsSingletonFilter::ColMapColors  )  const [inline]
 

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

Definition at line 216 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_Map* EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrixRowMap  )  const [inline]
 

Returns pointer to Epetra_Map describing the reduced system row distribution.

Definition at line 219 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_Map* EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrixColMap  )  const [inline]
 

Returns pointer to Epetra_Map describing the reduced system column distribution.

Definition at line 222 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_Map* EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrixDomainMap  )  const [inline]
 

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

Definition at line 225 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_Map* EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrixRangeMap  )  const [inline]
 

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

Definition at line 228 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_CrsMatrix* EpetraExt::LinearProblem_CrsSingletonFilter::FullCrsMatrix  )  const [inline, protected]
 

Definition at line 236 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

const Epetra_Map& EpetraExt::LinearProblem_CrsSingletonFilter::FullMatrixRowMap  )  const [inline, protected]
 

Definition at line 238 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

const Epetra_Map& EpetraExt::LinearProblem_CrsSingletonFilter::FullMatrixColMap  )  const [inline, protected]
 

Definition at line 239 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

const Epetra_Map& EpetraExt::LinearProblem_CrsSingletonFilter::FullMatrixDomainMap  )  const [inline, protected]
 

Definition at line 240 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

const Epetra_Map& EpetraExt::LinearProblem_CrsSingletonFilter::FullMatrixRangeMap  )  const [inline, protected]
 

Definition at line 241 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

void EpetraExt::LinearProblem_CrsSingletonFilter::InitializeDefaults  )  [protected]
 

Definition at line 170 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

int EpetraExt::LinearProblem_CrsSingletonFilter::ComputeEliminateMaps  )  [protected]
 

int EpetraExt::LinearProblem_CrsSingletonFilter::Setup Epetra_LinearProblem Problem  )  [protected]
 

int EpetraExt::LinearProblem_CrsSingletonFilter::InitFullMatrixAccess  )  [protected]
 

Definition at line 717 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

int EpetraExt::LinearProblem_CrsSingletonFilter::GetRow int  Row,
int &  NumIndices,
int *&  Indices
[protected]
 

Definition at line 730 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

int EpetraExt::LinearProblem_CrsSingletonFilter::GetRowGCIDs int  Row,
int &  NumIndices,
double *&  Values,
int *&  GlobalIndices
[protected]
 

Definition at line 758 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

int EpetraExt::LinearProblem_CrsSingletonFilter::GetRow int  Row,
int &  NumIndices,
double *&  Values,
int *&  Indices
[protected]
 

Definition at line 743 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

int EpetraExt::LinearProblem_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 769 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.

int EpetraExt::LinearProblem_CrsSingletonFilter::ConstructRedistributeExporter Epetra_Map SourceMap,
Epetra_Map TargetMap,
Epetra_Export *&  RedistributeExporter,
Epetra_Map *&  RedistributeMap
[protected]
 

Definition at line 631 of file EpetraExt_CrsSingletonFilter_LinearProblem.cpp.


Member Data Documentation

Epetra_LinearProblem* EpetraExt::LinearProblem_CrsSingletonFilter::FullProblem_ [protected]
 

Definition at line 259 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_LinearProblem* EpetraExt::LinearProblem_CrsSingletonFilter::ReducedProblem_ [protected]
 

Definition at line 260 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_RowMatrix* EpetraExt::LinearProblem_CrsSingletonFilter::FullMatrix_ [protected]
 

Definition at line 261 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_CrsMatrix* EpetraExt::LinearProblem_CrsSingletonFilter::FullCrsMatrix_ [protected]
 

Definition at line 262 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_CrsMatrix* EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrix_ [protected]
 

Definition at line 263 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_MultiVector* EpetraExt::LinearProblem_CrsSingletonFilter::ReducedRHS_ [protected]
 

Definition at line 264 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_MultiVector* EpetraExt::LinearProblem_CrsSingletonFilter::ReducedLHS_ [protected]
 

Definition at line 265 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_Map* EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrixRowMap_ [protected]
 

Definition at line 267 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_Map* EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrixColMap_ [protected]
 

Definition at line 268 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_Map* EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrixDomainMap_ [protected]
 

Definition at line 269 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_Map* EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrixRangeMap_ [protected]
 

Definition at line 270 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_Map* EpetraExt::LinearProblem_CrsSingletonFilter::OrigReducedMatrixDomainMap_ [protected]
 

Definition at line 271 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_Import* EpetraExt::LinearProblem_CrsSingletonFilter::Full2ReducedRHSImporter_ [protected]
 

Definition at line 272 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_Import* EpetraExt::LinearProblem_CrsSingletonFilter::Full2ReducedLHSImporter_ [protected]
 

Definition at line 273 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_Export* EpetraExt::LinearProblem_CrsSingletonFilter::RedistributeDomainExporter_ [protected]
 

Definition at line 274 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

int* EpetraExt::LinearProblem_CrsSingletonFilter::ColSingletonRowLIDs_ [protected]
 

Definition at line 276 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

int* EpetraExt::LinearProblem_CrsSingletonFilter::ColSingletonColLIDs_ [protected]
 

Definition at line 277 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

int* EpetraExt::LinearProblem_CrsSingletonFilter::ColSingletonPivotLIDs_ [protected]
 

Definition at line 278 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

double* EpetraExt::LinearProblem_CrsSingletonFilter::ColSingletonPivots_ [protected]
 

Definition at line 279 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

int EpetraExt::LinearProblem_CrsSingletonFilter::AbsoluteThreshold_ [protected]
 

Definition at line 282 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

double EpetraExt::LinearProblem_CrsSingletonFilter::RelativeThreshold_ [protected]
 

Definition at line 283 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

int EpetraExt::LinearProblem_CrsSingletonFilter::NumMyRowSingletons_ [protected]
 

Definition at line 285 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

int EpetraExt::LinearProblem_CrsSingletonFilter::NumMyColSingletons_ [protected]
 

Definition at line 286 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

int EpetraExt::LinearProblem_CrsSingletonFilter::NumGlobalRowSingletons_ [protected]
 

Definition at line 287 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

int EpetraExt::LinearProblem_CrsSingletonFilter::NumGlobalColSingletons_ [protected]
 

Definition at line 288 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

double EpetraExt::LinearProblem_CrsSingletonFilter::RatioOfDimensions_ [protected]
 

Definition at line 289 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

double EpetraExt::LinearProblem_CrsSingletonFilter::RatioOfNonzeros_ [protected]
 

Definition at line 290 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

bool EpetraExt::LinearProblem_CrsSingletonFilter::HaveReducedProblem_ [protected]
 

Definition at line 292 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

bool EpetraExt::LinearProblem_CrsSingletonFilter::UserDefinedEliminateMaps_ [protected]
 

Definition at line 293 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

bool EpetraExt::LinearProblem_CrsSingletonFilter::AnalysisDone_ [protected]
 

Definition at line 294 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

bool EpetraExt::LinearProblem_CrsSingletonFilter::SymmetricElimination_ [protected]
 

Definition at line 295 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_MultiVector* EpetraExt::LinearProblem_CrsSingletonFilter::tempExportX_ [protected]
 

Definition at line 297 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_MultiVector* EpetraExt::LinearProblem_CrsSingletonFilter::tempX_ [protected]
 

Definition at line 298 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_MultiVector* EpetraExt::LinearProblem_CrsSingletonFilter::tempB_ [protected]
 

Definition at line 299 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_MultiVector* EpetraExt::LinearProblem_CrsSingletonFilter::RedistributeReducedLHS_ [protected]
 

Definition at line 300 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

int* EpetraExt::LinearProblem_CrsSingletonFilter::Indices_ [protected]
 

Definition at line 301 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_SerialDenseVector EpetraExt::LinearProblem_CrsSingletonFilter::Values_ [protected]
 

Definition at line 302 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_MapColoring* EpetraExt::LinearProblem_CrsSingletonFilter::RowMapColors_ [protected]
 

Definition at line 304 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

Epetra_MapColoring* EpetraExt::LinearProblem_CrsSingletonFilter::ColMapColors_ [protected]
 

Definition at line 305 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

bool EpetraExt::LinearProblem_CrsSingletonFilter::FullMatrixIsCrsMatrix_ [protected]
 

Definition at line 306 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

int EpetraExt::LinearProblem_CrsSingletonFilter::MaxNumMyEntries_ [protected]
 

Definition at line 307 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.

bool EpetraExt::LinearProblem_CrsSingletonFilter::verbose_ [protected]
 

Definition at line 309 of file EpetraExt_CrsSingletonFilter_LinearProblem.h.


The documentation for this class was generated from the following files:
Generated on Thu Sep 18 12:32:06 2008 for EpetraExt Package Browser (Single Doxygen Collection) by doxygen 1.3.9.1