Tpetra parallel linear algebra  Version of the Day
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.hpp"
52 
53 namespace Tpetra {
54 namespace Experimental {
55 namespace Classes {
56 
134 template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type,
136  class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type,
139  virtual public ::Tpetra::RowMatrix<Scalar, LO, GO, Node>,
140  virtual public ::Tpetra::DistObject<char, LO, GO, Node>
141 {
142 private:
145  using STS = Teuchos::ScalarTraits<Scalar>;
146 
147 protected:
149  typedef char packet_type;
150 
151 public:
153 
154 
156  using scalar_type = Scalar;
157 
165 
167  typedef LO local_ordinal_type;
174  typedef Node node_type;
175 
177  typedef typename Node::device_type device_type;
179  typedef typename device_type::execution_space execution_space;
181  typedef typename device_type::memory_space memory_space;
182 
184  typedef ::Tpetra::Map<LO, GO, node_type> map_type;
186  typedef ::Tpetra::MultiVector<Scalar, LO, GO, node_type> mv_type;
188  typedef ::Tpetra::CrsGraph<LO, GO, node_type> crs_graph_type;
189 
191  typedef Kokkos::View<impl_scalar_type**,
192  Kokkos::LayoutRight,
193  device_type,
194  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
197  typedef Kokkos::View<const impl_scalar_type**,
198  Kokkos::LayoutRight,
199  device_type,
200  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
203  typedef Kokkos::View<impl_scalar_type*,
204  Kokkos::LayoutRight,
205  device_type,
206  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
209  typedef Kokkos::View<const impl_scalar_type*,
210  Kokkos::LayoutRight,
211  device_type,
212  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
214 
216 
218 
220  BlockCrsMatrix ();
221 
231  BlockCrsMatrix (const crs_graph_type& graph, const LO blockSize);
232 
240  BlockCrsMatrix (const crs_graph_type& graph,
241  const map_type& domainPointMap,
242  const map_type& rangePointMap,
243  const LO blockSize);
244 
246  virtual ~BlockCrsMatrix () {}
247 
249 
251 
253  Teuchos::RCP<const map_type> getDomainMap () const;
254 
256  Teuchos::RCP<const map_type> getRangeMap () const;
257 
259  Teuchos::RCP<const map_type> getRowMap () const;
260 
262  Teuchos::RCP<const map_type> getColMap () const;
263 
266 
268  size_t getNodeNumRows() const;
269 
270  size_t getNodeMaxNumRowEntries() const;
271 
281  void
282  apply (const mv_type& X,
283  mv_type& Y,
284  Teuchos::ETransp mode = Teuchos::NO_TRANS,
285  Scalar alpha = Teuchos::ScalarTraits<Scalar>::one (),
286  Scalar beta = Teuchos::ScalarTraits<Scalar>::zero ()) const;
287 
290  bool hasTransposeApply () const {
291  // FIXME (mfh 04 May 2014) Transpose and conjugate transpose modes
292  // are not implemented yet. Fill in applyBlockTrans() to fix this.
293  return false;
294  }
295 
297  void setAllToScalar (const Scalar& alpha);
298 
300 
302 
304  std::string description () const;
305 
329  void
330  describe (Teuchos::FancyOStream& out,
331  const Teuchos::EVerbosityLevel verbLevel) const;
332 
334 
336 
338  LO getBlockSize () const { return blockSize_; }
339 
341  virtual Teuchos::RCP<const ::Tpetra::RowGraph<LO,GO,Node> > getGraph () const;
342 
343  const crs_graph_type & getCrsGraph () const { return graph_; }
344 
349  void
352  Teuchos::ETransp mode = Teuchos::NO_TRANS,
353  const Scalar alpha = Teuchos::ScalarTraits<Scalar>::one (),
354  const Scalar beta = Teuchos::ScalarTraits<Scalar>::zero ());
355 
359  void
361  const ::Tpetra::MultiVector<Scalar,LO,GO,Node> &B,
362  const ::Tpetra::MultiVector<Scalar,LO,GO,Node> &D,
363  const Scalar& dampingFactor,
364  const ESweepDirection direction,
365  const int numSweeps,
366  const bool zeroInitialGuess) const;
367 
371  void
373  const ::Tpetra::MultiVector<Scalar,LO,GO,Node>& B,
374  const ::Tpetra::MultiVector<Scalar,LO,GO,Node>& D,
375  const Teuchos::ArrayView<LO>& rowIndices,
376  const Scalar& dampingFactor,
377  const ESweepDirection direction,
378  const int numSweeps,
379  const bool zeroInitialGuess) const;
380 
400  void
403  const Kokkos::View<impl_scalar_type***, device_type,
404  Kokkos::MemoryUnmanaged>& D_inv,
405  const Scalar& omega,
406  const ESweepDirection direction) const;
407 
434  LO
435  replaceLocalValues (const LO localRowInd,
436  const LO colInds[],
437  const Scalar vals[],
438  const LO numColInds) const;
439 
466  LO
467  sumIntoLocalValues (const LO localRowInd,
468  const LO colInds[],
469  const Scalar vals[],
470  const LO numColInds) const;
471 
501  LO
502  getLocalRowView (const LO localRowInd,
503  const LO*& colInds,
504  Scalar*& vals,
505  LO& numInds) const;
506 
508  void
509  getLocalRowView (LO LocalRow,
510  Teuchos::ArrayView<const LO> &indices,
511  Teuchos::ArrayView<const Scalar> &values) const;
512 
514  void
515  getLocalRowCopy (LO LocalRow,
516  const Teuchos::ArrayView<LO> &Indices,
517  const Teuchos::ArrayView<Scalar> &Values,
518  size_t &NumEntries) const;
519 
521  getLocalBlock (const LO localRowInd, const LO localColInd) const;
522 
546  LO
547  getLocalRowOffsets (const LO localRowInd,
548  ptrdiff_t offsets[],
549  const LO colInds[],
550  const LO numColInds) const;
551 
557  LO
558  replaceLocalValuesByOffsets (const LO localRowInd,
559  const ptrdiff_t offsets[],
560  const Scalar vals[],
561  const LO numOffsets) const;
562 
568  LO
569  sumIntoLocalValuesByOffsets (const LO localRowInd,
570  const ptrdiff_t offsets[],
571  const Scalar vals[],
572  const LO numOffsets) const;
573 
580  size_t getNumEntriesInLocalRow (const LO localRowInd) const;
581 
598  bool localError () const {
599  return *localError_;
600  }
601 
616  std::string errorMessages () const {
617  return (*errs_).is_null () ? std::string ("") : (*errs_)->str ();
618  }
619 
651  void
652  getLocalDiagOffsets (const Kokkos::View<size_t*, device_type,
653  Kokkos::MemoryUnmanaged>& offsets) const;
654 
660  void TPETRA_DEPRECATED
661  getLocalDiagOffsets (Teuchos::ArrayRCP<size_t>& offsets) const;
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 
713  void TPETRA_DEPRECATED
715  const Teuchos::ArrayView<const size_t>& offsets) const;
716 
717 protected:
719  LO
720  absMaxLocalValues (const LO localRowInd,
721  const LO colInds[],
722  const Scalar vals[],
723  const LO numColInds) const;
724 
726  LO
727  absMaxLocalValuesByOffsets (const LO localRowInd,
728  const ptrdiff_t offsets[],
729  const Scalar vals[],
730  const LO numOffsets) const;
731 
737 
738 
739  virtual bool checkSizes (const ::Tpetra::SrcDistObject& source);
740 
741  virtual void
742  copyAndPermute (const ::Tpetra::SrcDistObject& source,
743  size_t numSameIDs,
744  const Teuchos::ArrayView<const LO>& permuteToLIDs,
745  const Teuchos::ArrayView<const LO>& permuteFromLIDs);
746 
747  virtual void
748  packAndPrepare (const ::Tpetra::SrcDistObject& source,
749  const Teuchos::ArrayView<const LO>& exportLIDs,
750  Teuchos::Array<packet_type>& exports,
751  const Teuchos::ArrayView<size_t>& numPacketsPerLID,
752  size_t& constantNumPackets,
753  ::Tpetra::Distributor& distor);
754 
755  virtual void
756  unpackAndCombine (const Teuchos::ArrayView<const LO> &importLIDs,
757  const Teuchos::ArrayView<const packet_type> &imports,
758  const Teuchos::ArrayView<size_t> &numPacketsPerLID,
759  size_t constantNumPackets,
760  ::Tpetra::Distributor& distor,
761  ::Tpetra::CombineMode CM);
763 
764 private:
766  crs_graph_type graph_;
767  Teuchos::RCP<crs_graph_type> graphRCP_;
776  map_type rowMeshMap_;
783  map_type domainPointMap_;
790  map_type rangePointMap_;
792  LO blockSize_;
793 
807  typename crs_graph_type::local_graph_type::row_map_type::HostMirror ptrHost_;
808 
814  typename crs_graph_type::local_graph_type::entries_type::HostMirror indHost_;
815 
821  typename Kokkos::DualView<impl_scalar_type*, device_type> val_;
822 
844  Teuchos::RCP<Teuchos::RCP<BMV> > X_colMap_;
848  Teuchos::RCP<Teuchos::RCP<BMV> > Y_rowMap_;
849 
857  Teuchos::RCP<Teuchos::RCP<typename crs_graph_type::import_type> > pointImporter_;
858 
860  LO offsetPerBlock_;
861 
873  Teuchos::RCP<bool> localError_;
874 
882  Teuchos::RCP<Teuchos::RCP<std::ostringstream> > errs_;
883 
885  std::ostream& markLocalErrorAndGetStream ();
886 
887  // //! Clear the local error state and stream.
888  // void clearLocalErrorStateAndStream ();
889 
890 public:
892 
893 
895  template<class MemorySpace>
896  void modify ()
897  {
898  // It's legit to use a memory space, execution space, or
899  // Kokkos::Device specialization as the template parameter of
900  // Kokkos::DualView::modify. That method just extracts the
901  // memory_space typedef of its template parameter anyway.
902  // However, insisting on a memory space avoids unnecessary
903  // instantiations.
904  val_.template modify<typename MemorySpace::memory_space> ();
905  }
906 
908  template<class MemorySpace>
909  bool need_sync () const
910  {
911  // It's legit to use a memory space, execution space, or
912  // Kokkos::Device specialization as the template parameter of
913  // Kokkos::DualView::need_sync. That method just extracts the
914  // memory_space typedef of its template parameter anyway.
915  // However, insisting on a memory space avoids unnecessary
916  // instantiations.
917  return val_.template need_sync<typename MemorySpace::memory_space> ();
918  }
919 
921  template<class MemorySpace>
922  void sync ()
923  {
924  // It's legit to use a memory space, execution space, or
925  // Kokkos::Device specialization as the template parameter of
926  // Kokkos::DualView::sync. That method just extracts the
927  // memory_space typedef of its template parameter anyway.
928  // However, insisting on a memory space avoids unnecessary
929  // instantiations.
930  val_.template sync<typename MemorySpace::memory_space> ();
931  }
932 
947  template<class MemorySpace>
948  auto getValues () -> decltype (val_.template view<typename MemorySpace::memory_space> ())
949  {
950  return val_.template view<typename MemorySpace::memory_space> ();
951  }
952 
954 
955 private:
956 
966  void
967  applyBlockTrans (const BlockMultiVector<Scalar, LO, GO, Node>& X,
969  const Teuchos::ETransp mode,
970  const Scalar alpha,
971  const Scalar beta);
972 
980  void
981  applyBlockNoTrans (const BlockMultiVector<Scalar, LO, GO, Node>& X,
983  const Scalar alpha,
984  const Scalar beta);
985 
993  void
994  localApplyBlockNoTrans (const BlockMultiVector<Scalar, LO, GO, Node>& X,
996  const Scalar alpha,
997  const Scalar beta);
998 
1038  LO
1039  findRelOffsetOfColumnIndex (const LO localRowIndex,
1040  const LO colIndexToFind,
1041  const LO hint = 0) const;
1042 
1045  LO offsetPerBlock () const;
1046 
1048  getConstLocalBlockFromInput (const impl_scalar_type* val, const size_t pointOffset) const;
1049 
1051  getNonConstLocalBlockFromInput (impl_scalar_type* val, const size_t pointOffset) const;
1052 
1054  getConstLocalBlockFromAbsOffset (const size_t absBlockOffset) const;
1055 
1057  getNonConstLocalBlockFromAbsOffset (const size_t absBlockOffset) const;
1058 
1063  getConstLocalBlockFromRelOffset (const LO lclMeshRow,
1064  const size_t relMeshOffset) const;
1065 
1066 public:
1068  virtual Teuchos::RCP<const Teuchos::Comm<int> > getComm() const;
1069 
1071  virtual Teuchos::RCP<Node> getNode() const;
1072 
1074  virtual global_size_t getGlobalNumCols() const;
1075 
1076  virtual size_t getNodeNumCols() const;
1077 
1078  virtual GO getIndexBase() const;
1079 
1081  virtual global_size_t getGlobalNumEntries() const;
1082 
1084  virtual size_t getNodeNumEntries() const;
1085 
1095  virtual size_t getNumEntriesInGlobalRow (GO globalRow) const;
1096 
1104  virtual global_size_t TPETRA_DEPRECATED getGlobalNumDiags() const;
1105 
1113  virtual size_t TPETRA_DEPRECATED getNodeNumDiags() const;
1114 
1116  virtual size_t getGlobalMaxNumRowEntries() const;
1117 
1119  virtual bool hasColMap() const;
1120 
1131  virtual bool TPETRA_DEPRECATED isLowerTriangular () const;
1132 
1143  virtual bool TPETRA_DEPRECATED isUpperTriangular () const;
1144 
1154  virtual bool isLocallyIndexed() const;
1155 
1165  virtual bool isGloballyIndexed() const;
1166 
1168  virtual bool isFillComplete() const;
1169 
1171  virtual bool supportsRowViews() const;
1172 
1174 
1176 
1197  virtual void
1198  getGlobalRowCopy (GO GlobalRow,
1199  const Teuchos::ArrayView<GO> &Indices,
1200  const Teuchos::ArrayView<Scalar> &Values,
1201  size_t& NumEntries) const;
1202 
1227  virtual void
1228  getGlobalRowView (GO GlobalRow,
1229  Teuchos::ArrayView<const GO>& indices,
1230  Teuchos::ArrayView<const Scalar>& values) const;
1231 
1243  virtual void getLocalDiagCopy (::Tpetra::Vector<Scalar,LO,GO,Node>& diag) const;
1244 
1246 
1248 
1254  virtual void leftScale (const ::Tpetra::Vector<Scalar, LO, GO, Node>& x);
1255 
1261  virtual void rightScale (const ::Tpetra::Vector<Scalar, LO, GO, Node>& x);
1262 
1271  virtual typename ::Tpetra::RowMatrix<Scalar, LO, GO, Node>::mag_type
1272  getFrobeniusNorm () const;
1274 };
1275 
1276 } // namespace Classes
1277 } // namespace Experimental
1278 } // namespace Tpetra
1279 
1280 #endif // TPETRA_EXPERIMENTAL_BLOCKCRSMATRIX_DECL_HPP
Tpetra::Experimental::Classes::BlockCrsMatrix::describe
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
Print a description of this object to the given output stream.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3285
Tpetra::Experimental::Classes::BlockCrsMatrix::getNodeNumCols
virtual size_t getNodeNumCols() const
The number of columns needed to apply the forward operator on this node.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3563
Tpetra::ESweepDirection
ESweepDirection
Sweep direction for Gauss-Seidel or Successive Over-Relaxation (SOR).
Definition: Tpetra_ConfigDefs.hpp:245
Tpetra::Experimental::Classes::BlockCrsMatrix::replaceLocalValues
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,...
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:981
Tpetra::Experimental::Classes::BlockCrsMatrix::sumIntoLocalValuesByOffsets
LO sumIntoLocalValuesByOffsets(const LO localRowInd, const ptrdiff_t offsets[], const Scalar vals[], const LO numOffsets) const
Like sumIntoLocalValues, but avoids computing row offsets.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:1676
Tpetra::Experimental::Classes::BlockCrsMatrix::isFillComplete
virtual bool isFillComplete() const
Whether fillComplete() has been called.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3671
Tpetra::Experimental::Classes::BlockCrsMatrix::getGlobalNumEntries
virtual global_size_t getGlobalNumEntries() const
The global number of stored (structurally nonzero) entries.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3579
Tpetra::Experimental::Classes::BlockCrsMatrix::absMaxLocalValues
LO absMaxLocalValues(const LO localRowInd, const LO colInds[], const Scalar vals[], const LO numColInds) const
Like sumIntoLocalValues, but for the ABSMAX combine mode.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:1371
Tpetra::Experimental::Classes::BlockCrsMatrix::const_little_block_type
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.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:201
Tpetra::Experimental::Classes::BlockCrsMatrix::getGlobalMaxNumRowEntries
virtual size_t getGlobalMaxNumRowEntries() const
The maximum number of entries across all rows/columns on all nodes.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3621
Tpetra::Experimental::Classes::BlockCrsMatrix::description
std::string description() const
One-line description of this object.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3261
Tpetra::Experimental::Classes::BlockCrsMatrix::const_little_vec_type
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.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:213
Tpetra::Experimental::Classes::BlockCrsMatrix::getGraph
virtual Teuchos::RCP< const ::Tpetra::RowGraph< LO, GO, Node > > getGraph() const
Get the (mesh) graph.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3807
Tpetra::Experimental::Classes::BlockCrsMatrix::sync
void sync()
Sync the matrix's values to the given memory space.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:922
Tpetra::Experimental::Classes::BlockCrsMatrix::need_sync
bool need_sync() const
Whether the matrix's values need sync'ing to the given memory space.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:909
Tpetra::Experimental::Classes::BlockCrsMatrix::getRowMap
Teuchos::RCP< const map_type > getRowMap() const
get the (mesh) map for the rows of this block matrix.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:809
Tpetra::Experimental::Classes::BlockCrsMatrix::getLocalRowView
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.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:1491
Tpetra::Experimental::Classes::BlockCrsMatrix::replaceLocalValuesByOffsets
LO replaceLocalValuesByOffsets(const LO localRowInd, const ptrdiff_t offsets[], const Scalar vals[], const LO numOffsets) const
Like replaceLocalValues, but avoids computing row offsets.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:1598
Tpetra::Details::DefaultTypes::node_type
::Kokkos::Compat::KokkosDeviceWrapperNode< execution_space > node_type
Default value of Node template parameter.
Definition: Tpetra_Details_DefaultTypes.hpp:105
Tpetra::Experimental::Classes::BlockCrsMatrix::apply
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.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:849
Tpetra::Experimental::Classes::BlockCrsMatrix::getBlockSize
LO getBlockSize() const
The number of degrees of freedom per mesh point.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:338
Tpetra::Experimental::Classes::BlockCrsMatrix::modify
void modify()
Mark the matrix's values as modified in the given memory space.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:896
Tpetra::Experimental::Classes::BlockCrsMatrix::getValues
auto getValues() -> decltype(val_.template view< typename MemorySpace::memory_space >())
Get the host or device View of the matrix's values (val_).
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:948
Tpetra::Experimental::Classes::BlockCrsMatrix::errorMessages
std::string errorMessages() const
The current stream of error messages.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:616
Tpetra::Experimental::Classes::BlockCrsMatrix::getGlobalRowView
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.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3702
Tpetra::Experimental::Classes::BlockCrsMatrix::isUpperTriangular
virtual bool TPETRA_DEPRECATED isUpperTriangular() const
Whether the matrix's graph is locally upper triangular.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3646
Tpetra::Experimental::Classes::BlockCrsMatrix::little_vec_type
Kokkos::View< impl_scalar_type *, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > little_vec_type
The type used to access nonconst vector blocks.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:207
Tpetra::Experimental::Classes::BlockCrsMatrix::leftScale
virtual void leftScale(const ::Tpetra::Vector< Scalar, LO, GO, Node > &x)
Scale the RowMatrix on the left with the given Vector x.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3785
Tpetra::Experimental::Classes::BlockCrsMatrix::getLocalRowCopy
void getLocalRowCopy(LO LocalRow, const Teuchos::ArrayView< LO > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const
Not implemented.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:1539
Tpetra::Experimental::Classes::BlockCrsMatrix::getIndexBase
virtual GO getIndexBase() const
The index base for global indices in this matrix.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3571
Tpetra::Experimental::Classes::BlockCrsMatrix::getNodeMaxNumRowEntries
size_t getNodeMaxNumRowEntries() const
Maximum number of entries in any row of the matrix, on this process.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:841
Tpetra::Experimental::Classes::BlockCrsMatrix::packet_type
char packet_type
Implementation detail; tells.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:149
Tpetra::Experimental::Classes::BlockCrsMatrix
Sparse matrix whose entries are small dense square blocks, all of the same dimensions.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:138
Tpetra::Experimental::Classes::BlockCrsMatrix::localGaussSeidel
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.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:1105
Tpetra::Experimental::Classes::BlockCrsMatrix::setAllToScalar
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to alpha.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:942
Tpetra::Experimental::Classes::BlockCrsMatrix::getComm
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which this matrix is distributed.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3538
Tpetra::Experimental::Classes::BlockCrsMatrix::getColMap
Teuchos::RCP< const map_type > getColMap() const
get the (mesh) map for the columns of this block matrix.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:817
Tpetra::Experimental::Classes::BlockMultiVector::impl_scalar_type
mv_type::impl_scalar_type impl_scalar_type
The implementation type of entries in the matrix.
Definition: Tpetra_Experimental_BlockMultiVector_decl.hpp:173
Tpetra::Experimental::Classes::BlockCrsMatrix::scalar_type
Scalar scalar_type
The type of entries in the matrix (that is, of each entry in each block).
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:156
Tpetra::Experimental::Classes::BlockCrsMatrix::BlockCrsMatrix
BlockCrsMatrix()
Default constructor: Makes an empty block matrix.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:658
Tpetra::Experimental::Classes::BlockCrsMatrix::little_block_type
Kokkos::View< impl_scalar_type **, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > little_block_type
The type used to access nonconst matrix blocks.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:195
Tpetra_CrsMatrix_decl.hpp
Declaration of the Tpetra::CrsMatrix class.
Tpetra::Experimental::Classes::BlockCrsMatrix::impl_scalar_type
typename BMV::impl_scalar_type impl_scalar_type
The implementation type of entries in the matrix.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:164
Tpetra::Experimental::Classes::BlockCrsMatrix::getNodeNumDiags
virtual size_t TPETRA_DEPRECATED getNodeNumDiags() const
Number of diagonal entries in the matrix's graph, on the calling process.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3612
Tpetra::Classes::DistObject
Base class for distributed Tpetra objects that support data redistribution.
Definition: Tpetra_DistObject_decl.hpp:349
Tpetra::Experimental::Classes::BlockCrsMatrix::rightScale
virtual void rightScale(const ::Tpetra::Vector< Scalar, LO, GO, Node > &x)
Scale the RowMatrix on the right with the given Vector x.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3796
Tpetra::Experimental::Classes::BlockCrsMatrix::getNodeNumRows
size_t getNodeNumRows() const
get the local number of block rows
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:833
Tpetra::Experimental::Classes::BlockCrsMatrix::node_type
Node node_type
The Node type.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:174
Tpetra::Experimental::Classes::BlockCrsMatrix::hasColMap
virtual bool hasColMap() const
Whether this matrix has a well-defined column map.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3629
Tpetra::Distributor
Sets up and executes a communication plan for a Tpetra DistObject.
Definition: Tpetra_Distributor.hpp:188
Tpetra::Experimental::Classes::BlockCrsMatrix::getFrobeniusNorm
virtual typename ::Tpetra::RowMatrix< Scalar, LO, GO, Node >::mag_type getFrobeniusNorm() const
The Frobenius norm of the matrix.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3815
Tpetra::Experimental::Classes::BlockCrsMatrix::mv_type
::Tpetra::MultiVector< Scalar, LO, GO, node_type > mv_type
The implementation of MultiVector that this class uses.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:186
Tpetra::Experimental::Classes::BlockCrsMatrix::~BlockCrsMatrix
virtual ~BlockCrsMatrix()
Destructor (declared virtual for memory safety).
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:246
Tpetra::Classes::MultiVector
One or more distributed dense vectors.
Definition: Tpetra_MultiVector_decl.hpp:389
Tpetra::Details::DefaultTypes::local_ordinal_type
int local_ordinal_type
Default value of Scalar template parameter.
Definition: Tpetra_Details_DefaultTypes.hpp:72
Tpetra::Experimental::Classes::BlockCrsMatrix::execution_space
device_type::execution_space execution_space
The Kokkos execution space that this class uses.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:179
Tpetra::Experimental::Classes::BlockCrsMatrix::supportsRowViews
virtual bool supportsRowViews() const
Whether this object implements getLocalRowView() and getGlobalRowView().
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3679
Tpetra::Experimental::Classes::BlockCrsMatrix::getGlobalNumDiags
virtual global_size_t TPETRA_DEPRECATED getGlobalNumDiags() const
Number of diagonal entries in the matrix's graph, over all processes in the matrix's communicator.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3603
Tpetra::Experimental::Classes::BlockCrsMatrix::gaussSeidelCopy
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.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:1216
Tpetra::Experimental::Classes::BlockCrsMatrix::localError
bool localError() const
Whether this object had an error on the calling process.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:598
Tpetra::Experimental::Classes::BlockCrsMatrix::getNodeNumEntries
virtual size_t getNodeNumEntries() const
The local number of stored (structurally nonzero) entries.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3587
Tpetra::Experimental::Classes::BlockCrsMatrix::getGlobalRowCopy
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's entries.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3688
Tpetra::Experimental::Classes::BlockCrsMatrix::getNumEntriesInGlobalRow
virtual size_t getNumEntriesInGlobalRow(GO globalRow) const
The current number of entries on the calling process in the specified global row.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3595
Tpetra::Experimental::Classes::BlockCrsMatrix::getDomainMap
Teuchos::RCP< const map_type > getDomainMap() const
Get the (point) domain Map of this matrix.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:791
Tpetra::Experimental::Classes::BlockCrsMatrix::global_ordinal_type
GO global_ordinal_type
The type of global indices.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:169
Tpetra::Experimental::Classes::BlockCrsMatrix::map_type
::Tpetra::Map< LO, GO, node_type > map_type
The implementation of Map that this class uses.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:184
Tpetra::Experimental::Classes::BlockCrsMatrix::isLocallyIndexed
virtual bool isLocallyIndexed() const
Whether matrix indices are locally indexed.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3655
Tpetra::Experimental::Classes::BlockCrsMatrix::isLowerTriangular
virtual bool TPETRA_DEPRECATED isLowerTriangular() const
Whether the matrix's graph is locally lower triangular.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3637
Tpetra::Experimental::Classes::BlockCrsMatrix::device_type
Node::device_type device_type
The Kokkos::Device specialization that this class uses.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:177
Tpetra::Experimental::Classes::BlockMultiVector
MultiVector for multiple degrees of freedom per mesh point.
Definition: Tpetra_Experimental_BlockMultiVector_decl.hpp:147
Tpetra::Classes::Map
A parallel distribution of indices over processes.
Definition: Tpetra_Map_decl.hpp:247
Tpetra::Experimental::Classes::BlockCrsMatrix::absMaxLocalValuesByOffsets
LO absMaxLocalValuesByOffsets(const LO localRowInd, const ptrdiff_t offsets[], const Scalar vals[], const LO numOffsets) const
Like sumIntoLocalValuesByOffsets, but for the ABSMAX combine mode.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:1637
Tpetra::Experimental::Classes::BlockCrsMatrix::getLocalDiagOffsets
void getLocalDiagOffsets(const Kokkos::View< size_t *, device_type, Kokkos::MemoryUnmanaged > &offsets) const
Get offsets of the diagonal entries in the matrix.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:1058
Tpetra::Experimental::Classes::BlockCrsMatrix::applyBlock
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.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:914
Tpetra::Experimental::Classes::BlockCrsMatrix::getLocalDiagCopy
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:...
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:1284
Tpetra::Experimental::Classes::BlockCrsMatrix::local_ordinal_type
LO local_ordinal_type
The type of local indices.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:167
Tpetra::Experimental::Classes::BlockCrsMatrix::reorderedGaussSeidelCopy
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.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:1234
Tpetra::global_size_t
size_t global_size_t
Global size_t object.
Definition: Tpetra_ConfigDefs.hpp:109
Tpetra::Experimental::Classes::BlockCrsMatrix::hasTransposeApply
bool hasTransposeApply() const
Whether it is valid to apply the transpose or conjugate transpose of this matrix.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:290
Tpetra
Namespace Tpetra contains the class and methods constituting the Tpetra library.
Tpetra::Experimental::Classes::BlockCrsMatrix::sumIntoLocalValues
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,...
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:1414
Tpetra::Experimental::Classes::BlockCrsMatrix::getGlobalNumRows
global_size_t getGlobalNumRows() const
get the global number of block rows
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:825
Tpetra::Experimental::Classes::BlockCrsMatrix::getRangeMap
Teuchos::RCP< const map_type > getRangeMap() const
Get the (point) range Map of this matrix.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:800
Tpetra::Experimental::Classes::BlockCrsMatrix::getGlobalNumCols
virtual global_size_t getGlobalNumCols() const
The global number of columns of this matrix.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3555
Tpetra::Experimental::Classes::BlockCrsMatrix::getLocalRowOffsets
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.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:1565
Tpetra::Experimental::Classes::BlockCrsMatrix::memory_space
device_type::memory_space memory_space
The Kokkos memory space that this class uses.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:181
Tpetra::Experimental::Classes::BlockCrsMatrix::crs_graph_type
::Tpetra::CrsGraph< LO, GO, node_type > crs_graph_type
The implementation of CrsGraph that this class uses.
Definition: Tpetra_Experimental_BlockCrsMatrix_decl.hpp:188
Tpetra::Classes::CrsGraph
A distributed graph accessed by rows (adjacency lists) and stored sparsely.
Definition: Tpetra_CrsGraph_decl.hpp:259
Tpetra::Experimental::Classes::BlockCrsMatrix::isGloballyIndexed
virtual bool isGloballyIndexed() const
Whether matrix indices are globally indexed.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3663
Tpetra::Experimental::Classes::BlockCrsMatrix::getNumEntriesInLocalRow
size_t getNumEntriesInLocalRow(const LO localRowInd) const
Return the number of entries in the given row on the calling process.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:1717
Tpetra::Experimental::Classes::BlockCrsMatrix::getNode
virtual Teuchos::RCP< Node > getNode() const
The Kokkos Node instance.
Definition: Tpetra_Experimental_BlockCrsMatrix_def.hpp:3546
Tpetra::CombineMode
CombineMode
Rule for combining data in an Import or Export.
Definition: Tpetra_CombineMode.hpp:94
Tpetra::Classes::RowMatrix
A read-only, row-oriented interface to a sparse matrix.
Definition: Tpetra_RowMatrix_decl.hpp:85
Tpetra::Classes::Vector
A distributed dense vector.
Definition: Tpetra_Vector_decl.hpp:82