Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tpetra_BlockCrsMatrix_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Tpetra: Templated Linear Algebra Services Package
5 // Copyright (2008) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 // @HEADER
41 
42 #ifndef TPETRA_BLOCKCRSMATRIX_DECL_HPP
43 #define TPETRA_BLOCKCRSMATRIX_DECL_HPP
44 
47 
48 #include "Tpetra_CrsGraph.hpp"
49 #include "Tpetra_RowMatrix.hpp"
50 #include "Tpetra_BlockMultiVector_decl.hpp"
52 
53 namespace Tpetra {
54 
132 template<class Scalar,
133  class LO,
134  class GO,
135  class Node>
137  virtual public ::Tpetra::RowMatrix<Scalar, LO, GO, Node>,
138  virtual public ::Tpetra::DistObject<char, LO, GO, Node>
139 {
140 private:
143  using STS = Teuchos::ScalarTraits<Scalar>;
144 
145 protected:
147  typedef char packet_type;
148 
149 public:
151 
152 
154  using scalar_type = Scalar;
155 
163 
165  typedef LO local_ordinal_type;
172  typedef Node node_type;
173 
175  typedef typename Node::device_type device_type;
177  typedef typename device_type::execution_space execution_space;
179  typedef typename device_type::memory_space memory_space;
180 
182  typedef ::Tpetra::Map<LO, GO, node_type> map_type;
184  typedef ::Tpetra::MultiVector<Scalar, LO, GO, node_type> mv_type;
186  typedef ::Tpetra::CrsGraph<LO, GO, node_type> crs_graph_type;
187 
189  typedef Kokkos::View<impl_scalar_type**,
190  Kokkos::LayoutRight,
191  device_type,
192  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
195  typedef Kokkos::View<const impl_scalar_type**,
196  Kokkos::LayoutRight,
197  device_type,
198  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
201  typedef Kokkos::View<impl_scalar_type*,
202  Kokkos::LayoutRight,
203  device_type,
204  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
207  typedef Kokkos::View<const impl_scalar_type*,
208  Kokkos::LayoutRight,
209  device_type,
210  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
212 
214 
216 
218  BlockCrsMatrix ();
219 
229  BlockCrsMatrix (const crs_graph_type& graph, const LO blockSize);
230 
238  BlockCrsMatrix (const crs_graph_type& graph,
239  const map_type& domainPointMap,
240  const map_type& rangePointMap,
241  const LO blockSize);
242 
244  virtual ~BlockCrsMatrix () {}
245 
247 
249 
251  Teuchos::RCP<const map_type> getDomainMap () const;
252 
254  Teuchos::RCP<const map_type> getRangeMap () const;
255 
257  Teuchos::RCP<const map_type> getRowMap () const;
258 
260  Teuchos::RCP<const map_type> getColMap () const;
261 
264 
266  size_t getNodeNumRows() const;
267 
268  size_t getNodeMaxNumRowEntries() const;
269 
279  void
280  apply (const mv_type& X,
281  mv_type& Y,
282  Teuchos::ETransp mode = Teuchos::NO_TRANS,
283  Scalar alpha = Teuchos::ScalarTraits<Scalar>::one (),
284  Scalar beta = Teuchos::ScalarTraits<Scalar>::zero ()) const;
285 
288  bool hasTransposeApply () const {
289  // FIXME (mfh 04 May 2014) Transpose and conjugate transpose modes
290  // are not implemented yet. Fill in applyBlockTrans() to fix this.
291  return false;
292  }
293 
295  void setAllToScalar (const Scalar& alpha);
296 
298 
300 
302  std::string description () const;
303 
327  void
328  describe (Teuchos::FancyOStream& out,
329  const Teuchos::EVerbosityLevel verbLevel) const;
330 
332 
334 
336  LO getBlockSize () const { return blockSize_; }
337 
339  virtual Teuchos::RCP<const ::Tpetra::RowGraph<LO,GO,Node> > getGraph () const;
340 
341  const crs_graph_type & getCrsGraph () const { return graph_; }
342 
347  void
348  applyBlock (const BlockMultiVector<Scalar, LO, GO, Node>& X,
349  BlockMultiVector<Scalar, LO, GO, Node>& Y,
350  Teuchos::ETransp mode = Teuchos::NO_TRANS,
351  const Scalar alpha = Teuchos::ScalarTraits<Scalar>::one (),
352  const Scalar beta = Teuchos::ScalarTraits<Scalar>::zero ());
353 
357  void
358  gaussSeidelCopy (MultiVector<Scalar,LO,GO,Node> &X,
359  const ::Tpetra::MultiVector<Scalar,LO,GO,Node> &B,
360  const ::Tpetra::MultiVector<Scalar,LO,GO,Node> &D,
361  const Scalar& dampingFactor,
362  const ESweepDirection direction,
363  const int numSweeps,
364  const bool zeroInitialGuess) const;
365 
369  void
371  const ::Tpetra::MultiVector<Scalar,LO,GO,Node>& B,
372  const ::Tpetra::MultiVector<Scalar,LO,GO,Node>& D,
373  const Teuchos::ArrayView<LO>& rowIndices,
374  const Scalar& dampingFactor,
375  const ESweepDirection direction,
376  const int numSweeps,
377  const bool zeroInitialGuess) const;
378 
398  void
399  localGaussSeidel (const BlockMultiVector<Scalar, LO, GO, Node>& Residual,
400  BlockMultiVector<Scalar, LO, GO, Node>& Solution,
401  const Kokkos::View<impl_scalar_type***, device_type,
402  Kokkos::MemoryUnmanaged>& D_inv,
403  const Scalar& omega,
404  const ESweepDirection direction) const;
405 
432  LO
433  replaceLocalValues (const LO localRowInd,
434  const LO colInds[],
435  const Scalar vals[],
436  const LO numColInds) const;
437 
464  LO
465  sumIntoLocalValues (const LO localRowInd,
466  const LO colInds[],
467  const Scalar vals[],
468  const LO numColInds) const;
469 
499  LO
500  getLocalRowView (const LO localRowInd,
501  const LO*& colInds,
502  Scalar*& vals,
503  LO& numInds) const;
504 
506  void
507  getLocalRowView (LO LocalRow,
508  Teuchos::ArrayView<const LO> &indices,
509  Teuchos::ArrayView<const Scalar> &values) const;
510 
512  void
513  getLocalRowCopy (LO LocalRow,
514  const Teuchos::ArrayView<LO> &Indices,
515  const Teuchos::ArrayView<Scalar> &Values,
516  size_t &NumEntries) const;
517 
519  getLocalBlock (const LO localRowInd, const LO localColInd) const;
520 
544  LO
545  getLocalRowOffsets (const LO localRowInd,
546  ptrdiff_t offsets[],
547  const LO colInds[],
548  const LO numColInds) const;
549 
555  LO
556  replaceLocalValuesByOffsets (const LO localRowInd,
557  const ptrdiff_t offsets[],
558  const Scalar vals[],
559  const LO numOffsets) const;
560 
566  LO
567  sumIntoLocalValuesByOffsets (const LO localRowInd,
568  const ptrdiff_t offsets[],
569  const Scalar vals[],
570  const LO numOffsets) const;
571 
578  size_t getNumEntriesInLocalRow (const LO localRowInd) const;
579 
596  bool localError () const {
597  return *localError_;
598  }
599 
614  std::string errorMessages () const {
615  return (*errs_).is_null () ? std::string ("") : (*errs_)->str ();
616  }
617 
649  void
650  getLocalDiagOffsets (const Kokkos::View<size_t*, device_type,
651  Kokkos::MemoryUnmanaged>& offsets) const;
652 
653 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
654  void TPETRA_DEPRECATED
660  getLocalDiagOffsets (Teuchos::ArrayRCP<size_t>& offsets) const;
661 #endif // TPETRA_ENABLE_DEPRECATED_CODE
662 
676  void
677  getLocalDiagCopy (const Kokkos::View<impl_scalar_type***, device_type,
678  Kokkos::MemoryUnmanaged>& diag,
679  const Kokkos::View<const size_t*, device_type,
680  Kokkos::MemoryUnmanaged>& offsets) const;
681 
695  void
696  getLocalDiagCopy (const Kokkos::View<impl_scalar_type***, device_type,
697  Kokkos::MemoryUnmanaged>& diag,
698  const Teuchos::ArrayView<const size_t>& offsets) const;
699 
700 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
701  void TPETRA_DEPRECATED
716  const Teuchos::ArrayView<const size_t>& offsets) const;
717 #endif // TPETRA_ENABLE_DEPRECATED_CODE
718 
719 protected:
721  LO
722  absMaxLocalValues (const LO localRowInd,
723  const LO colInds[],
724  const Scalar vals[],
725  const LO numColInds) const;
726 
728  LO
729  absMaxLocalValuesByOffsets (const LO localRowInd,
730  const ptrdiff_t offsets[],
731  const Scalar vals[],
732  const LO numOffsets) const;
733 
739 
740 
745  using buffer_device_type = typename DistObject<Scalar, LO, GO,
747 
748  virtual bool checkSizes (const ::Tpetra::SrcDistObject& source);
749 
750  virtual void
751 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
752  copyAndPermuteNew
753 #else // TPETRA_ENABLE_DEPRECATED_CODE
754  copyAndPermute
755 #endif // TPETRA_ENABLE_DEPRECATED_CODE
756  (const SrcDistObject& sourceObj,
757  const size_t numSameIDs,
758  const Kokkos::DualView<const local_ordinal_type*,
759  buffer_device_type>& permuteToLIDs,
760  const Kokkos::DualView<const local_ordinal_type*,
761  buffer_device_type>& permuteFromLIDs);
762 
763  virtual void
764 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
765  packAndPrepareNew
766 #else // TPETRA_ENABLE_DEPRECATED_CODE
767  packAndPrepare
768 #endif // TPETRA_ENABLE_DEPRECATED_CODE
769  (const SrcDistObject& sourceObj,
770  const Kokkos::DualView<const local_ordinal_type*,
771  buffer_device_type>& exportLIDs,
772  Kokkos::DualView<packet_type*,
773  buffer_device_type>& exports,
774  Kokkos::DualView<size_t*,
775  buffer_device_type> numPacketsPerLID,
776  size_t& constantNumPackets,
777  Distributor& /* distor */);
778 
779  virtual void
780 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
781  unpackAndCombineNew
782 #else // TPETRA_ENABLE_DEPRECATED_CODE
783  unpackAndCombine
784 #endif // TPETRA_ENABLE_DEPRECATED_CODE
785  (const Kokkos::DualView<const local_ordinal_type*,
786  buffer_device_type>& importLIDs,
787  Kokkos::DualView<packet_type*,
788  buffer_device_type> imports,
789  Kokkos::DualView<size_t*,
790  buffer_device_type> numPacketsPerLID,
791  const size_t constantNumPackets,
792  Distributor& /* distor */,
793  const CombineMode combineMode);
795 
796 private:
798  crs_graph_type graph_;
799  Teuchos::RCP<crs_graph_type> graphRCP_;
808  map_type rowMeshMap_;
815  map_type domainPointMap_;
822  map_type rangePointMap_;
824  LO blockSize_;
825 
839  typename crs_graph_type::local_graph_type::row_map_type::HostMirror ptrHost_;
840 
846  typename crs_graph_type::local_graph_type::entries_type::HostMirror indHost_;
847 
853  typename Kokkos::DualView<impl_scalar_type*, device_type> val_;
854 
876  Teuchos::RCP<Teuchos::RCP<BMV> > X_colMap_;
880  Teuchos::RCP<Teuchos::RCP<BMV> > Y_rowMap_;
881 
889  Teuchos::RCP<Teuchos::RCP<typename crs_graph_type::import_type> > pointImporter_;
890 
892  LO offsetPerBlock_;
893 
905  Teuchos::RCP<bool> localError_;
906 
914  Teuchos::RCP<Teuchos::RCP<std::ostringstream> > errs_;
915 
917  std::ostream& markLocalErrorAndGetStream ();
918 
919  // //! Clear the local error state and stream.
920  // void clearLocalErrorStateAndStream ();
921 
922  template<class Device>
923  struct is_cuda {
924 #if defined(KOKKOS_ENABLE_CUDA)
925  // CudaHostPinnedSpace::execution_space ==
926  // HostSpace::execution_space. That's OK; it's host memory, that
927  // just happens to be Cuda accessible. But what if somebody gives
928  // us Device<Cuda, CudaHostPinnedSpace>? It looks like they mean
929  // to run on device then, so we should sync to device.
930  static constexpr bool value =
931  std::is_same<typename Device::execution_space, Kokkos::Cuda>::value;
932 #else
933  static constexpr bool value = false;
934 #endif // defined(KOKKOS_ENABLE_CUDA)
935  };
936 
937 public:
939 
940 
942  inline void modify_host()
943  {
944  val_.modify_host();
945  }
946 
948  inline void modify_device()
949  {
950  val_.modify_device();
951  }
952 
954  template<class MemorySpace>
955  void modify ()
956  {
957  if (is_cuda<MemorySpace>::value) {
958  this->modify_device ();
959  }
960  else {
961  this->modify_host ();
962  }
963  }
964 
966  inline bool need_sync_host() const
967  {
968  return val_.need_sync_host();
969  }
970 
972  inline bool need_sync_device() const
973  {
974  return val_.need_sync_device();
975  }
976 
978  template<class MemorySpace>
979  bool need_sync () const
980  {
981  if (is_cuda<MemorySpace>::value) {
982  return this->need_sync_device ();
983  }
984  else {
985  return this->need_sync_host ();
986  }
987  }
988 
990  inline void sync_host()
991  {
992  val_.sync_host();
993  }
994 
996  inline void sync_device()
997  {
998  val_.sync_device();
999  }
1000 
1002  template<class MemorySpace>
1003  void sync ()
1004  {
1005  if (is_cuda<MemorySpace>::value) {
1006  this->sync_device ();
1007  }
1008  else {
1009  this->sync_host ();
1010  }
1011  }
1012 
1013  // \brief Get the host view of the matrix's values
1014  typename Kokkos::DualView<impl_scalar_type*, device_type>::t_host getValuesHost () const {
1015  return val_.view_host();
1016  }
1017 
1018  // \brief Get the device view of the matrix's values
1019  typename Kokkos::DualView<impl_scalar_type*, device_type>::t_dev getValuesDevice () const {
1020  return val_.view_device();
1021  }
1022 
1040  template<class MemorySpace>
1041  typename std::conditional<is_cuda<MemorySpace>::value,
1042  typename Kokkos::DualView<impl_scalar_type*, device_type>::t_dev,
1043  typename Kokkos::DualView<impl_scalar_type*, device_type>::t_host>::type
1045  {
1046  // Unlike std::conditional, if_c has a select method.
1047  return Kokkos::Impl::if_c<
1048  is_cuda<MemorySpace>::value,
1049  typename Kokkos::DualView<impl_scalar_type*, device_type>::t_dev,
1050  typename Kokkos::DualView<impl_scalar_type*, device_type>::t_host
1051  >::select (this->getValuesDevice (), this->getValuesHost ());
1052  }
1053 
1055 
1056 private:
1057 
1067  void
1068  applyBlockTrans (const BlockMultiVector<Scalar, LO, GO, Node>& X,
1070  const Teuchos::ETransp mode,
1071  const Scalar alpha,
1072  const Scalar beta);
1073 
1081  void
1082  applyBlockNoTrans (const BlockMultiVector<Scalar, LO, GO, Node>& X,
1084  const Scalar alpha,
1085  const Scalar beta);
1086 
1094  void
1095  localApplyBlockNoTrans (const BlockMultiVector<Scalar, LO, GO, Node>& X,
1097  const Scalar alpha,
1098  const Scalar beta);
1099 
1139  LO
1140  findRelOffsetOfColumnIndex (const LO localRowIndex,
1141  const LO colIndexToFind,
1142  const LO hint = 0) const;
1143 
1146  LO offsetPerBlock () const;
1147 
1149  getConstLocalBlockFromInput (const impl_scalar_type* val, const size_t pointOffset) const;
1150 
1152  getNonConstLocalBlockFromInput (impl_scalar_type* val, const size_t pointOffset) const;
1153 
1155  getConstLocalBlockFromAbsOffset (const size_t absBlockOffset) const;
1156 
1158  getNonConstLocalBlockFromAbsOffset (const size_t absBlockOffset) const;
1159 
1164  getConstLocalBlockFromRelOffset (const LO lclMeshRow,
1165  const size_t relMeshOffset) const;
1166 
1167 public:
1169  virtual Teuchos::RCP<const Teuchos::Comm<int> > getComm() const;
1170 
1171 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
1172  virtual TPETRA_DEPRECATED Teuchos::RCP<Node> getNode() const;
1174 #endif // TPETRA_ENABLE_DEPRECATED_CODE
1175 
1177  virtual global_size_t getGlobalNumCols() const;
1178 
1179  virtual size_t getNodeNumCols() const;
1180 
1181  virtual GO getIndexBase() const;
1182 
1184  virtual global_size_t getGlobalNumEntries() const;
1185 
1187  virtual size_t getNodeNumEntries() const;
1188 
1198  virtual size_t getNumEntriesInGlobalRow (GO globalRow) const;
1199 
1202  virtual size_t getGlobalMaxNumRowEntries () const;
1203 
1205  virtual bool hasColMap () const;
1206 
1216  virtual bool isLocallyIndexed () const;
1217 
1227  virtual bool isGloballyIndexed () const;
1228 
1230  virtual bool isFillComplete () const;
1231 
1233  virtual bool supportsRowViews () const;
1234 
1235 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
1236  virtual global_size_t TPETRA_DEPRECATED getGlobalNumDiags() const;
1244 
1252  virtual size_t TPETRA_DEPRECATED getNodeNumDiags() const;
1253 
1264  virtual bool TPETRA_DEPRECATED isLowerTriangular () const;
1265 
1276  virtual bool TPETRA_DEPRECATED isUpperTriangular () const;
1277 #endif // TPETRA_ENABLE_DEPRECATED_CODE
1278 
1280 
1282 
1303  virtual void
1304  getGlobalRowCopy (GO GlobalRow,
1305  const Teuchos::ArrayView<GO> &Indices,
1306  const Teuchos::ArrayView<Scalar> &Values,
1307  size_t& NumEntries) const;
1308 
1333  virtual void
1334  getGlobalRowView (GO GlobalRow,
1335  Teuchos::ArrayView<const GO>& indices,
1336  Teuchos::ArrayView<const Scalar>& values) const;
1337 
1349  virtual void getLocalDiagCopy (::Tpetra::Vector<Scalar,LO,GO,Node>& diag) const;
1350 
1352 
1354 
1360  virtual void leftScale (const ::Tpetra::Vector<Scalar, LO, GO, Node>& x);
1361 
1367  virtual void rightScale (const ::Tpetra::Vector<Scalar, LO, GO, Node>& x);
1368 
1377  virtual typename ::Tpetra::RowMatrix<Scalar, LO, GO, Node>::mag_type
1378  getFrobeniusNorm () const;
1380 };
1381 
1382 } // namespace Tpetra
1383 
1384 #endif // TPETRA_BLOCKCRSMATRIX_DECL_HPP
virtual size_t getNumEntriesInGlobalRow(GO globalRow) const
The current number of entries on the calling process in the specified global row. ...
::Tpetra::MultiVector< Scalar, LO, GO, node_type > mv_type
The implementation of MultiVector that this class uses.
typename BMV::impl_scalar_type impl_scalar_type
The implementation type of entries in the matrix.
void getLocalDiagCopy(const Kokkos::View< impl_scalar_type ***, device_type, Kokkos::MemoryUnmanaged > &diag, const Kokkos::View< const size_t *, device_type, Kokkos::MemoryUnmanaged > &offsets) const
Variant of getLocalDiagCopy() that uses precomputed offsets and puts diagonal blocks in a 3-D Kokkos:...
LO replaceLocalValues(const LO localRowInd, const LO colInds[], const Scalar vals[], const LO numColInds) const
Replace values at the given (mesh, i.e., block) column indices, in the given (mesh, i.e., block) row.
void gaussSeidelCopy(MultiVector< Scalar, LO, GO, Node > &X, const ::Tpetra::MultiVector< Scalar, LO, GO, Node > &B, const ::Tpetra::MultiVector< Scalar, LO, GO, Node > &D, const Scalar &dampingFactor, const ESweepDirection direction, const int numSweeps, const bool zeroInitialGuess) const
Version of gaussSeidel(), with fewer requirements on X.
LO replaceLocalValuesByOffsets(const LO localRowInd, const ptrdiff_t offsets[], const Scalar vals[], const LO numOffsets) const
Like replaceLocalValues, but avoids computing row offsets.
Sparse matrix whose entries are small dense square blocks, all of the same dimensions.
size_t getNodeNumRows() const
get the local number of block rows
void modify_device()
Mark the matrix&#39;s valueas as modified in device space.
typename DistObject< Scalar, LO, GO, Node >::buffer_device_type buffer_device_type
Kokkos::Device specialization for communication buffers.
virtual bool isGloballyIndexed() const
Whether matrix indices are globally indexed.
size_t getNumEntriesInLocalRow(const LO localRowInd) const
Return the number of entries in the given row on the calling process.
LO local_ordinal_type
The type of local indices.
virtual void getGlobalRowView(GO GlobalRow, Teuchos::ArrayView< const GO > &indices, Teuchos::ArrayView< const Scalar > &values) const
Get a constant, nonpersisting, globally indexed view of the given row of the matrix.
virtual typename::Tpetra::RowMatrix< Scalar, LO, GO, Node >::mag_type getFrobeniusNorm() const
The Frobenius norm of the matrix.
virtual global_size_t getGlobalNumCols() const
The global number of columns of this matrix.
::Tpetra::CrsGraph< LO, GO, node_type > crs_graph_type
The implementation of CrsGraph that this class uses.
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to alpha.
virtual bool isLocallyIndexed() const
Whether matrix indices are locally indexed.
bool localError() const
Whether this object had an error on the calling process.
Declaration of the Tpetra::CrsMatrix class.
void sync_device()
Sync the matrix&#39;s values to device space.
Scalar scalar_type
The type of entries in the matrix (that is, of each entry in each block).
virtual void leftScale(const ::Tpetra::Vector< Scalar, LO, GO, Node > &x)
Scale the RowMatrix on the left with the given Vector x.
virtual bool hasColMap() const
Whether this matrix has a well-defined column Map.
Kokkos::View< impl_scalar_type *, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > little_vec_type
The type used to access nonconst vector blocks.
virtual size_t getNodeNumEntries() const
The local number of stored (structurally nonzero) entries.
::Tpetra::Map< LO, GO, node_type > map_type
The implementation of Map that this class uses.
virtual size_t getNodeNumCols() const
The number of columns needed to apply the forward operator on this node.
virtual ~BlockCrsMatrix()
Destructor (declared virtual for memory safety).
MultiVector for multiple degrees of freedom per mesh point.
GO global_ordinal_type
The type of global indices.
size_t global_size_t
Global size_t object.
BlockCrsMatrix()
Default constructor: Makes an empty block matrix.
virtual bool isFillComplete() const
Whether fillComplete() has been called.
void apply(const mv_type &X, mv_type &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
For this matrix A, compute Y := beta * Y + alpha * Op(A) * X.
LO absMaxLocalValues(const LO localRowInd, const LO colInds[], const Scalar vals[], const LO numColInds) const
Like sumIntoLocalValues, but for the ABSMAX combine mode.
virtual size_t getGlobalMaxNumRowEntries() const
The maximum number of entries in any row over all processes in the matrix&#39;s communicator.
device_type::memory_space memory_space
The Kokkos memory space that this class uses.
size_t getNodeMaxNumRowEntries() const
Maximum number of entries in any row of the matrix, on this process.
std::string errorMessages() const
The current stream of error messages.
LO sumIntoLocalValuesByOffsets(const LO localRowInd, const ptrdiff_t offsets[], const Scalar vals[], const LO numOffsets) const
Like sumIntoLocalValues, but avoids computing row offsets.
virtual void rightScale(const ::Tpetra::Vector< Scalar, LO, GO, Node > &x)
Scale the RowMatrix on the right with the given Vector x.
ESweepDirection
Sweep direction for Gauss-Seidel or Successive Over-Relaxation (SOR).
Teuchos::RCP< const map_type > getRangeMap() const
Get the (point) range Map of this matrix.
bool need_sync_device() const
Whether the matrix&#39;s values need sync&#39;ing to device space.
Teuchos::RCP< const map_type > getRowMap() const
get the (mesh) map for the rows of this block matrix.
std::conditional< is_cuda< MemorySpace >::value, typename Kokkos::DualView< impl_scalar_type *, device_type >::t_dev, typename Kokkos::DualView< impl_scalar_type *, device_type >::t_host >::type getValues()
Get the host or device View of the matrix&#39;s values (val_).
void modify_host()
Mark the matrix&#39;s valueas as modified in host space.
virtual global_size_t getGlobalNumEntries() const
The global number of stored (structurally nonzero) entries.
Sets up and executes a communication plan for a Tpetra DistObject.
void applyBlock(const BlockMultiVector< Scalar, LO, GO, Node > &X, BlockMultiVector< Scalar, LO, GO, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, const Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), const Scalar beta=Teuchos::ScalarTraits< Scalar >::zero())
Version of apply() that takes BlockMultiVector input and output.
CombineMode
Rule for combining data in an Import or Export.
global_size_t getGlobalNumRows() const
get the global number of block rows
CrsGraphType::global_ordinal_type getGlobalNumDiags(const CrsGraphType &G)
Number of populated diagonal entries in the given sparse graph, over all processes in the graph&#39;s (MP...
virtual void getGlobalRowCopy(GO GlobalRow, const Teuchos::ArrayView< GO > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const
Get a copy of the given global row&#39;s entries.
LO getLocalRowOffsets(const LO localRowInd, ptrdiff_t offsets[], const LO colInds[], const LO numColInds) const
Get relative offsets corresponding to the given rows, given by local row index.
void sync_host()
Sync the matrix&#39;s values to host space.
Abstract base class for objects that can be the source of an Import or Export operation.
device_type::execution_space execution_space
The Kokkos execution space that this class uses.
Kokkos::View< const impl_scalar_type *, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > const_little_vec_type
The type used to access const vector blocks.
Kokkos::View< impl_scalar_type **, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > little_block_type
The type used to access nonconst matrix blocks.
virtual bool supportsRowViews() const
Whether this object implements getLocalRowView() and getGlobalRowView().
LO absMaxLocalValuesByOffsets(const LO localRowInd, const ptrdiff_t offsets[], const Scalar vals[], const LO numOffsets) const
Like sumIntoLocalValuesByOffsets, but for the ABSMAX combine mode.
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which this matrix is distributed.
LO getLocalRowView(const LO localRowInd, const LO *&colInds, Scalar *&vals, LO &numInds) const
Get a view of the (mesh, i.e., block) row, using local (mesh, i.e., block) indices.
void getLocalDiagOffsets(const Kokkos::View< size_t *, device_type, Kokkos::MemoryUnmanaged > &offsets) const
Get offsets of the diagonal entries in the matrix.
char packet_type
Implementation detail; tells.
void sync()
Sync the matrix&#39;s values to the given memory space.
void localGaussSeidel(const BlockMultiVector< Scalar, LO, GO, Node > &Residual, BlockMultiVector< Scalar, LO, GO, Node > &Solution, const Kokkos::View< impl_scalar_type ***, device_type, Kokkos::MemoryUnmanaged > &D_inv, const Scalar &omega, const ESweepDirection direction) const
Local Gauss-Seidel solve, given a factorized diagonal.
A read-only, row-oriented interface to a sparse matrix.
void modify()
Mark the matrix&#39;s values as modified in the given memory space.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
Print a description of this object to the given output stream.
A distributed dense vector.
LO getBlockSize() const
The number of degrees of freedom per mesh point.
bool need_sync() const
Whether the matrix&#39;s values need sync&#39;ing to the given memory space.
Teuchos::RCP< const map_type > getDomainMap() const
Get the (point) domain Map of this matrix.
bool hasTransposeApply() const
Whether it is valid to apply the transpose or conjugate transpose of this matrix. ...
void getLocalRowCopy(LO LocalRow, const Teuchos::ArrayView< LO > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const
Not implemented.
Kokkos::View< const impl_scalar_type **, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > const_little_block_type
The type used to access const matrix blocks.
Node::device_type device_type
The Kokkos::Device specialization that this class uses.
Teuchos::RCP< const map_type > getColMap() const
get the (mesh) map for the columns of this block matrix.
virtual Teuchos::RCP< const ::Tpetra::RowGraph< LO, GO, Node > > getGraph() const
Get the (mesh) graph.
Base class for distributed Tpetra objects that support data redistribution.
typename mv_type::impl_scalar_type impl_scalar_type
The implementation type of entries in the object.
void reorderedGaussSeidelCopy(::Tpetra::MultiVector< Scalar, LO, GO, Node > &X, const ::Tpetra::MultiVector< Scalar, LO, GO, Node > &B, const ::Tpetra::MultiVector< Scalar, LO, GO, Node > &D, const Teuchos::ArrayView< LO > &rowIndices, const Scalar &dampingFactor, const ESweepDirection direction, const int numSweeps, const bool zeroInitialGuess) const
Version of reorderedGaussSeidel(), with fewer requirements on X.
LO sumIntoLocalValues(const LO localRowInd, const LO colInds[], const Scalar vals[], const LO numColInds) const
Sum into values at the given (mesh, i.e., block) column indices, in the given (mesh, i.e., block) row.
std::string description() const
One-line description of this object.
bool need_sync_host() const
Whether the matrix&#39;s values need sync&#39;ing to host space.
virtual GO getIndexBase() const
The index base for global indices in this matrix.