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_TPETRABLOCKCRSMATRIX_HPP
00047 #define XPETRA_TPETRABLOCKCRSMATRIX_HPP
00048
00049
00050
00051 #include "Xpetra_TpetraConfigDefs.hpp"
00052
00053 #include "Tpetra_Experimental_BlockCrsMatrix.hpp"
00054 #include "Tpetra_CrsMatrix.hpp"
00055
00056 #include "Xpetra_CrsMatrix.hpp"
00057 #include "Xpetra_TpetraMap.hpp"
00058 #include "Xpetra_TpetraMultiVector.hpp"
00059 #include "Xpetra_TpetraVector.hpp"
00060 #include "Xpetra_TpetraCrsGraph.hpp"
00061
00062 #include "Xpetra_Exceptions.hpp"
00063
00064
00065 namespace Xpetra {
00066
00067 template <class Scalar, class LocalOrdinal = int, class GlobalOrdinal = LocalOrdinal, class Node = KokkosClassic::DefaultNode::DefaultNodeType>
00068 class TpetraBlockCrsMatrix
00069 : public CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
00070 {
00071
00072
00073 typedef TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> TpetraBlockCrsMatrixClass;
00074 typedef TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TpetraVectorClass;
00075 typedef TpetraImport<LocalOrdinal,GlobalOrdinal,Node> TpetraImportClass;
00076 typedef TpetraExport<LocalOrdinal,GlobalOrdinal,Node> TpetraExportClass;
00077
00078 public:
00079
00081
00082
00084 TpetraBlockCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
00085 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00086
00088 TpetraBlockCrsMatrix(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)
00089 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00090
00092 TpetraBlockCrsMatrix(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)
00093 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00094
00096 TpetraBlockCrsMatrix(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)
00097 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00098
00100 TpetraBlockCrsMatrix(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node> > &graph, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
00101 : mtx_(Teuchos::rcp(new Tpetra::Experimental::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(graph), params)))
00102 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00103
00105 TpetraBlockCrsMatrix(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node> > &graph, const LocalOrdinal blockSize)
00106 : mtx_(Teuchos::rcp(new Tpetra::Experimental::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(graph), blockSize))) { }
00107
00108
00109
00110
00112 TpetraBlockCrsMatrix(const Teuchos::RCP<const Tpetra::Experimental::BlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& sourceMatrix,
00113 const Import<LocalOrdinal,GlobalOrdinal,Node> & importer,
00114 const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
00115 const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
00116 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null)
00117 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00118
00120 TpetraBlockCrsMatrix(const Teuchos::RCP<const Tpetra::Experimental::BlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& sourceMatrix,
00121 const Export<LocalOrdinal,GlobalOrdinal,Node> & exporter,
00122 const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
00123 const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
00124 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null)
00125 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00126
00127
00129 virtual ~TpetraBlockCrsMatrix() { }
00130
00132
00134
00135
00137 void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
00138 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00139
00141 void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
00142 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00143
00145 void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
00146 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00147
00149 void replaceLocalValues (LocalOrdinal localRow,const ArrayView<const LocalOrdinal> &cols,const ArrayView<const Scalar> &vals)
00150 {
00151 XPETRA_MONITOR("TpetraBlockCrsMatrix::replaceLocalValues");
00152 mtx_->replaceLocalValues(localRow,cols.getRawPtr(),vals.getRawPtr(),cols.size());
00153 }
00154
00156 void setAllToScalar(const Scalar &alpha) { XPETRA_MONITOR("TpetraBlockCrsMatrix::setAllToScalar"); mtx_->setAllToScalar(alpha); }
00157
00159 void scale(const Scalar &alpha)
00160 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00161
00163
00164 void allocateAllValues(size_t numNonZeros,ArrayRCP<size_t> & rowptr, ArrayRCP<LocalOrdinal> & colind, ArrayRCP<Scalar> & values)
00165 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00166
00168 void setAllValues(const ArrayRCP<size_t> & rowptr, const ArrayRCP<LocalOrdinal> & colind, const ArrayRCP<Scalar> & values)
00169 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00170
00172 void getAllValues(ArrayRCP<const size_t>& rowptr, ArrayRCP<const LocalOrdinal>& colind, ArrayRCP<const Scalar>& values) const
00173 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00174
00176
00178
00179
00181 void resumeFill(const RCP< ParameterList > ¶ms=null) { }
00182
00184 void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > ¶ms=null) { }
00185
00187 void fillComplete(const RCP< ParameterList > ¶ms=null) { }
00188
00189
00191 void replaceDomainMapAndImporter(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >& newDomainMap, Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > & newImporter)
00192 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00193
00195 void expertStaticFillComplete(const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & domainMap,
00196 const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & rangeMap,
00197 const RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > &importer=Teuchos::null,
00198 const RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > &exporter=Teuchos::null,
00199 const RCP<ParameterList> ¶ms=Teuchos::null)
00200 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00201
00203
00205
00206
00208 const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getRowMap"); return toXpetra(mtx_->getRowMap()); }
00209
00211 const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getColMap"); return toXpetra(mtx_->getColMap()); }
00212
00214 RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node> > getCrsGraph() const
00215 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00216
00218 global_size_t getGlobalNumRows() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalNumRows"); return mtx_->getGlobalNumRows(); }
00219
00221 global_size_t getGlobalNumCols() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalNumCols"); return mtx_->getGlobalNumCols(); }
00222
00224 size_t getNodeNumRows() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getNodeNumRows"); return mtx_->getNodeNumRows(); }
00225
00227 size_t getNodeNumCols() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getNodeNumCols"); return mtx_->getNodeNumCols(); }
00228
00230 global_size_t getGlobalNumEntries() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalNumEntries"); return mtx_->getGlobalNumEntries(); }
00231
00233 size_t getNodeNumEntries() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getNodeNumEntries"); return mtx_->getNodeNumEntries(); }
00234
00236 size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getNumEntriesInLocalRow"); return mtx_->getNumEntriesInLocalRow(localRow); }
00237
00239 global_size_t getGlobalNumDiags() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalNumDiags"); return mtx_->getGlobalNumDiags(); }
00240
00242 size_t getNodeNumDiags() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getNodeNumDiags"); return mtx_->getNodeNumDiags(); }
00243
00245 size_t getGlobalMaxNumRowEntries() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalMaxNumRowEntries"); return mtx_->getGlobalMaxNumRowEntries(); }
00246
00248 size_t getNodeMaxNumRowEntries() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getNodeMaxNumRowEntries"); return mtx_->getNodeMaxNumRowEntries(); }
00249
00251 bool isLocallyIndexed() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::isLocallyIndexed"); return mtx_->isLocallyIndexed(); }
00252
00254 bool isGloballyIndexed() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::isGloballyIndexed"); return mtx_->isGloballyIndexed(); }
00255
00257 bool isFillComplete() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::isFillComplete"); return mtx_->isFillComplete(); }
00258
00260 bool isFillActive() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::isFillActive"); return false; }
00261
00263 typename ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getFrobeniusNorm"); return mtx_->getFrobeniusNorm(); }
00264
00266 bool supportsRowViews() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::supportsRowViews"); return mtx_->supportsRowViews(); }
00267
00269 void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getLocalRowCopy"); mtx_->getLocalRowCopy(LocalRow, Indices, Values, NumEntries); }
00270
00272 void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalRowView"); mtx_->getGlobalRowView(GlobalRow, indices, values); }
00273
00275 void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalRowCopy"); mtx_->getGlobalRowCopy(GlobalRow, indices, values, numEntries); }
00276
00278 void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getLocalRowView"); mtx_->getLocalRowView(LocalRow, indices, values); }
00279
00281
00283
00284
00286 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 { XPETRA_MONITOR("TpetraBlockCrsMatrix::apply"); mtx_->apply(toTpetra(X), toTpetra(Y), mode, alpha, beta); }
00287
00289 const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getDomainMap"); return toXpetra(mtx_->getDomainMap()); }
00290
00292 const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getRangeMap"); return toXpetra(mtx_->getRangeMap()); }
00293
00295
00297
00298
00300 std::string description() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::description"); return mtx_->description(); }
00301
00303 void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const { XPETRA_MONITOR("TpetraBlockCrsMatrix::describe"); mtx_->describe(out, verbLevel); }
00304
00306
00308 TpetraBlockCrsMatrix(const TpetraBlockCrsMatrix& matrix)
00309 : mtx_ (matrix.mtx_->template clone<Node> (matrix.mtx_->getNode ())) {}
00310
00312 void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const {
00313 XPETRA_MONITOR("TpetraBlockCrsMatrix::getLocalDiagCopy");
00314 XPETRA_DYNAMIC_CAST(TpetraVectorClass, diag, tDiag, "Xpetra::TpetraBlockCrsMatrix.getLocalDiagCopy() only accept Xpetra::TpetraVector as input arguments.")
00315
00316 mtx_->getLocalDiagCopy(*tDiag.getTpetra_Vector());
00317 }
00318
00320 void getLocalDiagOffsets(Teuchos::ArrayRCP<size_t> &offsets) const {
00321 XPETRA_MONITOR("TpetraBlockCrsMatrix::getLocalDiagOffsets");
00322 mtx_->getLocalDiagOffsets(offsets);
00323 }
00324
00326 void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag, const Teuchos::ArrayView<const size_t> &offsets) const
00327 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00328
00330
00331
00333 Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getMap"); return rcp( new TpetraMap< LocalOrdinal, GlobalOrdinal, Node >(mtx_->getMap()) ); }
00334
00336 void doImport(const DistObject<char, LocalOrdinal, GlobalOrdinal, Node> &source,
00337 const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
00338 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00339
00341 void doExport(const DistObject<char, LocalOrdinal, GlobalOrdinal, Node> &dest,
00342 const Import< LocalOrdinal, GlobalOrdinal, Node >& importer, CombineMode CM)
00343 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00344
00346 void doImport(const DistObject<char, LocalOrdinal, GlobalOrdinal, Node> &source,
00347 const Export< LocalOrdinal, GlobalOrdinal, Node >& exporter, CombineMode CM)
00348 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00349
00351 void doExport(const DistObject<char, LocalOrdinal, GlobalOrdinal, Node> &dest,
00352 const Export< LocalOrdinal, GlobalOrdinal, Node >& exporter, CombineMode CM)
00353 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00354
00355 void removeEmptyProcessesInPlace (const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >& newMap)
00356 {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
00357
00358
00359
00360 template<class Node2>
00361 RCP<TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node2> > clone(const RCP<Node2> &node2) const {
00362 return RCP<TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node2> >(new TpetraBlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node2>(mtx_->clone(node2)));
00363 }
00364
00366
00367
00369 bool hasMatrix() const { return !mtx_.is_null();}
00370
00372 TpetraBlockCrsMatrix(const Teuchos::RCP<Tpetra::Experimental::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &mtx) : mtx_(mtx) { }
00373
00375 RCP<const Tpetra::Experimental::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > getTpetra_BlockCrsMatrix() const { return mtx_; }
00376
00378 RCP<Tpetra::Experimental::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > getTpetra_BlockCrsMatrixNonConst() const { return mtx_; }
00379
00381
00382 private:
00383
00384 RCP< Tpetra::Experimental::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > mtx_;
00385
00386 };
00387
00388
00389
00390
00391
00392
00393
00394
00395
00396
00397 }
00398
00399 #define XPETRA_TPETRABLOCKCRSMATRIX_SHORT
00400 #endif // XPETRA_TPETRABLOCKCRSMATRIX_HPP