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
00031 #ifndef IFPACK_OVERLAPPINGROWMATRIX_H
00032 #define IFPACK_OVERLAPPINGROWMATRIX_H
00033
00034 #include "Ifpack_ConfigDefs.h"
00035 #include "Epetra_RowMatrix.h"
00036 #include "Epetra_CombineMode.h"
00037 #include "Teuchos_RefCountPtr.hpp"
00038 #include "Epetra_Import.h"
00039 #ifdef IFPACK_NODE_AWARE_CODE
00040 #include "Epetra_IntVector.h"
00041 #endif
00042
00043 class Epetra_Map;
00044 class Epetra_BlockMap;
00045 class Epetra_CrsMatrix;
00046 class Epetra_Comm;
00047
00049
00050 class Ifpack_OverlappingRowMatrix : public virtual Epetra_RowMatrix {
00051
00052 public:
00053
00055 # ifdef IFPACK_NODE_AWARE_CODE
00056 Ifpack_OverlappingRowMatrix(const Teuchos::RefCountPtr<const Epetra_RowMatrix>& Matrix_in,
00057 int OverlapLevel_in, int myNodeID);
00058 # endif
00059 Ifpack_OverlappingRowMatrix(const Teuchos::RefCountPtr<const Epetra_RowMatrix>& Matrix_in,
00060 int OverlapLevel_in);
00061
00062 # ifdef IFPACK_NODE_AWARE_CODE
00063 ~Ifpack_OverlappingRowMatrix();
00064 # else
00065 ~Ifpack_OverlappingRowMatrix() {};
00066 # endif
00067
00068
00070
00072
00080 virtual int NumMyRowEntries(int MyRow, int & NumEntries) const;
00081
00083 virtual int MaxNumEntries() const
00084 {
00085 return(MaxNumEntries_);
00086 }
00087
00089
00103 virtual int ExtractMyRowCopy(int MyRow, int Length, int & NumEntries, double *Values, int * Indices) const;
00104 #ifdef IFPACK_NODE_AWARE_CODE
00105 virtual int ExtractGlobalRowCopy(int MyRow, int Length, int & NumEntries, double* Values, int* Indices) const;
00106 #endif
00107
00109
00115 virtual int ExtractDiagonalCopy(Epetra_Vector & Diagonal) const;
00117
00119
00121
00131 virtual int Multiply(bool TransA, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
00132
00134 virtual int Solve(bool Upper, bool Trans, bool UnitDiagonal, const Epetra_MultiVector& X,
00135 Epetra_MultiVector& Y) const
00136 {
00137 IFPACK_RETURN(-1);
00138 }
00139
00140 virtual int Apply(const Epetra_MultiVector& X,
00141 Epetra_MultiVector& Y) const;
00142
00143 virtual int ApplyInverse(const Epetra_MultiVector& X,
00144 Epetra_MultiVector& Y) const;
00146 virtual int InvRowSums(Epetra_Vector& x) const
00147 {
00148 IFPACK_RETURN(-1);
00149 }
00150
00152 virtual int LeftScale(const Epetra_Vector& x)
00153 {
00154 IFPACK_RETURN(-1);
00155 }
00156
00158 virtual int InvColSums(Epetra_Vector& x) const
00159 {
00160 IFPACK_RETURN(-1);
00161 }
00162
00163
00165 virtual int RightScale(const Epetra_Vector& x)
00166 {
00167 IFPACK_RETURN(-1);
00168 }
00169
00171
00173
00175 virtual bool Filled() const
00176 {
00177 return(true);
00178 }
00179
00181
00182
00183
00184 virtual double NormInf() const
00185 {
00186 return(A().NormInf());
00187 }
00188
00190
00191
00192
00193 virtual double NormOne() const
00194 {
00195 IFPACK_RETURN(A().NormOne());
00196 }
00197
00199 virtual int NumGlobalNonzeros() const
00200 {
00201 return(NumGlobalNonzeros_);
00202 }
00203
00205 virtual int NumGlobalRows() const
00206 {
00207 return(A().NumGlobalRows());
00208 }
00209
00211 virtual int NumGlobalCols() const
00212 {
00213 return(A().NumGlobalCols());
00214 }
00215
00217 virtual int NumGlobalDiagonals() const
00218 {
00219 return(A().NumGlobalDiagonals());
00220 }
00221
00223 virtual int NumMyNonzeros() const
00224 {
00225 return(NumMyNonzeros_);
00226 }
00227
00229 virtual int NumMyRows() const
00230 {
00231 return(NumMyRows_);
00232 }
00233
00235 virtual int NumMyCols() const
00236 {
00237 return(NumMyCols_);
00238 }
00239
00241 virtual int NumMyDiagonals() const
00242 {
00243 return(NumMyDiagonals_);
00244 }
00245
00247 virtual bool LowerTriangular() const
00248 {
00249 return(A().LowerTriangular());
00250 }
00251
00253 virtual bool UpperTriangular() const
00254 {
00255 return(A().UpperTriangular());
00256 }
00257
00259 virtual const Epetra_Map & RowMatrixRowMap() const
00260 {
00261 return(*Map_);
00262 }
00263
00265 virtual const Epetra_Map & RowMatrixColMap() const
00266 {
00267 # ifdef IFPACK_NODE_AWARE_CODE
00268 return(*colMap_);
00269 # else
00270 return(*Map_);
00271 # endif
00272 }
00273
00275 virtual const Epetra_Import * RowMatrixImporter() const
00276 {
00277 return(&*Importer_);
00278 }
00280
00281
00282
00284 int SetOwnership(bool ownership)
00285 {
00286 IFPACK_RETURN(-1);
00287 }
00288
00290 int SetUseTranspose(bool UseTranspose_in)
00291 {
00292 UseTranspose_ = UseTranspose_in;
00293 return(0);
00294 }
00295
00297 bool UseTranspose() const
00298 {
00299 return(UseTranspose_);
00300 }
00301
00303 bool HasNormInf() const
00304 {
00305 return(A().HasNormInf());
00306 }
00307
00309 const Epetra_Comm & Comm() const
00310 {
00311 return(A().Comm());
00312 }
00313
00315 const Epetra_Map & OperatorDomainMap() const
00316 {
00317 return(*Map_);
00318 }
00319
00321 const Epetra_Map & OperatorRangeMap() const
00322 {
00323 return(*Map_);
00324 }
00326
00327 const Epetra_BlockMap& Map() const;
00328
00329 const char* Label() const{
00330 return(Label_.c_str());
00331 };
00332
00333 int OverlapLevel() const
00334 {
00335 return(OverlapLevel_);
00336 }
00337
00338 int ImportMultiVector(const Epetra_MultiVector& X,
00339 Epetra_MultiVector& OvX,
00340 Epetra_CombineMode CM = Insert);
00341
00342 int ExportMultiVector(const Epetra_MultiVector& OvX,
00343 Epetra_MultiVector& X,
00344 Epetra_CombineMode CM = Add);
00345 #ifdef IFPACK_NODE_AWARE_CODE
00346 inline const Epetra_RowMatrix& A() const
00347 {
00348 return(*Matrix_);
00349 }
00350
00351 inline Epetra_CrsMatrix& B() const
00352 {
00353 return(*ExtMatrix_);
00354 }
00355 #endif
00356
00357 private:
00358 #ifndef IFPACK_NODE_AWARE_CODE
00359 inline const Epetra_RowMatrix& A() const
00360 {
00361 return(*Matrix_);
00362 }
00363
00364 inline Epetra_RowMatrix& B() const;
00365
00366 #endif
00367
00368 int NumMyRows_;
00369 int NumMyCols_;
00370 int NumMyDiagonals_;
00371 int NumMyNonzeros_;
00372
00373 int NumGlobalNonzeros_;
00374 int MaxNumEntries_;
00375
00376 int NumMyRowsA_;
00377 int NumMyRowsB_;
00378
00379 bool UseTranspose_;
00380
00381 Teuchos::RefCountPtr<const Epetra_Map> Map_;
00382 # ifdef IFPACK_NODE_AWARE_CODE
00383
00384 const Epetra_Map *colMap_;
00385
00386 # endif
00387 Teuchos::RefCountPtr<const Epetra_Import> Importer_;
00388
00389 Teuchos::RefCountPtr<const Epetra_RowMatrix> Matrix_;
00390 Teuchos::RefCountPtr<Epetra_CrsMatrix> ExtMatrix_;
00391 Teuchos::RefCountPtr<Epetra_Map> ExtMap_;
00392 Teuchos::RefCountPtr<Epetra_Import> ExtImporter_;
00393
00394 int OverlapLevel_;
00395 string Label_;
00396
00397 };
00398
00399 #endif // IFPACK_OVERLAPPINGROWMATRIX_H