All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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"
57 #include "Xpetra_TpetraMap.hpp"
58 #include "Xpetra_TpetraMultiVector.hpp"
59 #include "Xpetra_TpetraVector.hpp"
60 #include "Xpetra_TpetraCrsGraph.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,
85  const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null);
86 
87 
90  const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc,
92  const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null);
93 
94 
98  size_t maxNumEntriesPerRow,
100  const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null);
101 
102 
106  const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc,
107  ProfileType pftype=DynamicProfile,
108  const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null);
109 
110 
113  const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null);
114 
115 
118  const LocalOrdinal blockSize);
119 
120 
122  TpetraBlockCrsMatrix(const Teuchos::RCP<const Tpetra::BlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& sourceMatrix,
124  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
125  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
126  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
127 
128 
130  TpetraBlockCrsMatrix(const Teuchos::RCP<const Tpetra::BlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& sourceMatrix,
132  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
133  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
134  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
135 
136 
138  TpetraBlockCrsMatrix(const Teuchos::RCP<const Tpetra::BlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& sourceMatrix,
139  const Import<LocalOrdinal,GlobalOrdinal,Node> & RowImporter,
140  const Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > DomainImporter,
141  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
142  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
144 
145 
147  TpetraBlockCrsMatrix(const Teuchos::RCP<const Tpetra::BlockCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& sourceMatrix,
148  const Export<LocalOrdinal,GlobalOrdinal,Node> & RowExporter,
149  const Teuchos::RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > DomainExporter,
150  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
151  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
153 
155  virtual ~TpetraBlockCrsMatrix();
156 
157 
158 
159 
161 
162 
164  void insertGlobalValues(GlobalOrdinal globalRow,
166  const ArrayView< const Scalar > &vals);
167 
168 
170  void insertLocalValues(LocalOrdinal localRow,
171  const ArrayView< const LocalOrdinal > &cols,
172  const ArrayView< const Scalar > &vals);
173 
175  void replaceGlobalValues(GlobalOrdinal globalRow,
177  const ArrayView< const Scalar > &vals);
178 
179 
181  void replaceLocalValues (LocalOrdinal localRow,
182  const ArrayView<const LocalOrdinal> &cols,
183  const ArrayView<const Scalar> &vals);
184 
186  void setAllToScalar(const Scalar &alpha);
187 
189  void scale(const Scalar &alpha);
190 
192  //** \warning This is an expert-only routine and should not be called from user code. (not implemented)
193  void allocateAllValues(size_t numNonZeros,ArrayRCP<size_t> & rowptr,
194  ArrayRCP<LocalOrdinal> & colind,
195  ArrayRCP<Scalar> & values);
196 
198  void setAllValues(const ArrayRCP<size_t> & rowptr,
199  const ArrayRCP<LocalOrdinal> & colind,
200  const ArrayRCP<Scalar> & values);
201 
203  void getAllValues(ArrayRCP<const size_t>& rowptr,
205  ArrayRCP<const Scalar>& values) const;
206 
207 
209 
211  void resumeFill(const RCP< ParameterList > &params=null);
212 
214  void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap,
215  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null);
216 
218  void fillComplete(const RCP< ParameterList > &params=null);
219 
220 
224 
227  const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & rangeMap,
228  const RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > &importer=Teuchos::null,
229  const RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > &exporter=Teuchos::null,
230  const RCP<ParameterList> &params=Teuchos::null);
231 
232 
234 
237 
240 
243 
246 
249 
251  size_t getNodeNumRows() const;
252 
254  size_t getNodeNumCols() const;
255 
258 
260  size_t getNodeNumEntries() const;
261 
263  size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const;
264 
266  size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const;
267 
269  size_t getGlobalMaxNumRowEntries() const;
270 
272  size_t getNodeMaxNumRowEntries() const;
273 
275  bool isLocallyIndexed() const;
276 
278  bool isGloballyIndexed() const;
279 
281  bool isFillComplete() const;
282 
284  bool isFillActive() const;
285 
288 
290  bool supportsRowViews() const;
291 
293  void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const;
294 
296  void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const;
297 
299  void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const;
300 
302  void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const;
303 
305  virtual bool haveGlobalConstants() const;
306 
307 
309 
310 
313 
316 
319 
320 
322 
323 
325  std::string description() const;
326 
327 
330 
331 
333  void setObjectLabel( const std::string &objectLabel );
334 
335 
336 #ifdef XPETRA_ENABLE_DEPRECATED_CODE
337  // This probably never compiled but also never got called...
339  // We're leaving this in the decl file for now as part of the ETI work
340  // just to maintain status-quo but it'll go away anyhow with deprecations
341  // soon.
343 #endif // XPETRA_ENABLE_DEPRECATED_CODE
344 
345 
348 
349 
352  const Teuchos::ArrayView<const size_t> &offsets) const;
353 
354 
356  void getLocalDiagOffsets(Teuchos::ArrayRCP<size_t> &offsets) const;
357 
358 
360 
363 
366 
368 
371 
375 
379 
383 
387 
389 
390 
391 #ifdef XPETRA_ENABLE_DEPRECATED_CODE
392  template<class Node2>
394  XPETRA_DEPRECATED
395  clone(const RCP<Node2> &node2) const
396  {
399  }
400 #endif // XPETRA_ENABLE_DEPRECATED_CODE
401 
403 
405  bool hasMatrix() const;
406 
408  TpetraBlockCrsMatrix(const Teuchos::RCP<Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &mtx);
409 
412 
415 
416 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
417 #ifdef HAVE_XPETRA_TPETRA
418  //using local_matrix_type = typename Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_type;
420 
421  local_matrix_type getLocalMatrix () const;
422 
423  void setAllValues (const typename local_matrix_type::row_map_type& ptr,
424  const typename local_matrix_type::StaticCrsGraphType::entries_type::non_const_type& ind,
425  const typename local_matrix_type::values_type& val);
426 #endif // HAVE_XPETRA_TPETRA
427 #endif // HAVE_XPETRA_KOKKOS_REFACTOR
428 
429  private:
430 
432 
433  }; // TpetraBlockCrsMatrix class
434 
435 } // Xpetra namespace
436 
437 
438 #endif // XPETRA_TPETRABLOCKCRSMATRIX_DECL_HPP
439 
440 
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)
RCP< Map< LocalOrdinal, GlobalOrdinal, Node2 > > clone(const Map< LocalOrdinal, GlobalOrdinal, Node1 > &map, const RCP< Node2 > &node2)
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...
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.
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 (not implemented)
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs (not implemented)
RCP< Tpetra::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_BlockCrsMatrixNonConst() const
Get the underlying Tpetra matrix.
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale operator with given vector values.
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.
TpetraBlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraBlockCrsMatrixClass
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 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 setObjectLabel(const std::string &objectLabel)
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.
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this matrix. ...
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< const Tpetra::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_BlockCrsMatrix() 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.
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)
RCP< Tpetra::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > mtx_
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.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries in the (locally owned) global row.
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.
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 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)
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this operator. This will be null until fillComplete() i...
void replaceDiag(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag)
Replace the diagonal entries of the matrix.
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)
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
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 doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this 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.
void resumeFill(const RCP< ParameterList > &params=null)