Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tpetra_Experimental_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_EXPERIMENTAL_BLOCKCRSMATRIX_DECL_HPP
43 #define TPETRA_EXPERIMENTAL_BLOCKCRSMATRIX_DECL_HPP
44 
47 
48 #include "Tpetra_CrsGraph.hpp"
49 #include "Tpetra_RowMatrix.hpp"
50 #include "Tpetra_Experimental_BlockMultiVector_decl.hpp"
52 
53 namespace Tpetra {
54 namespace Experimental {
55 
133 template<class Scalar,
134  class LO,
135  class GO,
136  class Node>
138  virtual public ::Tpetra::RowMatrix<Scalar, LO, GO, Node>,
139  virtual public ::Tpetra::DistObject<char, LO, GO, Node>
140 {
141 private:
144  using STS = Teuchos::ScalarTraits<Scalar>;
145 
146 protected:
148  typedef char packet_type;
149 
150 public:
152 
153 
155  using scalar_type = Scalar;
156 
164 
166  typedef LO local_ordinal_type;
173  typedef Node node_type;
174 
176  typedef typename Node::device_type device_type;
178  typedef typename device_type::execution_space execution_space;
180  typedef typename device_type::memory_space memory_space;
181 
183  typedef ::Tpetra::Map<LO, GO, node_type> map_type;
185  typedef ::Tpetra::MultiVector<Scalar, LO, GO, node_type> mv_type;
187  typedef ::Tpetra::CrsGraph<LO, GO, node_type> crs_graph_type;
188 
190  typedef Kokkos::View<impl_scalar_type**,
191  Kokkos::LayoutRight,
192  device_type,
193  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
196  typedef Kokkos::View<const impl_scalar_type**,
197  Kokkos::LayoutRight,
198  device_type,
199  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
202  typedef Kokkos::View<impl_scalar_type*,
203  Kokkos::LayoutRight,
204  device_type,
205  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
208  typedef Kokkos::View<const impl_scalar_type*,
209  Kokkos::LayoutRight,
210  device_type,
211  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
213 
215 
217 
219  BlockCrsMatrix ();
220 
230  BlockCrsMatrix (const crs_graph_type& graph, const LO blockSize);
231 
239  BlockCrsMatrix (const crs_graph_type& graph,
240  const map_type& domainPointMap,
241  const map_type& rangePointMap,
242  const LO blockSize);
243 
245  virtual ~BlockCrsMatrix () {}
246 
248 
250 
252  Teuchos::RCP<const map_type> getDomainMap () const;
253 
255  Teuchos::RCP<const map_type> getRangeMap () const;
256 
258  Teuchos::RCP<const map_type> getRowMap () const;
259 
261  Teuchos::RCP<const map_type> getColMap () const;
262 
265 
267  size_t getNodeNumRows() const;
268 
269  size_t getNodeMaxNumRowEntries() const;
270 
280  void
281  apply (const mv_type& X,
282  mv_type& Y,
283  Teuchos::ETransp mode = Teuchos::NO_TRANS,
284  Scalar alpha = Teuchos::ScalarTraits<Scalar>::one (),
285  Scalar beta = Teuchos::ScalarTraits<Scalar>::zero ()) const;
286 
289  bool hasTransposeApply () const {
290  // FIXME (mfh 04 May 2014) Transpose and conjugate transpose modes
291  // are not implemented yet. Fill in applyBlockTrans() to fix this.
292  return false;
293  }
294 
296  void setAllToScalar (const Scalar& alpha);
297 
299 
301 
303  std::string description () const;
304 
328  void
329  describe (Teuchos::FancyOStream& out,
330  const Teuchos::EVerbosityLevel verbLevel) const;
331 
333 
335 
337  LO getBlockSize () const { return blockSize_; }
338 
340  virtual Teuchos::RCP<const ::Tpetra::RowGraph<LO,GO,Node> > getGraph () const;
341 
342  const crs_graph_type & getCrsGraph () const { return graph_; }
343 
348  void
349  applyBlock (const BlockMultiVector<Scalar, LO, GO, Node>& X,
350  BlockMultiVector<Scalar, LO, GO, Node>& Y,
351  Teuchos::ETransp mode = Teuchos::NO_TRANS,
352  const Scalar alpha = Teuchos::ScalarTraits<Scalar>::one (),
353  const Scalar beta = Teuchos::ScalarTraits<Scalar>::zero ());
354 
358  void
360  const ::Tpetra::MultiVector<Scalar,LO,GO,Node> &B,
361  const ::Tpetra::MultiVector<Scalar,LO,GO,Node> &D,
362  const Scalar& dampingFactor,
363  const ESweepDirection direction,
364  const int numSweeps,
365  const bool zeroInitialGuess) const;
366 
370  void
372  const ::Tpetra::MultiVector<Scalar,LO,GO,Node>& B,
373  const ::Tpetra::MultiVector<Scalar,LO,GO,Node>& D,
374  const Teuchos::ArrayView<LO>& rowIndices,
375  const Scalar& dampingFactor,
376  const ESweepDirection direction,
377  const int numSweeps,
378  const bool zeroInitialGuess) const;
379 
399  void
400  localGaussSeidel (const BlockMultiVector<Scalar, LO, GO, Node>& Residual,
401  BlockMultiVector<Scalar, LO, GO, Node>& Solution,
402  const Kokkos::View<impl_scalar_type***, device_type,
403  Kokkos::MemoryUnmanaged>& D_inv,
404  const Scalar& omega,
405  const ESweepDirection direction) const;
406 
433  LO
434  replaceLocalValues (const LO localRowInd,
435  const LO colInds[],
436  const Scalar vals[],
437  const LO numColInds) const;
438 
465  LO
466  sumIntoLocalValues (const LO localRowInd,
467  const LO colInds[],
468  const Scalar vals[],
469  const LO numColInds) const;
470 
500  LO
501  getLocalRowView (const LO localRowInd,
502  const LO*& colInds,
503  Scalar*& vals,
504  LO& numInds) const;
505 
507  void
508  getLocalRowView (LO LocalRow,
509  Teuchos::ArrayView<const LO> &indices,
510  Teuchos::ArrayView<const Scalar> &values) const;
511 
513  void
514  getLocalRowCopy (LO LocalRow,
515  const Teuchos::ArrayView<LO> &Indices,
516  const Teuchos::ArrayView<Scalar> &Values,
517  size_t &NumEntries) const;
518 
520  getLocalBlock (const LO localRowInd, const LO localColInd) const;
521 
545  LO
546  getLocalRowOffsets (const LO localRowInd,
547  ptrdiff_t offsets[],
548  const LO colInds[],
549  const LO numColInds) const;
550 
556  LO
557  replaceLocalValuesByOffsets (const LO localRowInd,
558  const ptrdiff_t offsets[],
559  const Scalar vals[],
560  const LO numOffsets) const;
561 
567  LO
568  sumIntoLocalValuesByOffsets (const LO localRowInd,
569  const ptrdiff_t offsets[],
570  const Scalar vals[],
571  const LO numOffsets) const;
572 
579  size_t getNumEntriesInLocalRow (const LO localRowInd) const;
580 
597  bool localError () const {
598  return *localError_;
599  }
600 
615  std::string errorMessages () const {
616  return (*errs_).is_null () ? std::string ("") : (*errs_)->str ();
617  }
618 
650  void
651  getLocalDiagOffsets (const Kokkos::View<size_t*, device_type,
652  Kokkos::MemoryUnmanaged>& offsets) const;
653 
654 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
655  void TPETRA_DEPRECATED
661  getLocalDiagOffsets (Teuchos::ArrayRCP<size_t>& offsets) const;
662 #endif // TPETRA_ENABLE_DEPRECATED_CODE
663 
677  void
678  getLocalDiagCopy (const Kokkos::View<impl_scalar_type***, device_type,
679  Kokkos::MemoryUnmanaged>& diag,
680  const Kokkos::View<const size_t*, device_type,
681  Kokkos::MemoryUnmanaged>& offsets) const;
682 
696  void
697  getLocalDiagCopy (const Kokkos::View<impl_scalar_type***, device_type,
698  Kokkos::MemoryUnmanaged>& diag,
699  const Teuchos::ArrayView<const size_t>& offsets) const;
700 
701 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
702  void TPETRA_DEPRECATED
717  const Teuchos::ArrayView<const size_t>& offsets) const;
718 #endif // TPETRA_ENABLE_DEPRECATED_CODE
719 
720 protected:
722  LO
723  absMaxLocalValues (const LO localRowInd,
724  const LO colInds[],
725  const Scalar vals[],
726  const LO numColInds) const;
727 
729  LO
730  absMaxLocalValuesByOffsets (const LO localRowInd,
731  const ptrdiff_t offsets[],
732  const Scalar vals[],
733  const LO numOffsets) const;
734 
740 
741 
746  using buffer_device_type = typename DistObject<Scalar, LO, GO,
748 
749  virtual bool checkSizes (const ::Tpetra::SrcDistObject& source);
750 
751  virtual void
752 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
753  copyAndPermuteNew
754 #else // TPETRA_ENABLE_DEPRECATED_CODE
755  copyAndPermute
756 #endif // TPETRA_ENABLE_DEPRECATED_CODE
757  (const SrcDistObject& sourceObj,
758  const size_t numSameIDs,
759  const Kokkos::DualView<const local_ordinal_type*,
760  buffer_device_type>& permuteToLIDs,
761  const Kokkos::DualView<const local_ordinal_type*,
762  buffer_device_type>& permuteFromLIDs);
763 
764  virtual void
765 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
766  packAndPrepareNew
767 #else // TPETRA_ENABLE_DEPRECATED_CODE
768  packAndPrepare
769 #endif // TPETRA_ENABLE_DEPRECATED_CODE
770  (const SrcDistObject& sourceObj,
771  const Kokkos::DualView<const local_ordinal_type*,
772  buffer_device_type>& exportLIDs,
773  Kokkos::DualView<packet_type*,
774  buffer_device_type>& exports,
775  Kokkos::DualView<size_t*,
776  buffer_device_type> numPacketsPerLID,
777  size_t& constantNumPackets,
778  Distributor& /* distor */);
779 
780  virtual void
781 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
782  unpackAndCombineNew
783 #else // TPETRA_ENABLE_DEPRECATED_CODE
784  unpackAndCombine
785 #endif // TPETRA_ENABLE_DEPRECATED_CODE
786  (const Kokkos::DualView<const local_ordinal_type*,
787  buffer_device_type>& importLIDs,
788  Kokkos::DualView<packet_type*,
789  buffer_device_type> imports,
790  Kokkos::DualView<size_t*,
791  buffer_device_type> numPacketsPerLID,
792  const size_t constantNumPackets,
793  Distributor& /* distor */,
794  const CombineMode combineMode);
796 
797 private:
799  crs_graph_type graph_;
800  Teuchos::RCP<crs_graph_type> graphRCP_;
809  map_type rowMeshMap_;
816  map_type domainPointMap_;
823  map_type rangePointMap_;
825  LO blockSize_;
826 
840  typename crs_graph_type::local_graph_type::row_map_type::HostMirror ptrHost_;
841 
847  typename crs_graph_type::local_graph_type::entries_type::HostMirror indHost_;
848 
854  typename Kokkos::DualView<impl_scalar_type*, device_type> val_;
855 
877  Teuchos::RCP<Teuchos::RCP<BMV> > X_colMap_;
881  Teuchos::RCP<Teuchos::RCP<BMV> > Y_rowMap_;
882 
890  Teuchos::RCP<Teuchos::RCP<typename crs_graph_type::import_type> > pointImporter_;
891 
893  LO offsetPerBlock_;
894 
906  Teuchos::RCP<bool> localError_;
907 
915  Teuchos::RCP<Teuchos::RCP<std::ostringstream> > errs_;
916 
918  std::ostream& markLocalErrorAndGetStream ();
919 
920  // //! Clear the local error state and stream.
921  // void clearLocalErrorStateAndStream ();
922 
923  template<class Device>
924  struct is_cuda {
925 #if defined(KOKKOS_ENABLE_CUDA)
926  // CudaHostPinnedSpace::execution_space ==
927  // HostSpace::execution_space. That's OK; it's host memory, that
928  // just happens to be Cuda accessible. But what if somebody gives
929  // us Device<Cuda, CudaHostPinnedSpace>? It looks like they mean
930  // to run on device then, so we should sync to device.
931  static constexpr bool value =
932  std::is_same<typename Device::execution_space, Kokkos::Cuda>::value;
933 #else
934  static constexpr bool value = false;
935 #endif // defined(KOKKOS_ENABLE_CUDA)
936  };
937 
938 public:
940 
941 
943  inline void modify_host()
944  {
945  val_.modify_host();
946  }
947 
949  inline void modify_device()
950  {
951  val_.modify_device();
952  }
953 
955  template<class MemorySpace>
956  void modify ()
957  {
958  if (is_cuda<MemorySpace>::value) {
959  this->modify_device ();
960  }
961  else {
962  this->modify_host ();
963  }
964  }
965 
967  inline bool need_sync_host() const
968  {
969  return val_.need_sync_host();
970  }
971 
973  inline bool need_sync_device() const
974  {
975  return val_.need_sync_device();
976  }
977 
979  template<class MemorySpace>
980  bool need_sync () const
981  {
982  if (is_cuda<MemorySpace>::value) {
983  return this->need_sync_device ();
984  }
985  else {
986  return this->need_sync_host ();
987  }
988  }
989 
991  inline void sync_host()
992  {
993  val_.sync_host();
994  }
995 
997  inline void sync_device()
998  {
999  val_.sync_device();
1000  }
1001 
1003  template<class MemorySpace>
1004  void sync ()
1005  {
1006  if (is_cuda<MemorySpace>::value) {
1007  this->sync_device ();
1008  }
1009  else {
1010  this->sync_host ();
1011  }
1012  }
1013 
1014  // \brief Get the host view of the matrix's values
1015  typename Kokkos::DualView<impl_scalar_type*, device_type>::t_host getValuesHost () const {
1016  return val_.view_host();
1017  }
1018 
1019  // \brief Get the device view of the matrix's values
1020  typename Kokkos::DualView<impl_scalar_type*, device_type>::t_dev getValuesDevice () const {
1021  return val_.view_device();
1022  }
1023 
1041  template<class MemorySpace>
1042  typename std::conditional<is_cuda<MemorySpace>::value,
1043  typename Kokkos::DualView<impl_scalar_type*, device_type>::t_dev,
1044  typename Kokkos::DualView<impl_scalar_type*, device_type>::t_host>::type
1046  {
1047  // Unlike std::conditional, if_c has a select method.
1048  return Kokkos::Impl::if_c<
1049  is_cuda<MemorySpace>::value,
1050  typename Kokkos::DualView<impl_scalar_type*, device_type>::t_dev,
1051  typename Kokkos::DualView<impl_scalar_type*, device_type>::t_host
1052  >::select (this->getValuesDevice (), this->getValuesHost ());
1053  }
1054 
1056 
1057 private:
1058 
1068  void
1069  applyBlockTrans (const BlockMultiVector<Scalar, LO, GO, Node>& X,
1071  const Teuchos::ETransp mode,
1072  const Scalar alpha,
1073  const Scalar beta);
1074 
1082  void
1083  applyBlockNoTrans (const BlockMultiVector<Scalar, LO, GO, Node>& X,
1085  const Scalar alpha,
1086  const Scalar beta);
1087 
1095  void
1096  localApplyBlockNoTrans (const BlockMultiVector<Scalar, LO, GO, Node>& X,
1098  const Scalar alpha,
1099  const Scalar beta);
1100 
1140  LO
1141  findRelOffsetOfColumnIndex (const LO localRowIndex,
1142  const LO colIndexToFind,
1143  const LO hint = 0) const;
1144 
1147  LO offsetPerBlock () const;
1148 
1150  getConstLocalBlockFromInput (const impl_scalar_type* val, const size_t pointOffset) const;
1151 
1153  getNonConstLocalBlockFromInput (impl_scalar_type* val, const size_t pointOffset) const;
1154 
1156  getConstLocalBlockFromAbsOffset (const size_t absBlockOffset) const;
1157 
1159  getNonConstLocalBlockFromAbsOffset (const size_t absBlockOffset) const;
1160 
1165  getConstLocalBlockFromRelOffset (const LO lclMeshRow,
1166  const size_t relMeshOffset) const;
1167 
1168 public:
1170  virtual Teuchos::RCP<const Teuchos::Comm<int> > getComm() const;
1171 
1172 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
1173  virtual TPETRA_DEPRECATED Teuchos::RCP<Node> getNode() const;
1175 #endif // TPETRA_ENABLE_DEPRECATED_CODE
1176 
1178  virtual global_size_t getGlobalNumCols() const;
1179 
1180  virtual size_t getNodeNumCols() const;
1181 
1182  virtual GO getIndexBase() const;
1183 
1185  virtual global_size_t getGlobalNumEntries() const;
1186 
1188  virtual size_t getNodeNumEntries() const;
1189 
1199  virtual size_t getNumEntriesInGlobalRow (GO globalRow) const;
1200 
1203  virtual size_t getGlobalMaxNumRowEntries () const;
1204 
1206  virtual bool hasColMap () const;
1207 
1217  virtual bool isLocallyIndexed () const;
1218 
1228  virtual bool isGloballyIndexed () const;
1229 
1231  virtual bool isFillComplete () const;
1232 
1234  virtual bool supportsRowViews () const;
1235 
1236 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
1237  virtual global_size_t TPETRA_DEPRECATED getGlobalNumDiags() const;
1245 
1253  virtual size_t TPETRA_DEPRECATED getNodeNumDiags() const;
1254 
1265  virtual bool TPETRA_DEPRECATED isLowerTriangular () const;
1266 
1277  virtual bool TPETRA_DEPRECATED isUpperTriangular () const;
1278 #endif // TPETRA_ENABLE_DEPRECATED_CODE
1279 
1281 
1283 
1304  virtual void
1305  getGlobalRowCopy (GO GlobalRow,
1306  const Teuchos::ArrayView<GO> &Indices,
1307  const Teuchos::ArrayView<Scalar> &Values,
1308  size_t& NumEntries) const;
1309 
1334  virtual void
1335  getGlobalRowView (GO GlobalRow,
1336  Teuchos::ArrayView<const GO>& indices,
1337  Teuchos::ArrayView<const Scalar>& values) const;
1338 
1350  virtual void getLocalDiagCopy (::Tpetra::Vector<Scalar,LO,GO,Node>& diag) const;
1351 
1353 
1355 
1361  virtual void leftScale (const ::Tpetra::Vector<Scalar, LO, GO, Node>& x);
1362 
1368  virtual void rightScale (const ::Tpetra::Vector<Scalar, LO, GO, Node>& x);
1369 
1378  virtual typename ::Tpetra::RowMatrix<Scalar, LO, GO, Node>::mag_type
1379  getFrobeniusNorm () const;
1381 };
1382 
1383 } // namespace Experimental
1384 } // namespace Tpetra
1385 
1386 #endif // TPETRA_EXPERIMENTAL_BLOCKCRSMATRIX_DECL_HPP
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 void leftScale(const ::Tpetra::Vector< Scalar, LO, GO, Node > &x)
Scale the RowMatrix on the left with the given Vector x.
virtual bool supportsRowViews() const
Whether this object implements getLocalRowView() and getGlobalRowView().
Teuchos::RCP< const map_type > getColMap() const
get the (mesh) map for the columns of this block matrix.
LO getBlockSize() const
The number of degrees of freedom per mesh point.
Kokkos::View< impl_scalar_type *, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > little_vec_type
The type used to access nonconst vector blocks.
Teuchos::RCP< const map_type > getRowMap() const
get the (mesh) map for the rows of this block matrix.
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.
LO replaceLocalValuesByOffsets(const LO localRowInd, const ptrdiff_t offsets[], const Scalar vals[], const LO numOffsets) const
Like replaceLocalValues, but avoids computing row offsets.
size_t getNodeNumRows() const
get the local number of block rows
Kokkos::View< impl_scalar_type **, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > little_block_type
The type used to access nonconst matrix blocks.
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.
Node::device_type device_type
The Kokkos::Device specialization that this class uses.
virtual void rightScale(const ::Tpetra::Vector< Scalar, LO, GO, Node > &x)
Scale the RowMatrix on the right with the given Vector x.
bool need_sync() const
Whether the matrix&#39;s values need sync&#39;ing to 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.
MultiVector for multiple degrees of freedom per mesh point.
Declaration of the Tpetra::CrsMatrix class.
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.
global_size_t getGlobalNumRows() const
get the global number of block rows
::Tpetra::MultiVector< Scalar, LO, GO, node_type > mv_type
The implementation of MultiVector that this class uses.
void modify_device()
Mark the matrix&#39;s valueas as modified in device space.
void modify_host()
Mark the matrix&#39;s valueas as modified in host space.
LO absMaxLocalValues(const LO localRowInd, const LO colInds[], const Scalar vals[], const LO numColInds) const
Like sumIntoLocalValues, but for the ABSMAX combine mode.
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.
virtual bool hasColMap() const
Whether this matrix has a well-defined column Map.
Teuchos::RCP< const map_type > getDomainMap() const
Get the (point) domain Map of this matrix.
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to alpha.
virtual Teuchos::RCP< const ::Tpetra::RowGraph< LO, GO, Node > > getGraph() const
Get the (mesh) graph.
typename DistObject< Scalar, LO, GO, Node >::buffer_device_type buffer_device_type
Kokkos::Device specialization for communication buffers.
::Tpetra::CrsGraph< LO, GO, node_type > crs_graph_type
The implementation of CrsGraph that this class uses.
virtual GO getIndexBase() const
The index base for global indices in this matrix.
virtual bool isGloballyIndexed() const
Whether matrix indices are globally indexed.
::Tpetra::Map< LO, GO, node_type > map_type
The implementation of Map that this class uses.
size_t global_size_t
Global size_t object.
device_type::execution_space execution_space
The Kokkos execution space that this class uses.
std::string description() const
One-line description of this object.
BlockCrsMatrix()
Default constructor: Makes an empty block matrix.
virtual global_size_t getGlobalNumCols() const
The global number of columns of this matrix.
ESweepDirection
Sweep direction for Gauss-Seidel or Successive Over-Relaxation (SOR).
virtual bool isFillComplete() const
Whether fillComplete() has been called.
Sets up and executes a communication plan for a Tpetra DistObject.
CombineMode
Rule for combining data in an Import or Export.
virtual size_t getNodeNumCols() const
The number of columns needed to apply the forward operator on this node.
virtual size_t getNodeNumEntries() const
The local number of stored (structurally nonzero) entries.
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 global_size_t getGlobalNumEntries() const
The global number of stored (structurally nonzero) entries.
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.
void getLocalRowCopy(LO LocalRow, const Teuchos::ArrayView< LO > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const
Not implemented.
Abstract base class for objects that can be the source of an Import or Export operation.
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.
Scalar scalar_type
The type of entries in the matrix (that is, of each entry in each block).
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.
std::string errorMessages() const
The current stream of error messages.
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.
void modify()
Mark the matrix&#39;s values as modified in the given memory space.
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.
virtual size_t getGlobalMaxNumRowEntries() const
The maximum number of entries in any row over all processes in the matrix&#39;s communicator.
void sync_device()
Sync the matrix&#39;s values to device space.
typename BMV::impl_scalar_type impl_scalar_type
The implementation type of entries in the matrix.
A read-only, row-oriented interface to a sparse matrix.
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.
void sync_host()
Sync the matrix&#39;s values to host space.
LO sumIntoLocalValuesByOffsets(const LO localRowInd, const ptrdiff_t offsets[], const Scalar vals[], const LO numOffsets) const
Like sumIntoLocalValues, but avoids computing row offsets.
A distributed dense vector.
Sparse matrix whose entries are small dense square blocks, all of the same dimensions.
virtual typename::Tpetra::RowMatrix< Scalar, LO, GO, Node >::mag_type getFrobeniusNorm() const
The Frobenius norm of the matrix.
void getLocalDiagOffsets(const Kokkos::View< size_t *, device_type, Kokkos::MemoryUnmanaged > &offsets) const
Get offsets of the diagonal entries in the matrix.
device_type::memory_space memory_space
The Kokkos memory space that this class uses.
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.
typename mv_type::impl_scalar_type impl_scalar_type
The implementation type of entries in the object.
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which this matrix is distributed.
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.
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_).
bool localError() const
Whether this object had an error on the calling process.
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.
virtual ~BlockCrsMatrix()
Destructor (declared virtual for memory safety).
void sync()
Sync the matrix&#39;s values to the given memory space.
Base class for distributed Tpetra objects that support data redistribution.
size_t getNodeMaxNumRowEntries() const
Maximum number of entries in any row of the matrix, on this process.
bool hasTransposeApply() const
Whether it is valid to apply the transpose or conjugate transpose of this matrix. ...
virtual size_t getNumEntriesInGlobalRow(GO globalRow) const
The current number of entries on the calling process in the specified global row. ...
size_t getNumEntriesInLocalRow(const LO localRowInd) const
Return the number of entries in the given row on the calling process.
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.
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:...
virtual bool isLocallyIndexed() const
Whether matrix indices are locally indexed.
bool need_sync_host() const
Whether the matrix&#39;s values need sync&#39;ing to host space.