Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Xpetra_TpetraCrsMatrix_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Xpetra: A linear algebra interface package
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef XPETRA_TPETRACRSMATRIX_DEF_HPP
47 #define XPETRA_TPETRACRSMATRIX_DEF_HPP
48 
50 #include "Tpetra_Details_residual.hpp"
51 
52 namespace Xpetra {
53 
54  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
55  TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, const Teuchos::RCP< Teuchos::ParameterList > &params)
56  : mtx_ (Teuchos::rcp (new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> (toTpetra(rowMap), maxNumEntriesPerRow, Tpetra::StaticProfile, params))) { }
57 
58  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
59  TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, const Teuchos::RCP< Teuchos::ParameterList > &params)
60  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> (toTpetra(rowMap), NumEntriesPerRowToAlloc(), Tpetra::StaticProfile, params))) { }
61 
62  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
63  TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, const Teuchos::RCP< Teuchos::ParameterList > &params)
64  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), maxNumEntriesPerRow, Tpetra::StaticProfile, params))) { }
65 
66  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
67  TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::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, const Teuchos::RCP< Teuchos::ParameterList > &params)
68  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(rowMap), toTpetra(colMap), NumEntriesPerRowToAlloc(), Tpetra::StaticProfile, params))) { }
69 
70  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
71  TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::TpetraCrsMatrix(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph, const Teuchos::RCP< Teuchos::ParameterList > &params)
72  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(graph), params))) { }
73 
74 
75 
76  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
79  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
80  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
81  const Teuchos::RCP<Teuchos::ParameterList>& params)
82  {
83  typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
84  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
85  RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
86 
87  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
88  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
89  mtx_=Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(importer),myDomainMap,myRangeMap,params);
90  bool restrictComm=false;
91  if(!params.is_null()) restrictComm = params->get("Restrict Communicator",restrictComm);
92  if(restrictComm && mtx_->getRowMap().is_null()) mtx_=Teuchos::null;
93 
94  }
95 
96  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
99  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
100  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
101  const Teuchos::RCP<Teuchos::ParameterList>& params)
102  {
103  typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
104  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
105  RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
106 
107  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
108  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
109  mtx_=Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(exporter),myDomainMap,myRangeMap,params);
110 
111  }
112 
113  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
115  const Import<LocalOrdinal,GlobalOrdinal,Node> & RowImporter,
116  const Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > DomainImporter,
117  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
118  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
119  const Teuchos::RCP<Teuchos::ParameterList>& params)
120  {
121  typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
122  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
123  RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
124 
125  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
126  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
127 
128  mtx_=Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(RowImporter),toTpetra(*DomainImporter),myDomainMap,myRangeMap,params);
129  bool restrictComm=false;
130  if(!params.is_null()) restrictComm = params->get("Restrict Communicator",restrictComm);
131  if(restrictComm && mtx_->getRowMap().is_null()) mtx_=Teuchos::null;
132  }
133 
134  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
136  const Export<LocalOrdinal,GlobalOrdinal,Node> & RowExporter,
137  const Teuchos::RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > DomainExporter,
138  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
139  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
140  const Teuchos::RCP<Teuchos::ParameterList>& params)
141  {
142  typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
143  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
144  RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
145 
146  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
147  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
148 
149  mtx_=Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(RowExporter),toTpetra(*DomainExporter),myDomainMap,myRangeMap,params);
150  }
151 
153 
154 
155 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
156 #ifdef HAVE_XPETRA_TPETRA
157  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
159  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& colMap,
160  const local_matrix_type& lclMatrix,
161  const Teuchos::RCP<Teuchos::ParameterList>& params)
162  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), lclMatrix, params))) { }
163 
164  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
166  const local_matrix_type& lclMatrix,
167  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rowMap,
168  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& colMap,
169  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
170  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
171  const Teuchos::RCP<Teuchos::ParameterList>& params)
172  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(lclMatrix, toTpetra(rowMap), toTpetra(colMap), toTpetra(domainMap), toTpetra(rangeMap), params))) { }
173 #endif
174 #endif
175 
176  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
178 
179  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
180  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals) { XPETRA_MONITOR("TpetraCrsMatrix::insertGlobalValues"); mtx_->insertGlobalValues(globalRow, cols, vals); }
181 
182  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
183  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals) { XPETRA_MONITOR("TpetraCrsMatrix::insertLocalValues"); mtx_->insertLocalValues(localRow, cols, vals); }
184 
185  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
186  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals) { XPETRA_MONITOR("TpetraCrsMatrix::replaceGlobalValues"); mtx_->replaceGlobalValues(globalRow, cols, vals); }
187 
188  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
189  void
191  const ArrayView<const LocalOrdinal> &cols,
192  const ArrayView<const Scalar> &vals)
193  {
194  XPETRA_MONITOR("TpetraCrsMatrix::replaceLocalValues");
195  typedef typename ArrayView<const LocalOrdinal>::size_type size_type;
196  const LocalOrdinal numValid =
197  mtx_->replaceLocalValues (localRow, cols, vals);
198  TEUCHOS_TEST_FOR_EXCEPTION(
199  static_cast<size_type> (numValid) != cols.size (), std::runtime_error,
200  "replaceLocalValues returned " << numValid << " != cols.size() = " <<
201  cols.size () << ".");
202  }
203 
204  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
205  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::setAllToScalar(const Scalar &alpha) { XPETRA_MONITOR("TpetraCrsMatrix::setAllToScalar"); mtx_->setAllToScalar(alpha); }
206 
207  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
208  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::scale(const Scalar &alpha) { XPETRA_MONITOR("TpetraCrsMatrix::scale"); mtx_->scale(alpha); }
209 
210  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
211  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::allocateAllValues(size_t numNonZeros,ArrayRCP<size_t> & rowptr, ArrayRCP<LocalOrdinal> & colind, ArrayRCP<Scalar> & values)
212  { XPETRA_MONITOR("TpetraCrsMatrix::allocateAllValues"); rowptr.resize(getNodeNumRows()+1); colind.resize(numNonZeros); values.resize(numNonZeros);}
213 
214  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
215  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::setAllValues(const ArrayRCP<size_t> & rowptr, const ArrayRCP<LocalOrdinal> & colind, const ArrayRCP<Scalar> & values)
216  { XPETRA_MONITOR("TpetraCrsMatrix::setAllValues"); mtx_->setAllValues(rowptr,colind,values); }
217 
218  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
219  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getAllValues(ArrayRCP<const size_t>& rowptr, ArrayRCP<const LocalOrdinal>& colind, ArrayRCP<const Scalar>& values) const
220  { XPETRA_MONITOR("TpetraCrsMatrix::getAllValues"); mtx_->getAllValues(rowptr,colind,values); }
221 
222  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
224  { return mtx_->haveGlobalConstants();}
225 
226  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
227  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::resumeFill(const RCP< ParameterList > &params) { XPETRA_MONITOR("TpetraCrsMatrix::resumeFill"); mtx_->resumeFill(params); }
228 
229  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
230  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params) { XPETRA_MONITOR("TpetraCrsMatrix::fillComplete"); mtx_->fillComplete(toTpetra(domainMap), toTpetra(rangeMap), params); }
231 
232  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
233  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::fillComplete(const RCP< ParameterList > &params) { XPETRA_MONITOR("TpetraCrsMatrix::fillComplete"); mtx_->fillComplete(params); }
234 
235 
236  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
238  XPETRA_MONITOR("TpetraCrsMatrix::replaceDomainMapAndImporter");
239  XPETRA_DYNAMIC_CAST( const TpetraImportClass , *newImporter, tImporter, "Xpetra::TpetraCrsMatrix::replaceDomainMapAndImporter only accepts Xpetra::TpetraImport.");
240  RCP<const Tpetra::Import<LocalOrdinal,GlobalOrdinal,Node> > myImport = tImporter.getTpetra_Import();
241  mtx_->replaceDomainMapAndImporter( toTpetra(newDomainMap),myImport);
242  }
243 
244  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
246  const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & rangeMap,
247  const RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > &importer,
248  const RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > &exporter,
249  const RCP<ParameterList> &params) {
250  XPETRA_MONITOR("TpetraCrsMatrix::expertStaticFillComplete");
251  RCP<const Tpetra::Import<LocalOrdinal,GlobalOrdinal,Node> > myImport;
252  RCP<const Tpetra::Export<LocalOrdinal,GlobalOrdinal,Node> > myExport;
253 
254  if(importer!=Teuchos::null) {
255  XPETRA_DYNAMIC_CAST( const TpetraImportClass , *importer, tImporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraImport.");
256  myImport = tImporter.getTpetra_Import();
257  }
258  if(exporter!=Teuchos::null) {
259  XPETRA_DYNAMIC_CAST( const TpetraExportClass , *exporter, tExporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraExport.");
260  myExport = tExporter.getTpetra_Export();
261  }
262 
263  mtx_->expertStaticFillComplete(toTpetra(domainMap),toTpetra(rangeMap),myImport,myExport,params);
264  }
265 
266  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
267  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getRowMap() const { XPETRA_MONITOR("TpetraCrsMatrix::getRowMap"); return toXpetra(mtx_->getRowMap()); }
268 
269  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
270  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getColMap() const { XPETRA_MONITOR("TpetraCrsMatrix::getColMap"); return toXpetra(mtx_->getColMap()); }
271 
272  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
273  RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getCrsGraph() const { XPETRA_MONITOR("TpetraCrsMatrix::getCrsGraph"); return toXpetra(mtx_->getCrsGraph()); }
274 
275  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
276  global_size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getGlobalNumRows() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumRows"); return mtx_->getGlobalNumRows(); }
277 
278  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
279  global_size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getGlobalNumCols() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumCols"); return mtx_->getGlobalNumCols(); }
280 
281  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
282  size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNodeNumRows() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeNumRows"); return mtx_->getNodeNumRows(); }
283 
284  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
285  size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNodeNumCols() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeNumCols"); return mtx_->getNodeNumCols(); }
286 
287  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
288  global_size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getGlobalNumEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumEntries"); return mtx_->getGlobalNumEntries(); }
289 
290  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
291  size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNodeNumEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeNumEntries"); return mtx_->getNodeNumEntries(); }
292 
293  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
294  size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNumEntriesInLocalRow(LocalOrdinal localRow) const { XPETRA_MONITOR("TpetraCrsMatrix::getNumEntriesInLocalRow"); return mtx_->getNumEntriesInLocalRow(localRow); }
295 
296  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
297  size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const { XPETRA_MONITOR("TpetraCrsMatrix::getNumEntriesInGlobalRow"); return mtx_->getNumEntriesInGlobalRow(globalRow); }
298 
299  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
300  size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getGlobalMaxNumRowEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalMaxNumRowEntries"); return mtx_->getGlobalMaxNumRowEntries(); }
301 
302  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
303  size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNodeMaxNumRowEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeMaxNumRowEntries"); return mtx_->getNodeMaxNumRowEntries(); }
304 
305  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
306  bool TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isLocallyIndexed() const { XPETRA_MONITOR("TpetraCrsMatrix::isLocallyIndexed"); return mtx_->isLocallyIndexed(); }
307 
308  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
309  bool TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isGloballyIndexed() const { XPETRA_MONITOR("TpetraCrsMatrix::isGloballyIndexed"); return mtx_->isGloballyIndexed(); }
310 
311  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
312  bool TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isFillComplete() const { XPETRA_MONITOR("TpetraCrsMatrix::isFillComplete"); return mtx_->isFillComplete(); }
313 
314  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
315  bool TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isFillActive() const { XPETRA_MONITOR("TpetraCrsMatrix::isFillActive"); return mtx_->isFillActive(); }
316 
317  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
318  typename ScalarTraits< Scalar >::magnitudeType TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getFrobeniusNorm() const { XPETRA_MONITOR("TpetraCrsMatrix::getFrobeniusNorm"); return mtx_->getFrobeniusNorm(); }
319 
320  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
321  bool TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::supportsRowViews() const { XPETRA_MONITOR("TpetraCrsMatrix::supportsRowViews"); return mtx_->supportsRowViews(); }
322 
323  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
324  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::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); }
325 
326  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
327  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowView"); mtx_->getGlobalRowView(GlobalRow, indices, values); }
328 
329  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
330  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::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); }
331 
332  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
333  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const { XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowView"); mtx_->getLocalRowView(LocalRow, indices, values); }
334 
335  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
336  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::apply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta) const { XPETRA_MONITOR("TpetraCrsMatrix::apply"); mtx_->apply(toTpetra(X), toTpetra(Y), mode, alpha, beta); }
337 
338  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
339  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getDomainMap() const { XPETRA_MONITOR("TpetraCrsMatrix::getDomainMap"); return toXpetra(mtx_->getDomainMap()); }
340 
341  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
342  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getRangeMap() const { XPETRA_MONITOR("TpetraCrsMatrix::getRangeMap"); return toXpetra(mtx_->getRangeMap()); }
343 
344  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
345  std::string TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::description() const { XPETRA_MONITOR("TpetraCrsMatrix::description"); return mtx_->description(); }
346 
347  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
348  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const { XPETRA_MONITOR("TpetraCrsMatrix::describe"); mtx_->describe(out, verbLevel); }
349 
350  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
352  XPETRA_MONITOR("TpetraCrsMatrix::setObjectLabel");
353  Teuchos::LabeledObject::setObjectLabel(objectLabel);
354  mtx_->setObjectLabel(objectLabel);
355  }
356 
357 
358 
359  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
361  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>(*(matrix.mtx_),Teuchos::Copy))) {}
362 
363  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
365  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
366  XPETRA_DYNAMIC_CAST(TpetraVectorClass, diag, tDiag, "Xpetra::TpetraCrsMatrix.getLocalDiagCopy() only accept Xpetra::TpetraVector as input arguments.");
367  mtx_->getLocalDiagCopy(*tDiag.getTpetra_Vector());
368  }
369 
370  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
372  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagOffsets");
373  mtx_->getLocalDiagOffsets(offsets);
374  }
375 
376  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
378  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
379  mtx_->getLocalDiagCopy(*(toTpetra(diag)), offsets);
380  }
381 
382  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
384  XPETRA_MONITOR("TpetraCrsMatrix::replaceDiag");
385  Tpetra::replaceDiagonalCrsMatrix(*mtx_, *(toTpetra(diag)));
386  }
387 
388  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
390  XPETRA_MONITOR("TpetraCrsMatrix::leftScale");
391  mtx_->leftScale(*(toTpetra(x)));
392  }
393 
394  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
396  XPETRA_MONITOR("TpetraCrsMatrix::rightScale");
397  mtx_->rightScale(*(toTpetra(x)));
398  }
399 
400  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
401  Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getMap() const { XPETRA_MONITOR("TpetraCrsMatrix::getMap"); return rcp( new TpetraMap< LocalOrdinal, GlobalOrdinal, Node >(mtx_->getMap()) ); }
402 
403  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
406  XPETRA_MONITOR("TpetraCrsMatrix::doImport");
407 
408  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
409  RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSource.getTpetra_CrsMatrix();
410  //mtx_->doImport(toTpetraCrsMatrix(source), *tImporter.getTpetra_Import(), toTpetra(CM));
411  mtx_->doImport(*v, toTpetra(importer), toTpetra(CM));
412  }
413 
415  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
418  XPETRA_MONITOR("TpetraCrsMatrix::doExport");
419 
420  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
421  RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tDest.getTpetra_CrsMatrix();
422  mtx_->doExport(*v, toTpetra(importer), toTpetra(CM));
423 
424  }
425 
426  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
429  XPETRA_MONITOR("TpetraCrsMatrix::doImport");
430 
431  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
432  RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSource.getTpetra_CrsMatrix();
433  mtx_->doImport(*v, toTpetra(exporter), toTpetra(CM));
434 
435  }
436 
437  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
440  XPETRA_MONITOR("TpetraCrsMatrix::doExport");
441 
442  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
443  RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tDest.getTpetra_CrsMatrix();
444  mtx_->doExport(*v, toTpetra(exporter), toTpetra(CM));
445 
446  }
447 
448  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
450  XPETRA_MONITOR("TpetraCrsMatrix::removeEmptyProcessesInPlace");
451  mtx_->removeEmptyProcessesInPlace(toTpetra(newMap));
452  }
453 
454  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
456  return ! mtx_.is_null ();
457  }
458 
459  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
460  TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::TpetraCrsMatrix(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &mtx) : mtx_(mtx) { }
461 
462  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
463  RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getTpetra_CrsMatrix() const { return mtx_; }
464 
466  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
467  RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getTpetra_CrsMatrixNonConst() const { return mtx_; } //TODO: remove
468 
470  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
474  Tpetra::Details::residual(*mtx_,toTpetra(X),toTpetra(B),toTpetra(R));
475  }
476 
477 
480 // End of TpetrCrsMatrix class definition //
483 
484 } // Xpetra namespace
485 
486 #endif // XPETRA_TPETRACRSMATRIX_DEF_HPP
global_size_t getGlobalNumRows() const
Number of global elements in the row map of this matrix.
std::string description() const
A simple one-line description of this object.
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalarThis.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this operator. This will be null until fillComplete() i...
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Get a copy of the diagonal entries owned by this node, with local row indices.
void residual(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &R) const
Compute a residual R = B - (*this) * X.
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries in the (locally owned) global row.
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
Computes the sparse matrix-multivector multiplication.
void setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)
Sets the 1D pointer arrays of the graph.
void allocateAllValues(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind, ArrayRCP< Scalar > &values)
Allocates and returns ArrayRCPs of the Crs arrays — This is an Xpetra-only routine.
void resumeFill(const RCP< ParameterList > &params=null)
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
void setObjectLabel(const std::string &objectLabel)
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this matrix. ...
void getAllValues(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind, ArrayRCP< const Scalar > &values) const
Gets the 1D pointer arrays of the graph.
void replaceDomainMapAndImporter(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newDomainMap, Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &newImporter)
Replaces the current domainMap and importer with the user-specified objects.
bool isFillComplete() const
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called.
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calli...
void expertStaticFillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &importer=Teuchos::null, const RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > &exporter=Teuchos::null, const RCP< ParameterList > &params=Teuchos::null)
Expert static fill complete.
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of global indices in a specified row of the matrix.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale operator with given vector values.
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
size_t global_size_t
Global size_t object.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying fixed number of entries for each row.
virtual ~TpetraCrsMatrix()
Destructor.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Map associated with the range of this operator, which must be compatible with Y...
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrix() const
Get the underlying Tpetra matrix.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
global_size_t getGlobalNumCols() const
Number of global columns in the matrix.
bool supportsRowViews() const
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class.
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrixNonConst() const
Get the underlying Tpetra matrix.
bool hasMatrix() const
Does this have an underlying matrix.
void replaceDiag(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag)
Replace the diagonal entries of the matrix.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
CombineMode
Xpetra::Combine Mode enumerable type.
#define XPETRA_MONITOR(funcName)
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale operator with given vector values.
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of local indices in a specified row of the matrix.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this matrix.
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null)
Signal that data entry is complete, specifying domain and range maps.
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const
Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage...
size_t getNodeMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
bool isFillActive() const
Returns true if the matrix is in edit mode.
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.