All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Defines
Xpetra_TpetraCrsMatrix.hpp
Go to the documentation of this file.
00001 // @HEADER
00002 //
00003 // ***********************************************************************
00004 //
00005 //             Xpetra: A linear algebra interface package
00006 //                  Copyright 2012 Sandia Corporation
00007 //
00008 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00009 // the U.S. Government retains certain rights in this software.
00010 //
00011 // Redistribution and use in source and binary forms, with or without
00012 // modification, are permitted provided that the following conditions are
00013 // met:
00014 //
00015 // 1. Redistributions of source code must retain the above copyright
00016 // notice, this list of conditions and the following disclaimer.
00017 //
00018 // 2. Redistributions in binary form must reproduce the above copyright
00019 // notice, this list of conditions and the following disclaimer in the
00020 // documentation and/or other materials provided with the distribution.
00021 //
00022 // 3. Neither the name of the Corporation nor the names of the
00023 // contributors may be used to endorse or promote products derived from
00024 // this software without specific prior written permission.
00025 //
00026 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00027 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00028 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00029 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00030 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00031 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00032 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00033 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00034 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00035 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00036 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00037 //
00038 // Questions? Contact
00039 //                    Jonathan Hu       (jhu@sandia.gov)
00040 //                    Andrey Prokopenko (aprokop@sandia.gov)
00041 //                    Ray Tuminaro      (rstumin@sandia.gov)
00042 //
00043 // ***********************************************************************
00044 //
00045 // @HEADER
00046 #ifndef XPETRA_TPETRACRSMATRIX_HPP
00047 #define XPETRA_TPETRACRSMATRIX_HPP
00048 
00049 /* this file is automatically generated - do not edit (see scripts/tpetra.py) */
00050 
00051 // FIXME (mfh 03 Sep 2014) The above is probably not true anymore.
00052 // Furthermore, I don't think anyone maintains the scripts.
00053 // Feel free to correct this comment if I'm wrong.
00054 
00055 #include "Xpetra_TpetraConfigDefs.hpp"
00056 
00057 #include "Tpetra_CrsMatrix.hpp"
00058 
00059 #include "Xpetra_CrsMatrix.hpp"
00060 #include "Xpetra_TpetraMap.hpp"
00061 #include "Xpetra_TpetraMultiVector.hpp"
00062 #include "Xpetra_TpetraVector.hpp"
00063 #include "Xpetra_TpetraCrsGraph.hpp"
00064 //#include "Xpetra_TpetraRowMatrix.hpp"
00065 #include "Xpetra_Exceptions.hpp"
00066 
00067 namespace Xpetra {
00068 
00069   // TODO: move that elsewhere
00070   // template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
00071   // const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> toTpetraCrsMatrix(const Xpetra::DistObject<char, LocalOrdinal, GlobalOrdinal, Node> &);
00072   //
00073 
00074   template<class Scalar = CrsMatrix<>::scalar_type,
00075            class LocalOrdinal =
00076              typename CrsMatrix<Scalar>::local_ordinal_type,
00077            class GlobalOrdinal =
00078              typename CrsMatrix<Scalar, LocalOrdinal>::global_ordinal_type,
00079            class Node =
00080              typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
00081   class TpetraCrsMatrix
00082     : public CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> //, public TpetraRowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
00083   {
00084 
00085     // The following typedef are used by the XPETRA_DYNAMIC_CAST() macro.
00086     typedef TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> TpetraCrsMatrixClass;
00087     typedef TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TpetraVectorClass;
00088     typedef TpetraImport<LocalOrdinal,GlobalOrdinal,Node> TpetraImportClass;
00089     typedef TpetraExport<LocalOrdinal,GlobalOrdinal,Node> TpetraExportClass;
00090 
00091   public:
00092 
00094 
00095 
00097     TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
00098       : mtx_ (Teuchos::rcp (new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> (toTpetra(rowMap), maxNumEntriesPerRow, toTpetra(pftype), params))) {  }
00099 
00101     TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
00102       : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> (toTpetra(rowMap), NumEntriesPerRowToAlloc, toTpetra(pftype), params))) {  }
00103 
00105     TpetraCrsMatrix(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 > &params=Teuchos::null)
00106       : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), maxNumEntriesPerRow, toTpetra(pftype), params))) {  }
00107 
00109     TpetraCrsMatrix(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 > &params=Teuchos::null)
00110       : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(rowMap), toTpetra(colMap), NumEntriesPerRowToAlloc, toTpetra(pftype), params))) {  }
00111 
00113     TpetraCrsMatrix(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
00114       : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(graph), params))) {  }
00115 
00116 
00117 
00119     TpetraCrsMatrix(const Teuchos::RCP<const CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& sourceMatrix,
00120                     const Import<LocalOrdinal,GlobalOrdinal,Node> & importer,
00121                     const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
00122                     const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
00123                     const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null)
00124     {
00125       typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
00126       XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
00127       RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
00128 
00129       RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
00130       RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap  = rangeMap!=Teuchos::null  ? toTpetra(rangeMap)  : Teuchos::null;
00131       mtx_=Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(importer),myDomainMap,myRangeMap,params);
00132       bool restrictComm=false;
00133       if(!params.is_null()) restrictComm = params->get("Restrict Communicator",restrictComm);
00134       if(restrictComm && mtx_->getRowMap().is_null()) mtx_=Teuchos::null;
00135 
00136     }
00137 
00139     TpetraCrsMatrix(const Teuchos::RCP<const CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& sourceMatrix,
00140                     const Export<LocalOrdinal,GlobalOrdinal,Node> & exporter,
00141                     const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
00142                     const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
00143                     const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null)
00144     {
00145       typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
00146       XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
00147       RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
00148 
00149       RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
00150       RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap  = rangeMap!=Teuchos::null  ? toTpetra(rangeMap)  : Teuchos::null;
00151       mtx_=Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(exporter),myDomainMap,myRangeMap,params);
00152 
00153     }
00154 
00156     virtual ~TpetraCrsMatrix() {  }
00157 
00159 
00161 
00162 
00164     void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals) { XPETRA_MONITOR("TpetraCrsMatrix::insertGlobalValues"); mtx_->insertGlobalValues(globalRow, cols, vals); }
00165 
00167     void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals) { XPETRA_MONITOR("TpetraCrsMatrix::insertLocalValues"); mtx_->insertLocalValues(localRow, cols, vals); }
00168 
00170     void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals) { XPETRA_MONITOR("TpetraCrsMatrix::replaceGlobalValues"); mtx_->replaceGlobalValues(globalRow, cols, vals); }
00171 
00173     void
00174     replaceLocalValues (LocalOrdinal localRow,
00175                         const ArrayView<const LocalOrdinal> &cols,
00176                         const ArrayView<const Scalar> &vals)
00177     {
00178       XPETRA_MONITOR("TpetraCrsMatrix::replaceLocalValues");
00179       typedef typename ArrayView<const LocalOrdinal>::size_type size_type;
00180       const LocalOrdinal numValid =
00181         mtx_->replaceLocalValues (localRow, cols, vals);
00182       TEUCHOS_TEST_FOR_EXCEPTION(
00183         static_cast<size_type> (numValid) != cols.size (), std::runtime_error,
00184         "replaceLocalValues returned " << numValid << " != cols.size() = " <<
00185         cols.size () << ".");
00186     }
00187 
00189     void setAllToScalar(const Scalar &alpha) { XPETRA_MONITOR("TpetraCrsMatrix::setAllToScalar"); mtx_->setAllToScalar(alpha); }
00190 
00192     void scale(const Scalar &alpha) { XPETRA_MONITOR("TpetraCrsMatrix::scale"); mtx_->scale(alpha); }
00193 
00195     //** \warning This is an expert-only routine and should not be called from user code. */
00196     void allocateAllValues(size_t numNonZeros,ArrayRCP<size_t> & rowptr, ArrayRCP<LocalOrdinal> & colind, ArrayRCP<Scalar> & values)
00197     { XPETRA_MONITOR("TpetraCrsMatrix::allocateAllValues"); rowptr.resize(getNodeNumRows()+1); colind.resize(numNonZeros); values.resize(numNonZeros);}
00198 
00200     void setAllValues(const ArrayRCP<size_t> & rowptr, const ArrayRCP<LocalOrdinal> & colind, const ArrayRCP<Scalar> & values)
00201     { XPETRA_MONITOR("TpetraCrsMatrix::setAllValues"); mtx_->setAllValues(rowptr,colind,values); }
00202 
00204     void getAllValues(ArrayRCP<const size_t>& rowptr, ArrayRCP<const LocalOrdinal>& colind, ArrayRCP<const Scalar>& values) const
00205     { XPETRA_MONITOR("TpetraCrsMatrix::getAllValues"); mtx_->getAllValues(rowptr,colind,values); }
00206 
00208 
00210 
00211 
00213     void resumeFill(const RCP< ParameterList > &params=null) { XPETRA_MONITOR("TpetraCrsMatrix::resumeFill"); mtx_->resumeFill(params); }
00214 
00216     void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null) { XPETRA_MONITOR("TpetraCrsMatrix::fillComplete"); mtx_->fillComplete(toTpetra(domainMap), toTpetra(rangeMap), params); }
00217 
00219     void fillComplete(const RCP< ParameterList > &params=null) { XPETRA_MONITOR("TpetraCrsMatrix::fillComplete"); mtx_->fillComplete(params); }
00220 
00221 
00223     void replaceDomainMapAndImporter(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >& newDomainMap, Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> >  & newImporter) {
00224       XPETRA_MONITOR("TpetraCrsMatrix::replaceDomainMapAndImporter");
00225       XPETRA_DYNAMIC_CAST( const TpetraImportClass , *newImporter, tImporter, "Xpetra::TpetraCrsMatrix::replaceDomainMapAndImporter only accepts Xpetra::TpetraImport.");
00226       RCP<const Tpetra::Import<LocalOrdinal,GlobalOrdinal,Node> > myImport = tImporter.getTpetra_Import();
00227             mtx_->replaceDomainMapAndImporter( toTpetra(newDomainMap),myImport);
00228     }
00229 
00231     void expertStaticFillComplete(const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & domainMap,
00232                                   const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & rangeMap,
00233                                   const RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > &importer=Teuchos::null,
00234                                   const RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > &exporter=Teuchos::null,
00235                                   const RCP<ParameterList> &params=Teuchos::null) {
00236       XPETRA_MONITOR("TpetraCrsMatrix::expertStaticFillComplete");
00237       RCP<const Tpetra::Import<LocalOrdinal,GlobalOrdinal,Node> > myImport;
00238       RCP<const Tpetra::Export<LocalOrdinal,GlobalOrdinal,Node> > myExport;
00239 
00240       if(importer!=Teuchos::null) {
00241         XPETRA_DYNAMIC_CAST( const TpetraImportClass , *importer, tImporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraImport.");
00242         myImport = tImporter.getTpetra_Import();
00243       }
00244       if(exporter!=Teuchos::null) {
00245         XPETRA_DYNAMIC_CAST( const TpetraExportClass , *exporter, tExporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraExport.");
00246         myExport = tExporter.getTpetra_Export();
00247       }
00248 
00249       mtx_->expertStaticFillComplete(toTpetra(domainMap),toTpetra(rangeMap),myImport,myExport,params);
00250     }
00251 
00253 
00255 
00256 
00258     const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  getRowMap() const { XPETRA_MONITOR("TpetraCrsMatrix::getRowMap"); return toXpetra(mtx_->getRowMap()); }
00259 
00261     const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  getColMap() const { XPETRA_MONITOR("TpetraCrsMatrix::getColMap"); return toXpetra(mtx_->getColMap()); }
00262 
00264     RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const { XPETRA_MONITOR("TpetraCrsMatrix::getCrsGraph"); return toXpetra(mtx_->getCrsGraph()); }
00265 
00267     global_size_t getGlobalNumRows() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumRows"); return mtx_->getGlobalNumRows(); }
00268 
00270     global_size_t getGlobalNumCols() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumCols"); return mtx_->getGlobalNumCols(); }
00271 
00273     size_t getNodeNumRows() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeNumRows"); return mtx_->getNodeNumRows(); }
00274 
00276     size_t getNodeNumCols() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeNumCols"); return mtx_->getNodeNumCols(); }
00277 
00279     global_size_t getGlobalNumEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumEntries"); return mtx_->getGlobalNumEntries(); }
00280 
00282     size_t getNodeNumEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeNumEntries"); return mtx_->getNodeNumEntries(); }
00283 
00285     size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const { XPETRA_MONITOR("TpetraCrsMatrix::getNumEntriesInLocalRow"); return mtx_->getNumEntriesInLocalRow(localRow); }
00286 
00288     global_size_t getGlobalNumDiags() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumDiags"); return mtx_->getGlobalNumDiags(); }
00289 
00291     size_t getNodeNumDiags() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeNumDiags"); return mtx_->getNodeNumDiags(); }
00292 
00294     size_t getGlobalMaxNumRowEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalMaxNumRowEntries"); return mtx_->getGlobalMaxNumRowEntries(); }
00295 
00297     size_t getNodeMaxNumRowEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeMaxNumRowEntries"); return mtx_->getNodeMaxNumRowEntries(); }
00298 
00300     bool isLocallyIndexed() const { XPETRA_MONITOR("TpetraCrsMatrix::isLocallyIndexed"); return mtx_->isLocallyIndexed(); }
00301 
00303     bool isGloballyIndexed() const { XPETRA_MONITOR("TpetraCrsMatrix::isGloballyIndexed"); return mtx_->isGloballyIndexed(); }
00304 
00306     bool isFillComplete() const { XPETRA_MONITOR("TpetraCrsMatrix::isFillComplete"); return mtx_->isFillComplete(); }
00307 
00309     bool isFillActive() const { XPETRA_MONITOR("TpetraCrsMatrix::isFillActive"); return mtx_->isFillActive(); }
00310 
00312     typename ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const { XPETRA_MONITOR("TpetraCrsMatrix::getFrobeniusNorm"); return mtx_->getFrobeniusNorm(); }
00313 
00315     bool supportsRowViews() const { XPETRA_MONITOR("TpetraCrsMatrix::supportsRowViews"); return mtx_->supportsRowViews(); }
00316 
00318     void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const { XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowCopy"); mtx_->getLocalRowCopy(LocalRow, Indices, Values, NumEntries); }
00319 
00321     void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowView"); mtx_->getGlobalRowView(GlobalRow, indices, values); }
00322 
00324     void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowCopy"); mtx_->getGlobalRowCopy(GlobalRow, indices, values, numEntries); }
00325 
00327     void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const { XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowView"); mtx_->getLocalRowView(LocalRow, indices, values); }
00328 
00330 
00332 
00333 
00335     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("TpetraCrsMatrix::apply"); mtx_->apply(toTpetra(X), toTpetra(Y), mode, alpha, beta); }
00336 
00338     const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  getDomainMap() const { XPETRA_MONITOR("TpetraCrsMatrix::getDomainMap"); return toXpetra(mtx_->getDomainMap()); }
00339 
00341     const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  getRangeMap() const { XPETRA_MONITOR("TpetraCrsMatrix::getRangeMap"); return toXpetra(mtx_->getRangeMap()); }
00342 
00344 
00346 
00347 
00349     std::string description() const { XPETRA_MONITOR("TpetraCrsMatrix::description"); return mtx_->description(); }
00350 
00352     void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const { XPETRA_MONITOR("TpetraCrsMatrix::describe"); mtx_->describe(out, verbLevel); }
00353 
00355 
00357     TpetraCrsMatrix(const TpetraCrsMatrix& matrix)
00358       : mtx_ (matrix.mtx_->template clone<Node> (matrix.mtx_->getNode ())) {}
00359 
00361     void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const {
00362       XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
00363       XPETRA_DYNAMIC_CAST(TpetraVectorClass, diag, tDiag, "Xpetra::TpetraCrsMatrix.getLocalDiagCopy() only accept Xpetra::TpetraVector as input arguments.");
00364       mtx_->getLocalDiagCopy(*tDiag.getTpetra_Vector());
00365       // mtx_->getLocalDiagCopy(toTpetra(diag));
00366     }
00367 
00369     void getLocalDiagOffsets(Teuchos::ArrayRCP<size_t> &offsets) const {
00370       XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagOffsets");
00371       mtx_->getLocalDiagOffsets(offsets);
00372     }
00373 
00375     void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag, const Teuchos::ArrayView<const size_t> &offsets) const {
00376       XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
00377       //XPETRA_DYNAMIC_CAST(TpetraVectorClass, diag, tDiag, "Xpetra::TpetraCrsMatrix.getLocalDiagCopy() only accept Xpetra::TpetraVector as input arguments.");
00378       //mtx_->getLocalDiagCopy(*tDiag.getTpetra_Vector(), offsets);
00379       mtx_->getLocalDiagCopy(*(toTpetra(diag)), offsets);
00380     }
00381 
00383     //{@
00384 
00386     Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const { XPETRA_MONITOR("TpetraCrsMatrix::getMap"); return rcp( new TpetraMap< LocalOrdinal, GlobalOrdinal, Node >(mtx_->getMap()) ); }
00387 
00389     void doImport(const DistObject<char, LocalOrdinal, GlobalOrdinal, Node> &source,
00390                   const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM) {
00391       XPETRA_MONITOR("TpetraCrsMatrix::doImport");
00392 
00393       XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
00394       RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSource.getTpetra_CrsMatrix();
00395       //mtx_->doImport(toTpetraCrsMatrix(source), *tImporter.getTpetra_Import(), toTpetra(CM));
00396       mtx_->doImport(*v, toTpetra(importer), toTpetra(CM));
00397     }
00398 
00400     void doExport(const DistObject<char, LocalOrdinal, GlobalOrdinal, Node> &dest,
00401                   const Import< LocalOrdinal, GlobalOrdinal, Node >& importer, CombineMode CM) {
00402       XPETRA_MONITOR("TpetraCrsMatrix::doExport");
00403 
00404       XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
00405       RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tDest.getTpetra_CrsMatrix();
00406       mtx_->doExport(*v, toTpetra(importer), toTpetra(CM));
00407 
00408     }
00409 
00411     void doImport(const DistObject<char, LocalOrdinal, GlobalOrdinal, Node> &source,
00412                   const Export< LocalOrdinal, GlobalOrdinal, Node >& exporter, CombineMode CM) {
00413       XPETRA_MONITOR("TpetraCrsMatrix::doImport");
00414 
00415       XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
00416       RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSource.getTpetra_CrsMatrix();
00417       mtx_->doImport(*v, toTpetra(exporter), toTpetra(CM));
00418 
00419     }
00420 
00422     void doExport(const DistObject<char, LocalOrdinal, GlobalOrdinal, Node> &dest,
00423                   const Export< LocalOrdinal, GlobalOrdinal, Node >& exporter, CombineMode CM) {
00424       XPETRA_MONITOR("TpetraCrsMatrix::doExport");
00425 
00426       XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
00427       RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tDest.getTpetra_CrsMatrix();
00428       mtx_->doExport(*v, toTpetra(exporter), toTpetra(CM));
00429 
00430     }
00431 
00432     void removeEmptyProcessesInPlace (const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >& newMap) {
00433       XPETRA_MONITOR("TpetraCrsMatrix::removeEmptyProcessesInPlace");
00434       mtx_->removeEmptyProcessesInPlace(toTpetra(newMap));
00435     }
00436 
00437     // @}
00438 
00439     template<class Node2>
00440     RCP<CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node2> >
00441     clone (const RCP<Node2> &node2) const
00442     {
00443       return RCP<CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node2> > (new TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node2> (mtx_->clone (node2)));
00444     }
00445 
00447 
00448 
00450     bool hasMatrix() const {
00451       return ! mtx_.is_null ();
00452     }
00453 
00455     TpetraCrsMatrix(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &mtx) : mtx_(mtx) {  }
00456 
00458     RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > getTpetra_CrsMatrix() const { return mtx_; }
00459 
00461     RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > getTpetra_CrsMatrixNonConst() const { return mtx_; } //TODO: remove
00462 
00464 
00465   private:
00466     RCP< Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > mtx_;
00467   }; // TpetraCrsMatrix class
00468 
00469 } // Xpetra namespace
00470 
00471 #define XPETRA_TPETRACRSMATRIX_SHORT
00472 #endif // XPETRA_TPETRACRSMATRIX_HPP
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Defines