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) {
 
   43    std::ostringstream os;
 
   44    os << *prefix << 
"Allocate Kokkos::View " << name
 
   45       << 
": " << size << std::endl;
 
   46    std::cerr << os.str();
 
   48  using Kokkos::view_alloc;
 
   49  using Kokkos::WithoutInitializing;
 
   50  return ViewType(view_alloc(name, WithoutInitializing), size);
 
   53template <
class ViewType>
 
   56    const std::string& name,
 
   59    const std::string* 
const prefix) {
 
   61    std::ostringstream os;
 
   62    os << *prefix << 
"Allocate & initialize Kokkos::View " 
   63       << name << 
": " << size << std::endl;
 
   64    std::cerr << os.str();
 
   66  return ViewType(name, size);
 
   69template <
class OutViewType, 
class InViewType>
 
   70void assign_to_view(OutViewType& out,
 
   72                    const char viewName[],
 
   74                    const std::string* 
const prefix) {
 
   76    std::ostringstream os;
 
   77    os << *prefix << 
"Assign to Kokkos::View " << viewName
 
   78       << 
": Old size: " << out.extent(0)
 
   79       << 
", New size: " << in.extent(0) << std::endl;
 
   80    std::cerr << os.str();
 
   85template <
class MemorySpace, 
class ViewType>
 
   86auto create_mirror_view(
 
   87    const MemorySpace& memSpace,
 
   90    const std::string* 
const prefix) -> 
decltype(Kokkos::create_mirror_view(memSpace, view)) {
 
   92    std::ostringstream os;
 
   93    os << *prefix << 
"Create mirror view: " 
   94       << 
"view.extent(0): " << view.extent(0) << std::endl;
 
   95    std::cerr << os.str();
 
   97  return Kokkos::create_mirror_view(memSpace, view);
 
  100enum class PadCrsAction {
 
  113template <
class RowPtr, 
class Indices, 
class Values, 
class Padding>
 
  115    const PadCrsAction 
action,
 
  122    const bool verbose) {
 
  123  using execution_space = 
typename Indices::t_dev::execution_space;
 
  124  using Kokkos::view_alloc;
 
  125  using Kokkos::WithoutInitializing;
 
  127  std::unique_ptr<std::string> 
prefix;
 
  131    std::ostringstream 
os;
 
  132    os << 
"Proc " << 
my_rank << 
": Tpetra::...::pad_crs_arrays: ";
 
  133    prefix = std::unique_ptr<std::string>(
new std::string(
os.str()));
 
  135    std::cerr << 
os.str();
 
  140    std::ostringstream 
os;
 
  156       << 
", values.extent(0): " << 
values_wdv.extent(0)
 
  160    std::cerr << 
os.str();
 
  165      std::ostringstream 
os;
 
  166      os << *
prefix << 
"Done; local matrix has no rows" << 
endl;
 
  167      std::cerr << 
os.str();
 
  177    std::ostringstream 
os;
 
  178    os << *
prefix << 
"Fill newAllocPerRow & compute increase" << 
endl;
 
  179    std::cerr << 
os.str();
 
  205        Kokkos::DefaultHostExecutionSpace, 
size_t>;
 
  206    Kokkos::parallel_reduce(
 
  207        "Tpetra::CrsGraph: Compute new allocation size per row",
 
  234      std::ostringstream 
os;
 
  239      std::cerr << 
os.str();
 
  250      typename Indices::t_dev::non_const_value_type;
 
  257        "Tpetra::CrsGraph column indices", 
newIndsSize, verbose,
 
  262    if (
action == PadCrsAction::INDICES_AND_VALUES) {
 
  271      std::ostringstream 
os;
 
  273      std::cerr << 
os.str();
 
  276    using range_type = Kokkos::RangePolicy<execution_space, size_t>;
 
  277    Kokkos::parallel_scan(
 
  278        "Tpetra::CrsGraph or CrsMatrix repack",
 
  279        range_type(
size_t(0), 
size_t(
lclNumRows + 1)),
 
  293              const Kokkos::pair<size_t, size_t> 
oldRange(
 
  295              const Kokkos::pair<size_t, size_t> 
newRange(
 
  303              if (
action == PadCrsAction::INDICES_AND_VALUES) {
 
  321      std::ostringstream 
os;
 
  341      std::cout << 
os.str();
 
  351    std::ostringstream 
os;
 
  362    std::ostringstream 
os;
 
  364    std::cerr << 
os.str();
 
 
  369template <
class Po
inters, 
class InOutIndices, 
class InIndices, 
class IndexMap>
 
  372    typename Pointers::value_type 
const row,
 
  378    std::function<
void(
size_t const, 
size_t const, 
size_t const)> 
cb) {
 
  385    return Teuchos::OrdinalTraits<size_t>::invalid();
 
  388  using offset_type  = 
typename std::decay<
decltype(
row_ptrs[0])>::type;
 
  389  using ordinal_type = 
typename std::decay<
decltype(
cur_indices[0])>::type;
 
  391  const offset_type start      = 
row_ptrs[row];
 
  416          return Teuchos::OrdinalTraits<size_t>::invalid();
 
  443          return Teuchos::OrdinalTraits<size_t>::invalid();
 
 
  464template <
class Po
inters, 
class Indices1, 
class Indices2, 
class IndexMap, 
class Callback>
 
  467    typename Pointers::value_type 
const row,
 
  478      typename std::remove_const<typename Indices1::value_type>::type;
 
  481  const size_t start = 
static_cast<size_t>(
row_ptrs[row]);
 
 
  518template <
class RowPtr, 
class Indices, 
class Padding>
 
  525    const bool verbose) {
 
  526  using impl::pad_crs_arrays;
 
 
  534template <
class RowPtr, 
class Indices, 
class Values, 
class Padding>
 
  542    const bool verbose) {
 
  543  using impl::pad_crs_arrays;
 
  594template <
class Po
inters, 
class InOutIndices, 
class InIndices>
 
  597    typename Pointers::value_type 
const row,
 
  602    std::function<
void(
const size_t, 
const size_t, 
const size_t)> 
cb =
 
  603        std::function<
void(
const size_t, 
const size_t, 
const size_t)>()) {
 
  604  static_assert(std::is_same<typename std::remove_const<typename InOutIndices::value_type>::type,
 
  605                             typename std::remove_const<typename InIndices::value_type>::type>::value,
 
  606                "Expected views to have same value type");
 
  609  using ordinal    = 
typename InOutIndices::value_type;
 
 
  616template <
class Po
inters, 
class InOutIndices, 
class InIndices>
 
  619    typename Pointers::value_type 
const row,
 
  624    std::function<
typename InOutIndices::value_type(
const typename InIndices::value_type)> 
map,
 
  625    std::function<
void(
const size_t, 
const size_t, 
const size_t)> 
cb =
 
  626        std::function<
void(
const size_t, 
const size_t, 
const size_t)>()) {
 
  661template <
class Po
inters, 
class Indices1, 
class Indices2, 
class Callback>
 
  664    typename Pointers::value_type 
const row,
 
  670  static_assert(std::is_same<typename std::remove_const<typename Indices1::value_type>::type,
 
  671                             typename std::remove_const<typename Indices2::value_type>::type>::value,
 
  672                "Expected views to have same value type");
 
  674  using ordinal = 
typename Indices2::value_type;
 
  675  auto numFound = impl::find_crs_indices(
 
 
  681template <
class Po
inters, 
class Indices1, 
class Indices2, 
class IndexMap, 
class Callback>
 
  684    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.