Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Epetra_CrsSingletonFilter.h
Go to the documentation of this file.
1 /*
2 //@HEADER
3 // ************************************************************************
4 //
5 // Epetra: Linear Algebra Services Package
6 // Copyright 2011 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39 //
40 // ************************************************************************
41 //@HEADER
42 */
43 
44 #ifndef EPETRA_CRSSINGLETONFILTER_H
45 #define EPETRA_CRSSINGLETONFILTER_H
46 
47 #include "Epetra_Object.h"
48 #include "Epetra_CrsMatrix.h"
49 #include "Epetra_MapColoring.h"
52 class Epetra_Map;
53 class Epetra_MultiVector;
54 class Epetra_Import;
55 class Epetra_Export;
56 class Epetra_IntVector;
57 
59 
116 
117  public:
118 
120 
123 
125  virtual ~Epetra_CrsSingletonFilter();
127 
129 
140 
142  bool SingletonsDetected() const {if (!AnalysisDone_) return(false); else return(RowMapColors_->MaxNumColors()>1);};
144 
146 
147 
154 
156 
162 
164 
166 
172  int ComputeFullSolution();
174 
176  int NumRowSingletons() const {return(NumGlobalRowSingletons_);};
178 
181 
183 
188  int NumSingletons() const {return(NumColSingletons()+NumRowSingletons());};
189 
191  double RatioOfDimensions() const {return(RatioOfDimensions_);};
192 
194  double RatioOfNonzeros() const {return(RatioOfNonzeros_);};
195 
197 
199 
202 
205 
208 
211 
214 
217 
220 
223 
226 
230 
231  protected:
232 
233 
234 
235  // This pointer will be zero if full matrix is not a CrsMatrix.
237 
238  const Epetra_Map & FullMatrixRowMap() const {return(FullMatrix()->RowMatrixRowMap());};
239  const Epetra_Map & FullMatrixColMap() const {return(FullMatrix()->RowMatrixColMap());};
240  const Epetra_Map & FullMatrixDomainMap() const {return((FullMatrix()->OperatorDomainMap()));};
241  const Epetra_Map & FullMatrixRangeMap() const {return((FullMatrix()->OperatorRangeMap()));};
242  void InitializeDefaults();
243  int ComputeEliminateMaps();
244  int Setup(Epetra_LinearProblem * Problem);
245  int InitFullMatrixAccess();
246  int GetRow(int Row, int & NumIndices, int * & Indices);
247  int GetRowGCIDs(int Row, int & NumIndices, double * & Values, int * & GlobalIndices);
248  int GetRow(int Row, int & NumIndices, double * & Values, int * & Indices);
249  int CreatePostSolveArrays(const Epetra_IntVector & RowIDs,
251  const Epetra_IntVector & ColProfiles,
252  const Epetra_IntVector & NewColProfiles,
253  const Epetra_IntVector & ColHasRowWithSingleton);
254 
255  int ConstructRedistributeExporter(Epetra_Map * SourceMap, Epetra_Map * TargetMap,
256  Epetra_Export * & RedistributeExporter,
257  Epetra_Map * & RedistributeMap);
258 
266 
275 
280 
281 
284 
291 
296 
301  int * Indices_;
303 
308 
309 
310  private:
314 };
315 #endif /* EPETRA_CRSSINGLETONFILTER_H */
Epetra_MultiVector: A class for constructing and using dense multi-vectors, vectors and matrices in p...
Epetra_Map: A class for partitioning vectors and matrices.
Definition: Epetra_Map.h:119
Epetra_MapColoring: A class for coloring Epetra_Map and Epetra_BlockMap objects.
Epetra_CrsSingletonFilter()
Epetra_CrsSingletonFilter default constructor.
int MaxNumColors() const
Returns maximum over all processors of the number of colors.
double RatioOfNonzeros() const
Returns ratio of reduced system to full system nonzero count, returns -1.0 if reduced problem not con...
int NumSingletons() const
Return total number of singletons detected, returns -1 if Analysis not performed yet.
Epetra_SerialDenseVector Values_
Epetra_IntVector: A class for constructing and using dense integer vectors on a parallel computer...
int UpdateReducedProblem(Epetra_LinearProblem *Problem)
Update a reduced linear problem using new values.
Epetra_MultiVector * RedistributeReducedLHS_
Epetra_LinearProblem * FullProblem_
virtual ~Epetra_CrsSingletonFilter()
Epetra_CrsSingletonFilter Destructor.
int ConstructRedistributeExporter(Epetra_Map *SourceMap, Epetra_Map *TargetMap, Epetra_Export *&RedistributeExporter, Epetra_Map *&RedistributeMap)
int Analyze(Epetra_RowMatrix *FullMatrix)
Analyze the input matrix, removing row/column pairs that have singletons.
const Epetra_Map & FullMatrixColMap() const
int ComputeFullSolution()
Compute a solution for the full problem using the solution of the reduced problem, put in LHS of FullProblem().
Epetra_Export: This class builds an export object for efficient exporting of off-processor elements...
Definition: Epetra_Export.h:62
Epetra_LinearProblem * ReducedProblem_
Epetra_Import: This class builds an import object for efficient importing of off-processor elements...
Definition: Epetra_Import.h:63
int NumRowSingletons() const
Return number of rows that contain a single entry, returns -1 if Analysis not performed yet...
const Epetra_Map & FullMatrixRangeMap() const
Epetra_SerialDenseVector: A class for constructing and using dense vectors.
Epetra_LinearProblem * ReducedProblem() const
Returns pointer to the derived reduced Epetra_LinearProblem.
Epetra_Map * ReducedMatrixDomainMap() const
Returns pointer to Epetra_Map describing the domain map for the reduced system.
int CreatePostSolveArrays(const Epetra_IntVector &RowIDs, const Epetra_MapColoring &RowMapColors, const Epetra_IntVector &ColProfiles, const Epetra_IntVector &NewColProfiles, const Epetra_IntVector &ColHasRowWithSingleton)
Epetra_Map * ReducedMatrixRowMap() const
Returns pointer to Epetra_Map describing the reduced system row distribution.
const Epetra_Map & FullMatrixDomainMap() const
Epetra_CrsSingletonFilter & operator=(const Epetra_CrsSingletonFilter &Problem)
int GetRow(int Row, int &NumIndices, int *&Indices)
Epetra_MapColoring * RowMapColors() const
Returns pointer to Epetra_MapColoring object: color 0 rows are part of reduced system.
Epetra_CrsSingletonFilter: A class for explicitly eliminating matrix rows and columns.
Epetra_Map * ReducedMatrixRangeMap() const
Returns pointer to Epetra_Map describing the range map for the reduced system.
Epetra_RowMatrix * FullMatrix() const
Returns pointer to Epetra_CrsMatrix from full problem.
bool SingletonsDetected() const
Returns true if singletons were detected in this matrix (must be called after Analyze() to be effecti...
Epetra_CrsMatrix * FullCrsMatrix() const
int NumColSingletons() const
Return number of columns that contain a single entry that are not associated with singleton row...
Epetra_Map * ReducedMatrixColMap() const
Returns pointer to Epetra_Map describing the reduced system column distribution.
double RatioOfDimensions() const
Returns ratio of reduced system to full system dimensions, returns -1.0 if reduced problem not constr...
int ConstructReducedProblem(Epetra_LinearProblem *Problem)
Return a reduced linear problem based on results of Analyze().
Epetra_CrsMatrix: A class for constructing and using real-valued double-precision sparse compressed r...
int Setup(Epetra_LinearProblem *Problem)
Epetra_LinearProblem * FullProblem() const
Returns pointer to the original unreduced Epetra_LinearProblem.
Epetra_CrsMatrix * ReducedMatrix() const
Returns pointer to Epetra_CrsMatrix from full problem.
Epetra_MapColoring * ColMapColors() const
Returns pointer to Epetra_MapColoring object: color 0 columns are part of reduced system...
Epetra_LinearProblem: The Epetra Linear Problem Class.
Epetra_RowMatrix: A pure virtual class for using real-valued double-precision row matrices...
int GetRowGCIDs(int Row, int &NumIndices, double *&Values, int *&GlobalIndices)
const Epetra_Map & FullMatrixRowMap() const