EpetraExt Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
EpetraExt_CrsSingletonFilter_LinearProblem.h
Go to the documentation of this file.
1 //@HEADER
2 // ***********************************************************************
3 //
4 // EpetraExt: Epetra Extended - Linear Algebra Services Package
5 // Copyright (2011) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 
42 #ifndef _EpetraExt_LINEARPROBLEM_CRSSINGLETONFILTER_H_
43 #define _EpetraExt_LINEARPROBLEM_CRSSINGLETONFILTER_H_
44 
45 #if defined(EpetraExt_SHOW_DEPRECATED_WARNINGS)
46 #ifdef __GNUC__
47 #warning "The EpetraExt package is deprecated"
48 #endif
49 #endif
50 
51 #include "Epetra_ConfigDefs.h"
52 #include "Epetra_Object.h"
53 #include "Epetra_CrsMatrix.h"
54 #include "Epetra_MapColoring.h"
56 
57 #include "EpetraExt_Transform.h"
58 #include "Teuchos_RCP.hpp"
59 
61 class Epetra_Map;
62 class Epetra_MultiVector;
63 class Epetra_Import;
64 class Epetra_Export;
65 class Epetra_IntVector;
66 
67 namespace EpetraExt {
68 
70 
126 class LinearProblem_CrsSingletonFilter : public SameTypeTransform<Epetra_LinearProblem> {
127 
128  public:
129 
131  LinearProblem_CrsSingletonFilter( bool verbose = false );
133 
137 
140 
142  bool analyze( OriginalTypeRef orig );
143 
146 
148  bool fwd();
149 
151  bool rvs();
152 
154 
165 
167  bool SingletonsDetected() const {if (!AnalysisDone_) return(false); else return(NumSingletons()>0);};
169 
171 
178 
180 
186 
188 
189 
195  int ComputeFullSolution();
197 
198  int NumRowSingletons() const {return(NumGlobalRowSingletons_);};
200 
203 
205 
210  int NumSingletons() const {return(NumColSingletons()+NumRowSingletons());};
211 
213  double RatioOfDimensions() const {return(RatioOfDimensions_);};
214 
216  double RatioOfNonzeros() const {return(RatioOfNonzeros_);};
217 
219 
220 
223 
226 
229 
232 
235 
238 
241 
244 
247 
251 
252  protected:
253 
254 
255 
256  // This pointer will be zero if full matrix is not a CrsMatrix.
258 
259  const Epetra_Map & FullMatrixRowMap() const {return(FullMatrix()->RowMatrixRowMap());};
260  const Epetra_Map & FullMatrixColMap() const {return(FullMatrix()->RowMatrixColMap());};
261  const Epetra_Map & FullMatrixDomainMap() const {return((FullMatrix()->OperatorDomainMap()));};
262  const Epetra_Map & FullMatrixRangeMap() const {return((FullMatrix()->OperatorRangeMap()));};
263  void InitializeDefaults();
264  int ComputeEliminateMaps();
265  int Setup(Epetra_LinearProblem * Problem);
266  int InitFullMatrixAccess();
267  int GetRow(int Row, int & NumIndices, int * & Indices);
268  int GetRow(int Row, int & NumIndices, double * & Values, int * & Indices);
269 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
270  int GetRowGCIDs(int Row, int & NumIndices, double * & Values, int * & GlobalIndices);
271 #endif
272 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
273  int GetRowGCIDs(int Row, int & NumIndices, double * & Values, long long * & GlobalIndices);
274 #endif
275  int CreatePostSolveArrays(const Epetra_IntVector & RowIDs,
277  const Epetra_IntVector & ColProfiles,
278  const Epetra_IntVector & NewColProfiles,
279  const Epetra_IntVector & ColHasRowWithSingleton);
280 
281  int ConstructRedistributeExporter(Epetra_Map * SourceMap, Epetra_Map * TargetMap,
282  Epetra_Export * & RedistributeExporter,
283  Epetra_Map * & RedistributeMap);
284 
292 
301 
306 
307 
310 
317 
322 
328 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
329  long long * Indices_LL_;
330 #endif
332 
337 
338  bool verbose_;
339 
340 
341  private:
344 
345  template<typename int_type>
347 
348  template<typename int_type>
350 
351  template<typename int_type>
352  int TConstructRedistributeExporter(Epetra_Map * SourceMap, Epetra_Map * TargetMap,
353  Epetra_Export * & RedistributeExporter,
354  Epetra_Map * & RedistributeMap);
355 };
356 
357 } //namespace EpetraExt
358 
359 #endif /* _EpetraExt_LINEARPROBLEM_CRSSINGLETONFILTER_H_ */
int GetRowGCIDs(int Row, int &NumIndices, double *&Values, int *&GlobalIndices)
Epetra_Map * ReducedMatrixDomainMap() const
Returns pointer to Epetra_Map describing the domain map for the reduced system.
Epetra_Map * ReducedMatrixRangeMap() const
Returns pointer to Epetra_Map describing the range map for the reduced system.
double RatioOfDimensions() const
Returns ratio of reduced system to full system dimensions, returns -1.0 if reduced problem not constr...
LinearProblem_CrsSingletonFilter(const LinearProblem_CrsSingletonFilter &)
Copy constructor (defined as private so it is unavailable to user).
NewTypeRef operator()(OriginalTypeRef orig)
Analysis of transform operation on original object and construction of new object.
bool fwd()
Forward transfer of data from orig object input in the operator() method call to the new object creat...
NewTypeRef construct()
Construction of new object as a result of the transform.
virtual ~LinearProblem_CrsSingletonFilter()
Epetra_CrsSingletonFilter Destructor.
bool rvs()
Reverse transfer of data from new object created in the operator() method call to the orig object inp...
T * get() const
bool SingletonsDetected() const
Returns true if singletons were detected in this matrix (must be called after Analyze() to be effecti...
int Analyze(Epetra_RowMatrix *FullMatrix)
Analyze the input matrix, removing row/column pairs that have singletons.
int Setup(Epetra_LinearProblem *Problem)
Epetra_CrsSingletonFilter: A class for explicitly eliminating matrix rows and columns.
Epetra_Map * ReducedMatrixRowMap() const
Returns pointer to Epetra_Map describing the reduced system row distribution.
Epetra_CrsMatrix * ReducedMatrix() const
Returns pointer to Epetra_CrsMatrix from full problem.
Epetra_MapColoring * RowMapColors() const
Returns pointer to Epetra_MapColoring object: color 0 rows are part of reduced system.
int ComputeFullSolution()
Compute a solution for the full problem using the solution of the reduced problem, put in LHS of FullProblem().
LinearProblem_CrsSingletonFilter(bool verbose=false)
Epetra_CrsSingletonFilter default constructor.
Epetra_LinearProblem * ReducedProblem() const
Returns pointer to the derived reduced Epetra_LinearProblem.
bool analyze(OriginalTypeRef orig)
Initial analysis phase of transform.
int NumRowSingletons() const
Return number of rows that contain a single entry, returns -1 if Analysis not performed yet...
double RatioOfNonzeros() const
Returns ratio of reduced system to full system nonzero count, returns -1.0 if reduced problem not con...
int ConstructRedistributeExporter(Epetra_Map *SourceMap, Epetra_Map *TargetMap, Epetra_Export *&RedistributeExporter, Epetra_Map *&RedistributeMap)
int CreatePostSolveArrays(const Epetra_IntVector &RowIDs, const Epetra_MapColoring &RowMapColors, const Epetra_IntVector &ColProfiles, const Epetra_IntVector &NewColProfiles, const Epetra_IntVector &ColHasRowWithSingleton)
Epetra_RowMatrix * FullMatrix() const
Returns pointer to Epetra_CrsMatrix from full problem.
int TConstructRedistributeExporter(Epetra_Map *SourceMap, Epetra_Map *TargetMap, Epetra_Export *&RedistributeExporter, Epetra_Map *&RedistributeMap)
int NumColSingletons() const
Return number of columns that contain a single entry that are not associated with singleton row...
int UpdateReducedProblem(Epetra_LinearProblem *Problem)
Update a reduced linear problem using new values.
Epetra_Map * ReducedMatrixColMap() const
Returns pointer to Epetra_Map describing the reduced system column distribution.
int ConstructReducedProblem(Epetra_LinearProblem *Problem)
Return a reduced linear problem based on results of Analyze().
Epetra_LinearProblem * FullProblem() const
Returns pointer to the original unreduced Epetra_LinearProblem.
int NumSingletons() const
Return total number of singletons detected, returns -1 if Analysis not performed yet.
Epetra_MapColoring * ColMapColors() const
Returns pointer to Epetra_MapColoring object: color 0 columns are part of reduced system...