EpetraExt Development
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#include "Epetra_ConfigDefs.h"
46#include "Epetra_Object.h"
47#include "Epetra_CrsMatrix.h"
48#include "Epetra_MapColoring.h"
49#include "Epetra_SerialDenseVector.h"
50
51#include "EpetraExt_Transform.h"
52#include "Teuchos_RCP.hpp"
53
55class Epetra_Map;
57class Epetra_Import;
58class Epetra_Export;
60
61namespace EpetraExt {
62
64
120class LinearProblem_CrsSingletonFilter : public SameTypeTransform<Epetra_LinearProblem> {
121
122 public:
123
125
126 LinearProblem_CrsSingletonFilter( bool verbose = false );
127
131
134
136 bool analyze( OriginalTypeRef orig );
137
140
142 bool fwd();
143
145 bool rvs();
146
148
159
161 bool SingletonsDetected() const {if (!AnalysisDone_) return(false); else return(NumSingletons()>0);};
163
165
172
174
180
182
183
191
192
194
197
199
205
207 double RatioOfDimensions() const {return(RatioOfDimensions_);};
208
210 double RatioOfNonzeros() const {return(RatioOfNonzeros_);};
211
213
214
217
220
223
226
229
232
235
238
241
245
246 protected:
247
248
249
250 // This pointer will be zero if full matrix is not a CrsMatrix.
252
253 const Epetra_Map & FullMatrixRowMap() const {return(FullMatrix()->RowMatrixRowMap());};
254 const Epetra_Map & FullMatrixColMap() const {return(FullMatrix()->RowMatrixColMap());};
255 const Epetra_Map & FullMatrixDomainMap() const {return((FullMatrix()->OperatorDomainMap()));};
256 const Epetra_Map & FullMatrixRangeMap() const {return((FullMatrix()->OperatorRangeMap()));};
257 void InitializeDefaults();
261 int GetRow(int Row, int & NumIndices, int * & Indices);
262 int GetRow(int Row, int & NumIndices, double * & Values, int * & Indices);
263#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
264 int GetRowGCIDs(int Row, int & NumIndices, double * & Values, int * & GlobalIndices);
265#endif
266#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
267 int GetRowGCIDs(int Row, int & NumIndices, double * & Values, long long * & GlobalIndices);
268#endif
269 int CreatePostSolveArrays(const Epetra_IntVector & RowIDs,
271 const Epetra_IntVector & ColProfiles,
272 const Epetra_IntVector & NewColProfiles,
273 const Epetra_IntVector & ColHasRowWithSingleton);
274
275 int ConstructRedistributeExporter(Epetra_Map * SourceMap, Epetra_Map * TargetMap,
276 Epetra_Export * & RedistributeExporter,
277 Epetra_Map * & RedistributeMap);
278
280 Teuchos::RCP<Epetra_LinearProblem> ReducedProblem_;
283 Teuchos::RCP<Epetra_CrsMatrix> ReducedMatrix_;
286
295
300
301
304
311
316
322#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
323 long long * Indices_LL_;
324#endif
326
331
333
334
335 private:
338
339 template<typename int_type>
340 int TConstructReducedProblem(Epetra_LinearProblem * Problem);
341
342 template<typename int_type>
343 int TUpdateReducedProblem(Epetra_LinearProblem * Problem);
344
345 template<typename int_type>
346 int TConstructRedistributeExporter(Epetra_Map * SourceMap, Epetra_Map * TargetMap,
347 Epetra_Export * & RedistributeExporter,
348 Epetra_Map * & RedistributeMap);
349};
350
351} //namespace EpetraExt
352
353#endif /* _EpetraExt_LINEARPROBLEM_CRSSINGLETONFILTER_H_ */
Epetra_CrsSingletonFilter: A class for explicitly eliminating matrix rows and columns.
int NumRowSingletons() const
Return number of rows that contain a single entry, returns -1 if Analysis not performed yet.
Epetra_RowMatrix * FullMatrix() const
Returns pointer to Epetra_CrsMatrix from full problem.
virtual ~LinearProblem_CrsSingletonFilter()
Epetra_CrsSingletonFilter Destructor.
int NumColSingletons() const
Return number of columns that contain a single entry that are not associated with singleton row,...
Epetra_LinearProblem * ReducedProblem() const
Returns pointer to the derived reduced Epetra_LinearProblem.
int GetRowGCIDs(int Row, int &NumIndices, double *&Values, int *&GlobalIndices)
bool rvs()
Reverse transfer of data from new object created in the operator() method call to the orig object inp...
int Setup(Epetra_LinearProblem *Problem)
Epetra_MapColoring * ColMapColors() const
Returns pointer to Epetra_MapColoring object: color 0 columns are part of reduced system.
Epetra_LinearProblem * FullProblem() const
Returns pointer to the original unreduced Epetra_LinearProblem.
int ConstructReducedProblem(Epetra_LinearProblem *Problem)
Return a reduced linear problem based on results of Analyze().
int CreatePostSolveArrays(const Epetra_IntVector &RowIDs, const Epetra_MapColoring &RowMapColors, const Epetra_IntVector &ColProfiles, const Epetra_IntVector &NewColProfiles, const Epetra_IntVector &ColHasRowWithSingleton)
NewTypeRef construct()
Construction of new object as a result of the transform.
Epetra_Map * ReducedMatrixRowMap() const
Returns pointer to Epetra_Map describing the reduced system row distribution.
double RatioOfNonzeros() const
Returns ratio of reduced system to full system nonzero count, returns -1.0 if reduced problem not con...
Epetra_Map * ReducedMatrixRangeMap() const
Returns pointer to Epetra_Map describing the range map for the reduced system.
int Analyze(Epetra_RowMatrix *FullMatrix)
Analyze the input matrix, removing row/column pairs that have singletons.
double RatioOfDimensions() const
Returns ratio of reduced system to full system dimensions, returns -1.0 if reduced problem not constr...
bool fwd()
Forward transfer of data from orig object input in the operator() method call to the new object creat...
Epetra_CrsMatrix * ReducedMatrix() const
Returns pointer to Epetra_CrsMatrix from full problem.
LinearProblem_CrsSingletonFilter(bool verbose=false)
Epetra_CrsSingletonFilter default constructor.
int ComputeFullSolution()
Compute a solution for the full problem using the solution of the reduced problem,...
int NumSingletons() const
Return total number of singletons detected, returns -1 if Analysis not performed yet.
Epetra_Map * ReducedMatrixDomainMap() const
Returns pointer to Epetra_Map describing the domain map for the reduced system.
Epetra_MapColoring * RowMapColors() const
Returns pointer to Epetra_MapColoring object: color 0 rows are part of reduced system.
int UpdateReducedProblem(Epetra_LinearProblem *Problem)
Update a reduced linear problem using new values.
bool SingletonsDetected() const
Returns true if singletons were detected in this matrix (must be called after Analyze() to be effecti...
int ConstructRedistributeExporter(Epetra_Map *SourceMap, Epetra_Map *TargetMap, Epetra_Export *&RedistributeExporter, Epetra_Map *&RedistributeMap)
Epetra_Map * ReducedMatrixColMap() const
Returns pointer to Epetra_Map describing the reduced system column distribution.
EpetraExt::BlockCrsMatrix: A class for constructing a distributed block matrix.