IFPACK Development
Loading...
Searching...
No Matches
Ifpack_OverlapSolveObject.h
1/*@HEADER
2// ***********************************************************************
3//
4// Ifpack: Object-Oriented Algebraic Preconditioner Package
5// Copyright (2002) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
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
43#ifndef IFPACK_SOLVEOBJECT_H
44#define IFPACK_SOLVEOBJECT_H
45
46#if defined(Ifpack_SHOW_DEPRECATED_WARNINGS)
47#ifdef __GNUC__
48#warning "The Ifpack package is deprecated"
49#endif
50#endif
51
52#include "Epetra_Operator.h"
53#include "Epetra_CrsMatrix.h"
54#include "Epetra_Vector.h"
55#include "Epetra_CombineMode.h"
56class Epetra_Flops;
57
59
61
62 public:
64
66
69
73
75
78
80 int SetLowerOperator (Epetra_CrsMatrix * L, bool UseLTrans) {L_ = L; UseLTrans_ = UseLTrans; return(0);};
81
83 int SetDiagonal (Epetra_Vector * D, bool UseDInv) {D_ = D; UseDInv_ = UseDInv; return(0);};
84
86 int SetUpperOperator (Epetra_CrsMatrix * U, bool UseUTrans) {U_ = U; UseUTrans_ = UseUTrans; return(0);};
88
90
91
93
103 int Solve(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
104
106
116 int Multiply(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
117
119
127 int Condest(bool Trans, double & ConditionNumberEstimate) const;
129
131
133
138 Epetra_CombineMode OverlapMode() const {return(OverlapMode_);}
139
142
144 int NumMyNonzeros() const {return(L().NumMyNonzeros()+U().NumMyNonzeros());};
145
147 const Epetra_CrsMatrix & L() const {return(*L_);};
148
150 const Epetra_Vector & D() const {return(*D_);};
151
153 const Epetra_CrsMatrix & U() const {return(*U_);};
155
157
159 const char * Label() const {return(Label_);};
160
162
171 int SetUseTranspose(bool UseTranspose) {UseTranspose_ = UseTranspose; return(0);};
172
174
187
189
204
206 double NormInf() const {return(0.0);};
207
209 bool HasNormInf() const {return(false);};
210
212 bool UseTranspose() const {return(UseTranspose_);};
213
215 const Epetra_Map & OperatorDomainMap() const {return(U_->OperatorDomainMap());};
216
218 const Epetra_Map & OperatorRangeMap() const{return(L_->OperatorRangeMap());};
219
221 const Epetra_Comm & Comm() const{return(Comm_);};
223
224 protected:
225 virtual int SetupXY(bool Trans,
226 const Epetra_MultiVector& Xin, const Epetra_MultiVector& Yin,
227 Epetra_MultiVector * & Xout, Epetra_MultiVector * & Yout) const=0;
228 private:
229 char * Label_;
230 Epetra_CrsMatrix * L_;
231 bool UseLTrans_;
232 Epetra_Vector * D_;
233 bool UseDInv_;
234 Epetra_CrsMatrix * U_;
235 bool UseUTrans_;
236 bool UseTranspose_;
237 const Epetra_Comm & Comm_;
238 mutable double Condest_;
239 Epetra_Flops * Counter_;
240 Epetra_CombineMode OverlapMode_;
241
242};
243#endif // IFPACK_SOLVEOBJECT_H
Epetra_CombineMode
const Epetra_Map & OperatorRangeMap() const
const Epetra_Map & OperatorDomainMap() const
Ifpack_OverlapSolveObject: Provides Overlapped Forward/back solve services for Ifpack.
const Epetra_CrsMatrix & L() const
Returns the address of the L factor associated with this factored matrix.
double NormInf() const
Returns 0.0 because this class cannot compute Inf-norm.
bool HasNormInf() const
Returns false because this class cannot compute an Inf-norm.
int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y.
const char * Label() const
Returns a character string describing the operator.
int Solve(bool Trans, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Ifpack_CrsIlut forward/back solve on a Epetra_MultiVector X in Y (works for E...
void SetOverlapMode(Epetra_CombineMode OverlapMode)
Generate Ifpack_OverlapGraph object using current settings.
Epetra_CombineMode OverlapMode() const
Returns the overlap mode used to combine terms that are redundantly computed.
int NumMyNonzeros() const
Returns the number of nonzero entries in the local graph.
const Epetra_Comm & Comm() const
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y.
int NumGlobalNonzeros() const
Returns the number of nonzero entries in the global graph.
int SetUseTranspose(bool UseTranspose)
If set true, transpose of this operator will be applied.
bool UseTranspose() const
Returns the current UseTranspose setting.
int Condest(bool Trans, double &ConditionNumberEstimate) const
Returns the maximum over all the condition number estimate for each local ILU set of factors.
const Epetra_Vector & D() const
Returns the address of the D factor associated with this factored matrix.
virtual ~Ifpack_OverlapSolveObject()
Ifpack_OverlapSolveObject Destructor.
int SetLowerOperator(Epetra_CrsMatrix *L, bool UseLTrans)
Define the operator to be used for the lower triangle.
const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
const Epetra_CrsMatrix & U() const
Returns the address of the L factor associated with this factored matrix.
int Multiply(bool Trans, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of multiplying U, D and L in that order on an Epetra_MultiVector X in Y.
const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
int SetUpperOperator(Epetra_CrsMatrix *U, bool UseUTrans)
Define the operator to be used for the upper triangle.
int SetDiagonal(Epetra_Vector *D, bool UseDInv)
Define the vector to be used for the diagonal.