Xpetra_TpetraBlockCrsMatrix.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_TPETRABLOCKCRSMATRIX_HPP
47 #define XPETRA_TPETRABLOCKCRSMATRIX_HPP
48 
49 /* this file is automatically generated - do not edit (see script/tpetra.py) */
50 
52 
53 #include "Tpetra_Experimental_BlockCrsMatrix.hpp"
54 #include "Tpetra_CrsMatrix.hpp"
55 
56 #include "Xpetra_CrsMatrix.hpp"
57 #include "Xpetra_TpetraMap.hpp"
59 #include "Xpetra_TpetraVector.hpp"
61 //#include "Xpetra_TpetraCrsMatrix.hpp"
62 #include "Xpetra_Exceptions.hpp"
63 
64 
65 namespace Xpetra {
66 
67  template <class Scalar, class LocalOrdinal = int, class GlobalOrdinal = LocalOrdinal, class Node = KokkosClassic::DefaultNode::DefaultNodeType>
69  : public CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>//, public TpetraRowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
70  {
71 
72  // The following typedef are used by the XPETRA_DYNAMIC_CAST() macro.
77 
78  public:
79 
81 
82 
84  TpetraBlockCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
85  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
86 
89  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
90 
93  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
94 
97  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
98 
101  : mtx_(Teuchos::rcp(new Tpetra::Experimental::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(graph), params)))
102  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
103 
105  TpetraBlockCrsMatrix(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node> > &graph, const LocalOrdinal blockSize)
106  : mtx_(Teuchos::rcp(new Tpetra::Experimental::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(graph), blockSize))) { }
107 
108 
109 
110 
112  TpetraBlockCrsMatrix(const Teuchos::RCP<const Tpetra::Experimental::BlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& sourceMatrix,
114  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
115  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
116  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null)
117  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
118 
120  TpetraBlockCrsMatrix(const Teuchos::RCP<const Tpetra::Experimental::BlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& sourceMatrix,
122  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
123  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
124  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null)
125  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
126 
127 
129  virtual ~TpetraBlockCrsMatrix() { }
130 
132 
134 
135 
137  void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
138  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
139 
141  void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
142  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
143 
145  void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
146  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
147 
149  void replaceLocalValues (LocalOrdinal localRow,const ArrayView<const LocalOrdinal> &cols,const ArrayView<const Scalar> &vals)
150  {
151  XPETRA_MONITOR("TpetraBlockCrsMatrix::replaceLocalValues");
152  mtx_->replaceLocalValues(localRow,cols.getRawPtr(),vals.getRawPtr(),cols.size());
153  }
154 
156  void setAllToScalar(const Scalar &alpha) { XPETRA_MONITOR("TpetraBlockCrsMatrix::setAllToScalar"); mtx_->setAllToScalar(alpha); }
157 
159  void scale(const Scalar &alpha)
160  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
161 
163  //** \warning This is an expert-only routine and should not be called from user code. (not implemented)
164  void allocateAllValues(size_t numNonZeros,ArrayRCP<size_t> & rowptr, ArrayRCP<LocalOrdinal> & colind, ArrayRCP<Scalar> & values)
165  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
166 
168  void setAllValues(const ArrayRCP<size_t> & rowptr, const ArrayRCP<LocalOrdinal> & colind, const ArrayRCP<Scalar> & values)
169  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
170 
173  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
174 
176 
178 
179 
181  void resumeFill(const RCP< ParameterList > &params=null) { /*noop*/ }
182 
184  void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null) { /*noop*/ }
185 
187  void fillComplete(const RCP< ParameterList > &params=null) { /*noop*/ }
188 
189 
192  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
193 
196  const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & rangeMap,
197  const RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > &importer=Teuchos::null,
198  const RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > &exporter=Teuchos::null,
199  const RCP<ParameterList> &params=Teuchos::null)
200  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
201 
203 
205 
206 
208  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getRowMap"); return toXpetra(mtx_->getRowMap()); }
209 
211  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getColMap"); return toXpetra(mtx_->getColMap()); }
212 
215  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
216 
218  global_size_t getGlobalNumRows() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalNumRows"); return mtx_->getGlobalNumRows(); }
219 
221  global_size_t getGlobalNumCols() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalNumCols"); return mtx_->getGlobalNumCols(); }
222 
224  size_t getNodeNumRows() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getNodeNumRows"); return mtx_->getNodeNumRows(); }
225 
227  size_t getNodeNumCols() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getNodeNumCols"); return mtx_->getNodeNumCols(); }
228 
230  global_size_t getGlobalNumEntries() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalNumEntries"); return mtx_->getGlobalNumEntries(); }
231 
233  size_t getNodeNumEntries() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getNodeNumEntries"); return mtx_->getNodeNumEntries(); }
234 
236  size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getNumEntriesInLocalRow"); return mtx_->getNumEntriesInLocalRow(localRow); }
237 
239  global_size_t getGlobalNumDiags() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalNumDiags"); return mtx_->getGlobalNumDiags(); }
240 
242  size_t getNodeNumDiags() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getNodeNumDiags"); return mtx_->getNodeNumDiags(); }
243 
245  size_t getGlobalMaxNumRowEntries() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalMaxNumRowEntries"); return mtx_->getGlobalMaxNumRowEntries(); }
246 
248  size_t getNodeMaxNumRowEntries() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getNodeMaxNumRowEntries"); return mtx_->getNodeMaxNumRowEntries(); }
249 
251  bool isLocallyIndexed() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::isLocallyIndexed"); return mtx_->isLocallyIndexed(); }
252 
254  bool isGloballyIndexed() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::isGloballyIndexed"); return mtx_->isGloballyIndexed(); }
255 
257  bool isFillComplete() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::isFillComplete"); return mtx_->isFillComplete(); }
258 
260  bool isFillActive() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::isFillActive"); return false; }
261 
263  typename ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getFrobeniusNorm"); return mtx_->getFrobeniusNorm(); }
264 
266  bool supportsRowViews() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::supportsRowViews"); return mtx_->supportsRowViews(); }
267 
269  void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getLocalRowCopy"); mtx_->getLocalRowCopy(LocalRow, Indices, Values, NumEntries); }
270 
272  void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalRowView"); mtx_->getGlobalRowView(GlobalRow, indices, values); }
273 
275  void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalRowCopy"); mtx_->getGlobalRowCopy(GlobalRow, indices, values, numEntries); }
276 
278  void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getLocalRowView"); mtx_->getLocalRowView(LocalRow, indices, values); }
279 
281 
283 
284 
287 
289  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getDomainMap"); return toXpetra(mtx_->getDomainMap()); }
290 
292  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::getRangeMap"); return toXpetra(mtx_->getRangeMap()); }
293 
295 
297 
298 
300  std::string description() const { XPETRA_MONITOR("TpetraBlockCrsMatrix::description"); return mtx_->description(); }
301 
303  void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const { XPETRA_MONITOR("TpetraBlockCrsMatrix::describe"); mtx_->describe(out, verbLevel); }
304 
306 
309  : mtx_ (matrix.mtx_->template clone<Node> (matrix.mtx_->getNode ())) {}
310 
313  XPETRA_MONITOR("TpetraBlockCrsMatrix::getLocalDiagCopy");
314  XPETRA_DYNAMIC_CAST(TpetraVectorClass, diag, tDiag, "Xpetra::TpetraBlockCrsMatrix.getLocalDiagCopy() only accept Xpetra::TpetraVector as input arguments.")
315 
316  mtx_->getLocalDiagCopy(*tDiag.getTpetra_Vector());
317  }
318 
321  XPETRA_MONITOR("TpetraBlockCrsMatrix::getLocalDiagOffsets");
322  mtx_->getLocalDiagOffsets(offsets);
323  }
324 
327  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
328 
330  //{@
331 
334 
338  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
339 
343  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
344 
348  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
349 
353  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
354 
356  {throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented");}
357 
358  // @}
359 
360  template<class Node2>
363  }
364 
366 
367 
369  bool hasMatrix() const { return !mtx_.is_null();}
370 
372  TpetraBlockCrsMatrix(const Teuchos::RCP<Tpetra::Experimental::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &mtx) : mtx_(mtx) { }
373 
376 
379 
380 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
382 
383  local_matrix_type getLocalMatrix () const {
384  throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix does not support getLocalMatrix due to missing Kokkos::CrsMatrix in Tpetra's experimental implementation");
385  local_matrix_type ret;
386  return ret; // make compiler happy
387  }
388 #endif
389 
390 
391  private:
392 
394 
395  }; // TpetraBlockCrsMatrix class
396 
397  // TODO: move that elsewhere
398  // template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class LocalMatOps>
399  // const Tpetra::DistObject<char, LocalOrdinal, GlobalOrdinal, Node> toTpetraBlockCrsMatrix(const DistObject<char, LocalOrdinal, GlobalOrdinal, Node> &mtx) {
400  // typedef TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps> TpetraBlockCrsMatrixClass;
401  // XPETRA_DYNAMIC_CAST(const TpetraBlockCrsMatrixClass, mtx, tMtx, "toTpetra");
402  // return *tMtx.getTpetra_BlockCrsMatrix();
403  // }
404  //
405 
406 } // Xpetra namespace
407 
408 #define XPETRA_TPETRABLOCKCRSMATRIX_SHORT
409 #endif // XPETRA_TPETRABLOCKCRSMATRIX_HPP
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this (not implemented)
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 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...
Teuchos::RCP< NodeType > getNode()
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalarThis.
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.
global_size_t getGlobalNumDiags() const
Returns the number of global diagonal entries, based on global row/column index comparisons.
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs (not implemented)
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs (not implemented)
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
RCP< Tpetra::Experimental::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_BlockCrsMatrixNonConst() const
Get the underlying Tpetra matrix.
TpetraBlockCrsMatrix(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying a previously constructed graph ( not implemented )
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
global_size_t getGlobalNumCols() const
Number of global columns in the matrix.
bool hasMatrix() const
Does this have an underlying matrix.
Xpetra namespace
TpetraBlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraBlockCrsMatrixClass
size_type size() const
global_size_t getGlobalNumRows() const
Number of global elements in the row map of this matrix.
RCP< Tpetra::Experimental::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > mtx_
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local 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.
RCP< const CrsGraph< int, GlobalOrdinal > > toXpetra(const Epetra_CrsGraph &g)
TpetraBlockCrsMatrix(const Teuchos::RCP< Tpetra::Experimental::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &mtx)
TpetraBlockCrsMatrix constructor to wrap a Tpetra::BlockCrsMatrix object.
TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraVectorClass
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...
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)
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 getGlobalNumEntries() const
Returns the global number of entries in this matrix.
TpetraBlockCrsMatrix(const Teuchos::RCP< const Tpetra::Experimental::BlockCrsMatrix< 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 (not implemented(.
size_t getNodeNumDiags() const
Returns the number of local diagonal entries, based on global row/column index comparisons.
TpetraBlockCrsMatrix(const Teuchos::RCP< const Tpetra::Experimental::BlockCrsMatrix< 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 ( not implemented )
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this matrix. ...
TpetraBlockCrsMatrix(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 (not implemented) ...
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
size_t getNodeMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
TpetraBlockCrsMatrix(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph, const LocalOrdinal blockSize)
Constructor specifying a previously constructed graph & blocksize.
TpetraBlockCrsMatrix(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 (not implemented) ...
void apply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const
Computes the sparse matrix-multivector multiplication.
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Get a copy of the diagonal entries owned by this node, with local row idices {.
std::string description() const
A simple one-line description of this object.
void setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)
Sets the 1D pointer arrays of the graph (not impelmented)
bool isFillComplete() const
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called.
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.
RCP< TpetraBlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node2 > > clone(const RCP< Node2 > &node2) const
T * getRawPtr() const
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
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.
size_t global_size_t
Global size_t object.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
static const EVerbosityLevel verbLevel_default
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
RCP< const Tpetra::Experimental::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_BlockCrsMatrix() const
Get the underlying Tpetra matrix.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
bool supportsRowViews() const
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class.
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
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 getAllValues(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind, ArrayRCP< const Scalar > &values) const
Gets the 1D pointer arrays of the graph (not implemented)
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...
TpetraBlockCrsMatrix(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 (not implemented)
CombineMode
Xpetra::Combine Mode enumerable type.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this matrix.
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs (not implemented)
#define XPETRA_MONITOR(funcName)
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter).
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 isFillActive() const
Returns true if the matrix is in edit mode.
void fillComplete(const RCP< ParameterList > &params=null)
Signal that data entry is complete.
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
TpetraBlockCrsMatrix(const TpetraBlockCrsMatrix &matrix)
Deep copy constructor.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
TpetraBlockCrsMatrix(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 (not implemented) ...
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 resumeFill(const RCP< ParameterList > &params=null)
bool is_null() const