Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
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
60namespace Tpetra {
61namespace Details {
62
63namespace impl {
64
65template<class ViewType>
66ViewType
67make_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
84template<class ViewType>
85ViewType
86make_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
101template<class OutViewType, class InViewType>
102void
103assign_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
119template<class MemorySpace, class ViewType>
120auto 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
136enum class PadCrsAction {
137 INDICES_ONLY,
138 INDICES_AND_VALUES
139};
140
149template<class RowPtr, class Indices, class Values, class Padding>
150void
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 execution_space = typename Indices::t_dev::execution_space;
162 using Kokkos::view_alloc;
163 using Kokkos::WithoutInitializing;
164 using std::endl;
165 std::unique_ptr<std::string> prefix;
166
167 const size_t maxNumToPrint = verbose ?
169 if (verbose) {
170 std::ostringstream os;
171 os << "Proc " << my_rank << ": Tpetra::...::pad_crs_arrays: ";
172 prefix = std::unique_ptr<std::string>(new std::string(os.str()));
173 os << "Start" << endl;
174 std::cerr << os.str();
175 }
176 Kokkos::HostSpace hostSpace;
177
178 if (verbose) {
179 std::ostringstream os;
180 os << *prefix << "On input: ";
181 auto row_ptr_beg_h =
182 Kokkos::create_mirror_view(hostSpace, row_ptr_beg);
183 // DEEP_COPY REVIEW - NOT TESTED
184 Kokkos::deep_copy(row_ptr_beg_h, row_ptr_beg);
185 verbosePrintArray(os, row_ptr_beg_h, "row_ptr_beg before scan",
186 maxNumToPrint);
187 os << ", ";
188 auto row_ptr_end_h =
189 Kokkos::create_mirror_view(hostSpace, row_ptr_end);
190 // DEEP_COPY REVIEW - NOT TESTED
191 Kokkos::deep_copy(row_ptr_end_h, row_ptr_end);
192 verbosePrintArray(os, row_ptr_end_h, "row_ptr_end before scan",
193 maxNumToPrint);
194 os << ", indices.extent(0): " << indices_wdv.extent(0)
195 << ", values.extent(0): " << values_wdv.extent(0)
196 << ", padding: ";
197 padding.print(os);
198 os << endl;
199 std::cerr << os.str();
200 }
201
202 if (row_ptr_beg.size() == 0) {
203 if (verbose) {
204 std::ostringstream os;
205 os << *prefix << "Done; local matrix has no rows" << endl;
206 std::cerr << os.str();
207 }
208 return; // nothing to do
209 }
210
211 const size_t lclNumRows(row_ptr_beg.size() - 1);
212 RowPtr newAllocPerRow =
213 make_uninitialized_view<RowPtr>("newAllocPerRow", lclNumRows,
214 verbose, prefix.get());
215 if (verbose) {
216 std::ostringstream os;
217 os << *prefix << "Fill newAllocPerRow & compute increase" << endl;
218 std::cerr << os.str();
219 }
220 size_t increase = 0;
221 {
222 // Must do on host because padding uses std::map
223 execution_space exec_space_instance = execution_space();
224 auto row_ptr_end_h = create_mirror_view(
225 hostSpace, row_ptr_end, verbose, prefix.get());
226 // DEEP_COPY REVIEW - DEVICE-TO-HOSTMIRROR
227 Kokkos::deep_copy(exec_space_instance, row_ptr_end_h, row_ptr_end);
228 auto row_ptr_beg_h = create_mirror_view(
229 hostSpace, row_ptr_beg, verbose, prefix.get());
230 // DEEP_COPY REVIEW - DEVICE-TO-HOSTMIRROR
231 Kokkos::deep_copy(exec_space_instance, row_ptr_beg_h, row_ptr_beg);
232
233 // lbv 03/15/23: The execution space deep_copy does an asynchronous
234 // copy so we really want to fence that space before touching the
235 // host data as it is not guarenteed to have arrived by the time we
236 // start the parallel_reduce below which might use a different
237 // execution space, see:
238 // https://kokkos.github.io/kokkos-core-wiki/API/core/view/deep_copy.html#semantics
239 exec_space_instance.fence();
240
241 auto newAllocPerRow_h = create_mirror_view(
242 hostSpace, newAllocPerRow, verbose, prefix.get());
243 using host_range_type = Kokkos::RangePolicy<
244 Kokkos::DefaultHostExecutionSpace, size_t>;
245 Kokkos::parallel_reduce
246 ("Tpetra::CrsGraph: Compute new allocation size per row",
247 host_range_type(0, lclNumRows),
248 [&] (const size_t lclRowInd, size_t& lclIncrease) {
249 const size_t start = row_ptr_beg_h[lclRowInd];
250 const size_t end = row_ptr_beg_h[lclRowInd+1];
251 TEUCHOS_ASSERT( end >= start );
252 const size_t oldAllocSize = end - start;
253 const size_t oldNumEnt = row_ptr_end_h[lclRowInd] - start;
254 TEUCHOS_ASSERT( oldNumEnt <= oldAllocSize );
255
256 // This is not a pack routine. Do not shrink! Shrinking now
257 // to fit the number of entries would ignore users' hint for
258 // the max number of entries in each row. Also, CrsPadding
259 // only counts entries and ignores any available free space.
260
261 auto result = padding.get_result(lclRowInd);
262 const size_t newNumEnt = oldNumEnt + result.numInSrcNotInTgt;
263 if (newNumEnt > oldAllocSize) {
264 lclIncrease += (newNumEnt - oldAllocSize);
265 newAllocPerRow_h[lclRowInd] = newNumEnt;
266 }
267 else {
268 newAllocPerRow_h[lclRowInd] = oldAllocSize;
269 }
270 }, increase);
271
272 if (verbose) {
273 std::ostringstream os;
274 os << *prefix << "increase: " << increase << ", ";
275 verbosePrintArray(os, newAllocPerRow_h, "newAllocPerRow",
276 maxNumToPrint);
277 os << endl;
278 std::cerr << os.str();
279 }
280
281 if (increase == 0) {
282 return;
283 }
284 // DEEP_COPY REVIEW - HOSTMIRROR-TO-DEVICE
285 Kokkos::deep_copy(execution_space(), newAllocPerRow, newAllocPerRow_h);
286 }
287
288 using inds_value_type =
289 typename Indices::t_dev::non_const_value_type;
290 using vals_value_type = typename Values::t_dev::non_const_value_type;
291
292 {
293 auto indices_old = indices_wdv.getDeviceView(Access::ReadOnly);
294 const size_t newIndsSize = size_t(indices_old.size()) + increase;
295 auto indices_new = make_uninitialized_view<typename Indices::t_dev>(
296 "Tpetra::CrsGraph column indices", newIndsSize, verbose,
297 prefix.get());
298
299 typename Values::t_dev values_new;
300 auto values_old = values_wdv.getDeviceView(Access::ReadOnly);
301 if (action == PadCrsAction::INDICES_AND_VALUES) {
302 const size_t newValsSize = newIndsSize;
303 // NOTE (mfh 10 Feb 2020) If we don't initialize values_new here,
304 // then the CrsMatrix tests fail.
305 values_new = make_initialized_view<typename Values::t_dev>(
306 "Tpetra::CrsMatrix values", newValsSize, verbose, prefix.get());
307 }
308
309 if (verbose) {
310 std::ostringstream os;
311 os << *prefix << "Repack" << endl;
312 std::cerr << os.str();
313 }
314
315 using range_type = Kokkos::RangePolicy<execution_space, size_t>;
316 Kokkos::parallel_scan(
317 "Tpetra::CrsGraph or CrsMatrix repack",
318 range_type(size_t(0), size_t(lclNumRows+1)),
319 KOKKOS_LAMBDA (const size_t lclRow, size_t& newRowBeg,
320 const bool finalPass)
321 {
322 // row_ptr_beg has lclNumRows + 1 entries.
323 // row_ptr_end has lclNumRows entries.
324 // newAllocPerRow has lclNumRows entries.
325 const size_t row_beg = row_ptr_beg[lclRow];
326 const size_t row_end =
327 lclRow < lclNumRows ? row_ptr_end[lclRow] : row_beg;
328 const size_t numEnt = row_end - row_beg;
329 const size_t newRowAllocSize =
330 lclRow < lclNumRows ? newAllocPerRow[lclRow] : size_t(0);
331 if (finalPass) {
332 if (lclRow < lclNumRows) {
333 const Kokkos::pair<size_t, size_t> oldRange(
334 row_beg, row_beg + numEnt);
335 const Kokkos::pair<size_t, size_t> newRange(
336 newRowBeg, newRowBeg + numEnt);
337 auto oldColInds = Kokkos::subview(indices_old, oldRange);
338 auto newColInds = Kokkos::subview(indices_new, newRange);
339 // memcpy works fine on device; the next step is to
340 // introduce two-level parallelism and use team copy.
341 memcpy(newColInds.data(), oldColInds.data(),
342 numEnt * sizeof(inds_value_type));
343 if (action == PadCrsAction::INDICES_AND_VALUES) {
344 auto oldVals =
345 Kokkos::subview(values_old, oldRange);
346 auto newVals = Kokkos::subview(values_new, newRange);
347 memcpy(newVals.data(), oldVals.data(),
348 numEnt * sizeof(vals_value_type));
349 }
350 }
351 // It's the final pass, so we can modify these arrays.
352 row_ptr_beg[lclRow] = newRowBeg;
353 if (lclRow < lclNumRows) {
354 row_ptr_end[lclRow] = newRowBeg + numEnt;
355 }
356 }
357 newRowBeg += newRowAllocSize;
358 });
359
360 if (verbose)
361 {
362 std::ostringstream os;
363
364 os << *prefix;
365 auto row_ptr_beg_h =
366 Kokkos::create_mirror_view(hostSpace, row_ptr_beg);
367 // DEEP_COPY REVIEW - NOT TESTED
368 Kokkos::deep_copy(row_ptr_beg_h, row_ptr_beg);
369 verbosePrintArray(os, row_ptr_beg_h, "row_ptr_beg after scan",
370 maxNumToPrint);
371 os << endl;
372
373 os << *prefix;
374 auto row_ptr_end_h =
375 Kokkos::create_mirror_view(hostSpace, row_ptr_end);
376 // DEEP_COPY REVIEW - NOT TESTED
377 Kokkos::deep_copy(row_ptr_end_h, row_ptr_end);
378 verbosePrintArray(os, row_ptr_end_h, "row_ptr_end after scan",
379 maxNumToPrint);
380 os << endl;
381
382 std::cout << os.str();
383 }
384
385 indices_wdv = Indices(indices_new);
386 values_wdv = Values(values_new);
387 }
388
389 if (verbose) {
390 auto indices_h = indices_wdv.getHostView(Access::ReadOnly);
391 auto values_h = values_wdv.getHostView(Access::ReadOnly);
392 std::ostringstream os;
393 os << "On output: ";
394 verbosePrintArray(os, indices_h, "indices", maxNumToPrint);
395 os << ", ";
396 verbosePrintArray(os, values_h, "values", maxNumToPrint);
397 os << ", padding: ";
398 padding.print(os);
399 os << endl;
400 }
401
402 if (verbose) {
403 std::ostringstream os;
404 os << *prefix << "Done" << endl;
405 std::cerr << os.str();
406 }
407}
408
410template <class Pointers, class InOutIndices, class InIndices, class IndexMap>
411size_t
413 typename Pointers::value_type const row,
414 Pointers const& row_ptrs,
415 InOutIndices& cur_indices,
416 size_t& num_assigned,
417 InIndices const& new_indices,
418 IndexMap&& map,
419 std::function<void(size_t const, size_t const, size_t const)> cb)
420{
421 if (new_indices.size() == 0) {
422 return 0;
423 }
424
425 if (cur_indices.size() == 0) {
426 // No room to insert new indices
427 return Teuchos::OrdinalTraits<size_t>::invalid();
428 }
429
430 using offset_type = typename std::decay<decltype (row_ptrs[0])>::type;
431 using ordinal_type = typename std::decay<decltype (cur_indices[0])>::type;
432
433 const offset_type start = row_ptrs[row];
434 offset_type end = start + static_cast<offset_type> (num_assigned);
435 const size_t num_avail = (row_ptrs[row + 1] < end) ? size_t (0) :
436 row_ptrs[row + 1] - end;
437 const size_t num_new_indices = static_cast<size_t> (new_indices.size ());
438 size_t num_inserted = 0;
439
440 size_t numIndicesLookup = num_assigned + num_new_indices;
441
442 // Threshold determined from test/Utils/insertCrsIndicesThreshold.cpp
443 const size_t useLookUpTableThreshold = 400;
444
445 if (numIndicesLookup <= useLookUpTableThreshold || num_new_indices == 1) {
446 // For rows with few nonzeros, can use a serial search to find duplicates
447 // Or if inserting only one index, serial search is as fast as anything else
448 for (size_t k = 0; k < num_new_indices; ++k) {
449 const ordinal_type idx = std::forward<IndexMap>(map)(new_indices[k]);
450 offset_type row_offset = start;
451 for (; row_offset < end; ++row_offset) {
452 if (idx == cur_indices[row_offset]) {
453 break;
454 }
455 }
456
457 if (row_offset == end) {
458 if (num_inserted >= num_avail) { // not enough room
459 return Teuchos::OrdinalTraits<size_t>::invalid();
460 }
461 // This index is not yet in indices
462 cur_indices[end++] = idx;
463 num_inserted++;
464 }
465 if (cb) {
466 cb(k, start, row_offset - start);
467 }
468 }
469 }
470 else {
471 // For rows with many nonzeros, use a lookup table to find duplicates
472 std::unordered_map<ordinal_type, offset_type> idxLookup(numIndicesLookup);
473
474 // Put existing indices into the lookup table
475 for (size_t k = 0; k < num_assigned; k++) {
476 idxLookup[cur_indices[start+k]] = start+k;
477 }
478
479 // Check for new indices in table; insert if not there yet
480 for (size_t k = 0; k < num_new_indices; k++) {
481 const ordinal_type idx = std::forward<IndexMap>(map)(new_indices[k]);
482 offset_type row_offset;
483
484 auto it = idxLookup.find(idx);
485 if (it == idxLookup.end()) {
486 if (num_inserted >= num_avail) { // not enough room
487 return Teuchos::OrdinalTraits<size_t>::invalid();
488 }
489 // index not found; insert it
490 row_offset = end;
491 cur_indices[end++] = idx;
492 idxLookup[idx] = row_offset;
493 num_inserted++;
494 }
495 else {
496 // index found; note its position
497 row_offset = it->second;
498 }
499 if (cb) {
500 cb(k, start, row_offset - start);
501 }
502 }
503 }
504 num_assigned += num_inserted;
505 return num_inserted;
506}
507
509template <class Pointers, class Indices1, class Indices2, class IndexMap, class Callback>
510size_t
512 typename Pointers::value_type const row,
513 Pointers const& row_ptrs,
514 const size_t curNumEntries,
515 Indices1 const& cur_indices,
516 Indices2 const& new_indices,
517 IndexMap&& map,
518 Callback&& cb)
519{
520 if (new_indices.size() == 0)
521 return 0;
522
523 using ordinal =
524 typename std::remove_const<typename Indices1::value_type>::type;
525 auto invalid_ordinal = Teuchos::OrdinalTraits<ordinal>::invalid();
526
527 const size_t start = static_cast<size_t> (row_ptrs[row]);
528 const size_t end = start + curNumEntries;
529 size_t num_found = 0;
530 for (size_t k = 0; k < new_indices.size(); k++)
531 {
532 auto row_offset = start;
533 auto idx = std::forward<IndexMap>(map)(new_indices[k]);
534 if (idx == invalid_ordinal)
535 continue;
536 for (; row_offset < end; row_offset++)
537 {
538 if (idx == cur_indices[row_offset])
539 {
540 std::forward<Callback>(cb)(k, start, row_offset - start);
541 num_found++;
542 }
543 }
544 }
545 return num_found;
546}
547
548} // namespace impl
549
550
568template<class RowPtr, class Indices, class Padding>
569void
571 const RowPtr& rowPtrBeg,
572 const RowPtr& rowPtrEnd,
573 Indices& indices_wdv,
574 const Padding& padding,
575 const int my_rank,
576 const bool verbose)
577{
579 // send empty values array
580 Indices values_null;
581 pad_crs_arrays<RowPtr, Indices, Indices, Padding>(
582 impl::PadCrsAction::INDICES_ONLY, rowPtrBeg, rowPtrEnd,
583 indices_wdv, values_null, padding, my_rank, verbose);
584}
585
586template<class RowPtr, class Indices, class Values, class Padding>
587void
589 const RowPtr& rowPtrBeg,
590 const RowPtr& rowPtrEnd,
591 Indices& indices_wdv,
592 Values& values_wdv,
593 const Padding& padding,
594 const int my_rank,
595 const bool verbose)
596{
597 using impl::pad_crs_arrays;
598 pad_crs_arrays<RowPtr, Indices, Values, Padding>(
599 impl::PadCrsAction::INDICES_AND_VALUES, rowPtrBeg, rowPtrEnd,
600 indices_wdv, values_wdv, padding, my_rank, verbose);
601}
602
648template <class Pointers, class InOutIndices, class InIndices>
649size_t
651 typename Pointers::value_type const row,
652 Pointers const& rowPtrs,
653 InOutIndices& curIndices,
654 size_t& numAssigned,
655 InIndices const& newIndices,
656 std::function<void(const size_t, const size_t, const size_t)> cb =
657 std::function<void(const size_t, const size_t, const size_t)>())
658{
659 static_assert(std::is_same<typename std::remove_const<typename InOutIndices::value_type>::type,
660 typename std::remove_const<typename InIndices::value_type>::type>::value,
661 "Expected views to have same value type");
662
663 // Provide a unit map for the more general insert_indices
664 using ordinal = typename InOutIndices::value_type;
665 auto numInserted = impl::insert_crs_indices(row, rowPtrs, curIndices,
666 numAssigned, newIndices, [](ordinal const idx) { return idx; }, cb);
667 return numInserted;
668}
669
670template <class Pointers, class InOutIndices, class InIndices>
671size_t
673 typename Pointers::value_type const row,
674 Pointers const& rowPtrs,
675 InOutIndices& curIndices,
676 size_t& numAssigned,
677 InIndices const& newIndices,
678 std::function<typename InOutIndices::value_type(const typename InIndices::value_type)> map,
679 std::function<void(const size_t, const size_t, const size_t)> cb =
680 std::function<void(const size_t, const size_t, const size_t)>())
681{
682 auto numInserted = impl::insert_crs_indices(row, rowPtrs, curIndices,
683 numAssigned, newIndices, map, cb);
684 return numInserted;
685}
686
687
717template <class Pointers, class Indices1, class Indices2, class Callback>
718size_t
720 typename Pointers::value_type const row,
721 Pointers const& rowPtrs,
722 const size_t curNumEntries,
723 Indices1 const& curIndices,
724 Indices2 const& newIndices,
725 Callback&& cb)
726{
727 static_assert(std::is_same<typename std::remove_const<typename Indices1::value_type>::type,
728 typename std::remove_const<typename Indices2::value_type>::type>::value,
729 "Expected views to have same value type");
730 // Provide a unit map for the more general find_crs_indices
731 using ordinal = typename Indices2::value_type;
732 auto numFound = impl::find_crs_indices(row, rowPtrs, curNumEntries, curIndices, newIndices,
733 [=](ordinal ind){ return ind; }, cb);
734 return numFound;
735}
736
737template <class Pointers, class Indices1, class Indices2, class IndexMap, class Callback>
738size_t
740 typename Pointers::value_type const row,
741 Pointers const& rowPtrs,
742 const size_t curNumEntries,
743 Indices1 const& curIndices,
744 Indices2 const& newIndices,
745 IndexMap&& map,
746 Callback&& cb)
747{
748 return impl::find_crs_indices(row, rowPtrs, curNumEntries, curIndices, newIndices, map, cb);
749}
750
751} // namespace Details
752} // namespace Tpetra
753
754#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.
Nonmember function that computes a residual Computes R = B - A * X.
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.