Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tpetra_CrsMatrix_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_CRSMATRIX_DECL_HPP
43 #define TPETRA_CRSMATRIX_DECL_HPP
44 
52 
53 #include "Tpetra_CrsMatrix_fwd.hpp"
54 #include "Tpetra_RowMatrix_decl.hpp"
55 #include "Tpetra_Exceptions.hpp"
56 #include "Tpetra_DistObject.hpp"
57 #include "Tpetra_CrsGraph.hpp"
58 #include "Tpetra_Vector.hpp"
60 
61 // localMultiply is templated on DomainScalar and RangeScalar, so we
62 // have to include this header file here, rather than in the _def
63 // header file, so that we can get KokkosSparse::spmv.
64 #include "KokkosSparse.hpp"
65 // localGaussSeidel and reorderedLocalGaussSeidel are templated on
66 // DomainScalar and RangeScalar, so we have to include this header
67 // file here, rather than in the _def header file, so that we can get
68 // the interfaces to the corresponding local computational kernels.
69 #include "KokkosSparse_sor_sequential_impl.hpp"
70 
71 namespace Tpetra {
72 
73  // Forward declaration for CrsMatrix::swap() test
74  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> class crsMatrix_Swap_Tester;
75 
127  template<class CrsMatrixType>
128  Teuchos::RCP<CrsMatrixType>
129  importAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
130  const Import<typename CrsMatrixType::local_ordinal_type,
131  typename CrsMatrixType::global_ordinal_type,
132  typename CrsMatrixType::node_type>& importer,
133  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
134  typename CrsMatrixType::global_ordinal_type,
135  typename CrsMatrixType::node_type> >& domainMap = Teuchos::null,
136  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
137  typename CrsMatrixType::global_ordinal_type,
138  typename CrsMatrixType::node_type> >& rangeMap = Teuchos::null,
139  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
140 
194  template<class CrsMatrixType>
195  Teuchos::RCP<CrsMatrixType>
196  importAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
197  const Import<typename CrsMatrixType::local_ordinal_type,
198  typename CrsMatrixType::global_ordinal_type,
199  typename CrsMatrixType::node_type>& rowImporter,
200  const Import<typename CrsMatrixType::local_ordinal_type,
201  typename CrsMatrixType::global_ordinal_type,
202  typename CrsMatrixType::node_type>& domainImporter,
203  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
204  typename CrsMatrixType::global_ordinal_type,
205  typename CrsMatrixType::node_type> >& domainMap,
206  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
207  typename CrsMatrixType::global_ordinal_type,
208  typename CrsMatrixType::node_type> >& rangeMap,
209  const Teuchos::RCP<Teuchos::ParameterList>& params);
210 
244  template<class CrsMatrixType>
245  Teuchos::RCP<CrsMatrixType>
246  exportAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
247  const Export<typename CrsMatrixType::local_ordinal_type,
248  typename CrsMatrixType::global_ordinal_type,
249  typename CrsMatrixType::node_type>& exporter,
250  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
251  typename CrsMatrixType::global_ordinal_type,
252  typename CrsMatrixType::node_type> >& domainMap = Teuchos::null,
253  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
254  typename CrsMatrixType::global_ordinal_type,
255  typename CrsMatrixType::node_type> >& rangeMap = Teuchos::null,
256  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
257 
291  template<class CrsMatrixType>
292  Teuchos::RCP<CrsMatrixType>
293  exportAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
294  const Export<typename CrsMatrixType::local_ordinal_type,
295  typename CrsMatrixType::global_ordinal_type,
296  typename CrsMatrixType::node_type>& rowExporter,
297  const Export<typename CrsMatrixType::local_ordinal_type,
298  typename CrsMatrixType::global_ordinal_type,
299  typename CrsMatrixType::node_type>& domainExporter,
300  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
301  typename CrsMatrixType::global_ordinal_type,
302  typename CrsMatrixType::node_type> >& domainMap,
303  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
304  typename CrsMatrixType::global_ordinal_type,
305  typename CrsMatrixType::node_type> >& rangeMap,
306  const Teuchos::RCP<Teuchos::ParameterList>& params);
307 
421  template <class Scalar,
422  class LocalOrdinal,
423  class GlobalOrdinal,
424  class Node>
425  class CrsMatrix :
426  public RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>,
427  public DistObject<char, LocalOrdinal, GlobalOrdinal, Node>
428 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
429  , public ::Tpetra::Details::HasDeprecatedMethods2630_WarningThisClassIsNotForUsers
430 #endif // TPETRA_ENABLE_DEPRECATED_CODE
431  {
432  public:
434 
435 
437  using scalar_type = Scalar;
447  using impl_scalar_type = typename Kokkos::ArithTraits<Scalar>::val_type;
449  using local_ordinal_type = LocalOrdinal;
451  using global_ordinal_type = GlobalOrdinal;
453  using device_type = typename Node::device_type;
455  using execution_space = typename device_type::execution_space;
456 
461  using node_type = Node;
462 
468  using mag_type = typename Kokkos::ArithTraits<impl_scalar_type>::mag_type;
469 
472 
475 
478 
481 
484 
487  using local_matrix_type =
488  KokkosSparse::CrsMatrix<impl_scalar_type,
491  void,
492  typename local_graph_type::size_type>;
493 
494 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
495  typedef typename local_matrix_type::row_map_type t_RowPtrs TPETRA_DEPRECATED;
498  typedef typename local_matrix_type::row_map_type::non_const_type t_RowPtrsNC TPETRA_DEPRECATED;
500  typedef typename local_graph_type::entries_type::non_const_type t_LocalOrdinal_1D TPETRA_DEPRECATED;
502  typedef typename local_matrix_type::values_type t_ValuesType TPETRA_DEPRECATED;
504  typedef local_matrix_type k_local_matrix_type TPETRA_DEPRECATED;
505 #endif // TPETRA_ENABLE_DEPRECATED_CODE
506 
508 
510 
528  CrsMatrix (const Teuchos::RCP<const map_type>& rowMap,
529  const size_t maxNumEntriesPerRow,
530  const ProfileType pftype = TPETRA_DEFAULT_PROFILE_TYPE,
531  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
532 
550  CrsMatrix (const Teuchos::RCP<const map_type>& rowMap,
551  const Teuchos::ArrayView<const size_t>& numEntPerRowToAlloc,
552  const ProfileType pftype = TPETRA_DEFAULT_PROFILE_TYPE,
553  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
554 
555 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
556  TPETRA_DEPRECATED
558  CrsMatrix (const Teuchos::RCP<const map_type>& rowMap,
559  const Teuchos::ArrayRCP<const size_t>& numEntPerRowToAlloc,
560  const ProfileType pftype = DynamicProfile,
561  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
562 #endif // TPETRA_ENABLE_DEPRECATED_CODE
563 
586  CrsMatrix (const Teuchos::RCP<const map_type>& rowMap,
587  const Teuchos::RCP<const map_type>& colMap,
588  const size_t maxNumEntPerRow,
589  const ProfileType pftype = TPETRA_DEFAULT_PROFILE_TYPE,
590  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
591 
614  CrsMatrix (const Teuchos::RCP<const map_type>& rowMap,
615  const Teuchos::RCP<const map_type>& colMap,
616  const Teuchos::ArrayView<const size_t>& numEntPerRowToAlloc,
617  const ProfileType pftype = TPETRA_DEFAULT_PROFILE_TYPE,
618  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
619 
620 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
621  TPETRA_DEPRECATED
623  CrsMatrix (const Teuchos::RCP<const map_type>& rowMap,
624  const Teuchos::RCP<const map_type>& colMap,
625  const Teuchos::ArrayRCP<const size_t>& numEntPerRowToAlloc,
626  const ProfileType pftype = DynamicProfile,
627  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
628 #endif // TPETRA_ENABLE_DEPRECATED_CODE
629 
654  explicit CrsMatrix (const Teuchos::RCP<const crs_graph_type>& graph,
655  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
656 
685  explicit CrsMatrix (const Teuchos::RCP<const crs_graph_type>& graph,
686  const typename local_matrix_type::values_type& values,
687  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
688 
712  CrsMatrix (const Teuchos::RCP<const map_type>& rowMap,
713  const Teuchos::RCP<const map_type>& colMap,
714  const typename local_matrix_type::row_map_type& rowPointers,
715  const typename local_graph_type::entries_type::non_const_type& columnIndices,
716  const typename local_matrix_type::values_type& values,
717  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
718 
742  CrsMatrix (const Teuchos::RCP<const map_type>& rowMap,
743  const Teuchos::RCP<const map_type>& colMap,
744  const Teuchos::ArrayRCP<size_t>& rowPointers,
745  const Teuchos::ArrayRCP<LocalOrdinal>& columnIndices,
746  const Teuchos::ArrayRCP<Scalar>& values,
747  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
748 
769  CrsMatrix (const Teuchos::RCP<const map_type>& rowMap,
770  const Teuchos::RCP<const map_type>& colMap,
771  const local_matrix_type& lclMatrix,
772  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
773 
800  CrsMatrix (const local_matrix_type& lclMatrix,
801  const Teuchos::RCP<const map_type>& rowMap,
802  const Teuchos::RCP<const map_type>& colMap,
803  const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
804  const Teuchos::RCP<const map_type>& rangeMap = Teuchos::null,
805  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
806 
807 
809  // This function in 'Copy' mode is only guaranteed to work correctly for matrices
810  // which are fillComplete.
812  const Teuchos::DataAccess copyOrView);
813 
814  // This friend declaration makes the clone() method work.
815  template <class S2, class LO2, class GO2, class N2>
816  friend class CrsMatrix;
817 
818 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
819  template <class Node2>
844  Teuchos::RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node2> > TPETRA_DEPRECATED
845  clone (const Teuchos::RCP<Node2>& node2,
846  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) const
847  {
848  using Teuchos::Array;
849  using Teuchos::ArrayRCP;
850  using Teuchos::ArrayView;
851  using Teuchos::null;
852  using Teuchos::ParameterList;
853  using Teuchos::RCP;
854  using Teuchos::rcp;
855  using Teuchos::sublist;
858  const char tfecfFuncName[] = "clone";
859 
860  // Get parameter values. Set them initially to their default values.
861  bool fillCompleteClone = true;
862  bool useLocalIndices = this->hasColMap ();
863  ProfileType pftype = StaticProfile;
864  if (! params.is_null ()) {
865  fillCompleteClone = params->get ("fillComplete clone", fillCompleteClone);
866  useLocalIndices = params->get ("Locally indexed clone", useLocalIndices);
867 
868  bool staticProfileClone = true;
869  staticProfileClone = params->get ("Static profile clone", staticProfileClone);
870  pftype = staticProfileClone ? StaticProfile : DynamicProfile;
871  }
872 
873  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(
874  ! this->hasColMap () && useLocalIndices, std::runtime_error,
875  ": You requested that the returned clone have local indices, but the "
876  "the source matrix does not have a column Map yet.");
877 
878  RCP<const Map2> clonedRowMap = this->getRowMap ()->template clone<Node2> (node2);
879 
880  // Get an upper bound on the number of entries per row.
881  RCP<CrsMatrix2> clonedMatrix;
882  ArrayRCP<const size_t> numEntriesPerRow;
883  size_t numEntriesForAll = 0;
884  bool boundSameForAllLocalRows = false;
885  staticGraph_->getNumEntriesPerLocalRowUpperBound (numEntriesPerRow,
886  numEntriesForAll,
887  boundSameForAllLocalRows);
888  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(
889  numEntriesForAll != 0 &&
890  static_cast<size_t> (numEntriesPerRow.size ()) != 0,
891  std::logic_error, ": getNumEntriesPerLocalRowUpperBound returned a "
892  "nonzero numEntriesForAll = " << numEntriesForAll << " , as well as a "
893  "numEntriesPerRow array of nonzero length " << numEntriesPerRow.size ()
894  << ". This should never happen. Please report this bug to the Tpetra "
895  "developers.");
896  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(
897  numEntriesForAll != 0 && ! boundSameForAllLocalRows,
898  std::logic_error, ": getNumEntriesPerLocalRowUpperBound returned a "
899  "nonzero numEntriesForAll = " << numEntriesForAll << " , but claims "
900  "(via its third output value) that the upper bound is not the same for "
901  "all rows. This should never happen. Please report this bug to the "
902  "Tpetra developers.");
903  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(
904  numEntriesPerRow.size () != 0 && boundSameForAllLocalRows,
905  std::logic_error, ": getNumEntriesPerLocalRowUpperBound returned a "
906  "numEntriesPerRow array of nonzero length " << numEntriesPerRow.size ()
907  << ", but claims (via its third output value) that the upper bound is "
908  "not the same for all rows. This should never happen. Please report "
909  "this bug to the Tpetra developers.");
910 
911  RCP<ParameterList> matParams =
912  params.is_null () ? null : sublist (params,"CrsMatrix");
913  if (useLocalIndices) {
914  RCP<const Map2> clonedColMap =
915  this->getColMap ()->template clone<Node2> (node2);
916  if (numEntriesPerRow.is_null ()) {
917  clonedMatrix = rcp (new CrsMatrix2 (clonedRowMap, clonedColMap,
918  numEntriesForAll, pftype,
919  matParams));
920  }
921  else {
922  clonedMatrix = rcp (new CrsMatrix2 (clonedRowMap, clonedColMap,
923  numEntriesPerRow (), pftype,
924  matParams));
925  }
926  }
927  else {
928  if (numEntriesPerRow.is_null ()) {
929  clonedMatrix = rcp (new CrsMatrix2 (clonedRowMap, numEntriesForAll,
930  pftype, matParams));
931  }
932  else {
933  clonedMatrix = rcp (new CrsMatrix2 (clonedRowMap,
934  numEntriesPerRow (),
935  pftype, matParams));
936  }
937  }
938  // done with these
939  numEntriesPerRow = Teuchos::null;
940  numEntriesForAll = 0;
941 
942  if (useLocalIndices) {
943  clonedMatrix->allocateValues (LocalIndices,
944  CrsMatrix2::GraphNotYetAllocated);
945  if (this->isLocallyIndexed ()) {
946  ArrayView<const LocalOrdinal> linds;
947  ArrayView<const Scalar> vals;
948  for (LocalOrdinal lrow = clonedRowMap->getMinLocalIndex ();
949  lrow <= clonedRowMap->getMaxLocalIndex ();
950  ++lrow) {
951  this->getLocalRowView (lrow, linds, vals);
952  if (linds.size ()) {
953  clonedMatrix->insertLocalValues (lrow, linds, vals);
954  }
955  }
956  }
957  else { // this->isGloballyIndexed()
958  Array<LocalOrdinal> linds;
959  Array<Scalar> vals;
960  for (LocalOrdinal lrow = clonedRowMap->getMinLocalIndex ();
961  lrow <= clonedRowMap->getMaxLocalIndex ();
962  ++lrow) {
963  size_t theNumEntries = this->getNumEntriesInLocalRow (lrow);
964  if (theNumEntries > static_cast<size_t> (linds.size ())) {
965  linds.resize (theNumEntries);
966  }
967  if (theNumEntries > static_cast<size_t> (vals.size ())) {
968  vals.resize (theNumEntries);
969  }
970  this->getLocalRowCopy (clonedRowMap->getGlobalElement (lrow),
971  linds (), vals (), theNumEntries);
972  if (theNumEntries != 0) {
973  clonedMatrix->insertLocalValues (lrow, linds (0, theNumEntries),
974  vals (0, theNumEntries));
975  }
976  }
977  }
978  }
979  else { // useGlobalIndices
980  clonedMatrix->allocateValues (GlobalIndices,
981  CrsMatrix2::GraphNotYetAllocated);
982  if (this->isGloballyIndexed ()) {
983  ArrayView<const GlobalOrdinal> ginds;
984  ArrayView<const Scalar> vals;
985  for (GlobalOrdinal grow = clonedRowMap->getMinGlobalIndex ();
986  grow <= clonedRowMap->getMaxGlobalIndex ();
987  ++grow) {
988  this->getGlobalRowView (grow, ginds, vals);
989  if (ginds.size () > 0) {
990  clonedMatrix->insertGlobalValues (grow, ginds, vals);
991  }
992  }
993  }
994  else { // this->isLocallyIndexed()
995  Array<GlobalOrdinal> ginds;
996  Array<Scalar> vals;
997  for (GlobalOrdinal grow = clonedRowMap->getMinGlobalIndex ();
998  grow <= clonedRowMap->getMaxGlobalIndex ();
999  ++grow) {
1000  size_t theNumEntries = this->getNumEntriesInGlobalRow (grow);
1001  if (theNumEntries > static_cast<size_t> (ginds.size ())) {
1002  ginds.resize (theNumEntries);
1003  }
1004  if (theNumEntries > static_cast<size_t> (vals.size ())) {
1005  vals.resize (theNumEntries);
1006  }
1007  this->getGlobalRowCopy (grow, ginds (), vals (), theNumEntries);
1008  if (theNumEntries != 0) {
1009  clonedMatrix->insertGlobalValues (grow, ginds (0, theNumEntries),
1010  vals (0, theNumEntries));
1011  }
1012  }
1013  }
1014  }
1015 
1016  if (fillCompleteClone) {
1017  RCP<const Map2> clonedRangeMap;
1018  RCP<const Map2> clonedDomainMap;
1019  try {
1020  if (! this->getRangeMap ().is_null () &&
1021  this->getRangeMap () != clonedRowMap) {
1022  clonedRangeMap = this->getRangeMap ()->template clone<Node2> (node2);
1023  }
1024  else {
1025  clonedRangeMap = clonedRowMap;
1026  }
1027  if (! this->getDomainMap ().is_null () &&
1028  this->getDomainMap () != clonedRowMap) {
1029  clonedDomainMap = this->getDomainMap ()->template clone<Node2> (node2);
1030  }
1031  else {
1032  clonedDomainMap = clonedRowMap;
1033  }
1034  }
1035  catch (std::exception &e) {
1036  const bool caughtExceptionOnClone = true;
1037  TEUCHOS_TEST_FOR_EXCEPTION
1038  (caughtExceptionOnClone, std::runtime_error,
1039  Teuchos::typeName (*this) << "::clone: Caught the following "
1040  "exception while cloning range and domain Maps on a clone of "
1041  "type " << Teuchos::typeName (*clonedMatrix) << ": " << e.what ());
1042  }
1043 
1044  RCP<ParameterList> fillparams =
1045  params.is_null () ? Teuchos::null : sublist (params, "fillComplete");
1046  try {
1047  clonedMatrix->fillComplete (clonedDomainMap, clonedRangeMap,
1048  fillparams);
1049  }
1050  catch (std::exception &e) {
1051  const bool caughtExceptionOnClone = true;
1052  TEUCHOS_TEST_FOR_EXCEPTION(
1053  caughtExceptionOnClone, std::runtime_error,
1054  Teuchos::typeName (*this) << "::clone: Caught the following "
1055  "exception while calling fillComplete() on a clone of type "
1056  << Teuchos::typeName (*clonedMatrix) << ": " << e.what ());
1057  }
1058  }
1059  return clonedMatrix;
1060  }
1061 #endif
1062 
1071  CrsMatrix (const CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>&) = delete;
1072 
1080  CrsMatrix (CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>&&) = delete;
1081 
1090  CrsMatrix&
1091  operator= (const CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>&) = delete;
1092 
1100  CrsMatrix&
1101  operator= (CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>&&) = delete;
1102 
1112  virtual ~CrsMatrix () = default;
1113 
1114  public:
1116 
1118 
1142  //
1188  void
1189  insertGlobalValues (const GlobalOrdinal globalRow,
1190  const Teuchos::ArrayView<const GlobalOrdinal>& cols,
1191  const Teuchos::ArrayView<const Scalar>& vals);
1192 
1207  void
1208  insertGlobalValues (const GlobalOrdinal globalRow,
1209  const LocalOrdinal numEnt,
1210  const Scalar vals[],
1211  const GlobalOrdinal inds[]);
1212 
1253  void
1254  insertLocalValues (const LocalOrdinal localRow,
1255  const Teuchos::ArrayView<const LocalOrdinal> &cols,
1256  const Teuchos::ArrayView<const Scalar> &vals);
1257 
1272  void
1273  insertLocalValues (const LocalOrdinal localRow,
1274  const LocalOrdinal numEnt,
1275  const Scalar vals[],
1276  const LocalOrdinal cols[]);
1277 
1278  private:
1289  LocalOrdinal
1290  replaceGlobalValuesImpl (impl_scalar_type rowVals[],
1291  const crs_graph_type& graph,
1292  const RowInfo& rowInfo,
1293  const GlobalOrdinal inds[],
1294  const impl_scalar_type newVals[],
1295  const LocalOrdinal numElts) const;
1296 
1297  public:
1334  template<class GlobalIndicesViewType,
1335  class ImplScalarViewType>
1336  LocalOrdinal
1337  replaceGlobalValues (const GlobalOrdinal globalRow,
1338  const typename UnmanagedView<GlobalIndicesViewType>::type& inputInds,
1339  const typename UnmanagedView<ImplScalarViewType>::type& inputVals) const
1340  {
1341  // We use static_assert here to check the template parameters,
1342  // rather than std::enable_if (e.g., on the return value, to
1343  // enable compilation only if the template parameters match the
1344  // desired attributes). This turns obscure link errors into
1345  // clear compilation errors. It also makes the return value a
1346  // lot easier to see.
1347  static_assert (Kokkos::is_view<GlobalIndicesViewType>::value,
1348  "First template parameter GlobalIndicesViewType must be "
1349  "a Kokkos::View.");
1350  static_assert (Kokkos::is_view<ImplScalarViewType>::value,
1351  "Second template parameter ImplScalarViewType must be a "
1352  "Kokkos::View.");
1353  static_assert (static_cast<int> (GlobalIndicesViewType::rank) == 1,
1354  "First template parameter GlobalIndicesViewType must "
1355  "have rank 1.");
1356  static_assert (static_cast<int> (ImplScalarViewType::rank) == 1,
1357  "Second template parameter ImplScalarViewType must have "
1358  "rank 1.");
1359  static_assert (std::is_same<
1360  typename GlobalIndicesViewType::non_const_value_type,
1361  global_ordinal_type>::value,
1362  "First template parameter GlobalIndicesViewType must "
1363  "contain values of type global_ordinal_type.");
1364  static_assert (std::is_same<
1365  typename ImplScalarViewType::non_const_value_type,
1366  impl_scalar_type>::value,
1367  "Second template parameter ImplScalarViewType must "
1368  "contain values of type impl_scalar_type.");
1369  typedef LocalOrdinal LO;
1370  const LO numInputEnt = inputInds.extent (0);
1371  if (static_cast<LO> (inputVals.extent (0)) != numInputEnt) {
1372  return Teuchos::OrdinalTraits<LO>::invalid ();
1373  }
1374  const Scalar* const inVals =
1375  reinterpret_cast<const Scalar*> (inputVals.data ());
1376  return this->replaceGlobalValues (globalRow, numInputEnt, inVals,
1377  inputInds.data ());
1378  }
1379 
1383  LocalOrdinal
1384  replaceGlobalValues (const GlobalOrdinal globalRow,
1385  const Teuchos::ArrayView<const GlobalOrdinal>& cols,
1386  const Teuchos::ArrayView<const Scalar>& vals) const;
1387 
1403  LocalOrdinal
1404  replaceGlobalValues (const GlobalOrdinal globalRow,
1405  const LocalOrdinal numEnt,
1406  const Scalar vals[],
1407  const GlobalOrdinal cols[]) const;
1408 
1409  private:
1420  LocalOrdinal
1421  replaceLocalValuesImpl (impl_scalar_type rowVals[],
1422  const crs_graph_type& graph,
1423  const RowInfo& rowInfo,
1424  const LocalOrdinal inds[],
1425  const impl_scalar_type newVals[],
1426  const LocalOrdinal numElts) const;
1427 
1428  public:
1464  template<class LocalIndicesViewType,
1465  class ImplScalarViewType>
1466  LocalOrdinal
1467  replaceLocalValues (const LocalOrdinal localRow,
1468  const typename UnmanagedView<LocalIndicesViewType>::type& inputInds,
1469  const typename UnmanagedView<ImplScalarViewType>::type& inputVals) const
1470  {
1471  // We use static_assert here to check the template parameters,
1472  // rather than std::enable_if (e.g., on the return value, to
1473  // enable compilation only if the template parameters match the
1474  // desired attributes). This turns obscure link errors into
1475  // clear compilation errors. It also makes the return value a
1476  // lot easier to see.
1477  static_assert (Kokkos::is_view<LocalIndicesViewType>::value,
1478  "First template parameter LocalIndicesViewType must be "
1479  "a Kokkos::View.");
1480  static_assert (Kokkos::is_view<ImplScalarViewType>::value,
1481  "Second template parameter ImplScalarViewType must be a "
1482  "Kokkos::View.");
1483  static_assert (static_cast<int> (LocalIndicesViewType::rank) == 1,
1484  "First template parameter LocalIndicesViewType must "
1485  "have rank 1.");
1486  static_assert (static_cast<int> (ImplScalarViewType::rank) == 1,
1487  "Second template parameter ImplScalarViewType must have "
1488  "rank 1.");
1489  static_assert (std::is_same<
1490  typename LocalIndicesViewType::non_const_value_type,
1491  local_ordinal_type>::value,
1492  "First template parameter LocalIndicesViewType must "
1493  "contain values of type local_ordinal_type.");
1494  static_assert (std::is_same<
1495  typename ImplScalarViewType::non_const_value_type,
1496  impl_scalar_type>::value,
1497  "Second template parameter ImplScalarViewType must "
1498  "contain values of type impl_scalar_type.");
1499 
1500  typedef LocalOrdinal LO;
1501  const LO numInputEnt = inputInds.extent (0);
1502  if (numInputEnt != inputVals.extent (0)) {
1503  return Teuchos::OrdinalTraits<LO>::invalid ();
1504  }
1505  const Scalar* const inVals =
1506  reinterpret_cast<const Scalar*> (inputVals.data ());
1507  return this->replaceLocalValues (localRow, numInputEnt,
1508  inVals, inputInds.data ());
1509  }
1510 
1514  LocalOrdinal
1515  replaceLocalValues (const LocalOrdinal localRow,
1516  const Teuchos::ArrayView<const LocalOrdinal>& cols,
1517  const Teuchos::ArrayView<const Scalar>& vals) const;
1518 
1536  LocalOrdinal
1537  replaceLocalValues (const LocalOrdinal localRow,
1538  const LocalOrdinal numEnt,
1539  const Scalar inputVals[],
1540  const LocalOrdinal inputCols[]) const;
1541 
1542  private:
1547  static const bool useAtomicUpdatesByDefault =
1548 #ifdef KOKKOS_ENABLE_SERIAL
1549  ! std::is_same<execution_space, Kokkos::Serial>::value;
1550 #else
1551  true;
1552 #endif // KOKKOS_ENABLE_SERIAL
1553 
1577  LocalOrdinal
1578  sumIntoGlobalValuesImpl (impl_scalar_type rowVals[],
1579  const crs_graph_type& graph,
1580  const RowInfo& rowInfo,
1581  const GlobalOrdinal inds[],
1582  const impl_scalar_type newVals[],
1583  const LocalOrdinal numElts,
1584  const bool atomic = useAtomicUpdatesByDefault) const;
1585 
1586  public:
1623  LocalOrdinal
1624  sumIntoGlobalValues (const GlobalOrdinal globalRow,
1625  const Teuchos::ArrayView<const GlobalOrdinal>& cols,
1626  const Teuchos::ArrayView<const Scalar>& vals,
1627  const bool atomic = useAtomicUpdatesByDefault);
1628 
1651  LocalOrdinal
1652  sumIntoGlobalValues (const GlobalOrdinal globalRow,
1653  const LocalOrdinal numEnt,
1654  const Scalar vals[],
1655  const GlobalOrdinal cols[],
1656  const bool atomic = useAtomicUpdatesByDefault);
1657 
1658  private:
1671  LocalOrdinal
1672  sumIntoLocalValuesImpl (impl_scalar_type rowVals[],
1673  const crs_graph_type& graph,
1674  const RowInfo& rowInfo,
1675  const LocalOrdinal inds[],
1676  const impl_scalar_type newVals[],
1677  const LocalOrdinal numElts,
1678  const bool atomic = useAtomicUpdatesByDefault) const;
1679 
1680  public:
1717  template<class LocalIndicesViewType,
1718  class ImplScalarViewType>
1719  LocalOrdinal
1720  sumIntoLocalValues (const LocalOrdinal localRow,
1721  const typename UnmanagedView<LocalIndicesViewType>::type& inputInds,
1722  const typename UnmanagedView<ImplScalarViewType>::type& inputVals,
1723  const bool atomic = useAtomicUpdatesByDefault) const
1724  {
1725  // We use static_assert here to check the template parameters,
1726  // rather than std::enable_if (e.g., on the return value, to
1727  // enable compilation only if the template parameters match the
1728  // desired attributes). This turns obscure link errors into
1729  // clear compilation errors. It also makes the return value a
1730  // lot easier to see.
1731  static_assert (Kokkos::is_view<LocalIndicesViewType>::value,
1732  "First template parameter LocalIndicesViewType must be "
1733  "a Kokkos::View.");
1734  static_assert (Kokkos::is_view<ImplScalarViewType>::value,
1735  "Second template parameter ImplScalarViewType must be a "
1736  "Kokkos::View.");
1737  static_assert (static_cast<int> (LocalIndicesViewType::rank) == 1,
1738  "First template parameter LocalIndicesViewType must "
1739  "have rank 1.");
1740  static_assert (static_cast<int> (ImplScalarViewType::rank) == 1,
1741  "Second template parameter ImplScalarViewType must have "
1742  "rank 1.");
1743  static_assert (std::is_same<
1744  typename LocalIndicesViewType::non_const_value_type,
1745  local_ordinal_type>::value,
1746  "First template parameter LocalIndicesViewType must "
1747  "contain values of type local_ordinal_type.");
1748  static_assert (std::is_same<
1749  typename ImplScalarViewType::non_const_value_type,
1750  impl_scalar_type>::value,
1751  "Second template parameter ImplScalarViewType must "
1752  "contain values of type impl_scalar_type.");
1753  typedef LocalOrdinal LO;
1754  const LO numInputEnt = inputInds.extent (0);
1755  if (static_cast<LO> (inputVals.extent (0)) != numInputEnt) {
1756  return Teuchos::OrdinalTraits<LO>::invalid ();
1757  }
1758  return this->sumIntoLocalValues (localRow,
1759  numInputEnt,
1760  reinterpret_cast<const Scalar*> (inputVals.data ()),
1761  inputInds.data (),
1762  atomic);
1763  }
1764 
1794  LocalOrdinal
1795  sumIntoLocalValues (const LocalOrdinal localRow,
1796  const Teuchos::ArrayView<const LocalOrdinal>& cols,
1797  const Teuchos::ArrayView<const Scalar>& vals,
1798  const bool atomic = useAtomicUpdatesByDefault) const;
1799 
1821  LocalOrdinal
1822  sumIntoLocalValues (const LocalOrdinal localRow,
1823  const LocalOrdinal numEnt,
1824  const Scalar vals[],
1825  const LocalOrdinal cols[],
1826  const bool atomic = useAtomicUpdatesByDefault) const;
1827 
1828  private:
1859  LocalOrdinal
1860  transformLocalValues (impl_scalar_type rowVals[],
1861  const crs_graph_type& graph,
1862  const RowInfo& rowInfo,
1863  const LocalOrdinal inds[],
1864  const impl_scalar_type newVals[],
1865  const LocalOrdinal numElts,
1866  std::function<impl_scalar_type (const impl_scalar_type&, const impl_scalar_type&) > f,
1867  const bool atomic = useAtomicUpdatesByDefault) const;
1868 
1899  LocalOrdinal
1900  transformGlobalValues (impl_scalar_type rowVals[],
1901  const crs_graph_type& graph,
1902  const RowInfo& rowInfo,
1903  const GlobalOrdinal inds[],
1904  const impl_scalar_type newVals[],
1905  const LocalOrdinal numElts,
1906  std::function<impl_scalar_type (const impl_scalar_type&, const impl_scalar_type&) > f,
1907  const bool atomic = useAtomicUpdatesByDefault) const;
1908 
1935  LocalOrdinal
1936  transformLocalValues (const LocalOrdinal lclRow,
1937  const LocalOrdinal numInputEnt,
1938  const impl_scalar_type inputVals[],
1939  const LocalOrdinal inputCols[],
1940  std::function<impl_scalar_type (const impl_scalar_type&, const impl_scalar_type&) > f,
1941  const bool atomic = useAtomicUpdatesByDefault) const;
1942 
1969  LocalOrdinal
1970  transformGlobalValues (const GlobalOrdinal gblRow,
1971  const LocalOrdinal numInputEnt,
1972  const impl_scalar_type inputVals[],
1973  const GlobalOrdinal inputCols[],
1974  std::function<impl_scalar_type (const impl_scalar_type&, const impl_scalar_type&) > f,
1975  const bool atomic = useAtomicUpdatesByDefault) const;
1976 
1977  public:
2021  template<class LocalIndicesViewType,
2022  class ImplScalarViewType,
2023  class BinaryFunction>
2024  LocalOrdinal
2025  transformLocalValues (const LocalOrdinal lclRow,
2026  const typename UnmanagedView<LocalIndicesViewType>::type& inputInds,
2027  const typename UnmanagedView<ImplScalarViewType>::type& inputVals,
2028  BinaryFunction f,
2029  const bool atomic = useAtomicUpdatesByDefault) const
2030  {
2031  // We use static_assert here to check the template parameters,
2032  // rather than std::enable_if (e.g., on the return value, to
2033  // enable compilation only if the template parameters match the
2034  // desired attributes). This turns obscure link errors into
2035  // clear compilation errors. It also makes the return value a
2036  // lot easier to see.
2037  static_assert (Kokkos::is_view<LocalIndicesViewType>::value,
2038  "First template parameter LocalIndicesViewType must be "
2039  "a Kokkos::View.");
2040  static_assert (Kokkos::is_view<ImplScalarViewType>::value,
2041  "Second template parameter ImplScalarViewType must be a "
2042  "Kokkos::View.");
2043  static_assert (static_cast<int> (LocalIndicesViewType::rank) == 1,
2044  "First template parameter LocalIndicesViewType must "
2045  "have rank 1.");
2046  static_assert (static_cast<int> (ImplScalarViewType::rank) == 1,
2047  "Second template parameter ImplScalarViewType must have "
2048  "rank 1.");
2049  static_assert (std::is_same<
2050  typename LocalIndicesViewType::non_const_value_type,
2051  local_ordinal_type>::value,
2052  "First template parameter LocalIndicesViewType must "
2053  "contain values of type local_ordinal_type.");
2054  static_assert (std::is_same<
2055  typename ImplScalarViewType::non_const_value_type,
2056  impl_scalar_type>::value,
2057  "Second template parameter ImplScalarViewType must "
2058  "contain values of type impl_scalar_type.");
2059  typedef LocalOrdinal LO;
2060  const LO numInputEnt = inputInds.extent (0);
2061  if (static_cast<LO> (inputVals.extent (0)) != numInputEnt) {
2062  return Teuchos::OrdinalTraits<LO>::invalid ();
2063  }
2064  return this->transformLocalValues (lclRow,
2065  numInputEnt,
2066  inputVals.data (),
2067  inputInds.data (),
2068  f,
2069  atomic);
2070  }
2071 
2113  template<class BinaryFunction, class InputMemorySpace>
2114  LocalOrdinal
2115  transformGlobalValues (const GlobalOrdinal gblRow,
2116  const Kokkos::View<const GlobalOrdinal*,
2117  InputMemorySpace,
2118  Kokkos::MemoryUnmanaged>& inputInds,
2119  const Kokkos::View<const impl_scalar_type*,
2120  InputMemorySpace,
2121  Kokkos::MemoryUnmanaged>& inputVals,
2122  BinaryFunction f,
2123  const bool atomic = useAtomicUpdatesByDefault) const
2124  {
2125  typedef LocalOrdinal LO;
2126  const LO numInputEnt = inputInds.extent (0);
2127  if (static_cast<LO> (inputVals.extent (0)) != numInputEnt) {
2128  return Teuchos::OrdinalTraits<LO>::invalid ();
2129  }
2130  return this->transformGlobalValues (gblRow,
2131  numInputEnt,
2132  inputVals.data (),
2133  inputInds.data (),
2134  f,
2135  atomic);
2136  }
2137 
2139  void setAllToScalar (const Scalar& alpha);
2140 
2142  void scale (const Scalar& alpha);
2143 
2167  void
2168  setAllValues (const typename local_matrix_type::row_map_type& ptr,
2169  const typename local_graph_type::entries_type::non_const_type& ind,
2170  const typename local_matrix_type::values_type& val);
2171 
2195  void
2196  setAllValues (const Teuchos::ArrayRCP<size_t>& ptr,
2197  const Teuchos::ArrayRCP<LocalOrdinal>& ind,
2198  const Teuchos::ArrayRCP<Scalar>& val);
2199 
2200  void
2201  getAllValues (Teuchos::ArrayRCP<const size_t>& rowPointers,
2202  Teuchos::ArrayRCP<const LocalOrdinal>& columnIndices,
2203  Teuchos::ArrayRCP<const Scalar>& values) const;
2204 
2206 
2208 
2237  void globalAssemble();
2238 
2252  void resumeFill (const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
2253 
2311  void
2312  fillComplete (const Teuchos::RCP<const map_type>& domainMap,
2313  const Teuchos::RCP<const map_type>& rangeMap,
2314  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
2315 
2342  void
2343  fillComplete (const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
2344 
2371  void
2372  expertStaticFillComplete (const Teuchos::RCP<const map_type>& domainMap,
2373  const Teuchos::RCP<const map_type>& rangeMap,
2374  const Teuchos::RCP<const import_type>& importer = Teuchos::null,
2375  const Teuchos::RCP<const export_type>& exporter = Teuchos::null,
2376  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
2377 
2387  void
2388  replaceColMap (const Teuchos::RCP<const map_type>& newColMap);
2389 
2471  void
2472  reindexColumns (crs_graph_type* const graph,
2473  const Teuchos::RCP<const map_type>& newColMap,
2474  const Teuchos::RCP<const import_type>& newImport = Teuchos::null,
2475  const bool sortEachRow = true);
2476 
2489  void
2490  replaceDomainMapAndImporter (const Teuchos::RCP<const map_type>& newDomainMap,
2491  Teuchos::RCP<const import_type>& newImporter);
2492 
2506  virtual void
2507  removeEmptyProcessesInPlace (const Teuchos::RCP<const map_type>& newMap) override;
2508 
2510 
2512 
2514  Teuchos::RCP<const Teuchos::Comm<int> > getComm() const override;
2515 
2516 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
2517  TPETRA_DEPRECATED Teuchos::RCP<node_type> getNode () const override;
2519 #endif // TPETRA_ENABLE_DEPRECATED_CODE
2520 
2522  Teuchos::RCP<const map_type> getRowMap () const override;
2523 
2525  Teuchos::RCP<const map_type> getColMap () const override;
2526 
2528  Teuchos::RCP<const RowGraph<LocalOrdinal, GlobalOrdinal, Node> >
2529  getGraph () const override;
2530 
2532  Teuchos::RCP<const crs_graph_type> getCrsGraph () const;
2533 
2534  private:
2545  const crs_graph_type& getCrsGraphRef () const;
2546 
2547  public:
2558 
2578  global_size_t getGlobalNumRows() const override;
2579 
2585  global_size_t getGlobalNumCols() const override;
2586 
2593  size_t getNodeNumRows() const override;
2594 
2598  size_t getNodeNumCols() const override;
2599 
2601  GlobalOrdinal getIndexBase() const override;
2602 
2604  global_size_t getGlobalNumEntries() const override;
2605 
2607  size_t getNodeNumEntries() const override;
2608 
2615  size_t getNumEntriesInGlobalRow (GlobalOrdinal globalRow) const override;
2616 
2623  size_t getNumEntriesInLocalRow (local_ordinal_type localRow) const override;
2624 
2632  size_t getGlobalMaxNumRowEntries () const override;
2633 
2641  size_t getNodeMaxNumRowEntries () const override;
2642 
2644  bool hasColMap () const override;
2645 
2646 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
2647  global_size_t TPETRA_DEPRECATED getGlobalNumDiags () const override;
2655 
2663  size_t TPETRA_DEPRECATED getNodeNumDiags () const override;
2664 
2676  global_size_t getGlobalNumDiagsImpl () const override;
2677 
2689  size_t getNodeNumDiagsImpl () const override;
2690 
2701  bool TPETRA_DEPRECATED isLowerTriangular () const override;
2702 
2713  bool TPETRA_DEPRECATED isUpperTriangular () const override;
2714 
2726  bool isLowerTriangularImpl () const override;
2727 
2739  bool isUpperTriangularImpl () const override;
2740 #endif // TPETRA_ENABLE_DEPRECATED_CODE
2741 
2762  bool isLocallyIndexed() const override;
2763 
2784  bool isGloballyIndexed() const override;
2785 
2808  bool isFillComplete() const override;
2809 
2832  bool isFillActive() const;
2833 
2835 
2841  bool isStorageOptimized () const;
2842 
2844  ProfileType getProfileType () const;
2845 
2847  bool isStaticGraph () const;
2848 
2860  mag_type getFrobeniusNorm () const override;
2861 
2864  virtual bool supportsRowViews () const override;
2865 
2914  void
2915  getGlobalRowCopy (GlobalOrdinal GlobalRow,
2916  const Teuchos::ArrayView<GlobalOrdinal>& Indices,
2917  const Teuchos::ArrayView<Scalar>& Values,
2918  size_t& NumEntries) const override;
2919 
2935  void
2936  getLocalRowCopy (LocalOrdinal localRow,
2937  const Teuchos::ArrayView<LocalOrdinal>& colInds,
2938  const Teuchos::ArrayView<Scalar>& vals,
2939  size_t& numEntries) const override;
2940 
2953  void
2954  getGlobalRowView (GlobalOrdinal GlobalRow,
2955  Teuchos::ArrayView<const GlobalOrdinal>& indices,
2956  Teuchos::ArrayView<const Scalar>& values) const override;
2957 
2970  void
2971  getLocalRowView (LocalOrdinal LocalRow,
2972  Teuchos::ArrayView<const LocalOrdinal>& indices,
2973  Teuchos::ArrayView<const Scalar>& values) const override;
2974 
2999  LocalOrdinal
3000  getLocalRowViewRaw (const LocalOrdinal lclRow,
3001  LocalOrdinal& numEnt,
3002  const LocalOrdinal*& lclColInds,
3003  const Scalar*& vals) const override;
3004 
3028  LocalOrdinal
3029  getLocalRowView (const LocalOrdinal lclRow,
3030  LocalOrdinal& numEnt,
3031  const impl_scalar_type*& val,
3032  const LocalOrdinal*& ind) const;
3033 
3041  template<class OutputScalarType>
3042  typename std::enable_if<! std::is_same<OutputScalarType, impl_scalar_type>::value &&
3043  std::is_convertible<impl_scalar_type, OutputScalarType>::value,
3044  LocalOrdinal>::type
3045  getLocalRowView (const LocalOrdinal lclRow,
3046  LocalOrdinal& numEnt,
3047  const OutputScalarType*& val,
3048  const LocalOrdinal*& ind) const
3049  {
3050  const impl_scalar_type* valTmp = NULL;
3051  const LocalOrdinal err = this->getLocalRowView (lclRow, numEnt, valTmp, ind);
3052  // Cast is legitimate because impl_scalar_type is convertible to
3053  // OutputScalarType.
3054  val = reinterpret_cast<const OutputScalarType*> (valTmp);
3055  return err;
3056  }
3057 
3064  void
3066 
3110  void getLocalDiagOffsets (Teuchos::ArrayRCP<size_t>& offsets) const;
3111 
3133  void
3135  const Kokkos::View<const size_t*, device_type,
3136  Kokkos::MemoryUnmanaged>& offsets) const;
3137 
3160  void
3162  const Teuchos::ArrayView<const size_t>& offsets) const;
3163 
3168  void
3170 
3175  void
3177 
3179 
3181 
3230  template <class DomainScalar, class RangeScalar>
3231  void
3234  Teuchos::ETransp mode,
3235  RangeScalar alpha,
3236  RangeScalar beta) const
3237  {
3238  using Teuchos::NO_TRANS;
3239  // Just like Scalar and impl_scalar_type may differ in CrsMatrix,
3240  // RangeScalar and its corresponding impl_scalar_type may differ in
3241  // MultiVector.
3242  typedef typename MultiVector<RangeScalar, LocalOrdinal, GlobalOrdinal,
3243  Node>::impl_scalar_type range_impl_scalar_type;
3244 #ifdef HAVE_TPETRA_DEBUG
3245  const char tfecfFuncName[] = "localMultiply: ";
3246 #endif // HAVE_TPETRA_DEBUG
3247 
3248  const range_impl_scalar_type theAlpha = static_cast<range_impl_scalar_type> (alpha);
3249  const range_impl_scalar_type theBeta = static_cast<range_impl_scalar_type> (beta);
3250  const bool conjugate = (mode == Teuchos::CONJ_TRANS);
3251  const bool transpose = (mode != Teuchos::NO_TRANS);
3252  auto X_lcl = X.template getLocalView<device_type> ();
3253  auto Y_lcl = Y.template getLocalView<device_type> ();
3254 
3255 #ifdef HAVE_TPETRA_DEBUG
3256  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC
3257  (X.getNumVectors () != Y.getNumVectors (), std::runtime_error,
3258  "X.getNumVectors() = " << X.getNumVectors () << " != Y.getNumVectors() = "
3259  << Y.getNumVectors () << ".");
3260  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC
3261  (! transpose && X.getLocalLength () != getColMap ()->getNodeNumElements (),
3262  std::runtime_error, "NO_TRANS case: X has the wrong number of local rows. "
3263  "X.getLocalLength() = " << X.getLocalLength () << " != getColMap()->"
3264  "getNodeNumElements() = " << getColMap ()->getNodeNumElements () << ".");
3265  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC
3266  (! transpose && Y.getLocalLength () != getRowMap ()->getNodeNumElements (),
3267  std::runtime_error, "NO_TRANS case: Y has the wrong number of local rows. "
3268  "Y.getLocalLength() = " << Y.getLocalLength () << " != getRowMap()->"
3269  "getNodeNumElements() = " << getRowMap ()->getNodeNumElements () << ".");
3270  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC
3271  (transpose && X.getLocalLength () != getRowMap ()->getNodeNumElements (),
3272  std::runtime_error, "TRANS or CONJ_TRANS case: X has the wrong number of "
3273  "local rows. X.getLocalLength() = " << X.getLocalLength () << " != "
3274  "getRowMap()->getNodeNumElements() = "
3275  << getRowMap ()->getNodeNumElements () << ".");
3276  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC
3277  (transpose && Y.getLocalLength () != getColMap ()->getNodeNumElements (),
3278  std::runtime_error, "TRANS or CONJ_TRANS case: X has the wrong number of "
3279  "local rows. Y.getLocalLength() = " << Y.getLocalLength () << " != "
3280  "getColMap()->getNodeNumElements() = "
3281  << getColMap ()->getNodeNumElements () << ".");
3282  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC
3283  (! isFillComplete (), std::runtime_error, "The matrix is not fill "
3284  "complete. You must call fillComplete() (possibly with domain and range "
3285  "Map arguments) without an intervening resumeFill() call before you may "
3286  "call this method.");
3287  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC
3288  (! X.isConstantStride () || ! Y.isConstantStride (), std::runtime_error,
3289  "X and Y must be constant stride.");
3290  // If the two pointers are NULL, then they don't alias one
3291  // another, even though they are equal.
3292  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(
3293  X_lcl.data () == Y_lcl.data () &&
3294  X_lcl.data () != NULL,
3295  std::runtime_error, "X and Y may not alias one another.");
3296 #endif // HAVE_TPETRA_DEBUG
3297 
3298  // Y = alpha*op(M) + beta*Y
3299  if (transpose) {
3300  KokkosSparse::spmv (conjugate ? KokkosSparse::ConjugateTranspose : KokkosSparse::Transpose,
3301  theAlpha,
3302  lclMatrix_,
3303  X.template getLocalView<device_type> (),
3304  theBeta,
3305  Y.template getLocalView<device_type> ());
3306  }
3307  else {
3308  KokkosSparse::spmv (KokkosSparse::NoTranspose,
3309  theAlpha,
3310  lclMatrix_,
3311  X.template getLocalView<device_type> (),
3312  theBeta,
3313  Y.template getLocalView<device_type> ());
3314  }
3315  }
3316 
3317  public:
3348  void
3351  const Teuchos::ETransp mode = Teuchos::NO_TRANS,
3352  const Scalar& alpha = Teuchos::ScalarTraits<Scalar>::one (),
3353  const Scalar& beta = Teuchos::ScalarTraits<Scalar>::zero ()) const;
3354 
3379  template <class DomainScalar, class RangeScalar>
3380  void
3384  const RangeScalar& dampingFactor,
3385  const ESweepDirection direction) const
3386  {
3387  typedef LocalOrdinal LO;
3388  typedef GlobalOrdinal GO;
3391  typedef typename Node::device_type::memory_space dev_mem_space;
3392  typedef typename MMV::dual_view_type::t_host::device_type host_mem_space;
3393  typedef typename Graph::local_graph_type k_local_graph_type;
3394  typedef typename k_local_graph_type::size_type offset_type;
3395  const char prefix[] = "Tpetra::CrsMatrix::localGaussSeidel: ";
3396 
3397  TEUCHOS_TEST_FOR_EXCEPTION
3398  (! this->isFillComplete (), std::runtime_error,
3399  prefix << "The matrix is not fill complete.");
3400  const size_t lclNumRows = this->getNodeNumRows ();
3401  const size_t numVecs = B.getNumVectors ();
3402  TEUCHOS_TEST_FOR_EXCEPTION
3403  (X.getNumVectors () != numVecs, std::invalid_argument,
3404  prefix << "B.getNumVectors() = " << numVecs << " != "
3405  "X.getNumVectors() = " << X.getNumVectors () << ".");
3406  TEUCHOS_TEST_FOR_EXCEPTION
3407  (B.getLocalLength () != lclNumRows, std::invalid_argument,
3408  prefix << "B.getLocalLength() = " << B.getLocalLength ()
3409  << " != this->getNodeNumRows() = " << lclNumRows << ".");
3410 
3411  // mfh 28 Aug 2017: The current local Gauss-Seidel kernel only
3412  // runs on host. (See comments below.) Thus, we need to access
3413  // the host versions of these data.
3414  const_cast<DMV&> (B).sync_host ();
3415  X.sync_host ();
3416  X.modify_host ();
3417  const_cast<MMV&> (D).sync_host ();
3418 
3419  auto B_lcl = B.template getLocalView<host_mem_space> ();
3420  auto X_lcl = X.template getLocalView<host_mem_space> ();
3421  auto D_lcl = D.template getLocalView<host_mem_space> ();
3422 
3423  offset_type B_stride[8], X_stride[8], D_stride[8];
3424  B_lcl.stride (B_stride);
3425  X_lcl.stride (X_stride);
3426  D_lcl.stride (D_stride);
3427 
3428  local_matrix_type lclMatrix = this->getLocalMatrix ();
3429  k_local_graph_type lclGraph = lclMatrix.graph;
3430  typename local_matrix_type::row_map_type ptr = lclGraph.row_map;
3431  typename local_matrix_type::index_type ind = lclGraph.entries;
3432  typename local_matrix_type::values_type val = lclMatrix.values;
3433  const offset_type* const ptrRaw = ptr.data ();
3434  const LO* const indRaw = ind.data ();
3435  const impl_scalar_type* const valRaw = val.data ();
3436 
3437  const std::string dir ((direction == Forward) ? "F" : "B");
3438  // NOTE (mfh 28 Aug 2017) This assumes UVM. We can't get around
3439  // that on GPUs without using a GPU-based sparse triangular
3440  // solve to implement Gauss-Seidel. This exists in cuSPARSE,
3441  // but we would need to implement a wrapper with a fall-back
3442  // algorithm for unsupported Scalar and LO types.
3443  KokkosSparse::Impl::Sequential::gaussSeidel (static_cast<LO> (lclNumRows),
3444  static_cast<LO> (numVecs),
3445  ptrRaw, indRaw, valRaw,
3446  B_lcl.data (), B_stride[1],
3447  X_lcl.data (), X_stride[1],
3448  D_lcl.data (),
3449  static_cast<impl_scalar_type> (dampingFactor),
3450  dir.c_str ());
3451  const_cast<DMV&> (B).template sync<dev_mem_space> ();
3452  X.template sync<dev_mem_space> ();
3453  const_cast<MMV&> (D).template sync<dev_mem_space> ();
3454  }
3455 
3482  template <class DomainScalar, class RangeScalar>
3483  void
3487  const Teuchos::ArrayView<LocalOrdinal>& rowIndices,
3488  const RangeScalar& dampingFactor,
3489  const ESweepDirection direction) const
3490  {
3491  typedef LocalOrdinal LO;
3492  typedef GlobalOrdinal GO;
3495  typedef typename Node::device_type::memory_space dev_mem_space;
3496  typedef typename MMV::dual_view_type::t_host::device_type host_mem_space;
3497  typedef typename Graph::local_graph_type k_local_graph_type;
3498  typedef typename k_local_graph_type::size_type offset_type;
3499  const char prefix[] = "Tpetra::CrsMatrix::reorderedLocalGaussSeidel: ";
3500 
3501  TEUCHOS_TEST_FOR_EXCEPTION
3502  (! this->isFillComplete (), std::runtime_error,
3503  prefix << "The matrix is not fill complete.");
3504  const size_t lclNumRows = this->getNodeNumRows ();
3505  const size_t numVecs = B.getNumVectors ();
3506  TEUCHOS_TEST_FOR_EXCEPTION
3507  (X.getNumVectors () != numVecs, std::invalid_argument,
3508  prefix << "B.getNumVectors() = " << numVecs << " != "
3509  "X.getNumVectors() = " << X.getNumVectors () << ".");
3510  TEUCHOS_TEST_FOR_EXCEPTION
3511  (B.getLocalLength () != lclNumRows, std::invalid_argument,
3512  prefix << "B.getLocalLength() = " << B.getLocalLength ()
3513  << " != this->getNodeNumRows() = " << lclNumRows << ".");
3514  TEUCHOS_TEST_FOR_EXCEPTION
3515  (static_cast<size_t> (rowIndices.size ()) < lclNumRows,
3516  std::invalid_argument, prefix << "rowIndices.size() = "
3517  << rowIndices.size () << " < this->getNodeNumRows() = "
3518  << lclNumRows << ".");
3519 
3520  // mfh 28 Aug 2017: The current local Gauss-Seidel kernel only
3521  // runs on host. (See comments below.) Thus, we need to access
3522  // the host versions of these data.
3523  const_cast<DMV&> (B).sync_host ();
3524  X.sync_host ();
3525  X.modify_host ();
3526  const_cast<MMV&> (D).sync_host ();
3527 
3528  auto B_lcl = B.template getLocalView<host_mem_space> ();
3529  auto X_lcl = X.template getLocalView<host_mem_space> ();
3530  auto D_lcl = D.template getLocalView<host_mem_space> ();
3531 
3532  offset_type B_stride[8], X_stride[8], D_stride[8];
3533  B_lcl.stride (B_stride);
3534  X_lcl.stride (X_stride);
3535  D_lcl.stride (D_stride);
3536 
3537  local_matrix_type lclMatrix = this->getLocalMatrix ();
3538  typename Graph::local_graph_type lclGraph = lclMatrix.graph;
3539  typename local_matrix_type::index_type ind = lclGraph.entries;
3540  typename local_matrix_type::row_map_type ptr = lclGraph.row_map;
3541  typename local_matrix_type::values_type val = lclMatrix.values;
3542  const offset_type* const ptrRaw = ptr.data ();
3543  const LO* const indRaw = ind.data ();
3544  const impl_scalar_type* const valRaw = val.data ();
3545 
3546  const std::string dir = (direction == Forward) ? "F" : "B";
3547  // NOTE (mfh 28 Aug 2017) This assumes UVM. We can't get around
3548  // that on GPUs without using a GPU-based sparse triangular
3549  // solve to implement Gauss-Seidel, and also handling the
3550  // permutations correctly.
3551  KokkosSparse::Impl::Sequential::reorderedGaussSeidel (static_cast<LO> (lclNumRows),
3552  static_cast<LO> (numVecs),
3553  ptrRaw, indRaw, valRaw,
3554  B_lcl.data (),
3555  B_stride[1],
3556  X_lcl.data (),
3557  X_stride[1],
3558  D_lcl.data (),
3559  rowIndices.getRawPtr (),
3560  static_cast<LO> (lclNumRows),
3561  static_cast<impl_scalar_type> (dampingFactor),
3562  dir.c_str ());
3563  const_cast<DMV&> (B).template sync<dev_mem_space> ();
3564  X.template sync<dev_mem_space> ();
3565  const_cast<MMV&> (D).template sync<dev_mem_space> ();
3566  }
3567 
3570  template <class T>
3571  Teuchos::RCP<CrsMatrix<T, LocalOrdinal, GlobalOrdinal, Node> >
3572  convert () const;
3573 
3575 
3577 
3588  void
3591  Teuchos::ETransp mode = Teuchos::NO_TRANS,
3592  Scalar alpha = Teuchos::ScalarTraits<Scalar>::one(),
3593  Scalar beta = Teuchos::ScalarTraits<Scalar>::zero()) const override;
3594 
3597  bool hasTransposeApply () const override;
3598 
3605  Teuchos::RCP<const map_type> getDomainMap () const override;
3606 
3613  Teuchos::RCP<const map_type> getRangeMap () const override;
3614 
3616 
3618 
3683  void
3687  const Scalar& dampingFactor,
3688  const ESweepDirection direction,
3689  const int numSweeps) const;
3690 
3757  void
3761  const Teuchos::ArrayView<LocalOrdinal>& rowIndices,
3762  const Scalar& dampingFactor,
3763  const ESweepDirection direction,
3764  const int numSweeps) const;
3765 
3794  void
3798  const Scalar& dampingFactor,
3799  const ESweepDirection direction,
3800  const int numSweeps,
3801  const bool zeroInitialGuess) const;
3802 
3832  void
3836  const Teuchos::ArrayView<LocalOrdinal>& rowIndices,
3837  const Scalar& dampingFactor,
3838  const ESweepDirection direction,
3839  const int numSweeps,
3840  const bool zeroInitialGuess) const;
3841 
3852  virtual Teuchos::RCP<RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
3853  add (const Scalar& alpha,
3855  const Scalar& beta,
3856  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >& domainMap,
3857  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >& rangeMap,
3858  const Teuchos::RCP<Teuchos::ParameterList>& params) const override;
3859 
3861 
3863 
3865  std::string description () const override;
3866 
3869  void
3870  describe (Teuchos::FancyOStream& out,
3871  const Teuchos::EVerbosityLevel verbLevel =
3872  Teuchos::Describable::verbLevel_default) const override;
3873 
3875 
3877 
3882  typedef typename DistObject<Scalar, LocalOrdinal, GlobalOrdinal,
3884 
3885  virtual bool
3886  checkSizes (const SrcDistObject& source) override;
3887 
3888  private:
3889  void
3890  copyAndPermuteImpl (const RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& source,
3891  const size_t numSameIDs,
3892  const LocalOrdinal permuteToLIDs[],
3893  const LocalOrdinal permuteFromLIDs[],
3894  const size_t numPermutes);
3895  protected:
3896  virtual void
3897 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
3898  copyAndPermuteNew
3899 #else // TPETRA_ENABLE_DEPRECATED_CODE
3900  copyAndPermute
3901 #endif // TPETRA_ENABLE_DEPRECATED_CODE
3902  (const SrcDistObject& source,
3903  const size_t numSameIDs,
3904  const Kokkos::DualView<
3905  const local_ordinal_type*,
3906  buffer_device_type>& permuteToLIDs,
3907  const Kokkos::DualView<
3908  const local_ordinal_type*,
3909  buffer_device_type>& permuteFromLIDs) override;
3910 
3911  virtual void
3912 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
3913  packAndPrepareNew
3914 #else // TPETRA_ENABLE_DEPRECATED_CODE
3915  packAndPrepare
3916 #endif // TPETRA_ENABLE_DEPRECATED_CODE
3917  (const SrcDistObject& source,
3918  const Kokkos::DualView<
3919  const local_ordinal_type*,
3920  buffer_device_type>& exportLIDs,
3921  Kokkos::DualView<char*, buffer_device_type>& exports,
3922  Kokkos::DualView<size_t*, buffer_device_type> numPacketsPerLID,
3923  size_t& constantNumPackets,
3924  Distributor& distor) override;
3925 
3926  private:
3929  void
3930  unpackAndCombineImpl (const Kokkos::DualView<const local_ordinal_type*,
3931  buffer_device_type>& importLIDs,
3932  const Kokkos::DualView<const char*,
3933  buffer_device_type>& imports,
3934  const Kokkos::DualView<const size_t*,
3935  buffer_device_type>& numPacketsPerLID,
3936  const size_t constantNumPackets,
3937  Distributor& distor,
3938  const CombineMode combineMode,
3939  const bool atomic = useAtomicUpdatesByDefault);
3940 
3943  void
3944  unpackAndCombineImplNonStatic (const Kokkos::DualView<const local_ordinal_type*,
3945  buffer_device_type>& importLIDs,
3946  const Kokkos::DualView<const char*,
3947  buffer_device_type>& imports,
3948  const Kokkos::DualView<const size_t*,
3949  buffer_device_type>& numPacketsPerLID,
3950  const size_t constantNumPackets,
3951  Distributor& distor,
3952  const CombineMode combineMode);
3953 
3954  public:
3964  void
3965 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
3966  unpackAndCombineNew
3967 #else // TPETRA_ENABLE_DEPRECATED_CODE
3969 #endif // TPETRA_ENABLE_DEPRECATED_CODE
3970  (const Kokkos::DualView<const local_ordinal_type*, buffer_device_type>& importLIDs,
3971  Kokkos::DualView<char*, buffer_device_type> imports,
3972  Kokkos::DualView<size_t*, buffer_device_type> numPacketsPerLID,
3973  const size_t constantNumPackets,
3974  Distributor& distor,
3975  const CombineMode CM) override;
3976 
4084  void
4085  packNew (const Kokkos::DualView<const local_ordinal_type*, buffer_device_type>& exportLIDs,
4086  Kokkos::DualView<char*, buffer_device_type>& exports,
4087  const Kokkos::DualView<size_t*, buffer_device_type>& numPacketsPerLID,
4088  size_t& constantNumPackets,
4089  Distributor& dist) const;
4090 
4091  private:
4098  void
4099  packNonStaticNew (const Kokkos::DualView<const local_ordinal_type*, buffer_device_type>& exportLIDs,
4100  Kokkos::DualView<char*, buffer_device_type>& exports,
4101  const Kokkos::DualView<size_t*, buffer_device_type>& numPacketsPerLID,
4102  size_t& constantNumPackets,
4103  Distributor& distor) const;
4104 
4134  size_t
4135  packRow (char exports[],
4136  const size_t offset,
4137  const size_t numEnt,
4138  const GlobalOrdinal gidsIn[],
4139  const impl_scalar_type valsIn[],
4140  const size_t numBytesPerValue) const;
4141 
4165  bool
4166  packRowStatic (char* const numEntOut,
4167  char* const valOut,
4168  char* const indOut,
4169  const size_t numEnt,
4170  const LocalOrdinal lclRow) const;
4171 
4197  size_t
4198  unpackRow (GlobalOrdinal gidsOut[],
4199  impl_scalar_type valsOut[],
4200  const char imports[],
4201  const size_t offset,
4202  const size_t numBytes,
4203  const size_t numEnt,
4204  const size_t numBytesPerValue);
4205 
4214  void
4215  allocatePackSpaceNew (Kokkos::DualView<char*, buffer_device_type>& exports,
4216  size_t& totalNumEntries,
4217  const Kokkos::DualView<const local_ordinal_type*,
4218  buffer_device_type>& exportLIDs) const;
4220 
4221  public:
4223  typename local_matrix_type::values_type getLocalValuesView () const {
4224  return k_values1D_;
4225  }
4226 
4227  private:
4228  // Friend declaration for nonmember function.
4229  template<class CrsMatrixType>
4230  friend Teuchos::RCP<CrsMatrixType>
4231  Tpetra::importAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
4232  const Import<typename CrsMatrixType::local_ordinal_type,
4233  typename CrsMatrixType::global_ordinal_type,
4234  typename CrsMatrixType::node_type>& importer,
4235  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4236  typename CrsMatrixType::global_ordinal_type,
4237  typename CrsMatrixType::node_type> >& domainMap,
4238  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4239  typename CrsMatrixType::global_ordinal_type,
4240  typename CrsMatrixType::node_type> >& rangeMap,
4241  const Teuchos::RCP<Teuchos::ParameterList>& params);
4242 
4243  // Friend declaration for nonmember function.
4244  template<class CrsMatrixType>
4245  friend Teuchos::RCP<CrsMatrixType>
4246  Tpetra::importAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
4247  const Import<typename CrsMatrixType::local_ordinal_type,
4248  typename CrsMatrixType::global_ordinal_type,
4249  typename CrsMatrixType::node_type>& rowImporter,
4250  const Import<typename CrsMatrixType::local_ordinal_type,
4251  typename CrsMatrixType::global_ordinal_type,
4252  typename CrsMatrixType::node_type>& domainImporter,
4253  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4254  typename CrsMatrixType::global_ordinal_type,
4255  typename CrsMatrixType::node_type> >& domainMap,
4256  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4257  typename CrsMatrixType::global_ordinal_type,
4258  typename CrsMatrixType::node_type> >& rangeMap,
4259  const Teuchos::RCP<Teuchos::ParameterList>& params);
4260 
4261 
4262  // Friend declaration for nonmember function.
4263  template<class CrsMatrixType>
4264  friend Teuchos::RCP<CrsMatrixType>
4265  Tpetra::exportAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
4266  const Export<typename CrsMatrixType::local_ordinal_type,
4267  typename CrsMatrixType::global_ordinal_type,
4268  typename CrsMatrixType::node_type>& exporter,
4269  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4270  typename CrsMatrixType::global_ordinal_type,
4271  typename CrsMatrixType::node_type> >& domainMap,
4272  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4273  typename CrsMatrixType::global_ordinal_type,
4274  typename CrsMatrixType::node_type> >& rangeMap,
4275  const Teuchos::RCP<Teuchos::ParameterList>& params);
4276 
4277  // Friend declaration for nonmember function.
4278  template<class CrsMatrixType>
4279  friend Teuchos::RCP<CrsMatrixType>
4280  Tpetra::exportAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
4281  const Export<typename CrsMatrixType::local_ordinal_type,
4282  typename CrsMatrixType::global_ordinal_type,
4283  typename CrsMatrixType::node_type>& rowExporter,
4284  const Export<typename CrsMatrixType::local_ordinal_type,
4285  typename CrsMatrixType::global_ordinal_type,
4286  typename CrsMatrixType::node_type>& domainExporter,
4287  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4288  typename CrsMatrixType::global_ordinal_type,
4289  typename CrsMatrixType::node_type> >& domainMap,
4290  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4291  typename CrsMatrixType::global_ordinal_type,
4292  typename CrsMatrixType::node_type> >& rangeMap,
4293  const Teuchos::RCP<Teuchos::ParameterList>& params);
4294 
4295  public:
4311  void
4313  const import_type& importer,
4314  const Teuchos::RCP<const map_type>& domainMap,
4315  const Teuchos::RCP<const map_type>& rangeMap,
4316  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) const;
4317 
4333  void
4335  const import_type& rowImporter,
4336  const import_type& domainImporter,
4337  const Teuchos::RCP<const map_type>& domainMap,
4338  const Teuchos::RCP<const map_type>& rangeMap,
4339  const Teuchos::RCP<Teuchos::ParameterList>& params) const;
4340 
4341 
4357  void
4359  const export_type& exporter,
4360  const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
4361  const Teuchos::RCP<const map_type>& rangeMap = Teuchos::null,
4362  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) const;
4363 
4379  void
4381  const export_type& rowExporter,
4382  const export_type& domainExporter,
4383  const Teuchos::RCP<const map_type>& domainMap,
4384  const Teuchos::RCP<const map_type>& rangeMap,
4385  const Teuchos::RCP<Teuchos::ParameterList>& params) const;
4386 
4387 
4388  private:
4409  void
4410  transferAndFillComplete (Teuchos::RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >& destMatrix,
4411  const ::Tpetra::Details::Transfer<LocalOrdinal, GlobalOrdinal, Node>& rowTransfer,
4412  const Teuchos::RCP<const ::Tpetra::Details::Transfer<LocalOrdinal, GlobalOrdinal, Node> > & domainTransfer,
4413  const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
4414  const Teuchos::RCP<const map_type>& rangeMap = Teuchos::null,
4415  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) const;
4416 
4428  void
4429  insertGlobalValuesImpl (crs_graph_type& graph,
4430  RowInfo& rowInfo,
4431  const GlobalOrdinal gblColInds[],
4432  const impl_scalar_type vals[],
4433  const size_t numInputEnt);
4434 
4444  void
4445  insertGlobalValuesFiltered (const GlobalOrdinal globalRow,
4446  const Teuchos::ArrayView<const GlobalOrdinal>& indices,
4447  const Teuchos::ArrayView<const Scalar>& values);
4448 
4460  void
4461  combineGlobalValues (const GlobalOrdinal globalRowIndex,
4462  const Teuchos::ArrayView<const GlobalOrdinal>& columnIndices,
4463  const Teuchos::ArrayView<const Scalar>& values,
4464  const Tpetra::CombineMode combineMode);
4465 
4483  LocalOrdinal
4484  combineGlobalValuesRaw (const LocalOrdinal lclRow,
4485  const LocalOrdinal numEnt,
4486  const impl_scalar_type vals[],
4487  const GlobalOrdinal cols[],
4488  const Tpetra::CombineMode combineMode);
4489 
4501  template<class BinaryFunction>
4502  LocalOrdinal
4503  transformGlobalValues (const GlobalOrdinal globalRow,
4504  const Teuchos::ArrayView<const GlobalOrdinal>& indices,
4505  const Teuchos::ArrayView<const Scalar>& values,
4506  BinaryFunction f,
4507  const bool atomic = useAtomicUpdatesByDefault) const
4508  {
4509  typedef impl_scalar_type IST;
4510  typedef LocalOrdinal LO;
4511  typedef GlobalOrdinal GO;
4512 
4513  const LO numInputEnt = static_cast<LO> (indices.size ());
4514  if (static_cast<LO> (values.size ()) != numInputEnt) {
4515  return Teuchos::OrdinalTraits<LO>::invalid ();
4516  }
4517 
4518  const GO* const inputCols = indices.getRawPtr ();
4519  const IST* const inputVals =
4520  reinterpret_cast<const IST*> (values.getRawPtr ());
4521  return this->transformGlobalValues (globalRow, numInputEnt, inputVals,
4522  inputCols, f, atomic);
4523  }
4524 
4531  void
4532  insertNonownedGlobalValues (const GlobalOrdinal globalRow,
4533  const Teuchos::ArrayView<const GlobalOrdinal>& indices,
4534  const Teuchos::ArrayView<const Scalar>& values);
4535 
4578  void
4579  insertIndicesAndValues (crs_graph_type& graph,
4580  RowInfo& rowInfo,
4581  const typename crs_graph_type::SLocalGlobalViews& newInds,
4582  const Teuchos::ArrayView<impl_scalar_type>& oldRowVals,
4583  const Teuchos::ArrayView<const impl_scalar_type>& newRowVals,
4584  const ELocalGlobal lg,
4585  const ELocalGlobal I);
4586 
4588  typedef DistObject<char, LocalOrdinal, GlobalOrdinal, Node> dist_object_type;
4589 
4590  protected:
4591  // useful typedefs
4592  typedef Teuchos::OrdinalTraits<LocalOrdinal> OTL;
4593  typedef Kokkos::Details::ArithTraits<impl_scalar_type> STS;
4594  typedef Kokkos::Details::ArithTraits<mag_type> STM;
4595  typedef MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> MV;
4596  typedef Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> V;
4597  typedef crs_graph_type Graph;
4598 
4599  // Enums
4600  enum GraphAllocationStatus {
4601  GraphAlreadyAllocated,
4602  GraphNotYetAllocated
4603  };
4604 
4605  private:
4609  Teuchos::ArrayRCP<Teuchos::Array<impl_scalar_type> >
4610  allocateValues2D ();
4611 
4612  protected:
4629  void allocateValues (ELocalGlobal lg, GraphAllocationStatus gas);
4630 
4640  size_t
4642  const RowInfo& rowInfo);
4643 
4658  void
4659  sortAndMergeIndicesAndValues (const bool sorted,
4660  const bool merged);
4661 
4669  void clearGlobalConstants();
4670 
4678  public:
4680  void computeGlobalConstants();
4682  bool haveGlobalConstants() const;
4683  protected:
4696  mutable Teuchos::RCP<MV> importMV_;
4697 
4710  mutable Teuchos::RCP<MV> exportMV_;
4711 
4731  Teuchos::RCP<MV>
4732  getColumnMapMultiVector (const MV& X_domainMap,
4733  const bool force = false) const;
4734 
4756  Teuchos::RCP<MV>
4757  getRowMapMultiVector (const MV& Y_rangeMap,
4758  const bool force = false) const;
4759 
4761  void
4762  applyNonTranspose (const MV& X_in,
4763  MV& Y_in,
4764  Scalar alpha,
4765  Scalar beta) const;
4766 
4768  void
4769  applyTranspose (const MV& X_in,
4770  MV& Y_in,
4771  const Teuchos::ETransp mode,
4772  Scalar alpha,
4773  Scalar beta) const;
4774 
4775  // matrix data accessors
4776 
4795  LocalOrdinal
4796  getViewRawConst (const impl_scalar_type*& vals,
4797  LocalOrdinal& numEnt,
4798  const RowInfo& rowinfo) const;
4799 
4818  LocalOrdinal
4819  getViewRaw (impl_scalar_type*& vals,
4820  LocalOrdinal& numEnt,
4821  const RowInfo& rowinfo) const;
4822 
4830  Teuchos::ArrayView<const impl_scalar_type> getView (RowInfo rowinfo) const;
4831 
4843  Teuchos::ArrayView<impl_scalar_type>
4844  getViewNonConst (const RowInfo& rowinfo) const;
4845 
4846  private:
4854  Kokkos::View<const impl_scalar_type*, execution_space, Kokkos::MemoryUnmanaged>
4855  getRowView (const RowInfo& rowInfo) const;
4856 
4868  Kokkos::View<impl_scalar_type*, execution_space, Kokkos::MemoryUnmanaged>
4869  getRowViewNonConst (const RowInfo& rowInfo) const;
4870 
4871 
4872  protected:
4873 
4874  // Friend the tester for CrsMatrix::swap
4875  friend class Tpetra::crsMatrix_Swap_Tester<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
4876 
4881 
4882 
4883  protected:
4884 
4890  void fillLocalMatrix (const Teuchos::RCP<Teuchos::ParameterList>& params);
4891 
4897  void fillLocalGraphAndMatrix (const Teuchos::RCP<Teuchos::ParameterList>& params);
4898 
4900  void checkInternalState () const;
4901 
4913 
4914  Teuchos::RCP<const Graph> staticGraph_;
4915  Teuchos::RCP< Graph> myGraph_;
4917 
4920 
4933 
4934  typename local_matrix_type::values_type k_values1D_;
4935  Teuchos::ArrayRCP<Teuchos::Array<impl_scalar_type> > values2D_;
4937 
4948 
4951 
4979  std::map<GlobalOrdinal, std::pair<Teuchos::Array<GlobalOrdinal>,
4980  Teuchos::Array<Scalar> > > nonlocals_;
4981 
4988 
4989  public:
4990  // FIXME (mfh 24 Feb 2014) Is it _really_ necessary to make this a
4991  // public inner class of CrsMatrix? It looks like it doesn't
4992  // depend on any implementation details of CrsMatrix at all. It
4993  // should really be declared and defined outside of CrsMatrix.
4994  template<class ViewType, class OffsetViewType>
4995  struct pack_functor {
4996  typedef typename ViewType::execution_space execution_space;
4997  ViewType src_;
4998  ViewType dst_;
4999  OffsetViewType src_offset_;
5000  OffsetViewType dst_offset_;
5001  typedef typename OffsetViewType::non_const_value_type scalar_index_type;
5002 
5003  pack_functor (ViewType dst, ViewType src,
5004  OffsetViewType dst_offset, OffsetViewType src_offset) :
5005  src_ (src),
5006  dst_ (dst),
5007  src_offset_ (src_offset),
5008  dst_offset_ (dst_offset)
5009  {}
5010 
5011  KOKKOS_INLINE_FUNCTION
5012  void operator () (const LocalOrdinal row) const {
5013  scalar_index_type srcPos = src_offset_(row);
5014  const scalar_index_type dstEnd = dst_offset_(row+1);
5015  scalar_index_type dstPos = dst_offset_(row);
5016  for ( ; dstPos < dstEnd; ++dstPos, ++srcPos) {
5017  dst_(dstPos) = src_(srcPos);
5018  }
5019  }
5020  };
5021  }; // class CrsMatrix
5022 
5031  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
5032  Teuchos::RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
5034  size_t maxNumEntriesPerRow = 0,
5035  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null)
5036  {
5038  return Teuchos::rcp (new matrix_type (map, maxNumEntriesPerRow,
5039  DynamicProfile, params));
5040  }
5041 
5042  template<class CrsMatrixType>
5043  Teuchos::RCP<CrsMatrixType>
5044  importAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
5045  const Import<typename CrsMatrixType::local_ordinal_type,
5046  typename CrsMatrixType::global_ordinal_type,
5047  typename CrsMatrixType::node_type>& importer,
5048  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
5049  typename CrsMatrixType::global_ordinal_type,
5050  typename CrsMatrixType::node_type> >& domainMap,
5051  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
5052  typename CrsMatrixType::global_ordinal_type,
5053  typename CrsMatrixType::node_type> >& rangeMap,
5054  const Teuchos::RCP<Teuchos::ParameterList>& params)
5055  {
5056  Teuchos::RCP<CrsMatrixType> destMatrix;
5057  sourceMatrix->importAndFillComplete (destMatrix, importer, domainMap, rangeMap, params);
5058  return destMatrix;
5059  }
5060 
5061  template<class CrsMatrixType>
5062  Teuchos::RCP<CrsMatrixType>
5063  importAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
5064  const Import<typename CrsMatrixType::local_ordinal_type,
5065  typename CrsMatrixType::global_ordinal_type,
5066  typename CrsMatrixType::node_type>& rowImporter,
5067  const Import<typename CrsMatrixType::local_ordinal_type,
5068  typename CrsMatrixType::global_ordinal_type,
5069  typename CrsMatrixType::node_type>& domainImporter,
5070  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
5071  typename CrsMatrixType::global_ordinal_type,
5072  typename CrsMatrixType::node_type> >& domainMap,
5073  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
5074  typename CrsMatrixType::global_ordinal_type,
5075  typename CrsMatrixType::node_type> >& rangeMap,
5076  const Teuchos::RCP<Teuchos::ParameterList>& params)
5077  {
5078  Teuchos::RCP<CrsMatrixType> destMatrix;
5079  sourceMatrix->importAndFillComplete (destMatrix, rowImporter, domainImporter, domainMap, rangeMap, params);
5080  return destMatrix;
5081  }
5082 
5083  template<class CrsMatrixType>
5084  Teuchos::RCP<CrsMatrixType>
5085  exportAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
5086  const Export<typename CrsMatrixType::local_ordinal_type,
5087  typename CrsMatrixType::global_ordinal_type,
5088  typename CrsMatrixType::node_type>& exporter,
5089  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
5090  typename CrsMatrixType::global_ordinal_type,
5091  typename CrsMatrixType::node_type> >& domainMap,
5092  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
5093  typename CrsMatrixType::global_ordinal_type,
5094  typename CrsMatrixType::node_type> >& rangeMap,
5095  const Teuchos::RCP<Teuchos::ParameterList>& params)
5096  {
5097  Teuchos::RCP<CrsMatrixType> destMatrix;
5098  sourceMatrix->exportAndFillComplete (destMatrix, exporter, domainMap, rangeMap, params);
5099  return destMatrix;
5100  }
5101 
5102  template<class CrsMatrixType>
5103  Teuchos::RCP<CrsMatrixType>
5104  exportAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
5105  const Export<typename CrsMatrixType::local_ordinal_type,
5106  typename CrsMatrixType::global_ordinal_type,
5107  typename CrsMatrixType::node_type>& rowExporter,
5108  const Export<typename CrsMatrixType::local_ordinal_type,
5109  typename CrsMatrixType::global_ordinal_type,
5110  typename CrsMatrixType::node_type>& domainExporter,
5111  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
5112  typename CrsMatrixType::global_ordinal_type,
5113  typename CrsMatrixType::node_type> >& domainMap,
5114  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
5115  typename CrsMatrixType::global_ordinal_type,
5116  typename CrsMatrixType::node_type> >& rangeMap,
5117  const Teuchos::RCP<Teuchos::ParameterList>& params)
5118  {
5119  Teuchos::RCP<CrsMatrixType> destMatrix;
5120  sourceMatrix->exportAndFillComplete (destMatrix, rowExporter, domainExporter, domainMap, rangeMap, params);
5121  return destMatrix;
5122  }
5123 } // namespace Tpetra
5124 
5132 #endif // TPETRA_CRSMATRIX_DECL_HPP
Communication plan for data redistribution from a uniquely-owned to a (possibly) multiply-owned distr...
bool hasColMap() const override
Whether the matrix has a well-defined column Map.
global_size_t getGlobalNumCols() const override
The number of global columns in the matrix.
void checkInternalState() const
Check that this object&#39;s state is sane; throw if it&#39;s not.
Sparse matrix that presents a row-oriented interface that lets users read or modify entries...
size_t getNodeNumCols() const override
The number of columns connected to the locally owned rows of this matrix.
Teuchos::RCP< CrsMatrixType > exportAndFillCompleteCrsMatrix(const Teuchos::RCP< const CrsMatrixType > &sourceMatrix, const Export< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &exporter, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Nonmember CrsMatrix constructor that fuses Export and fillComplete().
CrsGraph< LocalOrdinal, GlobalOrdinal, Node > crs_graph_type
The CrsGraph specialization suitable for this CrsMatrix specialization.
void replaceColMap(const Teuchos::RCP< const map_type > &newColMap)
Replace the matrix&#39;s column Map with the given Map.
virtual bool supportsRowViews() const override
Return true if getLocalRowView() and getGlobalRowView() are valid for this object.
LocalOrdinal getViewRaw(impl_scalar_type *&vals, LocalOrdinal &numEnt, const RowInfo &rowinfo) const
Nonconst pointer to all entries (including extra space) in the given row.
void replaceDomainMapAndImporter(const Teuchos::RCP< const map_type > &newDomainMap, Teuchos::RCP< const import_type > &newImporter)
Replace the current domain Map and Import with the given objects.
std::enable_if<!std::is_same< OutputScalarType, impl_scalar_type >::value &&std::is_convertible< impl_scalar_type, OutputScalarType >::value, LocalOrdinal >::type getLocalRowView(const LocalOrdinal lclRow, LocalOrdinal &numEnt, const OutputScalarType *&val, const LocalOrdinal *&ind) const
Get a constant, nonpersisting view of a row of this matrix, using local row and column indices...
LocalOrdinal replaceLocalValues(const LocalOrdinal localRow, const typename UnmanagedView< LocalIndicesViewType >::type &inputInds, const typename UnmanagedView< ImplScalarViewType >::type &inputVals) const
Replace one or more entries&#39; values, using local row and column indices.
Teuchos::RCP< const map_type > getRangeMap() const override
The range Map of this matrix.
global_size_t getGlobalNumRows() const override
Number of global elements in the row map of this matrix.
std::map< GlobalOrdinal, std::pair< Teuchos::Array< GlobalOrdinal >, Teuchos::Array< Scalar > > > nonlocals_
Nonlocal data added using insertGlobalValues().
void sortAndMergeIndicesAndValues(const bool sorted, const bool merged)
Sort and merge duplicate local column indices in all rows on the calling process, along with their co...
typename device_type::execution_space execution_space
The Kokkos execution space.
void getLocalRowCopy(LocalOrdinal localRow, const Teuchos::ArrayView< LocalOrdinal > &colInds, const Teuchos::ArrayView< Scalar > &vals, size_t &numEntries) const override
Fill given arrays with a deep copy of the locally owned entries of the matrix in a given row...
typename Kokkos::ArithTraits< impl_scalar_type >::mag_type mag_type
Type of a norm result.
size_t getNumVectors() const
Number of columns in the multivector.
size_t getLocalLength() const
Local number of rows on the calling process.
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to alpha.
bool isConstantStride() const
Whether this multivector has constant stride between columns.
One or more distributed dense vectors.
size_t getNodeNumEntries() const override
The local number of entries in this matrix.
void scale(const Scalar &alpha)
Scale the matrix&#39;s values: this := alpha*this.
size_t getNodeNumRows() const override
The number of matrix rows owned by the calling process.
void fillLocalGraphAndMatrix(const Teuchos::RCP< Teuchos::ParameterList > &params)
Fill data into the local graph and matrix.
void resumeFill(const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Resume operations that may change the values or structure of the matrix.
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x) override
Scale the matrix on the left with the given Vector.
local_matrix_type lclMatrix_
The local sparse matrix.
GlobalOrdinal global_ordinal_type
The type of each global index in the matrix.
Declaration and generic definition of traits class that tells Tpetra::CrsMatrix how to pack and unpac...
virtual ~CrsMatrix()=default
Destructor (virtual for memory safety of derived classes).
void unpackAndCombine(const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< char *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, Distributor &distor, const CombineMode CM) override
Unpack the imported column indices and values, and combine into matrix.
void reindexColumns(crs_graph_type *const graph, const Teuchos::RCP< const map_type > &newColMap, const Teuchos::RCP< const import_type > &newImport=Teuchos::null, const bool sortEachRow=true)
Reindex the column indices in place, and replace the column Map. Optionally, replace the Import objec...
Import< LocalOrdinal, GlobalOrdinal, Node > import_type
The Import specialization suitable for this CrsMatrix specialization.
void getGlobalRowCopy(GlobalOrdinal GlobalRow, const Teuchos::ArrayView< GlobalOrdinal > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const override
Fill given arrays with a deep copy of the locally owned entries of the matrix in a given row...
Teuchos::RCP< const crs_graph_type > getCrsGraph() const
This matrix&#39;s graph, as a CrsGraph.
size_t getGlobalMaxNumRowEntries() const override
Maximum number of entries in any row of the matrix, over all processes in the matrix&#39;s communicator...
void applyNonTranspose(const MV &X_in, MV &Y_in, Scalar alpha, Scalar beta) const
Special case of apply() for mode == Teuchos::NO_TRANS.
Teuchos::RCP< CrsMatrix< T, LocalOrdinal, GlobalOrdinal, Node > > convert() const
Return another CrsMatrix with the same entries, but converted to a different Scalar type T...
Scalar scalar_type
The type of each entry in the matrix.
Allocation information for a locally owned row in a CrsGraph or CrsMatrix.
LocalOrdinal sumIntoLocalValues(const LocalOrdinal localRow, const typename UnmanagedView< LocalIndicesViewType >::type &inputInds, const typename UnmanagedView< ImplScalarViewType >::type &inputVals, const bool atomic=useAtomicUpdatesByDefault) const
Sum into one or more sparse matrix entries, using local row and column indices.
void swap(CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &matrix)
Swaps the data from *this with the data and maps from crsMatrix.
void fillLocalMatrix(const Teuchos::RCP< Teuchos::ParameterList > &params)
Fill data into the local matrix.
void globalAssemble()
Communicate nonlocal contributions to other processes.
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const override
Get a copy of the diagonal entries of the matrix.
void gaussSeidel(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &D, const Scalar &dampingFactor, const ESweepDirection direction, const int numSweeps) const
&quot;Hybrid&quot; Jacobi + (Gauss-Seidel or SOR) on .
Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > createCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t maxNumEntriesPerRow=0, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Non-member function to create an empty CrsMatrix given a row map and a non-zero profile.
bool isFillActive() const
Whether the matrix is not fill complete.
Teuchos::RCP< MV > importMV_
Column Map MultiVector used in apply() and gaussSeidel().
typename Kokkos::ArithTraits< Scalar >::val_type impl_scalar_type
The type used internally in place of Scalar.
bool isStaticGraph() const
Indicates that the graph is static, so that new entries cannot be added to this matrix.
size_t global_size_t
Global size_t object.
bool hasTransposeApply() const override
Whether apply() allows applying the transpose or conjugate transpose.
void clearGlobalConstants()
Clear matrix properties that require collectives.
void gaussSeidelCopy(MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &D, const Scalar &dampingFactor, const ESweepDirection direction, const int numSweeps, const bool zeroInitialGuess) const
Version of gaussSeidel(), with fewer requirements on X.
LocalOrdinal getLocalRowViewRaw(const LocalOrdinal lclRow, LocalOrdinal &numEnt, const LocalOrdinal *&lclColInds, const Scalar *&vals) const override
Get a constant, nonpersisting, locally indexed view of the given row of the matrix, using &quot;raw&quot; pointers instead of Teuchos::ArrayView.
LocalOrdinal transformLocalValues(const LocalOrdinal lclRow, const typename UnmanagedView< LocalIndicesViewType >::type &inputInds, const typename UnmanagedView< ImplScalarViewType >::type &inputVals, BinaryFunction f, const bool atomic=useAtomicUpdatesByDefault) const
Transform CrsMatrix entries in place, using local indices to select the entries in the row to transfo...
void allocateValues(ELocalGlobal lg, GraphAllocationStatus gas)
Allocate values (and optionally indices) using the Node.
void exportAndFillComplete(Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &destMatrix, const export_type &exporter, const Teuchos::RCP< const map_type > &domainMap=Teuchos::null, const Teuchos::RCP< const map_type > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) const
Export from this to the given destination matrix, and make the result fill complete.
bool isFillComplete() const override
Whether the matrix is fill complete.
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const override
The communicator over which the matrix is distributed.
void sync_host()
Synchronize to Host.
void importAndFillComplete(Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &destMatrix, const import_type &importer, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) const
Import from this to the given destination matrix, and make the result fill complete.
void packNew(const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< char *, buffer_device_type > &exports, const Kokkos::DualView< size_t *, buffer_device_type > &numPacketsPerLID, size_t &constantNumPackets, Distributor &dist) const
Pack this object&#39;s data for an Import or Export.
void reorderedGaussSeidelCopy(MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &D, const Teuchos::ArrayView< LocalOrdinal > &rowIndices, const Scalar &dampingFactor, const ESweepDirection direction, const int numSweeps, const bool zeroInitialGuess) const
Version of reorderedGaussSeidel(), with fewer requirements on X.
ESweepDirection
Sweep direction for Gauss-Seidel or Successive Over-Relaxation (SOR).
GlobalOrdinal getIndexBase() const override
The index base for global indices for this matrix.
Communication plan for data redistribution from a (possibly) multiply-owned to a uniquely-owned distr...
typename crs_graph_type::local_graph_type local_graph_type
The part of the sparse matrix&#39;s graph on each MPI process.
Kokkos::StaticCrsGraph< local_ordinal_type, Kokkos::LayoutLeft, execution_space > local_graph_type
The type of the part of the sparse graph on each MPI process.
Node node_type
This class&#39; Kokkos Node type.
Teuchos::RCP< MV > getColumnMapMultiVector(const MV &X_domainMap, const bool force=false) const
Create a (or fetch a cached) column Map MultiVector.
Sets up and executes a communication plan for a Tpetra DistObject.
mag_type frobNorm_
Cached Frobenius norm of the (global) matrix.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override
Print this object with the given verbosity level to the given output stream.
Teuchos::RCP< const RowGraph< LocalOrdinal, GlobalOrdinal, Node > > getGraph() const override
This matrix&#39;s graph, as a RowGraph.
CombineMode
Rule for combining data in an Import or Export.
DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node >::buffer_device_type buffer_device_type
Kokkos::Device specialization for communication buffers.
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
CrsMatrix & operator=(const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &)=delete
Copy assignment (forbidden).
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...
bool isGloballyIndexed() const override
Whether the matrix is globally indexed on the calling process.
Forward declaration of Tpetra::CrsMatrix.
LocalOrdinal sumIntoGlobalValues(const GlobalOrdinal globalRow, const Teuchos::ArrayView< const GlobalOrdinal > &cols, const Teuchos::ArrayView< const Scalar > &vals, const bool atomic=useAtomicUpdatesByDefault)
Sum into one or more sparse matrix entries, using global indices.
void localGaussSeidel(const MultiVector< DomainScalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< RangeScalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &D, const RangeScalar &dampingFactor, const ESweepDirection direction) const
Gauss-Seidel or SOR on .
bool fillComplete_
Whether the matrix is fill complete.
virtual Teuchos::RCP< RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > add(const Scalar &alpha, const RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params) const override
Implementation of RowMatrix::add: return alpha*A + beta*this.
Abstract base class for objects that can be the source of an Import or Export operation.
local_matrix_type getLocalMatrix() const
The local sparse matrix.
typename Node::device_type device_type
The Kokkos device type.
size_t getNumEntriesInLocalRow(local_ordinal_type localRow) const override
Number of entries in the sparse matrix in the given local row, on the calling (MPI) process...
Teuchos::RCP< MV > getRowMapMultiVector(const MV &Y_rangeMap, const bool force=false) const
Create a (or fetch a cached) row Map MultiVector.
std::string description() const override
A one-line description of this object.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
void apply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const override
Compute a sparse matrix-MultiVector multiply.
void computeGlobalConstants()
Compute matrix properties that require collectives.
void fillComplete(const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Tell the matrix that you are done changing its structure or values, and that you are ready to do comp...
void insertGlobalValues(const GlobalOrdinal globalRow, const Teuchos::ArrayView< const GlobalOrdinal > &cols, const Teuchos::ArrayView< const Scalar > &vals)
Insert one or more entries into the matrix, using global column indices.
void modify_host()
Mark data as modified on the host side.
local_matrix_type::values_type getLocalValuesView() const
Get the Kokkos local values.
ProfileType getProfileType() const
Returns true if the matrix was allocated with static data structures.
void getGlobalRowView(GlobalOrdinal GlobalRow, Teuchos::ArrayView< const GlobalOrdinal > &indices, Teuchos::ArrayView< const Scalar > &values) const override
Get a constant, nonpersisting view of a row of this matrix, using global row and column indices...
LocalOrdinal transformGlobalValues(const GlobalOrdinal gblRow, const Kokkos::View< const GlobalOrdinal *, InputMemorySpace, Kokkos::MemoryUnmanaged > &inputInds, const Kokkos::View< const impl_scalar_type *, InputMemorySpace, Kokkos::MemoryUnmanaged > &inputVals, BinaryFunction f, const bool atomic=useAtomicUpdatesByDefault) const
Transform CrsMatrix entries in place, using global indices to select the entries in the row to transf...
virtual void removeEmptyProcessesInPlace(const Teuchos::RCP< const map_type > &newMap) override
Remove processes owning zero rows from the Maps and their communicator.
virtual bool checkSizes(const SrcDistObject &source) override
Compare the source and target (this) objects for compatibility.
void insertLocalValues(const LocalOrdinal localRow, const Teuchos::ArrayView< const LocalOrdinal > &cols, const Teuchos::ArrayView< const Scalar > &vals)
Insert one or more entries into the matrix, using local column indices.
A distributed graph accessed by rows (adjacency lists) and stored sparsely.
A parallel distribution of indices over processes.
Teuchos::RCP< const map_type > getDomainMap() const override
The domain Map of this matrix.
Teuchos::ArrayView< const impl_scalar_type > getView(RowInfo rowinfo) const
Constant view of all entries (including extra space) in the given row.
KokkosSparse::CrsMatrix< impl_scalar_type, local_ordinal_type, execution_space, void, typename local_graph_type::size_type > local_matrix_type
The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix on each MPI pro...
void setAllValues(const typename local_matrix_type::row_map_type &ptr, const typename local_graph_type::entries_type::non_const_type &ind, const typename local_matrix_type::values_type &val)
Set the local matrix using three (compressed sparse row) arrays.
void localMultiply(const MultiVector< DomainScalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< RangeScalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode, RangeScalar alpha, RangeScalar beta) const
Compute a sparse matrix-MultiVector product local to each process.
A read-only, row-oriented interface to a sparse matrix.
void getLocalRowView(LocalOrdinal LocalRow, Teuchos::ArrayView< const LocalOrdinal > &indices, Teuchos::ArrayView< const Scalar > &values) const override
Get a constant, nonpersisting view of a row of this matrix, using local row and column indices...
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x) override
Scale the matrix on the right with the given Vector.
size_t getNodeMaxNumRowEntries() const override
Maximum number of entries in any row of the matrix, on this process.
A distributed dense vector.
Teuchos::RCP< CrsMatrixType > importAndFillCompleteCrsMatrix(const Teuchos::RCP< const CrsMatrixType > &sourceMatrix, const Import< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &importer, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Nonmember CrsMatrix constructor that fuses Import and fillComplete().
void applyTranspose(const MV &X_in, MV &Y_in, const Teuchos::ETransp mode, Scalar alpha, Scalar beta) const
Special case of apply() for mode != Teuchos::NO_TRANS.
void expertStaticFillComplete(const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< const import_type > &importer=Teuchos::null, const Teuchos::RCP< const export_type > &exporter=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Perform a fillComplete on a matrix that already has data.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const override
Number of entries in the sparse matrix in the given global row, on the calling (MPI) process...
size_t mergeRowIndicesAndValues(crs_graph_type &graph, const RowInfo &rowInfo)
Merge duplicate row indices in the given row, along with their corresponding values.
LocalOrdinal replaceGlobalValues(const GlobalOrdinal globalRow, const typename UnmanagedView< GlobalIndicesViewType >::type &inputInds, const typename UnmanagedView< ImplScalarViewType >::type &inputVals) const
Replace one or more entries&#39; values, using global indices.
bool isLocallyIndexed() const override
Whether the matrix is locally indexed on the calling process.
Teuchos::RCP< const map_type > getColMap() const override
The Map that describes the column distribution in this matrix.
void reorderedGaussSeidel(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &D, const Teuchos::ArrayView< LocalOrdinal > &rowIndices, const Scalar &dampingFactor, const ESweepDirection direction, const int numSweeps) const
Reordered &quot;Hybrid&quot; Jacobi + (Gauss-Seidel or SOR) on .
Teuchos::ArrayView< impl_scalar_type > getViewNonConst(const RowInfo &rowinfo) const
Nonconst view of all entries (including extra space) in the given row.
Base class for distributed Tpetra objects that support data redistribution.
LocalOrdinal local_ordinal_type
The type of each local index in the matrix.
EStorageStatus
Status of the graph&#39;s or matrix&#39;s storage, when not in a fill-complete state.
mag_type getFrobeniusNorm() const override
Compute and return the Frobenius norm of the matrix.
::Tpetra::Details::EStorageStatus storageStatus_
Status of the matrix&#39;s storage, when not in a fill-complete state.
void localApply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, const Teuchos::ETransp mode=Teuchos::NO_TRANS, const Scalar &alpha=Teuchos::ScalarTraits< Scalar >::one(), const Scalar &beta=Teuchos::ScalarTraits< Scalar >::zero()) const
Compute the local part of a sparse matrix-(Multi)Vector multiply.
void reorderedLocalGaussSeidel(const MultiVector< DomainScalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< RangeScalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &D, const Teuchos::ArrayView< LocalOrdinal > &rowIndices, const RangeScalar &dampingFactor, const ESweepDirection direction) const
Reordered Gauss-Seidel or SOR on .
LocalOrdinal getViewRawConst(const impl_scalar_type *&vals, LocalOrdinal &numEnt, const RowInfo &rowinfo) const
Const pointer to all entries (including extra space) in the given row.
bool isStorageOptimized() const
Returns true if storage has been optimized.
Export< LocalOrdinal, GlobalOrdinal, Node > export_type
The Export specialization suitable for this CrsMatrix specialization.
Teuchos::RCP< MV > exportMV_
Row Map MultiVector used in apply().
Teuchos::RCP< const map_type > getRowMap() const override
The Map that describes the row distribution in this matrix.
global_size_t getGlobalNumEntries() const override
The global number of entries in this matrix.