Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Ifpack_CrsRick.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_CRSRICK_H_
44 #define _IFPACK_CRSRICK_H_
45 
46 #include "Ifpack_ScalingType.h"
47 #include "Ifpack_IlukGraph.h"
48 #include "Epetra_CompObject.h"
49 #include "Epetra_Operator.h"
50 #include "Epetra_CrsMatrix.h"
51 #include "Epetra_Object.h"
52 class Epetra_Comm;
53 class Epetra_Map;
54 class Epetra_CrsMatrix;
55 class Epetra_Vector;
56 class Epetra_MultiVector;
57 
58 namespace Teuchos {
59  class ParameterList;
60 }
61 
63 
205 class Ifpack_CrsRick: public Epetra_Object, public Epetra_CompObject, public virtual Epetra_Operator {
206 
207  // Give ostream << function some access to private and protected data/functions.
208 
209  friend std::ostream& operator << (std::ostream& os, const Ifpack_CrsRick& A);
210 
211  public:
213 
221 
223  Ifpack_CrsRick(const Ifpack_CrsRick & Matrix);
224 
226  virtual ~Ifpack_CrsRick();
227 
229 
231  int InitValues();
232 
234  bool ValuesInitialized() const {return(ValuesInitialized_);};
235 
237  void SetRelaxValue( double RelaxValue) {RelaxValue_ = RelaxValue; return;}
238 
240  void SetAbsoluteThreshold( double Athresh) {Athresh_ = Athresh; return;}
241 
243  void SetRelativeThreshold( double Rthresh) {Rthresh_ = Rthresh; return;}
244 
246  void SetOverlapMode( Epetra_CombineMode OverlapMode) {OverlapMode_ = OverlapMode; return;}
247 
249  /* This method is only available if the Teuchos package is enabled.
250  This method recognizes four parameter names: relax_value,
251  absolute_threshold, relative_threshold and overlap_mode. These names
252  are case insensitive. The ParameterEntry must have type double for all
253  except overlap_mode, which must have type Epetra_CombineMode.
254  */
255  int SetParameters(const Teuchos::ParameterList& parameterlist,
256  bool cerr_warning_if_unused=false);
257 
259 
267  int Factor();
268 
270  bool Factored() const {return(Factored_);};
271 
272 
273  // Mathematical functions.
274 
275 
277 
287  int Solve(bool Trans, const Epetra_Vector& x, Epetra_Vector& y) const;
288 
290 
300  int Solve(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
301 
303 
313  int Multiply(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
314 
316 
324  int Condest(bool Trans, double & ConditionNumberEstimate) const;
325  // Attribute access functions
326 
327 
329  double GetRelaxValue() {return RelaxValue_;}
330 
332  double GetAbsoluteThreshold() {return Athresh_;}
333 
335  double GetRelativeThreshold() {return Rthresh_;}
336 
339 
341  int NumGlobalRows() const {return(Graph().NumGlobalRows());};
342 
344  int NumGlobalCols() const {return(Graph().NumGlobalCols());};
345 
347  int NumGlobalNonzeros() const {return(Graph().NumGlobalNonzeros());};
348 
350  virtual int NumGlobalDiagonals() const {return(Graph().NumGlobalDiagonals());};
351 
353  int NumMyRows() const {return(Graph().NumMyRows());};
354 
356  int NumMyCols() const {return(Graph().NumMyCols());};
357 
359  int NumMyNonzeros() const {return(Graph().NumMyNonzeros());};
360 
362  virtual int NumMyDiagonals() const {return(Graph().NumMyDiagonals());};
363 
365  int IndexBase() const {return(Graph().IndexBase());};
366 
368  const Ifpack_IlukGraph & Graph() const {return(Graph_);};
369 
371  const Epetra_Vector & D() const {return(*D_);};
372 
374  const Epetra_CrsMatrix & U() const {return(*U_);};
375 
377 
379  const char * Label() const {return(Epetra_Object::Label());};
380 
382 
392 
394 
405  int Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const {
406  return(Multiply(Ifpack_CrsRick::UseTranspose(), X, Y));};
407 
409 
423  return(Solve(Ifpack_CrsRick::UseTranspose(), X, Y));};
424 
426  double NormInf() const {return(0.0);};
427 
429  bool HasNormInf() const {return(false);};
430 
432  bool UseTranspose() const {return(UseTranspose_);};
433 
435  const Epetra_Map & OperatorDomainMap() const {return(A_.DomainMap());};
436 
438  const Epetra_Map & OperatorRangeMap() const{return(A_.RangeMap());};
440 
441  protected:
442  void SetFactored(bool Flag) {Factored_ = Flag;};
443  void SetValuesInitialized(bool Flag) {ValuesInitialized_ = Flag;};
444  bool Allocated() const {return(Allocated_);};
445  int SetAllocated(bool Flag) {Allocated_ = Flag; return(0);};
446 
447  private:
448 
449 
450  int Allocate();
451 
457 
458 
461  bool Factored_;
462  double RelaxValue_;
463  double Athresh_;
464  double Rthresh_;
465  mutable double Condest_;
466 
470 
471 
472 };
473 
475 std::ostream& operator << (std::ostream& os, const Ifpack_CrsRick& A);
476 
477 #endif /* _IFPACK_CRSRICK_H_ */
void SetRelaxValue(double RelaxValue)
Set RIC(k) relaxation parameter.
const Epetra_Map & RangeMap() const
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.
const Epetra_Vector & D() const
Returns the address of the D factor associated with this factored matrix.
Epetra_CrsMatrix * U_
int IndexBase() const
Returns the index base for row and column indices for this graph.
const Epetra_CrsMatrix & U() const
Returns the address of the U factor associated with this factored matrix.
bool ValuesInitialized() const
If values have been initialized, this query returns true, otherwise it returns false.
bool Allocated() const
Epetra_MultiVector * OverlapY_
int SetAllocated(bool Flag)
Epetra_MultiVector * OverlapX_
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 NumMyRows() const
Returns the number of local matrix rows.
int NumGlobalNonzeros() const
Returns the number of nonzero entries in the global graph.
const Ifpack_IlukGraph & Graph() const
Returns the address of the Ifpack_IlukGraph associated with this factored matrix. ...
int InitValues()
Initialize L and U with values from user matrix A.
Epetra_CombineMode GetOverlapMode()
Get overlap mode type.
const Epetra_CrsMatrix & A_
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 Ifpack_IlukGraph & Graph_
int SetParameters(const Teuchos::ParameterList &parameterlist, bool cerr_warning_if_unused=false)
Set parameters using a Teuchos::ParameterList object.
int NumGlobalRows() const
Returns the number of global matrix rows.
int SetUseTranspose(bool UseTranspose)
If set true, transpose of this operator will be applied.
Ifpack_ScalingType enumerable type.
double GetRelaxValue()
Get RIC(k) relaxation parameter.
bool HasNormInf() const
Returns false because this class cannot compute an Inf-norm.
virtual int NumGlobalDiagonals() const
Returns the number of diagonal entries found in the global input graph.
Epetra_Vector * D_
adjacency_list< vecS, vecS, undirectedS, no_property, property< edge_weight_t, double > > Graph
Ifpack_CrsRick: A class for constructing and using an incomplete Cholesky (IC) factorization of a giv...
Epetra_CombineMode OverlapMode_
std::ostream & operator<<(std::ostream &os, const Ifpack_Container &obj)
const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
int NumMyNonzeros() const
Returns the number of nonzero entries in the local graph.
Ifpack_CrsRick(const Epetra_CrsMatrix &A, const Ifpack_IlukGraph &Graph)
Ifpack_CrsRick constuctor with variable number of indices per row.
virtual int NumMyDiagonals() const
Returns the number of diagonal entries found in the local input graph.
bool UseTranspose() const
Returns the current UseTranspose setting.
int NumGlobalCols() const
Returns the number of global matrix columns.
virtual ~Ifpack_CrsRick()
Ifpack_CrsRick Destructor.
void SetFactored(bool Flag)
void SetAbsoluteThreshold(double Athresh)
Set absolute threshold value.
int NumMyCols() const
Returns the number of local matrix columns.
int Solve(bool Trans, const Epetra_Vector &x, Epetra_Vector &y) const
Returns the result of a Ifpack_CrsRick forward/back solve on a Epetra_Vector x in y...
const Epetra_Map & DomainMap() const
void SetValuesInitialized(bool Flag)
Epetra_CombineMode
const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
int Factor()
Compute IC factor L using the specified graph, diagonal perturbation thresholds and relaxation parame...
friend std::ostream & operator<<(std::ostream &os, const Ifpack_CrsRick &A)
&lt;&lt; operator will work for Ifpack_CrsRick.
const char * Label() const
Returns a character string describing the operator.
int Multiply(bool Trans, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of multiplying U, D and U^T in that order on an Epetra_MultiVector X in Y...
int Condest(bool Trans, double &ConditionNumberEstimate) const
Returns the maximum over all the condition number estimate for each local IC set of factors...
bool Factored() const
If factor is completed, this query returns true, otherwise it returns false.
double GetRelativeThreshold()
Get relative threshold value.
double NormInf() const
Returns 0.0 because this class cannot compute Inf-norm.
void SetRelativeThreshold(double Rthresh)
Set relative threshold value.
double GetAbsoluteThreshold()
Get absolute threshold value.