Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Ifpack_CrsRiluk.h
Go to the documentation of this file.
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_CRSRILUK_H_
44 #define _IFPACK_CRSRILUK_H_
45 
46 #include "Ifpack_ConfigDefs.h"
47 #include "Ifpack_ScalingType.h"
48 #include "Ifpack_IlukGraph.h"
49 #include "Epetra_ConfigDefs.h"
50 #include "Epetra_CompObject.h"
51 #include "Epetra_Operator.h"
52 #include "Epetra_CrsMatrix.h"
53 #include "Epetra_Object.h"
54 #include "Epetra_MultiVector.h"
55 #include "Epetra_Vector.h"
56 #include "Epetra_Map.h"
57 class Epetra_Comm;
58 class Epetra_VbrMatrix;
59 class Epetra_RowMatrix;
60 
61 #include "Teuchos_RefCountPtr.hpp"
62 
63 namespace Teuchos {
64  class ParameterList;
65 }
66 
68 
210 class Ifpack_CrsRiluk: public Epetra_Object, public Epetra_CompObject, public virtual Epetra_Operator {
211 
212  // Give ostream << function some access to private and protected data/functions.
213 
214  friend std::ostream& operator << (std::ostream& os, const Ifpack_CrsRiluk& A);
215 
216  public:
218 
223  Ifpack_CrsRiluk(const Ifpack_IlukGraph & Graph_in);
224 
226  Ifpack_CrsRiluk(const Ifpack_CrsRiluk & Matrix);
227 
229  virtual ~Ifpack_CrsRiluk();
230 
232 
238  int InitValues(const Epetra_CrsMatrix &A);
239 
241 
247 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES // FIXME LONG LONG
248  int InitValues(const Epetra_VbrMatrix &A);
249 #endif
250 
252  bool ValuesInitialized() const {return(ValuesInitialized_);};
253 
255  void SetRelaxValue( double RelaxValue) {RelaxValue_ = RelaxValue; return;}
256 
258  void SetAbsoluteThreshold( double Athresh) {Athresh_ = Athresh; return;}
259 
261  void SetRelativeThreshold( double Rthresh) {Rthresh_ = Rthresh; return;}
262 
264  void SetOverlapMode( Epetra_CombineMode OverlapMode) {OverlapMode_ = OverlapMode; return;}
265 
267  /* This method is only available if the Teuchos package is enabled.
268  This method recognizes four parameter names: relax_value,
269  absolute_threshold, relative_threshold and overlap_mode. These names are
270  case insensitive, and in each case except overlap_mode, the ParameterEntry
271  must have type double. For overlap_mode, the ParameterEntry must have
272  type Epetra_CombineMode.
273  */
274  int SetParameters(const Teuchos::ParameterList& parameterlist,
275  bool cerr_warning_if_unused=false);
276 
278 
286  int Factor();
287 
289  bool Factored() const {return(Factored_);};
290 
291 
292  // Mathematical functions.
293 
294 
296 
306  int Solve(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
307 
309 
319  int Multiply(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
320 
322 
330  int Condest(bool Trans, double & ConditionNumberEstimate) const;
331  // Attribute access functions
332 
334  double GetRelaxValue() {return RelaxValue_;}
335 
337  double GetAbsoluteThreshold() {return Athresh_;}
338 
340  double GetRelativeThreshold() {return Rthresh_;}
341 
344 
345 
346 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
347  int NumGlobalRows() const {return(Graph().NumGlobalRows());};
349 
351  int NumGlobalCols() const {return(Graph().NumGlobalCols());};
352 
354  int NumGlobalNonzeros() const {return(L().NumGlobalNonzeros()+U().NumGlobalNonzeros());};
355 
357  virtual int NumGlobalBlockDiagonals() const {return(Graph().NumGlobalBlockDiagonals());};
358 #endif
359 
360  long long NumGlobalRows64() const {return(Graph().NumGlobalRows64());};
361  long long NumGlobalCols64() const {return(Graph().NumGlobalCols64());};
362  long long NumGlobalNonzeros64() const {return(L().NumGlobalNonzeros64()+U().NumGlobalNonzeros64());};
363  virtual long long NumGlobalBlockDiagonals64() const {return(Graph().NumGlobalBlockDiagonals64());};
364 
366  int NumMyRows() const {return(Graph().NumMyRows());};
367 
369  int NumMyCols() const {return(Graph().NumMyCols());};
370 
372  int NumMyNonzeros() const {return(L().NumMyNonzeros()+U().NumMyNonzeros());};
373 
375  virtual int NumMyBlockDiagonals() const {return(Graph().NumMyBlockDiagonals());};
376 
378  virtual int NumMyDiagonals() const {return(NumMyDiagonals_);};
379 
381 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
382  int IndexBase() const {return(Graph().IndexBase());};
383 #endif
384  long long IndexBase64() const {return(Graph().IndexBase64());};
385 
387  const Ifpack_IlukGraph & Graph() const {return(Graph_);};
388 
390  const Epetra_CrsMatrix & L() const {return(*L_);};
391 
393  const Epetra_Vector & D() const {return(*D_);};
394 
396  const Epetra_CrsMatrix & U() const {return(*U_);};
397 
399 
401  const char * Label() const {return(Epetra_Object::Label());};
402 
404 
413  int SetUseTranspose(bool UseTranspose_in) {UseTranspose_ = UseTranspose_in; return(0);};
414 
416 
427  int Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const {
428  return(Multiply(Ifpack_CrsRiluk::UseTranspose(), X, Y));};
429 
431 
445  return(Solve(Ifpack_CrsRiluk::UseTranspose(), X, Y));};
446 
448  double NormInf() const {return(0.0);};
449 
451  bool HasNormInf() const {return(false);};
452 
454  bool UseTranspose() const {return(UseTranspose_);};
455 
457  const Epetra_Map & OperatorDomainMap() const {return(U_->OperatorDomainMap());};
458 
460  const Epetra_Map & OperatorRangeMap() const{return(L_->OperatorRangeMap());};
461 
463  const Epetra_Comm & Comm() const{return(Comm_);};
465 
466  protected:
467  void SetFactored(bool Flag) {Factored_ = Flag;};
468  void SetValuesInitialized(bool Flag) {ValuesInitialized_ = Flag;};
469  bool Allocated() const {return(Allocated_);};
470  int SetAllocated(bool Flag) {Allocated_ = Flag; return(0);};
471  int BlockGraph2PointGraph(const Epetra_CrsGraph & BG, Epetra_CrsGraph & PG, bool Upper);
472 
473  private:
474 
475 
476  int AllocateCrs();
477  int AllocateVbr();
478  int InitAllValues(const Epetra_RowMatrix & A, int MaxNumEntries);
479  int BlockMap2PointMap(const Epetra_BlockMap & BlockMap, Teuchos::RefCountPtr<Epetra_Map>* PointMap);
480  int GenerateXY(bool Trans,
481  const Epetra_MultiVector& Xin, const Epetra_MultiVector& Yin,
482  Teuchos::RefCountPtr<Epetra_MultiVector>* Xout,
483  Teuchos::RefCountPtr<Epetra_MultiVector>* Yout) const;
488  Teuchos::RefCountPtr<Epetra_Map> IlukRowMap_;
489  Teuchos::RefCountPtr<Epetra_Map> IlukDomainMap_;
490  Teuchos::RefCountPtr<Epetra_Map> IlukRangeMap_;
491  Teuchos::RefCountPtr<const Epetra_Map> U_DomainMap_;
492  Teuchos::RefCountPtr<const Epetra_Map> L_RangeMap_;
494  Teuchos::RefCountPtr<Epetra_CrsMatrix> L_;
495  Teuchos::RefCountPtr<Epetra_CrsMatrix> U_;
496  Teuchos::RefCountPtr<Epetra_CrsGraph> L_Graph_;
497  Teuchos::RefCountPtr<Epetra_CrsGraph> U_Graph_;
498  Teuchos::RefCountPtr<Epetra_Vector> D_;
500 
504  bool Factored_;
505  double RelaxValue_;
506  double Athresh_;
507  double Rthresh_;
508  mutable double Condest_;
509 
510  mutable Teuchos::RefCountPtr<Epetra_MultiVector> OverlapX_;
511  mutable Teuchos::RefCountPtr<Epetra_MultiVector> OverlapY_;
512  mutable Teuchos::RefCountPtr<Epetra_MultiVector> VbrX_;
513  mutable Teuchos::RefCountPtr<Epetra_MultiVector> VbrY_;
515 
516 
517 };
518 
520 std::ostream& operator << (std::ostream& os, const Ifpack_CrsRiluk& A);
521 
522 #endif /* _IFPACK_CRSRILUK_H_ */
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 GenerateXY(bool Trans, const Epetra_MultiVector &Xin, const Epetra_MultiVector &Yin, Teuchos::RefCountPtr< Epetra_MultiVector > *Xout, Teuchos::RefCountPtr< Epetra_MultiVector > *Yout) const
Teuchos::RefCountPtr< const Epetra_Map > L_RangeMap_
double GetRelaxValue()
Get RILU(k) relaxation parameter.
int NumMyNonzeros() const
Returns the number of nonzero entries in the local graph.
int BlockMap2PointMap(const Epetra_BlockMap &BlockMap, Teuchos::RefCountPtr< Epetra_Map > *PointMap)
Ifpack_CrsRiluk: A class for constructing and using an incomplete lower/upper (ILU) factorization of ...
int Condest(bool Trans, double &ConditionNumberEstimate) const
Returns the maximum over all the condition number estimate for each local ILU set of factors...
void SetFactored(bool Flag)
long long NumGlobalRows64() const
Epetra_CombineMode GetOverlapMode()
Get overlap mode type.
double NormInf() const
Returns 0.0 because this class cannot compute Inf-norm.
Teuchos::RefCountPtr< Epetra_CrsGraph > U_Graph_
Teuchos::RefCountPtr< const Epetra_Map > U_DomainMap_
const Epetra_CrsMatrix & L() const
Returns the address of the L factor associated with this factored matrix.
Teuchos::RefCountPtr< Epetra_MultiVector > OverlapY_
const char * Label() const
Returns a character string describing the operator.
const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
void SetOverlapMode(Epetra_CombineMode OverlapMode)
Set overlap mode type.
Ifpack_IlukGraph: A class for constructing level filled graphs for use with ILU(k) class precondition...
int NumGlobalCols() const
Returns the number of global matrix columns.
int BlockGraph2PointGraph(const Epetra_CrsGraph &BG, Epetra_CrsGraph &PG, bool Upper)
void SetAbsoluteThreshold(double Athresh)
Set absolute threshold value.
const Epetra_Comm & Comm() const
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
int SetAllocated(bool Flag)
virtual long long NumGlobalBlockDiagonals64() const
void SetValuesInitialized(bool Flag)
void SetRelativeThreshold(double Rthresh)
Set relative threshold value.
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...
Teuchos::RefCountPtr< Epetra_CrsGraph > L_Graph_
Teuchos::RefCountPtr< Epetra_Map > IlukDomainMap_
long long NumGlobalCols64() const
Ifpack_ScalingType enumerable type.
Teuchos::RefCountPtr< Epetra_CrsMatrix > U_
int InitAllValues(const Epetra_RowMatrix &A, int MaxNumEntries)
const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
int SetParameters(const Teuchos::ParameterList &parameterlist, bool cerr_warning_if_unused=false)
Set parameters using a Teuchos::ParameterList object.
const Ifpack_IlukGraph & Graph_
bool Allocated() const
Teuchos::RefCountPtr< Epetra_Map > IlukRowMap_
int NumMyCols() const
Returns the number of local matrix columns.
Teuchos::RefCountPtr< Epetra_MultiVector > OverlapX_
long long IndexBase64() const
int NumGlobalRows() const
Returns the number of global matrix rows.
int InitValues(const Epetra_CrsMatrix &A)
Initialize L and U with values from user matrix A.
const Epetra_CrsMatrix & U() const
Returns the address of the L factor associated with this factored matrix.
double GetAbsoluteThreshold()
Get absolute threshold value.
std::ostream & operator<<(std::ostream &os, const Ifpack_Container &obj)
virtual int NumGlobalBlockDiagonals() const
Returns the number of diagonal entries found in the global input 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.
Teuchos::RefCountPtr< Epetra_CrsMatrix > L_
int Solve(bool Trans, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Ifpack_CrsRiluk forward/back solve on a Epetra_MultiVector X in Y (works for ...
int IndexBase() const
Returns the index base for row and column indices for this graph.
friend std::ostream & operator<<(std::ostream &os, const Ifpack_CrsRiluk &A)
&lt;&lt; operator will work for Ifpack_CrsRiluk.
int NumGlobalNonzeros() const
Returns the number of nonzero entries in the global graph.
bool HasNormInf() const
Returns false because this class cannot compute an Inf-norm.
const Epetra_Comm & Comm_
bool Factored() const
If factor is completed, this query returns true, otherwise it returns false.
int NumMyRows() const
Returns the number of local matrix rows.
bool ValuesInitialized() const
If values have been initialized, this query returns true, otherwise it returns false.
void SetRelaxValue(double RelaxValue)
Set RILU(k) relaxation parameter.
const Epetra_Vector & D() const
Returns the address of the D factor associated with this factored matrix.
bool UseTranspose() const
Returns the current UseTranspose setting.
Teuchos::RefCountPtr< Epetra_MultiVector > VbrY_
Epetra_CombineMode
long long NumGlobalNonzeros64() const
double GetRelativeThreshold()
Get relative threshold value.
int SetUseTranspose(bool UseTranspose_in)
If set true, transpose of this operator will be applied.
Teuchos::RefCountPtr< Epetra_Vector > D_
Teuchos::RefCountPtr< Epetra_MultiVector > VbrX_
virtual int NumMyBlockDiagonals() const
Returns the number of diagonal entries found in the local input graph.
virtual ~Ifpack_CrsRiluk()
Ifpack_CrsRiluk Destructor.
virtual int NumMyDiagonals() const
Returns the number of nonzero diagonal values found in matrix.
Epetra_CombineMode OverlapMode_
Teuchos::RefCountPtr< Epetra_Map > IlukRangeMap_
int Factor()
Compute ILU factors L and U using the specified graph, diagonal perturbation thresholds and relaxatio...
const Ifpack_IlukGraph & Graph() const
returns the address of the Ifpack_IlukGraph associated with this factored matrix. ...
Ifpack_CrsRiluk(const Ifpack_IlukGraph &Graph_in)
Ifpack_CrsRiluk constuctor with variable number of indices per row.