Tpetra parallel linear algebra  Version of the Day
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 
124  template<class CrsMatrixType>
125  Teuchos::RCP<CrsMatrixType>
126  importAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
127  const Import<typename CrsMatrixType::local_ordinal_type,
128  typename CrsMatrixType::global_ordinal_type,
129  typename CrsMatrixType::node_type>& importer,
130  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
131  typename CrsMatrixType::global_ordinal_type,
132  typename CrsMatrixType::node_type> >& domainMap = Teuchos::null,
133  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
134  typename CrsMatrixType::global_ordinal_type,
135  typename CrsMatrixType::node_type> >& rangeMap = Teuchos::null,
136  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
137 
191  template<class CrsMatrixType>
192  Teuchos::RCP<CrsMatrixType>
193  importAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
194  const Import<typename CrsMatrixType::local_ordinal_type,
195  typename CrsMatrixType::global_ordinal_type,
196  typename CrsMatrixType::node_type>& rowImporter,
197  const Import<typename CrsMatrixType::local_ordinal_type,
198  typename CrsMatrixType::global_ordinal_type,
199  typename CrsMatrixType::node_type>& domainImporter,
200  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
201  typename CrsMatrixType::global_ordinal_type,
202  typename CrsMatrixType::node_type> >& domainMap,
203  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
204  typename CrsMatrixType::global_ordinal_type,
205  typename CrsMatrixType::node_type> >& rangeMap,
206  const Teuchos::RCP<Teuchos::ParameterList>& params);
207 
241  template<class CrsMatrixType>
242  Teuchos::RCP<CrsMatrixType>
243  exportAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
244  const Export<typename CrsMatrixType::local_ordinal_type,
245  typename CrsMatrixType::global_ordinal_type,
246  typename CrsMatrixType::node_type>& exporter,
247  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
248  typename CrsMatrixType::global_ordinal_type,
249  typename CrsMatrixType::node_type> >& domainMap = Teuchos::null,
250  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
251  typename CrsMatrixType::global_ordinal_type,
252  typename CrsMatrixType::node_type> >& rangeMap = Teuchos::null,
253  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
254 
288  template<class CrsMatrixType>
289  Teuchos::RCP<CrsMatrixType>
290  exportAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
291  const Export<typename CrsMatrixType::local_ordinal_type,
292  typename CrsMatrixType::global_ordinal_type,
293  typename CrsMatrixType::node_type>& rowExporter,
294  const Export<typename CrsMatrixType::local_ordinal_type,
295  typename CrsMatrixType::global_ordinal_type,
296  typename CrsMatrixType::node_type>& domainExporter,
297  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
298  typename CrsMatrixType::global_ordinal_type,
299  typename CrsMatrixType::node_type> >& domainMap,
300  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
301  typename CrsMatrixType::global_ordinal_type,
302  typename CrsMatrixType::node_type> >& rangeMap,
303  const Teuchos::RCP<Teuchos::ParameterList>& params);
304 
305 namespace Classes {
306 
420  template <class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type,
422  class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type,
424  class CrsMatrix :
425  public RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>,
426  public DistObject<char, LocalOrdinal, GlobalOrdinal, Node>,
428  {
429  public:
431 
432 
435  typedef Scalar scalar_type;
445  typedef typename Kokkos::Details::ArithTraits<Scalar>::val_type impl_scalar_type;
447  typedef LocalOrdinal local_ordinal_type;
449  typedef GlobalOrdinal global_ordinal_type;
451  typedef Node node_type;
452 
454  typedef typename Node::device_type device_type;
456  typedef typename device_type::execution_space execution_space;
457 
463  typedef typename Kokkos::Details::ArithTraits<impl_scalar_type>::mag_type mag_type;
464 
467 
470 
473 
476 
479 
482  typedef KokkosSparse::CrsMatrix<impl_scalar_type, LocalOrdinal, execution_space, void,
483  typename local_graph_type::size_type> local_matrix_type;
484 
486  typedef typename local_matrix_type::row_map_type t_RowPtrs TPETRA_DEPRECATED;
488  typedef typename local_matrix_type::row_map_type::non_const_type t_RowPtrsNC TPETRA_DEPRECATED;
490  typedef typename local_graph_type::entries_type::non_const_type t_LocalOrdinal_1D TPETRA_DEPRECATED;
492  typedef typename local_matrix_type::values_type t_ValuesType TPETRA_DEPRECATED;
493 
495  typedef local_matrix_type k_local_matrix_type TPETRA_DEPRECATED;
496 
498 
500 
518  CrsMatrix (const Teuchos::RCP<const map_type>& rowMap,
519  size_t maxNumEntriesPerRow,
520  ProfileType pftype = DynamicProfile,
521  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
522 
540  CrsMatrix (const Teuchos::RCP<const map_type>& rowMap,
541  const Teuchos::ArrayRCP<const size_t>& NumEntriesPerRowToAlloc,
542  ProfileType pftype = DynamicProfile,
543  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
544 
567  CrsMatrix (const Teuchos::RCP<const map_type>& rowMap,
568  const Teuchos::RCP<const map_type>& colMap,
569  size_t maxNumEntriesPerRow,
570  ProfileType pftype = DynamicProfile,
571  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
572 
595  CrsMatrix (const Teuchos::RCP<const map_type>& rowMap,
596  const Teuchos::RCP<const map_type>& colMap,
597  const Teuchos::ArrayRCP<const size_t>& NumEntriesPerRowToAlloc,
598  ProfileType pftype = DynamicProfile,
599  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
600 
625  explicit CrsMatrix (const Teuchos::RCP<const crs_graph_type>& graph,
626  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
627 
651  CrsMatrix (const Teuchos::RCP<const map_type>& rowMap,
652  const Teuchos::RCP<const map_type>& colMap,
653  const typename local_matrix_type::row_map_type& rowPointers,
654  const typename local_graph_type::entries_type::non_const_type& columnIndices,
655  const typename local_matrix_type::values_type& values,
656  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
657 
681  CrsMatrix (const Teuchos::RCP<const map_type>& rowMap,
682  const Teuchos::RCP<const map_type>& colMap,
683  const Teuchos::ArrayRCP<size_t>& rowPointers,
684  const Teuchos::ArrayRCP<LocalOrdinal>& columnIndices,
685  const Teuchos::ArrayRCP<Scalar>& values,
686  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
687 
708  CrsMatrix (const Teuchos::RCP<const map_type>& rowMap,
709  const Teuchos::RCP<const map_type>& colMap,
710  const local_matrix_type& lclMatrix,
711  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
712 
739  CrsMatrix (const local_matrix_type& lclMatrix,
740  const Teuchos::RCP<const map_type>& rowMap,
741  const Teuchos::RCP<const map_type>& colMap,
742  const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
743  const Teuchos::RCP<const map_type>& rangeMap = Teuchos::null,
744  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
745 
746  // This friend declaration makes the clone() method work.
747  template <class S2, class LO2, class GO2, class N2>
748  friend class CrsMatrix;
749 
774  template <class Node2>
775  Teuchos::RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node2> >
776  clone (const Teuchos::RCP<Node2>& node2,
777  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) const
778  {
779  using Teuchos::Array;
780  using Teuchos::ArrayRCP;
781  using Teuchos::ArrayView;
782  using Teuchos::null;
783  using Teuchos::ParameterList;
784  using Teuchos::RCP;
785  using Teuchos::rcp;
786  using Teuchos::sublist;
789  const char tfecfFuncName[] = "clone";
790 
791  // Get parameter values. Set them initially to their default values.
792  bool fillCompleteClone = true;
793  bool useLocalIndices = this->hasColMap ();
794  ProfileType pftype = StaticProfile;
795  if (! params.is_null ()) {
796  fillCompleteClone = params->get ("fillComplete clone", fillCompleteClone);
797  useLocalIndices = params->get ("Locally indexed clone", useLocalIndices);
798 
799  bool staticProfileClone = true;
800  staticProfileClone = params->get ("Static profile clone", staticProfileClone);
801  pftype = staticProfileClone ? StaticProfile : DynamicProfile;
802  }
803 
804  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(
805  ! this->hasColMap () && useLocalIndices, std::runtime_error,
806  ": You requested that the returned clone have local indices, but the "
807  "the source matrix does not have a column Map yet.");
808 
809  RCP<const Map2> clonedRowMap = this->getRowMap ()->template clone<Node2> (node2);
810 
811  // Get an upper bound on the number of entries per row.
812  RCP<CrsMatrix2> clonedMatrix;
813  ArrayRCP<const size_t> numEntriesPerRow;
814  size_t numEntriesForAll = 0;
815  bool boundSameForAllLocalRows = false;
816  staticGraph_->getNumEntriesPerLocalRowUpperBound (numEntriesPerRow,
817  numEntriesForAll,
818  boundSameForAllLocalRows);
819  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(
820  numEntriesForAll != 0 &&
821  static_cast<size_t> (numEntriesPerRow.size ()) != 0,
822  std::logic_error, ": getNumEntriesPerLocalRowUpperBound returned a "
823  "nonzero numEntriesForAll = " << numEntriesForAll << " , as well as a "
824  "numEntriesPerRow array of nonzero length " << numEntriesPerRow.size ()
825  << ". This should never happen. Please report this bug to the Tpetra "
826  "developers.");
827  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(
828  numEntriesForAll != 0 && ! boundSameForAllLocalRows,
829  std::logic_error, ": getNumEntriesPerLocalRowUpperBound returned a "
830  "nonzero numEntriesForAll = " << numEntriesForAll << " , but claims "
831  "(via its third output value) that the upper bound is not the same for "
832  "all rows. This should never happen. Please report this bug to the "
833  "Tpetra developers.");
834  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(
835  numEntriesPerRow.size () != 0 && boundSameForAllLocalRows,
836  std::logic_error, ": getNumEntriesPerLocalRowUpperBound returned a "
837  "numEntriesPerRow array of nonzero length " << numEntriesPerRow.size ()
838  << ", but claims (via its third output value) that the upper bound is "
839  "not the same for all rows. This should never happen. Please report "
840  "this bug to the Tpetra developers.");
841 
842  RCP<ParameterList> matParams =
843  params.is_null () ? null : sublist (params,"CrsMatrix");
844  if (useLocalIndices) {
845  RCP<const Map2> clonedColMap =
846  this->getColMap ()->template clone<Node2> (node2);
847  if (numEntriesPerRow.is_null ()) {
848  clonedMatrix = rcp (new CrsMatrix2 (clonedRowMap, clonedColMap,
849  numEntriesForAll, pftype,
850  matParams));
851  }
852  else {
853  clonedMatrix = rcp (new CrsMatrix2 (clonedRowMap, clonedColMap,
854  numEntriesPerRow, pftype,
855  matParams));
856  }
857  }
858  else {
859  if (numEntriesPerRow.is_null ()) {
860  clonedMatrix = rcp (new CrsMatrix2 (clonedRowMap, numEntriesForAll,
861  pftype, matParams));
862  }
863  else {
864  clonedMatrix = rcp (new CrsMatrix2 (clonedRowMap, numEntriesPerRow,
865  pftype, matParams));
866  }
867  }
868  // done with these
869  numEntriesPerRow = Teuchos::null;
870  numEntriesForAll = 0;
871 
872  if (useLocalIndices) {
873  clonedMatrix->allocateValues (LocalIndices,
874  CrsMatrix2::GraphNotYetAllocated);
875  if (this->isLocallyIndexed ()) {
876  ArrayView<const LocalOrdinal> linds;
877  ArrayView<const Scalar> vals;
878  for (LocalOrdinal lrow = clonedRowMap->getMinLocalIndex ();
879  lrow <= clonedRowMap->getMaxLocalIndex ();
880  ++lrow) {
881  this->getLocalRowView (lrow, linds, vals);
882  if (linds.size ()) {
883  clonedMatrix->insertLocalValues (lrow, linds, vals);
884  }
885  }
886  }
887  else { // this->isGloballyIndexed()
888  Array<LocalOrdinal> linds;
889  Array<Scalar> vals;
890  for (LocalOrdinal lrow = clonedRowMap->getMinLocalIndex ();
891  lrow <= clonedRowMap->getMaxLocalIndex ();
892  ++lrow) {
893  size_t theNumEntries = this->getNumEntriesInLocalRow (lrow);
894  if (theNumEntries > static_cast<size_t> (linds.size ())) {
895  linds.resize (theNumEntries);
896  }
897  if (theNumEntries > static_cast<size_t> (vals.size ())) {
898  vals.resize (theNumEntries);
899  }
900  this->getLocalRowCopy (clonedRowMap->getGlobalElement (lrow),
901  linds (), vals (), theNumEntries);
902  if (theNumEntries != 0) {
903  clonedMatrix->insertLocalValues (lrow, linds (0, theNumEntries),
904  vals (0, theNumEntries));
905  }
906  }
907  }
908  }
909  else { // useGlobalIndices
910  clonedMatrix->allocateValues (GlobalIndices,
911  CrsMatrix2::GraphNotYetAllocated);
912  if (this->isGloballyIndexed ()) {
913  ArrayView<const GlobalOrdinal> ginds;
914  ArrayView<const Scalar> vals;
915  for (GlobalOrdinal grow = clonedRowMap->getMinGlobalIndex ();
916  grow <= clonedRowMap->getMaxGlobalIndex ();
917  ++grow) {
918  this->getGlobalRowView (grow, ginds, vals);
919  if (ginds.size () > 0) {
920  clonedMatrix->insertGlobalValues (grow, ginds, vals);
921  }
922  }
923  }
924  else { // this->isLocallyIndexed()
925  Array<GlobalOrdinal> ginds;
926  Array<Scalar> vals;
927  for (GlobalOrdinal grow = clonedRowMap->getMinGlobalIndex ();
928  grow <= clonedRowMap->getMaxGlobalIndex ();
929  ++grow) {
930  size_t theNumEntries = this->getNumEntriesInGlobalRow (grow);
931  if (theNumEntries > static_cast<size_t> (ginds.size ())) {
932  ginds.resize (theNumEntries);
933  }
934  if (theNumEntries > static_cast<size_t> (vals.size ())) {
935  vals.resize (theNumEntries);
936  }
937  this->getGlobalRowCopy (grow, ginds (), vals (), theNumEntries);
938  if (theNumEntries != 0) {
939  clonedMatrix->insertGlobalValues (grow, ginds (0, theNumEntries),
940  vals (0, theNumEntries));
941  }
942  }
943  }
944  }
945 
946  if (fillCompleteClone) {
947  RCP<const Map2> clonedRangeMap;
948  RCP<const Map2> clonedDomainMap;
949  try {
950  if (! this->getRangeMap ().is_null () &&
951  this->getRangeMap () != clonedRowMap) {
952  clonedRangeMap = this->getRangeMap ()->template clone<Node2> (node2);
953  }
954  else {
955  clonedRangeMap = clonedRowMap;
956  }
957  if (! this->getDomainMap ().is_null () &&
958  this->getDomainMap () != clonedRowMap) {
959  clonedDomainMap = this->getDomainMap ()->template clone<Node2> (node2);
960  }
961  else {
962  clonedDomainMap = clonedRowMap;
963  }
964  }
965  catch (std::exception &e) {
966  const bool caughtExceptionOnClone = true;
967  TEUCHOS_TEST_FOR_EXCEPTION
968  (caughtExceptionOnClone, std::runtime_error,
969  Teuchos::typeName (*this) << "::clone: Caught the following "
970  "exception while cloning range and domain Maps on a clone of "
971  "type " << Teuchos::typeName (*clonedMatrix) << ": " << e.what ());
972  }
973 
974  RCP<ParameterList> fillparams =
975  params.is_null () ? Teuchos::null : sublist (params, "fillComplete");
976  try {
977  clonedMatrix->fillComplete (clonedDomainMap, clonedRangeMap,
978  fillparams);
979  }
980  catch (std::exception &e) {
981  const bool caughtExceptionOnClone = true;
982  TEUCHOS_TEST_FOR_EXCEPTION(
983  caughtExceptionOnClone, std::runtime_error,
984  Teuchos::typeName (*this) << "::clone: Caught the following "
985  "exception while calling fillComplete() on a clone of type "
986  << Teuchos::typeName (*clonedMatrix) << ": " << e.what ());
987  }
988  }
989  return clonedMatrix;
990  }
991 
993  virtual ~CrsMatrix ();
994 
996 
998 
1022  //
1068  void
1069  insertGlobalValues (const GlobalOrdinal globalRow,
1070  const Teuchos::ArrayView<const GlobalOrdinal>& cols,
1071  const Teuchos::ArrayView<const Scalar>& vals);
1072 
1087  void
1088  insertGlobalValues (const GlobalOrdinal globalRow,
1089  const LocalOrdinal numEnt,
1090  const Scalar vals[],
1091  const GlobalOrdinal inds[]);
1092 
1133  void
1134  insertLocalValues (const LocalOrdinal localRow,
1135  const Teuchos::ArrayView<const LocalOrdinal> &cols,
1136  const Teuchos::ArrayView<const Scalar> &vals);
1137 
1152  void
1153  insertLocalValues (const LocalOrdinal localRow,
1154  const LocalOrdinal numEnt,
1155  const Scalar vals[],
1156  const LocalOrdinal cols[]);
1157 
1158  private:
1169  LocalOrdinal
1170  replaceGlobalValuesImpl (impl_scalar_type rowVals[],
1171  const crs_graph_type& graph,
1172  const RowInfo& rowInfo,
1173  const GlobalOrdinal inds[],
1174  const impl_scalar_type newVals[],
1175  const LocalOrdinal numElts) const;
1176 
1177  public:
1214  template<class GlobalIndicesViewType,
1215  class ImplScalarViewType>
1216  LocalOrdinal
1217  replaceGlobalValues (const GlobalOrdinal globalRow,
1218  const typename UnmanagedView<GlobalIndicesViewType>::type& inputInds,
1219  const typename UnmanagedView<ImplScalarViewType>::type& inputVals) const
1220  {
1221  // We use static_assert here to check the template parameters,
1222  // rather than std::enable_if (e.g., on the return value, to
1223  // enable compilation only if the template parameters match the
1224  // desired attributes). This turns obscure link errors into
1225  // clear compilation errors. It also makes the return value a
1226  // lot easier to see.
1227  static_assert (Kokkos::is_view<GlobalIndicesViewType>::value,
1228  "First template parameter GlobalIndicesViewType must be "
1229  "a Kokkos::View.");
1230  static_assert (Kokkos::is_view<ImplScalarViewType>::value,
1231  "Second template parameter ImplScalarViewType must be a "
1232  "Kokkos::View.");
1233  static_assert (static_cast<int> (GlobalIndicesViewType::rank) == 1,
1234  "First template parameter GlobalIndicesViewType must "
1235  "have rank 1.");
1236  static_assert (static_cast<int> (ImplScalarViewType::rank) == 1,
1237  "Second template parameter ImplScalarViewType must have "
1238  "rank 1.");
1239  static_assert (std::is_same<
1240  typename GlobalIndicesViewType::non_const_value_type,
1241  global_ordinal_type>::value,
1242  "First template parameter GlobalIndicesViewType must "
1243  "contain values of type global_ordinal_type.");
1244  static_assert (std::is_same<
1245  typename ImplScalarViewType::non_const_value_type,
1246  impl_scalar_type>::value,
1247  "Second template parameter ImplScalarViewType must "
1248  "contain values of type impl_scalar_type.");
1249  typedef LocalOrdinal LO;
1250  const LO numInputEnt = inputInds.extent (0);
1251  if (static_cast<LO> (inputVals.extent (0)) != numInputEnt) {
1252  return Teuchos::OrdinalTraits<LO>::invalid ();
1253  }
1254  const Scalar* const inVals =
1255  reinterpret_cast<const Scalar*> (inputVals.data ());
1256  return this->replaceGlobalValues (globalRow, numInputEnt, inVals,
1257  inputInds.data ());
1258  }
1259 
1263  LocalOrdinal
1264  replaceGlobalValues (const GlobalOrdinal globalRow,
1265  const Teuchos::ArrayView<const GlobalOrdinal>& cols,
1266  const Teuchos::ArrayView<const Scalar>& vals) const;
1267 
1283  LocalOrdinal
1284  replaceGlobalValues (const GlobalOrdinal globalRow,
1285  const LocalOrdinal numEnt,
1286  const Scalar vals[],
1287  const GlobalOrdinal cols[]) const;
1288 
1289  private:
1300  LocalOrdinal
1301  replaceLocalValuesImpl (impl_scalar_type rowVals[],
1302  const crs_graph_type& graph,
1303  const RowInfo& rowInfo,
1304  const LocalOrdinal inds[],
1305  const impl_scalar_type newVals[],
1306  const LocalOrdinal numElts) const;
1307 
1308  public:
1344  template<class LocalIndicesViewType,
1345  class ImplScalarViewType>
1346  LocalOrdinal
1347  replaceLocalValues (const LocalOrdinal localRow,
1348  const typename UnmanagedView<LocalIndicesViewType>::type& inputInds,
1349  const typename UnmanagedView<ImplScalarViewType>::type& inputVals) const
1350  {
1351  // We use static_assert here to check the template parameters,
1352  // rather than std::enable_if (e.g., on the return value, to
1353  // enable compilation only if the template parameters match the
1354  // desired attributes). This turns obscure link errors into
1355  // clear compilation errors. It also makes the return value a
1356  // lot easier to see.
1357  static_assert (Kokkos::is_view<LocalIndicesViewType>::value,
1358  "First template parameter LocalIndicesViewType must be "
1359  "a Kokkos::View.");
1360  static_assert (Kokkos::is_view<ImplScalarViewType>::value,
1361  "Second template parameter ImplScalarViewType must be a "
1362  "Kokkos::View.");
1363  static_assert (static_cast<int> (LocalIndicesViewType::rank) == 1,
1364  "First template parameter LocalIndicesViewType must "
1365  "have rank 1.");
1366  static_assert (static_cast<int> (ImplScalarViewType::rank) == 1,
1367  "Second template parameter ImplScalarViewType must have "
1368  "rank 1.");
1369  static_assert (std::is_same<
1370  typename LocalIndicesViewType::non_const_value_type,
1371  local_ordinal_type>::value,
1372  "First template parameter LocalIndicesViewType must "
1373  "contain values of type local_ordinal_type.");
1374  static_assert (std::is_same<
1375  typename ImplScalarViewType::non_const_value_type,
1376  impl_scalar_type>::value,
1377  "Second template parameter ImplScalarViewType must "
1378  "contain values of type impl_scalar_type.");
1379 
1380  typedef LocalOrdinal LO;
1381  const LO numInputEnt = inputInds.extent (0);
1382  if (numInputEnt != inputVals.extent (0)) {
1383  return Teuchos::OrdinalTraits<LO>::invalid ();
1384  }
1385  const Scalar* const inVals =
1386  reinterpret_cast<const Scalar*> (inputVals.data ());
1387  return this->replaceLocalValues (localRow, numInputEnt,
1388  inVals, inputInds.data ());
1389  }
1390 
1394  LocalOrdinal
1395  replaceLocalValues (const LocalOrdinal localRow,
1396  const Teuchos::ArrayView<const LocalOrdinal>& cols,
1397  const Teuchos::ArrayView<const Scalar>& vals) const;
1398 
1416  LocalOrdinal
1417  replaceLocalValues (const LocalOrdinal localRow,
1418  const LocalOrdinal numEnt,
1419  const Scalar inputVals[],
1420  const LocalOrdinal inputCols[]) const;
1421 
1422  private:
1427  static const bool useAtomicUpdatesByDefault =
1428 #ifdef KOKKOS_ENABLE_SERIAL
1429  ! std::is_same<execution_space, Kokkos::Serial>::value;
1430 #else
1431  true;
1432 #endif // KOKKOS_ENABLE_SERIAL
1433 
1457  LocalOrdinal
1458  sumIntoGlobalValuesImpl (impl_scalar_type rowVals[],
1459  const crs_graph_type& graph,
1460  const RowInfo& rowInfo,
1461  const GlobalOrdinal inds[],
1462  const impl_scalar_type newVals[],
1463  const LocalOrdinal numElts,
1464  const bool atomic = useAtomicUpdatesByDefault) const;
1465 
1466  public:
1503  LocalOrdinal
1504  sumIntoGlobalValues (const GlobalOrdinal globalRow,
1505  const Teuchos::ArrayView<const GlobalOrdinal>& cols,
1506  const Teuchos::ArrayView<const Scalar>& vals,
1507  const bool atomic = useAtomicUpdatesByDefault);
1508 
1531  LocalOrdinal
1532  sumIntoGlobalValues (const GlobalOrdinal globalRow,
1533  const LocalOrdinal numEnt,
1534  const Scalar vals[],
1535  const GlobalOrdinal cols[],
1536  const bool atomic = useAtomicUpdatesByDefault);
1537 
1538  private:
1551  LocalOrdinal
1552  sumIntoLocalValuesImpl (impl_scalar_type rowVals[],
1553  const crs_graph_type& graph,
1554  const RowInfo& rowInfo,
1555  const LocalOrdinal inds[],
1556  const impl_scalar_type newVals[],
1557  const LocalOrdinal numElts,
1558  const bool atomic = useAtomicUpdatesByDefault) const;
1559 
1560  public:
1597  template<class LocalIndicesViewType,
1598  class ImplScalarViewType>
1599  LocalOrdinal
1600  sumIntoLocalValues (const LocalOrdinal localRow,
1601  const typename UnmanagedView<LocalIndicesViewType>::type& inputInds,
1602  const typename UnmanagedView<ImplScalarViewType>::type& inputVals,
1603  const bool atomic = useAtomicUpdatesByDefault) const
1604  {
1605  // We use static_assert here to check the template parameters,
1606  // rather than std::enable_if (e.g., on the return value, to
1607  // enable compilation only if the template parameters match the
1608  // desired attributes). This turns obscure link errors into
1609  // clear compilation errors. It also makes the return value a
1610  // lot easier to see.
1611  static_assert (Kokkos::is_view<LocalIndicesViewType>::value,
1612  "First template parameter LocalIndicesViewType must be "
1613  "a Kokkos::View.");
1614  static_assert (Kokkos::is_view<ImplScalarViewType>::value,
1615  "Second template parameter ImplScalarViewType must be a "
1616  "Kokkos::View.");
1617  static_assert (static_cast<int> (LocalIndicesViewType::rank) == 1,
1618  "First template parameter LocalIndicesViewType must "
1619  "have rank 1.");
1620  static_assert (static_cast<int> (ImplScalarViewType::rank) == 1,
1621  "Second template parameter ImplScalarViewType must have "
1622  "rank 1.");
1623  static_assert (std::is_same<
1624  typename LocalIndicesViewType::non_const_value_type,
1625  local_ordinal_type>::value,
1626  "First template parameter LocalIndicesViewType must "
1627  "contain values of type local_ordinal_type.");
1628  static_assert (std::is_same<
1629  typename ImplScalarViewType::non_const_value_type,
1630  impl_scalar_type>::value,
1631  "Second template parameter ImplScalarViewType must "
1632  "contain values of type impl_scalar_type.");
1633  typedef LocalOrdinal LO;
1634  const LO numInputEnt = inputInds.extent (0);
1635  if (static_cast<LO> (inputVals.extent (0)) != numInputEnt) {
1636  return Teuchos::OrdinalTraits<LO>::invalid ();
1637  }
1638  return this->sumIntoLocalValues (localRow,
1639  numInputEnt,
1640  reinterpret_cast<const Scalar*> (inputVals.data ()),
1641  inputInds.data (),
1642  atomic);
1643  }
1644 
1674  LocalOrdinal
1675  sumIntoLocalValues (const LocalOrdinal localRow,
1676  const Teuchos::ArrayView<const LocalOrdinal>& cols,
1677  const Teuchos::ArrayView<const Scalar>& vals,
1678  const bool atomic = useAtomicUpdatesByDefault) const;
1679 
1701  LocalOrdinal
1702  sumIntoLocalValues (const LocalOrdinal localRow,
1703  const LocalOrdinal numEnt,
1704  const Scalar vals[],
1705  const LocalOrdinal cols[],
1706  const bool atomic = useAtomicUpdatesByDefault) const;
1707 
1708  private:
1739  LocalOrdinal
1740  transformLocalValues (impl_scalar_type rowVals[],
1741  const crs_graph_type& graph,
1742  const RowInfo& rowInfo,
1743  const LocalOrdinal inds[],
1744  const impl_scalar_type newVals[],
1745  const LocalOrdinal numElts,
1746  std::function<impl_scalar_type (const impl_scalar_type&, const impl_scalar_type&) > f,
1747  const bool atomic = useAtomicUpdatesByDefault) const;
1748 
1779  LocalOrdinal
1780  transformGlobalValues (impl_scalar_type rowVals[],
1781  const crs_graph_type& graph,
1782  const RowInfo& rowInfo,
1783  const GlobalOrdinal inds[],
1784  const impl_scalar_type newVals[],
1785  const LocalOrdinal numElts,
1786  std::function<impl_scalar_type (const impl_scalar_type&, const impl_scalar_type&) > f,
1787  const bool atomic = useAtomicUpdatesByDefault) const;
1788 
1815  LocalOrdinal
1816  transformLocalValues (const LocalOrdinal lclRow,
1817  const LocalOrdinal numInputEnt,
1818  const impl_scalar_type inputVals[],
1819  const LocalOrdinal inputCols[],
1820  std::function<impl_scalar_type (const impl_scalar_type&, const impl_scalar_type&) > f,
1821  const bool atomic = useAtomicUpdatesByDefault) const;
1822 
1849  LocalOrdinal
1850  transformGlobalValues (const GlobalOrdinal gblRow,
1851  const LocalOrdinal numInputEnt,
1852  const impl_scalar_type inputVals[],
1853  const GlobalOrdinal inputCols[],
1854  std::function<impl_scalar_type (const impl_scalar_type&, const impl_scalar_type&) > f,
1855  const bool atomic = useAtomicUpdatesByDefault) const;
1856 
1857  public:
1901  template<class LocalIndicesViewType,
1902  class ImplScalarViewType,
1903  class BinaryFunction>
1904  LocalOrdinal
1905  transformLocalValues (const LocalOrdinal lclRow,
1906  const typename UnmanagedView<LocalIndicesViewType>::type& inputInds,
1907  const typename UnmanagedView<ImplScalarViewType>::type& inputVals,
1908  BinaryFunction f,
1909  const bool atomic = useAtomicUpdatesByDefault) const
1910  {
1911  // We use static_assert here to check the template parameters,
1912  // rather than std::enable_if (e.g., on the return value, to
1913  // enable compilation only if the template parameters match the
1914  // desired attributes). This turns obscure link errors into
1915  // clear compilation errors. It also makes the return value a
1916  // lot easier to see.
1917  static_assert (Kokkos::is_view<LocalIndicesViewType>::value,
1918  "First template parameter LocalIndicesViewType must be "
1919  "a Kokkos::View.");
1920  static_assert (Kokkos::is_view<ImplScalarViewType>::value,
1921  "Second template parameter ImplScalarViewType must be a "
1922  "Kokkos::View.");
1923  static_assert (static_cast<int> (LocalIndicesViewType::rank) == 1,
1924  "First template parameter LocalIndicesViewType must "
1925  "have rank 1.");
1926  static_assert (static_cast<int> (ImplScalarViewType::rank) == 1,
1927  "Second template parameter ImplScalarViewType must have "
1928  "rank 1.");
1929  static_assert (std::is_same<
1930  typename LocalIndicesViewType::non_const_value_type,
1931  local_ordinal_type>::value,
1932  "First template parameter LocalIndicesViewType must "
1933  "contain values of type local_ordinal_type.");
1934  static_assert (std::is_same<
1935  typename ImplScalarViewType::non_const_value_type,
1936  impl_scalar_type>::value,
1937  "Second template parameter ImplScalarViewType must "
1938  "contain values of type impl_scalar_type.");
1939  typedef LocalOrdinal LO;
1940  const LO numInputEnt = inputInds.extent (0);
1941  if (static_cast<LO> (inputVals.extent (0)) != numInputEnt) {
1942  return Teuchos::OrdinalTraits<LO>::invalid ();
1943  }
1944  return this->transformLocalValues (lclRow,
1945  numInputEnt,
1946  inputVals.data (),
1947  inputInds.data (),
1948  f,
1949  atomic);
1950  }
1951 
1993  template<class BinaryFunction, class InputMemorySpace>
1994  LocalOrdinal
1995  transformGlobalValues (const GlobalOrdinal gblRow,
1996  const Kokkos::View<const GlobalOrdinal*,
1997  InputMemorySpace,
1998  Kokkos::MemoryUnmanaged>& inputInds,
1999  const Kokkos::View<const impl_scalar_type*,
2000  InputMemorySpace,
2001  Kokkos::MemoryUnmanaged>& inputVals,
2002  BinaryFunction f,
2003  const bool atomic = useAtomicUpdatesByDefault) const
2004  {
2005  typedef LocalOrdinal LO;
2006  const LO numInputEnt = inputInds.extent (0);
2007  if (static_cast<LO> (inputVals.extent (0)) != numInputEnt) {
2008  return Teuchos::OrdinalTraits<LO>::invalid ();
2009  }
2010  return this->transformGlobalValues (gblRow,
2011  numInputEnt,
2012  inputVals.data (),
2013  inputInds.data (),
2014  f,
2015  atomic);
2016  }
2017 
2019  void setAllToScalar (const Scalar& alpha);
2020 
2022  void scale (const Scalar& alpha);
2023 
2047  void
2048  setAllValues (const typename local_matrix_type::row_map_type& ptr,
2049  const typename local_graph_type::entries_type::non_const_type& ind,
2050  const typename local_matrix_type::values_type& val);
2051 
2075  void
2076  setAllValues (const Teuchos::ArrayRCP<size_t>& ptr,
2077  const Teuchos::ArrayRCP<LocalOrdinal>& ind,
2078  const Teuchos::ArrayRCP<Scalar>& val);
2079 
2080  void
2081  getAllValues (Teuchos::ArrayRCP<const size_t>& rowPointers,
2082  Teuchos::ArrayRCP<const LocalOrdinal>& columnIndices,
2083  Teuchos::ArrayRCP<const Scalar>& values) const;
2084 
2086 
2088 
2117  void globalAssemble();
2118 
2132  void resumeFill (const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
2133 
2191  void
2192  fillComplete (const Teuchos::RCP<const map_type>& domainMap,
2193  const Teuchos::RCP<const map_type>& rangeMap,
2194  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
2195 
2222  void
2223  fillComplete (const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
2224 
2251  void
2252  expertStaticFillComplete (const Teuchos::RCP<const map_type>& domainMap,
2253  const Teuchos::RCP<const map_type>& rangeMap,
2254  const Teuchos::RCP<const import_type>& importer = Teuchos::null,
2255  const Teuchos::RCP<const export_type>& exporter = Teuchos::null,
2256  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
2257 
2267  void
2268  replaceColMap (const Teuchos::RCP<const map_type>& newColMap);
2269 
2351  void
2352  reindexColumns (crs_graph_type* const graph,
2353  const Teuchos::RCP<const map_type>& newColMap,
2354  const Teuchos::RCP<const import_type>& newImport = Teuchos::null,
2355  const bool sortEachRow = true);
2356 
2369  void
2370  replaceDomainMapAndImporter (const Teuchos::RCP<const map_type>& newDomainMap,
2371  Teuchos::RCP<const import_type>& newImporter);
2372 
2386  virtual void
2387  removeEmptyProcessesInPlace (const Teuchos::RCP<const map_type>& newMap) override;
2388 
2390 
2392 
2394  Teuchos::RCP<const Teuchos::Comm<int> > getComm() const override;
2395 
2397  Teuchos::RCP<node_type> getNode () const override;
2398 
2400  Teuchos::RCP<const map_type> getRowMap () const override;
2401 
2403  Teuchos::RCP<const map_type> getColMap () const override;
2404 
2406  Teuchos::RCP<const RowGraph<LocalOrdinal, GlobalOrdinal, Node> >
2407  getGraph () const override;
2408 
2410  Teuchos::RCP<const crs_graph_type> getCrsGraph () const;
2411 
2412  private:
2423  const crs_graph_type& getCrsGraphRef () const;
2424 
2425  public:
2436 
2456  global_size_t getGlobalNumRows() const override;
2457 
2463  global_size_t getGlobalNumCols() const override;
2464 
2471  size_t getNodeNumRows() const override;
2472 
2476  size_t getNodeNumCols() const override;
2477 
2479  GlobalOrdinal getIndexBase() const override;
2480 
2482  global_size_t getGlobalNumEntries() const override;
2483 
2485  size_t getNodeNumEntries() const override;
2486 
2493  size_t getNumEntriesInGlobalRow (GlobalOrdinal globalRow) const override;
2494 
2501  size_t getNumEntriesInLocalRow (LocalOrdinal localRow) const override;
2502 
2514  global_size_t getGlobalNumDiagsImpl () const override;
2515 
2524 
2536  size_t getNodeNumDiagsImpl () const override;
2537 
2545  size_t TPETRA_DEPRECATED getNodeNumDiags () const override;
2546 
2554  size_t getGlobalMaxNumRowEntries () const override;
2555 
2563  size_t getNodeMaxNumRowEntries () const override;
2564 
2566  bool hasColMap () const override;
2567 
2579  bool isLowerTriangularImpl () const override;
2580 
2591  bool TPETRA_DEPRECATED isLowerTriangular () const override;
2592 
2604  bool isUpperTriangularImpl () const override;
2605 
2616  bool TPETRA_DEPRECATED isUpperTriangular () const override;
2617 
2638  bool isLocallyIndexed() const override;
2639 
2659  bool isGloballyIndexed() const override;
2660 
2683  bool isFillComplete() const override;
2684 
2707  bool isFillActive() const;
2708 
2710 
2716  bool isStorageOptimized () const;
2717 
2719  ProfileType getProfileType () const;
2720 
2722  bool isStaticGraph () const;
2723 
2735  mag_type getFrobeniusNorm () const override;
2736 
2739  virtual bool supportsRowViews () const override;
2740 
2789  void
2790  getGlobalRowCopy (GlobalOrdinal GlobalRow,
2791  const Teuchos::ArrayView<GlobalOrdinal>& Indices,
2792  const Teuchos::ArrayView<Scalar>& Values,
2793  size_t& NumEntries) const override;
2794 
2812  void
2813  getLocalRowCopy (LocalOrdinal localRow,
2814  const Teuchos::ArrayView<LocalOrdinal>& colInds,
2815  const Teuchos::ArrayView<Scalar>& vals,
2816  size_t& numEntries) const override;
2817 
2831  void
2832  getGlobalRowView (GlobalOrdinal GlobalRow,
2833  Teuchos::ArrayView<const GlobalOrdinal>& indices,
2834  Teuchos::ArrayView<const Scalar>& values) const override;
2835 
2849  void
2850  getLocalRowView (LocalOrdinal LocalRow,
2851  Teuchos::ArrayView<const LocalOrdinal>& indices,
2852  Teuchos::ArrayView<const Scalar>& values) const override;
2853 
2880  LocalOrdinal
2881  getLocalRowViewRaw (const LocalOrdinal lclRow,
2882  LocalOrdinal& numEnt,
2883  const LocalOrdinal*& lclColInds,
2884  const Scalar*& vals) const override;
2885 
2911  LocalOrdinal
2912  getLocalRowView (const LocalOrdinal lclRow,
2913  LocalOrdinal& numEnt,
2914  const impl_scalar_type*& val,
2915  const LocalOrdinal*& ind) const;
2916 
2924  template<class OutputScalarType>
2925  typename std::enable_if<! std::is_same<OutputScalarType, impl_scalar_type>::value &&
2926  std::is_convertible<impl_scalar_type, OutputScalarType>::value,
2927  LocalOrdinal>::type
2928  getLocalRowView (const LocalOrdinal lclRow,
2929  LocalOrdinal& numEnt,
2930  const OutputScalarType*& val,
2931  const LocalOrdinal*& ind) const
2932  {
2933  const impl_scalar_type* valTmp = NULL;
2934  const LocalOrdinal err = this->getLocalRowView (lclRow, numEnt, valTmp, ind);
2935  // Cast is legitimate because impl_scalar_type is convertible to
2936  // OutputScalarType.
2937  val = reinterpret_cast<const OutputScalarType*> (valTmp);
2938  return err;
2939  }
2940 
2947  void
2949 
2993  void getLocalDiagOffsets (Teuchos::ArrayRCP<size_t>& offsets) const;
2994 
3016  void
3018  const Kokkos::View<const size_t*, device_type,
3019  Kokkos::MemoryUnmanaged>& offsets) const;
3020 
3043  void
3045  const Teuchos::ArrayView<const size_t>& offsets) const;
3046 
3051  void
3053 
3058  void
3060 
3062 
3064 
3113  template <class DomainScalar, class RangeScalar>
3114  void
3117  Teuchos::ETransp mode,
3118  RangeScalar alpha,
3119  RangeScalar beta) const
3120  {
3121  using Teuchos::NO_TRANS;
3122  // Just like Scalar and impl_scalar_type may differ in CrsMatrix,
3123  // RangeScalar and its corresponding impl_scalar_type may differ in
3124  // MultiVector.
3125  typedef typename MultiVector<RangeScalar, LocalOrdinal, GlobalOrdinal,
3126  Node>::impl_scalar_type range_impl_scalar_type;
3127 #ifdef HAVE_TPETRA_DEBUG
3128  const char tfecfFuncName[] = "localMultiply: ";
3129 #endif // HAVE_TPETRA_DEBUG
3130 
3131  const range_impl_scalar_type theAlpha = static_cast<range_impl_scalar_type> (alpha);
3132  const range_impl_scalar_type theBeta = static_cast<range_impl_scalar_type> (beta);
3133  const bool conjugate = (mode == Teuchos::CONJ_TRANS);
3134  const bool transpose = (mode != Teuchos::NO_TRANS);
3135  auto X_lcl = X.template getLocalView<device_type> ();
3136  auto Y_lcl = Y.template getLocalView<device_type> ();
3137 
3138 #ifdef HAVE_TPETRA_DEBUG
3139  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC
3140  (X.getNumVectors () != Y.getNumVectors (), std::runtime_error,
3141  "X.getNumVectors() = " << X.getNumVectors () << " != Y.getNumVectors() = "
3142  << Y.getNumVectors () << ".");
3143  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC
3144  (! transpose && X.getLocalLength () != getColMap ()->getNodeNumElements (),
3145  std::runtime_error, "NO_TRANS case: X has the wrong number of local rows. "
3146  "X.getLocalLength() = " << X.getLocalLength () << " != getColMap()->"
3147  "getNodeNumElements() = " << getColMap ()->getNodeNumElements () << ".");
3148  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC
3149  (! transpose && Y.getLocalLength () != getRowMap ()->getNodeNumElements (),
3150  std::runtime_error, "NO_TRANS case: Y has the wrong number of local rows. "
3151  "Y.getLocalLength() = " << Y.getLocalLength () << " != getRowMap()->"
3152  "getNodeNumElements() = " << getRowMap ()->getNodeNumElements () << ".");
3153  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC
3154  (transpose && X.getLocalLength () != getRowMap ()->getNodeNumElements (),
3155  std::runtime_error, "TRANS or CONJ_TRANS case: X has the wrong number of "
3156  "local rows. X.getLocalLength() = " << X.getLocalLength () << " != "
3157  "getRowMap()->getNodeNumElements() = "
3158  << getRowMap ()->getNodeNumElements () << ".");
3159  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC
3160  (transpose && Y.getLocalLength () != getColMap ()->getNodeNumElements (),
3161  std::runtime_error, "TRANS or CONJ_TRANS case: X has the wrong number of "
3162  "local rows. Y.getLocalLength() = " << Y.getLocalLength () << " != "
3163  "getColMap()->getNodeNumElements() = "
3164  << getColMap ()->getNodeNumElements () << ".");
3165  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC
3166  (! isFillComplete (), std::runtime_error, "The matrix is not fill "
3167  "complete. You must call fillComplete() (possibly with domain and range "
3168  "Map arguments) without an intervening resumeFill() call before you may "
3169  "call this method.");
3170  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC
3171  (! X.isConstantStride () || ! Y.isConstantStride (), std::runtime_error,
3172  "X and Y must be constant stride.");
3173  // If the two pointers are NULL, then they don't alias one
3174  // another, even though they are equal.
3175  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(
3176  X_lcl.data () == Y_lcl.data () &&
3177  X_lcl.data () != NULL,
3178  std::runtime_error, "X and Y may not alias one another.");
3179 #endif // HAVE_TPETRA_DEBUG
3180 
3181  // Y = alpha*op(M) + beta*Y
3182  if (transpose) {
3183  KokkosSparse::spmv (conjugate ? KokkosSparse::ConjugateTranspose : KokkosSparse::Transpose,
3184  theAlpha,
3185  lclMatrix_,
3186  X.template getLocalView<device_type> (),
3187  theBeta,
3188  Y.template getLocalView<device_type> ());
3189  }
3190  else {
3191  KokkosSparse::spmv (KokkosSparse::NoTranspose,
3192  theAlpha,
3193  lclMatrix_,
3194  X.template getLocalView<device_type> (),
3195  theBeta,
3196  Y.template getLocalView<device_type> ());
3197  }
3198  }
3199 
3200  private:
3201 
3232  void
3235  const Teuchos::ETransp mode = Teuchos::NO_TRANS,
3236  const Scalar& alpha = Teuchos::ScalarTraits<Scalar>::one (),
3237  const Scalar& beta = Teuchos::ScalarTraits<Scalar>::zero ()) const;
3238 
3239  public:
3240 
3265  template <class DomainScalar, class RangeScalar>
3266  void
3270  const RangeScalar& dampingFactor,
3271  const ESweepDirection direction) const
3272  {
3273  typedef LocalOrdinal LO;
3274  typedef GlobalOrdinal GO;
3277  typedef typename Node::device_type::memory_space dev_mem_space;
3278  typedef Kokkos::HostSpace host_mem_space;
3279  typedef typename Graph::local_graph_type k_local_graph_type;
3280  typedef typename k_local_graph_type::size_type offset_type;
3281  const char prefix[] = "Tpetra::CrsMatrix::localGaussSeidel: ";
3282 
3283  TEUCHOS_TEST_FOR_EXCEPTION
3284  (! this->isFillComplete (), std::runtime_error,
3285  prefix << "The matrix is not fill complete.");
3286  const size_t lclNumRows = this->getNodeNumRows ();
3287  const size_t numVecs = B.getNumVectors ();
3288  TEUCHOS_TEST_FOR_EXCEPTION
3289  (X.getNumVectors () != numVecs, std::invalid_argument,
3290  prefix << "B.getNumVectors() = " << numVecs << " != "
3291  "X.getNumVectors() = " << X.getNumVectors () << ".");
3292  TEUCHOS_TEST_FOR_EXCEPTION
3293  (B.getLocalLength () != lclNumRows, std::invalid_argument,
3294  prefix << "B.getLocalLength() = " << B.getLocalLength ()
3295  << " != this->getNodeNumRows() = " << lclNumRows << ".");
3296 
3297  // mfh 28 Aug 2017: The current local Gauss-Seidel kernel only
3298  // runs on host. (See comments below.) Thus, we need to access
3299  // the host versions of these data.
3300  const_cast<DMV&> (B).template sync<host_mem_space> ();
3301  X.template sync<host_mem_space> ();
3302  X.template modify<host_mem_space> ();
3303  const_cast<MMV&> (D).template sync<host_mem_space> ();
3304 
3305  auto B_lcl = B.template getLocalView<host_mem_space> ();
3306  auto X_lcl = X.template getLocalView<host_mem_space> ();
3307  auto D_lcl = D.template getLocalView<host_mem_space> ();
3308 
3309  offset_type B_stride[8], X_stride[8], D_stride[8];
3310  B_lcl.stride (B_stride);
3311  X_lcl.stride (X_stride);
3312  D_lcl.stride (D_stride);
3313 
3314  local_matrix_type lclMatrix = this->getLocalMatrix ();
3315  k_local_graph_type lclGraph = lclMatrix.graph;
3316  typename local_matrix_type::row_map_type ptr = lclGraph.row_map;
3317  typename local_matrix_type::index_type ind = lclGraph.entries;
3318  typename local_matrix_type::values_type val = lclMatrix.values;
3319  const offset_type* const ptrRaw = ptr.data ();
3320  const LO* const indRaw = ind.data ();
3321  const impl_scalar_type* const valRaw = val.data ();
3322 
3323  const std::string dir ((direction == Forward) ? "F" : "B");
3324  // NOTE (mfh 28 Aug 2017) This assumes UVM. We can't get around
3325  // that on GPUs without using a GPU-based sparse triangular
3326  // solve to implement Gauss-Seidel. This exists in cuSPARSE,
3327  // but we would need to implement a wrapper with a fall-back
3328  // algorithm for unsupported Scalar and LO types.
3329  KokkosSparse::Impl::Sequential::gaussSeidel (static_cast<LO> (lclNumRows),
3330  static_cast<LO> (numVecs),
3331  ptrRaw, indRaw, valRaw,
3332  B_lcl.data (), B_stride[1],
3333  X_lcl.data (), X_stride[1],
3334  D_lcl.data (),
3335  static_cast<impl_scalar_type> (dampingFactor),
3336  dir.c_str ());
3337  const_cast<DMV&> (B).template sync<dev_mem_space> ();
3338  X.template sync<dev_mem_space> ();
3339  const_cast<MMV&> (D).template sync<dev_mem_space> ();
3340  }
3341 
3368  template <class DomainScalar, class RangeScalar>
3369  void
3373  const Teuchos::ArrayView<LocalOrdinal>& rowIndices,
3374  const RangeScalar& dampingFactor,
3375  const ESweepDirection direction) const
3376  {
3377  typedef LocalOrdinal LO;
3378  typedef GlobalOrdinal GO;
3381  typedef typename Node::device_type::memory_space dev_mem_space;
3382  typedef Kokkos::HostSpace host_mem_space;
3383  typedef typename Graph::local_graph_type k_local_graph_type;
3384  typedef typename k_local_graph_type::size_type offset_type;
3385  const char prefix[] = "Tpetra::CrsMatrix::reorderedLocalGaussSeidel: ";
3386 
3387  TEUCHOS_TEST_FOR_EXCEPTION
3388  (! this->isFillComplete (), std::runtime_error,
3389  prefix << "The matrix is not fill complete.");
3390  const size_t lclNumRows = this->getNodeNumRows ();
3391  const size_t numVecs = B.getNumVectors ();
3392  TEUCHOS_TEST_FOR_EXCEPTION
3393  (X.getNumVectors () != numVecs, std::invalid_argument,
3394  prefix << "B.getNumVectors() = " << numVecs << " != "
3395  "X.getNumVectors() = " << X.getNumVectors () << ".");
3396  TEUCHOS_TEST_FOR_EXCEPTION
3397  (B.getLocalLength () != lclNumRows, std::invalid_argument,
3398  prefix << "B.getLocalLength() = " << B.getLocalLength ()
3399  << " != this->getNodeNumRows() = " << lclNumRows << ".");
3400  TEUCHOS_TEST_FOR_EXCEPTION
3401  (static_cast<size_t> (rowIndices.size ()) < lclNumRows,
3402  std::invalid_argument, prefix << "rowIndices.size() = "
3403  << rowIndices.size () << " < this->getNodeNumRows() = "
3404  << lclNumRows << ".");
3405 
3406  // mfh 28 Aug 2017: The current local Gauss-Seidel kernel only
3407  // runs on host. (See comments below.) Thus, we need to access
3408  // the host versions of these data.
3409  const_cast<DMV&> (B).template sync<host_mem_space> ();
3410  X.template sync<host_mem_space> ();
3411  X.template modify<host_mem_space> ();
3412  const_cast<MMV&> (D).template sync<host_mem_space> ();
3413 
3414  auto B_lcl = B.template getLocalView<host_mem_space> ();
3415  auto X_lcl = X.template getLocalView<host_mem_space> ();
3416  auto D_lcl = D.template getLocalView<host_mem_space> ();
3417 
3418  offset_type B_stride[8], X_stride[8], D_stride[8];
3419  B_lcl.stride (B_stride);
3420  X_lcl.stride (X_stride);
3421  D_lcl.stride (D_stride);
3422 
3423  local_matrix_type lclMatrix = this->getLocalMatrix ();
3424  typename Graph::local_graph_type lclGraph = lclMatrix.graph;
3425  typename local_matrix_type::index_type ind = lclGraph.entries;
3426  typename local_matrix_type::row_map_type ptr = lclGraph.row_map;
3427  typename local_matrix_type::values_type val = lclMatrix.values;
3428  const offset_type* const ptrRaw = ptr.data ();
3429  const LO* const indRaw = ind.data ();
3430  const impl_scalar_type* const valRaw = val.data ();
3431 
3432  const std::string dir = (direction == Forward) ? "F" : "B";
3433  // NOTE (mfh 28 Aug 2017) This assumes UVM. We can't get around
3434  // that on GPUs without using a GPU-based sparse triangular
3435  // solve to implement Gauss-Seidel, and also handling the
3436  // permutations correctly.
3437  KokkosSparse::Impl::Sequential::reorderedGaussSeidel (static_cast<LO> (lclNumRows),
3438  static_cast<LO> (numVecs),
3439  ptrRaw, indRaw, valRaw,
3440  B_lcl.data (),
3441  B_stride[1],
3442  X_lcl.data (),
3443  X_stride[1],
3444  D_lcl.data (),
3445  rowIndices.getRawPtr (),
3446  static_cast<LO> (lclNumRows),
3447  static_cast<impl_scalar_type> (dampingFactor),
3448  dir.c_str ());
3449  const_cast<DMV&> (B).template sync<dev_mem_space> ();
3450  X.template sync<dev_mem_space> ();
3451  const_cast<MMV&> (D).template sync<dev_mem_space> ();
3452  }
3453 
3456  template <class T>
3457  Teuchos::RCP<CrsMatrix<T, LocalOrdinal, GlobalOrdinal, Node> >
3458  convert () const;
3459 
3461 
3463 
3474  void
3477  Teuchos::ETransp mode = Teuchos::NO_TRANS,
3478  Scalar alpha = Teuchos::ScalarTraits<Scalar>::one(),
3479  Scalar beta = Teuchos::ScalarTraits<Scalar>::zero()) const override;
3480 
3483  bool hasTransposeApply () const override;
3484 
3491  Teuchos::RCP<const map_type> getDomainMap () const override;
3492 
3499  Teuchos::RCP<const map_type> getRangeMap () const override;
3500 
3502 
3504 
3569  void
3573  const Scalar& dampingFactor,
3574  const ESweepDirection direction,
3575  const int numSweeps) const;
3576 
3643  void
3647  const Teuchos::ArrayView<LocalOrdinal>& rowIndices,
3648  const Scalar& dampingFactor,
3649  const ESweepDirection direction,
3650  const int numSweeps) const;
3651 
3680  void
3684  const Scalar& dampingFactor,
3685  const ESweepDirection direction,
3686  const int numSweeps,
3687  const bool zeroInitialGuess) const;
3688 
3718  void
3722  const Teuchos::ArrayView<LocalOrdinal>& rowIndices,
3723  const Scalar& dampingFactor,
3724  const ESweepDirection direction,
3725  const int numSweeps,
3726  const bool zeroInitialGuess) const;
3727 
3738  virtual Teuchos::RCP<RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
3739  add (const Scalar& alpha,
3741  const Scalar& beta,
3742  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >& domainMap,
3743  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >& rangeMap,
3744  const Teuchos::RCP<Teuchos::ParameterList>& params) const override;
3745 
3747 
3749 
3751  std::string description () const override;
3752 
3755  void
3756  describe (Teuchos::FancyOStream& out,
3757  const Teuchos::EVerbosityLevel verbLevel =
3758  Teuchos::Describable::verbLevel_default) const override;
3759 
3761 
3763 
3768  typedef typename DistObject<Scalar, LocalOrdinal, GlobalOrdinal,
3770 
3771  virtual bool
3772  checkSizes (const SrcDistObject& source) override;
3773 
3787  virtual bool
3788  useNewInterface () override;
3789 
3790  private:
3791 
3792  void
3793  copyAndPermuteImpl (const RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& source,
3794  const size_t numSameIDs,
3795  const LocalOrdinal permuteToLIDs[],
3796  const LocalOrdinal permuteFromLIDs[],
3797  const size_t numPermutes);
3798 
3799  protected:
3800  virtual void
3801  copyAndPermuteNew (const SrcDistObject& source,
3802  const size_t numSameIDs,
3803  const Kokkos::DualView<const local_ordinal_type*, device_type>& permuteToLIDs,
3804  const Kokkos::DualView<const local_ordinal_type*, device_type>& permuteFromLIDs) override;
3805 
3806  virtual void
3807  packAndPrepareNew (const SrcDistObject& source,
3808  const Kokkos::DualView<const local_ordinal_type*, device_type>& exportLIDs,
3809  Kokkos::DualView<char*, buffer_device_type>& exports,
3810  const Kokkos::DualView<size_t*, buffer_device_type>& numPacketsPerLID,
3811  size_t& constantNumPackets,
3812  Distributor& distor) override;
3813 
3814  private:
3817  void
3818  unpackAndCombineNewImpl (const Kokkos::DualView<const LocalOrdinal*, device_type>& importLIDs,
3819  const Kokkos::DualView<const char*, buffer_device_type>& imports,
3820  const Kokkos::DualView<const size_t*, buffer_device_type>& numPacketsPerLID,
3821  const size_t constantNumPackets,
3822  Distributor& distor,
3823  const CombineMode combineMode,
3824  const bool atomic = useAtomicUpdatesByDefault);
3827  void
3828  unpackAndCombineNewImplNonStatic (const Kokkos::DualView<const LocalOrdinal*, device_type>& importLIDs,
3829  const Kokkos::DualView<const char*, buffer_device_type>& imports,
3830  const Kokkos::DualView<const size_t*, buffer_device_type>& numPacketsPerLID,
3831  const size_t constantNumPackets,
3832  Distributor& distor,
3833  const CombineMode combineMode);
3834 
3835  public:
3846  void
3847  unpackAndCombineNew (const Kokkos::DualView<const local_ordinal_type*, device_type>& importLIDs,
3848  const Kokkos::DualView<const char*, buffer_device_type>& imports,
3849  const Kokkos::DualView<const size_t*, buffer_device_type>& numPacketsPerLID,
3850  const size_t constantNumPackets,
3851  Distributor& distor,
3852  const CombineMode CM) override;
3853 
3961  void
3962  packNew (const Kokkos::DualView<const local_ordinal_type*, device_type>& exportLIDs,
3963  Kokkos::DualView<char*, buffer_device_type>& exports,
3964  const Kokkos::DualView<size_t*, buffer_device_type>& numPacketsPerLID,
3965  size_t& constantNumPackets,
3966  Distributor& dist) const;
3967 
3968  private:
3975  void
3976  packNonStaticNew (const Kokkos::DualView<const local_ordinal_type*, device_type>& exportLIDs,
3977  Kokkos::DualView<char*, buffer_device_type>& exports,
3978  const Kokkos::DualView<size_t*, buffer_device_type>& numPacketsPerLID,
3979  size_t& constantNumPackets,
3980  Distributor& distor) const;
3981 
4011  size_t
4012  packRow (char exports[],
4013  const size_t offset,
4014  const size_t numEnt,
4015  const GlobalOrdinal gidsIn[],
4016  const impl_scalar_type valsIn[],
4017  const size_t numBytesPerValue) const;
4018 
4042  bool
4043  packRowStatic (char* const numEntOut,
4044  char* const valOut,
4045  char* const indOut,
4046  const size_t numEnt,
4047  const LocalOrdinal lclRow) const;
4048 
4074  size_t
4075  unpackRow (GlobalOrdinal gidsOut[],
4076  impl_scalar_type valsOut[],
4077  const char imports[],
4078  const size_t offset,
4079  const size_t numBytes,
4080  const size_t numEnt,
4081  const size_t numBytesPerValue);
4082 
4091  void
4092  allocatePackSpaceNew (Kokkos::DualView<char*, buffer_device_type>& exports,
4093  size_t& totalNumEntries,
4094  const Kokkos::DualView<const local_ordinal_type*, device_type>& exportLIDs) const;
4096 
4097  public:
4099  typename local_matrix_type::values_type getLocalValuesView () const {
4100  return k_values1D_;
4101  }
4102 
4103  private:
4104  // Friend declaration for nonmember function.
4105  template<class CrsMatrixType>
4106  friend Teuchos::RCP<CrsMatrixType>
4107  Tpetra::importAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
4108  const Import<typename CrsMatrixType::local_ordinal_type,
4109  typename CrsMatrixType::global_ordinal_type,
4110  typename CrsMatrixType::node_type>& importer,
4111  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4112  typename CrsMatrixType::global_ordinal_type,
4113  typename CrsMatrixType::node_type> >& domainMap,
4114  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4115  typename CrsMatrixType::global_ordinal_type,
4116  typename CrsMatrixType::node_type> >& rangeMap,
4117  const Teuchos::RCP<Teuchos::ParameterList>& params);
4118 
4119  // Friend declaration for nonmember function.
4120  template<class CrsMatrixType>
4121  friend Teuchos::RCP<CrsMatrixType>
4122  Tpetra::importAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
4123  const Import<typename CrsMatrixType::local_ordinal_type,
4124  typename CrsMatrixType::global_ordinal_type,
4125  typename CrsMatrixType::node_type>& rowImporter,
4126  const Import<typename CrsMatrixType::local_ordinal_type,
4127  typename CrsMatrixType::global_ordinal_type,
4128  typename CrsMatrixType::node_type>& domainImporter,
4129  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4130  typename CrsMatrixType::global_ordinal_type,
4131  typename CrsMatrixType::node_type> >& domainMap,
4132  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4133  typename CrsMatrixType::global_ordinal_type,
4134  typename CrsMatrixType::node_type> >& rangeMap,
4135  const Teuchos::RCP<Teuchos::ParameterList>& params);
4136 
4137 
4138  // Friend declaration for nonmember function.
4139  template<class CrsMatrixType>
4140  friend Teuchos::RCP<CrsMatrixType>
4141  Tpetra::exportAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
4142  const Export<typename CrsMatrixType::local_ordinal_type,
4143  typename CrsMatrixType::global_ordinal_type,
4144  typename CrsMatrixType::node_type>& exporter,
4145  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4146  typename CrsMatrixType::global_ordinal_type,
4147  typename CrsMatrixType::node_type> >& domainMap,
4148  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4149  typename CrsMatrixType::global_ordinal_type,
4150  typename CrsMatrixType::node_type> >& rangeMap,
4151  const Teuchos::RCP<Teuchos::ParameterList>& params);
4152 
4153  // Friend declaration for nonmember function.
4154  template<class CrsMatrixType>
4155  friend Teuchos::RCP<CrsMatrixType>
4156  Tpetra::exportAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
4157  const Export<typename CrsMatrixType::local_ordinal_type,
4158  typename CrsMatrixType::global_ordinal_type,
4159  typename CrsMatrixType::node_type>& rowExporter,
4160  const Export<typename CrsMatrixType::local_ordinal_type,
4161  typename CrsMatrixType::global_ordinal_type,
4162  typename CrsMatrixType::node_type>& domainExporter,
4163  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4164  typename CrsMatrixType::global_ordinal_type,
4165  typename CrsMatrixType::node_type> >& domainMap,
4166  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4167  typename CrsMatrixType::global_ordinal_type,
4168  typename CrsMatrixType::node_type> >& rangeMap,
4169  const Teuchos::RCP<Teuchos::ParameterList>& params);
4170 
4171  public:
4187  void
4189  const import_type& importer,
4190  const Teuchos::RCP<const map_type>& domainMap,
4191  const Teuchos::RCP<const map_type>& rangeMap,
4192  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) const;
4193 
4209  void
4211  const import_type& rowImporter,
4212  const import_type& domainImporter,
4213  const Teuchos::RCP<const map_type>& domainMap,
4214  const Teuchos::RCP<const map_type>& rangeMap,
4215  const Teuchos::RCP<Teuchos::ParameterList>& params) const;
4216 
4217 
4233  void
4235  const export_type& exporter,
4236  const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
4237  const Teuchos::RCP<const map_type>& rangeMap = Teuchos::null,
4238  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) const;
4239 
4255  void
4257  const export_type& rowExporter,
4258  const export_type& domainExporter,
4259  const Teuchos::RCP<const map_type>& domainMap,
4260  const Teuchos::RCP<const map_type>& rangeMap,
4261  const Teuchos::RCP<Teuchos::ParameterList>& params) const;
4262 
4263 
4264  private:
4285  void
4286  transferAndFillComplete (Teuchos::RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >& destMatrix,
4287  const ::Tpetra::Details::Transfer<LocalOrdinal, GlobalOrdinal, Node>& rowTransfer,
4288  const Teuchos::RCP<const ::Tpetra::Details::Transfer<LocalOrdinal, GlobalOrdinal, Node> > & domainTransfer,
4289  const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
4290  const Teuchos::RCP<const map_type>& rangeMap = Teuchos::null,
4291  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) const;
4292  // We forbid copy construction by declaring this method private
4293  // and not implementing it.
4295 
4296  // We forbid assignment (operator=) by declaring this method
4297  // private and not implementing it.
4300 
4312  void
4313  insertGlobalValuesImpl (crs_graph_type& graph,
4314  RowInfo& rowInfo,
4315  const GlobalOrdinal gblColInds[],
4316  const impl_scalar_type vals[],
4317  const size_t numInputEnt);
4318 
4328  void
4329  insertGlobalValuesFiltered (const GlobalOrdinal globalRow,
4330  const Teuchos::ArrayView<const GlobalOrdinal>& indices,
4331  const Teuchos::ArrayView<const Scalar>& values);
4332 
4344  void
4345  combineGlobalValues (const GlobalOrdinal globalRowIndex,
4346  const Teuchos::ArrayView<const GlobalOrdinal>& columnIndices,
4347  const Teuchos::ArrayView<const Scalar>& values,
4348  const Tpetra::CombineMode combineMode);
4349 
4367  LocalOrdinal
4368  combineGlobalValuesRaw (const LocalOrdinal lclRow,
4369  const LocalOrdinal numEnt,
4370  const impl_scalar_type vals[],
4371  const GlobalOrdinal cols[],
4372  const Tpetra::CombineMode combineMode);
4373 
4385  template<class BinaryFunction>
4386  LocalOrdinal
4387  transformGlobalValues (const GlobalOrdinal globalRow,
4388  const Teuchos::ArrayView<const GlobalOrdinal>& indices,
4389  const Teuchos::ArrayView<const Scalar>& values,
4390  BinaryFunction f,
4391  const bool atomic = useAtomicUpdatesByDefault) const
4392  {
4393  typedef impl_scalar_type IST;
4394  typedef LocalOrdinal LO;
4395  typedef GlobalOrdinal GO;
4396 
4397  const LO numInputEnt = static_cast<LO> (indices.size ());
4398  if (static_cast<LO> (values.size ()) != numInputEnt) {
4399  return Teuchos::OrdinalTraits<LO>::invalid ();
4400  }
4401 
4402  const GO* const inputCols = indices.getRawPtr ();
4403  const IST* const inputVals =
4404  reinterpret_cast<const IST*> (values.getRawPtr ());
4405  return this->transformGlobalValues (globalRow, numInputEnt, inputVals,
4406  inputCols, f, atomic);
4407  }
4408 
4415  void
4416  insertNonownedGlobalValues (const GlobalOrdinal globalRow,
4417  const Teuchos::ArrayView<const GlobalOrdinal>& indices,
4418  const Teuchos::ArrayView<const Scalar>& values);
4419 
4462  void
4463  insertIndicesAndValues (crs_graph_type& graph,
4464  RowInfo& rowInfo,
4465  const typename crs_graph_type::SLocalGlobalViews& newInds,
4466  const Teuchos::ArrayView<impl_scalar_type>& oldRowVals,
4467  const Teuchos::ArrayView<const impl_scalar_type>& newRowVals,
4468  const ELocalGlobal lg,
4469  const ELocalGlobal I);
4470 
4472  typedef DistObject<char, LocalOrdinal, GlobalOrdinal, Node> dist_object_type;
4473 
4474  protected:
4475  // useful typedefs
4476  typedef Teuchos::OrdinalTraits<LocalOrdinal> OTL;
4477  typedef Kokkos::Details::ArithTraits<impl_scalar_type> STS;
4478  typedef Kokkos::Details::ArithTraits<mag_type> STM;
4479  typedef MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> MV;
4480  typedef Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> V;
4481  typedef crs_graph_type Graph;
4482 
4483  // Enums
4484  enum GraphAllocationStatus {
4485  GraphAlreadyAllocated,
4486  GraphNotYetAllocated
4487  };
4488 
4489  private:
4493  Teuchos::ArrayRCP<Teuchos::Array<impl_scalar_type> >
4494  allocateValues2D ();
4495 
4496  protected:
4513  void allocateValues (ELocalGlobal lg, GraphAllocationStatus gas);
4514 
4524  size_t
4526  const RowInfo& rowInfo);
4527 
4542  void
4543  sortAndMergeIndicesAndValues (const bool sorted,
4544  const bool merged);
4545 
4553  void clearGlobalConstants();
4554 
4563  void computeGlobalConstants();
4564 
4565  public:
4567  bool haveGlobalConstants() const;
4568  protected:
4581  mutable Teuchos::RCP<MV> importMV_;
4582 
4595  mutable Teuchos::RCP<MV> exportMV_;
4596 
4616  Teuchos::RCP<MV>
4617  getColumnMapMultiVector (const MV& X_domainMap,
4618  const bool force = false) const;
4619 
4641  Teuchos::RCP<MV>
4642  getRowMapMultiVector (const MV& Y_rangeMap,
4643  const bool force = false) const;
4644 
4646  void
4647  applyNonTranspose (const MV& X_in,
4648  MV& Y_in,
4649  Scalar alpha,
4650  Scalar beta) const;
4651 
4653  void
4654  applyTranspose (const MV& X_in,
4655  MV& Y_in,
4656  const Teuchos::ETransp mode,
4657  Scalar alpha,
4658  Scalar beta) const;
4659 
4660  // matrix data accessors
4661 
4680  LocalOrdinal
4681  getViewRawConst (const impl_scalar_type*& vals,
4682  LocalOrdinal& numEnt,
4683  const RowInfo& rowinfo) const;
4684 
4703  LocalOrdinal
4704  getViewRaw (impl_scalar_type*& vals,
4705  LocalOrdinal& numEnt,
4706  const RowInfo& rowinfo) const;
4707 
4715  Teuchos::ArrayView<const impl_scalar_type> getView (RowInfo rowinfo) const;
4716 
4728  Teuchos::ArrayView<impl_scalar_type>
4729  getViewNonConst (const RowInfo& rowinfo) const;
4730 
4731  private:
4739  Kokkos::View<const impl_scalar_type*, execution_space, Kokkos::MemoryUnmanaged>
4740  getRowView (const RowInfo& rowInfo) const;
4741 
4753  Kokkos::View<impl_scalar_type*, execution_space, Kokkos::MemoryUnmanaged>
4754  getRowViewNonConst (const RowInfo& rowInfo) const;
4755 
4756  protected:
4757 
4763  void fillLocalMatrix (const Teuchos::RCP<Teuchos::ParameterList>& params);
4764 
4770  void fillLocalGraphAndMatrix (const Teuchos::RCP<Teuchos::ParameterList>& params);
4771 
4773  void checkInternalState () const;
4774 
4786 
4787  Teuchos::RCP<const Graph> staticGraph_;
4788  Teuchos::RCP< Graph> myGraph_;
4790 
4793 
4806 
4807  typename local_matrix_type::values_type k_values1D_;
4808  Teuchos::ArrayRCP<Teuchos::Array<impl_scalar_type> > values2D_;
4810 
4821 
4824 
4852  std::map<GlobalOrdinal, std::pair<Teuchos::Array<GlobalOrdinal>,
4853  Teuchos::Array<Scalar> > > nonlocals_;
4854 
4861 
4862  public:
4863  // FIXME (mfh 24 Feb 2014) Is it _really_ necessary to make this a
4864  // public inner class of CrsMatrix? It looks like it doesn't
4865  // depend on any implementation details of CrsMatrix at all. It
4866  // should really be declared and defined outside of CrsMatrix.
4867  template<class ViewType, class OffsetViewType>
4868  struct pack_functor {
4869  typedef typename ViewType::execution_space execution_space;
4870  ViewType src_;
4871  ViewType dst_;
4872  OffsetViewType src_offset_;
4873  OffsetViewType dst_offset_;
4874  typedef typename OffsetViewType::non_const_value_type scalar_index_type;
4875 
4876  pack_functor (ViewType dst, ViewType src,
4877  OffsetViewType dst_offset, OffsetViewType src_offset) :
4878  src_ (src),
4879  dst_ (dst),
4880  src_offset_ (src_offset),
4881  dst_offset_ (dst_offset)
4882  {}
4883 
4884  KOKKOS_INLINE_FUNCTION
4885  void operator () (const LocalOrdinal row) const {
4886  scalar_index_type srcPos = src_offset_(row);
4887  const scalar_index_type dstEnd = dst_offset_(row+1);
4888  scalar_index_type dstPos = dst_offset_(row);
4889  for ( ; dstPos < dstEnd; ++dstPos, ++srcPos) {
4890  dst_(dstPos) = src_(srcPos);
4891  }
4892  }
4893  };
4894  }; // class CrsMatrix
4895 } // namespace Classes
4896 
4905  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
4906  Teuchos::RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
4908  size_t maxNumEntriesPerRow = 0,
4909  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null)
4910  {
4912  return Teuchos::rcp (new matrix_type (map, maxNumEntriesPerRow,
4913  DynamicProfile, params));
4914  }
4915 
4916  template<class CrsMatrixType>
4917  Teuchos::RCP<CrsMatrixType>
4918  importAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
4919  const Import<typename CrsMatrixType::local_ordinal_type,
4920  typename CrsMatrixType::global_ordinal_type,
4921  typename CrsMatrixType::node_type>& importer,
4922  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4923  typename CrsMatrixType::global_ordinal_type,
4924  typename CrsMatrixType::node_type> >& domainMap,
4925  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4926  typename CrsMatrixType::global_ordinal_type,
4927  typename CrsMatrixType::node_type> >& rangeMap,
4928  const Teuchos::RCP<Teuchos::ParameterList>& params)
4929  {
4930  Teuchos::RCP<CrsMatrixType> destMatrix;
4931  sourceMatrix->importAndFillComplete (destMatrix, importer, domainMap, rangeMap, params);
4932  return destMatrix;
4933  }
4934 
4935  template<class CrsMatrixType>
4936  Teuchos::RCP<CrsMatrixType>
4937  importAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
4938  const Import<typename CrsMatrixType::local_ordinal_type,
4939  typename CrsMatrixType::global_ordinal_type,
4940  typename CrsMatrixType::node_type>& rowImporter,
4941  const Import<typename CrsMatrixType::local_ordinal_type,
4942  typename CrsMatrixType::global_ordinal_type,
4943  typename CrsMatrixType::node_type>& domainImporter,
4944  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4945  typename CrsMatrixType::global_ordinal_type,
4946  typename CrsMatrixType::node_type> >& domainMap,
4947  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4948  typename CrsMatrixType::global_ordinal_type,
4949  typename CrsMatrixType::node_type> >& rangeMap,
4950  const Teuchos::RCP<Teuchos::ParameterList>& params)
4951  {
4952  Teuchos::RCP<CrsMatrixType> destMatrix;
4953  sourceMatrix->importAndFillComplete (destMatrix, rowImporter, domainImporter, domainMap, rangeMap, params);
4954  return destMatrix;
4955  }
4956 
4957  template<class CrsMatrixType>
4958  Teuchos::RCP<CrsMatrixType>
4959  exportAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
4960  const Export<typename CrsMatrixType::local_ordinal_type,
4961  typename CrsMatrixType::global_ordinal_type,
4962  typename CrsMatrixType::node_type>& exporter,
4963  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4964  typename CrsMatrixType::global_ordinal_type,
4965  typename CrsMatrixType::node_type> >& domainMap,
4966  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4967  typename CrsMatrixType::global_ordinal_type,
4968  typename CrsMatrixType::node_type> >& rangeMap,
4969  const Teuchos::RCP<Teuchos::ParameterList>& params)
4970  {
4971  Teuchos::RCP<CrsMatrixType> destMatrix;
4972  sourceMatrix->exportAndFillComplete (destMatrix, exporter, domainMap, rangeMap, params);
4973  return destMatrix;
4974  }
4975 
4976  template<class CrsMatrixType>
4977  Teuchos::RCP<CrsMatrixType>
4978  exportAndFillCompleteCrsMatrix (const Teuchos::RCP<const CrsMatrixType>& sourceMatrix,
4979  const Export<typename CrsMatrixType::local_ordinal_type,
4980  typename CrsMatrixType::global_ordinal_type,
4981  typename CrsMatrixType::node_type>& rowExporter,
4982  const Export<typename CrsMatrixType::local_ordinal_type,
4983  typename CrsMatrixType::global_ordinal_type,
4984  typename CrsMatrixType::node_type>& domainExporter,
4985  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4986  typename CrsMatrixType::global_ordinal_type,
4987  typename CrsMatrixType::node_type> >& domainMap,
4988  const Teuchos::RCP<const Map<typename CrsMatrixType::local_ordinal_type,
4989  typename CrsMatrixType::global_ordinal_type,
4990  typename CrsMatrixType::node_type> >& rangeMap,
4991  const Teuchos::RCP<Teuchos::ParameterList>& params)
4992  {
4993  Teuchos::RCP<CrsMatrixType> destMatrix;
4994  sourceMatrix->exportAndFillComplete (destMatrix, rowExporter, domainExporter, domainMap, rangeMap, params);
4995  return destMatrix;
4996  }
4997 } // namespace Tpetra
4998 
5006 #endif // TPETRA_CRSMATRIX_DECL_HPP
Tpetra::Classes::MultiVector::isConstantStride
bool isConstantStride() const
Whether this multivector has constant stride between columns.
Definition: Tpetra_MultiVector_def.hpp:741
Tpetra::Classes::CrsMatrix::getNodeNumCols
size_t getNodeNumCols() const override
The number of columns connected to the locally owned rows of this matrix.
Definition: Tpetra_CrsMatrix_def.hpp:725
Tpetra::Classes::CrsMatrix::setAllToScalar
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to alpha.
Definition: Tpetra_CrsMatrix_def.hpp:3711
Tpetra::Classes::CrsMatrix::add
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.
Definition: Tpetra_CrsMatrix_def.hpp:7888
Tpetra::Classes::CrsMatrix::map_type
Map< LocalOrdinal, GlobalOrdinal, Node > map_type
The Map specialization suitable for this CrsMatrix specialization.
Definition: Tpetra_CrsMatrix_decl.hpp:466
Tpetra::Classes::CrsMatrix::CrsMatrix
friend class CrsMatrix
Alias for Tpetra::Classes::CrsMatrix.
Definition: Tpetra_CrsMatrix_decl.hpp:748
Tpetra::Classes::CrsMatrix::rightScale
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x) override
Scale the matrix on the right with the given Vector.
Definition: Tpetra_CrsMatrix_def.hpp:4109
Tpetra::Classes::CrsMatrix::crs_graph_type
CrsGraph< LocalOrdinal, GlobalOrdinal, Node > crs_graph_type
The CrsGraph specialization suitable for this CrsMatrix specialization.
Definition: Tpetra_CrsMatrix_decl.hpp:475
Tpetra::ProfileType
ProfileType
Definition: Tpetra_ConfigDefs.hpp:130
Tpetra::ESweepDirection
ESweepDirection
Sweep direction for Gauss-Seidel or Successive Over-Relaxation (SOR).
Definition: Tpetra_ConfigDefs.hpp:245
Tpetra::Classes::CrsMatrix::clone
Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node2 > > clone(const Teuchos::RCP< Node2 > &node2, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) const
Create a deep copy of this CrsMatrix, where the copy may have a different Node type.
Definition: Tpetra_CrsMatrix_decl.hpp:776
Tpetra::StaticProfile
Definition: Tpetra_ConfigDefs.hpp:131
Tpetra::Classes::CrsMatrix::getLocalValuesView
local_matrix_type::values_type getLocalValuesView() const
Get the Kokkos local values.
Definition: Tpetra_CrsMatrix_decl.hpp:4099
Tpetra::Classes::CrsMatrix::getNodeMaxNumRowEntries
size_t getNodeMaxNumRowEntries() const override
Maximum number of entries in any row of the matrix, on this process.
Definition: Tpetra_CrsMatrix_def.hpp:785
Tpetra::Classes::CrsMatrix::applyTranspose
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.
Definition: Tpetra_CrsMatrix_def.hpp:5153
Tpetra::Details::HasDeprecatedMethods2630_WarningThisClassIsNotForUsers
Mix-in to avoid spurious deprecation warnings due to #2630.
Definition: Tpetra_CrsGraph_decl.hpp:180
Tpetra::Classes::CrsMatrix::getRowMap
Teuchos::RCP< const map_type > getRowMap() const override
The Map that describes the row distribution in this matrix.
Definition: Tpetra_CrsMatrix_def.hpp:799
Tpetra::Classes::CrsMatrix::getNodeNumDiags
size_t TPETRA_DEPRECATED getNodeNumDiags() const override
Number of diagonal entries in the matrix's graph, on the calling process.
Definition: Tpetra_CrsMatrix_def.hpp:757
Tpetra::Classes::CrsMatrix::sumIntoLocalValues
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.
Definition: Tpetra_CrsMatrix_decl.hpp:1600
Tpetra::Classes::MultiVector::getLocalLength
size_t getLocalLength() const
Local number of rows on the calling process.
Definition: Tpetra_MultiVector_def.hpp:748
Tpetra::Classes::CrsMatrix::expertStaticFillComplete
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.
Definition: Tpetra_CrsMatrix_def.hpp:4802
Tpetra::Classes::CrsMatrix::isFillActive
bool isFillActive() const
Whether the matrix is not fill complete.
Definition: Tpetra_CrsMatrix_def.hpp:655
Tpetra::Classes::CrsMatrix::getGlobalRowView
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.
Definition: Tpetra_CrsMatrix_def.hpp:3605
Tpetra::Classes::CrsMatrix::gaussSeidelCopy
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.
Definition: Tpetra_CrsMatrix_def.hpp:5620
Tpetra::Classes::CrsMatrix::leftScale
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x) override
Scale the matrix on the left with the given Vector.
Definition: Tpetra_CrsMatrix_def.hpp:4032
Tpetra::Classes::CrsMatrix::scale
void scale(const Scalar &alpha)
Scale the matrix's values: this := alpha*this.
Definition: Tpetra_CrsMatrix_def.hpp:3667
Tpetra::RowInfo
Allocation information for a locally owned row in a CrsGraph or CrsMatrix.
Definition: Tpetra_CrsGraph_decl.hpp:112
Tpetra::Classes::CrsMatrix::getColMap
Teuchos::RCP< const map_type > getColMap() const override
The Map that describes the column distribution in this matrix.
Definition: Tpetra_CrsMatrix_def.hpp:806
Tpetra::Classes::CrsMatrix::insertLocalValues
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.
Definition: Tpetra_CrsMatrix_def.hpp:1799
Tpetra::Classes::CrsMatrix::getViewRaw
LocalOrdinal getViewRaw(impl_scalar_type *&vals, LocalOrdinal &numEnt, const RowInfo &rowinfo) const
Nonconst pointer to all entries (including extra space) in the given row.
Definition: Tpetra_CrsMatrix_def.hpp:3175
Tpetra::Details::DefaultTypes::node_type
::Kokkos::Compat::KokkosDeviceWrapperNode< execution_space > node_type
Default value of Node template parameter.
Definition: Tpetra_Details_DefaultTypes.hpp:105
Tpetra::Classes::CrsMatrix::getGlobalNumRows
global_size_t getGlobalNumRows() const override
Number of global elements in the row map of this matrix.
Definition: Tpetra_CrsMatrix_def.hpp:704
Tpetra::Classes::CrsMatrix::isUpperTriangular
bool TPETRA_DEPRECATED isUpperTriangular() const override
Whether the matrix is locally upper triangular.
Definition: Tpetra_CrsMatrix_def.hpp:891
Tpetra::Classes::CrsMatrix::applyNonTranspose
void applyNonTranspose(const MV &X_in, MV &Y_in, Scalar alpha, Scalar beta) const
Special case of apply() for mode == Teuchos::NO_TRANS.
Definition: Tpetra_CrsMatrix_def.hpp:4985
Tpetra::Classes::CrsMatrix::localMultiply
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.
Definition: Tpetra_CrsMatrix_decl.hpp:3115
Tpetra::Classes::CrsMatrix::reorderedGaussSeidelCopy
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.
Definition: Tpetra_CrsMatrix_def.hpp:5635
Tpetra::Classes::CrsMatrix::getCrsGraph
Teuchos::RCP< const crs_graph_type > getCrsGraph() const
This matrix's graph, as a CrsGraph.
Definition: Tpetra_CrsMatrix_def.hpp:837
Tpetra::Classes::CrsMatrix::replaceColMap
void replaceColMap(const Teuchos::RCP< const map_type > &newColMap)
Replace the matrix's column Map with the given Map.
Definition: Tpetra_CrsMatrix_def.hpp:4252
Tpetra::Classes::CrsMatrix::isFillComplete
bool isFillComplete() const override
Whether the matrix is fill complete.
Definition: Tpetra_CrsMatrix_def.hpp:648
Tpetra::Classes::CrsMatrix::getNodeNumEntries
size_t getNodeNumEntries() const override
The local number of entries in this matrix.
Definition: Tpetra_CrsMatrix_def.hpp:697
Tpetra::Classes::CrsMatrix::storageStatus_
::Tpetra::Details::EStorageStatus storageStatus_
Status of the matrix's storage, when not in a fill-complete state.
Definition: Tpetra_CrsMatrix_decl.hpp:4820
Tpetra::Classes::CrsMatrix::fillComplete_
bool fillComplete_
Whether the matrix is fill complete.
Definition: Tpetra_CrsMatrix_decl.hpp:4823
Tpetra::Classes::CrsMatrix::lclMatrix_
local_matrix_type lclMatrix_
The local sparse matrix.
Definition: Tpetra_CrsMatrix_decl.hpp:4792
Tpetra::Classes::CrsMatrix::hasColMap
bool hasColMap() const override
Whether the matrix has a well-defined column Map.
Definition: Tpetra_CrsMatrix_def.hpp:683
Tpetra::importAndFillCompleteCrsMatrix
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().
Definition: Tpetra_CrsMatrix_decl.hpp:4918
Tpetra_Details_PackTraits.hpp
Declaration and generic definition of traits class that tells Tpetra::CrsMatrix how to pack and unpac...
Tpetra::Classes::CrsMatrix::execution_space
device_type::execution_space execution_space
The Kokkos execution space.
Definition: Tpetra_CrsMatrix_decl.hpp:456
Tpetra::Classes::CrsMatrix::reindexColumns
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...
Definition: Tpetra_CrsMatrix_def.hpp:4270
Tpetra::Classes::CrsMatrix::packNew
void packNew(const Kokkos::DualView< const local_ordinal_type *, 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's data for an Import or Export.
Definition: Tpetra_CrsMatrix_def.hpp:7080
Tpetra::Classes::CrsMatrix::exportAndFillComplete
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.
Definition: Tpetra_CrsMatrix_def.hpp:8859
Tpetra::Classes::CrsMatrix::getViewNonConst
Teuchos::ArrayView< impl_scalar_type > getViewNonConst(const RowInfo &rowinfo) const
Nonconst view of all entries (including extra space) in the given row.
Definition: Tpetra_CrsMatrix_def.hpp:3276
Tpetra::Classes::CrsMatrix::description
std::string description() const override
A one-line description of this object.
Definition: Tpetra_CrsMatrix_def.hpp:6037
Tpetra::Classes::CrsMatrix::getGlobalNumCols
global_size_t getGlobalNumCols() const override
The number of global columns in the matrix.
Definition: Tpetra_CrsMatrix_def.hpp:711
Tpetra::Classes::CrsMatrix::getNodeNumDiagsImpl
size_t getNodeNumDiagsImpl() const override
DO NOT CALL THIS METHOD; THIS IS NOT FOR USERS.
Definition: Tpetra_CrsMatrix_def.hpp:748
Tpetra::Classes::CrsMatrix::isLocallyIndexed
bool isLocallyIndexed() const override
Whether the matrix is locally indexed on the calling process.
Definition: Tpetra_CrsMatrix_def.hpp:669
Tpetra::Classes::CrsMatrix::computeGlobalConstants
void computeGlobalConstants()
Compute matrix properties that require collectives.
Definition: Tpetra_CrsMatrix_def.hpp:4555
Tpetra::Classes::CrsMatrix::mag_type
Kokkos::Details::ArithTraits< impl_scalar_type >::mag_type mag_type
Type of a norm result.
Definition: Tpetra_CrsMatrix_decl.hpp:463
Tpetra::Classes::CrsMatrix::isStaticGraph
bool isStaticGraph() const
Indicates that the graph is static, so that new entries cannot be added to this matrix.
Definition: Tpetra_CrsMatrix_def.hpp:898
Tpetra::Classes::CrsMatrix::mergeRowIndicesAndValues
size_t mergeRowIndicesAndValues(crs_graph_type &graph, const RowInfo &rowInfo)
Merge duplicate row indices in the given row, along with their corresponding values.
Definition: Tpetra_CrsMatrix_def.hpp:4870
Tpetra::Classes::CrsMatrix::getNodeNumRows
size_t getNodeNumRows() const override
The number of matrix rows owned by the calling process.
Definition: Tpetra_CrsMatrix_def.hpp:718
Tpetra::Classes::CrsMatrix::replaceGlobalValues
LocalOrdinal replaceGlobalValues(const GlobalOrdinal globalRow, const typename UnmanagedView< GlobalIndicesViewType >::type &inputInds, const typename UnmanagedView< ImplScalarViewType >::type &inputVals) const
Replace one or more entries' values, using global indices.
Definition: Tpetra_CrsMatrix_decl.hpp:1217
Tpetra::Classes::CrsMatrix::getView
Teuchos::ArrayView< const impl_scalar_type > getView(RowInfo rowinfo) const
Constant view of all entries (including extra space) in the given row.
Definition: Tpetra_CrsMatrix_def.hpp:3094
Tpetra::Classes::CrsMatrix::getNumEntriesInGlobalRow
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const override
Number of entries in the sparse matrix in the given globa row, on the calling (MPI) process.
Definition: Tpetra_CrsMatrix_def.hpp:764
Tpetra::Classes::CrsMatrix::getRangeMap
Teuchos::RCP< const map_type > getRangeMap() const override
The range Map of this matrix.
Definition: Tpetra_CrsMatrix_def.hpp:820
Tpetra::Classes::CrsMatrix::useNewInterface
virtual bool useNewInterface() override
Whether the subclass implements the "old" or "new" (Kokkos-friendly) interface.
Definition: Tpetra_CrsMatrix_def.hpp:6316
Tpetra::Classes::CrsMatrix::isLowerTriangular
bool TPETRA_DEPRECATED isLowerTriangular() const override
Whether the matrix is locally lower triangular.
Definition: Tpetra_CrsMatrix_def.hpp:875
Tpetra::Classes::CrsMatrix::resumeFill
void resumeFill(const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Resume operations that may change the values or structure of the matrix.
Definition: Tpetra_CrsMatrix_def.hpp:4543
Tpetra::Classes::CrsMatrix::getFrobeniusNorm
mag_type getFrobeniusNorm() const override
Compute and return the Frobenius norm of the matrix.
Definition: Tpetra_CrsMatrix_def.hpp:4186
Tpetra::Classes::CrsMatrix::getRowMapMultiVector
Teuchos::RCP< MV > getRowMapMultiVector(const MV &Y_rangeMap, const bool force=false) const
Create a (or fetch a cached) row Map MultiVector.
Definition: Tpetra_CrsMatrix_def.hpp:7819
Tpetra::Classes::CrsMatrix::getGlobalNumDiags
global_size_t TPETRA_DEPRECATED getGlobalNumDiags() const override
Number of diagonal entries in the matrix's graph, over all processes in the matrix's communicator.
Definition: Tpetra_CrsMatrix_def.hpp:741
Tpetra::Classes::CrsMatrix::TPETRA_DEPRECATED
local_matrix_type::row_map_type::non_const_type t_RowPtrsNC TPETRA_DEPRECATED
DEPRECATED; use local_matrix_type::row_map_type::non_const_type instead.
Definition: Tpetra_CrsMatrix_decl.hpp:488
Tpetra::Classes::CrsMatrix::clearGlobalConstants
void clearGlobalConstants()
Clear matrix properties that require collectives.
Definition: Tpetra_CrsMatrix_def.hpp:4576
Tpetra::Classes::CrsMatrix::getDomainMap
Teuchos::RCP< const map_type > getDomainMap() const override
The domain Map of this matrix.
Definition: Tpetra_CrsMatrix_def.hpp:813
Tpetra::DynamicProfile
Definition: Tpetra_ConfigDefs.hpp:132
Tpetra::Classes::CrsMatrix::removeEmptyProcessesInPlace
virtual void removeEmptyProcessesInPlace(const Teuchos::RCP< const map_type > &newMap) override
Remove processes owning zero rows from the Maps and their communicator.
Definition: Tpetra_CrsMatrix_def.hpp:7865
Tpetra::Classes::CrsMatrix::unpackAndCombineNew
void unpackAndCombineNew(const Kokkos::DualView< const local_ordinal_type *, device_type > &importLIDs, const Kokkos::DualView< const char *, buffer_device_type > &imports, const Kokkos::DualView< const 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; implements "new" DistObject i...
Definition: Tpetra_CrsMatrix_def.hpp:7393
Tpetra::Classes::CrsMatrix::checkInternalState
void checkInternalState() const
Check that this object's state is sane; throw if it's not.
Definition: Tpetra_CrsMatrix_def.hpp:5980
Tpetra::Classes::CrsMatrix::TPETRA_DEPRECATED
local_matrix_type::row_map_type t_RowPtrs TPETRA_DEPRECATED
DEPRECATED; use local_matrix_type::row_map_type instead.
Definition: Tpetra_CrsMatrix_decl.hpp:486
Tpetra::Classes::CrsMatrix
Sparse matrix that presents a row-oriented interface that lets users read or modify entries.
Definition: Tpetra_CrsMatrix_decl.hpp:424
Tpetra::Classes::DistObject
Base class for distributed Tpetra objects that support data redistribution.
Definition: Tpetra_DistObject_decl.hpp:349
Tpetra::Classes::CrsMatrix::haveGlobalConstants
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
Definition: Tpetra_CrsMatrix_def.hpp:4569
Tpetra::Classes::CrsMatrix::replaceDomainMapAndImporter
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.
Definition: Tpetra_CrsMatrix_def.hpp:4303
Tpetra::Classes::CrsMatrix::setAllValues
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.
Definition: Tpetra_CrsMatrix_def.hpp:3748
Tpetra::Classes::CrsMatrix::transformLocalValues
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...
Definition: Tpetra_CrsMatrix_decl.hpp:1905
Tpetra::Classes::CrsMatrix::local_graph_type
crs_graph_type::local_graph_type local_graph_type
The part of the sparse matrix's graph on each MPI process.
Definition: Tpetra_CrsMatrix_decl.hpp:478
Tpetra::Classes::CrsMatrix::~CrsMatrix
virtual ~CrsMatrix()
Destructor.
Definition: Tpetra_CrsMatrix_def.hpp:621
Tpetra_CrsMatrix_fwd.hpp
Forward declaration of Tpetra::CrsMatrix.
Tpetra::Classes::CrsMatrix::getGlobalRowCopy
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,...
Definition: Tpetra_CrsMatrix_def.hpp:3395
Tpetra::Classes::CrsMatrix::fillLocalMatrix
void fillLocalMatrix(const Teuchos::RCP< Teuchos::ParameterList > &params)
Fill data into the local matrix.
Definition: Tpetra_CrsMatrix_def.hpp:1543
Tpetra::Classes::CrsMatrix::gaussSeidel
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
"Hybrid" Jacobi + (Gauss-Seidel or SOR) on .
Definition: Tpetra_CrsMatrix_def.hpp:5338
Tpetra::Classes::CrsMatrix::apply
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.
Definition: Tpetra_CrsMatrix_def.hpp:5301
Tpetra::Classes::CrsMatrix::getLocalMatrix
local_matrix_type getLocalMatrix() const
The local sparse matrix.
Definition: Tpetra_CrsMatrix_decl.hpp:2435
Tpetra::Classes::CrsMatrix::hasTransposeApply
bool hasTransposeApply() const override
Whether apply() allows applying the transpose or conjugate transpose.
Definition: Tpetra_CrsMatrix_def.hpp:905
Tpetra::Classes::CrsMatrix::global_ordinal_type
GlobalOrdinal global_ordinal_type
This class' third template parameter; the type of global indices.
Definition: Tpetra_CrsMatrix_decl.hpp:449
Tpetra::Classes::CrsMatrix::isLowerTriangularImpl
bool isLowerTriangularImpl() const override
DO NOT CALL THIS METHOD; THIS IS NOT FOR USERS.
Definition: Tpetra_CrsMatrix_def.hpp:866
Tpetra::Classes::CrsMatrix::getComm
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const override
The communicator over which the matrix is distributed.
Definition: Tpetra_CrsMatrix_def.hpp:627
Tpetra::Classes::CrsMatrix::fillLocalGraphAndMatrix
void fillLocalGraphAndMatrix(const Teuchos::RCP< Teuchos::ParameterList > &params)
Fill data into the local graph and matrix.
Definition: Tpetra_CrsMatrix_def.hpp:1102
Tpetra::Classes::CrsMatrix::getLocalDiagOffsets
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
Definition: Tpetra_CrsMatrix_def.hpp:3837
Tpetra::Classes::Import
Communication plan for data redistribution from a uniquely-owned to a (possibly) multiply-owned distr...
Definition: Tpetra_Import_decl.hpp:115
Tpetra::Classes::CrsMatrix::globalAssemble
void globalAssemble()
Communicate nonlocal contributions to other processes.
Definition: Tpetra_CrsMatrix_def.hpp:4345
Tpetra::Distributor
Sets up and executes a communication plan for a Tpetra DistObject.
Definition: Tpetra_Distributor.hpp:188
Tpetra::Classes::CrsMatrix::TPETRA_DEPRECATED
local_matrix_type k_local_matrix_type TPETRA_DEPRECATED
DEPRECATED; use local_matrix_type instead.
Definition: Tpetra_CrsMatrix_decl.hpp:495
Tpetra::Classes::CrsMatrix::getIndexBase
GlobalOrdinal getIndexBase() const override
The index base for global indices for this matrix.
Definition: Tpetra_CrsMatrix_def.hpp:792
Tpetra::exportAndFillCompleteCrsMatrix
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().
Definition: Tpetra_CrsMatrix_decl.hpp:4959
Tpetra::Classes::CrsMatrix::impl_scalar_type
Kokkos::Details::ArithTraits< Scalar >::val_type impl_scalar_type
The type used internally in place of Scalar.
Definition: Tpetra_CrsMatrix_decl.hpp:445
Tpetra::Classes::CrsMatrix::getLocalRowCopy
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,...
Definition: Tpetra_CrsMatrix_def.hpp:3284
Tpetra::Classes::CrsMatrix::reorderedGaussSeidel
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 "Hybrid" Jacobi + (Gauss-Seidel or SOR) on .
Definition: Tpetra_CrsMatrix_def.hpp:5351
Tpetra::createCrsMatrix
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.
Definition: Tpetra_CrsMatrix_decl.hpp:4907
Tpetra::Classes::CrsMatrix< SC, LO, GO, NT >::buffer_device_type
DistObject< SC, LO, GO, NT >::buffer_device_type buffer_device_type
Kokkos::Device specialization for communication buffers.
Definition: Tpetra_CrsMatrix_decl.hpp:3769
Tpetra::Classes::CrsMatrix::export_type
Export< LocalOrdinal, GlobalOrdinal, Node > export_type
The Export specialization suitable for this CrsMatrix specialization.
Definition: Tpetra_CrsMatrix_decl.hpp:472
Tpetra::Classes::CrsMatrix::importAndFillComplete
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.
Definition: Tpetra_CrsMatrix_def.hpp:8834
Tpetra::Classes::MultiVector
One or more distributed dense vectors.
Definition: Tpetra_MultiVector_decl.hpp:389
Tpetra::Classes::CrsMatrix::getNode
Teuchos::RCP< node_type > getNode() const override
The Kokkos Node instance.
Definition: Tpetra_CrsMatrix_def.hpp:634
Tpetra::Details::DefaultTypes::local_ordinal_type
int local_ordinal_type
Default value of Scalar template parameter.
Definition: Tpetra_Details_DefaultTypes.hpp:72
Tpetra::Classes::CrsMatrix::getGraph
Teuchos::RCP< const RowGraph< LocalOrdinal, GlobalOrdinal, Node > > getGraph() const override
This matrix's graph, as a RowGraph.
Definition: Tpetra_CrsMatrix_def.hpp:827
Tpetra::Classes::CrsMatrix::fillComplete
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...
Definition: Tpetra_CrsMatrix_def.hpp:4617
Tpetra::Classes::CrsMatrix::local_matrix_type
KokkosSparse::CrsMatrix< impl_scalar_type, LocalOrdinal, 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...
Definition: Tpetra_CrsMatrix_decl.hpp:483
Tpetra::Classes::CrsMatrix::isUpperTriangularImpl
bool isUpperTriangularImpl() const override
DO NOT CALL THIS METHOD; THIS IS NOT FOR USERS.
Definition: Tpetra_CrsMatrix_def.hpp:882
Tpetra::Classes::MultiVector::getNumVectors
size_t getNumVectors() const
Number of columns in the multivector.
Definition: Tpetra_MultiVector_def.hpp:1739
Tpetra::Classes::CrsMatrix::getLocalRowView
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,...
Definition: Tpetra_CrsMatrix_decl.hpp:2928
Tpetra::Classes::CrsMatrix::getNumEntriesInLocalRow
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const override
Number of entries in the sparse matrix in the given local row, on the calling (MPI) process.
Definition: Tpetra_CrsMatrix_def.hpp:771
Tpetra::Export
Classes::Export< LocalOrdinal, GlobalOrdinal, Node > Export
Alias for Tpetra::Classes::Export.
Definition: Tpetra_Export_fwd.hpp:71
Tpetra::Classes::CrsMatrix::checkSizes
virtual bool checkSizes(const SrcDistObject &source) override
Compare the source and target (this) objects for compatibility.
Definition: Tpetra_CrsMatrix_def.hpp:6298
Tpetra::Classes::CrsMatrix::isStorageOptimized
bool isStorageOptimized() const
Returns true if storage has been optimized.
Definition: Tpetra_CrsMatrix_def.hpp:662
Tpetra::Classes::CrsMatrix::isGloballyIndexed
bool isGloballyIndexed() const override
Whether the matrix is globally indexed on the calling process.
Definition: Tpetra_CrsMatrix_def.hpp:676
Tpetra::Classes::CrsMatrix::device_type
Node::device_type device_type
The Kokkos device type.
Definition: Tpetra_CrsMatrix_decl.hpp:454
Tpetra::Classes::Map
A parallel distribution of indices over processes.
Definition: Tpetra_Map_decl.hpp:247
Tpetra::Classes::CrsMatrix::nonlocals_
std::map< GlobalOrdinal, std::pair< Teuchos::Array< GlobalOrdinal >, Teuchos::Array< Scalar > > > nonlocals_
Nonlocal data added using insertGlobalValues().
Definition: Tpetra_CrsMatrix_decl.hpp:4853
Tpetra::CrsMatrix
Classes::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrix
Alias for Tpetra::Classes::CrsMatrix.
Definition: Tpetra_CrsMatrix_fwd.hpp:72
Tpetra::Map
Classes::Map< LocalOrdinal, GlobalOrdinal, Node > Map
Alias for Tpetra::Classes::Map.
Definition: Tpetra_Map_fwd.hpp:71
Tpetra::Classes::CrsMatrix::getLocalDiagCopy
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const override
Get a copy of the diagonal entries of the matrix.
Definition: Tpetra_CrsMatrix_def.hpp:3879
Tpetra::Classes::CrsMatrix::getProfileType
ProfileType getProfileType() const
Returns true if the matrix was allocated with static data structures.
Definition: Tpetra_CrsMatrix_def.hpp:641
Tpetra::Classes::CrsMatrix::insertGlobalValues
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.
Definition: Tpetra_CrsMatrix_def.hpp:2049
Tpetra::Classes::CrsMatrix::sumIntoGlobalValues
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.
Definition: Tpetra_CrsMatrix_def.hpp:2606
Tpetra::Classes::CrsMatrix::TPETRA_DEPRECATED
local_matrix_type::values_type t_ValuesType TPETRA_DEPRECATED
DEPRECATED; use local_matrix_type::values_type instead.
Definition: Tpetra_CrsMatrix_decl.hpp:492
Tpetra::Classes::CrsMatrix::getGlobalNumEntries
global_size_t getGlobalNumEntries() const override
The global number of entries in this matrix.
Definition: Tpetra_CrsMatrix_def.hpp:690
Tpetra::Classes::CrsMatrix::supportsRowViews
virtual bool supportsRowViews() const override
Return true if getLocalRowView() and getGlobalRowView() are valid for this object.
Definition: Tpetra_CrsMatrix_def.hpp:912
Tpetra::Classes::CrsGraph::local_graph_type
Kokkos::StaticCrsGraph< LocalOrdinal, Kokkos::LayoutLeft, execution_space > local_graph_type
The type of the part of the sparse graph on each MPI process.
Definition: Tpetra_CrsGraph_decl.hpp:292
Tpetra::Import
Classes::Import< LocalOrdinal, GlobalOrdinal, Node > Import
Alias for Tpetra::Classes::Import.
Definition: Tpetra_Import_fwd.hpp:71
Tpetra::Classes::CrsMatrix::describe
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.
Definition: Tpetra_CrsMatrix_def.hpp:6063
Tpetra::Classes::CrsMatrix::node_type
Node node_type
This class' fourth template parameter; the Kokkos device type.
Definition: Tpetra_CrsMatrix_decl.hpp:451
Tpetra::global_size_t
size_t global_size_t
Global size_t object.
Definition: Tpetra_ConfigDefs.hpp:109
Tpetra::Classes::CrsMatrix::getGlobalNumDiagsImpl
global_size_t getGlobalNumDiagsImpl() const override
DO NOT CALL THIS METHOD; THIS IS NOT FOR USERS.
Definition: Tpetra_CrsMatrix_def.hpp:732
Tpetra::Classes::CrsMatrix::TPETRA_DEPRECATED
local_graph_type::entries_type::non_const_type t_LocalOrdinal_1D TPETRA_DEPRECATED
DEPRECATED; use local_graph_type::entries_type::non_const_type instead.
Definition: Tpetra_CrsMatrix_decl.hpp:490
Tpetra
Namespace Tpetra contains the class and methods constituting the Tpetra library.
Tpetra::Classes::CrsMatrix::allocateValues
void allocateValues(ELocalGlobal lg, GraphAllocationStatus gas)
Allocate values (and optionally indices) using the Node.
Definition: Tpetra_CrsMatrix_def.hpp:954
Tpetra::Classes::CrsMatrix::frobNorm_
mag_type frobNorm_
Cached Frobenius norm of the (global) matrix.
Definition: Tpetra_CrsMatrix_decl.hpp:4860
Tpetra::Classes::Export
Communication plan for data redistribution from a (possibly) multiply-owned to a uniquely-owned distr...
Definition: Tpetra_Export_decl.hpp:124
Tpetra::SrcDistObject
Abstract base class for objects that can be the source of an Import or Export operation.
Definition: Tpetra_SrcDistObject.hpp:89
Tpetra::Classes::CrsMatrix::localGaussSeidel
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 .
Definition: Tpetra_CrsMatrix_decl.hpp:3267
Tpetra::Classes::CrsMatrix::importMV_
Teuchos::RCP< MV > importMV_
Column Map MultiVector used in apply() and gaussSeidel().
Definition: Tpetra_CrsMatrix_decl.hpp:4581
Tpetra::Classes::CrsMatrix::getColumnMapMultiVector
Teuchos::RCP< MV > getColumnMapMultiVector(const MV &X_domainMap, const bool force=false) const
Create a (or fetch a cached) column Map MultiVector.
Definition: Tpetra_CrsMatrix_def.hpp:7762
Tpetra::Classes::CrsMatrix::getViewRawConst
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.
Definition: Tpetra_CrsMatrix_def.hpp:3135
Tpetra::Classes::CrsMatrix::local_ordinal_type
LocalOrdinal local_ordinal_type
This class' second template parameter; the type of local indices.
Definition: Tpetra_CrsMatrix_decl.hpp:447
Tpetra::Classes::CrsMatrix::replaceLocalValues
LocalOrdinal replaceLocalValues(const LocalOrdinal localRow, const typename UnmanagedView< LocalIndicesViewType >::type &inputInds, const typename UnmanagedView< ImplScalarViewType >::type &inputVals) const
Replace one or more entries' values, using local row and column indices.
Definition: Tpetra_CrsMatrix_decl.hpp:1347
Tpetra::Classes::CrsMatrix::exportMV_
Teuchos::RCP< MV > exportMV_
Row Map MultiVector used in apply().
Definition: Tpetra_CrsMatrix_decl.hpp:4595
Tpetra::Classes::CrsMatrix::getLocalRowViewRaw
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,...
Definition: Tpetra_CrsMatrix_def.hpp:3590
Tpetra::Classes::CrsGraph
A distributed graph accessed by rows (adjacency lists) and stored sparsely.
Definition: Tpetra_CrsGraph_decl.hpp:259
Tpetra::Classes::CrsMatrix::scalar_type
Scalar scalar_type
This class' first template parameter; the type of each entry in the matrix.
Definition: Tpetra_CrsMatrix_decl.hpp:435
Tpetra::Classes::CrsMatrix::getGlobalMaxNumRowEntries
size_t getGlobalMaxNumRowEntries() const override
Maximum number of entries in any row of the matrix, over all processes in the matrix's communicator.
Definition: Tpetra_CrsMatrix_def.hpp:778
Tpetra::Classes::CrsMatrix::reorderedLocalGaussSeidel
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 .
Definition: Tpetra_CrsMatrix_decl.hpp:3370
Tpetra::Classes::CrsMatrix::convert
Teuchos::RCP< CrsMatrix< T, LocalOrdinal, GlobalOrdinal, Node > > convert() const
Return another CrsMatrix with the same entries, but converted to a different Scalar type T.
Definition: Tpetra_CrsMatrix_def.hpp:5945
Tpetra::Classes::CrsMatrix::transformGlobalValues
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...
Definition: Tpetra_CrsMatrix_decl.hpp:1995
Tpetra::Details::EStorageStatus
EStorageStatus
Status of the graph's or matrix's storage, when not in a fill-complete state.
Definition: Tpetra_CrsGraph_decl.hpp:164
Tpetra::Classes::CrsMatrix::sortAndMergeIndicesAndValues
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...
Definition: Tpetra_CrsMatrix_def.hpp:4929
Tpetra::Classes::CrsMatrix::import_type
Import< LocalOrdinal, GlobalOrdinal, Node > import_type
The Import specialization suitable for this CrsMatrix specialization.
Definition: Tpetra_CrsMatrix_decl.hpp:469
Tpetra::CombineMode
CombineMode
Rule for combining data in an Import or Export.
Definition: Tpetra_CombineMode.hpp:94
Tpetra::Classes::RowMatrix
A read-only, row-oriented interface to a sparse matrix.
Definition: Tpetra_RowMatrix_decl.hpp:85
Tpetra::Classes::CrsMatrix::getLocalRowView
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.
Definition: Tpetra_CrsMatrix_def.hpp:3499
Tpetra::Classes::Vector
A distributed dense vector.
Definition: Tpetra_Vector_decl.hpp:82