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_TPETRAVBRMATRIX_HPP
00047 #define XPETRA_TPETRAVBRMATRIX_HPP
00048
00049 #include "Xpetra_ConfigDefs.hpp"
00050
00051 #ifndef HAVE_XPETRA_TPETRA
00052 #error This file should be included only if HAVE_XPETRA_TPETRA is defined.
00053 #endif
00054
00055 #include "Xpetra_VbrMatrix.hpp"
00056
00057 #include <Tpetra_VbrMatrix.hpp>
00058
00059 namespace Xpetra {
00060
00061 template <class Scalar,
00062 class LocalOrdinal = int,
00063 class GlobalOrdinal = LocalOrdinal,
00064 class Node = Kokkos::DefaultNode::DefaultNodeType,
00065 class LocalMatOps = typename Kokkos::DefaultKernels<Scalar,LocalOrdinal,Node>::BlockSparseOps >
00066 class TpetraVbrMatrix :
00067 public VbrMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node, LocalOrdinal> {
00068 public:
00069
00071
00072
00073 TpetraVbrMatrix(const Teuchos::RCP<const Tpetra::VbrMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps> > &mtx) : mtx_(mtx) { }
00074
00076 virtual ~TpetraVbrMatrix();
00077
00079
00081
00082
00084
00086 const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & getDomainMap() const { return mtx_->getDomainMap(); }
00087
00089
00091 const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & getRangeMap() const { return mtx_->getRangeMap(); }
00092
00094
00099 void apply(const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &X,
00100 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y,
00101 Teuchos::ETransp trans = Teuchos::NO_TRANS,
00102 Scalar alpha = Teuchos::ScalarTraits<Scalar>::one(),
00103 Scalar beta = Teuchos::ScalarTraits<Scalar>::zero()) const { mtx_->apply(X,Y,trans,alpha,beta); }
00104
00106
00109 void applyInverse(const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> & Y,
00110 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &X,
00111 Teuchos::ETransp trans) const { mtx_->applyInverse(Y,X,trans); }
00112
00114 bool hasTransposeApply() const { return mtx_->hasTransposeApply(); }
00115
00117
00119
00120
00122 const Teuchos::RCP<const BlockMap<LocalOrdinal,GlobalOrdinal,Node> > & getBlockRowMap() const { return mtx_->getBlockRowMap(); }
00123
00125 const Teuchos::RCP<const BlockMap<LocalOrdinal,GlobalOrdinal,Node> > & getBlockColMap() const { return mtx_->getBlockColMap(); }
00126
00128 const Teuchos::RCP<const BlockMap<LocalOrdinal,GlobalOrdinal,Node> > & getBlockDomainMap() const { return mtx_->getBlockDomainMap(); }
00129
00131 const Teuchos::RCP<const BlockMap<LocalOrdinal,GlobalOrdinal,Node> > & getBlockRangeMap() const { return mtx_->getBlockRangeMap(); }
00132
00134 const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & getPointRowMap() const { return mtx_->getPointRowMap(); }
00135
00137 const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & getPointColMap() const { return mtx_->getPointColMap(); }
00138
00140 bool isFillComplete() const { return mtx_->isFillComplete(); }
00142
00144
00145
00147
00150 void putScalar(Scalar s) { mtx_->putScalar(s); }
00151
00153
00162 void setGlobalBlockEntry(GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, const Teuchos::SerialDenseMatrix<GlobalOrdinal,Scalar>& blockEntry) { mtx_->setGlobalBlockEntry(globalBlockRow, globalBlockCol, blockEntry); }
00163
00165
00172 void setLocalBlockEntry(LocalOrdinal localBlockRow, LocalOrdinal localBlockCol, const Teuchos::SerialDenseMatrix<LocalOrdinal,Scalar>& blockEntry) { mtx_->setLocalBlockEntry(localBlockRow, localBlockCol, blockEntry); }
00173
00175
00184 void sumIntoGlobalBlockEntry(GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, const Teuchos::SerialDenseMatrix<GlobalOrdinal,Scalar>& blockEntry) { mtx_->sumIntoGlobalBlockEntry(globalBlockRow, globalBlockCol, blockEntry); }
00185
00187
00194 void sumIntoLocalBlockEntry(LocalOrdinal localBlockRow, LocalOrdinal localBlockCol, const Teuchos::SerialDenseMatrix<LocalOrdinal,Scalar>& blockEntry) { mtx_->sumIntoLocalBlockEntry(localBlockRow, localBlockCol, blockEntry); }
00195
00197
00206 void setGlobalBlockEntry(GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, LocalOrdinal blkRowSize, LocalOrdinal blkColSize, LocalOrdinal LDA, const Teuchos::ArrayView<const Scalar>& blockEntry) { mtx_->setGlobalBlockEntry(globalBlockRow, globalBlockCol, blkRowSize, blkColSize, LDA, blockEntry); }
00207
00209
00216 void setLocalBlockEntry(LocalOrdinal localBlockRow, LocalOrdinal localBlockCol, LocalOrdinal blkRowSize, LocalOrdinal blkColSize, LocalOrdinal LDA, const Teuchos::ArrayView<const Scalar>& blockEntry) { mtx_->setLocalBlockEntry(localBlockRow, localBlockCol, blkRowSize, blkColSize, LDA, blockEntry); }
00217
00219
00228 void sumIntoGlobalBlockEntry(GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, LocalOrdinal blkRowSize, LocalOrdinal blkColSize, LocalOrdinal LDA, const Teuchos::ArrayView<const Scalar>& blockEntry) { mtx_->sumIntoGlobalBlockEntry(globalBlockRow, globalBlockCol, blkRowSize, blkColSize, LDA, blockEntry); }
00229
00231
00238 void sumIntoLocalBlockEntry(LocalOrdinal localBlockRow, LocalOrdinal localBlockCol, LocalOrdinal blkRowSize, LocalOrdinal blkColSize, LocalOrdinal LDA, const Teuchos::ArrayView<const Scalar>& blockEntry) { mtx_->sumIntoLocalBlockEntry(localBlockRow, localBlockCol, blkRowSize, blkColSize, LDA, blockEntry); }
00239
00241
00243
00244
00246
00249 void fillComplete(const Teuchos::RCP<const BlockMap<LocalOrdinal,GlobalOrdinal,Node> >& blockDomainMap, const Teuchos::RCP<const BlockMap<LocalOrdinal,GlobalOrdinal,Node> >& blockRangeMap) { mtx_->fillComplete(blockDomainMap, blockRangeMap); }
00250
00252
00255 void fillComplete() { mtx_->fillComplete(); }
00257
00259
00260
00262
00270 void getGlobalBlockEntryView(GlobalOrdinal globalBlockRow,
00271 GlobalOrdinal globalBlockCol,
00272 LocalOrdinal& numPtRows,
00273 LocalOrdinal& numPtCols,
00274 Teuchos::ArrayRCP<const Scalar>& blockEntry) const { mtx_->getGlobalBlockEntryView(globalBlockRow, globalBlockCol, numPtRows, numPtCols, blockEntry); }
00275
00277
00287 void getGlobalBlockEntryViewNonConst(GlobalOrdinal globalBlockRow,
00288 GlobalOrdinal globalBlockCol,
00289 LocalOrdinal& numPtRows,
00290 LocalOrdinal& numPtCols,
00291 Teuchos::ArrayRCP<Scalar>& blockEntry) { mtx_->getGlobalBlockEntryViewNonConst(globalBlockRow, globalBlockCol, numPtRows, numPtCols, blockEntry); }
00292
00294
00304 void getLocalBlockEntryView(LocalOrdinal localBlockRow,
00305 LocalOrdinal localBlockCol,
00306 LocalOrdinal& numPtRows,
00307 LocalOrdinal& numPtCols,
00308 Teuchos::ArrayRCP<const Scalar>& blockEntry) const { mtx_->getLocalBlockEntryView(localBlockRow, localBlockCol, numPtRows, numPtCols, blockEntry); }
00309
00311
00327 void getLocalBlockEntryViewNonConst(LocalOrdinal localBlockRow,
00328 LocalOrdinal localBlockCol,
00329 LocalOrdinal& numPtRows,
00330 LocalOrdinal& numPtCols,
00331 Teuchos::ArrayRCP<Scalar>& blockEntry) { mtx_->getLocalBlockEntryViewNonConst(localBlockRow, localBlockCol, numPtRows, numPtCols, blockEntry); }
00332
00334
00336
00337 std::string description() const { return mtx_->description(); }
00338
00341 void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const { mtx_->describable(); }
00343
00344 RCP< const Tpetra::VbrMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps> > getTpetra_VbrMatrix() const { return mtx_; }
00345
00346 private:
00347
00348 const RCP< const Tpetra::VbrMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps> > mtx_;
00349
00350 };
00351
00352 }
00353
00354 #define XPETRA_TPETRAVBRMATRIX_SHORT
00355 #endif //XPETRA_VBRMATRIX_DECL_HPP