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_CRSRILUK_H_
00031 #define _IFPACK_CRSRILUK_H_
00032
00033 #include "Ifpack_ConfigDefs.h"
00034 #include "Ifpack_ScalingType.h"
00035 #include "Ifpack_IlukGraph.h"
00036 #include "Epetra_CompObject.h"
00037 #include "Epetra_Operator.h"
00038 #include "Epetra_CrsMatrix.h"
00039 #include "Epetra_Object.h"
00040 class Epetra_Comm;
00041 class Epetra_Map;
00042 class Epetra_CrsMatrix;
00043 class Epetra_VbrMatrix;
00044 class Epetra_RowMatrix;
00045 class Epetra_Vector;
00046 class Epetra_MultiVector;
00047
00048 #ifdef HAVE_IFPACK_TEUCHOS
00049 namespace Teuchos {
00050 class ParameterList;
00051 }
00052 #endif
00053
00055
00197 class Ifpack_CrsRiluk: public Epetra_Object, public Epetra_CompObject, public virtual Epetra_Operator {
00198
00199
00200
00201 friend ostream& operator << (ostream& os, const Ifpack_CrsRiluk& A);
00202
00203 public:
00205
00210 Ifpack_CrsRiluk(const Ifpack_IlukGraph & Graph);
00211
00213 Ifpack_CrsRiluk(const Ifpack_CrsRiluk & Matrix);
00214
00216 virtual ~Ifpack_CrsRiluk();
00217
00219
00225 int InitValues(const Epetra_CrsMatrix &A);
00226
00228
00234 int InitValues(const Epetra_VbrMatrix &A);
00235
00237 bool ValuesInitialized() const {return(ValuesInitialized_);};
00238
00240 void SetRelaxValue( double RelaxValue) {RelaxValue_ = RelaxValue; return;}
00241
00243 void SetAbsoluteThreshold( double Athresh) {Athresh_ = Athresh; return;}
00244
00246 void SetRelativeThreshold( double Rthresh) {Rthresh_ = Rthresh; return;}
00247
00249 void SetOverlapMode( Epetra_CombineMode OverlapMode) {OverlapMode_ = OverlapMode; return;}
00250
00251 #ifdef HAVE_IFPACK_TEUCHOS
00252
00253
00254
00255
00256
00257
00258
00259
00260 int SetParameters(const Teuchos::ParameterList& parameterlist,
00261 bool cerr_warning_if_unused=false);
00262 #endif
00263
00265
00273 int Factor();
00274
00276 bool Factored() const {return(Factored_);};
00277
00278
00279
00280
00281
00283
00293 int Solve(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
00294
00296
00306 int Multiply(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
00307
00309
00317 int Condest(bool Trans, double & ConditionNumberEstimate) const;
00318
00319
00321 double GetRelaxValue() {return RelaxValue_;}
00322
00324 double GetAbsoluteThreshold() {return Athresh_;}
00325
00327 double GetRelativeThreshold() {return Rthresh_;}
00328
00330 Epetra_CombineMode GetOverlapMode() {return OverlapMode_;}
00331
00332
00334 int NumGlobalRows() const {return(Graph().NumGlobalRows());};
00335
00337 int NumGlobalCols() const {return(Graph().NumGlobalCols());};
00338
00340 int NumGlobalNonzeros() const {return(L().NumGlobalNonzeros()+U().NumGlobalNonzeros());};
00341
00343 virtual int NumGlobalBlockDiagonals() const {return(Graph().NumGlobalBlockDiagonals());};
00344
00346 int NumMyRows() const {return(Graph().NumMyRows());};
00347
00349 int NumMyCols() const {return(Graph().NumMyCols());};
00350
00352 int NumMyNonzeros() const {return(L().NumMyNonzeros()+U().NumMyNonzeros());};
00353
00355 virtual int NumMyBlockDiagonals() const {return(Graph().NumMyBlockDiagonals());};
00356
00358 virtual int NumMyDiagonals() const {return(NumMyDiagonals_);};
00359
00361 int IndexBase() const {return(Graph().IndexBase());};
00362
00364 const Ifpack_IlukGraph & Graph() const {return(Graph_);};
00365
00367 const Epetra_CrsMatrix & L() const {return(*L_);};
00368
00370 const Epetra_Vector & D() const {return(*D_);};
00371
00373 const Epetra_CrsMatrix & U() const {return(*U_);};
00374
00376
00378 const char * Label() const {return(Epetra_Object::Label());};
00379
00381
00390 int SetUseTranspose(bool UseTranspose) {UseTranspose_ = UseTranspose; return(0);};
00391
00393
00404 int Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const {
00405 return(Multiply(Ifpack_CrsRiluk::UseTranspose(), X, Y));};
00406
00408
00421 int ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const {
00422 return(Solve(Ifpack_CrsRiluk::UseTranspose(), X, Y));};
00423
00425 double NormInf() const {return(0.0);};
00426
00428 bool HasNormInf() const {return(false);};
00429
00431 bool UseTranspose() const {return(UseTranspose_);};
00432
00434 const Epetra_Map & OperatorDomainMap() const {return(U_->OperatorDomainMap());};
00435
00437 const Epetra_Map & OperatorRangeMap() const{return(L_->OperatorRangeMap());};
00438
00440 const Epetra_Comm & Comm() const{return(Comm_);};
00442
00443 protected:
00444 void SetFactored(bool Flag) {Factored_ = Flag;};
00445 void SetValuesInitialized(bool Flag) {ValuesInitialized_ = Flag;};
00446 bool Allocated() const {return(Allocated_);};
00447 int SetAllocated(bool Flag) {Allocated_ = Flag; return(0);};
00448 int BlockGraph2PointGraph(const Epetra_CrsGraph & BG, Epetra_CrsGraph & PG, bool Upper);
00449
00450 private:
00451
00452
00453 int AllocateCrs();
00454 int AllocateVbr();
00455 int InitAllValues(const Epetra_RowMatrix & A, int MaxNumEntries);
00456 int BlockMap2PointMap(const Epetra_BlockMap & BlockMap, Epetra_Map * & PointMap);
00457 int GenerateXY(bool Trans,
00458 const Epetra_MultiVector& Xin, const Epetra_MultiVector& Yin,
00459 Epetra_MultiVector * & Xout, Epetra_MultiVector * & Yout) const;
00460 bool UserMatrixIsVbr_;
00461 bool UserMatrixIsCrs_;
00462 bool IsOverlapped_;
00463 const Ifpack_IlukGraph & Graph_;
00464 Epetra_Map * IlukRowMap_;
00465 Epetra_Map * IlukDomainMap_;
00466 Epetra_Map * IlukRangeMap_;
00467 const Epetra_Map * U_DomainMap_;
00468 const Epetra_Map * L_RangeMap_;
00469 const Epetra_Comm & Comm_;
00470 Epetra_CrsMatrix * L_;
00471 Epetra_CrsMatrix * U_;
00472 Epetra_CrsGraph * L_Graph_;
00473 Epetra_CrsGraph * U_Graph_;
00474 Epetra_Vector * D_;
00475 bool UseTranspose_;
00476
00477 int NumMyDiagonals_;
00478 bool Allocated_;
00479 bool ValuesInitialized_;
00480 bool Factored_;
00481 double RelaxValue_;
00482 double Athresh_;
00483 double Rthresh_;
00484 mutable double Condest_;
00485
00486 mutable Epetra_MultiVector * OverlapX_;
00487 mutable Epetra_MultiVector * OverlapY_;
00488 mutable Epetra_MultiVector * VbrX_;
00489 mutable Epetra_MultiVector * VbrY_;
00490 Epetra_CombineMode OverlapMode_;
00491
00492
00493 };
00494
00496 ostream& operator << (ostream& os, const Ifpack_CrsRiluk& A);
00497
00498 #endif