Xpetra  Version of the Day
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 
53 namespace 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, Tpetra::StaticProfile, 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(), Tpetra::StaticProfile, 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, Tpetra::StaticProfile, 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(), Tpetra::StaticProfile, params))) { }
70 
71  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
73  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(graph), params))) { }
74 
75 
76 
77  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
80  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
81  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
83  {
84  typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
85  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
86  RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
87 
88  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
89  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
90  mtx_=Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(importer),myDomainMap,myRangeMap,params);
91  bool restrictComm=false;
92  if(!params.is_null()) restrictComm = params->get("Restrict Communicator",restrictComm);
93  if(restrictComm && mtx_->getRowMap().is_null()) mtx_=Teuchos::null;
94 
95  }
96 
97  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
100  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
101  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
103  {
104  typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
105  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
106  RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
107 
108  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
109  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
110  mtx_=Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(exporter),myDomainMap,myRangeMap,params);
111 
112  }
113 
114  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
116  const Import<LocalOrdinal,GlobalOrdinal,Node> & RowImporter,
117  const Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > DomainImporter,
118  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
119  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
121  {
122  typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
123  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
124  RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
125 
126  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
127  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
128 
129  mtx_=Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(RowImporter),toTpetra(*DomainImporter),myDomainMap,myRangeMap,params);
130  bool restrictComm=false;
131  if(!params.is_null()) restrictComm = params->get("Restrict Communicator",restrictComm);
132  if(restrictComm && mtx_->getRowMap().is_null()) mtx_=Teuchos::null;
133  }
134 
135  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
137  const Export<LocalOrdinal,GlobalOrdinal,Node> & RowExporter,
138  const Teuchos::RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > DomainExporter,
139  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
140  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
142  {
143  typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
144  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
145  RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
146 
147  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
148  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
149 
150  mtx_=Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(RowExporter),toTpetra(*DomainExporter),myDomainMap,myRangeMap,params);
151  }
152 
154 
155 
156 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
157 #ifdef HAVE_XPETRA_TPETRA
158  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
161  const local_matrix_type& lclMatrix,
163  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), lclMatrix, params))) { }
164 
165  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
167  const local_matrix_type& lclMatrix,
168  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rowMap,
169  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& colMap,
170  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
171  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
173  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(lclMatrix, toTpetra(rowMap), toTpetra(colMap), toTpetra(domainMap), toTpetra(rangeMap), params))) { }
174 
175  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
177  const local_matrix_type& lclMatrix,
178  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rowMap,
179  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& colMap,
180  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
181  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
182  const Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> >& importer,
183  const Teuchos::RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> >& exporter,
185  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(lclMatrix, toTpetra(rowMap), toTpetra(colMap), toTpetra(domainMap), toTpetra(rangeMap), toTpetra(importer), toTpetra(exporter), params))) { }
186 #endif
187 #endif
188 
189  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
191 
192  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
193  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); }
194 
195  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
196  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); }
197 
198  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
199  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); }
200 
201  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
202  void
204  const ArrayView<const LocalOrdinal> &cols,
205  const ArrayView<const Scalar> &vals)
206  {
207  XPETRA_MONITOR("TpetraCrsMatrix::replaceLocalValues");
208  typedef typename ArrayView<const LocalOrdinal>::size_type size_type;
209  const LocalOrdinal numValid =
210  mtx_->replaceLocalValues (localRow, cols, vals);
212  static_cast<size_type> (numValid) != cols.size (), std::runtime_error,
213  "replaceLocalValues returned " << numValid << " != cols.size() = " <<
214  cols.size () << ".");
215  }
216 
217  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
218  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::setAllToScalar(const Scalar &alpha) { XPETRA_MONITOR("TpetraCrsMatrix::setAllToScalar"); mtx_->setAllToScalar(alpha); }
219 
220  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
221  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::scale(const Scalar &alpha) { XPETRA_MONITOR("TpetraCrsMatrix::scale"); mtx_->scale(alpha); }
222 
223  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
225  { XPETRA_MONITOR("TpetraCrsMatrix::allocateAllValues"); rowptr.resize(getNodeNumRows()+1); colind.resize(numNonZeros); values.resize(numNonZeros);}
226 
227  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
229  { XPETRA_MONITOR("TpetraCrsMatrix::setAllValues"); mtx_->setAllValues(rowptr,colind,values); }
230 
231  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
233  { XPETRA_MONITOR("TpetraCrsMatrix::getAllValues"); mtx_->getAllValues(rowptr,colind,values); }
234 
235  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
237  { XPETRA_MONITOR("TpetraCrsMatrix::getAllValues"); mtx_->getAllValues(values); }
238 
239  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
241  { return mtx_->haveGlobalConstants();}
242 
243  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
244  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::resumeFill(const RCP< ParameterList > &params) { XPETRA_MONITOR("TpetraCrsMatrix::resumeFill"); mtx_->resumeFill(params); }
245 
246  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
247  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); }
248 
249  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
250  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::fillComplete(const RCP< ParameterList > &params) { XPETRA_MONITOR("TpetraCrsMatrix::fillComplete"); mtx_->fillComplete(params); }
251 
252 
253  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
255  XPETRA_MONITOR("TpetraCrsMatrix::replaceDomainMapAndImporter");
256  XPETRA_DYNAMIC_CAST( const TpetraImportClass , *newImporter, tImporter, "Xpetra::TpetraCrsMatrix::replaceDomainMapAndImporter only accepts Xpetra::TpetraImport.");
257  RCP<const Tpetra::Import<LocalOrdinal,GlobalOrdinal,Node> > myImport = tImporter.getTpetra_Import();
258  mtx_->replaceDomainMapAndImporter( toTpetra(newDomainMap),myImport);
259  }
260 
261  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
263  const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & rangeMap,
264  const RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > &importer,
265  const RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > &exporter,
266  const RCP<ParameterList> &params) {
267  XPETRA_MONITOR("TpetraCrsMatrix::expertStaticFillComplete");
270 
271  if(importer!=Teuchos::null) {
272  XPETRA_DYNAMIC_CAST( const TpetraImportClass , *importer, tImporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraImport.");
273  myImport = tImporter.getTpetra_Import();
274  }
275  if(exporter!=Teuchos::null) {
276  XPETRA_DYNAMIC_CAST( const TpetraExportClass , *exporter, tExporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraExport.");
277  myExport = tExporter.getTpetra_Export();
278  }
279 
280  mtx_->expertStaticFillComplete(toTpetra(domainMap),toTpetra(rangeMap),myImport,myExport,params);
281  }
282 
283  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
285 
286  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
288 
289  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
291 
292  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
293  global_size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getGlobalNumRows() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumRows"); return mtx_->getGlobalNumRows(); }
294 
295  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
296  global_size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getGlobalNumCols() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumCols"); return mtx_->getGlobalNumCols(); }
297 
298  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
299  size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNodeNumRows() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeNumRows"); return mtx_->getNodeNumRows(); }
300 
301  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
302  size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNodeNumCols() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeNumCols"); return mtx_->getNodeNumCols(); }
303 
304  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
305  global_size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getGlobalNumEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumEntries"); return mtx_->getGlobalNumEntries(); }
306 
307  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
308  size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNodeNumEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeNumEntries"); return mtx_->getNodeNumEntries(); }
309 
310  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
311  size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNumEntriesInLocalRow(LocalOrdinal localRow) const { XPETRA_MONITOR("TpetraCrsMatrix::getNumEntriesInLocalRow"); return mtx_->getNumEntriesInLocalRow(localRow); }
312 
313  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
314  size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const { XPETRA_MONITOR("TpetraCrsMatrix::getNumEntriesInGlobalRow"); return mtx_->getNumEntriesInGlobalRow(globalRow); }
315 
316  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
317  size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getGlobalMaxNumRowEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalMaxNumRowEntries"); return mtx_->getGlobalMaxNumRowEntries(); }
318 
319  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
320  size_t TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNodeMaxNumRowEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeMaxNumRowEntries"); return mtx_->getNodeMaxNumRowEntries(); }
321 
322  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
323  bool TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isLocallyIndexed() const { XPETRA_MONITOR("TpetraCrsMatrix::isLocallyIndexed"); return mtx_->isLocallyIndexed(); }
324 
325  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
326  bool TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isGloballyIndexed() const { XPETRA_MONITOR("TpetraCrsMatrix::isGloballyIndexed"); return mtx_->isGloballyIndexed(); }
327 
328  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
329  bool TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isFillComplete() const { XPETRA_MONITOR("TpetraCrsMatrix::isFillComplete"); return mtx_->isFillComplete(); }
330 
331  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
332  bool TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isFillActive() const { XPETRA_MONITOR("TpetraCrsMatrix::isFillActive"); return mtx_->isFillActive(); }
333 
334  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
335  typename ScalarTraits< Scalar >::magnitudeType TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getFrobeniusNorm() const { XPETRA_MONITOR("TpetraCrsMatrix::getFrobeniusNorm"); return mtx_->getFrobeniusNorm(); }
336 
337  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
338  bool TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::supportsRowViews() const { XPETRA_MONITOR("TpetraCrsMatrix::supportsRowViews"); return mtx_->supportsRowViews(); }
339 
340  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
341  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const {
342  XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowCopy");
343  typename Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node>::nonconst_local_inds_host_view_type indices("indices",Indices.size());
344  typename Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::nonconst_values_host_view_type values("values", Values.size());
345 
346  mtx_->getLocalRowCopy(LocalRow, indices, values, NumEntries);
347  for (size_t i=0;i<NumEntries; ++i) {
348  Indices[i]=indices(i);
349  Values[i]=values(i);
350  }
351  }
352 
353  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
354  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const {
355  XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowCopy");
356  typename Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node>::nonconst_global_inds_host_view_type indices("indices",Indices.size());
357  typename Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::nonconst_values_host_view_type values("values", Values.size());
358 
359  mtx_->getGlobalRowCopy(GlobalRow, indices, values, NumEntries);
360  for (size_t i=0;i<NumEntries; ++i) {
361  Indices[i]=indices(i);
362  Values[i]=values(i);
363  }
364 
365  }
366 
367  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
369  XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowView");
370  typename Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node>::global_inds_host_view_type indices;
371  typename Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::values_host_view_type values;
372 
373  mtx_->getGlobalRowView(GlobalRow, indices, values);
374  Indices = ArrayView<const GlobalOrdinal> (indices.data(), indices.extent(0));
375  Values = ArrayView<const Scalar> (reinterpret_cast<const Scalar*>(values.data()), values.extent(0));
376  }
377 
378  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
380  XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowView");
381  typename Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node>::local_inds_host_view_type indices;
382  typename Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::values_host_view_type values;
383 
384  mtx_->getLocalRowView(LocalRow, indices, values);
385  Indices = ArrayView<const LocalOrdinal> (indices.data(), indices.extent(0));
386  Values = ArrayView<const Scalar> (reinterpret_cast<const Scalar*>(values.data()), values.extent(0));
387  }
388 
389  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
390  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); }
391 
392  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
393  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
394  {
395  XPETRA_MONITOR("TpetraCrsMatrix::apply(region)");
396  RCP<const Map< LocalOrdinal, GlobalOrdinal, Node >> regionInterfaceMap = regionInterfaceImporter->getTargetMap();
397  mtx_->localApply(toTpetra(X), toTpetra(Y), mode, alpha, beta);
398  if (sumInterfaceValues)
399  {
400  // preform communication to propagate local interface
401  // values to all the processor that share interfaces.
403  matvecInterfaceTmp->doImport(Y, *regionInterfaceImporter, INSERT);
404 
405  // sum all contributions to interface values
406  // on all ranks
407  ArrayRCP<Scalar> YData = Y.getDataNonConst(0);
408  ArrayRCP<Scalar> interfaceData = matvecInterfaceTmp->getDataNonConst(0);
409  for(LocalOrdinal interfaceIdx = 0; interfaceIdx < static_cast<LocalOrdinal>(interfaceData.size()); ++interfaceIdx) {
410  YData[regionInterfaceLIDs[interfaceIdx]] += interfaceData[interfaceIdx];
411  }
412  }
413  }
414 
415  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
417 
418  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
420 
421  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
422  std::string TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::description() const { XPETRA_MONITOR("TpetraCrsMatrix::description"); return mtx_->description(); }
423 
424  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
425  void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const { XPETRA_MONITOR("TpetraCrsMatrix::describe"); mtx_->describe(out, verbLevel); }
426 
427  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
429  XPETRA_MONITOR("TpetraCrsMatrix::setObjectLabel");
431  mtx_->setObjectLabel(objectLabel);
432  }
433 
434 
435 
436  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
438  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>(*(matrix.mtx_),Teuchos::Copy))) {}
439 
440  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
442  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
443  XPETRA_DYNAMIC_CAST(TpetraVectorClass, diag, tDiag, "Xpetra::TpetraCrsMatrix.getLocalDiagCopy() only accept Xpetra::TpetraVector as input arguments.");
444  mtx_->getLocalDiagCopy(*tDiag.getTpetra_Vector());
445  }
446 
447  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
449  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagOffsets");
450  mtx_->getLocalDiagOffsets(offsets);
451  }
452 
453  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
455  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
456  mtx_->getLocalDiagCopy(*(toTpetra(diag)), offsets);
457  }
458 
459  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
461  XPETRA_MONITOR("TpetraCrsMatrix::replaceDiag");
462  Tpetra::replaceDiagonalCrsMatrix(*mtx_, *(toTpetra(diag)));
463  }
464 
465  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
467  XPETRA_MONITOR("TpetraCrsMatrix::leftScale");
468  mtx_->leftScale(*(toTpetra(x)));
469  }
470 
471  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
473  XPETRA_MONITOR("TpetraCrsMatrix::rightScale");
474  mtx_->rightScale(*(toTpetra(x)));
475  }
476 
477  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
479 
480  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
483  XPETRA_MONITOR("TpetraCrsMatrix::doImport");
484 
485  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
487  //mtx_->doImport(toTpetraCrsMatrix(source), *tImporter.getTpetra_Import(), toTpetra(CM));
488  mtx_->doImport(*v, toTpetra(importer), toTpetra(CM));
489  }
490 
492  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
495  XPETRA_MONITOR("TpetraCrsMatrix::doExport");
496 
497  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
499  mtx_->doExport(*v, toTpetra(importer), toTpetra(CM));
500 
501  }
502 
503  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
506  XPETRA_MONITOR("TpetraCrsMatrix::doImport");
507 
508  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
510  mtx_->doImport(*v, toTpetra(exporter), toTpetra(CM));
511 
512  }
513 
514  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
517  XPETRA_MONITOR("TpetraCrsMatrix::doExport");
518 
519  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
521  mtx_->doExport(*v, toTpetra(exporter), toTpetra(CM));
522 
523  }
524 
525  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
527  XPETRA_MONITOR("TpetraCrsMatrix::removeEmptyProcessesInPlace");
528  mtx_->removeEmptyProcessesInPlace(toTpetra(newMap));
529  }
530 
531  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
533  return ! mtx_.is_null ();
534  }
535 
536  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
537  TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::TpetraCrsMatrix(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &mtx) : mtx_(mtx) { }
538 
539  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
541 
543  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
545 
547  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
549  const MultiVector & B,
550  MultiVector & R) const {
551  Tpetra::Details::residual(*mtx_,toTpetra(X),toTpetra(B),toTpetra(R));
552  }
553 
554 
557 // End of TpetrCrsMatrix class definition //
560 
561 } // Xpetra namespace
562 
563 #endif // XPETRA_TPETRACRSMATRIX_DEF_HPP
#define XPETRA_MONITOR(funcName)
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
void resize(const size_type n, const T &val=T())
size_type size() const
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.
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.
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
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.
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...
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.
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
bool hasMatrix() const
Does this have an underlying matrix.
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.
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)
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.
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.
size_t getNodeMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
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.
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.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this matrix.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Copy
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
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)
CombineMode
Xpetra::Combine Mode enumerable type.