Tpetra parallel linear algebra  Version of the Day
Tpetra_Details_crsUtils.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 // ************************************************************************
38 // @HEADER
39 
40 #ifndef TPETRA_DETAILS_CRSUTILS_HPP
41 #define TPETRA_DETAILS_CRSUTILS_HPP
42 #include <numeric>
43 #include <type_traits>
44 
45 #include "TpetraCore_config.h"
46 #include "Kokkos_Core.hpp"
48 #include "Tpetra_Details_CrsPadding.hpp"
49 #include "Tpetra_Details_WrappedDualView.hpp"
50 #include <iostream>
51 #include <memory>
52 #include <unordered_map>
53 
59 
60 namespace Tpetra {
61 namespace Details {
62 
63 namespace impl {
64 
65 template<class ViewType>
66 ViewType
67 make_uninitialized_view(
68  const std::string& name,
69  const size_t size,
70  const bool verbose,
71  const std::string* const prefix)
72 {
73  if (verbose) {
74  std::ostringstream os;
75  os << *prefix << "Allocate Kokkos::View " << name
76  << ": " << size << std::endl;
77  std::cerr << os.str();
78  }
79  using Kokkos::view_alloc;
80  using Kokkos::WithoutInitializing;
81  return ViewType(view_alloc(name, WithoutInitializing), size);
82 }
83 
84 template<class ViewType>
85 ViewType
86 make_initialized_view(
87  const std::string& name,
88  const size_t size,
89  const bool verbose,
90  const std::string* const prefix)
91 {
92  if (verbose) {
93  std::ostringstream os;
94  os << *prefix << "Allocate & initialize Kokkos::View "
95  << name << ": " << size << std::endl;
96  std::cerr << os.str();
97  }
98  return ViewType(name, size);
99 }
100 
101 template<class OutViewType, class InViewType>
102 void
103 assign_to_view(OutViewType& out,
104  const InViewType& in,
105  const char viewName[],
106  const bool verbose,
107  const std::string* const prefix)
108 {
109  if (verbose) {
110  std::ostringstream os;
111  os << *prefix << "Assign to Kokkos::View " << viewName
112  << ": Old size: " << out.extent(0)
113  << ", New size: " << in.extent(0) << std::endl;
114  std::cerr << os.str();
115  }
116  out = in;
117 }
118 
119 template<class MemorySpace, class ViewType>
120 auto create_mirror_view(
121  const MemorySpace& memSpace,
122  const ViewType& view,
123  const bool verbose,
124  const std::string* const prefix) ->
125  decltype(Kokkos::create_mirror_view(memSpace, view))
126 {
127  if (verbose) {
128  std::ostringstream os;
129  os << *prefix << "Create mirror view: "
130  << "view.extent(0): " << view.extent(0) << std::endl;
131  std::cerr << os.str();
132  }
133  return Kokkos::create_mirror_view(memSpace, view);
134 }
135 
136 enum class PadCrsAction {
137  INDICES_ONLY,
138  INDICES_AND_VALUES
139 };
140 
149 template<class RowPtr, class Indices, class Values, class Padding>
150 void
152  const PadCrsAction action,
153  const RowPtr& row_ptr_beg,
154  const RowPtr& row_ptr_end,
155  Indices& indices_wdv,
156  Values& values_wdv,
157  const Padding& padding,
158  const int my_rank,
159  const bool verbose)
160 {
161  using Kokkos::view_alloc;
162  using Kokkos::WithoutInitializing;
163  using std::endl;
164  std::unique_ptr<std::string> prefix;
165 
166  const size_t maxNumToPrint = verbose ?
168  if (verbose) {
169  std::ostringstream os;
170  os << "Proc " << my_rank << ": Tpetra::...::pad_crs_arrays: ";
171  prefix = std::unique_ptr<std::string>(new std::string(os.str()));
172  os << "Start" << endl;
173  std::cerr << os.str();
174  }
175  Kokkos::HostSpace hostSpace;
176 
177  if (verbose) {
178  std::ostringstream os;
179  os << *prefix << "On input: ";
180  auto row_ptr_beg_h =
181  Kokkos::create_mirror_view(hostSpace, row_ptr_beg);
182  Kokkos::deep_copy(row_ptr_beg_h, row_ptr_beg);
183  verbosePrintArray(os, row_ptr_beg_h, "row_ptr_beg before scan",
184  maxNumToPrint);
185  os << ", ";
186  auto row_ptr_end_h =
187  Kokkos::create_mirror_view(hostSpace, row_ptr_end);
188  Kokkos::deep_copy(row_ptr_end_h, row_ptr_end);
189  verbosePrintArray(os, row_ptr_end_h, "row_ptr_end before scan",
190  maxNumToPrint);
191  os << ", indices.extent(0): " << indices_wdv.extent(0)
192  << ", values.extent(0): " << values_wdv.extent(0)
193  << ", padding: ";
194  padding.print(os);
195  os << endl;
196  std::cerr << os.str();
197  }
198 
199  if (row_ptr_beg.size() == 0) {
200  if (verbose) {
201  std::ostringstream os;
202  os << *prefix << "Done; local matrix has no rows" << endl;
203  std::cerr << os.str();
204  }
205  return; // nothing to do
206  }
207 
208  const size_t lclNumRows(row_ptr_beg.size() - 1);
209  RowPtr newAllocPerRow =
210  make_uninitialized_view<RowPtr>("newAllocPerRow", lclNumRows,
211  verbose, prefix.get());
212  if (verbose) {
213  std::ostringstream os;
214  os << *prefix << "Fill newAllocPerRow & compute increase" << endl;
215  std::cerr << os.str();
216  }
217  size_t increase = 0;
218  {
219  // Must do on host because padding uses std::map
220  auto row_ptr_end_h = create_mirror_view(
221  hostSpace, row_ptr_end, verbose, prefix.get());
222  Kokkos::deep_copy(row_ptr_end_h, row_ptr_end);
223  auto row_ptr_beg_h = create_mirror_view(
224  hostSpace, row_ptr_beg, verbose, prefix.get());
225  Kokkos::deep_copy(row_ptr_beg_h, row_ptr_beg);
226 
227  auto newAllocPerRow_h = create_mirror_view(
228  hostSpace, newAllocPerRow, verbose, prefix.get());
229  using host_range_type = Kokkos::RangePolicy<
230  Kokkos::DefaultHostExecutionSpace, size_t>;
231  Kokkos::parallel_reduce
232  ("Tpetra::CrsGraph: Compute new allocation size per row",
233  host_range_type(0, lclNumRows),
234  [&] (const size_t lclRowInd, size_t& lclIncrease) {
235  const size_t start = row_ptr_beg_h[lclRowInd];
236  const size_t end = row_ptr_beg_h[lclRowInd+1];
237  TEUCHOS_ASSERT( end >= start );
238  const size_t oldAllocSize = end - start;
239  const size_t oldNumEnt = row_ptr_end_h[lclRowInd] - start;
240  TEUCHOS_ASSERT( oldNumEnt <= oldAllocSize );
241 
242  // This is not a pack routine. Do not shrink! Shrinking now
243  // to fit the number of entries would ignore users' hint for
244  // the max number of entries in each row. Also, CrsPadding
245  // only counts entries and ignores any available free space.
246 
247  auto result = padding.get_result(lclRowInd);
248  const size_t newNumEnt = oldNumEnt + result.numInSrcNotInTgt;
249  if (newNumEnt > oldAllocSize) {
250  lclIncrease += (newNumEnt - oldAllocSize);
251  newAllocPerRow_h[lclRowInd] = newNumEnt;
252  }
253  else {
254  newAllocPerRow_h[lclRowInd] = oldAllocSize;
255  }
256  }, increase);
257 
258  if (verbose) {
259  std::ostringstream os;
260  os << *prefix << "increase: " << increase << ", ";
261  verbosePrintArray(os, newAllocPerRow_h, "newAllocPerRow",
262  maxNumToPrint);
263  os << endl;
264  std::cerr << os.str();
265  }
266 
267  if (increase == 0) {
268  return;
269  }
270  Kokkos::deep_copy(newAllocPerRow, newAllocPerRow_h);
271  }
272 
273  using inds_value_type =
274  typename Indices::DeviceViewType::non_const_value_type;
275  using vals_value_type = typename Values::DeviceViewType::non_const_value_type;
276 
277  auto indices_old = indices_wdv.getDeviceView(Access::ReadOnly);
278  const size_t newIndsSize = size_t(indices_old.size()) + increase;
279  auto indices_new = make_uninitialized_view<typename Indices::DeviceViewType>(
280  "Tpetra::CrsGraph column indices", newIndsSize, verbose,
281  prefix.get());
282 
283  typename Values::DeviceViewType values_new;
284  auto values_old = values_wdv.getDeviceView(Access::ReadOnly);
285  if (action == PadCrsAction::INDICES_AND_VALUES) {
286  const size_t newValsSize = newIndsSize;
287  // NOTE (mfh 10 Feb 2020) If we don't initialize values_new here,
288  // then the CrsMatrix tests fail.
289  values_new = make_initialized_view<typename Values::DeviceViewType>(
290  "Tpetra::CrsMatrix values", newValsSize, verbose, prefix.get());
291  }
292 
293  if (verbose) {
294  std::ostringstream os;
295  os << *prefix << "Repack" << endl;
296  std::cerr << os.str();
297  }
298  using execution_space = typename Indices::DeviceViewType::execution_space;
299  using range_type = Kokkos::RangePolicy<execution_space, size_t>;
300  Kokkos::parallel_scan(
301  "Tpetra::CrsGraph or CrsMatrix repack",
302  range_type(size_t(0), size_t(lclNumRows+1)),
303  KOKKOS_LAMBDA (const size_t lclRow, size_t& newRowBeg,
304  const bool finalPass)
305  {
306  // row_ptr_beg has lclNumRows + 1 entries.
307  // row_ptr_end has lclNumRows entries.
308  // newAllocPerRow has lclNumRows entries.
309  const size_t row_beg = row_ptr_beg[lclRow];
310  const size_t row_end =
311  lclRow < lclNumRows ? row_ptr_end[lclRow] : row_beg;
312  const size_t numEnt = row_end - row_beg;
313  const size_t newRowAllocSize =
314  lclRow < lclNumRows ? newAllocPerRow[lclRow] : size_t(0);
315  if (finalPass) {
316  if (lclRow < lclNumRows) {
317  const Kokkos::pair<size_t, size_t> oldRange(
318  row_beg, row_beg + numEnt);
319  const Kokkos::pair<size_t, size_t> newRange(
320  newRowBeg, newRowBeg + numEnt);
321  auto oldColInds = Kokkos::subview(indices_old, oldRange);
322  auto newColInds = Kokkos::subview(indices_new, newRange);
323  // memcpy works fine on device; the next step is to
324  // introduce two-level parallelism and use team copy.
325  memcpy(newColInds.data(), oldColInds.data(),
326  numEnt * sizeof(inds_value_type));
327  if (action == PadCrsAction::INDICES_AND_VALUES) {
328  auto oldVals =
329  Kokkos::subview(values_old, oldRange);
330  auto newVals = Kokkos::subview(values_new, newRange);
331  memcpy(newVals.data(), oldVals.data(),
332  numEnt * sizeof(vals_value_type));
333  }
334  }
335  // It's the final pass, so we can modify these arrays.
336  row_ptr_beg[lclRow] = newRowBeg;
337  if (lclRow < lclNumRows) {
338  row_ptr_end[lclRow] = newRowBeg + numEnt;
339  }
340  }
341  newRowBeg += newRowAllocSize;
342  });
343 
344  if (verbose)
345  {
346  std::ostringstream os;
347 
348  os << *prefix;
349  auto row_ptr_beg_h =
350  Kokkos::create_mirror_view(hostSpace, row_ptr_beg);
351  Kokkos::deep_copy(row_ptr_beg_h, row_ptr_beg);
352  verbosePrintArray(os, row_ptr_beg_h, "row_ptr_beg after scan",
353  maxNumToPrint);
354  os << endl;
355 
356  os << *prefix;
357  auto row_ptr_end_h =
358  Kokkos::create_mirror_view(hostSpace, row_ptr_end);
359  Kokkos::deep_copy(row_ptr_end_h, row_ptr_end);
360  verbosePrintArray(os, row_ptr_end_h, "row_ptr_end after scan",
361  maxNumToPrint);
362  os << endl;
363 
364  std::cout << os.str();
365  }
366 
367  indices_wdv = Indices(indices_new);
368  values_wdv = Values(values_new);
369 
370  if (verbose) {
371  auto indices_h = indices_wdv.getHostView(Access::ReadOnly);
372  auto values_h = values_wdv.getHostView(Access::ReadOnly);
373  std::ostringstream os;
374  os << "On output: ";
375  verbosePrintArray(os, indices_h, "indices", maxNumToPrint);
376  os << ", ";
377  verbosePrintArray(os, values_h, "values", maxNumToPrint);
378  os << ", padding: ";
379  padding.print(os);
380  os << endl;
381  }
382 
383  if (verbose) {
384  std::ostringstream os;
385  os << *prefix << "Done" << endl;
386  std::cerr << os.str();
387  }
388 }
389 
391 template <class Pointers, class InOutIndices, class InIndices, class IndexMap>
392 size_t
394  typename Pointers::value_type const row,
395  Pointers const& row_ptrs,
396  InOutIndices& cur_indices,
397  size_t& num_assigned,
398  InIndices const& new_indices,
399  IndexMap&& map,
400  std::function<void(size_t const, size_t const, size_t const)> cb)
401 {
402  if (new_indices.size() == 0) {
403  return 0;
404  }
405 
406  if (cur_indices.size() == 0) {
407  // No room to insert new indices
408  return Teuchos::OrdinalTraits<size_t>::invalid();
409  }
410 
411  using offset_type = typename std::decay<decltype (row_ptrs[0])>::type;
412  using ordinal_type = typename std::decay<decltype (cur_indices[0])>::type;
413 
414  const offset_type start = row_ptrs[row];
415  offset_type end = start + static_cast<offset_type> (num_assigned);
416  const size_t num_avail = (row_ptrs[row + 1] < end) ? size_t (0) :
417  row_ptrs[row + 1] - end;
418  const size_t num_new_indices = static_cast<size_t> (new_indices.size ());
419  size_t num_inserted = 0;
420 
421  size_t numIndicesLookup = num_assigned + num_new_indices;
422 
423  // Threshold determined from test/Utils/insertCrsIndicesThreshold.cpp
424  const size_t useLookUpTableThreshold = 400;
425 
426  if (numIndicesLookup <= useLookUpTableThreshold || num_new_indices == 1) {
427  // For rows with few nonzeros, can use a serial search to find duplicates
428  // Or if inserting only one index, serial search is as fast as anything else
429  for (size_t k = 0; k < num_new_indices; ++k) {
430  const ordinal_type idx = std::forward<IndexMap>(map)(new_indices[k]);
431  offset_type row_offset = start;
432  for (; row_offset < end; ++row_offset) {
433  if (idx == cur_indices[row_offset]) {
434  break;
435  }
436  }
437 
438  if (row_offset == end) {
439  if (num_inserted >= num_avail) { // not enough room
440  return Teuchos::OrdinalTraits<size_t>::invalid();
441  }
442  // This index is not yet in indices
443  cur_indices[end++] = idx;
444  num_inserted++;
445  }
446  if (cb) {
447  cb(k, start, row_offset - start);
448  }
449  }
450  }
451  else {
452  // For rows with many nonzeros, use a lookup table to find duplicates
453  std::unordered_map<ordinal_type, offset_type> idxLookup(numIndicesLookup);
454 
455  // Put existing indices into the lookup table
456  for (size_t k = 0; k < num_assigned; k++) {
457  idxLookup[cur_indices[start+k]] = start+k;
458  }
459 
460  // Check for new indices in table; insert if not there yet
461  for (size_t k = 0; k < num_new_indices; k++) {
462  const ordinal_type idx = std::forward<IndexMap>(map)(new_indices[k]);
463  offset_type row_offset;
464 
465  auto it = idxLookup.find(idx);
466  if (it == idxLookup.end()) {
467  if (num_inserted >= num_avail) { // not enough room
468  return Teuchos::OrdinalTraits<size_t>::invalid();
469  }
470  // index not found; insert it
471  row_offset = end;
472  cur_indices[end++] = idx;
473  idxLookup[idx] = row_offset;
474  num_inserted++;
475  }
476  else {
477  // index found; note its position
478  row_offset = it->second;
479  }
480  if (cb) {
481  cb(k, start, row_offset - start);
482  }
483  }
484  }
485  num_assigned += num_inserted;
486  return num_inserted;
487 }
488 
490 template <class Pointers, class Indices1, class Indices2, class IndexMap, class Callback>
491 size_t
493  typename Pointers::value_type const row,
494  Pointers const& row_ptrs,
495  const size_t curNumEntries,
496  Indices1 const& cur_indices,
497  Indices2 const& new_indices,
498  IndexMap&& map,
499  Callback&& cb)
500 {
501  if (new_indices.size() == 0)
502  return 0;
503 
504  using ordinal =
505  typename std::remove_const<typename Indices1::value_type>::type;
506  auto invalid_ordinal = Teuchos::OrdinalTraits<ordinal>::invalid();
507 
508  const size_t start = static_cast<size_t> (row_ptrs[row]);
509  const size_t end = start + curNumEntries;
510  size_t num_found = 0;
511  for (size_t k = 0; k < new_indices.size(); k++)
512  {
513  auto row_offset = start;
514  auto idx = std::forward<IndexMap>(map)(new_indices[k]);
515  if (idx == invalid_ordinal)
516  continue;
517  for (; row_offset < end; row_offset++)
518  {
519  if (idx == cur_indices[row_offset])
520  {
521  std::forward<Callback>(cb)(k, start, row_offset - start);
522  num_found++;
523  }
524  }
525  }
526  return num_found;
527 }
528 
529 } // namespace impl
530 
531 
549 template<class RowPtr, class Indices, class Padding>
550 void
552  const RowPtr& rowPtrBeg,
553  const RowPtr& rowPtrEnd,
554  Indices& indices_wdv,
555  const Padding& padding,
556  const int my_rank,
557  const bool verbose)
558 {
559  using impl::pad_crs_arrays;
560  // send empty values array
561  Indices values_null;
562  pad_crs_arrays<RowPtr, Indices, Indices, Padding>(
563  impl::PadCrsAction::INDICES_ONLY, rowPtrBeg, rowPtrEnd,
564  indices_wdv, values_null, padding, my_rank, verbose);
565 }
566 
567 template<class RowPtr, class Indices, class Values, class Padding>
568 void
570  const RowPtr& rowPtrBeg,
571  const RowPtr& rowPtrEnd,
572  Indices& indices_wdv,
573  Values& values_wdv,
574  const Padding& padding,
575  const int my_rank,
576  const bool verbose)
577 {
578  using impl::pad_crs_arrays;
579  pad_crs_arrays<RowPtr, Indices, Values, Padding>(
580  impl::PadCrsAction::INDICES_AND_VALUES, rowPtrBeg, rowPtrEnd,
581  indices_wdv, values_wdv, padding, my_rank, verbose);
582 }
583 
629 template <class Pointers, class InOutIndices, class InIndices>
630 size_t
632  typename Pointers::value_type const row,
633  Pointers const& rowPtrs,
634  InOutIndices& curIndices,
635  size_t& numAssigned,
636  InIndices const& newIndices,
637  std::function<void(const size_t, const size_t, const size_t)> cb =
638  std::function<void(const size_t, const size_t, const size_t)>())
639 {
640  static_assert(std::is_same<typename std::remove_const<typename InOutIndices::value_type>::type,
641  typename std::remove_const<typename InIndices::value_type>::type>::value,
642  "Expected views to have same value type");
643 
644  // Provide a unit map for the more general insert_indices
645  using ordinal = typename InOutIndices::value_type;
646  auto numInserted = impl::insert_crs_indices(row, rowPtrs, curIndices,
647  numAssigned, newIndices, [](ordinal const idx) { return idx; }, cb);
648  return numInserted;
649 }
650 
651 template <class Pointers, class InOutIndices, class InIndices>
652 size_t
654  typename Pointers::value_type const row,
655  Pointers const& rowPtrs,
656  InOutIndices& curIndices,
657  size_t& numAssigned,
658  InIndices const& newIndices,
659  std::function<typename InOutIndices::value_type(const typename InIndices::value_type)> map,
660  std::function<void(const size_t, const size_t, const size_t)> cb =
661  std::function<void(const size_t, const size_t, const size_t)>())
662 {
663  auto numInserted = impl::insert_crs_indices(row, rowPtrs, curIndices,
664  numAssigned, newIndices, map, cb);
665  return numInserted;
666 }
667 
668 
698 template <class Pointers, class Indices1, class Indices2, class Callback>
699 size_t
701  typename Pointers::value_type const row,
702  Pointers const& rowPtrs,
703  const size_t curNumEntries,
704  Indices1 const& curIndices,
705  Indices2 const& newIndices,
706  Callback&& cb)
707 {
708  static_assert(std::is_same<typename std::remove_const<typename Indices1::value_type>::type,
709  typename std::remove_const<typename Indices2::value_type>::type>::value,
710  "Expected views to have same value type");
711  // Provide a unit map for the more general find_crs_indices
712  using ordinal = typename Indices2::value_type;
713  auto numFound = impl::find_crs_indices(row, rowPtrs, curNumEntries, curIndices, newIndices,
714  [=](ordinal ind){ return ind; }, cb);
715  return numFound;
716 }
717 
718 template <class Pointers, class Indices1, class Indices2, class IndexMap, class Callback>
719 size_t
721  typename Pointers::value_type const row,
722  Pointers const& rowPtrs,
723  const size_t curNumEntries,
724  Indices1 const& curIndices,
725  Indices2 const& newIndices,
726  IndexMap&& map,
727  Callback&& cb)
728 {
729  return impl::find_crs_indices(row, rowPtrs, curNumEntries, curIndices, newIndices, map, cb);
730 }
731 
732 } // namespace Details
733 } // namespace Tpetra
734 
735 #endif // TPETRA_DETAILS_CRSUTILS_HPP
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.
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.
void deep_copy(MultiVector< DS, DL, DG, DN > &dst, const MultiVector< SS, SL, SG, SN > &src)
Copy the contents of the MultiVector src into dst.