IFPACK  Development
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 #include "Epetra_Operator.h"
47 #include "Epetra_CrsMatrix.h"
48 #include "Epetra_Vector.h"
49 #include "Epetra_CombineMode.h"
50 class Epetra_Flops;
51 
53 
55 
56  public:
60 
63 
67 
69 
71  void SetOverlapMode( Epetra_CombineMode OverlapMode) {OverlapMode_ = OverlapMode; return;}
72 
74  int SetLowerOperator (Epetra_CrsMatrix * L, bool UseLTrans) {L_ = L; UseLTrans_ = UseLTrans; return(0);};
75 
77  int SetDiagonal (Epetra_Vector * D, bool UseDInv) {D_ = D; UseDInv_ = UseDInv; return(0);};
78 
80  int SetUpperOperator (Epetra_CrsMatrix * U, bool UseUTrans) {U_ = U; UseUTrans_ = UseUTrans; return(0);};
82 
84 
85 
87 
97  int Solve(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
98 
100 
110  int Multiply(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
111 
113 
121  int Condest(bool Trans, double & ConditionNumberEstimate) const;
123 
125 
127 
132  Epetra_CombineMode OverlapMode() const {return(OverlapMode_);}
133 
135  int NumGlobalNonzeros() const {return(L().NumGlobalNonzeros()+U().NumGlobalNonzeros());};
136 
138  int NumMyNonzeros() const {return(L().NumMyNonzeros()+U().NumMyNonzeros());};
139 
141  const Epetra_CrsMatrix & L() const {return(*L_);};
142 
144  const Epetra_Vector & D() const {return(*D_);};
145 
147  const Epetra_CrsMatrix & U() const {return(*U_);};
149 
151 
153  const char * Label() const {return(Label_);};
154 
156 
165  int SetUseTranspose(bool UseTranspose) {UseTranspose_ = UseTranspose; return(0);};
166 
168 
179  int Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const {
181 
183 
198 
200  double NormInf() const {return(0.0);};
201 
203  bool HasNormInf() const {return(false);};
204 
206  bool UseTranspose() const {return(UseTranspose_);};
207 
209  const Epetra_Map & OperatorDomainMap() const {return(U_->OperatorDomainMap());};
210 
212  const Epetra_Map & OperatorRangeMap() const{return(L_->OperatorRangeMap());};
213 
215  const Epetra_Comm & Comm() const{return(Comm_);};
217 
218  protected:
219  virtual int SetupXY(bool Trans,
220  const Epetra_MultiVector& Xin, const Epetra_MultiVector& Yin,
221  Epetra_MultiVector * & Xout, Epetra_MultiVector * & Yout) const=0;
222  private:
223  char * Label_;
224  Epetra_CrsMatrix * L_;
225  bool UseLTrans_;
226  Epetra_Vector * D_;
227  bool UseDInv_;
228  Epetra_CrsMatrix * U_;
229  bool UseUTrans_;
230  bool UseTranspose_;
231  const Epetra_Comm & Comm_;
232  mutable double Condest_;
233  Epetra_Flops * Counter_;
234  Epetra_CombineMode OverlapMode_;
235 
236 };
237 #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.
Ifpack_OverlapSolveObject(char *Label, const Epetra_Comm &Comm)
Constructor.
const Epetra_Comm & Comm() const
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
const Epetra_CrsMatrix & L() const
Returns the address of the L factor associated with this factored matrix.
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.
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.
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.
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 char * Label() const
Returns a character string describing the operator.
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 & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
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.
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.
const Epetra_Vector & D() const
Returns the address of the D factor associated with this factored matrix.