PlayaEpetraMatrix.hpp
Go to the documentation of this file.
00001 /* @HEADER@ */
00002 //   
00003  /* @HEADER@ */
00004 
00005 #ifndef PLAYA_EPETRAMATRIX_HPP
00006 #define PLAYA_EPETRAMATRIX_HPP
00007 
00008 #include "PlayaLoadableMatrix.hpp"
00009 #include "PlayaLinearOperatorDecl.hpp"
00010 #include "PlayaLinearOpWithSpacesDecl.hpp"
00011 #include "PlayaRowAccessibleOp.hpp"
00012 #include "PlayaPrintable.hpp"
00013 #include "PlayaILUFactorizableOp.hpp"
00014 #include "PlayaICCFactorizableOp.hpp"
00015 #include "Epetra_CrsMatrix.h"
00016 #include "PlayaLinearOpWithSpacesDecl.hpp"
00017 
00018 namespace Playa
00019 {
00020 using namespace Teuchos;
00021 
00022 
00023 
00024 /** Playa wrapper for epetra matrix */
00025 class EpetraMatrix : public LinearOpWithSpaces<double>,
00026                      public LoadableMatrix<double>,
00027                      public RowAccessibleOp<double>,
00028                      public Printable,
00029                      public ILUFactorizableOp<double>,
00030                      public ICCFactorizableOp<double>
00031 {
00032 public:
00033 
00034   /** Construct an empty EpetraMatrix structured according to the graph 
00035    * argument */
00036   EpetraMatrix(const Epetra_CrsGraph& graph,
00037     const VectorSpace<double>& domain,
00038     const VectorSpace<double>& range);
00039 
00040   /** Wrap an existing Epetra CRS Matrix */
00041   EpetraMatrix(const RCP<Epetra_CrsMatrix>& mat,
00042     const VectorSpace<double>& domain,
00043     const VectorSpace<double>& range);
00044 
00045   /** Apply the operator */
00046   virtual void apply(
00047     Teuchos::ETransp applyType,
00048     const Vector<double>& in,
00049     Vector<double> out) const ;
00050   
00051 
00052 
00053   /** \name LoadableMatrix interface functions */
00054   //@{
00055   /** Insert a set of elements in a row, adding to any previously
00056    * existing values. 
00057    * @param globalRowIndex the global index of the row to which these
00058    * elements belong.
00059    * @param nElemsToInsert the number of elements being inserted in this
00060    * step
00061    * @param globalColumnIndices array of column indices. Must 
00062    * be nElemsToInsert in length. 
00063    * @param elements array of element values. Must be nElemsToInsert in
00064    * length
00065    */
00066   virtual void addToRow(int globalRowIndex,
00067     int nElemsToInsert,
00068     const int* globalColumnIndices,
00069     const double* elementValues) ;
00070 
00071 
00072   /** 
00073    * Add to a batch of elements
00074    */
00075   virtual void addToElementBatch(int numRows, 
00076     int rowBlockSize,
00077     const int* globalRowIndices,
00078     int numColumnsPerRow,
00079     const int* globalColumnIndices,
00080     const double* values,
00081     const int* skipRow);
00082 
00083   /** Set all elements to zero, preserving the existing structure */
00084   virtual void zero() ;
00085   
00086   //@}
00087 
00088   /** \name incomplete factorization preconditioning interface */
00089   //@{
00090   /** create an incomplete factorization. 
00091    * @param fillLevels number of levels of fill on the local processor
00092    * @param overlapFill number of levels of fill on remote processors
00093    * @param relaxationValue fraction of dropped values to be added to the
00094    * diagonal
00095    * @param relativeThreshold relative diagonal perutrbation
00096    * @param absoluteThreshold absolute diagonal perturbation
00097    * @param leftOrRight whether this preconditioner is to be applied
00098    * from the left or right 
00099    * @param rtn newly created preconditioner, returned 
00100    * by reference argument.
00101    */
00102   virtual void getILUKPreconditioner(int fillLevels,
00103     int overlapFill,
00104     double relaxationValue,
00105     double relativeThreshold,
00106     double absoluteThreshold,
00107     LeftOrRight leftOrRight,
00108     Preconditioner<double>& rtn) const ;
00109   //@}
00110 
00111 
00112    /** \name incomplete factorization preconditioning interface */
00113   //@{
00114   /** create an incomplete factorization. 
00115    * @param fillLevels number of levels of fill on the local processor
00116    * @param overlapFill number of levels of fill on remote processors
00117    * @param relaxationValue fraction of dropped values to be added to the
00118    * diagonal
00119    * @param relativeThreshold relative diagonal perutrbation
00120    * @param absoluteThreshold absolute diagonal perturbation
00121    * @param rtn newly created preconditioner, returned 
00122    * by reference argument.
00123    */
00124   virtual void getICCPreconditioner(int fillLevels,
00125     int overlapFill,
00126     double relaxationValue,
00127     double relativeThreshold,
00128     double absoluteThreshold,
00129     Preconditioner<double>& rtn) const ;
00130 
00131   
00132   /** \name Row access interface */
00133   //@{
00134   /** Get the specified row as defined by RowAccessible  */
00135   void getRow(const int& row, 
00136     Teuchos::Array<int>& indices, 
00137     Teuchos::Array<double>& values) const;
00138   //@}
00139   
00140 
00141   /** \name Diagnostic output */
00142   //@{
00143   /** Print the matrix */
00144   virtual void print(std::ostream& os) const ;
00145   //@}
00146 
00147   /** */
00148   std::ostream& describe(
00149     std::ostream                         &out
00150     ,const Teuchos::EVerbosityLevel      verbLevel
00151     ,const std::string                   leadingIndent
00152     , const std::string                   indentSpacer
00153     ) const 
00154     {
00155       out << leadingIndent << indentSpacer << this->description() << std::endl;
00156       return out;
00157     }
00158   /** */
00159   std::string description() const ;
00160   //@}
00161 
00162   /** */
00163   static Epetra_CrsMatrix& getConcrete(const LinearOperator<double>& A);
00164 
00165   /** */
00166   static RCP<const Epetra_CrsMatrix> getConcretePtr(const LinearOperator<double>& A);
00167 
00168   /** 
00169    * Read-only access to the underlying crs matrix. Needed by Ifpack and ML.
00170    */
00171   const Epetra_CrsMatrix* crsMatrix() const ;
00172 
00173 private:
00174 
00175   Epetra_CrsMatrix* crsMatrix();
00176 
00177   RCP<Epetra_CrsMatrix> matrix_;
00178 
00179   const Epetra_Map& getRangeMap() const;
00180   const Epetra_Map& getDomainMap() const;
00181 };
00182 
00183   
00184 
00185 }
00186 
00187 #endif

Site Contact