EpetraExt
Development
Toggle main menu visibility
Loading...
Searching...
No Matches
src
transform
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
54
class
Epetra_LinearProblem
;
55
class
Epetra_Map
;
56
class
Epetra_MultiVector
;
57
class
Epetra_Import
;
58
class
Epetra_Export
;
59
class
Epetra_IntVector
;
60
61
namespace
EpetraExt
{
62
64
119
120
class
LinearProblem_CrsSingletonFilter
:
public
SameTypeTransform
<Epetra_LinearProblem> {
121
122
public
:
123
125
126
LinearProblem_CrsSingletonFilter
(
bool
verbose =
false
);
127
129
virtual
~LinearProblem_CrsSingletonFilter
();
131
133
NewTypeRef
operator()
(
OriginalTypeRef
orig );
134
136
bool
analyze
(
OriginalTypeRef
orig );
137
139
NewTypeRef
construct
();
140
142
bool
fwd
();
143
145
bool
rvs
();
146
148
158
int
Analyze
(
Epetra_RowMatrix
*
FullMatrix
);
159
161
bool
SingletonsDetected
()
const
{
if
(!
AnalysisDone_
)
return
(
false
);
else
return
(
NumSingletons
()>0);};
163
165
171
int
ConstructReducedProblem
(
Epetra_LinearProblem
* Problem);
172
174
179
int
UpdateReducedProblem
(
Epetra_LinearProblem
* Problem);
180
182
183
189
int
ComputeFullSolution
();
191
192
193
int
NumRowSingletons
()
const
{
return
(
NumGlobalRowSingletons_
);};
194
196
int
NumColSingletons
()
const
{
return
(
NumGlobalColSingletons_
);};
197
199
204
int
NumSingletons
()
const
{
return
(
NumColSingletons
()+
NumRowSingletons
());};
205
207
double
RatioOfDimensions
()
const
{
return
(
RatioOfDimensions_
);};
208
210
double
RatioOfNonzeros
()
const
{
return
(
RatioOfNonzeros_
);};
211
213
214
216
Epetra_LinearProblem
*
FullProblem
()
const
{
return
(
FullProblem_
);};
217
219
Epetra_LinearProblem
*
ReducedProblem
()
const
{
return
(
ReducedProblem_
.get());};
220
222
Epetra_RowMatrix
*
FullMatrix
()
const
{
return
(
FullMatrix_
);};
223
225
Epetra_CrsMatrix
*
ReducedMatrix
()
const
{
return
(
ReducedMatrix_
.get());};
226
228
Epetra_MapColoring
*
RowMapColors
()
const
{
return
(
RowMapColors_
);};
229
231
Epetra_MapColoring
*
ColMapColors
()
const
{
return
(
ColMapColors_
);};
232
234
Epetra_Map
*
ReducedMatrixRowMap
()
const
{
return
(
ReducedMatrixRowMap_
);};
235
237
Epetra_Map
*
ReducedMatrixColMap
()
const
{
return
(
ReducedMatrixColMap_
);};
238
240
Epetra_Map
*
ReducedMatrixDomainMap
()
const
{
return
(
ReducedMatrixDomainMap_
);};
241
243
Epetra_Map
*
ReducedMatrixRangeMap
()
const
{
return
(
ReducedMatrixRangeMap_
);};
245
246
protected
:
247
248
249
250
// This pointer will be zero if full matrix is not a CrsMatrix.
251
Epetra_CrsMatrix
*
FullCrsMatrix
()
const
{
return
(
FullCrsMatrix_
);};
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
();
258
int
ComputeEliminateMaps
();
259
int
Setup
(
Epetra_LinearProblem
* Problem);
260
int
InitFullMatrixAccess
();
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,
270
const
Epetra_MapColoring
&
RowMapColors
,
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
279
Epetra_LinearProblem
*
FullProblem_
;
280
Teuchos::RCP<Epetra_LinearProblem>
ReducedProblem_
;
281
Epetra_RowMatrix
*
FullMatrix_
;
282
Epetra_CrsMatrix
*
FullCrsMatrix_
;
283
Teuchos::RCP<Epetra_CrsMatrix>
ReducedMatrix_
;
284
Epetra_MultiVector
*
ReducedRHS_
;
285
Epetra_MultiVector
*
ReducedLHS_
;
286
287
Epetra_Map
*
ReducedMatrixRowMap_
;
288
Epetra_Map
*
ReducedMatrixColMap_
;
289
Epetra_Map
*
ReducedMatrixDomainMap_
;
290
Epetra_Map
*
ReducedMatrixRangeMap_
;
291
Epetra_Map
*
OrigReducedMatrixDomainMap_
;
292
Epetra_Import
*
Full2ReducedRHSImporter_
;
293
Epetra_Import
*
Full2ReducedLHSImporter_
;
294
Epetra_Export
*
RedistributeDomainExporter_
;
295
296
int
*
ColSingletonRowLIDs_
;
297
int
*
ColSingletonColLIDs_
;
298
int
*
ColSingletonPivotLIDs_
;
299
double
*
ColSingletonPivots_
;
300
301
302
int
AbsoluteThreshold_
;
303
double
RelativeThreshold_
;
304
305
int
NumMyRowSingletons_
;
306
int
NumMyColSingletons_
;
307
int
NumGlobalRowSingletons_
;
308
int
NumGlobalColSingletons_
;
309
double
RatioOfDimensions_
;
310
double
RatioOfNonzeros_
;
311
312
bool
HaveReducedProblem_
;
313
bool
UserDefinedEliminateMaps_
;
314
bool
AnalysisDone_
;
315
bool
SymmetricElimination_
;
316
317
Epetra_MultiVector
*
tempExportX_
;
318
Epetra_MultiVector
*
tempX_
;
319
Epetra_MultiVector
*
tempB_
;
320
Epetra_MultiVector
*
RedistributeReducedLHS_
;
321
int
*
Indices_int_
;
322
#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
323
long
long
*
Indices_LL_
;
324
#endif
325
Epetra_SerialDenseVector
Values_
;
326
327
Epetra_MapColoring
*
RowMapColors_
;
328
Epetra_MapColoring
*
ColMapColors_
;
329
bool
FullMatrixIsCrsMatrix_
;
330
int
MaxNumMyEntries_
;
331
332
bool
verbose_
;
333
334
335
private
:
337
LinearProblem_CrsSingletonFilter
(
const
LinearProblem_CrsSingletonFilter
&
/* Problem */
){};
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_ */
EpetraExt_Transform.h
EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrixDomainMap_
Epetra_Map * ReducedMatrixDomainMap_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:289
EpetraExt::LinearProblem_CrsSingletonFilter::FullMatrixIsCrsMatrix_
bool FullMatrixIsCrsMatrix_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:329
EpetraExt::LinearProblem_CrsSingletonFilter::ColSingletonColLIDs_
int * ColSingletonColLIDs_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:297
EpetraExt::LinearProblem_CrsSingletonFilter::ReducedRHS_
Epetra_MultiVector * ReducedRHS_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:284
EpetraExt::LinearProblem_CrsSingletonFilter::AnalysisDone_
bool AnalysisDone_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:314
EpetraExt::LinearProblem_CrsSingletonFilter::analyze
bool analyze(OriginalTypeRef orig)
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:109
EpetraExt::LinearProblem_CrsSingletonFilter::ColSingletonPivotLIDs_
int * ColSingletonPivotLIDs_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:298
EpetraExt::LinearProblem_CrsSingletonFilter::Values_
Epetra_SerialDenseVector Values_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:325
EpetraExt::LinearProblem_CrsSingletonFilter::NumRowSingletons
int NumRowSingletons() const
Return number of rows that contain a single entry, returns -1 if Analysis not performed yet.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:193
EpetraExt::LinearProblem_CrsSingletonFilter::FullMatrixDomainMap
const Epetra_Map & FullMatrixDomainMap() const
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:255
EpetraExt::LinearProblem_CrsSingletonFilter::RelativeThreshold_
double RelativeThreshold_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:303
EpetraExt::LinearProblem_CrsSingletonFilter::FullCrsMatrix
Epetra_CrsMatrix * FullCrsMatrix() const
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:251
EpetraExt::LinearProblem_CrsSingletonFilter::FullMatrixRangeMap
const Epetra_Map & FullMatrixRangeMap() const
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:256
EpetraExt::LinearProblem_CrsSingletonFilter::FullMatrix
Epetra_RowMatrix * FullMatrix() const
Returns pointer to Epetra_CrsMatrix from full problem.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:222
EpetraExt::LinearProblem_CrsSingletonFilter::ColMapColors_
Epetra_MapColoring * ColMapColors_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:328
EpetraExt::LinearProblem_CrsSingletonFilter::~LinearProblem_CrsSingletonFilter
virtual ~LinearProblem_CrsSingletonFilter()
Epetra_CrsSingletonFilter Destructor.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:67
EpetraExt::LinearProblem_CrsSingletonFilter::NumColSingletons
int NumColSingletons() const
Return number of columns that contain a single entry that are not associated with singleton row,...
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:196
EpetraExt::LinearProblem_CrsSingletonFilter::ReducedProblem
Epetra_LinearProblem * ReducedProblem() const
Returns pointer to the derived reduced Epetra_LinearProblem.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:219
EpetraExt::LinearProblem_CrsSingletonFilter::GetRowGCIDs
int GetRowGCIDs(int Row, int &NumIndices, double *&Values, int *&GlobalIndices)
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:861
EpetraExt::LinearProblem_CrsSingletonFilter::MaxNumMyEntries_
int MaxNumMyEntries_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:330
EpetraExt::LinearProblem_CrsSingletonFilter::rvs
bool rvs()
Reverse transfer of data from new object created in the operator() method call to the orig object inp...
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:196
EpetraExt::LinearProblem_CrsSingletonFilter::Setup
int Setup(Epetra_LinearProblem *Problem)
EpetraExt::LinearProblem_CrsSingletonFilter::OrigReducedMatrixDomainMap_
Epetra_Map * OrigReducedMatrixDomainMap_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:291
EpetraExt::LinearProblem_CrsSingletonFilter::ColMapColors
Epetra_MapColoring * ColMapColors() const
Returns pointer to Epetra_MapColoring object: color 0 columns are part of reduced system.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:231
EpetraExt::LinearProblem_CrsSingletonFilter::verbose_
bool verbose_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:332
EpetraExt::LinearProblem_CrsSingletonFilter::ComputeEliminateMaps
int ComputeEliminateMaps()
EpetraExt::LinearProblem_CrsSingletonFilter::GetRow
int GetRow(int Row, int &NumIndices, int *&Indices)
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:832
EpetraExt::LinearProblem_CrsSingletonFilter::HaveReducedProblem_
bool HaveReducedProblem_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:312
EpetraExt::LinearProblem_CrsSingletonFilter::Indices_LL_
long long * Indices_LL_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:323
EpetraExt::LinearProblem_CrsSingletonFilter::RedistributeReducedLHS_
Epetra_MultiVector * RedistributeReducedLHS_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:320
EpetraExt::LinearProblem_CrsSingletonFilter::NumGlobalColSingletons_
int NumGlobalColSingletons_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:308
EpetraExt::LinearProblem_CrsSingletonFilter::ReducedProblem_
Teuchos::RCP< Epetra_LinearProblem > ReducedProblem_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:280
EpetraExt::LinearProblem_CrsSingletonFilter::NumMyRowSingletons_
int NumMyRowSingletons_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:305
EpetraExt::LinearProblem_CrsSingletonFilter::FullMatrixColMap
const Epetra_Map & FullMatrixColMap() const
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:254
EpetraExt::LinearProblem_CrsSingletonFilter::InitFullMatrixAccess
int InitFullMatrixAccess()
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:815
EpetraExt::LinearProblem_CrsSingletonFilter::FullProblem
Epetra_LinearProblem * FullProblem() const
Returns pointer to the original unreduced Epetra_LinearProblem.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:216
EpetraExt::LinearProblem_CrsSingletonFilter::ConstructReducedProblem
int ConstructReducedProblem(Epetra_LinearProblem *Problem)
Return a reduced linear problem based on results of Analyze().
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:569
EpetraExt::LinearProblem_CrsSingletonFilter::tempB_
Epetra_MultiVector * tempB_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:319
EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrixRangeMap_
Epetra_Map * ReducedMatrixRangeMap_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:290
EpetraExt::LinearProblem_CrsSingletonFilter::ReducedLHS_
Epetra_MultiVector * ReducedLHS_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:285
EpetraExt::LinearProblem_CrsSingletonFilter::Indices_int_
int * Indices_int_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:321
EpetraExt::LinearProblem_CrsSingletonFilter::CreatePostSolveArrays
int CreatePostSolveArrays(const Epetra_IntVector &RowIDs, const Epetra_MapColoring &RowMapColors, const Epetra_IntVector &ColProfiles, const Epetra_IntVector &NewColProfiles, const Epetra_IntVector &ColHasRowWithSingleton)
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:885
EpetraExt::LinearProblem_CrsSingletonFilter::Full2ReducedLHSImporter_
Epetra_Import * Full2ReducedLHSImporter_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:293
EpetraExt::LinearProblem_CrsSingletonFilter::construct
NewTypeRef construct()
Construction of new object as a result of the transform.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:159
EpetraExt::LinearProblem_CrsSingletonFilter::UserDefinedEliminateMaps_
bool UserDefinedEliminateMaps_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:313
EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrixRowMap_
Epetra_Map * ReducedMatrixRowMap_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:287
EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrixRowMap
Epetra_Map * ReducedMatrixRowMap() const
Returns pointer to Epetra_Map describing the reduced system row distribution.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:234
EpetraExt::LinearProblem_CrsSingletonFilter::operator()
NewTypeRef operator()(OriginalTypeRef orig)
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:100
EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrix_
Teuchos::RCP< Epetra_CrsMatrix > ReducedMatrix_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:283
EpetraExt::LinearProblem_CrsSingletonFilter::FullCrsMatrix_
Epetra_CrsMatrix * FullCrsMatrix_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:282
EpetraExt::LinearProblem_CrsSingletonFilter::NumGlobalRowSingletons_
int NumGlobalRowSingletons_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:307
EpetraExt::LinearProblem_CrsSingletonFilter::RatioOfNonzeros
double RatioOfNonzeros() const
Returns ratio of reduced system to full system nonzero count, returns -1.0 if reduced problem not con...
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:210
EpetraExt::LinearProblem_CrsSingletonFilter::Full2ReducedRHSImporter_
Epetra_Import * Full2ReducedRHSImporter_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:292
EpetraExt::LinearProblem_CrsSingletonFilter::RowMapColors_
Epetra_MapColoring * RowMapColors_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:327
EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrixRangeMap
Epetra_Map * ReducedMatrixRangeMap() const
Returns pointer to Epetra_Map describing the range map for the reduced system.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:243
EpetraExt::LinearProblem_CrsSingletonFilter::Analyze
int Analyze(Epetra_RowMatrix *FullMatrix)
Analyze the input matrix, removing row/column pairs that have singletons.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:259
EpetraExt::LinearProblem_CrsSingletonFilter::RatioOfDimensions
double RatioOfDimensions() const
Returns ratio of reduced system to full system dimensions, returns -1.0 if reduced problem not constr...
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:207
EpetraExt::LinearProblem_CrsSingletonFilter::FullProblem_
Epetra_LinearProblem * FullProblem_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:279
EpetraExt::LinearProblem_CrsSingletonFilter::fwd
bool fwd()
Forward transfer of data from orig object input in the operator() method call to the new object creat...
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:187
EpetraExt::LinearProblem_CrsSingletonFilter::ColSingletonRowLIDs_
int * ColSingletonRowLIDs_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:296
EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrix
Epetra_CrsMatrix * ReducedMatrix() const
Returns pointer to Epetra_CrsMatrix from full problem.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:225
EpetraExt::LinearProblem_CrsSingletonFilter::SymmetricElimination_
bool SymmetricElimination_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:315
EpetraExt::LinearProblem_CrsSingletonFilter::LinearProblem_CrsSingletonFilter
LinearProblem_CrsSingletonFilter(bool verbose=false)
Epetra_CrsSingletonFilter default constructor.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:60
EpetraExt::LinearProblem_CrsSingletonFilter::NumMyColSingletons_
int NumMyColSingletons_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:306
EpetraExt::LinearProblem_CrsSingletonFilter::ComputeFullSolution
int ComputeFullSolution()
Compute a solution for the full problem using the solution of the reduced problem,...
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:768
EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrixColMap_
Epetra_Map * ReducedMatrixColMap_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:288
EpetraExt::LinearProblem_CrsSingletonFilter::RatioOfNonzeros_
double RatioOfNonzeros_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:310
EpetraExt::LinearProblem_CrsSingletonFilter::FullMatrix_
Epetra_RowMatrix * FullMatrix_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:281
EpetraExt::LinearProblem_CrsSingletonFilter::NumSingletons
int NumSingletons() const
Return total number of singletons detected, returns -1 if Analysis not performed yet.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:204
EpetraExt::LinearProblem_CrsSingletonFilter::RedistributeDomainExporter_
Epetra_Export * RedistributeDomainExporter_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:294
EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrixDomainMap
Epetra_Map * ReducedMatrixDomainMap() const
Returns pointer to Epetra_Map describing the domain map for the reduced system.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:240
EpetraExt::LinearProblem_CrsSingletonFilter::FullMatrixRowMap
const Epetra_Map & FullMatrixRowMap() const
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:253
EpetraExt::LinearProblem_CrsSingletonFilter::RowMapColors
Epetra_MapColoring * RowMapColors() const
Returns pointer to Epetra_MapColoring object: color 0 rows are part of reduced system.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:228
EpetraExt::LinearProblem_CrsSingletonFilter::UpdateReducedProblem
int UpdateReducedProblem(Epetra_LinearProblem *Problem)
Update a reduced linear problem using new values.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:692
EpetraExt::LinearProblem_CrsSingletonFilter::SingletonsDetected
bool SingletonsDetected() const
Returns true if singletons were detected in this matrix (must be called after Analyze() to be effecti...
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:161
EpetraExt::LinearProblem_CrsSingletonFilter::ColSingletonPivots_
double * ColSingletonPivots_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:299
EpetraExt::LinearProblem_CrsSingletonFilter::tempExportX_
Epetra_MultiVector * tempExportX_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:317
EpetraExt::LinearProblem_CrsSingletonFilter::ConstructRedistributeExporter
int ConstructRedistributeExporter(Epetra_Map *SourceMap, Epetra_Map *TargetMap, Epetra_Export *&RedistributeExporter, Epetra_Map *&RedistributeMap)
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:750
EpetraExt::LinearProblem_CrsSingletonFilter::ReducedMatrixColMap
Epetra_Map * ReducedMatrixColMap() const
Returns pointer to Epetra_Map describing the reduced system column distribution.
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:237
EpetraExt::LinearProblem_CrsSingletonFilter::RatioOfDimensions_
double RatioOfDimensions_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:309
EpetraExt::LinearProblem_CrsSingletonFilter::tempX_
Epetra_MultiVector * tempX_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:318
EpetraExt::LinearProblem_CrsSingletonFilter::AbsoluteThreshold_
int AbsoluteThreshold_
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.h:302
EpetraExt::LinearProblem_CrsSingletonFilter::InitializeDefaults
void InitializeDefaults()
Definition
EpetraExt_CrsSingletonFilter_LinearProblem.cpp:205
EpetraExt::SameTypeTransform
Definition
EpetraExt_Transform.h:271
EpetraExt::Transform< Epetra_LinearProblem, Epetra_LinearProblem >::OriginalTypeRef
Epetra_LinearProblem & OriginalTypeRef
Definition
EpetraExt_Transform.h:74
EpetraExt::Transform< Epetra_LinearProblem, Epetra_LinearProblem >::NewTypeRef
Epetra_LinearProblem & NewTypeRef
Definition
EpetraExt_Transform.h:79
Epetra_CrsMatrix
Epetra_Export
Epetra_Import
Epetra_IntVector
Epetra_LinearProblem
Epetra_MapColoring
Epetra_Map
Epetra_MultiVector
Epetra_RowMatrix
Epetra_SerialDenseVector
EpetraExt
EpetraExt::BlockCrsMatrix: A class for constructing a distributed block matrix.
Definition
EpetraExt_BlockCrsMatrix.cpp:46
Generated by
1.17.0