IFPACK Development
Loading...
Searching...
No Matches
Ifpack_CrsRick.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_CRSRICK_H_
44#define _IFPACK_CRSRICK_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 "Ifpack_ScalingType.h"
53#include "Ifpack_IlukGraph.h"
54#include "Epetra_CompObject.h"
55#include "Epetra_Operator.h"
56#include "Epetra_CrsMatrix.h"
57#include "Epetra_Object.h"
58class Epetra_Comm;
59class Epetra_Map;
61class Epetra_Vector;
63
64namespace Teuchos {
65 class ParameterList;
66}
67
69
211class Ifpack_CrsRick: public Epetra_Object, public Epetra_CompObject, public virtual Epetra_Operator {
212
213 // Give ostream << function some access to private and protected data/functions.
214
215 friend std::ostream& operator << (std::ostream& os, const Ifpack_CrsRick& A);
216
217 public:
219
226 Ifpack_CrsRick(const Epetra_CrsMatrix &A, const Ifpack_IlukGraph & Graph);
227
229 Ifpack_CrsRick(const Ifpack_CrsRick & Matrix);
230
232 virtual ~Ifpack_CrsRick();
233
235
237 int InitValues();
238
240 bool ValuesInitialized() const {return(ValuesInitialized_);};
241
243 void SetRelaxValue( double RelaxValue) {RelaxValue_ = RelaxValue; return;}
244
246 void SetAbsoluteThreshold( double Athresh) {Athresh_ = Athresh; return;}
247
249 void SetRelativeThreshold( double Rthresh) {Rthresh_ = Rthresh; return;}
250
252 void SetOverlapMode( Epetra_CombineMode OverlapMode) {OverlapMode_ = OverlapMode; return;}
253
255 /* This method is only available if the Teuchos package is enabled.
256 This method recognizes four parameter names: relax_value,
257 absolute_threshold, relative_threshold and overlap_mode. These names
258 are case insensitive. The ParameterEntry must have type double for all
259 except overlap_mode, which must have type Epetra_CombineMode.
260 */
261 int SetParameters(const Teuchos::ParameterList& parameterlist,
262 bool cerr_warning_if_unused=false);
263
265
273 int Factor();
274
276 bool Factored() const {return(Factored_);};
277
278
279 // Mathematical functions.
280
281
283
293 int Solve(bool Trans, const Epetra_Vector& x, Epetra_Vector& y) const;
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
333
335 double GetRelaxValue() {return RelaxValue_;}
336
338 double GetAbsoluteThreshold() {return Athresh_;}
339
341 double GetRelativeThreshold() {return Rthresh_;}
342
344 Epetra_CombineMode GetOverlapMode() {return OverlapMode_;}
345
347 int NumGlobalRows() const {return(Graph().NumGlobalRows());};
348
350 int NumGlobalCols() const {return(Graph().NumGlobalCols());};
351
353 int NumGlobalNonzeros() const {return(Graph().NumGlobalNonzeros());};
354
356 virtual int NumGlobalDiagonals() const {return(Graph().NumGlobalDiagonals());};
357
359 int NumMyRows() const {return(Graph().NumMyRows());};
360
362 int NumMyCols() const {return(Graph().NumMyCols());};
363
365 int NumMyNonzeros() const {return(Graph().NumMyNonzeros());};
366
368 virtual int NumMyDiagonals() const {return(Graph().NumMyDiagonals());};
369
371 int IndexBase() const {return(Graph().IndexBase());};
372
374 const Ifpack_IlukGraph & Graph() const {return(Graph_);};
375
377 const Epetra_Vector & D() const {return(*D_);};
378
380 const Epetra_CrsMatrix & U() const {return(*U_);};
381
383
385 const char * Label() const {return(Epetra_Object::Label());};
386
388
397 int SetUseTranspose(bool UseTranspose) {UseTranspose_ = UseTranspose; return(0);};
398
400
412 return(Multiply(Ifpack_CrsRick::UseTranspose(), X, Y));};
413
415
430
432 double NormInf() const {return(0.0);};
433
435 bool HasNormInf() const {return(false);};
436
438 bool UseTranspose() const {return(UseTranspose_);};
439
441 const Epetra_Map & OperatorDomainMap() const {return(A_.DomainMap());};
442
444 const Epetra_Map & OperatorRangeMap() const{return(A_.RangeMap());};
446
447 protected:
448 void SetFactored(bool Flag) {Factored_ = Flag;};
449 void SetValuesInitialized(bool Flag) {ValuesInitialized_ = Flag;};
450 bool Allocated() const {return(Allocated_);};
451 int SetAllocated(bool Flag) {Allocated_ = Flag; return(0);};
452
453 private:
454
455
456 int Allocate();
457
458 const Epetra_CrsMatrix &A_;
459 const Ifpack_IlukGraph & Graph_;
460 Epetra_CrsMatrix * U_;
461 Epetra_Vector * D_;
462 bool UseTranspose_;
463
464
465 bool Allocated_;
466 bool ValuesInitialized_;
467 bool Factored_;
468 double RelaxValue_;
469 double Athresh_;
470 double Rthresh_;
471 mutable double Condest_;
472
473 mutable Epetra_MultiVector * OverlapX_;
474 mutable Epetra_MultiVector * OverlapY_;
475 Epetra_CombineMode OverlapMode_;
476
477
478};
479
481std::ostream& operator << (std::ostream& os, const Ifpack_CrsRick& A);
482
483#endif /* _IFPACK_CRSRICK_H_ */
Epetra_CombineMode
Ifpack_ScalingType enumerable type.
const Epetra_Map & RangeMap() const
const Epetra_Map & DomainMap() const
virtual const char * Label() const
Ifpack_CrsRick: A class for constructing and using an incomplete Cholesky (IC) factorization of a giv...
Epetra_CombineMode GetOverlapMode()
Get overlap mode type.
bool UseTranspose() const
Returns the current UseTranspose setting.
int NumMyCols() const
Returns the number of local matrix columns.
virtual int NumGlobalDiagonals() const
Returns the number of diagonal entries found in the global input graph.
void SetAbsoluteThreshold(double Athresh)
Set absolute threshold value.
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_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
virtual int NumMyDiagonals() const
Returns the number of diagonal entries found in the local input graph.
const char * Label() const
Returns a character string describing the operator.
int SetParameters(const Teuchos::ParameterList &parameterlist, bool cerr_warning_if_unused=false)
Set parameters using a Teuchos::ParameterList object.
double GetAbsoluteThreshold()
Get absolute threshold value.
int NumGlobalCols() const
Returns the number of global matrix columns.
int Factor()
Compute IC factor L using the specified graph, diagonal perturbation thresholds and relaxation parame...
const Epetra_Vector & D() const
Returns the address of the D factor associated with this factored matrix.
double GetRelaxValue()
Get RIC(k) relaxation parameter.
double NormInf() const
Returns 0.0 because this class cannot compute Inf-norm.
int NumMyRows() const
Returns the number of local matrix rows.
int InitValues()
Initialize L and U with values from user matrix A.
int NumGlobalRows() const
Returns the number of global matrix rows.
const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
const Ifpack_IlukGraph & Graph() const
Returns the address of the Ifpack_IlukGraph associated with this factored matrix.
virtual ~Ifpack_CrsRick()
Ifpack_CrsRick Destructor.
bool Factored() const
If factor is completed, this query returns true, otherwise it returns false.
int SetUseTranspose(bool UseTranspose)
If set true, transpose of this operator will be applied.
const Epetra_CrsMatrix & U() const
Returns the address of the U factor associated with this factored matrix.
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.
double GetRelativeThreshold()
Get relative threshold value.
void SetRelaxValue(double RelaxValue)
Set RIC(k) relaxation parameter.
void SetRelativeThreshold(double Rthresh)
Set relative threshold value.
int NumGlobalNonzeros() const
Returns the number of nonzero entries in the global graph.
int Condest(bool Trans, double &ConditionNumberEstimate) const
Returns the maximum over all the condition number estimate for each local IC set of factors.
void SetOverlapMode(Epetra_CombineMode OverlapMode)
Set overlap mode type.
int NumMyNonzeros() const
Returns the number of nonzero entries in the local graph.
int IndexBase() const
Returns the index base for row and column indices for this graph.
friend std::ostream & operator<<(std::ostream &os, const Ifpack_CrsRick &A)
<< operator will work for Ifpack_CrsRick.
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 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.
bool HasNormInf() const
Returns false because this class cannot compute an Inf-norm.
bool ValuesInitialized() const
If values have been initialized, this query returns true, otherwise it returns false.
Ifpack_IlukGraph: A class for constructing level filled graphs for use with ILU(k) class precondition...