Xpetra Version of the Day
Loading...
Searching...
No Matches
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
49#include <Xpetra_MultiVectorFactory.hpp>
51#include "Tpetra_Details_residual.hpp"
52
53namespace Xpetra {
54
55 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
57 : mtx_ (Teuchos::rcp (new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> (toTpetra(rowMap), maxNumEntriesPerRow, params))) { }
58
59 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
61 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> (toTpetra(rowMap), NumEntriesPerRowToAlloc(), params))) { }
62
63 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
65 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), maxNumEntriesPerRow, params))) { }
66
67 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
69 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(rowMap), toTpetra(colMap), NumEntriesPerRowToAlloc(), params))) { }
70
71 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
74
75 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
78
79
80
81
82 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
88 {
89 typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
90 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
92
93 RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
94 RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
95 mtx_=Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(importer),myDomainMap,myRangeMap,params);
96 bool restrictComm=false;
97 if(!params.is_null()) restrictComm = params->get("Restrict Communicator",restrictComm);
98 if(restrictComm && mtx_->getRowMap().is_null()) mtx_=Teuchos::null;
99
100 }
101
102 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
105 const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
108 {
109 typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
110 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
111 RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
112
113 RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
114 RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
115 mtx_=Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(exporter),myDomainMap,myRangeMap,params);
116
117 }
118
119 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
121 const Import<LocalOrdinal,GlobalOrdinal,Node> & RowImporter,
122 const Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > DomainImporter,
123 const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
126 {
127 typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
128 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
129 RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
130
131 RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
132 RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
133
134 mtx_=Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(RowImporter),toTpetra(*DomainImporter),myDomainMap,myRangeMap,params);
135 bool restrictComm=false;
136 if(!params.is_null()) restrictComm = params->get("Restrict Communicator",restrictComm);
137 if(restrictComm && mtx_->getRowMap().is_null()) mtx_=Teuchos::null;
138 }
139
140 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
142 const Export<LocalOrdinal,GlobalOrdinal,Node> & RowExporter,
143 const Teuchos::RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > DomainExporter,
144 const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
147 {
148 typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
149 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
150 RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
151
152 RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
153 RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
154
155 mtx_=Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(RowExporter),toTpetra(*DomainExporter),myDomainMap,myRangeMap,params);
156 }
157
159
160
161#ifdef HAVE_XPETRA_TPETRA
162 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
165 const local_matrix_type& lclMatrix,
167 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), lclMatrix, params))) { }
168
169 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
171 const local_matrix_type& lclMatrix,
174 const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
177 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(lclMatrix, toTpetra(rowMap), toTpetra(colMap), toTpetra(domainMap), toTpetra(rangeMap), params))) { }
178
179 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
181 const local_matrix_type& lclMatrix,
184 const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
189 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(lclMatrix, toTpetra(rowMap), toTpetra(colMap), toTpetra(domainMap), toTpetra(rangeMap), toTpetra(importer), toTpetra(exporter), params))) { }
190#endif
191
192 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
194
195 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
196 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); }
197
198 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
199 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); }
200
201 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
202 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); }
203
204 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
205 void
208 const ArrayView<const Scalar> &vals)
209 {
210 XPETRA_MONITOR("TpetraCrsMatrix::replaceLocalValues");
211 typedef typename ArrayView<const LocalOrdinal>::size_type size_type;
212 const LocalOrdinal numValid =
213 mtx_->replaceLocalValues (localRow, cols, vals);
215 static_cast<size_type> (numValid) != cols.size (), std::runtime_error,
216 "replaceLocalValues returned " << numValid << " != cols.size() = " <<
217 cols.size () << ".");
218 }
219
220 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
221 void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::setAllToScalar(const Scalar &alpha) { XPETRA_MONITOR("TpetraCrsMatrix::setAllToScalar"); mtx_->setAllToScalar(alpha); }
222
223 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
224 void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::scale(const Scalar &alpha) { XPETRA_MONITOR("TpetraCrsMatrix::scale"); mtx_->scale(alpha); }
225
226 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
228 { XPETRA_MONITOR("TpetraCrsMatrix::allocateAllValues"); rowptr.resize(getLocalNumRows()+1); colind.resize(numNonZeros); values.resize(numNonZeros);}
229
230 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
232 { XPETRA_MONITOR("TpetraCrsMatrix::setAllValues"); mtx_->setAllValues(rowptr,colind,values); }
233
234 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
236 { XPETRA_MONITOR("TpetraCrsMatrix::getAllValues");
237 // TODO: Change getAllValues interface to return Kokkos views,
238 // TODO: or add getLocalRowPtrsHost, getLocalIndicesHost, etc., interfaces
239 // TODO: and use them in MueLu
240 rowptr = Kokkos::Compat::persistingView(mtx_->getLocalRowPtrsHost());
241 colind = Kokkos::Compat::persistingView(mtx_->getLocalIndicesHost());
242 values = Teuchos::arcp_reinterpret_cast<const Scalar>(
243 Kokkos::Compat::persistingView(mtx_->getLocalValuesHost(
244 Tpetra::Access::ReadOnly)));
245 }
246
247 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
249 { XPETRA_MONITOR("TpetraCrsMatrix::getAllValues");
250 // TODO: Change getAllValues interface to return Kokkos view,
251 // TODO: or add getLocalValuesDevice() interfaces
252 values = Teuchos::arcp_reinterpret_cast<Scalar>(
253 Kokkos::Compat::persistingView(mtx_->getLocalValuesDevice(
254 Tpetra::Access::ReadWrite)));
255 }
256
257 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
260
261 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
262 void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::resumeFill(const RCP< ParameterList > &params) { XPETRA_MONITOR("TpetraCrsMatrix::resumeFill"); mtx_->resumeFill(params); }
263
264 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
265 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); }
266
267 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
268 void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::fillComplete(const RCP< ParameterList > &params) { XPETRA_MONITOR("TpetraCrsMatrix::fillComplete"); mtx_->fillComplete(params); }
269
270
271 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
273 XPETRA_MONITOR("TpetraCrsMatrix::replaceDomainMapAndImporter");
274 XPETRA_DYNAMIC_CAST( const TpetraImportClass , *newImporter, tImporter, "Xpetra::TpetraCrsMatrix::replaceDomainMapAndImporter only accepts Xpetra::TpetraImport.");
275 RCP<const Tpetra::Import<LocalOrdinal,GlobalOrdinal,Node> > myImport = tImporter.getTpetra_Import();
276 mtx_->replaceDomainMapAndImporter( toTpetra(newDomainMap),myImport);
277 }
278
279 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
281 const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & rangeMap,
282 const RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > &importer,
283 const RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > &exporter,
284 const RCP<ParameterList> &params) {
285 XPETRA_MONITOR("TpetraCrsMatrix::expertStaticFillComplete");
288
289 if(importer!=Teuchos::null) {
290 XPETRA_DYNAMIC_CAST( const TpetraImportClass , *importer, tImporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraImport.");
291 myImport = tImporter.getTpetra_Import();
292 }
293 if(exporter!=Teuchos::null) {
294 XPETRA_DYNAMIC_CAST( const TpetraExportClass , *exporter, tExporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraExport.");
295 myExport = tExporter.getTpetra_Export();
296 }
297
298 mtx_->expertStaticFillComplete(toTpetra(domainMap),toTpetra(rangeMap),myImport,myExport,params);
299 }
300
301 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
303
304 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
306
307 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
309
310 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
311 global_size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getGlobalNumRows() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumRows"); return mtx_->getGlobalNumRows(); }
312
313 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
314 global_size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getGlobalNumCols() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumCols"); return mtx_->getGlobalNumCols(); }
315
316 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
317 size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getLocalNumRows() const { XPETRA_MONITOR("TpetraCrsMatrix::getLocalNumRows"); return mtx_->getLocalNumRows(); }
318
319 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
320 size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getLocalNumCols() const { XPETRA_MONITOR("TpetraCrsMatrix::getLocalNumCols"); return mtx_->getLocalNumCols(); }
321
322 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
323 global_size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getGlobalNumEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumEntries"); return mtx_->getGlobalNumEntries(); }
324
325 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
326 size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getLocalNumEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getLocalNumEntries"); return mtx_->getLocalNumEntries(); }
327
328 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
329 size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNumEntriesInLocalRow(LocalOrdinal localRow) const { XPETRA_MONITOR("TpetraCrsMatrix::getNumEntriesInLocalRow"); return mtx_->getNumEntriesInLocalRow(localRow); }
330
331 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
332 size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const { XPETRA_MONITOR("TpetraCrsMatrix::getNumEntriesInGlobalRow"); return mtx_->getNumEntriesInGlobalRow(globalRow); }
333
334 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
335 size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getGlobalMaxNumRowEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalMaxNumRowEntries"); return mtx_->getGlobalMaxNumRowEntries(); }
336
337 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
338 size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getLocalMaxNumRowEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getLocalMaxNumRowEntries"); return mtx_->getLocalMaxNumRowEntries(); }
339
340 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
341 bool TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isLocallyIndexed() const { XPETRA_MONITOR("TpetraCrsMatrix::isLocallyIndexed"); return mtx_->isLocallyIndexed(); }
342
343 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
344 bool TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isGloballyIndexed() const { XPETRA_MONITOR("TpetraCrsMatrix::isGloballyIndexed"); return mtx_->isGloballyIndexed(); }
345
346 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
347 bool TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isFillComplete() const { XPETRA_MONITOR("TpetraCrsMatrix::isFillComplete"); return mtx_->isFillComplete(); }
348
349 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
350 bool TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isFillActive() const { XPETRA_MONITOR("TpetraCrsMatrix::isFillActive"); return mtx_->isFillActive(); }
351
352 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
353 typename ScalarTraits< Scalar >::magnitudeType TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getFrobeniusNorm() const { XPETRA_MONITOR("TpetraCrsMatrix::getFrobeniusNorm"); return mtx_->getFrobeniusNorm(); }
354
355 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
356 bool TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::supportsRowViews() const { XPETRA_MONITOR("TpetraCrsMatrix::supportsRowViews"); return mtx_->supportsRowViews(); }
357
358 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
359 void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const {
360 XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowCopy");
361 typename Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node>::nonconst_local_inds_host_view_type indices("indices",Indices.size());
362 typename Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::nonconst_values_host_view_type values("values", Values.size());
363
364 mtx_->getLocalRowCopy(LocalRow, indices, values, NumEntries);
365 for (size_t i=0;i<NumEntries; ++i) {
366 Indices[i]=indices(i);
367 Values[i]=values(i);
368 }
369 }
370
371 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
372 void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const {
373 XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowCopy");
374 typename Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node>::nonconst_global_inds_host_view_type indices("indices",Indices.size());
375 typename Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::nonconst_values_host_view_type values("values", Values.size());
376
377 mtx_->getGlobalRowCopy(GlobalRow, indices, values, NumEntries);
378 for (size_t i=0;i<NumEntries; ++i) {
379 Indices[i]=indices(i);
380 Values[i]=values(i);
381 }
382
383 }
384
385 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
387 XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowView");
388 typename Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node>::global_inds_host_view_type indices;
389 typename Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::values_host_view_type values;
390
391 mtx_->getGlobalRowView(GlobalRow, indices, values);
392 Indices = ArrayView<const GlobalOrdinal> (indices.data(), indices.extent(0));
393 Values = ArrayView<const Scalar> (reinterpret_cast<const Scalar*>(values.data()), values.extent(0));
394 }
395
396 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
398 XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowView");
399 typename Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node>::local_inds_host_view_type indices;
400 typename Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::values_host_view_type values;
401
402 mtx_->getLocalRowView(LocalRow, indices, values);
403 Indices = ArrayView<const LocalOrdinal> (indices.data(), indices.extent(0));
404 Values = ArrayView<const Scalar> (reinterpret_cast<const Scalar*>(values.data()), values.extent(0));
405 }
406
407 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
408 void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta) const { XPETRA_MONITOR("TpetraCrsMatrix::apply"); mtx_->apply(toTpetra(X), toTpetra(Y), mode, alpha, beta); }
409
410 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
411 void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues, const RCP<Import<LocalOrdinal, GlobalOrdinal, Node> >& regionInterfaceImporter, const Teuchos::ArrayRCP<LocalOrdinal>& regionInterfaceLIDs ) const
412 {
413 XPETRA_MONITOR("TpetraCrsMatrix::apply(region)");
414 RCP<const Map< LocalOrdinal, GlobalOrdinal, Node >> regionInterfaceMap = regionInterfaceImporter->getTargetMap();
415 mtx_->localApply(toTpetra(X), toTpetra(Y), mode, alpha, beta);
416 if (sumInterfaceValues)
417 {
418 // preform communication to propagate local interface
419 // values to all the processor that share interfaces.
421 matvecInterfaceTmp->doImport(Y, *regionInterfaceImporter, INSERT);
422
423 // sum all contributions to interface values
424 // on all ranks
426 ArrayRCP<Scalar> interfaceData = matvecInterfaceTmp->getDataNonConst(0);
427 for(LocalOrdinal interfaceIdx = 0; interfaceIdx < static_cast<LocalOrdinal>(interfaceData.size()); ++interfaceIdx) {
428 YData[regionInterfaceLIDs[interfaceIdx]] += interfaceData[interfaceIdx];
429 }
430 }
431 }
432
433 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
435
436 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
438
439 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
440 std::string TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::description() const { XPETRA_MONITOR("TpetraCrsMatrix::description"); return mtx_->description(); }
441
442 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
443 void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const { XPETRA_MONITOR("TpetraCrsMatrix::describe"); mtx_->describe(out, verbLevel); }
444
445 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
447 XPETRA_MONITOR("TpetraCrsMatrix::setObjectLabel");
449 mtx_->setObjectLabel(objectLabel);
450 }
451
452
453
454 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
456 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>(*(matrix.mtx_),Teuchos::Copy))) {}
457
458 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
460 XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
461 XPETRA_DYNAMIC_CAST(TpetraVectorClass, diag, tDiag, "Xpetra::TpetraCrsMatrix.getLocalDiagCopy() only accept Xpetra::TpetraVector as input arguments.");
462 mtx_->getLocalDiagCopy(*tDiag.getTpetra_Vector());
463 }
464
465 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
467 XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagOffsets");
468 mtx_->getLocalDiagOffsets(offsets);
469 }
470
471 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
473 XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
474 mtx_->getLocalDiagCopy(*(toTpetra(diag)), offsets);
475 }
476
477 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
479 XPETRA_MONITOR("TpetraCrsMatrix::replaceDiag");
480 Tpetra::replaceDiagonalCrsMatrix(*mtx_, *(toTpetra(diag)));
481 }
482
483 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
485 XPETRA_MONITOR("TpetraCrsMatrix::leftScale");
486 mtx_->leftScale(*(toTpetra(x)));
487 }
488
489 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
491 XPETRA_MONITOR("TpetraCrsMatrix::rightScale");
492 mtx_->rightScale(*(toTpetra(x)));
493 }
494
495 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
497
498 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
501 XPETRA_MONITOR("TpetraCrsMatrix::doImport");
502
503 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
505 //mtx_->doImport(toTpetraCrsMatrix(source), *tImporter.getTpetra_Import(), toTpetra(CM));
506 mtx_->doImport(*v, toTpetra(importer), toTpetra(CM));
507 }
508
510 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
513 XPETRA_MONITOR("TpetraCrsMatrix::doExport");
514
515 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
517 mtx_->doExport(*v, toTpetra(importer), toTpetra(CM));
518
519 }
520
521 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
524 XPETRA_MONITOR("TpetraCrsMatrix::doImport");
525
526 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
528 mtx_->doImport(*v, toTpetra(exporter), toTpetra(CM));
529
530 }
531
532 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
535 XPETRA_MONITOR("TpetraCrsMatrix::doExport");
536
537 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
539 mtx_->doExport(*v, toTpetra(exporter), toTpetra(CM));
540
541 }
542
543 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
545 XPETRA_MONITOR("TpetraCrsMatrix::removeEmptyProcessesInPlace");
546 mtx_->removeEmptyProcessesInPlace(toTpetra(newMap));
547 }
548
549 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
553
554 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
555 TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::TpetraCrsMatrix(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &mtx) : mtx_(mtx) { }
556
557 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
559
561 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
563
565 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
567 const MultiVector & B,
568 MultiVector & R) const {
569 Tpetra::Details::residual(*mtx_,toTpetra(X),toTpetra(B),toTpetra(R));
570 }
571
572
575// End of TpetrCrsMatrix class definition //
578
579} // Xpetra namespace
580
581#endif // XPETRA_TPETRACRSMATRIX_DEF_HPP
Copy
#define XPETRA_MONITOR(funcName)
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
size_type size() const
void resize(const size_type n, const T &val=T())
size_type size() const
virtual void setObjectLabel(const std::string &objectLabel)
bool is_null() const
T * get() const
static Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Build(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Constructor specifying the number of non-zeros for all rows.
virtual Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)=0
View of the local values in a particular vector of this multivector.
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
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.
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
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.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
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.
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_type local_matrix_type
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
void getLocalDiagCopy(Vector &diag) const
Get a copy of the diagonal entries owned by this node, with local row indices.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries in the (locally owned) global row.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
void rightScale(const Vector &x)
Right scale operator with given vector values.
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...
global_size_t getGlobalNumCols() const
Number of global columns in the matrix.
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.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this matrix.
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 residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const
Compute a residual R = B - (*this) * X.
TpetraCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraCrsMatrixClass
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > mtx_
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrix
void resumeFill(const RCP< ParameterList > &params=null)
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
void leftScale(const Vector &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.
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
bool hasMatrix() const
Does this have an underlying matrix.
size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
std::string description() const
A simple one-line description of this object.
RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrix() const
Get the underlying Tpetra matrix.
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
size_t getLocalNumRows() const
Returns the number of matrix rows owned on the calling node.
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalarThis.
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
bool isFillActive() const
Returns true if the matrix is in edit mode.
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.
void setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)
Sets the 1D pointer arrays of the graph.
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)
size_t getLocalNumCols() const
Returns the number of columns connected to the locally owned rows of this matrix.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in 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 doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
global_size_t getGlobalNumRows() const
Number of global elements in the row map of this 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.
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
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.
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 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 replaceDiag(const Vector &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.
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
void setObjectLabel(const std::string &objectLabel)
void apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const
Computes the sparse matrix-multivector multiplication.
Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > Vector
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
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....
bool isFillComplete() const
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called.
size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
basic_FancyOStream< char > FancyOStream
Xpetra namespace
size_t global_size_t
Global size_t object.
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
CombineMode
Xpetra::Combine Mode enumerable type.