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>
35template<
class ViewType>
37make_uninitialized_view(
38 const std::string& name,
41 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)
63 std::ostringstream os;
64 os << *prefix <<
"Allocate & initialize Kokkos::View "
65 << name <<
": " << size << std::endl;
66 std::cerr << os.str();
68 return ViewType(name, size);
71template<
class OutViewType,
class InViewType>
73assign_to_view(OutViewType& out,
75 const char viewName[],
77 const std::string*
const prefix)
80 std::ostringstream os;
81 os << *prefix <<
"Assign to Kokkos::View " << viewName
82 <<
": Old size: " << out.extent(0)
83 <<
", New size: " << in.extent(0) << std::endl;
84 std::cerr << os.str();
89template<
class MemorySpace,
class ViewType>
90auto create_mirror_view(
91 const MemorySpace& memSpace,
94 const std::string*
const prefix) ->
95 decltype(Kokkos::create_mirror_view(memSpace, view))
98 std::ostringstream os;
99 os << *prefix <<
"Create mirror view: "
100 <<
"view.extent(0): " << view.extent(0) << std::endl;
101 std::cerr << os.str();
103 return Kokkos::create_mirror_view(memSpace, view);
106enum class PadCrsAction {
119template<
class RowPtr,
class Indices,
class Values,
class Padding>
122 const PadCrsAction
action,
131 using execution_space =
typename Indices::t_dev::execution_space;
132 using Kokkos::view_alloc;
133 using Kokkos::WithoutInitializing;
135 std::unique_ptr<std::string>
prefix;
140 std::ostringstream
os;
141 os <<
"Proc " <<
my_rank <<
": Tpetra::...::pad_crs_arrays: ";
142 prefix = std::unique_ptr<std::string>(
new std::string(
os.str()));
144 std::cerr <<
os.str();
149 std::ostringstream
os;
165 <<
", values.extent(0): " <<
values_wdv.extent(0)
169 std::cerr <<
os.str();
174 std::ostringstream
os;
175 os << *
prefix <<
"Done; local matrix has no rows" <<
endl;
176 std::cerr <<
os.str();
186 std::ostringstream
os;
187 os << *
prefix <<
"Fill newAllocPerRow & compute increase" <<
endl;
188 std::cerr <<
os.str();
214 Kokkos::DefaultHostExecutionSpace,
size_t>;
215 Kokkos::parallel_reduce
216 (
"Tpetra::CrsGraph: Compute new allocation size per row",
243 std::ostringstream
os;
248 std::cerr <<
os.str();
259 typename Indices::t_dev::non_const_value_type;
266 "Tpetra::CrsGraph column indices",
newIndsSize, verbose,
271 if (
action == PadCrsAction::INDICES_AND_VALUES) {
280 std::ostringstream
os;
282 std::cerr <<
os.str();
285 using range_type = Kokkos::RangePolicy<execution_space, size_t>;
286 Kokkos::parallel_scan(
287 "Tpetra::CrsGraph or CrsMatrix repack",
303 const Kokkos::pair<size_t, size_t>
oldRange(
305 const Kokkos::pair<size_t, size_t>
newRange(
313 if (
action == PadCrsAction::INDICES_AND_VALUES) {
332 std::ostringstream
os;
352 std::cout <<
os.str();
362 std::ostringstream
os;
373 std::ostringstream
os;
375 std::cerr <<
os.str();
380template <
class Po
inters,
class InOutIndices,
class InIndices,
class IndexMap>
383 typename Pointers::value_type
const row,
389 std::function<
void(
size_t const,
size_t const,
size_t const)>
cb)
397 return Teuchos::OrdinalTraits<size_t>::invalid();
400 using offset_type =
typename std::decay<
decltype (
row_ptrs[0])>::type;
401 using ordinal_type =
typename std::decay<
decltype (
cur_indices[0])>::type;
403 const offset_type start =
row_ptrs[row];
429 return Teuchos::OrdinalTraits<size_t>::invalid();
457 return Teuchos::OrdinalTraits<size_t>::invalid();
479template <
class Po
inters,
class Indices1,
class Indices2,
class IndexMap,
class Callback>
482 typename Pointers::value_type
const row,
494 typename std::remove_const<typename Indices1::value_type>::type;
497 const size_t start =
static_cast<size_t> (
row_ptrs[row]);
538template<
class RowPtr,
class Indices,
class Padding>
548 using impl::pad_crs_arrays;
556template<
class RowPtr,
class Indices,
class Values,
class Padding>
567 using impl::pad_crs_arrays;
618template <
class Po
inters,
class InOutIndices,
class InIndices>
621 typename Pointers::value_type
const row,
626 std::function<
void(
const size_t,
const size_t,
const size_t)>
cb =
627 std::function<
void(
const size_t,
const size_t,
const size_t)>())
629 static_assert(std::is_same<typename std::remove_const<typename InOutIndices::value_type>::type,
630 typename std::remove_const<typename InIndices::value_type>::type>::value,
631 "Expected views to have same value type");
634 using ordinal =
typename InOutIndices::value_type;
640template <
class Po
inters,
class InOutIndices,
class InIndices>
643 typename Pointers::value_type
const row,
648 std::function<
typename InOutIndices::value_type(
const typename InIndices::value_type)>
map,
649 std::function<
void(
const size_t,
const size_t,
const size_t)>
cb =
650 std::function<
void(
const size_t,
const size_t,
const size_t)>())
687template <
class Po
inters,
class Indices1,
class Indices2,
class Callback>
690 typename Pointers::value_type
const row,
697 static_assert(std::is_same<typename std::remove_const<typename Indices1::value_type>::type,
698 typename std::remove_const<typename Indices2::value_type>::type>::value,
699 "Expected views to have same value type");
701 using ordinal =
typename Indices2::value_type;
707template <
class Po
inters,
class Indices1,
class Indices2,
class IndexMap,
class Callback>
710 typename Pointers::value_type
const row,
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 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.