10#ifndef TPETRA_CRSMATRIX_DECL_HPP
11#define TPETRA_CRSMATRIX_DECL_HPP
18#include "KokkosSparse_Utils.hpp"
19#include "KokkosSparse_CrsMatrix.hpp"
20#include "Tpetra_Details_MatrixApplyHelper.hpp"
21#include "Tpetra_RowMatrix_decl.hpp"
22#include "Tpetra_Exceptions.hpp"
23#include "Tpetra_DistObject.hpp"
24#include "Tpetra_CrsGraph.hpp"
25#include "Tpetra_Vector.hpp"
27#include "Tpetra_Details_ExecutionSpacesUser.hpp"
28#include "Teuchos_DataAccess.hpp"
35template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
36class crsMatrix_Swap_Tester;
89template <
class CrsMatrixType>
90Teuchos::RCP<CrsMatrixType>
92 const Import<
typename CrsMatrixType::local_ordinal_type,
93 typename CrsMatrixType::global_ordinal_type,
94 typename CrsMatrixType::node_type>&
importer,
95 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
96 typename CrsMatrixType::global_ordinal_type,
97 typename CrsMatrixType::node_type> >& domainMap = Teuchos::null,
98 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
99 typename CrsMatrixType::global_ordinal_type,
100 typename CrsMatrixType::node_type> >&
rangeMap = Teuchos::null,
101 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
156template <
class CrsMatrixType>
157Teuchos::RCP<CrsMatrixType>
159 const Import<
typename CrsMatrixType::local_ordinal_type,
160 typename CrsMatrixType::global_ordinal_type,
162 const Import<
typename CrsMatrixType::local_ordinal_type,
163 typename CrsMatrixType::global_ordinal_type,
165 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
166 typename CrsMatrixType::global_ordinal_type,
167 typename CrsMatrixType::node_type> >& domainMap,
168 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
169 typename CrsMatrixType::global_ordinal_type,
170 typename CrsMatrixType::node_type> >&
rangeMap,
171 const Teuchos::RCP<Teuchos::ParameterList>&
params);
206template <
class CrsMatrixType>
207Teuchos::RCP<CrsMatrixType>
209 const Export<
typename CrsMatrixType::local_ordinal_type,
210 typename CrsMatrixType::global_ordinal_type,
211 typename CrsMatrixType::node_type>&
exporter,
212 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
213 typename CrsMatrixType::global_ordinal_type,
214 typename CrsMatrixType::node_type> >& domainMap = Teuchos::null,
215 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
216 typename CrsMatrixType::global_ordinal_type,
217 typename CrsMatrixType::node_type> >&
rangeMap = Teuchos::null,
218 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
253template <
class CrsMatrixType>
254Teuchos::RCP<CrsMatrixType>
256 const Export<
typename CrsMatrixType::local_ordinal_type,
257 typename CrsMatrixType::global_ordinal_type,
259 const Export<
typename CrsMatrixType::local_ordinal_type,
260 typename CrsMatrixType::global_ordinal_type,
262 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
263 typename CrsMatrixType::global_ordinal_type,
264 typename CrsMatrixType::node_type> >& domainMap,
265 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
266 typename CrsMatrixType::global_ordinal_type,
267 typename CrsMatrixType::node_type> >&
rangeMap,
268 const Teuchos::RCP<Teuchos::ParameterList>&
params);
273template <
class SC,
class LO,
class GO,
class NO>
398 public DistObject<char, LocalOrdinal, GlobalOrdinal, Node>,
456#if KOKKOS_VERSION >= 40799
457 using mag_type =
typename KokkosKernels::ArithTraits<impl_scalar_type>::mag_type;
459 using mag_type =
typename Kokkos::ArithTraits<impl_scalar_type>::mag_type;
476 typename local_graph_device_type::size_type>;
477#if KOKKOS_VERSION >= 40799
478 using local_matrix_host_type =
479 typename local_matrix_device_type::host_mirror_type;
481 using local_matrix_host_type =
482 typename local_matrix_device_type::HostMirror;
485 using row_ptrs_device_view_type =
486 typename row_matrix_type::row_ptrs_device_view_type;
487 using row_ptrs_host_view_type =
488 typename row_matrix_type::row_ptrs_host_view_type;
490 using local_inds_device_view_type =
491 typename row_matrix_type::local_inds_device_view_type;
492 using local_inds_host_view_type =
493 typename row_matrix_type::local_inds_host_view_type;
494 using nonconst_local_inds_host_view_type =
495 typename row_matrix_type::nonconst_local_inds_host_view_type;
497 using global_inds_device_view_type =
498 typename row_matrix_type::global_inds_device_view_type;
499 using global_inds_host_view_type =
500 typename row_matrix_type::global_inds_host_view_type;
501 using nonconst_global_inds_host_view_type =
502 typename row_matrix_type::nonconst_global_inds_host_view_type;
504 using values_device_view_type =
505 typename row_matrix_type::values_device_view_type;
506 using values_host_view_type =
507 typename row_matrix_type::values_host_view_type;
508 using nonconst_values_host_view_type =
509 typename row_matrix_type::nonconst_values_host_view_type;
546 CrsMatrix(
const Teuchos::RCP<const map_type>& rowMap,
548 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
562 CrsMatrix(
const Teuchos::RCP<const map_type>& rowMap,
564 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
584 CrsMatrix(
const Teuchos::RCP<const map_type>& rowMap,
585 const Teuchos::RCP<const map_type>& colMap,
587 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
607 CrsMatrix(
const Teuchos::RCP<const map_type>& rowMap,
608 const Teuchos::RCP<const map_type>& colMap,
610 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
641 const Teuchos::RCP<const crs_graph_type>&
graph,
642 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
668 explicit CrsMatrix(
const Teuchos::RCP<const crs_graph_type>&
graph,
669 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
699 explicit CrsMatrix(
const Teuchos::RCP<const crs_graph_type>&
graph,
700 const typename local_matrix_device_type::values_type&
values,
701 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
729 CrsMatrix(
const Teuchos::RCP<const map_type>& rowMap,
730 const Teuchos::RCP<const map_type>& colMap,
731 const typename local_graph_device_type::row_map_type&
rowPointers,
732 const typename local_graph_device_type::entries_type::non_const_type&
columnIndices,
733 const typename local_matrix_device_type::values_type&
values,
734 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
762 CrsMatrix(
const Teuchos::RCP<const map_type>& rowMap,
763 const Teuchos::RCP<const map_type>& colMap,
766 const Teuchos::ArrayRCP<Scalar>&
values,
767 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
790 CrsMatrix(
const Teuchos::RCP<const map_type>& rowMap,
791 const Teuchos::RCP<const map_type>& colMap,
793 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
826 const Teuchos::RCP<const map_type>& rowMap,
827 const Teuchos::RCP<const map_type>& colMap,
828 const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
829 const Teuchos::RCP<const map_type>&
rangeMap = Teuchos::null,
830 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
837 const Teuchos::RCP<const map_type>& rowMap,
838 const Teuchos::RCP<const map_type>& colMap,
839 const Teuchos::RCP<const map_type>& domainMap,
840 const Teuchos::RCP<const map_type>&
rangeMap,
841 const Teuchos::RCP<const import_type>&
importer,
842 const Teuchos::RCP<const export_type>&
exporter,
843 const Teuchos::RCP<Teuchos::ParameterList>&
params =
864 template <
class S2,
class LO2,
class GO2,
class N2>
868 template <
class S2,
class LO2,
class GO2,
class N2>
875 template <
class MatrixArray,
class MultiVectorArray>
877 const typename std::remove_pointer<typename MultiVectorArray::value_type>::type&
X,
879 typename std::remove_pointer<typename MatrixArray::value_type>::type::scalar_type
alpha,
880 typename std::remove_pointer<typename MatrixArray::value_type>::type::scalar_type
beta,
881 Teuchos::RCP<Teuchos::ParameterList>
params);
955 const Teuchos::ArrayView<const GlobalOrdinal>&
cols,
956 const Teuchos::ArrayView<const Scalar>&
vals);
1022 const Teuchos::ArrayView<const LocalOrdinal>&
cols,
1023 const Teuchos::ArrayView<const Scalar>&
vals,
1111 const Kokkos::View<const global_ordinal_type*, Kokkos::AnonymousSpace>&
inputInds,
1112 const Kokkos::View<const impl_scalar_type*, Kokkos::AnonymousSpace>&
inputVals);
1118 const Teuchos::ArrayView<const GlobalOrdinal>&
cols,
1119 const Teuchos::ArrayView<const Scalar>&
vals);
1199 const Kokkos::View<const local_ordinal_type*, Kokkos::AnonymousSpace>&
inputInds,
1200 const Kokkos::View<const impl_scalar_type*, Kokkos::AnonymousSpace>&
inputVals);
1207 const Teuchos::ArrayView<const LocalOrdinal>&
cols,
1208 const Teuchos::ArrayView<const Scalar>&
vals);
1238 static const bool useAtomicUpdatesByDefault =
1239#ifdef KOKKOS_ENABLE_SERIAL
1240 !std::is_same<execution_space, Kokkos::Serial>::value;
1276 const bool atomic = useAtomicUpdatesByDefault);
1317 const Teuchos::ArrayView<const GlobalOrdinal>&
cols,
1318 const Teuchos::ArrayView<const Scalar>&
vals,
1319 const bool atomic = useAtomicUpdatesByDefault);
1348 const bool atomic = useAtomicUpdatesByDefault);
1370 const bool atomic = useAtomicUpdatesByDefault);
1412 const Kokkos::View<const local_ordinal_type*, Kokkos::AnonymousSpace>&
inputInds,
1413 const Kokkos::View<const impl_scalar_type*, Kokkos::AnonymousSpace>&
inputVals,
1414 const bool atomic = useAtomicUpdatesByDefault);
1447 const Teuchos::ArrayView<const LocalOrdinal>&
cols,
1448 const Teuchos::ArrayView<const Scalar>&
vals,
1449 const bool atomic = useAtomicUpdatesByDefault);
1477 const bool atomic = useAtomicUpdatesByDefault);
1518 const bool atomic = useAtomicUpdatesByDefault);
1558 const bool atomic = useAtomicUpdatesByDefault);
1592 const bool atomic = useAtomicUpdatesByDefault);
1626 const bool atomic = useAtomicUpdatesByDefault);
1677 const typename UnmanagedView<LocalIndicesViewType>::type&
inputInds,
1678 const typename UnmanagedView<ImplScalarViewType>::type&
inputVals,
1680 const bool atomic = useAtomicUpdatesByDefault) {
1687 static_assert(Kokkos::is_view<LocalIndicesViewType>::value,
1688 "First template parameter LocalIndicesViewType must be "
1690 static_assert(Kokkos::is_view<ImplScalarViewType>::value,
1691 "Second template parameter ImplScalarViewType must be a "
1693 static_assert(
static_cast<int>(LocalIndicesViewType::rank) == 1,
1694 "First template parameter LocalIndicesViewType must "
1696 static_assert(
static_cast<int>(ImplScalarViewType::rank) == 1,
1697 "Second template parameter ImplScalarViewType must have "
1699 static_assert(std::is_same<
1700 typename LocalIndicesViewType::non_const_value_type,
1702 "First template parameter LocalIndicesViewType must "
1703 "contain values of type local_ordinal_type.");
1704 static_assert(std::is_same<
1705 typename ImplScalarViewType::non_const_value_type,
1707 "Second template parameter ImplScalarViewType must "
1708 "contain values of type impl_scalar_type.");
1712 return Teuchos::OrdinalTraits<LO>::invalid();
1714 return this->transformLocalValues(
lclRow,
1763 template <
class BinaryFunction,
class InputMemorySpace>
1773 const bool atomic = useAtomicUpdatesByDefault) {
1777 return Teuchos::OrdinalTraits<LO>::invalid();
1779 return this->transformGlobalValues(
gblRow,
1818 setAllValues(
const typename local_graph_device_type::row_map_type&
ptr,
1819 const typename local_graph_device_type::entries_type::non_const_type&
ind,
1820 const typename local_matrix_device_type::values_type&
val);
1871 const Teuchos::ArrayRCP<LocalOrdinal>&
ind,
1872 const Teuchos::ArrayRCP<Scalar>&
val);
1933 void resumeFill(
const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
1993 fillComplete(
const Teuchos::RCP<const map_type>& domainMap,
1994 const Teuchos::RCP<const map_type>&
rangeMap,
1995 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
2054 const Teuchos::RCP<const map_type>&
rangeMap,
2055 const Teuchos::RCP<const import_type>&
importer = Teuchos::null,
2056 const Teuchos::RCP<const export_type>&
exporter = Teuchos::null,
2057 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
2162 const Teuchos::RCP<const map_type>&
newColMap,
2163 const Teuchos::RCP<const import_type>&
newImport = Teuchos::null,
2245 Teuchos::RCP<const Teuchos::Comm<int> >
getComm()
const override;
2248 Teuchos::RCP<const map_type>
getRowMap()
const override;
2251 Teuchos::RCP<const map_type>
getColMap()
const override;
2254 Teuchos::RCP<const RowGraph<LocalOrdinal, GlobalOrdinal, Node> >
2258 Teuchos::RCP<const crs_graph_type>
getCrsGraph()
const;
2274#if __armclang_major__ == 22 && __armclang_minor__ == 1
2279#define TPETRA_DETAILS_ALWAYS_INLINE __attribute__((always_inline))
2281#define TPETRA_DETAILS_ALWAYS_INLINE
2297 local_matrix_host_type getLocalMatrixHost()
const;
2298#undef TPETRA_DETAILS_ALWAYS_INLINE
2528 using values_dualv_type =
2529 Kokkos::DualView<impl_scalar_type*, device_type>;
2530 using values_wdv_type =
2532 values_wdv_type valuesUnpacked_wdv;
2533 mutable values_wdv_type valuesPacked_wdv;
2586 nonconst_global_inds_host_view_type&
Indices,
2587 nonconst_values_host_view_type&
Values,
2606 nonconst_local_inds_host_view_type&
Indices,
2607 nonconst_values_host_view_type&
Values,
2625 global_inds_host_view_type& indices,
2626 values_host_view_type&
values)
const override;
2642 local_inds_host_view_type& indices,
2643 values_host_view_type&
values)
const override;
2731 Kokkos::MemoryUnmanaged>& offsets)
const;
2757 const Teuchos::ArrayView<const size_t>& offsets)
const;
2841 const Teuchos::ETransp
mode = Teuchos::NO_TRANS,
2842 const Scalar&
alpha = Teuchos::ScalarTraits<Scalar>::one(),
2843 const Scalar&
beta = Teuchos::ScalarTraits<Scalar>::zero())
const;
2848 Teuchos::RCP<CrsMatrix<T, LocalOrdinal, GlobalOrdinal, Node> >
2868 Teuchos::ETransp
mode = Teuchos::NO_TRANS,
2869 Scalar alpha = Teuchos::ScalarTraits<Scalar>::one(),
2870 Scalar beta = Teuchos::ScalarTraits<Scalar>::zero())
const override;
2882 Teuchos::RCP<const map_type>
getDomainMap()
const override;
2890 Teuchos::RCP<const map_type>
getRangeMap()
const override;
2906 virtual Teuchos::RCP<RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
2912 const Teuchos::RCP<Teuchos::ParameterList>&
params)
const override;
2925 const Teuchos::EVerbosityLevel
verbLevel =
2926 Teuchos::Describable::verbLevel_default)
const override;
2945 const bool verbose);
2949 copyAndPermuteStaticGraph(
2951 const size_t numSameIDs,
2954 const size_t numPermutes);
2957 copyAndPermuteNonStaticGraph(
2959 const size_t numSameIDs,
2960 const Kokkos::DualView<const local_ordinal_type*, buffer_device_type>&
permuteToLIDs_dv,
2961 const Kokkos::DualView<const local_ordinal_type*, buffer_device_type>&
permuteFromLIDs_dv,
2962 const size_t numPermutes);
2971 const size_t numSameIDs,
2972 const Kokkos::DualView<
2975 const Kokkos::DualView<
2982 const Kokkos::DualView<
2985 Kokkos::DualView<char*, buffer_device_type>& exports,
2998 unpackAndCombineImpl(
3001 Kokkos::DualView<char*, buffer_device_type> imports,
3005 const bool verbose);
3010 unpackAndCombineImplNonStatic(
3013 Kokkos::DualView<char*, buffer_device_type> imports,
3030 Kokkos::DualView<char*, buffer_device_type> imports,
3146 packNew(
const Kokkos::DualView<const local_ordinal_type*, buffer_device_type>&
exportLIDs,
3147 Kokkos::DualView<char*, buffer_device_type>& exports,
3159 packNonStaticNew(
const Kokkos::DualView<const local_ordinal_type*, buffer_device_type>&
exportLIDs,
3160 Kokkos::DualView<char*, buffer_device_type>& exports,
3194 packRow(
char exports[],
3259 const char imports[],
3274 allocatePackSpaceNew(Kokkos::DualView<char*, buffer_device_type>& exports,
3275 size_t& totalNumEntries,
3282 typename local_matrix_host_type::values_type::const_type
3284 return valuesPacked_wdv.getHostView(
s);
3288 typename local_matrix_host_type::values_type
3290 return valuesPacked_wdv.getHostView(
s);
3294 typename local_matrix_host_type::values_type
3296 return valuesPacked_wdv.getHostView(
s);
3300 typename local_matrix_device_type::values_type::const_type
3302 return valuesPacked_wdv.getDeviceView(
s);
3306 typename local_matrix_device_type::values_type
3308 return valuesPacked_wdv.getDeviceView(
s);
3312 typename local_matrix_device_type::values_type
3314 return valuesPacked_wdv.getDeviceView(
s);
3319 template <
class CrsMatrixType>
3320 friend Teuchos::RCP<CrsMatrixType>
3322 const Import<
typename CrsMatrixType::local_ordinal_type,
3323 typename CrsMatrixType::global_ordinal_type,
3324 typename CrsMatrixType::node_type>&
importer,
3325 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
3326 typename CrsMatrixType::global_ordinal_type,
3327 typename CrsMatrixType::node_type> >& domainMap,
3328 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
3329 typename CrsMatrixType::global_ordinal_type,
3330 typename CrsMatrixType::node_type> >&
rangeMap,
3331 const Teuchos::RCP<Teuchos::ParameterList>&
params);
3334 template <
class CrsMatrixType>
3335 friend Teuchos::RCP<CrsMatrixType>
3337 const Import<
typename CrsMatrixType::local_ordinal_type,
3338 typename CrsMatrixType::global_ordinal_type,
3340 const Import<
typename CrsMatrixType::local_ordinal_type,
3341 typename CrsMatrixType::global_ordinal_type,
3343 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
3344 typename CrsMatrixType::global_ordinal_type,
3345 typename CrsMatrixType::node_type> >& domainMap,
3346 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
3347 typename CrsMatrixType::global_ordinal_type,
3348 typename CrsMatrixType::node_type> >&
rangeMap,
3349 const Teuchos::RCP<Teuchos::ParameterList>&
params);
3352 template <
class CrsMatrixType>
3353 friend Teuchos::RCP<CrsMatrixType>
3355 const Export<
typename CrsMatrixType::local_ordinal_type,
3356 typename CrsMatrixType::global_ordinal_type,
3357 typename CrsMatrixType::node_type>&
exporter,
3358 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
3359 typename CrsMatrixType::global_ordinal_type,
3360 typename CrsMatrixType::node_type> >& domainMap,
3361 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
3362 typename CrsMatrixType::global_ordinal_type,
3363 typename CrsMatrixType::node_type> >&
rangeMap,
3364 const Teuchos::RCP<Teuchos::ParameterList>&
params);
3367 template <
class CrsMatrixType>
3368 friend Teuchos::RCP<CrsMatrixType>
3370 const Export<
typename CrsMatrixType::local_ordinal_type,
3371 typename CrsMatrixType::global_ordinal_type,
3373 const Export<
typename CrsMatrixType::local_ordinal_type,
3374 typename CrsMatrixType::global_ordinal_type,
3376 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
3377 typename CrsMatrixType::global_ordinal_type,
3378 typename CrsMatrixType::node_type> >& domainMap,
3379 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
3380 typename CrsMatrixType::global_ordinal_type,
3381 typename CrsMatrixType::node_type> >&
rangeMap,
3382 const Teuchos::RCP<Teuchos::ParameterList>&
params);
3403 const Teuchos::RCP<const map_type>& domainMap,
3404 const Teuchos::RCP<const map_type>&
rangeMap,
3405 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null)
const;
3426 const Teuchos::RCP<const map_type>& domainMap,
3427 const Teuchos::RCP<const map_type>&
rangeMap,
3428 const Teuchos::RCP<Teuchos::ParameterList>&
params)
const;
3448 const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
3449 const Teuchos::RCP<const map_type>&
rangeMap = Teuchos::null,
3450 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null)
const;
3471 const Teuchos::RCP<const map_type>& domainMap,
3472 const Teuchos::RCP<const map_type>&
rangeMap,
3473 const Teuchos::RCP<Teuchos::ParameterList>&
params)
const;
3498 const ::Tpetra::Details::Transfer<LocalOrdinal, GlobalOrdinal, Node>&
rowTransfer,
3499 const Teuchos::RCP<const ::Tpetra::Details::Transfer<LocalOrdinal, GlobalOrdinal, Node> >&
domainTransfer,
3500 const Teuchos::RCP<const map_type>& domainMap = Teuchos::null,
3501 const Teuchos::RCP<const map_type>&
rangeMap = Teuchos::null,
3502 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null)
const;
3534 insertGlobalValuesFiltered(
3536 const Teuchos::ArrayView<const GlobalOrdinal>& indices,
3537 const Teuchos::ArrayView<const Scalar>&
values,
3543 insertGlobalValuesFilteredChecked(
3545 const Teuchos::ArrayView<const GlobalOrdinal>& indices,
3546 const Teuchos::ArrayView<const Scalar>&
values,
3547 const char*
const prefix,
3549 const bool verbose);
3563 combineGlobalValues(
3565 const Teuchos::ArrayView<const GlobalOrdinal>&
columnIndices,
3566 const Teuchos::ArrayView<const Scalar>&
values,
3568 const char*
const prefix,
3570 const bool verbose);
3599 const char*
const prefix,
3601 const bool verbose);
3614 template <
class BinaryFunction>
3617 const Teuchos::ArrayView<const GlobalOrdinal>& indices,
3618 const Teuchos::ArrayView<const Scalar>&
values,
3620 const bool atomic = useAtomicUpdatesByDefault) {
3625 const LO
numInputEnt =
static_cast<LO
>(indices.size());
3627 return Teuchos::OrdinalTraits<LO>::invalid();
3630 const GO*
const inputCols = indices.getRawPtr();
3631 const IST*
const inputVals =
3632 reinterpret_cast<const IST*
>(values.getRawPtr());
3633 return this->transformGlobalValues(globalRow, numInputEnt, inputVals,
3634 inputCols, f, atomic);
3644 insertNonownedGlobalValues(
const GlobalOrdinal globalRow,
3645 const Teuchos::ArrayView<const GlobalOrdinal>& indices,
3646 const Teuchos::ArrayView<const Scalar>& values);
3693 const typename crs_graph_type::SLocalGlobalViews& newInds,
3694 const Teuchos::ArrayView<impl_scalar_type>& oldRowVals,
3695 const Teuchos::ArrayView<const impl_scalar_type>& newRowVals,
3696 const ELocalGlobal lg,
3697 const ELocalGlobal I);
3701 typedef Teuchos::OrdinalTraits<LocalOrdinal> OTL;
3702#if KOKKOS_VERSION >= 40799
3703 typedef KokkosKernels::ArithTraits<impl_scalar_type> STS;
3705 typedef Kokkos::ArithTraits<impl_scalar_type> STS;
3707#if KOKKOS_VERSION >= 40799
3708 typedef KokkosKernels::ArithTraits<mag_type> STM;
3710 typedef Kokkos::ArithTraits<mag_type> STM;
3712 typedef MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> MV;
3713 typedef Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> V;
3717 enum GraphAllocationStatus {
3718 GraphAlreadyAllocated,
3719 GraphNotYetAllocated
3742 const bool verbose);
3829 const bool force =
false)
const;
3854 const bool force =
false)
const;
3867 const Teuchos::ETransp
mode,
3875 typename values_dualv_type::t_host::const_type
3880 typename values_dualv_type::t_dev::const_type
3885 typename values_dualv_type::t_host
3890 typename values_dualv_type::t_dev
3898 using local_matrix_int_rowptrs_device_type =
3908 local_matrix_int_rowptrs_device_type,
3909 typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::device_view_type>;
3911 std::shared_ptr<ApplyHelper> getApplyHelper()
const {
3914 applyHelper = std::make_shared<ApplyHelper>(A_lcl.nnz(), A_lcl.graph.row_map);
3921 friend class Tpetra::crsMatrix_Swap_Tester<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
3925 template <
typename S,
typename LO,
typename GO,
typename NODE,
typename LOV>
3926 friend struct Tpetra::MMdetails::KernelWrappers;
3927 template <
typename S,
typename LO,
typename GO,
typename NODE,
typename LOV>
3928 friend struct Tpetra::MMdetails::KernelWrappers2;
3931 friend void Tpetra::MMdetails::import_and_extract_views<Scalar, LocalOrdinal, GlobalOrdinal, Node>(
3932 const CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
3933 Teuchos::RCP<
const Map<LocalOrdinal, GlobalOrdinal, Node> > targetMap,
3934 CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Aview,
3935 Teuchos::RCP<
const Import<LocalOrdinal, GlobalOrdinal, Node> > prototypeImporter,
3936 bool userAssertsThereAreNoRemotes,
3937 const std::string& label,
3938 const Teuchos::RCP<Teuchos::ParameterList>& params);
3943 void swap(CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& matrix);
3951 void fillLocalMatrix(
const Teuchos::RCP<Teuchos::ParameterList>& params);
3975 Teuchos::RCP<const Graph> staticGraph_;
3976 Teuchos::RCP<Graph> myGraph_;
3990 Details::STORAGE_1D_UNPACKED;
4022 std::map<GlobalOrdinal, std::pair<Teuchos::Array<GlobalOrdinal>,
4023 Teuchos::Array<Scalar> > >
4032 mutable std::shared_ptr<ApplyHelper> applyHelper;
4041 struct pack_functor {
4047 typedef typename DestOffsetViewType::non_const_value_type scalar_index_type;
4055 , src_offset_(src_offset)
4060 scalar_index_type
srcPos = src_offset_(row);
4061 const scalar_index_type
dstEnd = dst_offset_(row + 1);
4062 scalar_index_type
dstPos = dst_offset_(row);
4074template <
class Scalar,
4076 class GlobalOrdinal,
4078Teuchos::RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
4082 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null) {
4089template <
class CrsMatrixType>
4090Teuchos::RCP<CrsMatrixType>
4092 const Import<
typename CrsMatrixType::local_ordinal_type,
4093 typename CrsMatrixType::global_ordinal_type,
4094 typename CrsMatrixType::node_type>&
importer,
4095 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
4096 typename CrsMatrixType::global_ordinal_type,
4097 typename CrsMatrixType::node_type> >& domainMap,
4098 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
4099 typename CrsMatrixType::global_ordinal_type,
4100 typename CrsMatrixType::node_type> >&
rangeMap,
4101 const Teuchos::RCP<Teuchos::ParameterList>&
params) {
4107template <
class CrsMatrixType>
4108Teuchos::RCP<CrsMatrixType>
4110 const Import<
typename CrsMatrixType::local_ordinal_type,
4111 typename CrsMatrixType::global_ordinal_type,
4113 const Import<
typename CrsMatrixType::local_ordinal_type,
4114 typename CrsMatrixType::global_ordinal_type,
4116 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
4117 typename CrsMatrixType::global_ordinal_type,
4118 typename CrsMatrixType::node_type> >& domainMap,
4119 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
4120 typename CrsMatrixType::global_ordinal_type,
4121 typename CrsMatrixType::node_type> >&
rangeMap,
4122 const Teuchos::RCP<Teuchos::ParameterList>&
params) {
4128template <
class CrsMatrixType>
4129Teuchos::RCP<CrsMatrixType>
4131 const Export<
typename CrsMatrixType::local_ordinal_type,
4132 typename CrsMatrixType::global_ordinal_type,
4133 typename CrsMatrixType::node_type>&
exporter,
4134 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
4135 typename CrsMatrixType::global_ordinal_type,
4136 typename CrsMatrixType::node_type> >& domainMap,
4137 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
4138 typename CrsMatrixType::global_ordinal_type,
4139 typename CrsMatrixType::node_type> >&
rangeMap,
4140 const Teuchos::RCP<Teuchos::ParameterList>&
params) {
4146template <
class CrsMatrixType>
4147Teuchos::RCP<CrsMatrixType>
4149 const Export<
typename CrsMatrixType::local_ordinal_type,
4150 typename CrsMatrixType::global_ordinal_type,
4152 const Export<
typename CrsMatrixType::local_ordinal_type,
4153 typename CrsMatrixType::global_ordinal_type,
4155 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
4156 typename CrsMatrixType::global_ordinal_type,
4157 typename CrsMatrixType::node_type> >& domainMap,
4158 const Teuchos::RCP<
const Map<
typename CrsMatrixType::local_ordinal_type,
4159 typename CrsMatrixType::global_ordinal_type,
4160 typename CrsMatrixType::node_type> >&
rangeMap,
4161 const Teuchos::RCP<Teuchos::ParameterList>&
params) {
4173template <
class CrsMatrixType>
4175 typename Teuchos::ScalarTraits<typename CrsMatrixType::scalar_type>::magnitudeType
const&
threshold =
4176 Teuchos::ScalarTraits<typename CrsMatrixType::scalar_type>::magnitude(Teuchos::ScalarTraits<typename CrsMatrixType::scalar_type>::zero())) {
Forward declaration of some Tpetra Matrix Matrix objects.
Forward declaration of Tpetra::CrsMatrix.
Declaration and generic definition of traits class that tells Tpetra::CrsMatrix how to pack and unpac...
A distributed graph accessed by rows (adjacency lists) and stored sparsely.
typename local_graph_device_type::HostMirror local_graph_host_type
The type of the part of the sparse graph on each MPI process.
KokkosSparse::StaticCrsGraph< local_ordinal_type, Kokkos::LayoutLeft, device_type, void, size_t > local_graph_device_type
The type of the part of the sparse graph on each MPI process.
Struct that holds views of the contents of a CrsMatrix.
Sparse matrix that presents a row-oriented interface that lets users read or modify entries.
virtual void insertGlobalValuesImpl(crs_graph_type &graph, RowInfo &rowInfo, const GlobalOrdinal gblColInds[], const impl_scalar_type vals[], const size_t numInputEnt)
Common implementation detail of insertGlobalValues and insertGlobalValuesFiltered.
bool isGloballyIndexed() const override
Whether the matrix is globally indexed on the calling process.
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.
std::map< GlobalOrdinal, std::pair< Teuchos::Array< GlobalOrdinal >, Teuchos::Array< Scalar > > > nonlocals_
Nonlocal data added using insertGlobalValues().
void localApply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, const Teuchos::ETransp mode=Teuchos::NO_TRANS, const Scalar &alpha=Teuchos::ScalarTraits< Scalar >::one(), const Scalar &beta=Teuchos::ScalarTraits< Scalar >::zero()) const
Compute the local part of a sparse matrix-(Multi)Vector multiply.
void unpackAndCombine(const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< char *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, const CombineMode CM) override
Unpack the imported column indices and values, and combine into matrix.
CrsMatrix & operator=(const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &)=default
Copy assignment.
void replaceRangeMap(const Teuchos::RCP< const map_type > &newRangeMap)
Replace the current range Map with the given objects.
Details::EStorageStatus storageStatus_
Status of the matrix's storage, when not in a fill-complete state.
typename device_type::execution_space execution_space
The Kokkos execution space.
void applyNonTranspose(const MV &X_in, MV &Y_in, Scalar alpha, Scalar beta) const
Special case of apply() for mode == Teuchos::NO_TRANS.
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 > ¶ms=Teuchos::null) const
Import from this to the given destination matrix, and make the result fill complete.
CrsGraph< LocalOrdinal, GlobalOrdinal, Node > crs_graph_type
The CrsGraph specialization suitable for this CrsMatrix specialization.
Node node_type
This class' Kokkos Node type.
local_ordinal_type replaceGlobalValues(const global_ordinal_type globalRow, const Kokkos::View< const global_ordinal_type *, Kokkos::AnonymousSpace > &inputInds, const Kokkos::View< const impl_scalar_type *, Kokkos::AnonymousSpace > &inputVals)
Replace one or more entries' values, using global indices.
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 > &rowImporter, const Import< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &domainImporter, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &domainMap, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Nonmember CrsMatrix constructor that fuses Import and fillComplete().
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
size_t getGlobalMaxNumRowEntries() const override
Maximum number of entries in any row of the matrix, over all processes in the matrix's communicator.
void getGlobalRowCopy(GlobalOrdinal GlobalRow, nonconst_global_inds_host_view_type &Indices, nonconst_values_host_view_type &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,...
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const override
Number of entries in the sparse matrix in the given global row, on the calling (MPI) process.
void scale(const Scalar &alpha)
Scale the matrix's values: this := alpha*this.
GlobalOrdinal global_ordinal_type
The type of each global index in the matrix.
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...
void packNew(const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< char *, buffer_device_type > &exports, const Kokkos::DualView< size_t *, buffer_device_type > &numPacketsPerLID, size_t &constantNumPackets) const
Pack this object's data for an Import or Export.
size_t getLocalNumCols() const override
The number of columns connected to the locally owned rows of this matrix.
Teuchos::RCP< const map_type > getDomainMap() const override
The domain Map of this matrix.
bool hasColMap() const override
Whether the matrix has a well-defined column Map.
row_ptrs_host_view_type getLocalRowPtrsHost() const
Get a host view of the CRS packed row pointers.
mag_type getNormInf() const
Compute and return the infinity norm of the matrix.
Teuchos::RCP< CrsMatrix< T, LocalOrdinal, GlobalOrdinal, Node > > convert() const
Return another CrsMatrix with the same entries, but converted to a different Scalar type T.
values_dualv_type::t_dev getValuesViewDeviceNonConst(const RowInfo &rowinfo)
Get a non-const Device view of the locally owned values row myRow, such that rowinfo = getRowInfo(myR...
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 > ¶ms=Teuchos::null)
Perform a fillComplete on a matrix that already has data.
LocalOrdinal transformLocalValues(const LocalOrdinal lclRow, const typename UnmanagedView< LocalIndicesViewType >::type &inputInds, const typename UnmanagedView< ImplScalarViewType >::type &inputVals, BinaryFunction f, const bool atomic=useAtomicUpdatesByDefault)
Transform CrsMatrix entries in place, using local indices to select the entries in the row to transfo...
local_ordinal_type sumIntoLocalValues(const local_ordinal_type localRow, const Kokkos::View< const local_ordinal_type *, Kokkos::AnonymousSpace > &inputInds, const Kokkos::View< const impl_scalar_type *, Kokkos::AnonymousSpace > &inputVals, const bool atomic=useAtomicUpdatesByDefault)
Sum into one or more sparse matrix entries, using local row and column indices.
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 > ¶ms) const override
Implementation of RowMatrix::add: return alpha*A + beta*this.
CrsMatrix(CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &&)=default
Move constructor.
friend void batchedApply(const MatrixArray &Matrices, const typename std::remove_pointer< typename MultiVectorArray::value_type >::type &X, MultiVectorArray &Y, typename std::remove_pointer< typename MatrixArray::value_type >::type::scalar_type alpha, typename std::remove_pointer< typename MatrixArray::value_type >::type::scalar_type beta, Teuchos::RCP< Teuchos::ParameterList > params)
Does multiply matrix apply() calls with a single X vector.
DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node >::buffer_device_type buffer_device_type
Kokkos::Device specialization for communication buffers.
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.
size_t getNumEntriesInLocalRow(local_ordinal_type localRow) const override
Number of entries in the sparse matrix in the given local row, on the calling (MPI) process.
Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > createCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const size_t maxNumEntriesPerRow=0, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Create an empty CrsMatrix given a row map and a single integer upper bound on the number of stored en...
Teuchos::RCP< MV > exportMV_
Row Map MultiVector used in apply().
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const override
The communicator over which the matrix is distributed.
bool isFillActive() const
Whether the matrix is not fill complete.
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.
local_matrix_host_type::values_type getLocalValuesHost(Access::ReadWriteStruct s)
Get the Kokkos local values on host, read write.
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.
mag_type getNorm1(bool assumeSymmetric=false) const
Compute and return the 1-norm of the matrix.
virtual ~CrsMatrix()=default
Destructor (virtual for memory safety of derived classes).
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.
mag_type getFrobeniusNorm() const override
Compute and return the Frobenius norm of the matrix.
void insertLocalValues(const LocalOrdinal localRow, const Teuchos::ArrayView< const LocalOrdinal > &cols, const Teuchos::ArrayView< const Scalar > &vals, const CombineMode CM=ADD)
Insert one or more entries into the matrix, using local column indices.
global_size_t getGlobalNumCols() const override
The number of global columns in the matrix.
Teuchos::RCP< const map_type > getRangeMap() const override
The range Map of this matrix.
Teuchos::RCP< MV > importMV_
Column Map MultiVector used in apply().
void allocateValues(ELocalGlobal lg, GraphAllocationStatus gas, const bool verbose)
Allocate values (and optionally indices) using the Node.
size_t getLocalNumEntries() const override
The local number of entries in this matrix.
typename Node::device_type device_type
The Kokkos device type.
bool fillComplete_
Whether the matrix is fill complete.
virtual LocalOrdinal sumIntoGlobalValuesImpl(impl_scalar_type rowVals[], const crs_graph_type &graph, const RowInfo &rowInfo, const GlobalOrdinal inds[], const impl_scalar_type newVals[], const LocalOrdinal numElts, const bool atomic=useAtomicUpdatesByDefault)
Implementation detail of sumIntoGlobalValues.
void replaceDomainMap(const Teuchos::RCP< const map_type > &newDomainMap)
Replace the current domain Map with the given objects.
void removeCrsMatrixZeros(CrsMatrixType &matrix, typename Teuchos::ScalarTraits< typename CrsMatrixType::scalar_type >::magnitudeType const &threshold=Teuchos::ScalarTraits< typename CrsMatrixType::scalar_type >::magnitude(Teuchos::ScalarTraits< typename CrsMatrixType::scalar_type >::zero()))
Remove zero entries from a matrix.
std::string description() const override
A one-line description of this object.
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...
Teuchos::RCP< MV > getColumnMapMultiVector(const MV &X_domainMap, const bool force=false) const
Create a (or fetch a cached) column Map MultiVector.
KokkosSparse::CrsMatrix< impl_scalar_type, local_ordinal_type, device_type, void, typename local_graph_device_type::size_type > local_matrix_device_type
The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix on each MPI pro...
void replaceRangeMapAndExporter(const Teuchos::RCP< const map_type > &newRangeMap, Teuchos::RCP< const export_type > &newExporter)
Replace the current Range Map and Export with the given objects.
size_t getLocalMaxNumRowEntries() const override
Maximum number of entries in any row of the matrix, on this process.
void replaceColMap(const Teuchos::RCP< const map_type > &newColMap)
Replace the matrix's column Map with the given Map.
global_size_t getGlobalNumRows() const override
Number of global elements in the row map of this matrix.
local_matrix_host_type::values_type::const_type getLocalValuesHost(Access::ReadOnlyStruct s) const
Get the Kokkos local values on host, read only.
void globalAssemble()
Communicate nonlocal contributions to other processes.
void checkInternalState() const
Check that this object's state is sane; throw if it's not.
bool hasTransposeApply() const override
Whether apply() allows applying the transpose or conjugate transpose.
GlobalOrdinal getIndexBase() const override
The index base for global indices for this matrix.
LocalOrdinal local_ordinal_type
The type of each local index in the matrix.
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const override
Get a constant, nonpersisting view of a row of this matrix, using local row and column indices,...
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to alpha.
void fillLocalGraphAndMatrix(const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Fill data into the local graph and matrix.
row_ptrs_device_view_type getLocalRowPtrsDevice() const
Get a device view of the CRS packed row pointers.
Export< LocalOrdinal, GlobalOrdinal, Node > export_type
The Export specialization suitable for this CrsMatrix specialization.
local_inds_host_view_type getLocalIndicesHost() const
Get a host view of the CRS packed column indicies.
TPETRA_DETAILS_ALWAYS_INLINE local_matrix_device_type getLocalMatrixDevice() const
The local sparse matrix.
void getLocalRowView(LocalOrdinal LocalRow, local_inds_host_view_type &indices, values_host_view_type &values) const override
Get a constant view of a row of this matrix, using local row and column indices.
Teuchos::RCP< const map_type > getColMap() const override
The Map that describes the column distribution in this matrix.
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.
typename Kokkos::ArithTraits< impl_scalar_type >::mag_type mag_type
Type of a norm result.
void fillComplete(const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Tell the matrix that you are done changing its structure or values, and that you are ready to do comp...
typename crs_graph_type::local_graph_device_type local_graph_device_type
The part of the sparse matrix's graph on each MPI process.
void getGlobalRowView(GlobalOrdinal GlobalRow, global_inds_host_view_type &indices, values_host_view_type &values) const override
Get a constant, nonpersisting view of a row of this matrix, using global row and column indices.
void setAllValues(const typename local_graph_device_type::row_map_type &ptr, const typename local_graph_device_type::entries_type::non_const_type &ind, const typename local_matrix_device_type::values_type &val)
Set the local matrix using three (compressed sparse row) arrays.
Teuchos::RCP< const RowGraph< LocalOrdinal, GlobalOrdinal, Node > > getGraph() const override
This matrix's graph, as a RowGraph.
virtual void copyAndPermute(const SrcDistObject &source, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs, const CombineMode CM) override
virtual void removeEmptyProcessesInPlace(const Teuchos::RCP< const map_type > &newMap) override
Remove processes owning zero rows from the Maps and their communicator.
local_matrix_device_type::values_type getLocalValuesDevice(Access::ReadWriteStruct s)
Get the Kokkos local values on device, read write.
virtual LocalOrdinal sumIntoLocalValuesImpl(impl_scalar_type rowVals[], const crs_graph_type &graph, const RowInfo &rowInfo, const LocalOrdinal inds[], const impl_scalar_type newVals[], const LocalOrdinal numElts, const bool atomic=useAtomicUpdatesByDefault)
Implementation detail of sumIntoLocalValues.
void swap(CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &matrix)
Swaps the data from *this with the data and maps from crsMatrix.
bool isStaticGraph() const
Indicates that the graph is static, so that new entries cannot be added to this matrix.
global_size_t getGlobalNumEntries() const override
The global number of entries in this matrix.
virtual LocalOrdinal replaceLocalValuesImpl(impl_scalar_type rowVals[], const crs_graph_type &graph, const RowInfo &rowInfo, const LocalOrdinal inds[], const impl_scalar_type newVals[], const LocalOrdinal numElts)
Implementation detail of replaceLocalValues.
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 > ¶ms=Teuchos::null)
Nonmember CrsMatrix constructor that fuses Import and fillComplete().
size_t getLocalNumRows() const override
The number of matrix rows owned by the calling process.
bool isFillComplete() const override
Whether the matrix is fill complete.
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)
Transform CrsMatrix entries in place, using global indices to select the entries in the row to transf...
virtual bool checkSizes(const SrcDistObject &source) override
Compare the source and target (this) objects for compatibility.
Teuchos::RCP< const map_type > getRowMap() const override
The Map that describes the row distribution in this matrix.
typename device_type::memory_space memory_space
The Kokkos memory space.
local_matrix_device_type::values_type getLocalValuesDevice(Access::OverwriteAllStruct s)
Get the Kokkos local values on device, overwrite all.
local_ordinal_type replaceLocalValues(const local_ordinal_type localRow, const Kokkos::View< const local_ordinal_type *, Kokkos::AnonymousSpace > &inputInds, const Kokkos::View< const impl_scalar_type *, Kokkos::AnonymousSpace > &inputVals)
Replace one or more entries' values, using local row and column indices.
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 > ¶ms=Teuchos::null) const
Export from this to the given destination matrix, and make the result fill complete.
local_matrix_host_type::values_type getLocalValuesHost(Access::OverwriteAllStruct s)
Get the Kokkos local values on host, overwrite all.
values_dualv_type::t_host::const_type getValuesViewHost(const RowInfo &rowinfo) const
Get a const Host view of the locally owned values row myRow, such that rowinfo = getRowInfo(myRow).
CrsMatrix(const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &)=default
Copy constructor.
bool isLocallyIndexed() const override
Whether the matrix is locally indexed on the calling process.
typename row_matrix_type::impl_scalar_type impl_scalar_type
The type used internally in place of Scalar.
Teuchos::RCP< MV > getRowMapMultiVector(const MV &Y_rangeMap, const bool force=false) const
Create a (or fetch a cached) row Map MultiVector.
local_matrix_device_type::values_type::const_type getLocalValuesDevice(Access::ReadOnlyStruct s) const
Get the Kokkos local values on device, read only.
virtual LocalOrdinal replaceGlobalValuesImpl(impl_scalar_type rowVals[], const crs_graph_type &graph, const RowInfo &rowInfo, const GlobalOrdinal inds[], const impl_scalar_type newVals[], const LocalOrdinal numElts)
Implementation detail of replaceGlobalValues.
values_dualv_type::t_host getValuesViewHostNonConst(const RowInfo &rowinfo)
Get a non-const Host view of the locally owned values row myRow, such that rowinfo = getRowInfo(myRow...
void resumeFill(const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Resume operations that may change the values or structure of the matrix.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
void fillLocalMatrix(const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Fill data into the local matrix.
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x) override
Scale the matrix on the right with the given Vector.
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 > &rowExporter, const Export< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &domainExporter, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &domainMap, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Nonmember CrsMatrix constructor that fuses Export and fillComplete().
bool isStorageOptimized() const
Returns true if storage has been optimized.
Import< LocalOrdinal, GlobalOrdinal, Node > import_type
The Import specialization suitable for this CrsMatrix specialization.
void getLocalRowCopy(LocalOrdinal LocalRow, nonconst_local_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const override
Fill given arrays with a deep copy of the locally owned entries of the matrix in a given row,...
Teuchos::RCP< 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 > ¶ms=Teuchos::null)
Nonmember CrsMatrix constructor that fuses Export and fillComplete().
values_dualv_type::t_dev::const_type getValuesViewDevice(const RowInfo &rowinfo) const
Get a const Device view of the locally owned values row myRow, such that rowinfo = getRowInfo(myRow).
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x) override
Scale the matrix on the left with the given Vector.
CrsMatrix & operator=(CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &&)=default
Move assignment.
virtual bool supportsRowViews() const override
Return true if getLocalRowView() and getGlobalRowView() are valid for this object.
local_inds_device_view_type getLocalIndicesDevice() const
Get a device_view of the CRS packed column indicies.
static size_t mergeRowIndicesAndValues(size_t rowLen, local_ordinal_type *cols, impl_scalar_type *vals)
Merge duplicate row indices in the given row, along with their corresponding values.
Teuchos::RCP< const crs_graph_type > getCrsGraph() const
This matrix's graph, as a CrsGraph.
virtual LocalOrdinal getBlockSize() const override
The number of degrees of freedom per mesh point.
Keep track of how much more space a CrsGraph or CrsMatrix needs, when the graph or matrix is the targ...
A class can inherit from this if it wants to use Tpetra managed spaces.
A wrapper around Kokkos::DualView to safely manage data that might be replicated between host and dev...
Base class for distributed Tpetra objects that support data redistribution.
virtual void packAndPrepare(const SrcDistObject &source, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< packet_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, size_t &constantNumPackets)
Pack data and metadata for communication (sends).
virtual void unpackAndCombine(const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< packet_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, const CombineMode combineMode)
Perform any unpacking and combining after communication.
Communication plan for data redistribution from a (possibly) multiply-owned to a uniquely-owned distr...
Communication plan for data redistribution from a uniquely-owned to a (possibly) multiply-owned distr...
A parallel distribution of indices over processes.
One or more distributed dense vectors.
A read-only, row-oriented interface to a sparse matrix.
typename Kokkos::ArithTraits< Scalar >::val_type impl_scalar_type
The type used internally in place of Scalar.
Abstract base class for objects that can be the source of an Import or Export operation.
Implementation details of Tpetra.
EStorageStatus
Status of the graph's or matrix's storage, when not in a fill-complete state.
void residual(const Operator< SC, LO, GO, NO > &A, const MultiVector< SC, LO, GO, NO > &X, const MultiVector< SC, LO, GO, NO > &B, MultiVector< SC, LO, GO, NO > &R)
Computes R = B - A * X.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
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 > ¶ms=Teuchos::null)
Nonmember CrsMatrix constructor that fuses Export and fillComplete().
size_t global_size_t
Global size_t object.
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 > ¶ms=Teuchos::null)
Nonmember CrsMatrix constructor that fuses Import and fillComplete().
CombineMode
Rule for combining data in an Import or Export.
Allocation information for a locally owned row in a CrsGraph or CrsMatrix.