00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030 #ifndef _IFPACK_CRSRICK_H_
00031 #define _IFPACK_CRSRICK_H_
00032
00033 #include "Ifpack_ScalingType.h"
00034 #include "Ifpack_IlukGraph.h"
00035 #include "Epetra_CompObject.h"
00036 #include "Epetra_Operator.h"
00037 #include "Epetra_CrsMatrix.h"
00038 #include "Epetra_Object.h"
00039 class Epetra_Comm;
00040 class Epetra_Map;
00041 class Epetra_CrsMatrix;
00042 class Epetra_Vector;
00043 class Epetra_MultiVector;
00044
00045 namespace Teuchos {
00046 class ParameterList;
00047 }
00048
00050
00192 class Ifpack_CrsRick: public Epetra_Object, public Epetra_CompObject, public virtual Epetra_Operator {
00193
00194
00195
00196 friend ostream& operator << (ostream& os, const Ifpack_CrsRick& A);
00197
00198 public:
00200
00207 Ifpack_CrsRick(const Epetra_CrsMatrix &A, const Ifpack_IlukGraph & Graph);
00208
00210 Ifpack_CrsRick(const Ifpack_CrsRick & Matrix);
00211
00213 virtual ~Ifpack_CrsRick();
00214
00216
00218 int InitValues();
00219
00221 bool ValuesInitialized() const {return(ValuesInitialized_);};
00222
00224 void SetRelaxValue( double RelaxValue) {RelaxValue_ = RelaxValue; return;}
00225
00227 void SetAbsoluteThreshold( double Athresh) {Athresh_ = Athresh; return;}
00228
00230 void SetRelativeThreshold( double Rthresh) {Rthresh_ = Rthresh; return;}
00231
00233 void SetOverlapMode( Epetra_CombineMode OverlapMode) {OverlapMode_ = OverlapMode; return;}
00234
00236
00237
00238
00239
00240
00241
00242 int SetParameters(const Teuchos::ParameterList& parameterlist,
00243 bool cerr_warning_if_unused=false);
00244
00246
00254 int Factor();
00255
00257 bool Factored() const {return(Factored_);};
00258
00259
00260
00261
00262
00264
00274 int Solve(bool Trans, const Epetra_Vector& x, Epetra_Vector& y) const;
00275
00277
00287 int Solve(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
00288
00290
00300 int Multiply(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
00301
00303
00311 int Condest(bool Trans, double & ConditionNumberEstimate) const;
00312
00313
00314
00316 double GetRelaxValue() {return RelaxValue_;}
00317
00319 double GetAbsoluteThreshold() {return Athresh_;}
00320
00322 double GetRelativeThreshold() {return Rthresh_;}
00323
00325 Epetra_CombineMode GetOverlapMode() {return OverlapMode_;}
00326
00328 int NumGlobalRows() const {return(Graph().NumGlobalRows());};
00329
00331 int NumGlobalCols() const {return(Graph().NumGlobalCols());};
00332
00334 int NumGlobalNonzeros() const {return(Graph().NumGlobalNonzeros());};
00335
00337 virtual int NumGlobalDiagonals() const {return(Graph().NumGlobalDiagonals());};
00338
00340 int NumMyRows() const {return(Graph().NumMyRows());};
00341
00343 int NumMyCols() const {return(Graph().NumMyCols());};
00344
00346 int NumMyNonzeros() const {return(Graph().NumMyNonzeros());};
00347
00349 virtual int NumMyDiagonals() const {return(Graph().NumMyDiagonals());};
00350
00352 int IndexBase() const {return(Graph().IndexBase());};
00353
00355 const Ifpack_IlukGraph & Graph() const {return(Graph_);};
00356
00358 const Epetra_Vector & D() const {return(*D_);};
00359
00361 const Epetra_CrsMatrix & U() const {return(*U_);};
00362
00364
00366 char * Label() const {return(Epetra_Object::Label());};
00367
00369
00378 int SetUseTranspose(bool UseTranspose) {UseTranspose_ = UseTranspose; return(0);};
00379
00381
00392 int Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const {
00393 return(Multiply(Ifpack_CrsRick::UseTranspose(), X, Y));};
00394
00396
00409 int ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const {
00410 return(Solve(Ifpack_CrsRick::UseTranspose(), X, Y));};
00411
00413 double NormInf() const {return(0.0);};
00414
00416 bool HasNormInf() const {return(false);};
00417
00419 bool UseTranspose() const {return(UseTranspose_);};
00420
00422 const Epetra_Map & OperatorDomainMap() const {return(A_.DomainMap());};
00423
00425 const Epetra_Map & OperatorRangeMap() const{return(A_.RangeMap());};
00427
00428 protected:
00429 void SetFactored(bool Flag) {Factored_ = Flag;};
00430 void SetValuesInitialized(bool Flag) {ValuesInitialized_ = Flag;};
00431 bool Allocated() const {return(Allocated_);};
00432 int SetAllocated(bool Flag) {Allocated_ = Flag; return(0);};
00433
00434 private:
00435
00436
00437 int Allocate();
00438
00439 const Epetra_CrsMatrix &A_;
00440 const Ifpack_IlukGraph & Graph_;
00441 Epetra_CrsMatrix * U_;
00442 Epetra_Vector * D_;
00443 bool UseTranspose_;
00444
00445
00446 bool Allocated_;
00447 bool ValuesInitialized_;
00448 bool Factored_;
00449 double RelaxValue_;
00450 double Athresh_;
00451 double Rthresh_;
00452 mutable double Condest_;
00453
00454 mutable Epetra_MultiVector * OverlapX_;
00455 mutable Epetra_MultiVector * OverlapY_;
00456 Epetra_CombineMode OverlapMode_;
00457
00458
00459 };
00460
00462 ostream& operator << (ostream& os, const Ifpack_CrsRick& A);
00463
00464 #endif