Xpetra  Version of the Day
Xpetra_TpetraBlockCrsMatrix_decl.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_DECL_HPP
47 #define XPETRA_TPETRABLOCKCRSMATRIX_DECL_HPP
48 
49 /* this file is automatically generated - do not edit (see script/tpetra.py) */
50 
52 
53 #include "Tpetra_BlockCrsMatrix.hpp"
54 #include "Tpetra_CrsMatrix.hpp"
55 
56 #include "Xpetra_CrsMatrix.hpp"
61 #include "Xpetra_Exceptions.hpp"
62 
63 
64 namespace Xpetra {
65 
66  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node = KokkosClassic::DefaultNode::DefaultNodeType>
68  : public CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>//, public TpetraRowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
69  {
70 
71  // The following typedef are used by the XPETRA_DYNAMIC_CAST() macro.
76 
77  public:
78 
80 
83  size_t maxNumEntriesPerRow,
84  const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null);
85 
86 
89  const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc,
90  const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null);
91 
92 
96  size_t maxNumEntriesPerRow,
97  const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null);
98 
99 
103  const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc,
104  const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null);
105 
106 
109  const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null);
110 
111 
114  const LocalOrdinal blockSize);
115 
116 
118  TpetraBlockCrsMatrix(const Teuchos::RCP<const Tpetra::BlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& sourceMatrix,
120  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
121  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
122  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
123 
124 
126  TpetraBlockCrsMatrix(const Teuchos::RCP<const Tpetra::BlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& sourceMatrix,
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 
134  TpetraBlockCrsMatrix(const Teuchos::RCP<const Tpetra::BlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& sourceMatrix,
135  const Import<LocalOrdinal,GlobalOrdinal,Node> & RowImporter,
136  const Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > DomainImporter,
137  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
138  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
140 
141 
143  TpetraBlockCrsMatrix(const Teuchos::RCP<const Tpetra::BlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& sourceMatrix,
144  const Export<LocalOrdinal,GlobalOrdinal,Node> & RowExporter,
145  const Teuchos::RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > DomainExporter,
146  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
147  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
149 
151  virtual ~TpetraBlockCrsMatrix();
152 
153 
154 
155 
157 
158 
160  void insertGlobalValues(GlobalOrdinal globalRow,
162  const ArrayView< const Scalar > &vals);
163 
164 
166  void insertLocalValues(LocalOrdinal localRow,
167  const ArrayView< const LocalOrdinal > &cols,
168  const ArrayView< const Scalar > &vals);
169 
171  void replaceGlobalValues(GlobalOrdinal globalRow,
173  const ArrayView< const Scalar > &vals);
174 
175 
177  void replaceLocalValues (LocalOrdinal localRow,
178  const ArrayView<const LocalOrdinal> &cols,
179  const ArrayView<const Scalar> &vals);
180 
182  void setAllToScalar(const Scalar &alpha);
183 
185  void scale(const Scalar &alpha);
186 
188  //** \warning This is an expert-only routine and should not be called from user code. (not implemented)
189  void allocateAllValues(size_t numNonZeros,ArrayRCP<size_t> & rowptr,
190  ArrayRCP<LocalOrdinal> & colind,
191  ArrayRCP<Scalar> & values);
192 
194  void setAllValues(const ArrayRCP<size_t> & rowptr,
195  const ArrayRCP<LocalOrdinal> & colind,
196  const ArrayRCP<Scalar> & values);
197 
199  void getAllValues(ArrayRCP<const size_t>& rowptr,
201  ArrayRCP<const Scalar>& values) const;
202 
204  void getAllValues(ArrayRCP<Scalar>& values);
205 
207 
209  void resumeFill(const RCP< ParameterList > &params=null);
210 
212  void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap,
213  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null);
214 
216  void fillComplete(const RCP< ParameterList > &params=null);
217 
218 
222 
225  const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & rangeMap,
226  const RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > &importer=Teuchos::null,
227  const RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > &exporter=Teuchos::null,
228  const RCP<ParameterList> &params=Teuchos::null);
229 
230 
232 
235 
238 
241 
244 
247 
249  size_t getNodeNumRows() const;
250 
252  size_t getNodeNumCols() const;
253 
256 
258  size_t getNodeNumEntries() const;
259 
261  size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const;
262 
264  size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const;
265 
267  size_t getGlobalMaxNumRowEntries() const;
268 
270  size_t getNodeMaxNumRowEntries() const;
271 
273  bool isLocallyIndexed() const;
274 
276  bool isGloballyIndexed() const;
277 
279  bool isFillComplete() const;
280 
282  bool isFillActive() const;
283 
286 
288  bool supportsRowViews() const;
289 
291  void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const;
292 
294  void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const;
295 
297  void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const;
298 
300  void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const;
301 
303  virtual bool haveGlobalConstants() const;
304 
305 
307 
308 
311 
313  void apply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues, const RCP<Import<LocalOrdinal, GlobalOrdinal, Node> >& regionInterfaceImporter, const Teuchos::ArrayRCP<LocalOrdinal>& regionInterfaceLIDs) const;
314 
317 
320 
321 
323 
324 
326  std::string description() const;
327 
328 
331 
332 
334  void setObjectLabel( const std::string &objectLabel );
335 
336 
337 
338 
341 
342 
345  const Teuchos::ArrayView<const size_t> &offsets) const;
346 
347 
349  void getLocalDiagOffsets(Teuchos::ArrayRCP<size_t> &offsets) const;
350 
351 
353 
356 
359 
361 
364 
368 
372 
376 
380 
382 
383 
384 
386 
388  bool hasMatrix() const;
389 
391  TpetraBlockCrsMatrix(const Teuchos::RCP<Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &mtx);
392 
395 
398 
399 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
400 #ifdef HAVE_XPETRA_TPETRA
401  //using local_matrix_type = typename Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_type;
403 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
404  local_matrix_type getLocalMatrix () const {
405  return getLocalMatrixDevice();
406  }
407 #endif
408  local_matrix_type getLocalMatrixDevice () const;
409  typename local_matrix_type::HostMirror getLocalMatrixHost () const;
410 
411  void setAllValues (const typename local_matrix_type::row_map_type& ptr,
412  const typename local_matrix_type::StaticCrsGraphType::entries_type::non_const_type& ind,
413  const typename local_matrix_type::values_type& val);
414 #endif // HAVE_XPETRA_TPETRA
415 #endif // HAVE_XPETRA_KOKKOS_REFACTOR
416 
421  using STS = Teuchos::ScalarTraits<Scalar>;
422  R.update(STS::one(),B,STS::zero());
423  this->apply (X, R, Teuchos::NO_TRANS, -STS::one(), STS::one());
424  }
425 
426 
427 
428 
429  private:
430 
432 
433  }; // TpetraBlockCrsMatrix class
434 
435 } // Xpetra namespace
436 
437 
438 #endif // XPETRA_TPETRABLOCKCRSMATRIX_DECL_HPP
439 
440 
static const EVerbosityLevel verbLevel_default
virtual void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)=0
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs (not implemented)
TpetraBlockCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying fixed number of entries for each row (not implemented)
void setObjectLabel(const std::string &objectLabel)
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalarThis.
std::string description() const
A simple one-line description of this object.
virtual bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
size_t getNodeMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
RCP< Tpetra::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > mtx_
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 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...
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
bool supportsRowViews() const
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class.
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
global_size_t getGlobalNumRows() const
Number of global elements in the row map of this matrix.
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale operator with given vector values.
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)
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 hasMatrix() const
Does this have an underlying matrix.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
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)
void replaceDiag(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag)
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale operator with given vector values.
RCP< const Tpetra::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_BlockCrsMatrix() const
Get the underlying Tpetra matrix.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this (not implemented)
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this matrix.
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs (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.
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this operator. This will be null until fillComplete() i...
TpetraBlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraBlockCrsMatrixClass
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 getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
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.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
void resumeFill(const RCP< ParameterList > &params=null)
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this matrix.
global_size_t getGlobalNumCols() const
Number of global columns in the matrix.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries in the (locally owned) global row.
bool isFillComplete() const
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called.
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const
Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage.
void expertStaticFillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &importer=Teuchos::null, const RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > &exporter=Teuchos::null, const RCP< ParameterList > &params=Teuchos::null)
Expert static fill complete.
void 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 (not implemented)
bool isFillActive() const
Returns true if the matrix is in edit mode.
void residual(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &R) const
Compute a residual R = B - (*this) * X.
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
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.
RCP< Tpetra::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_BlockCrsMatrixNonConst() const
Get the underlying Tpetra matrix.
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 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.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)
TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraVectorClass
Xpetra namespace
size_t global_size_t
Global size_t object.
CombineMode
Xpetra::Combine Mode enumerable type.