10#ifndef TPETRA_DETAILS_CRSUTILS_HPP
11#define TPETRA_DETAILS_CRSUTILS_HPP
15#include "TpetraCore_config.h"
16#include "Kokkos_Core.hpp"
18#include "Tpetra_Details_CrsPadding.hpp"
19#include "Tpetra_Details_WrappedDualView.hpp"
22#include <unordered_map>
36template <
class ViewType>
38make_uninitialized_view(
39 const std::string& name,
42 const std::string*
const prefix) {
44 std::ostringstream os;
45 os << *prefix <<
"Allocate Kokkos::View " << name
46 <<
": " << size << std::endl;
47 std::cerr << os.str();
49 using Kokkos::view_alloc;
50 using Kokkos::WithoutInitializing;
51 return ViewType(view_alloc(name, WithoutInitializing), size);
54template <
class ViewType>
57 const std::string& name,
60 const std::string*
const prefix) {
62 std::ostringstream os;
63 os << *prefix <<
"Allocate & initialize Kokkos::View "
64 << name <<
": " << size << std::endl;
65 std::cerr << os.str();
67 return ViewType(name, size);
70template <
class OutViewType,
class InViewType>
71void assign_to_view(OutViewType& out,
73 const char viewName[],
75 const std::string*
const prefix) {
77 std::ostringstream os;
78 os << *prefix <<
"Assign to Kokkos::View " << viewName
79 <<
": Old size: " << out.extent(0)
80 <<
", New size: " << in.extent(0) << std::endl;
81 std::cerr << os.str();
86template <
class MemorySpace,
class ViewType>
87auto create_mirror_view(
88 const MemorySpace& memSpace,
91 const std::string*
const prefix) ->
decltype(Kokkos::create_mirror_view(memSpace, view)) {
93 std::ostringstream os;
94 os << *prefix <<
"Create mirror view: "
95 <<
"view.extent(0): " << view.extent(0) << std::endl;
96 std::cerr << os.str();
98 return Kokkos::create_mirror_view(memSpace, view);
101enum class PadCrsAction {
114template <
class RowPtr,
class Indices,
class Values,
class Padding>
116 const PadCrsAction
action,
123 const bool verbose) {
124 using execution_space =
typename Indices::t_dev::execution_space;
125 using Kokkos::view_alloc;
126 using Kokkos::WithoutInitializing;
128 std::unique_ptr<std::string>
prefix;
132 std::ostringstream
os;
133 os <<
"Proc " <<
my_rank <<
": Tpetra::...::pad_crs_arrays: ";
134 prefix = std::unique_ptr<std::string>(
new std::string(
os.str()));
136 std::cerr <<
os.str();
141 std::ostringstream
os;
157 <<
", values.extent(0): " <<
values_wdv.extent(0)
161 std::cerr <<
os.str();
166 std::ostringstream
os;
167 os << *
prefix <<
"Done; local matrix has no rows" <<
endl;
168 std::cerr <<
os.str();
178 std::ostringstream
os;
179 os << *
prefix <<
"Fill newAllocPerRow & compute increase" <<
endl;
180 std::cerr <<
os.str();
206 Kokkos::DefaultHostExecutionSpace,
size_t>;
207 Kokkos::parallel_reduce(
208 "Tpetra::CrsGraph: Compute new allocation size per row",
235 std::ostringstream
os;
240 std::cerr <<
os.str();
251 typename Indices::t_dev::non_const_value_type;
258 "Tpetra::CrsGraph column indices",
newIndsSize, verbose,
263 if (
action == PadCrsAction::INDICES_AND_VALUES) {
272 std::ostringstream
os;
274 std::cerr <<
os.str();
277 using range_type = Kokkos::RangePolicy<execution_space, size_t>;
278 Kokkos::parallel_scan(
279 "Tpetra::CrsGraph or CrsMatrix repack",
280 range_type(
size_t(0),
size_t(
lclNumRows + 1)),
294 const Kokkos::pair<size_t, size_t>
oldRange(
296 const Kokkos::pair<size_t, size_t>
newRange(
304 if (
action == PadCrsAction::INDICES_AND_VALUES) {
322 std::ostringstream
os;
342 std::cout <<
os.str();
352 std::ostringstream
os;
363 std::ostringstream
os;
365 std::cerr <<
os.str();
370template <
class Po
inters,
class InOutIndices,
class InIndices,
class IndexMap>
373 typename Pointers::value_type
const row,
379 std::function<
void(
size_t const,
size_t const,
size_t const)>
cb) {
386 return Teuchos::OrdinalTraits<size_t>::invalid();
389 using offset_type =
typename std::decay<
decltype(
row_ptrs[0])>::type;
390 using ordinal_type =
typename std::decay<
decltype(
cur_indices[0])>::type;
392 const offset_type start =
row_ptrs[row];
417 return Teuchos::OrdinalTraits<size_t>::invalid();
444 return Teuchos::OrdinalTraits<size_t>::invalid();
465template <
class Po
inters,
class Indices1,
class Indices2,
class IndexMap,
class Callback>
468 typename Pointers::value_type
const row,
479 typename std::remove_const<typename Indices1::value_type>::type;
482 const size_t start =
static_cast<size_t>(
row_ptrs[row]);
493 std::forward<Callback>(
cb)(
k, start,
off);
502template <
class Po
inters,
class Indices1,
class Indices2,
class IndexMap,
class Callback>
504 typename Pointers::value_type
const row,
514 using ordinal =
typename std::remove_const<typename Indices1::value_type>::type;
517 const size_t start =
static_cast<size_t>(
row_ptrs[row]);
532 std::forward<Callback>(
cb)(
k, start,
off);
558template <
class RowPtr,
class Indices,
class Padding>
565 const bool verbose) {
566 using impl::pad_crs_arrays;
574template <
class RowPtr,
class Indices,
class Values,
class Padding>
582 const bool verbose) {
583 using impl::pad_crs_arrays;
634template <
class Po
inters,
class InOutIndices,
class InIndices>
637 typename Pointers::value_type
const row,
642 std::function<
void(
const size_t,
const size_t,
const size_t)>
cb =
643 std::function<
void(
const size_t,
const size_t,
const size_t)>()) {
644 static_assert(std::is_same<typename std::remove_const<typename InOutIndices::value_type>::type,
645 typename std::remove_const<typename InIndices::value_type>::type>::value,
646 "Expected views to have same value type");
649 using ordinal =
typename InOutIndices::value_type;
656template <
class Po
inters,
class InOutIndices,
class InIndices>
659 typename Pointers::value_type
const row,
664 std::function<
typename InOutIndices::value_type(
const typename InIndices::value_type)>
map,
665 std::function<
void(
const size_t,
const size_t,
const size_t)>
cb =
666 std::function<
void(
const size_t,
const size_t,
const size_t)>()) {
701template <
class Po
inters,
class Indices1,
class Indices2,
class Callback>
704 typename Pointers::value_type
const row,
710 static_assert(std::is_same<typename std::remove_const<typename Indices1::value_type>::type,
711 typename std::remove_const<typename Indices2::value_type>::type>::value,
712 "Expected views to have same value type");
714 using ordinal =
typename Indices2::value_type;
715 auto numFound = impl::find_crs_indices(
721template <
class Po
inters,
class Indices1,
class Indices2,
class IndexMap,
class Callback>
724 typename Pointers::value_type
const row,
734template <
class Po
inters,
class Indices1,
class Indices2,
class IndexMap,
class Callback>
735size_t findCrsIndicesSorted(
736 typename Pointers::value_type
const row,
737 Pointers
const& rowPtrs,
738 const size_t curNumEntries,
739 Indices1
const& curIndices,
740 Indices2
const& newIndices,
743 return impl::find_crs_indices_sorted(row, rowPtrs, curNumEntries, curIndices, newIndices, map, cb);
Declaration of Tpetra::Details::Behavior, a class that describes Tpetra's behavior.
void pad_crs_arrays(const PadCrsAction action, const RowPtr &row_ptr_beg, const RowPtr &row_ptr_end, Indices &indices_wdv, Values &values_wdv, const Padding &padding, const int my_rank, const bool verbose)
Implementation of padCrsArrays.
size_t find_crs_indices(typename Pointers::value_type const row, Pointers const &row_ptrs, const size_t curNumEntries, Indices1 const &cur_indices, Indices2 const &new_indices, IndexMap &&map, Callback &&cb)
Implementation of findCrsIndices.
size_t find_crs_indices_sorted(typename Pointers::value_type const row, Pointers const &row_ptrs, const size_t curNumEntries, Indices1 const &cur_indices, Indices2 const &new_indices, IndexMap &&map, Callback &&cb)
Implementation of findCrsIndices.
size_t insert_crs_indices(typename Pointers::value_type const row, Pointers const &row_ptrs, InOutIndices &cur_indices, size_t &num_assigned, InIndices const &new_indices, IndexMap &&map, std::function< void(size_t const, size_t const, size_t const)> cb)
Implementation of insertCrsIndices.
Struct that holds views of the contents of a CrsMatrix.
static size_t verbosePrintCountThreshold()
Number of entries below which arrays, lists, etc. will be printed in debug mode.
Implementation details of Tpetra.
void padCrsArrays(const RowPtr &rowPtrBeg, const RowPtr &rowPtrEnd, Indices &indices_wdv, const Padding &padding, const int my_rank, const bool verbose)
Determine if the row pointers and indices arrays need to be resized to accommodate new entries....
void verbosePrintArray(std::ostream &out, const ArrayType &x, const char name[], const size_t maxNumToPrint)
Print min(x.size(), maxNumToPrint) entries of x.
size_t insertCrsIndices(typename Pointers::value_type const row, Pointers const &rowPtrs, InOutIndices &curIndices, size_t &numAssigned, InIndices const &newIndices, std::function< void(const size_t, const size_t, const size_t)> cb=std::function< void(const size_t, const size_t, const size_t)>())
Insert new indices in to current list of indices.
size_t findCrsIndices(typename Pointers::value_type const row, Pointers const &rowPtrs, const size_t curNumEntries, Indices1 const &curIndices, Indices2 const &newIndices, Callback &&cb)
Finds offsets in to current list of indices.
Namespace Tpetra contains the class and methods constituting the Tpetra library.