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
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046 #ifndef XPETRA_EPETRACRSMATRIX_HPP
00047 #define XPETRA_EPETRACRSMATRIX_HPP
00048
00049
00050
00051 #include "Xpetra_EpetraConfigDefs.hpp"
00052
00053 #include "Xpetra_CrsMatrix.hpp"
00054
00055 #include <Epetra_CrsMatrix.h>
00056 #include <Epetra_Map.h>
00057
00058 #include "Xpetra_EpetraMap.hpp"
00059 #include "Xpetra_EpetraVector.hpp"
00060 #include "Xpetra_EpetraMultiVector.hpp"
00061 #include "Xpetra_EpetraCrsGraph.hpp"
00062
00063 #include "Xpetra_Utils.hpp"
00064 #include "Xpetra_Exceptions.hpp"
00065
00066 namespace Xpetra {
00067
00068 class EpetraCrsMatrix
00069 : public CrsMatrix<double, int, int>
00070 {
00071
00072 typedef double Scalar;
00073 typedef int LocalOrdinal;
00074 typedef int GlobalOrdinal;
00075 typedef Kokkos::DefaultNode::DefaultNodeType Node;
00076 typedef Kokkos::DefaultKernels<void,int,Node>::SparseOps LocalMatOps;
00077
00078 public:
00079
00081
00082
00084 EpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null);
00085
00087 EpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null);
00088
00090 EpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null);
00091
00093 EpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null);
00094
00096 EpetraCrsMatrix(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node, LocalMatOps > > &graph, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null);
00097
00099 virtual ~EpetraCrsMatrix() { }
00100
00102
00104
00105
00107 void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals);
00108
00110 void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals);
00111
00113 void setAllToScalar(const Scalar &alpha) { XPETRA_MONITOR("EpetraCrsMatrix::setAllToScalar"); mtx_->PutScalar(alpha); }
00114
00116 void scale(const Scalar &alpha) { XPETRA_MONITOR("EpetraCrsMatrix::scale"); mtx_->Scale(alpha); }
00117
00119
00121
00122
00124 void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > ¶ms=null);
00125
00127 void fillComplete(const RCP< ParameterList > ¶ms=null);
00128
00130
00132
00133
00135 const RCP< const Comm< int > > getComm() const { XPETRA_MONITOR("EpetraCrsMatrix::getComm"); return toXpetra(mtx_->Comm()); }
00136
00138 const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const { XPETRA_MONITOR("EpetraCrsMatrix::getRowMap"); return toXpetra(mtx_->RowMap()); }
00139
00141 const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const { XPETRA_MONITOR("EpetraCrsMatrix::getColMap"); return toXpetra(mtx_->ColMap()); }
00142
00144 RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node, LocalMatOps > > getCrsGraph() const { XPETRA_MONITOR("EpetraCrsMatrix::getCrsGraph"); return toXpetra(mtx_->Graph()); }
00145
00147 global_size_t getGlobalNumRows() const { XPETRA_MONITOR("EpetraCrsMatrix::getGlobalNumRows"); return mtx_->NumGlobalRows(); }
00148
00150 global_size_t getGlobalNumCols() const { XPETRA_MONITOR("EpetraCrsMatrix::getGlobalNumCols"); return mtx_->NumGlobalCols(); }
00151
00153 size_t getNodeNumRows() const { XPETRA_MONITOR("EpetraCrsMatrix::getNodeNumRows"); return mtx_->NumMyRows(); }
00154
00156 size_t getNodeNumCols() const { XPETRA_MONITOR("EpetraCrsMatrix::getNodeNumCols"); return mtx_->NumMyCols(); }
00157
00159 global_size_t getGlobalNumEntries() const { XPETRA_MONITOR("EpetraCrsMatrix::getGlobalNumEntries"); return mtx_->NumGlobalNonzeros(); }
00160
00162 size_t getNodeNumEntries() const { XPETRA_MONITOR("EpetraCrsMatrix::getNodeNumEntries"); return mtx_->NumMyNonzeros(); }
00163
00165 size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const { XPETRA_MONITOR("EpetraCrsMatrix::getNumEntriesInLocalRow"); return mtx_->NumMyEntries(localRow); }
00166
00168 global_size_t getGlobalNumDiags() const { XPETRA_MONITOR("EpetraCrsMatrix::getGlobalNumDiags"); return mtx_->NumGlobalDiagonals(); }
00169
00171 size_t getNodeNumDiags() const { XPETRA_MONITOR("EpetraCrsMatrix::getNodeNumDiags"); return mtx_->NumMyDiagonals(); }
00172
00174 size_t getGlobalMaxNumRowEntries() const { XPETRA_MONITOR("EpetraCrsMatrix::getGlobalMaxNumRowEntries"); return mtx_->GlobalMaxNumEntries(); }
00175
00177 size_t getNodeMaxNumRowEntries() const { XPETRA_MONITOR("EpetraCrsMatrix::getNodeMaxNumRowEntries"); return mtx_->MaxNumEntries(); }
00178
00180 bool isLocallyIndexed() const { XPETRA_MONITOR("EpetraCrsMatrix::isLocallyIndexed"); return mtx_->IndicesAreLocal(); }
00181
00183 bool isGloballyIndexed() const { XPETRA_MONITOR("EpetraCrsMatrix::isGloballyIndexed"); return mtx_->IndicesAreGlobal(); }
00184
00186 bool isFillComplete() const { XPETRA_MONITOR("EpetraCrsMatrix::isFillComplete"); return mtx_->Filled(); }
00187
00189 ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const { XPETRA_MONITOR("EpetraCrsMatrix::getFrobeniusNorm"); return mtx_->NormFrobenius(); }
00190
00192 void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const;
00193
00195 void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const;
00196
00198 void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const;
00199
00201 void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const { XPETRA_MONITOR("EpetraCrsMatrix::getLocalDiagCopy"); mtx_->ExtractDiagonalCopy(toEpetra(diag)); }
00202
00204
00206
00207
00209 void apply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const;
00210
00212 const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const { XPETRA_MONITOR("EpetraCrsMatrix::getDomainMap"); return toXpetra(mtx_->DomainMap()); }
00213
00215 const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const { XPETRA_MONITOR("EpetraCrsMatrix::getRangeMap"); return toXpetra(mtx_->RangeMap()); }
00216
00218
00220
00221
00223 std::string description() const;
00224
00226 void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const;
00227
00229
00231
00232
00234 const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const { XPETRA_MONITOR("EpetraCrsMatrix::getMap"); return toXpetra(mtx_->Map()); }
00235
00237 void doImport(const DistObject<char, LocalOrdinal, GlobalOrdinal, Node> &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM);
00238
00240 void doExport(const DistObject<char, LocalOrdinal, GlobalOrdinal, Node> &dest, const Import< LocalOrdinal, GlobalOrdinal, Node >& importer, CombineMode CM);
00241
00243 void doImport(const DistObject<char, LocalOrdinal, GlobalOrdinal, Node> &source, const Export< LocalOrdinal, GlobalOrdinal, Node >& exporter, CombineMode CM);
00244
00246 void doExport(const DistObject<char, LocalOrdinal, GlobalOrdinal, Node> &dest, const Export< LocalOrdinal, GlobalOrdinal, Node >& exporter, CombineMode CM);
00247
00249
00251
00252
00254 EpetraCrsMatrix(const Teuchos::RCP<Epetra_CrsMatrix > &mtx) : mtx_(mtx) { }
00255
00257 RCP<const Epetra_CrsMatrix> getEpetra_CrsMatrix() const { return mtx_; }
00258
00260 RCP<Epetra_CrsMatrix> getEpetra_CrsMatrixNonConst() const { return mtx_; }
00261
00263
00264 private:
00265
00266 RCP<Epetra_CrsMatrix> mtx_;
00267
00268 };
00269
00270 }
00271
00272 #define XPETRA_EPETRACRSMATRIX_SHORT
00273 #endif // XPETRA_EPETRACRSMATRIX_HPP