Xpetra_TpetraCrsMatrix.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_HPP
47 #define XPETRA_TPETRACRSMATRIX_HPP
48 
49 /* this file is automatically generated - do not edit (see scripts/tpetra.py) */
50 
51 // FIXME (mfh 03 Sep 2014) The above is probably not true anymore.
52 // Furthermore, I don't think anyone maintains the scripts.
53 // Feel free to correct this comment if I'm wrong.
54 
56 
57 #include "Tpetra_CrsMatrix.hpp"
58 
59 #include "Xpetra_CrsMatrix.hpp"
60 #include "Xpetra_TpetraMap.hpp"
62 #include "Xpetra_TpetraVector.hpp"
64 //#include "Xpetra_TpetraRowMatrix.hpp"
65 #include "Xpetra_Exceptions.hpp"
66 
67 namespace Xpetra {
68 
69  // TODO: move that elsewhere
70  // template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
71  // const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> toTpetraCrsMatrix(const Xpetra::DistObject<char, LocalOrdinal, GlobalOrdinal, Node> &);
72  //
73 
74  template<class Scalar = CrsMatrix<>::scalar_type,
75  class LocalOrdinal =
77  class GlobalOrdinal =
79  class Node =
82  : public CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> //, public TpetraRowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
83  {
84 
85  // The following typedef are used by the XPETRA_DYNAMIC_CAST() macro.
90 
91  // The following typedefs are used by the Kokkos interface
92 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
93 #ifdef HAVE_XPETRA_TPETRA
95 #endif
96 #endif
97 
98  public:
99 
101 
102 
104  TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
105  : mtx_ (Teuchos::rcp (new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> (toTpetra(rowMap), maxNumEntriesPerRow, toTpetra(pftype), params))) { }
106 
109  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> (toTpetra(rowMap), NumEntriesPerRowToAlloc, toTpetra(pftype), params))) { }
110 
113  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), maxNumEntriesPerRow, toTpetra(pftype), params))) { }
114 
117  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(rowMap), toTpetra(colMap), NumEntriesPerRowToAlloc, toTpetra(pftype), params))) { }
118 
121  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(graph), params))) { }
122 
123 
124 
128  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
129  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
130  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null)
131  {
132  typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
133  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
134  RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
135 
136  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
137  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
138  mtx_=Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(importer),myDomainMap,myRangeMap,params);
139  bool restrictComm=false;
140  if(!params.is_null()) restrictComm = params->get("Restrict Communicator",restrictComm);
141  if(restrictComm && mtx_->getRowMap().is_null()) mtx_=Teuchos::null;
142 
143  }
144 
148  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
149  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
150  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null)
151  {
152  typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
153  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
154  RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
155 
156  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
157  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
158  mtx_=Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(exporter),myDomainMap,myRangeMap,params);
159 
160  }
161 
164  const Import<LocalOrdinal,GlobalOrdinal,Node> & RowImporter,
165  const Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > DomainImporter,
166  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
167  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
169  {
170  typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
171  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
172  RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
173 
174  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
175  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
176 
177  mtx_=Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(RowImporter),toTpetra(*DomainImporter),myDomainMap,myRangeMap,params);
178  bool restrictComm=false;
179  if(!params.is_null()) restrictComm = params->get("Restrict Communicator",restrictComm);
180  if(restrictComm && mtx_->getRowMap().is_null()) mtx_=Teuchos::null;
181  }
182 
185  const Export<LocalOrdinal,GlobalOrdinal,Node> & RowExporter,
186  const Teuchos::RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > DomainExporter,
187  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
188  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
190  {
191  typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
192  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument.");//TODO: remove and use toTpetra()
193  RCP< const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > v = tSourceMatrix.getTpetra_CrsMatrix();
194 
195  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myDomainMap = domainMap!=Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
196  RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > myRangeMap = rangeMap!=Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
197 
198  mtx_=Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),toTpetra(RowExporter),toTpetra(*DomainExporter),myDomainMap,myRangeMap,params);
199  }
200 
201 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
202 #ifdef HAVE_XPETRA_TPETRA
225  const local_matrix_type& lclMatrix,
226  const Teuchos::RCP<Teuchos::ParameterList>& params = null)
227  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), lclMatrix, params))) { }
228 #endif
229 #endif
230 
232  virtual ~TpetraCrsMatrix() { }
233 
235 
237 
238 
240  void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals) { XPETRA_MONITOR("TpetraCrsMatrix::insertGlobalValues"); mtx_->insertGlobalValues(globalRow, cols, vals); }
241 
243  void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals) { XPETRA_MONITOR("TpetraCrsMatrix::insertLocalValues"); mtx_->insertLocalValues(localRow, cols, vals); }
244 
246  void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals) { XPETRA_MONITOR("TpetraCrsMatrix::replaceGlobalValues"); mtx_->replaceGlobalValues(globalRow, cols, vals); }
247 
249  void
250  replaceLocalValues (LocalOrdinal localRow,
251  const ArrayView<const LocalOrdinal> &cols,
252  const ArrayView<const Scalar> &vals)
253  {
254  XPETRA_MONITOR("TpetraCrsMatrix::replaceLocalValues");
255  typedef typename ArrayView<const LocalOrdinal>::size_type size_type;
256  const LocalOrdinal numValid =
257  mtx_->replaceLocalValues (localRow, cols, vals);
259  static_cast<size_type> (numValid) != cols.size (), std::runtime_error,
260  "replaceLocalValues returned " << numValid << " != cols.size() = " <<
261  cols.size () << ".");
262  }
263 
265  void setAllToScalar(const Scalar &alpha) { XPETRA_MONITOR("TpetraCrsMatrix::setAllToScalar"); mtx_->setAllToScalar(alpha); }
266 
268  void scale(const Scalar &alpha) { XPETRA_MONITOR("TpetraCrsMatrix::scale"); mtx_->scale(alpha); }
269 
271  //** \warning This is an expert-only routine and should not be called from user code. */
272  void allocateAllValues(size_t numNonZeros,ArrayRCP<size_t> & rowptr, ArrayRCP<LocalOrdinal> & colind, ArrayRCP<Scalar> & values)
273  { XPETRA_MONITOR("TpetraCrsMatrix::allocateAllValues"); rowptr.resize(getNodeNumRows()+1); colind.resize(numNonZeros); values.resize(numNonZeros);}
274 
276  void setAllValues(const ArrayRCP<size_t> & rowptr, const ArrayRCP<LocalOrdinal> & colind, const ArrayRCP<Scalar> & values)
277  { XPETRA_MONITOR("TpetraCrsMatrix::setAllValues"); mtx_->setAllValues(rowptr,colind,values); }
278 
281  { XPETRA_MONITOR("TpetraCrsMatrix::getAllValues"); mtx_->getAllValues(rowptr,colind,values); }
282 
283  bool haveGlobalConstants() const
284  { return mtx_->haveGlobalConstants();}
285 
287 
289 
290 
292  void resumeFill(const RCP< ParameterList > &params=null) { XPETRA_MONITOR("TpetraCrsMatrix::resumeFill"); mtx_->resumeFill(params); }
293 
295  void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null) { XPETRA_MONITOR("TpetraCrsMatrix::fillComplete"); mtx_->fillComplete(toTpetra(domainMap), toTpetra(rangeMap), params); }
296 
298  void fillComplete(const RCP< ParameterList > &params=null) { XPETRA_MONITOR("TpetraCrsMatrix::fillComplete"); mtx_->fillComplete(params); }
299 
300 
303  XPETRA_MONITOR("TpetraCrsMatrix::replaceDomainMapAndImporter");
304  XPETRA_DYNAMIC_CAST( const TpetraImportClass , *newImporter, tImporter, "Xpetra::TpetraCrsMatrix::replaceDomainMapAndImporter only accepts Xpetra::TpetraImport.");
305  RCP<const Tpetra::Import<LocalOrdinal,GlobalOrdinal,Node> > myImport = tImporter.getTpetra_Import();
306  mtx_->replaceDomainMapAndImporter( toTpetra(newDomainMap),myImport);
307  }
308 
311  const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & rangeMap,
312  const RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > &importer=Teuchos::null,
313  const RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > &exporter=Teuchos::null,
314  const RCP<ParameterList> &params=Teuchos::null) {
315  XPETRA_MONITOR("TpetraCrsMatrix::expertStaticFillComplete");
318 
319  if(importer!=Teuchos::null) {
320  XPETRA_DYNAMIC_CAST( const TpetraImportClass , *importer, tImporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraImport.");
321  myImport = tImporter.getTpetra_Import();
322  }
323  if(exporter!=Teuchos::null) {
324  XPETRA_DYNAMIC_CAST( const TpetraExportClass , *exporter, tExporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraExport.");
325  myExport = tExporter.getTpetra_Export();
326  }
327 
328  mtx_->expertStaticFillComplete(toTpetra(domainMap),toTpetra(rangeMap),myImport,myExport,params);
329  }
330 
332 
334 
335 
337  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const { XPETRA_MONITOR("TpetraCrsMatrix::getRowMap"); return toXpetra(mtx_->getRowMap()); }
338 
340  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const { XPETRA_MONITOR("TpetraCrsMatrix::getColMap"); return toXpetra(mtx_->getColMap()); }
341 
343  RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const { XPETRA_MONITOR("TpetraCrsMatrix::getCrsGraph"); return toXpetra(mtx_->getCrsGraph()); }
344 
346  global_size_t getGlobalNumRows() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumRows"); return mtx_->getGlobalNumRows(); }
347 
349  global_size_t getGlobalNumCols() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumCols"); return mtx_->getGlobalNumCols(); }
350 
352  size_t getNodeNumRows() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeNumRows"); return mtx_->getNodeNumRows(); }
353 
355  size_t getNodeNumCols() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeNumCols"); return mtx_->getNodeNumCols(); }
356 
358  global_size_t getGlobalNumEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumEntries"); return mtx_->getGlobalNumEntries(); }
359 
361  size_t getNodeNumEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeNumEntries"); return mtx_->getNodeNumEntries(); }
362 
364  size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const { XPETRA_MONITOR("TpetraCrsMatrix::getNumEntriesInLocalRow"); return mtx_->getNumEntriesInLocalRow(localRow); }
365 
367  global_size_t getGlobalNumDiags() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumDiags"); return mtx_->getGlobalNumDiags(); }
368 
370  size_t getNodeNumDiags() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeNumDiags"); return mtx_->getNodeNumDiags(); }
371 
373  size_t getGlobalMaxNumRowEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalMaxNumRowEntries"); return mtx_->getGlobalMaxNumRowEntries(); }
374 
376  size_t getNodeMaxNumRowEntries() const { XPETRA_MONITOR("TpetraCrsMatrix::getNodeMaxNumRowEntries"); return mtx_->getNodeMaxNumRowEntries(); }
377 
379  bool isLocallyIndexed() const { XPETRA_MONITOR("TpetraCrsMatrix::isLocallyIndexed"); return mtx_->isLocallyIndexed(); }
380 
382  bool isGloballyIndexed() const { XPETRA_MONITOR("TpetraCrsMatrix::isGloballyIndexed"); return mtx_->isGloballyIndexed(); }
383 
385  bool isFillComplete() const { XPETRA_MONITOR("TpetraCrsMatrix::isFillComplete"); return mtx_->isFillComplete(); }
386 
388  bool isFillActive() const { XPETRA_MONITOR("TpetraCrsMatrix::isFillActive"); return mtx_->isFillActive(); }
389 
391  typename ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const { XPETRA_MONITOR("TpetraCrsMatrix::getFrobeniusNorm"); return mtx_->getFrobeniusNorm(); }
392 
394  bool supportsRowViews() const { XPETRA_MONITOR("TpetraCrsMatrix::supportsRowViews"); return mtx_->supportsRowViews(); }
395 
397  void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const { XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowCopy"); mtx_->getLocalRowCopy(LocalRow, Indices, Values, NumEntries); }
398 
400  void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowView"); mtx_->getGlobalRowView(GlobalRow, indices, values); }
401 
403  void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const { XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowCopy"); mtx_->getGlobalRowCopy(GlobalRow, indices, values, numEntries); }
404 
406  void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const { XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowView"); mtx_->getLocalRowView(LocalRow, indices, values); }
407 
409 
411 
412 
415 
417  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const { XPETRA_MONITOR("TpetraCrsMatrix::getDomainMap"); return toXpetra(mtx_->getDomainMap()); }
418 
420  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const { XPETRA_MONITOR("TpetraCrsMatrix::getRangeMap"); return toXpetra(mtx_->getRangeMap()); }
421 
423 
425 
426 
428  std::string description() const { XPETRA_MONITOR("TpetraCrsMatrix::description"); return mtx_->description(); }
429 
431  void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const { XPETRA_MONITOR("TpetraCrsMatrix::describe"); mtx_->describe(out, verbLevel); }
432 
434 
436  TpetraCrsMatrix(const TpetraCrsMatrix& matrix)
437  : mtx_ (matrix.mtx_->template clone<Node> (matrix.mtx_->getNode ())) {}
438 
441  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
442  XPETRA_DYNAMIC_CAST(TpetraVectorClass, diag, tDiag, "Xpetra::TpetraCrsMatrix.getLocalDiagCopy() only accept Xpetra::TpetraVector as input arguments.");
443  mtx_->getLocalDiagCopy(*tDiag.getTpetra_Vector());
444  // mtx_->getLocalDiagCopy(toTpetra(diag));
445  }
446 
449  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagOffsets");
450  mtx_->getLocalDiagOffsets(offsets);
451  }
452 
455  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
456  //XPETRA_DYNAMIC_CAST(TpetraVectorClass, diag, tDiag, "Xpetra::TpetraCrsMatrix.getLocalDiagCopy() only accept Xpetra::TpetraVector as input arguments.");
457  //mtx_->getLocalDiagCopy(*tDiag.getTpetra_Vector(), offsets);
458  mtx_->getLocalDiagCopy(*(toTpetra(diag)), offsets);
459  }
460 
463  XPETRA_MONITOR("TpetraCrsMatrix::leftScale");
464  mtx_->leftScale(*(toTpetra(x)));
465  }
466 
469  XPETRA_MONITOR("TpetraCrsMatrix::rightScale");
470  mtx_->rightScale(*(toTpetra(x)));
471  }
472 
474  //{@
475 
478 
482  XPETRA_MONITOR("TpetraCrsMatrix::doImport");
483 
484  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
486  //mtx_->doImport(toTpetraCrsMatrix(source), *tImporter.getTpetra_Import(), toTpetra(CM));
487  mtx_->doImport(*v, toTpetra(importer), toTpetra(CM));
488  }
489 
493  XPETRA_MONITOR("TpetraCrsMatrix::doExport");
494 
495  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
497  mtx_->doExport(*v, toTpetra(importer), toTpetra(CM));
498 
499  }
500 
504  XPETRA_MONITOR("TpetraCrsMatrix::doImport");
505 
506  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
508  mtx_->doImport(*v, toTpetra(exporter), toTpetra(CM));
509 
510  }
511 
515  XPETRA_MONITOR("TpetraCrsMatrix::doExport");
516 
517  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments.");//TODO: remove and use toTpetra()
519  mtx_->doExport(*v, toTpetra(exporter), toTpetra(CM));
520 
521  }
522 
524  XPETRA_MONITOR("TpetraCrsMatrix::removeEmptyProcessesInPlace");
525  mtx_->removeEmptyProcessesInPlace(toTpetra(newMap));
526  }
527 
528  // @}
529 
530  template<class Node2>
532  clone (const RCP<Node2> &node2) const
533  {
535  }
536 
538 
539 
541  bool hasMatrix() const {
542  return ! mtx_.is_null ();
543  }
544 
546  TpetraCrsMatrix(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &mtx) : mtx_(mtx) { }
547 
550 
553 
554 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
555 #ifdef HAVE_XPETRA_TPETRA
556  local_matrix_type getLocalMatrix () const {
558  if(isFillComplete() == false)
559  throw std::runtime_error("Xpetra::EpetraCrsMatrix::getLocalMatrix: matrix must be filled and completed before you can access the data through the Kokkos interface!");
560 
561  return getTpetra_CrsMatrixNonConst()->getLocalMatrix();
562  }
563 #endif
564 #endif
565 
567 
568  private:
570  }; // TpetraCrsMatrix class
571 
572 #ifdef HAVE_XPETRA_EPETRA
573 
574 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
575  (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
576 
577  // specialization of TpetraCrsMatrix for GO=LO=int
578  template<class Scalar>
579  class TpetraCrsMatrix<Scalar,int,int,EpetraNode>
580  : public CrsMatrix<Scalar,int,int,EpetraNode> //, public TpetraRowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
581  {
582  typedef int LocalOrdinal;
583  typedef int GlobalOrdinal;
584  typedef EpetraNode Node;
585  // The following typedef are used by the XPETRA_DYNAMIC_CAST() macro.
590 
591  // The following typedefs are used by the Kokkos interface
592 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
593 #ifdef HAVE_XPETRA_TPETRA
595 #endif
596 #endif
597 
598  public:
599 
601 
602 
604  TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) {
606  }
607 
611  }
612 
614  TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) {
616  }
617 
621  }
622 
626  }
627 
628 
629 
633  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
634  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
635  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
637  }
638 
642  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
643  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
644  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
646  }
647 
650  const Import<LocalOrdinal,GlobalOrdinal,Node> & RowImporter,
651  const Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > DomainImporter,
652  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
653  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
655  {
657  }
658 
661  const Export<LocalOrdinal,GlobalOrdinal,Node> & RowExporter,
662  const Teuchos::RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > DomainExporter,
663  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
664  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
666  {
668  }
669 
670 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
671 #ifdef HAVE_XPETRA_TPETRA
694  const local_matrix_type& lclMatrix,
695  const Teuchos::RCP<Teuchos::ParameterList>& params = null) {
696  XPETRA_TPETRA_ETI_EXCEPTION( typeid(TpetraCrsMatrix<Scalar,int,int,EpetraNode>).name() , "TpetraCrsMatrix<int,int>", "int", typeid(EpetraNode).name() );
697  }
698 #endif
699 #endif
700 
702  virtual ~TpetraCrsMatrix() { }
703 
705 
707 
708 
710  void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals) { }
711 
713  void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals) { }
714 
716  void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals) { }
717 
719  void
720  replaceLocalValues (LocalOrdinal localRow,
721  const ArrayView<const LocalOrdinal> &cols,
722  const ArrayView<const Scalar> &vals)
723  { }
724 
726  void setAllToScalar(const Scalar &alpha) { }
727 
729  void scale(const Scalar &alpha) { }
730 
732  //** \warning This is an expert-only routine and should not be called from user code. */
733  void allocateAllValues(size_t numNonZeros,ArrayRCP<size_t> & rowptr, ArrayRCP<LocalOrdinal> & colind, ArrayRCP<Scalar> & values) { }
734 
736  void setAllValues(const ArrayRCP<size_t> & rowptr, const ArrayRCP<LocalOrdinal> & colind, const ArrayRCP<Scalar> & values) { }
737 
740 
741  bool haveGlobalConstants() const { return false;}
742 
743 
745 
747 
748 
750  void resumeFill(const RCP< ParameterList > &params=null) { }
751 
753  void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null) { }
754 
756  void fillComplete(const RCP< ParameterList > &params=null) { }
757 
758 
761 
764  const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & rangeMap,
765  const RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > &importer=Teuchos::null,
766  const RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > &exporter=Teuchos::null,
767  const RCP<ParameterList> &params=Teuchos::null) { }
768 
770 
772 
773 
775  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const { return Teuchos::null; }
776 
778  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const { return Teuchos::null; }
779 
782 
784  global_size_t getGlobalNumRows() const { return 0; }
785 
787  global_size_t getGlobalNumCols() const { return 0; }
788 
790  size_t getNodeNumRows() const { return 0; }
791 
793  size_t getNodeNumCols() const { return 0; }
794 
796  global_size_t getGlobalNumEntries() const { return 0; }
797 
799  size_t getNodeNumEntries() const { return 0; }
800 
802  size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const { return 0; }
803 
805  global_size_t getGlobalNumDiags() const { return 0; }
806 
808  size_t getNodeNumDiags() const { return 0; }
809 
811  size_t getGlobalMaxNumRowEntries() const { return 0; }
812 
814  size_t getNodeMaxNumRowEntries() const { return 0; }
815 
817  bool isLocallyIndexed() const { return false; }
818 
820  bool isGloballyIndexed() const { return false; }
821 
823  bool isFillComplete() const { return false; }
824 
826  bool isFillActive() const { return false; }
827 
830 
832  bool supportsRowViews() const { return false; }
833 
835  void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const { }
836 
838  void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const { }
839 
841  void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const { }
842 
844  void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const { }
845 
847 
849 
850 
853 
855  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const { return Teuchos::null; }
856 
858  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const { return Teuchos::null; }
859 
861 
863 
864 
866  std::string description() const { return std::string(""); }
867 
870 
872 
874  TpetraCrsMatrix(const TpetraCrsMatrix& matrix) {}
875 
878 
881 
884 
887 
890 
892  //{@
893 
896 
900 
904 
908 
912 
914 
915  // @}
916 
917  template<class Node2>
919  clone (const RCP<Node2> &node2) const { return Teuchos::null; }
920 
922 
923 
925  bool hasMatrix() const { return false; }
926 
928  TpetraCrsMatrix(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &mtx) {
930  }
931 
934 
937 
938 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
939 #ifdef HAVE_XPETRA_TPETRA
940  local_matrix_type getLocalMatrix () const {
942  if(isFillComplete() == false)
943  throw std::runtime_error("Xpetra::EpetraCrsMatrix::getLocalMatrix: matrix must be filled and completed before you can access the data through the Kokkos interface!");
944  local_matrix_type ret;
945  return ret;
946  }
947 #endif
948 #endif
949 
951  }; // TpetraCrsMatrix class (specialization for GO=int, NO=EpetraNode)
952 #endif
953 
954 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))) || \
955  (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))))
956 
957  // specialization of TpetraCrsMatrix for GO=long long, NO=EpetraNode
958  template<class Scalar>
959  class TpetraCrsMatrix<Scalar,int,long long,EpetraNode>
960  : public CrsMatrix<Scalar,int,long long,EpetraNode> //, public TpetraRowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
961  {
962  typedef int LocalOrdinal;
963  typedef long long GlobalOrdinal;
964  typedef EpetraNode Node;
965  // The following typedef are used by the XPETRA_DYNAMIC_CAST() macro.
970 
971  // The following typedefs are used by the Kokkos interface
972 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
973 #ifdef HAVE_XPETRA_TPETRA
975 #endif
976 #endif
977 
978  public:
979 
981 
982 
984  TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) {
986  }
987 
991  }
992 
994  TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) {
996  }
997 
1001  }
1002 
1006  }
1007 
1008 
1009 
1012  const Import<LocalOrdinal,GlobalOrdinal,Node> & importer,
1013  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
1014  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
1015  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
1017  }
1018 
1021  const Export<LocalOrdinal,GlobalOrdinal,Node> & exporter,
1022  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
1023  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
1024  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
1026  }
1027 
1030  const Import<LocalOrdinal,GlobalOrdinal,Node> & RowImporter,
1031  const Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > DomainImporter,
1032  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
1033  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
1035  {
1037  }
1038 
1041  const Export<LocalOrdinal,GlobalOrdinal,Node> & RowExporter,
1042  const Teuchos::RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > DomainExporter,
1043  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
1044  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
1046  {
1048  }
1049 
1050 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
1051 #ifdef HAVE_XPETRA_TPETRA
1073  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& colMap,
1074  const local_matrix_type& lclMatrix,
1075  const Teuchos::RCP<Teuchos::ParameterList>& params = null) {
1077  }
1078 #endif
1079 #endif
1080 
1082  virtual ~TpetraCrsMatrix() { }
1083 
1085 
1087 
1088 
1090  void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals) { }
1091 
1093  void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals) { }
1094 
1096  void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals) { }
1097 
1099  void
1100  replaceLocalValues (LocalOrdinal localRow,
1101  const ArrayView<const LocalOrdinal> &cols,
1102  const ArrayView<const Scalar> &vals)
1103  { }
1104 
1106  void setAllToScalar(const Scalar &alpha) { }
1107 
1109  void scale(const Scalar &alpha) { }
1110 
1112  //** \warning This is an expert-only routine and should not be called from user code. */
1113  void allocateAllValues(size_t numNonZeros,ArrayRCP<size_t> & rowptr, ArrayRCP<LocalOrdinal> & colind, ArrayRCP<Scalar> & values) { }
1114 
1116  void setAllValues(const ArrayRCP<size_t> & rowptr, const ArrayRCP<LocalOrdinal> & colind, const ArrayRCP<Scalar> & values) { }
1117 
1120 
1121  bool haveGlobalConstants() const { return false;}
1122 
1124 
1126 
1127 
1129  void resumeFill(const RCP< ParameterList > &params=null) { }
1130 
1132  void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null) { }
1133 
1135  void fillComplete(const RCP< ParameterList > &params=null) { }
1136 
1137 
1140 
1143  const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & rangeMap,
1144  const RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > &importer=Teuchos::null,
1145  const RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > &exporter=Teuchos::null,
1146  const RCP<ParameterList> &params=Teuchos::null) { }
1147 
1149 
1151 
1152 
1154  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const { return Teuchos::null; }
1155 
1157  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const { return Teuchos::null; }
1158 
1161 
1163  global_size_t getGlobalNumRows() const { return 0; }
1164 
1166  global_size_t getGlobalNumCols() const { return 0; }
1167 
1169  size_t getNodeNumRows() const { return 0; }
1170 
1172  size_t getNodeNumCols() const { return 0; }
1173 
1175  global_size_t getGlobalNumEntries() const { return 0; }
1176 
1178  size_t getNodeNumEntries() const { return 0; }
1179 
1181  size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const { return 0; }
1182 
1184  global_size_t getGlobalNumDiags() const { return 0; }
1185 
1187  size_t getNodeNumDiags() const { return 0; }
1188 
1190  size_t getGlobalMaxNumRowEntries() const { return 0; }
1191 
1193  size_t getNodeMaxNumRowEntries() const { return 0; }
1194 
1196  bool isLocallyIndexed() const { return false; }
1197 
1199  bool isGloballyIndexed() const { return false; }
1200 
1202  bool isFillComplete() const { return false; }
1203 
1205  bool isFillActive() const { return false; }
1206 
1209 
1211  bool supportsRowViews() const { return false; }
1212 
1214  void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const { }
1215 
1217  void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const { }
1218 
1220  void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const { }
1221 
1223  void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const { }
1224 
1226 
1228 
1229 
1232 
1234  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const { return Teuchos::null; }
1235 
1237  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const { return Teuchos::null; }
1238 
1240 
1242 
1243 
1245  std::string description() const { return std::string(""); }
1246 
1249 
1251 
1253  TpetraCrsMatrix(const TpetraCrsMatrix& matrix) {}
1254 
1257 
1260 
1263 
1266 
1269 
1271  //{@
1272 
1275 
1279 
1283 
1287 
1291 
1293 
1294  // @}
1295 
1296  template<class Node2>
1298  clone (const RCP<Node2> &node2) const { return Teuchos::null; }
1299 
1301 
1302 
1304  bool hasMatrix() const { return false; }
1305 
1307  TpetraCrsMatrix(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &mtx) {
1309  }
1310 
1313 
1316 
1317 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
1318 #ifdef HAVE_XPETRA_TPETRA
1319  local_matrix_type getLocalMatrix () const {
1321  if(isFillComplete() == false)
1322  throw std::runtime_error("Xpetra::EpetraCrsMatrix::getLocalMatrix: matrix must be filled and completed before you can access the data through the Kokkos interface!");
1323  local_matrix_type ret;
1324  return ret;
1325  }
1326 #endif
1327 #endif
1328 
1330  }; // TpetraCrsMatrix class (specialization for GO=long long, NO=EpetraNode)
1331 #endif
1332 
1333 #endif // HAVE_XPETRA_EPETRA
1334 
1335 } // Xpetra namespace
1336 
1337 #define XPETRA_TPETRACRSMATRIX_SHORT
1338 #endif // XPETRA_TPETRACRSMATRIX_HPP
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node2 > > clone(const RCP< Node2 > &node2) const
TpetraCrsMatrix(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying a previously constructed graph.
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor for a fused export.
TpetraCrsMatrix(const TpetraCrsMatrix &matrix)
Deep copy constructor.
bool isFillComplete() const
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called.
void fillComplete(const RCP< ParameterList > &params=null)
Signal that data entry is complete.
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag, const Teuchos::ArrayView< const size_t > &offsets) const
Get a copy of the diagonal entries owned by this node, with local row indices.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying column Map and fixed number of entries for each row.
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalarThis.
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 getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag, const Teuchos::ArrayView< const size_t > &offsets) const
Get a copy of the diagonal entries owned by this node, with local row indices.
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraVectorClass
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this matrix.
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.
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this matrix. ...
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)
bool isFillComplete() const
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called.
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
Teuchos::RCP< NodeType > getNode()
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.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying fixed number of entries for each row.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
TpetraCrsMatrix(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying a previously constructed graph.
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
global_size_t getGlobalNumDiags() const
Returns the number of global diagonal entries, based on global row/column index comparisons.
TpetraCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraCrsMatrixClass
TpetraCrsMatrix(const Teuchos::RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &mtx)
TpetraCrsMatrix constructor to wrap a Tpetra::CrsMatrix object.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
GlobalOrdinal global_ordinal_type
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.
RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrix() const
Get the underlying Tpetra matrix.
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.
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...
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying column Map and number of entries in each row.
T & get(const std::string &name, T def_value)
size_t getNodeNumDiags() const
Returns the number of local diagonal entries, based on global row/column index comparisons.
size_t getNodeNumDiags() const
Returns the number of local diagonal entries, based on global row/column index comparisons.
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > mtx_
RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node2 > > clone(const RCP< Node2 > &node2) const
void setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)
Sets the 1D pointer arrays of the graph.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
void fillComplete(const RCP< ParameterList > &params=null)
Signal that data entry is complete.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale operator with given vector values.
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 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 getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
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 replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
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.
Xpetra namespace
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
TpetraCrsMatrix(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying a previously constructed graph.
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying fixed number of entries for each row.
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
global_size_t getGlobalNumDiags() const
Returns the number of global diagonal entries, based on global row/column index comparisons.
size_type size() const
LocalOrdinal local_ordinal_type
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 setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)
Sets the 1D pointer arrays of the graph.
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrixNonConst() const
Get the underlying Tpetra matrix.
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
global_size_t getGlobalNumRows() const
Number of global elements in the row map of this matrix.
#define XPETRA_TPETRA_ETI_EXCEPTION(cl, obj, go, node)
void fillComplete(const RCP< ParameterList > &params=null)
Signal that data entry is complete.
size_t getNodeNumDiags() const
Returns the number of local diagonal entries, based on global row/column index comparisons.
size_t getNodeMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this matrix.
T * get() const
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...
bool supportsRowViews() const
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class.
std::string description() const
A simple one-line description of this object.
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale operator with given vector values.
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 is_null() const
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
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...
bool isFillComplete() const
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called.
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalarThis.
void resize(const size_type n, const T &val=T())
global_size_t getGlobalNumDiags() const
Returns the number of global diagonal entries, based on global row/column index comparisons.
global_size_t getGlobalNumCols() const
Number of global columns in the matrix.
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
void getAllValues(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind, ArrayRCP< const Scalar > &values) const
Gets the 1D pointer arrays of the graph.
TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraVectorClass
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying (possibly different) number of entries in each row.
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor for a fused import.
std::string description() const
A simple one-line description of this object.
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
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.
bool supportsRowViews() const
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class.
TpetraCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraCrsMatrixClass
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
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.
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
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.
virtual ~TpetraCrsMatrix()
Destructor.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Export< LocalOrdinal, GlobalOrdinal, Node > &RowExporter, const Teuchos::RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > DomainExporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
Constructor for a fused export.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying (possibly different) number of entries in each row.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
global_size_t getGlobalNumCols() const
Number of global columns in the matrix.
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrixNonConst() const
Get the underlying Tpetra matrix.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
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...
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this matrix. ...
RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrix() const
Get the underlying Tpetra matrix.
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
TpetraCrsMatrix(const Teuchos::RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &mtx)
TpetraCrsMatrix constructor to wrap a Tpetra::CrsMatrix object.
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter).
TpetraCrsMatrix(const Teuchos::RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &mtx)
TpetraCrsMatrix constructor to wrap a Tpetra::CrsMatrix object.
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.
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Import< LocalOrdinal, GlobalOrdinal, Node > &RowImporter, const Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > DomainImporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
Constructor for a fused import.
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 removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)
bool hasMatrix() const
Does this have an underlying matrix.
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor for a fused import.
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
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.
global_size_t getGlobalNumCols() const
Number of global columns in the matrix.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying column Map and number of entries in each row.
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.
bool hasMatrix() const
Does this have an underlying 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 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 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.
TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraVectorClass
void resumeFill(const RCP< ParameterList > &params=null)
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale operator with given vector values.
bool isFillActive() const
Returns true if the matrix is in edit mode.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
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.
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
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.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Export< LocalOrdinal, GlobalOrdinal, Node > &RowExporter, const Teuchos::RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > DomainExporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
Constructor for a fused export.
size_t global_size_t
Global size_t object.
RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node2 > > clone(const RCP< Node2 > &node2) const
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
static const EVerbosityLevel verbLevel_default
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...
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor for a fused export.
RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrix() const
Get the underlying Tpetra matrix.
size_t getNodeMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
static magnitudeType magnitude(T a)
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Export< LocalOrdinal, GlobalOrdinal, Node > &RowExporter, const Teuchos::RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > DomainExporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
Constructor for a fused export.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
global_size_t getGlobalNumRows() const
Number of global elements in the row map of this matrix.
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Get a copy of the diagonal entries owned by this node, with local row idices.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
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...
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...
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying (possibly different) number of entries in each row.
void setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)
Sets the 1D pointer arrays of the graph.
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying column Map and fixed number of entries for each row.
bool hasMatrix() const
Does this have an underlying matrix.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale operator with given vector values.
bool isFillActive() const
Returns true if the matrix is in edit mode.
size_t getNodeMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying column Map and fixed number of entries for each row.
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Import< LocalOrdinal, GlobalOrdinal, Node > &RowImporter, const Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > DomainImporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
Constructor for a fused import.
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 insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag, const Teuchos::ArrayView< const size_t > &offsets) const
Get a copy of the diagonal entries owned by this node, with local row indices.
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Get a copy of the diagonal entries owned by this node, with local row idices.
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 getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
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.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
CombineMode
Xpetra::Combine Mode enumerable type.
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...
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
TpetraCrsMatrix(const TpetraCrsMatrix &matrix)
Deep copy constructor.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this matrix.
#define XPETRA_MONITOR(funcName)
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale operator with given vector values.
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter).
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.
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Import< LocalOrdinal, GlobalOrdinal, Node > &RowImporter, const Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > DomainImporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
Constructor for a fused import.
global_size_t getGlobalNumRows() const
Number of global elements in the row map of this matrix.
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor for a fused export.
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
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.
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrixNonConst() const
Get the underlying Tpetra matrix.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale operator with given vector values.
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
std::string description() const
A simple one-line description of this object.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
TpetraCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraCrsMatrixClass
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.
bool supportsRowViews() const
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class.
TpetraCrsMatrix(const TpetraCrsMatrix &matrix)
Deep copy constructor.
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 insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying fixed number of entries for each row.
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalarThis.
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor for a fused import.
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...
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying column Map and number of entries in each row.
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter).
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Get a copy of the diagonal entries owned by this node, with local row idices.
bool isFillActive() const
Returns true if the matrix is in edit mode.