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 #ifdef HAVE_IFPACK_TEUCHOS
00046 namespace Teuchos {
00047 class ParameterList;
00048 }
00049 #endif
00050
00052
00194 class Ifpack_CrsRick: public Epetra_Object, public Epetra_CompObject, public virtual Epetra_Operator {
00195
00196
00197
00198 friend ostream& operator << (ostream& os, const Ifpack_CrsRick& A);
00199
00200 public:
00202
00209 Ifpack_CrsRick(const Epetra_CrsMatrix &A, const Ifpack_IlukGraph & Graph);
00210
00212 Ifpack_CrsRick(const Ifpack_CrsRick & Matrix);
00213
00215 virtual ~Ifpack_CrsRick();
00216
00218
00220 int InitValues();
00221
00223 bool ValuesInitialized() const {return(ValuesInitialized_);};
00224
00226 void SetRelaxValue( double RelaxValue) {RelaxValue_ = RelaxValue; return;}
00227
00229 void SetAbsoluteThreshold( double Athresh) {Athresh_ = Athresh; return;}
00230
00232 void SetRelativeThreshold( double Rthresh) {Rthresh_ = Rthresh; return;}
00233
00235 void SetOverlapMode( Epetra_CombineMode OverlapMode) {OverlapMode_ = OverlapMode; return;}
00236
00237 #ifdef HAVE_IFPACK_TEUCHOS
00238
00239
00240
00241
00242
00243
00244
00245 int SetParameters(const Teuchos::ParameterList& parameterlist,
00246 bool cerr_warning_if_unused=false);
00247 #endif
00248
00250
00258 int Factor();
00259
00261 bool Factored() const {return(Factored_);};
00262
00263
00264
00265
00266
00268
00278 int Solve(bool Trans, const Epetra_Vector& x, Epetra_Vector& y) const;
00279
00281
00291 int Solve(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
00292
00294
00304 int Multiply(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
00305
00307
00315 int Condest(bool Trans, double & ConditionNumberEstimate) const;
00316
00317
00318
00320 double GetRelaxValue() {return RelaxValue_;}
00321
00323 double GetAbsoluteThreshold() {return Athresh_;}
00324
00326 double GetRelativeThreshold() {return Rthresh_;}
00327
00329 Epetra_CombineMode GetOverlapMode() {return OverlapMode_;}
00330
00332 int NumGlobalRows() const {return(Graph().NumGlobalRows());};
00333
00335 int NumGlobalCols() const {return(Graph().NumGlobalCols());};
00336
00338 int NumGlobalNonzeros() const {return(Graph().NumGlobalNonzeros());};
00339
00341 virtual int NumGlobalDiagonals() const {return(Graph().NumGlobalDiagonals());};
00342
00344 int NumMyRows() const {return(Graph().NumMyRows());};
00345
00347 int NumMyCols() const {return(Graph().NumMyCols());};
00348
00350 int NumMyNonzeros() const {return(Graph().NumMyNonzeros());};
00351
00353 virtual int NumMyDiagonals() const {return(Graph().NumMyDiagonals());};
00354
00356 int IndexBase() const {return(Graph().IndexBase());};
00357
00359 const Ifpack_IlukGraph & Graph() const {return(Graph_);};
00360
00362 const Epetra_Vector & D() const {return(*D_);};
00363
00365 const Epetra_CrsMatrix & U() const {return(*U_);};
00366
00368
00370 char * Label() const {return(Epetra_Object::Label());};
00371
00373
00382 int SetUseTranspose(bool UseTranspose) {UseTranspose_ = UseTranspose; return(0);};
00383
00385
00396 int Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const {
00397 return(Multiply(Ifpack_CrsRick::UseTranspose(), X, Y));};
00398
00400
00413 int ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const {
00414 return(Solve(Ifpack_CrsRick::UseTranspose(), X, Y));};
00415
00417 double NormInf() const {return(0.0);};
00418
00420 bool HasNormInf() const {return(false);};
00421
00423 bool UseTranspose() const {return(UseTranspose_);};
00424
00426 const Epetra_Map & OperatorDomainMap() const {return(A_.DomainMap());};
00427
00429 const Epetra_Map & OperatorRangeMap() const{return(A_.RangeMap());};
00431
00432 protected:
00433 void SetFactored(bool Flag) {Factored_ = Flag;};
00434 void SetValuesInitialized(bool Flag) {ValuesInitialized_ = Flag;};
00435 bool Allocated() const {return(Allocated_);};
00436 int SetAllocated(bool Flag) {Allocated_ = Flag; return(0);};
00437
00438 private:
00439
00440
00441 int Allocate();
00442
00443 const Epetra_CrsMatrix &A_;
00444 const Ifpack_IlukGraph & Graph_;
00445 Epetra_CrsMatrix * U_;
00446 Epetra_Vector * D_;
00447 bool UseTranspose_;
00448
00449
00450 bool Allocated_;
00451 bool ValuesInitialized_;
00452 bool Factored_;
00453 double RelaxValue_;
00454 double Athresh_;
00455 double Rthresh_;
00456 mutable double Condest_;
00457
00458 mutable Epetra_MultiVector * OverlapX_;
00459 mutable Epetra_MultiVector * OverlapY_;
00460 Epetra_CombineMode OverlapMode_;
00461
00462
00463 };
00464
00466 ostream& operator << (ostream& os, const Ifpack_CrsRick& A);
00467
00468 #endif