MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_DroppingCommon.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_DROPPINGCOMMON_HPP
11#define MUELU_DROPPINGCOMMON_HPP
12
13// #define MUELU_COALESCE_DROP_DEBUG 1
14
15#include "Kokkos_Core.hpp"
16#include "KokkosKernels_ArithTraits.hpp"
17#include "Tpetra_Access.hpp"
18#include "Xpetra_Matrix.hpp"
19#include "Xpetra_VectorFactory.hpp"
20#include "MueLu_Utilities.hpp"
21
22namespace MueLu {
23
28enum DecisionType : char {
29 UNDECIDED = 0, // no decision has been taken yet, used for initialization
30 KEEP = 1, // keeep the entry
31 DROP = 2, // drop it
32 BOUNDARY = 3 // entry is a boundary
33};
34
35namespace Misc {
36
37template <class local_ordinal_type>
39 public:
41
42 KOKKOS_FORCEINLINE_FUNCTION
43 void operator()(local_ordinal_type rlid) const {
44 }
45};
46
51template <class local_matrix_type>
53 private:
54 using scalar_type = typename local_matrix_type::value_type;
55 using local_ordinal_type = typename local_matrix_type::ordinal_type;
56 using memory_space = typename local_matrix_type::memory_space;
57 using results_view = Kokkos::View<DecisionType*, memory_space>;
58 using boundary_nodes_view = Kokkos::View<const bool*, memory_space>;
59
60 local_matrix_type A;
63
64 public:
65 PointwiseDropBoundaryFunctor(local_matrix_type& A_, boundary_nodes_view boundaryNodes_, results_view& results_)
66 : A(A_)
67 , boundaryNodes(boundaryNodes_)
68 , results(results_) {}
69
70 KOKKOS_FORCEINLINE_FUNCTION
71 void operator()(local_ordinal_type rlid) const {
72 auto row = A.rowConst(rlid);
73 const size_t offset = A.graph.row_map(rlid);
74 const bool isBoundaryRow = boundaryNodes(rlid);
75 if (isBoundaryRow) {
76 for (local_ordinal_type k = 0; k < row.length; ++k) {
77 auto clid = row.colidx(k);
78 results(offset + k) = Kokkos::max(rlid == clid ? KEEP : DROP,
79 results(offset + k));
80 }
81 }
82 }
83};
84
89template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
91 private:
92 using local_matrix_type = typename Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_device_type;
93 using scalar_type = typename local_matrix_type::value_type;
94 using local_ordinal_type = typename local_matrix_type::ordinal_type;
95 using memory_space = typename local_matrix_type::memory_space;
96 using results_view = Kokkos::View<DecisionType*, memory_space>;
97 using boundary_nodes_view = Kokkos::View<bool*, memory_space>;
98
103
104 public:
105 PointwiseSymmetricDropBoundaryFunctor(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A_, boundary_nodes_view boundaryNodes_, results_view& results_)
106 : boundaryNodes(boundaryNodes_)
107 , results(results_) {
108 A = A_.getLocalMatrixDevice();
109 auto boundaryNodesColumn = typename boundary_nodes_view::non_const_type("boundaryNodesColumn", A_.getColMap()->getLocalNumElements());
110 auto boundaryNodesDomain = typename boundary_nodes_view::non_const_type("boundaryNodesDomain", A_.getDomainMap()->getLocalNumElements());
112 boundaryNodesCol = boundaryNodesColumn;
113 }
114
115 KOKKOS_FORCEINLINE_FUNCTION
117 auto row = A.rowConst(rlid);
118 const size_t offset = A.graph.row_map(rlid);
119 const bool isBoundaryRow = boundaryNodes(rlid);
120 for (local_ordinal_type k = 0; k < row.length; ++k) {
121 auto clid = row.colidx(k);
122 if (isBoundaryRow || boundaryNodesCol(clid))
123 results(offset + k) = Kokkos::max(rlid == clid ? KEEP : DROP,
124 results(offset + k));
125 }
126 }
127};
128
133template <class local_matrix_type>
135 private:
136 using scalar_type = typename local_matrix_type::value_type;
137 using local_ordinal_type = typename local_matrix_type::ordinal_type;
138 using memory_space = typename local_matrix_type::memory_space;
139 using results_view = Kokkos::View<DecisionType*, memory_space>;
140 using boundary_nodes_view = Kokkos::View<const bool*, memory_space>;
141 using block_indices_view_type = Kokkos::View<local_ordinal_type*, memory_space>;
142
143 local_matrix_type A;
147
148 public:
149 VectorDropBoundaryFunctor(local_matrix_type& A_, block_indices_view_type point_to_block_, boundary_nodes_view boundaryNodes_, results_view& results_)
150 : A(A_)
151 , point_to_block(point_to_block_)
152 , boundaryNodes(boundaryNodes_)
153 , results(results_) {}
154
155 KOKKOS_FORCEINLINE_FUNCTION
157 auto row = A.rowConst(rlid);
158 const size_t offset = A.graph.row_map(rlid);
159 const bool isBoundaryRow = boundaryNodes(point_to_block(rlid));
160 if (isBoundaryRow) {
161 for (local_ordinal_type k = 0; k < row.length; ++k) {
162 auto clid = row.colidx(k);
163 results(offset + k) = Kokkos::max(rlid == clid ? KEEP : DROP,
164 results(offset + k));
165 }
166 }
167 }
168};
169
174template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
176 private:
177 using local_matrix_type = typename Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_device_type;
178 using scalar_type = typename local_matrix_type::value_type;
179 using local_ordinal_type = typename local_matrix_type::ordinal_type;
180 using memory_space = typename local_matrix_type::memory_space;
181 using results_view = Kokkos::View<DecisionType*, memory_space>;
182 using boundary_nodes_view = Kokkos::View<bool*, memory_space>;
183 using block_indices_view_type = Kokkos::View<local_ordinal_type*, memory_space>;
184
191
192 public:
193 VectorSymmetricDropBoundaryFunctor(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A_, block_indices_view_type point_to_block_, block_indices_view_type ghosted_point_to_block_, boundary_nodes_view boundaryNodes_, results_view& results_)
194 : point_to_block(point_to_block_)
195 , ghosted_point_to_block(ghosted_point_to_block_)
196 , boundaryNodes(boundaryNodes_)
197 , results(results_) {
198 A = A_.getLocalMatrixDevice();
199 auto boundaryNodesColumn = typename boundary_nodes_view::non_const_type("boundaryNodesColumn", A_.getColMap()->getLocalNumElements());
200 auto boundaryNodesDomain = typename boundary_nodes_view::non_const_type("boundaryNodesDomain", A_.getDomainMap()->getLocalNumElements());
202 boundaryNodesCol = boundaryNodesColumn;
203 }
204
205 KOKKOS_FORCEINLINE_FUNCTION
207 auto row = A.rowConst(rlid);
208 const size_t offset = A.graph.row_map(rlid);
209 const bool isBoundaryRow = boundaryNodes(point_to_block(rlid));
210 for (local_ordinal_type k = 0; k < row.length; ++k) {
211 auto clid = row.colidx(k);
212 if (isBoundaryRow || boundaryNodesCol(ghosted_point_to_block(clid))) {
213 results(offset + k) = Kokkos::max(rlid == clid ? KEEP : DROP,
214 results(offset + k));
215 }
216 }
217 }
218};
219
224template <class local_matrix_type>
226 private:
227 using scalar_type = typename local_matrix_type::value_type;
228 using local_ordinal_type = typename local_matrix_type::ordinal_type;
229 using memory_space = typename local_matrix_type::memory_space;
230 using results_view = Kokkos::View<DecisionType*, memory_space>;
231
232 local_matrix_type A;
234
235 public:
236 KeepDiagonalFunctor(local_matrix_type& A_, results_view& results_)
237 : A(A_)
238 , results(results_) {}
239
240 KOKKOS_FORCEINLINE_FUNCTION
242 auto row = A.rowConst(rlid);
243 const size_t offset = A.graph.row_map(rlid);
244 for (local_ordinal_type k = 0; k < row.length; ++k) {
245 auto clid = row.colidx(k);
246 if ((rlid == clid) && (results(offset + k) != BOUNDARY)) {
247 results(offset + k) = KEEP;
248 break;
249 }
250 }
251 }
252};
253
258template <class local_matrix_type>
260 private:
261 using scalar_type = typename local_matrix_type::value_type;
262 using local_ordinal_type = typename local_matrix_type::ordinal_type;
263 using memory_space = typename local_matrix_type::memory_space;
264 using results_view = Kokkos::View<DecisionType*, memory_space>;
265
266 local_matrix_type A;
268
269 public:
270 DropOffRankFunctor(local_matrix_type& A_, results_view& results_)
271 : A(A_)
272 , results(results_) {}
273
274 KOKKOS_FORCEINLINE_FUNCTION
276 auto row = A.rowConst(rlid);
277 const size_t offset = A.graph.row_map(rlid);
278 for (local_ordinal_type k = 0; k < row.length; ++k) {
279 auto clid = row.colidx(k);
280 if (clid >= A.numRows()) {
281 results(offset + k) = Kokkos::max(DROP, results(offset + k));
282 }
283 }
284 }
285};
286
291template <class local_matrix_type>
293 private:
294 using scalar_type = typename local_matrix_type::value_type;
295 using local_ordinal_type = typename local_matrix_type::ordinal_type;
296 using memory_space = typename local_matrix_type::memory_space;
297 using results_view = Kokkos::View<DecisionType*, memory_space>;
298
299 using boundary_nodes_view = Kokkos::View<bool*, memory_space>;
300
301 local_matrix_type A;
304
305 public:
306 MarkSingletonFunctor(local_matrix_type& A_, boundary_nodes_view boundaryNodes_, results_view& results_)
307 : A(A_)
308 , boundaryNodes(boundaryNodes_)
309 , results(results_) {}
310
311 KOKKOS_FORCEINLINE_FUNCTION
313 auto row = A.rowConst(rlid);
314 const size_t offset = A.graph.row_map(rlid);
315 for (local_ordinal_type k = 0; k < row.length; ++k) {
316 auto clid = row.colidx(k);
317 if ((results(offset + k) == KEEP) && (rlid != clid))
318 return;
319 }
320 boundaryNodes(rlid) = true;
321 for (local_ordinal_type k = 0; k < row.length; ++k) {
322 auto clid = row.colidx(k);
323 if (rlid == clid)
324 results(offset + k) = KEEP;
325 else
326 results(offset + k) = BOUNDARY;
327 }
328 }
329};
330
335template <class local_matrix_type>
337 private:
338 using scalar_type = typename local_matrix_type::value_type;
339 using local_ordinal_type = typename local_matrix_type::ordinal_type;
340 using memory_space = typename local_matrix_type::memory_space;
341 using results_view = Kokkos::View<DecisionType*, memory_space>;
342 using block_indices_view_type = Kokkos::View<local_ordinal_type*, memory_space>;
343
344 using boundary_nodes_view = Kokkos::View<bool*, memory_space>;
345
346 local_matrix_type A;
350
351 public:
352 MarkSingletonVectorFunctor(local_matrix_type& A_, block_indices_view_type point_to_block_, boundary_nodes_view boundaryNodes_, results_view& results_)
353 : A(A_)
354 , point_to_block(point_to_block_)
355 , boundaryNodes(boundaryNodes_)
356 , results(results_) {}
357
358 KOKKOS_FORCEINLINE_FUNCTION
360 auto row = A.rowConst(rlid);
361 const size_t offset = A.graph.row_map(rlid);
362 for (local_ordinal_type k = 0; k < row.length; ++k) {
363 auto clid = row.colidx(k);
364 if ((results(offset + k) == KEEP) && (rlid != clid))
365 return;
366 }
367 auto brlid = point_to_block(rlid);
368 boundaryNodes(brlid) = true;
369 for (local_ordinal_type k = 0; k < row.length; ++k) {
370 auto clid = row.colidx(k);
371 if (rlid == clid)
372 results(offset + k) = KEEP;
373 else
374 results(offset + k) = BOUNDARY;
375 }
376 }
377};
378
383template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
385 private:
386 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
387 using local_matrix_type = typename matrix_type::local_matrix_device_type;
388
389 using scalar_type = typename local_matrix_type::value_type;
390 using local_ordinal_type = typename local_matrix_type::ordinal_type;
391 using memory_space = typename local_matrix_type::memory_space;
392 using results_view = Kokkos::View<DecisionType*, memory_space>;
393
394 using block_indices_type = Xpetra::MultiVector<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>;
395 using local_block_indices_view_type = typename block_indices_type::dual_view_type_const::t_dev;
396
399 Teuchos::RCP<block_indices_type> ghosted_point_to_blockMV;
402
403 public:
405 : A(A_.getLocalMatrixDevice())
406 , point_to_block(point_to_block_.getLocalViewDevice(Tpetra::Access::ReadOnly))
407 , results(results_) {
408 auto importer = A_.getCrsGraph()->getImporter();
409
410 if (!importer.is_null()) {
411 ghosted_point_to_blockMV = Xpetra::VectorFactory<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(), false);
412 ghosted_point_to_blockMV->doImport(point_to_block_, *importer, Xpetra::INSERT);
413 ghosted_point_to_block = ghosted_point_to_blockMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
414 } else
416 }
417
418 KOKKOS_FORCEINLINE_FUNCTION
420 auto row = A.rowConst(rlid);
421 const size_t offset = A.graph.row_map(rlid);
422 for (local_ordinal_type k = 0; k < row.length; ++k) {
423 auto clid = row.colidx(k);
424 if (point_to_block(rlid, 0) == ghosted_point_to_block(clid, 0)) {
425 results(offset + k) = Kokkos::max(KEEP, results(offset + k));
426 } else {
427 results(offset + k) = Kokkos::max(DROP, results(offset + k));
428 }
429 }
430 }
431};
432
437template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
439 private:
440 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
441 using local_matrix_type = typename matrix_type::local_matrix_device_type;
442
443 using scalar_type = typename local_matrix_type::value_type;
444 using local_ordinal_type = typename local_matrix_type::ordinal_type;
445 using memory_space = typename local_matrix_type::memory_space;
446 using results_view = Kokkos::View<DecisionType*, memory_space>;
447
448 using block_indices_type = Xpetra::Vector<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>;
449 using local_block_indices_view_type = typename block_indices_type::dual_view_type_const::t_dev;
450 using id_translation_type = Kokkos::View<local_ordinal_type*, memory_space>;
451 using importer_type = Xpetra::Import<LocalOrdinal, GlobalOrdinal, Node>;
452
459 Teuchos::RCP<block_indices_type> ghosted_point_to_blockMV;
460
461 public:
462 BlockDiagonalizeVectorFunctor(matrix_type& A_, block_indices_type& point_to_block_, const RCP<const importer_type>& importer, results_view& results_, id_translation_type row_translation_, id_translation_type col_translation_)
463 : A(A_.getLocalMatrixDevice())
464 , point_to_block(point_to_block_.getLocalViewDevice(Tpetra::Access::ReadOnly))
465 , results(results_)
466 , row_translation(row_translation_)
467 , col_translation(col_translation_) {
468 if (!importer.is_null()) {
469 ghosted_point_to_blockMV = Xpetra::VectorFactory<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(), false);
470 ghosted_point_to_blockMV->doImport(point_to_block_, *importer, Xpetra::INSERT);
471 ghosted_point_to_block = ghosted_point_to_blockMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
472 } else
474 }
475
476 KOKKOS_FORCEINLINE_FUNCTION
478 auto row = A.rowConst(rlid);
479 const size_t offset = A.graph.row_map(rlid);
480 auto brlid = row_translation(rlid);
481 for (local_ordinal_type k = 0; k < row.length; ++k) {
482 auto clid = row.colidx(k);
483 auto bclid = col_translation(clid);
484 if (point_to_block(brlid, 0) == ghosted_point_to_block(bclid, 0)) {
485 results(offset + k) = Kokkos::max(KEEP, results(offset + k));
486 } else {
487 results(offset + k) = Kokkos::max(DROP, results(offset + k));
488 }
489 }
490 }
491};
492
497template <class local_matrix_type>
499 private:
500 using scalar_type = typename local_matrix_type::value_type;
501 using local_ordinal_type = typename local_matrix_type::ordinal_type;
502 using memory_space = typename local_matrix_type::memory_space;
503 using results_view = Kokkos::View<DecisionType*, memory_space>;
504
505 using boundary_nodes_view = Kokkos::View<bool*, memory_space>;
506
507 local_matrix_type A;
509
510 public:
511 DebugFunctor(local_matrix_type& A_, results_view& results_)
512 : A(A_)
513 , results(results_) {}
514
515 KOKKOS_FORCEINLINE_FUNCTION
517 auto row = A.rowConst(rlid);
518 const size_t offset = A.graph.row_map(rlid);
519 for (local_ordinal_type k = 0; k < row.length; ++k) {
520 if (results(offset + k) == UNDECIDED) {
521 Kokkos::printf("No dropping decision was taken for entry (%d, %d)\n", rlid, row.colidx(k));
522 assert(false);
523 }
524 }
525 }
526};
527
532template <class local_matrix_type>
534 private:
535 using scalar_type = typename local_matrix_type::value_type;
536 using local_ordinal_type = typename local_matrix_type::ordinal_type;
537 using memory_space = typename local_matrix_type::memory_space;
538 using results_view = Kokkos::View<DecisionType*, memory_space>;
539
540 local_matrix_type A;
542
543 public:
544 SymmetrizeFunctor(local_matrix_type& A_, results_view& results_)
545 : A(A_)
546 , results(results_) {}
547
548 KOKKOS_FORCEINLINE_FUNCTION
550 auto row = A.rowConst(rlid);
551 const size_t offset = A.graph.row_map(rlid);
552 for (local_ordinal_type k = 0; k < row.length; ++k) {
553 if (results(offset + k) == KEEP) {
554 auto clid = row.colidx(k);
555 if (clid >= A.numRows())
556 continue;
557 auto row2 = A.rowConst(clid);
558 const size_t offset2 = A.graph.row_map(clid);
559 for (local_ordinal_type k2 = 0; k2 < row2.length; ++k2) {
560 auto clid2 = row2.colidx(k2);
561 if (clid2 == rlid) {
562 if (results(offset2 + k2) == DROP)
563 results(offset2 + k2) = KEEP;
564 break;
565 }
566 }
567 }
568 }
569 }
570};
571
572template <class view_type, class comparator_type>
573KOKKOS_INLINE_FUNCTION void serialHeapSort(view_type& v, comparator_type comparator) {
574 auto N = v.extent(0);
575 size_t start = N / 2;
576 size_t end = N;
577 while (end > 1) {
578 if (start > 0)
579 start = start - 1;
580 else {
581 end = end - 1;
582 auto temp = v(0);
583 v(0) = v(end);
584 v(end) = temp;
585 }
586 size_t root = start;
587 while (2 * root + 1 < end) {
588 size_t child = 2 * root + 1;
589 if ((child + 1 < end) and (comparator(v(child), v(child + 1))))
590 ++child;
591
592 if (comparator(v(root), v(child))) {
593 auto temp = v(root);
594 v(root) = v(child);
595 v(child) = temp;
596 root = child;
597 } else
598 break;
599 }
600 }
601}
602
605enum StrengthMeasure : int {
606 /*
607 \f[
608 \frac{|A_{ij}|^2}{|A_{ii}| |A_{jj}|} \le \theta^2
609 \f]
610 */
612 /*
613 \f[
614 \frac{-\operatorname{Re}A_{ij}}{| max_j -A_{ij}|} \le \theta
615 \f]
616 */
618
619 /*
620 \f[
621 \frac{-\operatorname{sign}(A_{ij}) |A_{ij}|^2}{|A_{ii}| |A_{jj}|} \le \theta^2
622 \f]
623 */
625
626 /*
627 \f[
628 |A_{ij}| \le \theta
629 \f]
630 */
633
634} // namespace Misc
635
636} // namespace MueLu
637
638#define MUELU_ETI_SLGN_SoC(CLASSNAME, SC, LO, GO, NO) \
639 template class CLASSNAME<SC, LO, GO, NO, MueLu::Misc::SmoothedAggregationMeasure>; \
640 template class CLASSNAME<SC, LO, GO, NO, MueLu::Misc::SignedRugeStuebenMeasure>; \
641 template class CLASSNAME<SC, LO, GO, NO, MueLu::Misc::SignedSmoothedAggregationMeasure>; \
642 template class CLASSNAME<SC, LO, GO, NO, MueLu::Misc::UnscaledMeasure>;
643
644#endif
Functor that drops all entries that are not on the block diagonal.
Teuchos::RCP< block_indices_type > ghosted_point_to_blockMV
typename matrix_type::local_matrix_device_type local_matrix_type
Xpetra::MultiVector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > block_indices_type
typename local_matrix_type::value_type scalar_type
local_block_indices_view_type point_to_block
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
local_block_indices_view_type ghosted_point_to_block
typename block_indices_type::dual_view_type_const::t_dev local_block_indices_view_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
typename local_matrix_type::memory_space memory_space
typename local_matrix_type::ordinal_type local_ordinal_type
Kokkos::View< DecisionType *, memory_space > results_view
BlockDiagonalizeFunctor(matrix_type &A_, block_indices_type &point_to_block_, results_view &results_)
Functor that drops all entries that are not on the block diagonal.
Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > block_indices_type
typename local_matrix_type::memory_space memory_space
Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > importer_type
typename block_indices_type::dual_view_type_const::t_dev local_block_indices_view_type
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
BlockDiagonalizeVectorFunctor(matrix_type &A_, block_indices_type &point_to_block_, const RCP< const importer_type > &importer, results_view &results_, id_translation_type row_translation_, id_translation_type col_translation_)
typename local_matrix_type::value_type scalar_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
Kokkos::View< DecisionType *, memory_space > results_view
typename local_matrix_type::ordinal_type local_ordinal_type
Teuchos::RCP< block_indices_type > ghosted_point_to_blockMV
Kokkos::View< local_ordinal_type *, memory_space > id_translation_type
typename matrix_type::local_matrix_device_type local_matrix_type
local_block_indices_view_type ghosted_point_to_block
Functor that checks that all entries have been marked.
Kokkos::View< DecisionType *, memory_space > results_view
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
DebugFunctor(local_matrix_type &A_, results_view &results_)
typename local_matrix_type::ordinal_type local_ordinal_type
Kokkos::View< bool *, memory_space > boundary_nodes_view
typename local_matrix_type::value_type scalar_type
typename local_matrix_type::memory_space memory_space
Functor that drops off-rank entries.
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
Kokkos::View< DecisionType *, memory_space > results_view
typename local_matrix_type::ordinal_type local_ordinal_type
typename local_matrix_type::memory_space memory_space
typename local_matrix_type::value_type scalar_type
DropOffRankFunctor(local_matrix_type &A_, results_view &results_)
Functor that marks diagonal as kept, unless the are already marked as boundary.
Kokkos::View< DecisionType *, memory_space > results_view
typename local_matrix_type::memory_space memory_space
KeepDiagonalFunctor(local_matrix_type &A_, results_view &results_)
typename local_matrix_type::ordinal_type local_ordinal_type
typename local_matrix_type::value_type scalar_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
Functor that marks singletons (all off-diagonal entries in a row are dropped) as boundary.
Kokkos::View< DecisionType *, memory_space > results_view
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
typename local_matrix_type::value_type scalar_type
MarkSingletonFunctor(local_matrix_type &A_, boundary_nodes_view boundaryNodes_, results_view &results_)
typename local_matrix_type::memory_space memory_space
typename local_matrix_type::ordinal_type local_ordinal_type
Kokkos::View< bool *, memory_space > boundary_nodes_view
Functor that marks singletons (all off-diagonal entries in a row are dropped) as boundary.
MarkSingletonVectorFunctor(local_matrix_type &A_, block_indices_view_type point_to_block_, boundary_nodes_view boundaryNodes_, results_view &results_)
Kokkos::View< DecisionType *, memory_space > results_view
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
typename local_matrix_type::value_type scalar_type
Kokkos::View< local_ordinal_type *, memory_space > block_indices_view_type
typename local_matrix_type::ordinal_type local_ordinal_type
Kokkos::View< bool *, memory_space > boundary_nodes_view
typename local_matrix_type::memory_space memory_space
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
Functor that drops boundary nodes for a blockSize == 1 problem.
typename local_matrix_type::ordinal_type local_ordinal_type
typename local_matrix_type::memory_space memory_space
Kokkos::View< const bool *, memory_space > boundary_nodes_view
typename local_matrix_type::value_type scalar_type
Kokkos::View< DecisionType *, memory_space > results_view
PointwiseDropBoundaryFunctor(local_matrix_type &A_, boundary_nodes_view boundaryNodes_, results_view &results_)
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
Functor that drops boundary nodes for a blockSize == 1 problem.
Kokkos::View< DecisionType *, memory_space > results_view
typename Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_device_type local_matrix_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
typename local_matrix_type::value_type scalar_type
typename local_matrix_type::ordinal_type local_ordinal_type
Kokkos::View< bool *, memory_space > boundary_nodes_view
typename local_matrix_type::memory_space memory_space
PointwiseSymmetricDropBoundaryFunctor(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A_, boundary_nodes_view boundaryNodes_, results_view &results_)
Functor that symmetrizes the dropping decisions.
typename local_matrix_type::value_type scalar_type
SymmetrizeFunctor(local_matrix_type &A_, results_view &results_)
Kokkos::View< DecisionType *, memory_space > results_view
typename local_matrix_type::memory_space memory_space
typename local_matrix_type::ordinal_type local_ordinal_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
Functor that drops boundary nodes for a blockSize > 1 problem.
VectorDropBoundaryFunctor(local_matrix_type &A_, block_indices_view_type point_to_block_, boundary_nodes_view boundaryNodes_, results_view &results_)
Kokkos::View< local_ordinal_type *, memory_space > block_indices_view_type
Kokkos::View< const bool *, memory_space > boundary_nodes_view
typename local_matrix_type::memory_space memory_space
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
typename local_matrix_type::ordinal_type local_ordinal_type
Kokkos::View< DecisionType *, memory_space > results_view
typename local_matrix_type::value_type scalar_type
Functor that drops boundary nodes for a blockSize > 1 problem.
typename local_matrix_type::value_type scalar_type
typename local_matrix_type::memory_space memory_space
typename Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_device_type local_matrix_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
Kokkos::View< bool *, memory_space > boundary_nodes_view
Kokkos::View< DecisionType *, memory_space > results_view
VectorSymmetricDropBoundaryFunctor(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A_, block_indices_view_type point_to_block_, block_indices_view_type ghosted_point_to_block_, boundary_nodes_view boundaryNodes_, results_view &results_)
Kokkos::View< local_ordinal_type *, memory_space > block_indices_view_type
typename local_matrix_type::ordinal_type local_ordinal_type
static void DetectDirichletColsAndDomains(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayRCP< bool > &dirichletRows, Teuchos::ArrayRCP< bool > dirichletCols, Teuchos::ArrayRCP< bool > dirichletDomain)
Detects Dirichlet columns & domains from a list of Dirichlet rows.
KOKKOS_INLINE_FUNCTION void serialHeapSort(view_type &v, comparator_type comparator)
Namespace for MueLu classes and methods.