MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_UtilitiesBase_def.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_UTILITIESBASE_DEF_HPP
11#define MUELU_UTILITIESBASE_DEF_HPP
12
13#include "MueLu_ConfigDefs.hpp"
14
16
17#include <Kokkos_Core.hpp>
18#include <KokkosSparse_CrsMatrix.hpp>
19#include <KokkosSparse_getDiagCopy.hpp>
20
21#include <Xpetra_BlockedVector.hpp>
22#include <Xpetra_BlockedMap.hpp>
23#include <Xpetra_BlockedMultiVector.hpp>
24#include <Xpetra_ExportFactory.hpp>
25
26#include <Xpetra_Import.hpp>
27#include <Xpetra_ImportFactory.hpp>
28#include <Xpetra_MatrixMatrix.hpp>
29#include <Xpetra_CrsGraph.hpp>
30#include <Xpetra_CrsGraphFactory.hpp>
31#include <Xpetra_CrsMatrixWrap.hpp>
32#include <Xpetra_StridedMap.hpp>
33
34#include "MueLu_Exceptions.hpp"
35#include "Xpetra_CrsMatrixFactory.hpp"
36
37#include <KokkosKernels_Handle.hpp>
38#include <KokkosGraph_RCM.hpp>
39
40namespace MueLu {
41
42template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
43RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
45 Crs2Op(RCP<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Op) {
46 if (Op.is_null())
47 return Teuchos::null;
48 return rcp(new CrsMatrixWrap(Op));
49}
50
51template <class Scalar,
52 class LocalOrdinal,
53 class GlobalOrdinal,
54 class Node>
55Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
56removeSmallEntries(Teuchos::RCP<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& A,
57 const typename Teuchos::ScalarTraits<Scalar>::magnitudeType threshold,
58 const bool keepDiagonal) {
59 using ATS = KokkosKernels::ArithTraits<Scalar>;
60 using impl_SC = typename ATS::val_type;
61 using impl_ATS = KokkosKernels::ArithTraits<impl_SC>;
62
63 RCP<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> crsA;
64 if (keepDiagonal) {
65 crsA = applyFilter_GID(
66 A, KOKKOS_LAMBDA(GlobalOrdinal rgid,
67 GlobalOrdinal cgid,
68 impl_SC val) {
69 return ((impl_ATS::magnitude(val) > threshold) || (rgid == cgid));
70 });
71
72 } else {
73 crsA = applyFilter_vals(
74 A, KOKKOS_LAMBDA(impl_SC val) {
75 return (impl_ATS::magnitude(val) > threshold);
76 });
77 }
78 return rcp(new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(crsA));
79}
80
81template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
82RCP<Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
84 GetThresholdedMatrix(const RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& Ain, const typename Teuchos::ScalarTraits<Scalar>::magnitudeType threshold, const bool keepDiagonal, const GlobalOrdinal expectedNNZperRow) {
85 auto crsWrap = rcp_dynamic_cast<CrsMatrixWrap>(Ain);
86 if (!crsWrap.is_null()) {
87 auto crsMat = crsWrap->getCrsMatrix();
88 auto filteredMat = removeSmallEntries(crsMat, threshold, keepDiagonal);
89 return rcp_static_cast<CrsMatrixWrap>(filteredMat);
90 }
91
92 RCP<const Map> rowmap = Ain->getRowMap();
93 RCP<const Map> colmap = Ain->getColMap();
94 RCP<CrsMatrixWrap> Aout = rcp(new CrsMatrixWrap(rowmap, expectedNNZperRow <= 0 ? Ain->getGlobalMaxNumRowEntries() : expectedNNZperRow));
95 // loop over local rows
96 for (size_t row = 0; row < Ain->getLocalNumRows(); row++) {
97 size_t nnz = Ain->getNumEntriesInLocalRow(row);
98
99 Teuchos::ArrayView<const LocalOrdinal> indices;
100 Teuchos::ArrayView<const Scalar> vals;
101 Ain->getLocalRowView(row, indices, vals);
102
103 TEUCHOS_TEST_FOR_EXCEPTION(Teuchos::as<size_t>(indices.size()) != nnz, Exceptions::RuntimeError, "MueLu::ThresholdAFilterFactory::Build: number of nonzeros not equal to number of indices? Error.");
104
105 Teuchos::ArrayRCP<GlobalOrdinal> indout(indices.size(), Teuchos::ScalarTraits<GlobalOrdinal>::zero());
106 Teuchos::ArrayRCP<Scalar> valout(indices.size(), Teuchos::ScalarTraits<Scalar>::zero());
107 size_t nNonzeros = 0;
108 if (keepDiagonal) {
109 GlobalOrdinal glbRow = rowmap->getGlobalElement(row);
110 LocalOrdinal lclColIdx = colmap->getLocalElement(glbRow);
111 for (size_t i = 0; i < (size_t)indices.size(); i++) {
112 if (Teuchos::ScalarTraits<Scalar>::magnitude(vals[i]) > Teuchos::ScalarTraits<Scalar>::magnitude(threshold) || indices[i] == lclColIdx) {
113 indout[nNonzeros] = colmap->getGlobalElement(indices[i]); // LID -> GID (column)
114 valout[nNonzeros] = vals[i];
115 nNonzeros++;
116 }
117 }
118 } else
119 for (size_t i = 0; i < (size_t)indices.size(); i++) {
120 if (Teuchos::ScalarTraits<Scalar>::magnitude(vals[i]) > Teuchos::ScalarTraits<Scalar>::magnitude(threshold)) {
121 indout[nNonzeros] = colmap->getGlobalElement(indices[i]); // LID -> GID (column)
122 valout[nNonzeros] = vals[i];
123 nNonzeros++;
124 }
125 }
126
127 indout.resize(nNonzeros);
128 valout.resize(nNonzeros);
129
130 Aout->insertGlobalValues(Ain->getRowMap()->getGlobalElement(row), indout.view(0, indout.size()), valout.view(0, valout.size()));
131 }
132 Aout->fillComplete(Ain->getDomainMap(), Ain->getRangeMap());
133
134 return Aout;
135}
136
137template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
138RCP<Xpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>>
140 GetThresholdedGraph(const RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& A, const Magnitude threshold, const GlobalOrdinal expectedNNZperRow) {
141 using STS = Teuchos::ScalarTraits<Scalar>;
142 RCP<CrsGraph> sparsityPattern = CrsGraphFactory::Build(A->getRowMap(), expectedNNZperRow <= 0 ? A->getGlobalMaxNumRowEntries() : expectedNNZperRow);
143
144 RCP<Vector> diag = GetMatrixOverlappedDiagonal(*A);
145 ArrayRCP<const Scalar> D = diag->getData(0);
146
147 for (size_t row = 0; row < A->getLocalNumRows(); row++) {
148 ArrayView<const LocalOrdinal> indices;
149 ArrayView<const Scalar> vals;
150 A->getLocalRowView(row, indices, vals);
151
152 GlobalOrdinal globalRow = A->getRowMap()->getGlobalElement(row);
153 LocalOrdinal col = A->getColMap()->getLocalElement(globalRow);
154
155 const Scalar Dk = STS::magnitude(D[col]) > 0.0 ? STS::magnitude(D[col]) : 1.0;
156 Array<GlobalOrdinal> indicesNew;
157
158 for (size_t i = 0; i < size_t(indices.size()); i++)
159 // keep diagonal per default
160 if (col == indices[i] || STS::magnitude(STS::squareroot(Dk) * vals[i] * STS::squareroot(Dk)) > STS::magnitude(threshold))
161 indicesNew.append(A->getColMap()->getGlobalElement(indices[i]));
162
163 sparsityPattern->insertGlobalIndices(globalRow, ArrayView<const GlobalOrdinal>(indicesNew.data(), indicesNew.length()));
164 }
165 sparsityPattern->fillComplete();
166
167 return sparsityPattern;
168}
169
170template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
171Teuchos::ArrayRCP<Scalar>
173 GetMatrixDiagonal_arcp(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A) {
174 size_t numRows = A.getRowMap()->getLocalNumElements();
175 Teuchos::ArrayRCP<Scalar> diag(numRows);
176 Teuchos::ArrayView<const LocalOrdinal> cols;
177 Teuchos::ArrayView<const Scalar> vals;
178 for (size_t i = 0; i < numRows; ++i) {
179 A.getLocalRowView(i, cols, vals);
181 for (; j < cols.size(); ++j) {
182 if (Teuchos::as<size_t>(cols[j]) == i) {
183 diag[i] = vals[j];
184 break;
185 }
186 }
187 if (j == cols.size()) {
188 // Diagonal entry is absent
189 diag[i] = Teuchos::ScalarTraits<Scalar>::zero();
190 }
191 }
192 return diag;
193}
194
195template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
196Teuchos::RCP<Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
198 GetMatrixDiagonal(const Matrix& A) {
199 const auto rowMap = A.getRowMap();
200 auto diag = Xpetra::VectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(rowMap, false);
201
202 const CrsMatrixWrap* crsOp = dynamic_cast<const CrsMatrixWrap*>(&A);
203 if ((crsOp != NULL) && (rowMap->lib() == Xpetra::UseTpetra)) {
204 using device_type = typename CrsGraph::device_type;
205 Kokkos::View<size_t*, device_type> offsets("offsets", rowMap->getLocalNumElements());
206 crsOp->getCrsGraph()->getLocalDiagOffsets(offsets);
207 crsOp->getCrsMatrix()->getLocalDiagCopy(*diag, offsets);
208 } else {
209 A.getLocalDiagCopy(*diag);
210 }
211
212 return diag;
213}
215// template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
216// RCP<Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
217// UtilitiesBase<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
218// GetMatrixDiagonalInverse(const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> & A, Magnitude tol, Scalar valReplacement) {
219// Teuchos::TimeMonitor MM = *Teuchos::TimeMonitor::getNewTimer("UtilitiesBase::GetMatrixDiagonalInverse");
221// RCP<const Map> rowMap = A.getRowMap();
222// RCP<Vector> diag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap,true);
223
224// A.getLocalDiagCopy(*diag);
225
226// RCP<Vector> inv = MueLu::UtilitiesBase<Scalar,LocalOrdinal,GlobalOrdinal,Node>::GetInverse(diag, tol, valReplacement);
227
228// return inv;
229// }
230
231template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
232Teuchos::RCP<Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
234 GetMatrixDiagonalInverse(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
235 typename Teuchos::ScalarTraits<Scalar>::magnitudeType tol,
236 Scalar valReplacement,
237 const bool doLumped) {
238 Teuchos::TimeMonitor MM = *Teuchos::TimeMonitor::getNewTimer("Utilities::GetMatrixDiagonalInverse");
239
240 RCP<const BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(rcpFromRef(A));
241 if (!bA.is_null()) {
242 RCP<const Map> rowMap = A.getRowMap();
243 RCP<Vector> diag = Xpetra::VectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(rowMap, false);
244 A.getLocalDiagCopy(*diag);
245 RCP<Vector> inv = MueLu::UtilitiesBase<Scalar, LocalOrdinal, GlobalOrdinal, Node>::GetInverse(diag, tol, valReplacement);
246 return inv;
247 }
248
249 // Some useful type definitions
250 using local_matrix_type = typename Matrix::local_matrix_device_type;
251 // using local_graph_type = typename local_matrix_type::staticcrsgraph_type;
252 using value_type = typename local_matrix_type::value_type;
253 using values_type = typename local_matrix_type::values_type;
254 using scalar_type = typename values_type::non_const_value_type;
255 using ordinal_type = typename local_matrix_type::ordinal_type;
256 using execution_space = typename local_matrix_type::execution_space;
257 // using memory_space = typename local_matrix_type::memory_space;
258 // Be careful with this one, if using KokkosKernels::ArithTraits<Scalar>
259 // you are likely to run into errors when handling std::complex<>
260 // a good way to work around that is to use the following:
261 // using KAT = KokkosKernels::ArithTraits<KokkosKernels::ArithTraits<Scalar>::val_type> >
262 // here we have: value_type = KokkosKernels::ArithTraits<Scalar>::val_type
263 using KAT = KokkosKernels::ArithTraits<value_type>;
265 // Get/Create distributed objects
266 RCP<const Map> rowMap = A.getRowMap();
267 RCP<Vector> diag = VectorFactory::Build(rowMap, false);
268
269 // Now generate local objects
270 local_matrix_type localMatrix = A.getLocalMatrixDevice();
271 auto diagVals = diag->getLocalViewDevice(Tpetra::Access::ReadWrite);
272
273 ordinal_type numRows = localMatrix.graph.numRows();
274
275 scalar_type valReplacement_dev = valReplacement;
276
277 // Note: 2019-11-21, LBV
278 // This could be implemented with a TeamPolicy over the rows
279 // and a TeamVectorRange over the entries in a row if performance
280 // becomes more important here.
281 if (!doLumped)
282 Kokkos::parallel_for(
283 "Utilities::GetMatrixDiagonalInverse",
284 Kokkos::RangePolicy<ordinal_type, execution_space>(0, numRows),
285 KOKKOS_LAMBDA(const ordinal_type rowIdx) {
286 bool foundDiagEntry = false;
287 auto myRow = localMatrix.rowConst(rowIdx);
288 for (ordinal_type entryIdx = 0; entryIdx < myRow.length; ++entryIdx) {
289 if (myRow.colidx(entryIdx) == rowIdx) {
290 foundDiagEntry = true;
291 if (KAT::magnitude(myRow.value(entryIdx)) > KAT::magnitude(tol)) {
292 diagVals(rowIdx, 0) = KAT::one() / myRow.value(entryIdx);
293 } else {
294 diagVals(rowIdx, 0) = valReplacement_dev;
295 }
296 break;
297 }
298 }
299
300 if (!foundDiagEntry) {
301 diagVals(rowIdx, 0) = KAT::zero();
302 }
303 });
304 else
305 Kokkos::parallel_for(
306 "Utilities::GetMatrixDiagonalInverse",
307 Kokkos::RangePolicy<ordinal_type, execution_space>(0, numRows),
308 KOKKOS_LAMBDA(const ordinal_type rowIdx) {
309 auto myRow = localMatrix.rowConst(rowIdx);
310 for (ordinal_type entryIdx = 0; entryIdx < myRow.length; ++entryIdx) {
311 diagVals(rowIdx, 0) += KAT::magnitude(myRow.value(entryIdx));
313 if (KAT::magnitude(diagVals(rowIdx, 0)) > KAT::magnitude(tol))
314 diagVals(rowIdx, 0) = KAT::one() / diagVals(rowIdx, 0);
315 else
316 diagVals(rowIdx, 0) = valReplacement_dev;
317 });
318
319 return diag;
320}
321
322template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
323Teuchos::RCP<Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
325 GetLumpedMatrixDiagonal(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> const& A, const bool doReciprocal,
326 Magnitude tol,
327 Scalar valReplacement,
328 const bool replaceSingleEntryRowWithZero,
329 const bool useAverageAbsDiagVal) {
330 typedef Teuchos::ScalarTraits<Scalar> TST;
331
332 RCP<Vector> diag = Teuchos::null;
333 const Scalar zero = TST::zero();
334 const Scalar one = TST::one();
335 const Scalar two = one + one;
336
337 Teuchos::RCP<const Matrix> rcpA = Teuchos::rcpFromRef(A);
338
339 RCP<const Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> bA =
340 Teuchos::rcp_dynamic_cast<const Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(rcpA);
341 if (bA == Teuchos::null) {
342 RCP<const Map> rowMap = rcpA->getRowMap();
343 diag = Xpetra::VectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(rowMap, false);
344
345 if (rowMap->lib() == Xpetra::UnderlyingLib::UseTpetra) {
346 Teuchos::TimeMonitor MM = *Teuchos::TimeMonitor::getNewTimer("UtilitiesBase::GetLumpedMatrixDiagonal (Kokkos implementation)");
347 // Implement using Kokkos
348 using local_vector_type = typename Vector::dual_view_type::t_dev_um;
349 using local_matrix_type = typename Matrix::local_matrix_device_type;
350 using execution_space = typename local_vector_type::execution_space;
351 // using rowmap_type = typename local_matrix_type::row_map_type;
352 // using entries_type = typename local_matrix_type::index_type;
353 using values_type = typename local_matrix_type::values_type;
354 using scalar_type = typename values_type::non_const_value_type;
355 using mag_type = typename KokkosKernels::ArithTraits<scalar_type>::mag_type;
356 using KAT_S = typename KokkosKernels::ArithTraits<scalar_type>;
357 using KAT_M = typename KokkosKernels::ArithTraits<mag_type>;
358 using size_type = typename local_matrix_type::non_const_size_type;
359
360 local_vector_type diag_dev = diag->getLocalViewDevice(Tpetra::Access::OverwriteAll);
361 local_matrix_type local_mat_dev = rcpA->getLocalMatrixDevice();
362 Kokkos::RangePolicy<execution_space, int> my_policy(0, static_cast<int>(diag_dev.extent(0)));
363 scalar_type valReplacement_dev = valReplacement;
364
365 if (doReciprocal) {
366 Kokkos::View<int*, execution_space> nnzPerRow("nnz per rows", diag_dev.extent(0));
367 Kokkos::View<scalar_type*, execution_space> regSum("regSum", diag_dev.extent(0));
368 Kokkos::View<mag_type, execution_space> avgAbsDiagVal_dev("avgAbsDiagVal");
369 Kokkos::View<int, execution_space> numDiagsEqualToOne_dev("numDiagsEqualToOne");
370
371 {
372 Teuchos::TimeMonitor MMM = *Teuchos::TimeMonitor::getNewTimer("GetLumpedMatrixDiagonal: parallel_for (doReciprocal)");
373 Kokkos::parallel_for(
374 "GetLumpedMatrixDiagonal", my_policy,
375 KOKKOS_LAMBDA(const int rowIdx) {
376 diag_dev(rowIdx, 0) = KAT_S::zero();
377 for (size_type entryIdx = local_mat_dev.graph.row_map(rowIdx);
378 entryIdx < local_mat_dev.graph.row_map(rowIdx + 1);
379 ++entryIdx) {
380 regSum(rowIdx) += local_mat_dev.values(entryIdx);
381 if (KAT_M::zero() < KAT_S::abs(local_mat_dev.values(entryIdx))) {
382 ++nnzPerRow(rowIdx);
383 }
384 diag_dev(rowIdx, 0) += KAT_S::abs(local_mat_dev.values(entryIdx));
385 if (rowIdx == local_mat_dev.graph.entries(entryIdx)) {
386 Kokkos::atomic_add(&avgAbsDiagVal_dev(), KAT_S::abs(local_mat_dev.values(entryIdx)));
387 }
388 }
389
390 if (nnzPerRow(rowIdx) == 1 && KAT_S::magnitude(diag_dev(rowIdx, 0)) == KAT_M::one()) {
391 Kokkos::atomic_add(&numDiagsEqualToOne_dev(), 1);
392 }
393 });
394 }
395 if (useAverageAbsDiagVal) {
396 Teuchos::TimeMonitor MMM = *Teuchos::TimeMonitor::getNewTimer("GetLumpedMatrixDiagonal: useAverageAbsDiagVal");
397 typename Kokkos::View<mag_type, execution_space>::host_mirror_type avgAbsDiagVal = Kokkos::create_mirror_view(avgAbsDiagVal_dev);
398 Kokkos::deep_copy(avgAbsDiagVal, avgAbsDiagVal_dev);
399 int numDiagsEqualToOne;
400 Kokkos::deep_copy(numDiagsEqualToOne, numDiagsEqualToOne_dev);
401
402 tol = TST::magnitude(100 * Teuchos::ScalarTraits<Scalar>::eps()) * (avgAbsDiagVal() - numDiagsEqualToOne) / (rowMap->getLocalNumElements() - numDiagsEqualToOne);
403 }
404
406 Teuchos::TimeMonitor MMM = *Teuchos::TimeMonitor::getNewTimer("ComputeLumpedDiagonalInverse: parallel_for (doReciprocal)");
407 Kokkos::parallel_for(
408 "ComputeLumpedDiagonalInverse", my_policy,
409 KOKKOS_LAMBDA(const int rowIdx) {
410 if (replaceSingleEntryRowWithZero && nnzPerRow(rowIdx) <= 1) {
411 diag_dev(rowIdx, 0) = KAT_S::zero();
412 } else if ((diag_dev(rowIdx, 0) != KAT_S::zero()) && (KAT_S::magnitude(diag_dev(rowIdx, 0)) < KAT_S::magnitude(2 * regSum(rowIdx)))) {
413 diag_dev(rowIdx, 0) = KAT_S::one() / KAT_S::magnitude(2 * regSum(rowIdx));
414 } else {
415 if (KAT_S::magnitude(diag_dev(rowIdx, 0)) > tol) {
416 diag_dev(rowIdx, 0) = KAT_S::one() / diag_dev(rowIdx, 0);
417 } else {
418 diag_dev(rowIdx, 0) = valReplacement_dev;
419 }
420 }
421 });
422 }
423
424 } else {
425 Teuchos::TimeMonitor MMM = *Teuchos::TimeMonitor::getNewTimer("GetLumpedMatrixDiagonal: parallel_for");
426 Kokkos::parallel_for(
427 "GetLumpedMatrixDiagonal", my_policy,
428 KOKKOS_LAMBDA(const int rowIdx) {
429 diag_dev(rowIdx, 0) = KAT_S::zero();
430 for (size_type entryIdx = local_mat_dev.graph.row_map(rowIdx);
431 entryIdx < local_mat_dev.graph.row_map(rowIdx + 1);
432 ++entryIdx) {
433 diag_dev(rowIdx, 0) += KAT_S::magnitude(local_mat_dev.values(entryIdx));
435 });
436 }
437 } else {
438 // Implement using Teuchos
439 Teuchos::TimeMonitor MMM = *Teuchos::TimeMonitor::getNewTimer("UtilitiesBase: GetLumpedMatrixDiagonal: (Teuchos implementation)");
440 ArrayRCP<Scalar> diagVals = diag->getDataNonConst(0);
441 Teuchos::Array<Scalar> regSum(diag->getLocalLength());
442 Teuchos::ArrayView<const LocalOrdinal> cols;
443 Teuchos::ArrayView<const Scalar> vals;
444
445 std::vector<int> nnzPerRow(rowMap->getLocalNumElements());
446
447 // FIXME 2021-10-22 JHU If this is called with doReciprocal=false, what should the correct behavior be? Currently,
448 // FIXME 2021-10-22 JHU the diagonal entry is set to be the sum of the absolute values of the row entries.
449
450 const Magnitude zeroMagn = TST::magnitude(zero);
451 Magnitude avgAbsDiagVal = TST::magnitude(zero);
452 int numDiagsEqualToOne = 0;
453 for (size_t i = 0; i < rowMap->getLocalNumElements(); ++i) {
454 nnzPerRow[i] = 0;
455 rcpA->getLocalRowView(i, cols, vals);
456 diagVals[i] = zero;
457 for (LocalOrdinal j = 0; j < cols.size(); ++j) {
458 regSum[i] += vals[j];
459 const Magnitude rowEntryMagn = TST::magnitude(vals[j]);
460 if (rowEntryMagn > zeroMagn)
461 nnzPerRow[i]++;
462 diagVals[i] += rowEntryMagn;
463 if (static_cast<size_t>(cols[j]) == i)
464 avgAbsDiagVal += rowEntryMagn;
465 }
466 if (nnzPerRow[i] == 1 && TST::magnitude(diagVals[i]) == 1.)
467 numDiagsEqualToOne++;
468 }
469 if (useAverageAbsDiagVal)
470 tol = TST::magnitude(100 * Teuchos::ScalarTraits<Scalar>::eps()) * (avgAbsDiagVal - numDiagsEqualToOne) / (rowMap->getLocalNumElements() - numDiagsEqualToOne);
471 if (doReciprocal) {
472 for (size_t i = 0; i < rowMap->getLocalNumElements(); ++i) {
473 if (replaceSingleEntryRowWithZero && nnzPerRow[i] <= static_cast<int>(1))
474 diagVals[i] = zero;
475 else if ((diagVals[i] != zero) && (TST::magnitude(diagVals[i]) < TST::magnitude(two * regSum[i])))
476 diagVals[i] = one / TST::magnitude((two * regSum[i]));
477 else {
478 if (TST::magnitude(diagVals[i]) > tol)
479 diagVals[i] = one / diagVals[i];
480 else {
481 diagVals[i] = valReplacement;
482 }
483 }
484 }
485 }
486 }
487 } else {
488 TEUCHOS_TEST_FOR_EXCEPTION(doReciprocal, Xpetra::Exceptions::RuntimeError,
489 "UtilitiesBase::GetLumpedMatrixDiagonal(): extracting reciprocal of diagonal of a blocked matrix is not supported");
490 diag = Xpetra::VectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(bA->getRangeMapExtractor()->getFullMap(), true);
491
492 for (size_t row = 0; row < bA->Rows(); ++row) {
493 for (size_t col = 0; col < bA->Cols(); ++col) {
494 if (!bA->getMatrix(row, col).is_null()) {
495 // if we are in Thyra mode, but the block (row,row) is again a blocked operator, we have to use (pseudo) Xpetra-style GIDs with offset!
496 bool bThyraMode = bA->getRangeMapExtractor()->getThyraMode() && (Teuchos::rcp_dynamic_cast<Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(bA->getMatrix(row, col)) == Teuchos::null);
497 RCP<Vector> ddtemp = bA->getRangeMapExtractor()->ExtractVector(diag, row, bThyraMode);
498 RCP<const Vector> dd = GetLumpedMatrixDiagonal(*(bA->getMatrix(row, col)));
499 ddtemp->update(Teuchos::as<Scalar>(1.0), *dd, Teuchos::as<Scalar>(1.0));
500 bA->getRangeMapExtractor()->InsertVector(ddtemp, row, diag, bThyraMode);
501 }
502 }
503 }
504 }
505
506 return diag;
507}
508
509template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
510Teuchos::RCP<Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
512 GetMatrixMaxMinusOffDiagonal(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A) {
513 // Get/Create distributed objects
514 RCP<const Map> rowMap = A.getRowMap();
515 auto diag = Xpetra::VectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(rowMap, false);
516
517 // Implement using Kokkos
518 using local_vector_type = typename Vector::dual_view_type::t_dev_um;
519 using local_matrix_type = typename Matrix::local_matrix_device_type;
520 using execution_space = typename local_vector_type::execution_space;
521 using values_type = typename local_matrix_type::values_type;
522 using scalar_type = typename values_type::non_const_value_type;
523 using KAT_S = typename KokkosKernels::ArithTraits<scalar_type>;
524
525 auto diag_dev = diag->getLocalViewDevice(Tpetra::Access::OverwriteAll);
526 auto local_mat_dev = A.getLocalMatrixDevice();
527 Kokkos::RangePolicy<execution_space, int> my_policy(0, static_cast<int>(diag_dev.extent(0)));
528
529 Kokkos::parallel_for(
530 "GetMatrixMaxMinusOffDiagonal", my_policy,
531 KOKKOS_LAMBDA(const LocalOrdinal rowIdx) {
532 auto mymax = KAT_S::zero();
533 auto row = local_mat_dev.rowConst(rowIdx);
534 for (LocalOrdinal entryIdx = 0; entryIdx < row.length; ++entryIdx) {
535 if (rowIdx != row.colidx(entryIdx)) {
536 if (KAT_S::real(mymax) < -KAT_S::real(row.value(entryIdx)))
537 mymax = -KAT_S::real(row.value(entryIdx));
538 }
539 }
540 diag_dev(rowIdx, 0) = mymax;
541 });
542
543 return diag;
544}
545
546template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
547Teuchos::RCP<Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
549 GetMatrixMaxMinusOffDiagonal(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A, const Xpetra::Vector<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>& BlockNumber) {
550 TEUCHOS_TEST_FOR_EXCEPTION(!A.getColMap()->isSameAs(*BlockNumber.getMap()), std::runtime_error, "GetMatrixMaxMinusOffDiagonal: BlockNumber must match's A's column map.");
551
552 // Get/Create distributed objects
553 RCP<const Map> rowMap = A.getRowMap();
554 auto diag = Xpetra::VectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(rowMap, false);
555
556 // Implement using Kokkos
557 using local_vector_type = typename Vector::dual_view_type::t_dev_um;
558 using local_matrix_type = typename Matrix::local_matrix_device_type;
559 using execution_space = typename local_vector_type::execution_space;
560 using values_type = typename local_matrix_type::values_type;
561 using scalar_type = typename values_type::non_const_value_type;
562 using KAT_S = typename KokkosKernels::ArithTraits<scalar_type>;
563
564 auto diag_dev = diag->getLocalViewDevice(Tpetra::Access::OverwriteAll);
565 auto local_mat_dev = A.getLocalMatrixDevice();
566 auto local_block_dev = BlockNumber.getLocalViewDevice(Tpetra::Access::ReadOnly);
567 Kokkos::RangePolicy<execution_space, int> my_policy(0, static_cast<int>(diag_dev.extent(0)));
568
569 Kokkos::parallel_for(
570 "GetMatrixMaxMinusOffDiagonal", my_policy,
571 KOKKOS_LAMBDA(const LocalOrdinal rowIdx) {
572 auto mymax = KAT_S::zero();
573 auto row = local_mat_dev.row(rowIdx);
574 for (LocalOrdinal entryIdx = 0; entryIdx < row.length; ++entryIdx) {
575 if ((rowIdx != row.colidx(entryIdx)) && (local_block_dev(rowIdx, 0) == local_block_dev(row.colidx(entryIdx), 0))) {
576 if (KAT_S::real(mymax) < -KAT_S::real(row.value(entryIdx)))
577 mymax = -KAT_S::real(row.value(entryIdx));
578 }
579 }
580 diag_dev(rowIdx, 0) = mymax;
581 });
582
583 return diag;
584}
585
586template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
587Teuchos::RCP<Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
589 GetInverse(Teuchos::RCP<const Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v, typename Teuchos::ScalarTraits<Scalar>::magnitudeType tol, Scalar valReplacement) {
590 RCP<Vector> ret = Xpetra::VectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(v->getMap(), true);
591
592 // check whether input vector "v" is a BlockedVector
593 RCP<const BlockedVector> bv = Teuchos::rcp_dynamic_cast<const BlockedVector>(v);
594 if (bv.is_null() == false) {
595 RCP<BlockedVector> bret = Teuchos::rcp_dynamic_cast<BlockedVector>(ret);
596 TEUCHOS_TEST_FOR_EXCEPTION(bret.is_null() == true, MueLu::Exceptions::RuntimeError, "MueLu::UtilitiesBase::GetInverse: return vector should be of type BlockedVector");
597 RCP<const BlockedMap> bmap = bv->getBlockedMap();
598 for (size_t r = 0; r < bmap->getNumMaps(); ++r) {
599 RCP<const MultiVector> submvec = bv->getMultiVector(r, bmap->getThyraMode());
600 RCP<const Vector> subvec = submvec->getVector(0);
601 RCP<Vector> subvecinf = MueLu::UtilitiesBase<Scalar, LocalOrdinal, GlobalOrdinal, Node>::GetInverse(subvec, tol, valReplacement);
602 bret->setMultiVector(r, subvecinf, bmap->getThyraMode());
603 }
604 return ret;
605 }
606
607 // v is an {Epetra,Tpetra}Vector: work with the underlying raw data
608 ArrayRCP<Scalar> retVals = ret->getDataNonConst(0);
609 ArrayRCP<const Scalar> inputVals = v->getData(0);
610 for (size_t i = 0; i < v->getMap()->getLocalNumElements(); ++i) {
611 if (Teuchos::ScalarTraits<Scalar>::magnitude(inputVals[i]) > tol)
612 retVals[i] = Teuchos::ScalarTraits<Scalar>::one() / inputVals[i];
613 else
614 retVals[i] = valReplacement;
615 }
616 return ret;
617}
618
619// template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
620// RCP<Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
621// UtilitiesBase<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
622// GetMatrixOverlappedDiagonal(const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> & A) {
623// RCP<const Map> rowMap = A.getRowMap(), colMap = A.getColMap();
624
625// // Undo block map (if we have one)
626// RCP<const BlockedMap> browMap = Teuchos::rcp_dynamic_cast<const BlockedMap>(rowMap);
627// if(!browMap.is_null()) rowMap = browMap->getMap();
628
629// RCP<Vector> localDiag = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap);
630// try {
631// const CrsMatrixWrap* crsOp = dynamic_cast<const CrsMatrixWrap*>(&A);
632// if (crsOp == NULL) {
633// throw Exceptions::RuntimeError("cast to CrsMatrixWrap failed");
634// }
635// Teuchos::ArrayRCP<size_t> offsets;
636// crsOp->getLocalDiagOffsets(offsets);
637// crsOp->getLocalDiagCopy(*localDiag,offsets());
638// }
639// catch (...) {
640// ArrayRCP<Scalar> localDiagVals = localDiag->getDataNonConst(0);
641// Teuchos::ArrayRCP<Scalar> diagVals = GetMatrixDiagonal(A);
642// for (LocalOrdinal i = 0; i < localDiagVals.size(); i++)
643// localDiagVals[i] = diagVals[i];
644// localDiagVals = diagVals = null;
645// }
646
647// RCP<Vector> diagonal = Xpetra::VectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(colMap);
648// RCP< const Xpetra::Import<LocalOrdinal,GlobalOrdinal,Node> > importer;
649// importer = A.getCrsGraph()->getImporter();
650// if (importer == Teuchos::null) {
651// importer = Xpetra::ImportFactory<LocalOrdinal,GlobalOrdinal,Node>::Build(rowMap, colMap);
652// }
653// diagonal->doImport(*localDiag, *(importer), Xpetra::INSERT);
654// return diagonal;
655// }
656
657template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
658RCP<Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
660 GetMatrixOverlappedDiagonal(const Matrix& A) {
661 RCP<const Map> rowMap = A.getRowMap(), colMap = A.getColMap();
662 RCP<Vector> localDiag = GetMatrixDiagonal(A);
663 RCP<const Import> importer = A.getCrsGraph()->getImporter();
664 if (importer.is_null() && !rowMap->isSameAs(*colMap)) {
665 importer = ImportFactory::Build(rowMap, colMap);
666 }
667 if (!importer.is_null()) {
668 RCP<Vector> diagonal = VectorFactory::Build(colMap, false);
669 diagonal->doImport(*localDiag, *(importer), Xpetra::INSERT);
670 return diagonal;
671 } else
672 return localDiag;
673}
674
675template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
676RCP<Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
678 GetMatrixOverlappedDeletedRowsum(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A) {
679 using STS = typename Teuchos::ScalarTraits<SC>;
680
681 // Undo block map (if we have one)
682 RCP<const Map> rowMap = A.getRowMap(), colMap = A.getColMap();
683 RCP<const BlockedMap> browMap = Teuchos::rcp_dynamic_cast<const BlockedMap>(rowMap);
684 if (!browMap.is_null()) rowMap = browMap->getMap();
685
686 RCP<Vector> local = Xpetra::VectorFactory<SC, LO, GO, Node>::Build(rowMap);
687 RCP<Vector> ghosted = Xpetra::VectorFactory<SC, LO, GO, Node>::Build(colMap, true);
688 ArrayRCP<SC> localVals = local->getDataNonConst(0);
689
690 for (LO row = 0; row < static_cast<LO>(A.getRowMap()->getLocalNumElements()); ++row) {
691 size_t nnz = A.getNumEntriesInLocalRow(row);
692 ArrayView<const LO> indices;
693 ArrayView<const SC> vals;
694 A.getLocalRowView(row, indices, vals);
695
696 SC si = STS::zero();
697
698 for (LO colID = 0; colID < static_cast<LO>(nnz); colID++) {
699 if (indices[colID] != row) {
700 si += vals[colID];
701 }
702 }
703 localVals[row] = si;
704 }
705
706 RCP<const Xpetra::Import<LO, GO, Node>> importer;
707 importer = A.getCrsGraph()->getImporter();
708 if (importer == Teuchos::null) {
709 importer = Xpetra::ImportFactory<LO, GO, Node>::Build(rowMap, colMap);
710 }
711 ghosted->doImport(*local, *(importer), Xpetra::INSERT);
712 return ghosted;
713}
714
715template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
716RCP<Xpetra::Vector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node>>
718 GetMatrixOverlappedAbsDeletedRowsum(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A) {
719 RCP<const Map> rowMap = A.getRowMap(), colMap = A.getColMap();
720 using STS = typename Teuchos::ScalarTraits<Scalar>;
721 using MTS = typename Teuchos::ScalarTraits<Magnitude>;
722 using MT = Magnitude;
723 using RealValuedVector = Xpetra::Vector<MT, LO, GO, Node>;
724
725 // Undo block map (if we have one)
726 RCP<const BlockedMap> browMap = Teuchos::rcp_dynamic_cast<const BlockedMap>(rowMap);
727 if (!browMap.is_null()) rowMap = browMap->getMap();
728
729 RCP<RealValuedVector> local = Xpetra::VectorFactory<MT, LO, GO, Node>::Build(rowMap);
730 RCP<RealValuedVector> ghosted = Xpetra::VectorFactory<MT, LO, GO, Node>::Build(colMap, true);
731 ArrayRCP<MT> localVals = local->getDataNonConst(0);
732
733 for (LO rowIdx = 0; rowIdx < static_cast<LO>(A.getRowMap()->getLocalNumElements()); ++rowIdx) {
734 size_t nnz = A.getNumEntriesInLocalRow(rowIdx);
735 ArrayView<const LO> indices;
736 ArrayView<const SC> vals;
737 A.getLocalRowView(rowIdx, indices, vals);
738
739 MT si = MTS::zero();
740
741 for (LO colID = 0; colID < static_cast<LO>(nnz); ++colID) {
742 if (indices[colID] != rowIdx) {
743 si += STS::magnitude(vals[colID]);
744 }
745 }
746 localVals[rowIdx] = si;
747 }
748
749 RCP<const Xpetra::Import<LO, GO, Node>> importer;
750 importer = A.getCrsGraph()->getImporter();
751 if (importer == Teuchos::null) {
752 importer = Xpetra::ImportFactory<LO, GO, Node>::Build(rowMap, colMap);
753 }
754 ghosted->doImport(*local, *(importer), Xpetra::INSERT);
755 return ghosted;
756}
757
758template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
761 CountNegativeDiagonalEntries(const Matrix& A) {
762 using local_matrix_type = typename Matrix::local_matrix_device_type;
763 using execution_space = typename local_matrix_type::execution_space;
764 using KAT_S = typename KokkosKernels::ArithTraits<typename local_matrix_type::value_type>;
765
766 auto local_mat_dev = A.getLocalMatrixDevice();
767 Kokkos::RangePolicy<execution_space, int> my_policy(0, static_cast<int>(local_mat_dev.numRows()));
768 GlobalOrdinal count_l = 0, count_g = 0;
769
770 Kokkos::parallel_reduce(
771 "CountNegativeDiagonalEntries", my_policy,
772 KOKKOS_LAMBDA(const LocalOrdinal rowIdx, GlobalOrdinal& sum) {
773 auto row = local_mat_dev.row(rowIdx);
774 for (LocalOrdinal entryIdx = 0; entryIdx < row.length; ++entryIdx) {
775 if (rowIdx == row.colidx(entryIdx) && KAT_S::real(row.value(entryIdx)) < 0)
776 sum++;
777 }
778 },
779 count_l);
780
781 MueLu_sumAll(A.getRowMap()->getComm(), count_l, count_g);
782 return count_g;
783}
784
785template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
786Teuchos::Array<typename Teuchos::ScalarTraits<Scalar>::magnitudeType>
788 ResidualNorm(const Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op, const MultiVector& X, const MultiVector& RHS) {
789 TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != RHS.getNumVectors(), Exceptions::RuntimeError, "Number of solution vectors != number of right-hand sides")
790 const size_t numVecs = X.getNumVectors();
791 RCP<MultiVector> RES = Residual(Op, X, RHS);
792 Teuchos::Array<Magnitude> norms(numVecs);
793 RES->norm2(norms);
794 return norms;
795}
796
797template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
798Teuchos::Array<typename Teuchos::ScalarTraits<Scalar>::magnitudeType>
800 ResidualNorm(const Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op, const MultiVector& X, const MultiVector& RHS, MultiVector& Resid) {
801 TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != RHS.getNumVectors(), Exceptions::RuntimeError, "Number of solution vectors != number of right-hand sides")
802 const size_t numVecs = X.getNumVectors();
803 Residual(Op, X, RHS, Resid);
804 Teuchos::Array<Magnitude> norms(numVecs);
805 Resid.norm2(norms);
806 return norms;
807}
808
809template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
810RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
812 Residual(const Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op, const MultiVector& X, const MultiVector& RHS) {
813 TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != RHS.getNumVectors(), Exceptions::RuntimeError, "Number of solution vectors != number of right-hand sides")
814 const size_t numVecs = X.getNumVectors();
815 // TODO Op.getRangeMap should return a BlockedMap if it is a BlockedCrsOperator
816 RCP<MultiVector> RES = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(RHS.getMap(), numVecs, false); // no need to initialize to zero
817 Op.residual(X, RHS, *RES);
818 return RES;
819}
820
821template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
823 Residual(const Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op, const MultiVector& X, const MultiVector& RHS, MultiVector& Resid) {
824 TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != RHS.getNumVectors(), Exceptions::RuntimeError, "Number of solution vectors != number of right-hand sides");
825 TEUCHOS_TEST_FOR_EXCEPTION(Resid.getNumVectors() != RHS.getNumVectors(), Exceptions::RuntimeError, "Number of residual vectors != number of right-hand sides");
826 Op.residual(X, RHS, Resid);
827}
828
829template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
830Scalar
832 PowerMethod(const Matrix& A, bool scaleByDiag,
833 LocalOrdinal niters, typename Teuchos::ScalarTraits<Scalar>::magnitudeType tolerance, bool verbose, unsigned int seed) {
834 TEUCHOS_TEST_FOR_EXCEPTION(!(A.getRangeMap()->isSameAs(*(A.getDomainMap()))), Exceptions::Incompatible,
835 "Utils::PowerMethod: operator must have domain and range maps that are equivalent.");
836
837 // power iteration
838 RCP<Vector> diagInvVec;
839 if (scaleByDiag) {
840 diagInvVec = GetMatrixDiagonalInverse(A);
841 }
842
843 Scalar lambda = PowerMethod(A, diagInvVec, niters, tolerance, verbose, seed);
844 return lambda;
845}
846
847template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
848Scalar
849UtilitiesBase<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
850 PowerMethod(const Matrix& A, const RCP<Vector>& diagInvVec,
851 LocalOrdinal niters, typename Teuchos::ScalarTraits<Scalar>::magnitudeType tolerance, bool verbose, unsigned int seed) {
852 TEUCHOS_TEST_FOR_EXCEPTION(!(A.getRangeMap()->isSameAs(*(A.getDomainMap()))), Exceptions::Incompatible,
853 "Utils::PowerMethod: operator must have domain and range maps that are equivalent.");
854
855 // Create three vectors, fill z with random numbers
856 RCP<Vector> q = Xpetra::VectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(A.getDomainMap(), true);
857 RCP<Vector> r = Xpetra::VectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(A.getRangeMap(), true);
858 RCP<Vector> z = Xpetra::VectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(A.getRangeMap(), false);
859
860 z->setSeed(seed); // seed random number generator
861 z->randomize(true); // use Xpetra implementation: -> same results for Epetra and Tpetra
862
863 Teuchos::Array<Magnitude> norms(1);
864
865 typedef Teuchos::ScalarTraits<Scalar> STS;
866
867 const Scalar zero = STS::zero(), one = STS::one();
868
869 Scalar lambda = zero;
870 Magnitude residual = STS::magnitude(zero);
871
872 // power iteration
873 for (int iter = 0; iter < niters; ++iter) {
874 z->norm2(norms); // Compute 2-norm of z
875 q->update(one / norms[0], *z, zero); // Set q = z / normz
876 A.apply(*q, *z); // Compute z = A*q
877 if (diagInvVec != Teuchos::null)
878 z->elementWiseMultiply(one, *diagInvVec, *z, zero);
879 lambda = q->dot(*z); // Approximate maximum eigenvalue: lamba = dot(q,z)
880
881 if (iter % 100 == 0 || iter + 1 == niters) {
882 r->update(1.0, *z, -lambda, *q, zero); // Compute A*q - lambda*q
883 r->norm2(norms);
884 residual = STS::magnitude(norms[0] / lambda);
885 if (verbose) {
886 std::cout << "Iter = " << iter
887 << " Lambda = " << lambda
888 << " Residual of A*q - lambda*q = " << residual
889 << std::endl;
890 }
891 }
892 if (residual < tolerance)
893 break;
894 }
895 return lambda;
896}
897
898template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
899RCP<Teuchos::FancyOStream>
901 MakeFancy(std::ostream& os) {
902 RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(os));
903 return fancy;
904}
905
906template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
907typename Teuchos::ScalarTraits<Scalar>::magnitudeType
909 Distance2(const Teuchos::Array<Teuchos::ArrayRCP<const Scalar>>& v, LocalOrdinal i0, LocalOrdinal i1) {
910 const size_t numVectors = v.size();
911
912 Scalar d = Teuchos::ScalarTraits<Scalar>::zero();
913 for (size_t j = 0; j < numVectors; j++) {
914 d += (v[j][i0] - v[j][i1]) * (v[j][i0] - v[j][i1]);
915 }
916 return Teuchos::ScalarTraits<Scalar>::magnitude(d);
917}
918
919template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
920typename Teuchos::ScalarTraits<Scalar>::magnitudeType
922 Distance2(const Teuchos::ArrayView<double>& weight, const Teuchos::Array<Teuchos::ArrayRCP<const Scalar>>& v, LocalOrdinal i0, LocalOrdinal i1) {
923 const size_t numVectors = v.size();
924 using MT = typename Teuchos::ScalarTraits<Scalar>::magnitudeType;
925
926 Scalar d = Teuchos::ScalarTraits<Scalar>::zero();
927 for (size_t j = 0; j < numVectors; j++) {
928 d += Teuchos::as<MT>(weight[j]) * (v[j][i0] - v[j][i1]) * (v[j][i0] - v[j][i1]);
929 }
930 return Teuchos::ScalarTraits<Scalar>::magnitude(d);
931}
932
933template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
934Teuchos::ArrayRCP<const bool>
936 DetectDirichletRows(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A, const typename Teuchos::ScalarTraits<Scalar>::magnitudeType& tol, bool count_twos_as_dirichlet) {
937 LocalOrdinal numRows = A.getLocalNumRows();
938 typedef Teuchos::ScalarTraits<Scalar> STS;
939 ArrayRCP<bool> boundaryNodes(numRows, true);
940 if (count_twos_as_dirichlet) {
941 for (LocalOrdinal row = 0; row < numRows; row++) {
942 ArrayView<const LocalOrdinal> indices;
943 ArrayView<const Scalar> vals;
944 A.getLocalRowView(row, indices, vals);
945 size_t nnz = A.getNumEntriesInLocalRow(row);
946 if (nnz > 2) {
947 size_t col;
948 for (col = 0; col < nnz; col++)
949 if ((indices[col] != row) && STS::magnitude(vals[col]) > tol) {
950 if (!boundaryNodes[row])
951 break;
952 boundaryNodes[row] = false;
953 }
954 if (col == nnz)
955 boundaryNodes[row] = true;
956 }
957 }
958 } else {
959 for (LocalOrdinal row = 0; row < numRows; row++) {
960 ArrayView<const LocalOrdinal> indices;
961 ArrayView<const Scalar> vals;
962 A.getLocalRowView(row, indices, vals);
963 size_t nnz = A.getNumEntriesInLocalRow(row);
964 if (nnz > 1)
965 for (size_t col = 0; col < nnz; col++)
966 if ((indices[col] != row) && STS::magnitude(vals[col]) > tol) {
967 boundaryNodes[row] = false;
968 break;
969 }
970 }
971 }
972 return boundaryNodes;
973}
974
975template <class CrsMatrix>
976KOKKOS_FORCEINLINE_FUNCTION bool isDirichletRow(typename CrsMatrix::ordinal_type rowId,
977 KokkosSparse::SparseRowViewConst<CrsMatrix>& row,
978 const typename KokkosKernels::ArithTraits<typename CrsMatrix::value_type>::magnitudeType& tol,
979 const bool count_twos_as_dirichlet) {
980 using ATS = KokkosKernels::ArithTraits<typename CrsMatrix::value_type>;
981
982 auto length = row.length;
983 bool boundaryNode = true;
984
985 if (count_twos_as_dirichlet) {
986 if (length > 2) {
987 decltype(length) colID = 0;
988 for (; colID < length; colID++)
989 if ((row.colidx(colID) != rowId) &&
990 (ATS::magnitude(row.value(colID)) > tol)) {
991 if (!boundaryNode)
992 break;
993 boundaryNode = false;
994 }
995 if (colID == length)
996 boundaryNode = true;
997 }
998 } else {
999 for (decltype(length) colID = 0; colID < length; colID++)
1000 if ((row.colidx(colID) != rowId) &&
1001 (ATS::magnitude(row.value(colID)) > tol)) {
1002 boundaryNode = false;
1003 break;
1004 }
1005 }
1006 return boundaryNode;
1007}
1008
1009template <class SC, class LO, class GO, class NO, class memory_space>
1010Kokkos::View<bool*, memory_space>
1011DetectDirichletRows_kokkos(const Xpetra::Matrix<SC, LO, GO, NO>& A,
1012 const typename Teuchos::ScalarTraits<SC>::magnitudeType& tol,
1013 const bool count_twos_as_dirichlet) {
1014 using impl_scalar_type = typename KokkosKernels::ArithTraits<SC>::val_type;
1015 using ATS = KokkosKernels::ArithTraits<impl_scalar_type>;
1016 using range_type = Kokkos::RangePolicy<LO, typename NO::execution_space>;
1017 using helpers = Xpetra::Helpers<SC, LO, GO, NO>;
1018
1019 Kokkos::View<bool*, typename NO::device_type::memory_space> boundaryNodes;
1020
1021 if (helpers::isTpetraBlockCrs(A)) {
1022 const Tpetra::BlockCrsMatrix<SC, LO, GO, NO>& Am = toTpetraBlock(A);
1023 auto b_graph = Am.getCrsGraph().getLocalGraphDevice();
1024 auto b_rowptr = Am.getCrsGraph().getLocalRowPtrsDevice();
1025 auto values = Am.getValuesDevice();
1026 LO numBlockRows = Am.getLocalNumRows();
1027 const LO stride = Am.getBlockSize() * Am.getBlockSize();
1028
1029 boundaryNodes = Kokkos::View<bool*, typename NO::device_type::memory_space>(Kokkos::ViewAllocateWithoutInitializing("boundaryNodes"), numBlockRows);
1030
1031 if (count_twos_as_dirichlet)
1032 throw Exceptions::RuntimeError("BlockCrs does not support counting twos as Dirichlet");
1033
1034 Kokkos::parallel_for(
1035 "MueLu:Utils::DetectDirichletRowsBlockCrs", range_type(0, numBlockRows),
1036 KOKKOS_LAMBDA(const LO row) {
1037 auto rowView = b_graph.rowConst(row);
1038 auto length = rowView.length;
1039 LO valstart = b_rowptr[row] * stride;
1040
1041 boundaryNodes(row) = true;
1042 decltype(length) colID = 0;
1043 for (; colID < length; colID++) {
1044 if (rowView.colidx(colID) != row) {
1045 LO current = valstart + colID * stride;
1046 for (LO k = 0; k < stride; k++) {
1047 if (ATS::magnitude(values[current + k]) > tol) {
1048 boundaryNodes(row) = false;
1049 break;
1050 }
1051 }
1052 }
1053 if (boundaryNodes(row) == false)
1054 break;
1055 }
1056 });
1057 } else {
1058 auto localMatrix = A.getLocalMatrixDevice();
1059 LO numRows = A.getLocalNumRows();
1060 boundaryNodes = Kokkos::View<bool*, typename NO::device_type::memory_space>(Kokkos::ViewAllocateWithoutInitializing("boundaryNodes"), numRows);
1061
1062 Kokkos::parallel_for(
1063 "MueLu:Utils::DetectDirichletRows", range_type(0, numRows),
1064 KOKKOS_LAMBDA(const LO row) {
1065 auto rowView = localMatrix.rowConst(row);
1066 boundaryNodes(row) = isDirichletRow(row, rowView, tol, count_twos_as_dirichlet);
1067 });
1068 }
1069 if constexpr (std::is_same<memory_space, typename NO::device_type::memory_space>::value)
1070 return boundaryNodes;
1071 else {
1072 Kokkos::View<bool*, memory_space> boundaryNodes2(Kokkos::ViewAllocateWithoutInitializing("boundaryNodes"), boundaryNodes.extent(0));
1073 Kokkos::deep_copy(boundaryNodes2, boundaryNodes);
1074 return boundaryNodes2;
1075 }
1076 // CAG: No idea why this is needed to avoid "warning: missing return statement at end of non-void function"
1077 Kokkos::View<bool*, memory_space> dummy("dummy", 0);
1078 return dummy;
1079}
1080
1081template <class SC, class LO, class GO, class NO>
1082Kokkos::View<bool*, typename NO::device_type::memory_space>
1084 DetectDirichletRows_kokkos(const Xpetra::Matrix<SC, LO, GO, NO>& A,
1085 const typename Teuchos::ScalarTraits<SC>::magnitudeType& tol,
1086 const bool count_twos_as_dirichlet) {
1087 return MueLu::DetectDirichletRows_kokkos<SC, LO, GO, NO, typename NO::device_type::memory_space>(A, tol, count_twos_as_dirichlet);
1088}
1089
1090template <class SC, class LO, class GO, class NO>
1091Kokkos::View<bool*, typename Kokkos::HostSpace>
1093 DetectDirichletRows_kokkos_host(const Xpetra::Matrix<SC, LO, GO, NO>& A,
1094 const typename Teuchos::ScalarTraits<SC>::magnitudeType& tol,
1095 const bool count_twos_as_dirichlet) {
1096 return MueLu::DetectDirichletRows_kokkos<SC, LO, GO, NO, typename Kokkos::HostSpace>(A, tol, count_twos_as_dirichlet);
1097}
1098
1099template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1100Teuchos::ArrayRCP<const bool>
1102 DetectDirichletRowsExt(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A, bool& bHasZeroDiagonal, const typename Teuchos::ScalarTraits<Scalar>::magnitudeType& tol) {
1103 // assume that there is no zero diagonal in matrix
1104 bHasZeroDiagonal = false;
1105
1106 Teuchos::RCP<Vector> diagVec = Xpetra::VectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(A.getRowMap());
1107 A.getLocalDiagCopy(*diagVec);
1108 Teuchos::ArrayRCP<const Scalar> diagVecData = diagVec->getData(0);
1109
1110 LocalOrdinal numRows = A.getLocalNumRows();
1111 typedef Teuchos::ScalarTraits<Scalar> STS;
1112 ArrayRCP<bool> boundaryNodes(numRows, false);
1113 for (LocalOrdinal row = 0; row < numRows; row++) {
1114 ArrayView<const LocalOrdinal> indices;
1115 ArrayView<const Scalar> vals;
1116 A.getLocalRowView(row, indices, vals);
1117 size_t nnz = 0; // collect nonzeros in row (excluding the diagonal)
1118 bool bHasDiag = false;
1119 for (decltype(indices.size()) col = 0; col < indices.size(); col++) {
1120 if (indices[col] != row) {
1121 if (STS::magnitude(vals[col] / STS::magnitude(sqrt(STS::magnitude(diagVecData[row]) * STS::magnitude(diagVecData[col])))) > tol) {
1122 nnz++;
1123 }
1124 } else
1125 bHasDiag = true; // found a diagonal entry
1126 }
1127 if (bHasDiag == false)
1128 bHasZeroDiagonal = true; // we found at least one row without a diagonal
1129 else if (nnz == 0)
1130 boundaryNodes[row] = true;
1131 }
1132 return boundaryNodes;
1133}
1134
1135template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1137 EnforceInitialCondition(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
1138 const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& RHS,
1139 Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& InitialGuess,
1140 const typename Teuchos::ScalarTraits<SC>::magnitudeType& tol,
1141 const bool count_twos_as_dirichlet) {
1142 using range_type = Kokkos::RangePolicy<LO, typename Node::execution_space>;
1143
1144 LocalOrdinal numRows = A.getLocalNumRows();
1145 LocalOrdinal numVectors = RHS.getNumVectors();
1146 TEUCHOS_ASSERT_EQUALITY(numVectors, Teuchos::as<LocalOrdinal>(InitialGuess.getNumVectors()));
1147#ifdef MUELU_DEBUG
1148 TEUCHOS_ASSERT(RHS.getMap()->isCompatible(InitialGuess.getMap()));
1149#endif
1150
1151 auto lclRHS = RHS.getLocalViewDevice(Tpetra::Access::ReadOnly);
1152 auto lclInitialGuess = InitialGuess.getLocalViewDevice(Tpetra::Access::ReadWrite);
1153 auto lclA = A.getLocalMatrixDevice();
1154
1155 Kokkos::parallel_for(
1156 "MueLu:Utils::EnforceInitialCondition", range_type(0, numRows),
1157 KOKKOS_LAMBDA(const LO i) {
1158 auto row = lclA.rowConst(i);
1159 if (isDirichletRow(i, row, tol, count_twos_as_dirichlet)) {
1160 for (LocalOrdinal j = 0; j < numVectors; ++j)
1161 lclInitialGuess(i, j) = lclRHS(i, j);
1162 }
1163 });
1164}
1165
1166template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1168 FindNonZeros(const Teuchos::ArrayRCP<const Scalar> vals,
1169 Teuchos::ArrayRCP<bool> nonzeros) {
1170 TEUCHOS_ASSERT(vals.size() == nonzeros.size());
1171 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType magnitudeType;
1172 const magnitudeType eps = 2.0 * Teuchos::ScalarTraits<magnitudeType>::eps();
1173 for (size_t i = 0; i < static_cast<size_t>(vals.size()); i++) {
1174 nonzeros[i] = (Teuchos::ScalarTraits<Scalar>::magnitude(vals[i]) > eps);
1175 }
1176}
1177
1178// Find Nonzeros in a device view
1179template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1181 FindNonZeros(const typename Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::dual_view_type::t_dev_const_um vals,
1182 Kokkos::View<bool*, typename Node::device_type> nonzeros) {
1183 using ATS = KokkosKernels::ArithTraits<Scalar>;
1184 using impl_ATS = KokkosKernels::ArithTraits<typename ATS::val_type>;
1185 using range_type = Kokkos::RangePolicy<LocalOrdinal, typename Node::execution_space>;
1186 TEUCHOS_ASSERT(vals.extent(0) == nonzeros.extent(0));
1187 const typename ATS::magnitudeType eps = 2.0 * impl_ATS::eps();
1188
1189 Kokkos::parallel_for(
1190 "MueLu:Maxwell1::FindNonZeros", range_type(0, vals.extent(0)),
1191 KOKKOS_LAMBDA(const size_t i) {
1192 nonzeros(i) = (impl_ATS::magnitude(vals(i, 0)) > eps);
1193 });
1194}
1195
1196template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1198 DetectDirichletColsAndDomains(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
1199 const Teuchos::ArrayRCP<bool>& dirichletRows,
1200 Teuchos::ArrayRCP<bool> dirichletCols,
1201 Teuchos::ArrayRCP<bool> dirichletDomain) {
1202 const Scalar one = Teuchos::ScalarTraits<Scalar>::one();
1203 RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> domMap = A.getDomainMap();
1204 RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> rowMap = A.getRowMap();
1205 RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> colMap = A.getColMap();
1206 TEUCHOS_ASSERT(static_cast<size_t>(dirichletRows.size()) == rowMap->getLocalNumElements());
1207 TEUCHOS_ASSERT(static_cast<size_t>(dirichletCols.size()) == colMap->getLocalNumElements());
1208 TEUCHOS_ASSERT(static_cast<size_t>(dirichletDomain.size()) == domMap->getLocalNumElements());
1209 RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> myColsToZero = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(colMap, 1, /*zeroOut=*/true);
1210 // Find all local column indices that are in Dirichlet rows, record in myColsToZero as 1.0
1211 for (size_t i = 0; i < (size_t)dirichletRows.size(); i++) {
1212 if (dirichletRows[i]) {
1213 ArrayView<const LocalOrdinal> indices;
1214 ArrayView<const Scalar> values;
1215 A.getLocalRowView(i, indices, values);
1216 for (size_t j = 0; j < static_cast<size_t>(indices.size()); j++)
1217 myColsToZero->replaceLocalValue(indices[j], 0, one);
1218 }
1219 }
1220
1221 RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> globalColsToZero;
1222 RCP<const Xpetra::Import<LocalOrdinal, GlobalOrdinal, Node>> importer = A.getCrsGraph()->getImporter();
1223 if (!importer.is_null()) {
1224 globalColsToZero = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(domMap, 1, /*zeroOut=*/true);
1225 // export to domain map
1226 globalColsToZero->doExport(*myColsToZero, *importer, Xpetra::ADD);
1227 // import to column map
1228 myColsToZero->doImport(*globalColsToZero, *importer, Xpetra::INSERT);
1229 } else
1230 globalColsToZero = myColsToZero;
1231
1232 FindNonZeros(globalColsToZero->getData(0), dirichletDomain);
1233 FindNonZeros(myColsToZero->getData(0), dirichletCols);
1234}
1235
1236template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1238 DetectDirichletColsAndDomains(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
1239 const Kokkos::View<bool*, typename Node::device_type>& dirichletRows,
1240 Kokkos::View<bool*, typename Node::device_type> dirichletCols,
1241 Kokkos::View<bool*, typename Node::device_type> dirichletDomain) {
1242 using ATS = KokkosKernels::ArithTraits<Scalar>;
1243 using impl_ATS = KokkosKernels::ArithTraits<typename ATS::val_type>;
1244 using range_type = Kokkos::RangePolicy<LocalOrdinal, typename Node::execution_space>;
1245 RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> domMap = A.getDomainMap();
1246 RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> rowMap = A.getRowMap();
1247 RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> colMap = A.getColMap();
1248 TEUCHOS_ASSERT(dirichletRows.extent(0) == rowMap->getLocalNumElements());
1249 TEUCHOS_ASSERT(dirichletCols.extent(0) == colMap->getLocalNumElements());
1250 TEUCHOS_ASSERT(dirichletDomain.extent(0) == domMap->getLocalNumElements());
1251 RCP<Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> myColsToZero = Xpetra::VectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(colMap, /*zeroOut=*/true);
1252 // Find all local column indices that are in Dirichlet rows, record in myColsToZero as 1.0
1253 auto myColsToZeroView = myColsToZero->getLocalViewDevice(Tpetra::Access::ReadWrite);
1254 auto localMatrix = A.getLocalMatrixDevice();
1255 Kokkos::parallel_for(
1256 "MueLu:Maxwell1::DetectDirichletCols", range_type(0, rowMap->getLocalNumElements()),
1257 KOKKOS_LAMBDA(const LocalOrdinal row) {
1258 if (dirichletRows(row)) {
1259 auto rowView = localMatrix.row(row);
1260 auto length = rowView.length;
1261
1262 for (decltype(length) colID = 0; colID < length; colID++)
1263 myColsToZeroView(rowView.colidx(colID), 0) = impl_ATS::one();
1264 }
1265 });
1266
1267 RCP<Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> globalColsToZero;
1268 RCP<const Xpetra::Import<LocalOrdinal, GlobalOrdinal, Node>> importer = A.getCrsGraph()->getImporter();
1269 if (!importer.is_null()) {
1270 globalColsToZero = Xpetra::VectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(domMap, /*zeroOut=*/true);
1271 // export to domain map
1272 globalColsToZero->doExport(*myColsToZero, *importer, Xpetra::ADD);
1273 // import to column map
1274 myColsToZero->doImport(*globalColsToZero, *importer, Xpetra::INSERT);
1275 } else
1276 globalColsToZero = myColsToZero;
1277 UtilitiesBase<Scalar, LocalOrdinal, GlobalOrdinal, Node>::FindNonZeros(globalColsToZero->getLocalViewDevice(Tpetra::Access::ReadOnly), dirichletDomain);
1278 UtilitiesBase<Scalar, LocalOrdinal, GlobalOrdinal, Node>::FindNonZeros(myColsToZero->getLocalViewDevice(Tpetra::Access::ReadOnly), dirichletCols);
1279}
1280
1281template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1283 ApplyRowSumCriterion(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A, const typename Teuchos::ScalarTraits<Scalar>::magnitudeType rowSumTol, Teuchos::ArrayRCP<bool>& dirichletRows) {
1284 typedef Teuchos::ScalarTraits<Scalar> STS;
1285 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType MT;
1286 typedef Teuchos::ScalarTraits<MT> MTS;
1287 RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> rowmap = A.getRowMap();
1288 for (LocalOrdinal row = 0; row < Teuchos::as<LocalOrdinal>(rowmap->getLocalNumElements()); ++row) {
1289 size_t nnz = A.getNumEntriesInLocalRow(row);
1290 ArrayView<const LocalOrdinal> indices;
1291 ArrayView<const Scalar> vals;
1292 A.getLocalRowView(row, indices, vals);
1293
1294 Scalar rowsum = STS::zero();
1295 Scalar diagval = STS::zero();
1296
1297 for (LocalOrdinal colID = 0; colID < Teuchos::as<LocalOrdinal>(nnz); colID++) {
1298 LocalOrdinal col = indices[colID];
1299 if (row == col)
1300 diagval = vals[colID];
1301 rowsum += vals[colID];
1302 }
1303 // printf("A(%d,:) row_sum(point) = %6.4e\n",row,rowsum);
1304 if (rowSumTol < MTS::one() && STS::magnitude(rowsum) > STS::magnitude(diagval) * rowSumTol) {
1305 // printf("Row %d triggers rowsum\n",(int)row);
1306 dirichletRows[row] = true;
1307 }
1308 }
1309}
1310
1311template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1312void UtilitiesBase<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
1313 ApplyRowSumCriterion(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A, const Xpetra::Vector<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>& BlockNumber, const typename Teuchos::ScalarTraits<Scalar>::magnitudeType rowSumTol, Teuchos::ArrayRCP<bool>& dirichletRows) {
1314 typedef Teuchos::ScalarTraits<Scalar> STS;
1315 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType MT;
1316 typedef Teuchos::ScalarTraits<MT> MTS;
1317 RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> rowmap = A.getRowMap();
1318
1319 TEUCHOS_TEST_FOR_EXCEPTION(!A.getColMap()->isSameAs(*BlockNumber.getMap()), std::runtime_error, "ApplyRowSumCriterion: BlockNumber must match's A's column map.");
1320
1321 Teuchos::ArrayRCP<const LocalOrdinal> block_id = BlockNumber.getData(0);
1322 for (LocalOrdinal row = 0; row < Teuchos::as<LocalOrdinal>(rowmap->getLocalNumElements()); ++row) {
1323 size_t nnz = A.getNumEntriesInLocalRow(row);
1324 ArrayView<const LocalOrdinal> indices;
1325 ArrayView<const Scalar> vals;
1326 A.getLocalRowView(row, indices, vals);
1327
1328 Scalar rowsum = STS::zero();
1329 Scalar diagval = STS::zero();
1330 for (LocalOrdinal colID = 0; colID < Teuchos::as<LocalOrdinal>(nnz); colID++) {
1331 LocalOrdinal col = indices[colID];
1332 if (row == col)
1333 diagval = vals[colID];
1334 if (block_id[row] == block_id[col])
1335 rowsum += vals[colID];
1336 }
1337
1338 // printf("A(%d,:) row_sum(block) = %6.4e\n",row,rowsum);
1339 if (rowSumTol < MTS::one() && STS::magnitude(rowsum) > STS::magnitude(diagval) * rowSumTol) {
1340 // printf("Row %d triggers rowsum\n",(int)row);
1341 dirichletRows[row] = true;
1342 }
1343 }
1344}
1345
1346// Applies rowsum criterion
1347template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class memory_space>
1348void ApplyRowSumCriterion(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
1349 const typename Teuchos::ScalarTraits<Scalar>::magnitudeType rowSumTol,
1350 Kokkos::View<bool*, memory_space>& dirichletRows) {
1351 typedef Teuchos::ScalarTraits<Scalar> STS;
1352 RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> rowmap = A.getRowMap();
1353
1354 auto dirichletRowsHost = Kokkos::create_mirror_view(dirichletRows);
1355 Kokkos::deep_copy(dirichletRowsHost, dirichletRows);
1356
1357 for (LocalOrdinal row = 0; row < Teuchos::as<LocalOrdinal>(rowmap->getLocalNumElements()); ++row) {
1358 size_t nnz = A.getNumEntriesInLocalRow(row);
1359 ArrayView<const LocalOrdinal> indices;
1360 ArrayView<const Scalar> vals;
1361 A.getLocalRowView(row, indices, vals);
1362
1363 Scalar rowsum = STS::zero();
1364 Scalar diagval = STS::zero();
1365 for (LocalOrdinal colID = 0; colID < Teuchos::as<LocalOrdinal>(nnz); colID++) {
1366 LocalOrdinal col = indices[colID];
1367 if (row == col)
1368 diagval = vals[colID];
1369 rowsum += vals[colID];
1370 }
1371 if (STS::real(rowsum) > STS::magnitude(diagval) * rowSumTol)
1372 dirichletRowsHost(row) = true;
1373 }
1374
1375 Kokkos::deep_copy(dirichletRows, dirichletRowsHost);
1376}
1377
1378template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1379void UtilitiesBase<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
1380 ApplyRowSumCriterion(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
1381 const typename Teuchos::ScalarTraits<Scalar>::magnitudeType rowSumTol,
1382 Kokkos::View<bool*, typename Node::device_type::memory_space>& dirichletRows) {
1383 MueLu::ApplyRowSumCriterion<Scalar, LocalOrdinal, GlobalOrdinal, Node, typename Node::device_type::memory_space>(A, rowSumTol, dirichletRows);
1384}
1385
1386template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1387void UtilitiesBase<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
1388 ApplyRowSumCriterionHost(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
1389 const typename Teuchos::ScalarTraits<Scalar>::magnitudeType rowSumTol,
1390 Kokkos::View<bool*, Kokkos::HostSpace>& dirichletRows) {
1391 MueLu::ApplyRowSumCriterion<Scalar, LocalOrdinal, GlobalOrdinal, Node, Kokkos::HostSpace>(A, rowSumTol, dirichletRows);
1392}
1393
1394// Applies rowsum criterion
1395template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class memory_space>
1396void ApplyRowSumCriterion(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
1397 const Xpetra::Vector<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>& BlockNumber,
1398 const typename Teuchos::ScalarTraits<Scalar>::magnitudeType rowSumTol,
1399 Kokkos::View<bool*, memory_space>& dirichletRows) {
1400 typedef Teuchos::ScalarTraits<Scalar> STS;
1401 RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> rowmap = A.getRowMap();
1402
1403 TEUCHOS_TEST_FOR_EXCEPTION(!A.getColMap()->isSameAs(*BlockNumber.getMap()), std::runtime_error, "ApplyRowSumCriterion: BlockNumber must match's A's column map.");
1404
1405 auto dirichletRowsHost = Kokkos::create_mirror_view(dirichletRows);
1406 Kokkos::deep_copy(dirichletRowsHost, dirichletRows);
1407
1408 Teuchos::ArrayRCP<const LocalOrdinal> block_id = BlockNumber.getData(0);
1409 for (LocalOrdinal row = 0; row < Teuchos::as<LocalOrdinal>(rowmap->getLocalNumElements()); ++row) {
1410 size_t nnz = A.getNumEntriesInLocalRow(row);
1411 ArrayView<const LocalOrdinal> indices;
1412 ArrayView<const Scalar> vals;
1413 A.getLocalRowView(row, indices, vals);
1414
1415 Scalar rowsum = STS::zero();
1416 Scalar diagval = STS::zero();
1417 for (LocalOrdinal colID = 0; colID < Teuchos::as<LocalOrdinal>(nnz); colID++) {
1418 LocalOrdinal col = indices[colID];
1419 if (row == col)
1420 diagval = vals[colID];
1421 if (block_id[row] == block_id[col])
1422 rowsum += vals[colID];
1423 }
1424 if (STS::real(rowsum) > STS::magnitude(diagval) * rowSumTol)
1425 dirichletRowsHost(row) = true;
1426 }
1427
1428 Kokkos::deep_copy(dirichletRows, dirichletRowsHost);
1429}
1430
1431template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1432void UtilitiesBase<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
1433 ApplyRowSumCriterion(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
1434 const Xpetra::Vector<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>& BlockNumber,
1435 const typename Teuchos::ScalarTraits<Scalar>::magnitudeType rowSumTol,
1436 Kokkos::View<bool*, typename Node::device_type::memory_space>& dirichletRows) {
1437 MueLu::ApplyRowSumCriterion<Scalar, LocalOrdinal, GlobalOrdinal, Node, typename Node::device_type::memory_space>(A, BlockNumber, rowSumTol, dirichletRows);
1438}
1439
1440template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1441void UtilitiesBase<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
1442 ApplyRowSumCriterionHost(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
1443 const Xpetra::Vector<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>& BlockNumber,
1444 const typename Teuchos::ScalarTraits<Scalar>::magnitudeType rowSumTol,
1445 Kokkos::View<bool*, Kokkos::HostSpace>& dirichletRows) {
1446 MueLu::ApplyRowSumCriterion<Scalar, LocalOrdinal, GlobalOrdinal, Node, Kokkos::HostSpace>(A, BlockNumber, rowSumTol, dirichletRows);
1447}
1448
1449template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1450Teuchos::ArrayRCP<const bool>
1452 DetectDirichletCols(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
1453 const Teuchos::ArrayRCP<const bool>& dirichletRows) {
1454 Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
1455 Scalar one = Teuchos::ScalarTraits<Scalar>::one();
1456 Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> domMap = A.getDomainMap();
1457 Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> colMap = A.getColMap();
1458 Teuchos::RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> myColsToZero = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(colMap, 1);
1459 myColsToZero->putScalar(zero);
1460 // Find all local column indices that are in Dirichlet rows, record in myColsToZero as 1.0
1461 for (size_t i = 0; i < (size_t)dirichletRows.size(); i++) {
1462 if (dirichletRows[i]) {
1463 Teuchos::ArrayView<const LocalOrdinal> indices;
1464 Teuchos::ArrayView<const Scalar> values;
1465 A.getLocalRowView(i, indices, values);
1466 for (size_t j = 0; j < static_cast<size_t>(indices.size()); j++)
1467 myColsToZero->replaceLocalValue(indices[j], 0, one);
1468 }
1469 }
1470
1471 Teuchos::RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> globalColsToZero = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(domMap, 1);
1472 globalColsToZero->putScalar(zero);
1473 Teuchos::RCP<Xpetra::Export<LocalOrdinal, GlobalOrdinal, Node>> exporter = Xpetra::ExportFactory<LocalOrdinal, GlobalOrdinal, Node>::Build(colMap, domMap);
1474 // export to domain map
1475 globalColsToZero->doExport(*myColsToZero, *exporter, Xpetra::ADD);
1476 // import to column map
1477 myColsToZero->doImport(*globalColsToZero, *exporter, Xpetra::INSERT);
1478 Teuchos::ArrayRCP<const Scalar> myCols = myColsToZero->getData(0);
1479 Teuchos::ArrayRCP<bool> dirichletCols(colMap->getLocalNumElements(), true);
1480 Magnitude eps = Teuchos::ScalarTraits<Magnitude>::eps();
1481 for (size_t i = 0; i < colMap->getLocalNumElements(); i++) {
1482 dirichletCols[i] = Teuchos::ScalarTraits<Scalar>::magnitude(myCols[i]) > 2.0 * eps;
1483 }
1484 return dirichletCols;
1485}
1486
1487template <class SC, class LO, class GO, class NO>
1488Kokkos::View<bool*, typename NO::device_type>
1490 DetectDirichletCols(const Xpetra::Matrix<SC, LO, GO, NO>& A,
1491 const Kokkos::View<const bool*, typename NO::device_type>& dirichletRows) {
1492 using ATS = KokkosKernels::ArithTraits<SC>;
1493 using impl_ATS = KokkosKernels::ArithTraits<typename ATS::val_type>;
1494 using range_type = Kokkos::RangePolicy<LO, typename NO::execution_space>;
1495
1496 SC zero = ATS::zero();
1497
1498 auto localMatrix = A.getLocalMatrixDevice();
1499 LO numRows = A.getLocalNumRows();
1500
1501 Teuchos::RCP<const Xpetra::Map<LO, GO, NO>> domMap = A.getDomainMap();
1502 Teuchos::RCP<const Xpetra::Map<LO, GO, NO>> colMap = A.getColMap();
1503 Teuchos::RCP<Xpetra::MultiVector<SC, LO, GO, NO>> myColsToZero = Xpetra::MultiVectorFactory<SC, LO, GO, NO>::Build(colMap, 1);
1504 myColsToZero->putScalar(zero);
1505 auto myColsToZeroView = myColsToZero->getLocalViewDevice(Tpetra::Access::ReadWrite);
1506 // Find all local column indices that are in Dirichlet rows, record in myColsToZero as 1.0
1507 Kokkos::parallel_for(
1508 "MueLu:Utils::DetectDirichletCols1", range_type(0, numRows),
1509 KOKKOS_LAMBDA(const LO row) {
1510 if (dirichletRows(row)) {
1511 auto rowView = localMatrix.row(row);
1512 auto length = rowView.length;
1513
1514 for (decltype(length) colID = 0; colID < length; colID++)
1515 myColsToZeroView(rowView.colidx(colID), 0) = impl_ATS::one();
1516 }
1517 });
1518
1519 Teuchos::RCP<Xpetra::MultiVector<SC, LO, GO, NO>> globalColsToZero = Xpetra::MultiVectorFactory<SC, LO, GO, NO>::Build(domMap, 1);
1520 globalColsToZero->putScalar(zero);
1521 Teuchos::RCP<Xpetra::Export<LO, GO, NO>> exporter = Xpetra::ExportFactory<LO, GO, NO>::Build(colMap, domMap);
1522 // export to domain map
1523 globalColsToZero->doExport(*myColsToZero, *exporter, Xpetra::ADD);
1524 // import to column map
1525 myColsToZero->doImport(*globalColsToZero, *exporter, Xpetra::INSERT);
1526
1527 auto myCols = myColsToZero->getLocalViewDevice(Tpetra::Access::ReadOnly);
1528 size_t numColEntries = colMap->getLocalNumElements();
1529 Kokkos::View<bool*, typename NO::device_type> dirichletCols(Kokkos::ViewAllocateWithoutInitializing("dirichletCols"), numColEntries);
1530 const typename ATS::magnitudeType eps = 2.0 * ATS::eps();
1531
1532 Kokkos::parallel_for(
1533 "MueLu:Utils::DetectDirichletCols2", range_type(0, numColEntries),
1534 KOKKOS_LAMBDA(const size_t i) {
1535 dirichletCols(i) = impl_ATS::magnitude(myCols(i, 0)) > eps;
1536 });
1537 return dirichletCols;
1538}
1539
1540template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1541Scalar
1543 Frobenius(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A, const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& B) {
1544 // We check only row maps. Column may be different. One would hope that they are the same, as we typically
1545 // calculate frobenius norm of the specified sparsity pattern with an updated matrix from the previous step,
1546 // but matrix addition, even when one is submatrix of the other, changes column map (though change may be as
1547 // simple as couple of elements swapped)
1548 TEUCHOS_TEST_FOR_EXCEPTION(!A.getRowMap()->isSameAs(*B.getRowMap()), Exceptions::Incompatible, "MueLu::CGSolver::Frobenius: row maps are incompatible");
1549 TEUCHOS_TEST_FOR_EXCEPTION(!A.isFillComplete() || !B.isFillComplete(), Exceptions::RuntimeError, "Matrices must be fill completed");
1550
1551 const Map& AColMap = *A.getColMap();
1552 const Map& BColMap = *B.getColMap();
1553
1554 Teuchos::ArrayView<const LocalOrdinal> indA, indB;
1555 Teuchos::ArrayView<const Scalar> valA, valB;
1556 size_t nnzA = 0, nnzB = 0;
1557
1558 // We use a simple algorithm
1559 // for each row we fill valBAll array with the values in the corresponding row of B
1560 // as such, it serves as both sorted array and as storage, so we don't need to do a
1561 // tricky problem: "find a value in the row of B corresponding to the specific GID"
1562 // Once we do that, we translate LID of entries of row of A to LID of B, and multiply
1563 // corresponding entries.
1564 // The algorithm should be reasonably cheap, as it does not sort anything, provided
1565 // that getLocalElement and getGlobalElement functions are reasonably effective. It
1566 // *is* possible that the costs are hidden in those functions, but if maps are close
1567 // to linear maps, we should be fine
1568 Teuchos::Array<Scalar> valBAll(BColMap.getLocalNumElements());
1569
1570 LocalOrdinal invalid = Teuchos::OrdinalTraits<LocalOrdinal>::invalid();
1571 Scalar zero = Teuchos::ScalarTraits<Scalar>::zero(), f = zero, gf;
1572 size_t numRows = A.getLocalNumRows();
1573 for (size_t i = 0; i < numRows; i++) {
1574 A.getLocalRowView(i, indA, valA);
1575 B.getLocalRowView(i, indB, valB);
1576 nnzA = indA.size();
1577 nnzB = indB.size();
1578
1579 // Set up array values
1580 for (size_t j = 0; j < nnzB; j++)
1581 valBAll[indB[j]] = valB[j];
1582
1583 for (size_t j = 0; j < nnzA; j++) {
1584 // The cost of the whole Frobenius dot product function depends on the
1585 // cost of the getLocalElement and getGlobalElement functions here.
1586 LocalOrdinal ind = BColMap.getLocalElement(AColMap.getGlobalElement(indA[j]));
1587 if (ind != invalid)
1588 f += valBAll[ind] * valA[j];
1589 }
1590
1591 // Clean up array values
1592 for (size_t j = 0; j < nnzB; j++)
1593 valBAll[indB[j]] = zero;
1594 }
1595
1596 MueLu_sumAll(AColMap.getComm(), f, gf);
1597
1598 return gf;
1599}
1600
1601template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1604 // Distribute the seeds evenly in [1,maxint-1]. This guarantees nothing
1605 // about where in random number stream we are, but avoids overflow situations
1606 // in parallel when multiplying by a PID. It would be better to use
1607 // a good parallel random number generator.
1608 double one = 1.0;
1609 int maxint = INT_MAX; //= 2^31-1 = 2147483647 for 32-bit integers
1610 int mySeed = Teuchos::as<int>((maxint - 1) * (one - (comm.getRank() + 1) / (comm.getSize() + one)));
1611 if (mySeed < 1 || mySeed == maxint) {
1612 std::ostringstream errStr;
1613 errStr << "Error detected with random seed = " << mySeed << ". It should be in the interval [1,2^31-2].";
1614 throw Exceptions::RuntimeError(errStr.str());
1615 }
1616 std::srand(mySeed);
1617 // For Tpetra, we could use Kokkos' random number generator here.
1618 Teuchos::ScalarTraits<Scalar>::seedrandom(mySeed);
1619 // Epetra
1620 // MultiVector::Random() -> Epetra_Util::RandomDouble() -> Epetra_Utils::RandomInt()
1621 // Its own random number generator, based on Seed_. Seed_ is initialized in Epetra_Util constructor with std::rand()
1622 // So our setting std::srand() affects that too
1623}
1624
1625template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1627 FindDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& A,
1628 std::vector<LocalOrdinal>& dirichletRows, bool count_twos_as_dirichlet) {
1629 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType MT;
1630 dirichletRows.resize(0);
1631 for (size_t i = 0; i < A->getLocalNumRows(); i++) {
1632 Teuchos::ArrayView<const LocalOrdinal> indices;
1633 Teuchos::ArrayView<const Scalar> values;
1634 A->getLocalRowView(i, indices, values);
1635 int nnz = 0;
1636 for (size_t j = 0; j < (size_t)indices.size(); j++) {
1637 if (Teuchos::ScalarTraits<Scalar>::magnitude(values[j]) > Teuchos::ScalarTraits<MT>::eps()) {
1638 nnz++;
1639 }
1640 }
1641 if (nnz == 1 || (count_twos_as_dirichlet && nnz == 2)) {
1642 dirichletRows.push_back(i);
1643 }
1644 }
1645}
1646
1647template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1649 ApplyOAZToMatrixRows(Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& A,
1650 const std::vector<LocalOrdinal>& dirichletRows) {
1651 RCP<const Map> Rmap = A->getRowMap();
1652 RCP<const Map> Cmap = A->getColMap();
1653 Scalar one = Teuchos::ScalarTraits<Scalar>::one();
1654 Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
1655
1656 for (size_t i = 0; i < dirichletRows.size(); i++) {
1657 GlobalOrdinal row_gid = Rmap->getGlobalElement(dirichletRows[i]);
1658
1659 Teuchos::ArrayView<const LocalOrdinal> indices;
1660 Teuchos::ArrayView<const Scalar> values;
1661 A->getLocalRowView(dirichletRows[i], indices, values);
1662 // NOTE: This won't work with fancy node types.
1663 Scalar* valuesNC = const_cast<Scalar*>(values.getRawPtr());
1664 for (size_t j = 0; j < (size_t)indices.size(); j++) {
1665 if (Cmap->getGlobalElement(indices[j]) == row_gid)
1666 valuesNC[j] = one;
1667 else
1668 valuesNC[j] = zero;
1669 }
1670 }
1671}
1672
1673template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1675 ApplyOAZToMatrixRows(Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& A,
1676 const Teuchos::ArrayRCP<const bool>& dirichletRows) {
1677 TEUCHOS_ASSERT(A->isFillComplete());
1678 RCP<const Map> domMap = A->getDomainMap();
1679 RCP<const Map> ranMap = A->getRangeMap();
1680 RCP<const Map> Rmap = A->getRowMap();
1681 RCP<const Map> Cmap = A->getColMap();
1682 TEUCHOS_ASSERT(static_cast<size_t>(dirichletRows.size()) == Rmap->getLocalNumElements());
1683 const Scalar one = Teuchos::ScalarTraits<Scalar>::one();
1684 const Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
1685 A->resumeFill();
1686 for (size_t i = 0; i < (size_t)dirichletRows.size(); i++) {
1687 if (dirichletRows[i]) {
1688 GlobalOrdinal row_gid = Rmap->getGlobalElement(i);
1689
1690 Teuchos::ArrayView<const LocalOrdinal> indices;
1691 Teuchos::ArrayView<const Scalar> values;
1692 A->getLocalRowView(i, indices, values);
1693
1694 Teuchos::ArrayRCP<Scalar> valuesNC(values.size());
1695 for (size_t j = 0; j < (size_t)indices.size(); j++) {
1696 if (Cmap->getGlobalElement(indices[j]) == row_gid)
1697 valuesNC[j] = one;
1698 else
1699 valuesNC[j] = zero;
1700 }
1701 A->replaceLocalValues(i, indices, valuesNC());
1702 }
1703 }
1704 A->fillComplete(domMap, ranMap);
1705}
1706
1707template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1709 ApplyOAZToMatrixRows(Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& A,
1710 const Kokkos::View<const bool*, typename Node::device_type>& dirichletRows) {
1711 TEUCHOS_ASSERT(A->isFillComplete());
1712 using ATS = KokkosKernels::ArithTraits<Scalar>;
1713 using impl_ATS = KokkosKernels::ArithTraits<typename ATS::val_type>;
1714 using range_type = Kokkos::RangePolicy<LocalOrdinal, typename Node::execution_space>;
1715
1716 RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> domMap = A->getDomainMap();
1717 RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> ranMap = A->getRangeMap();
1718 RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> Rmap = A->getRowMap();
1719 RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> Cmap = A->getColMap();
1720
1721 TEUCHOS_ASSERT(static_cast<size_t>(dirichletRows.size()) == Rmap->getLocalNumElements());
1722
1723 auto localMatrix = A->getLocalMatrixDevice();
1724 auto localRmap = Rmap->getLocalMap();
1725 auto localCmap = Cmap->getLocalMap();
1726
1727 Kokkos::parallel_for(
1728 "MueLu::Utils::ApplyOAZ", range_type(0, dirichletRows.extent(0)),
1729 KOKKOS_LAMBDA(const LocalOrdinal row) {
1730 if (dirichletRows(row)) {
1731 auto rowView = localMatrix.row(row);
1732 auto length = rowView.length;
1733 auto row_gid = localRmap.getGlobalElement(row);
1734 auto row_lid = localCmap.getLocalElement(row_gid);
1735
1736 for (decltype(length) colID = 0; colID < length; colID++)
1737 if (rowView.colidx(colID) == row_lid)
1738 rowView.value(colID) = impl_ATS::one();
1739 else
1740 rowView.value(colID) = impl_ATS::zero();
1741 }
1742 });
1743}
1744
1745template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1747 ZeroDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& A,
1748 const std::vector<LocalOrdinal>& dirichletRows,
1749 Scalar replaceWith) {
1750 for (size_t i = 0; i < dirichletRows.size(); i++) {
1751 Teuchos::ArrayView<const LocalOrdinal> indices;
1752 Teuchos::ArrayView<const Scalar> values;
1753 A->getLocalRowView(dirichletRows[i], indices, values);
1754 // NOTE: This won't work with fancy node types.
1755 Scalar* valuesNC = const_cast<Scalar*>(values.getRawPtr());
1756 for (size_t j = 0; j < (size_t)indices.size(); j++)
1757 valuesNC[j] = replaceWith;
1758 }
1759}
1760
1761template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1763 ZeroDirichletRows(Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& A,
1764 const Teuchos::ArrayRCP<const bool>& dirichletRows,
1765 Scalar replaceWith) {
1766 TEUCHOS_ASSERT(static_cast<size_t>(dirichletRows.size()) == A->getRowMap()->getLocalNumElements());
1767 for (size_t i = 0; i < (size_t)dirichletRows.size(); i++) {
1768 if (dirichletRows[i]) {
1769 Teuchos::ArrayView<const LocalOrdinal> indices;
1770 Teuchos::ArrayView<const Scalar> values;
1771 A->getLocalRowView(i, indices, values);
1772 // NOTE: This won't work with fancy node types.
1773 Scalar* valuesNC = const_cast<Scalar*>(values.getRawPtr());
1774 for (size_t j = 0; j < (size_t)indices.size(); j++)
1775 valuesNC[j] = replaceWith;
1776 }
1777 }
1778}
1779
1780template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1782 ZeroDirichletRows(Teuchos::RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& X,
1783 const Teuchos::ArrayRCP<const bool>& dirichletRows,
1784 Scalar replaceWith) {
1785 TEUCHOS_ASSERT(static_cast<size_t>(dirichletRows.size()) == X->getMap()->getLocalNumElements());
1786 for (size_t i = 0; i < (size_t)dirichletRows.size(); i++) {
1787 if (dirichletRows[i]) {
1788 for (size_t j = 0; j < X->getNumVectors(); j++)
1789 X->replaceLocalValue(i, j, replaceWith);
1790 }
1791 }
1792}
1793
1794template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1796 ZeroDirichletRows(RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& A,
1797 const Kokkos::View<const bool*, typename Node::device_type>& dirichletRows,
1798 Scalar replaceWith) {
1799 using ATS = KokkosKernels::ArithTraits<Scalar>;
1800 using range_type = Kokkos::RangePolicy<LocalOrdinal, typename Node::execution_space>;
1801
1802 typename ATS::val_type impl_replaceWith = replaceWith;
1803
1804 auto localMatrix = A->getLocalMatrixDevice();
1805 LocalOrdinal numRows = A->getLocalNumRows();
1806
1807 Kokkos::parallel_for(
1808 "MueLu:Utils::ZeroDirichletRows", range_type(0, numRows),
1809 KOKKOS_LAMBDA(const LocalOrdinal row) {
1810 if (dirichletRows(row)) {
1811 auto rowView = localMatrix.row(row);
1812 auto length = rowView.length;
1813 for (decltype(length) colID = 0; colID < length; colID++)
1814 rowView.value(colID) = impl_replaceWith;
1815 }
1816 });
1817}
1818
1819template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1820void UtilitiesBase<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
1821 ZeroDirichletRows(RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& X,
1822 const Kokkos::View<const bool*, typename Node::device_type>& dirichletRows,
1823 Scalar replaceWith) {
1824 using ATS = KokkosKernels::ArithTraits<Scalar>;
1825 using range_type = Kokkos::RangePolicy<LocalOrdinal, typename Node::execution_space>;
1826
1827 typename ATS::val_type impl_replaceWith = replaceWith;
1828
1829 auto myCols = X->getLocalViewDevice(Tpetra::Access::ReadWrite);
1830 size_t numVecs = X->getNumVectors();
1831 Kokkos::parallel_for(
1832 "MueLu:Utils::ZeroDirichletRows_MV", range_type(0, dirichletRows.size()),
1833 KOKKOS_LAMBDA(const size_t i) {
1834 if (dirichletRows(i)) {
1835 for (size_t j = 0; j < numVecs; j++)
1836 myCols(i, j) = impl_replaceWith;
1837 }
1838 });
1839}
1840
1841template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1843 ZeroDirichletCols(Teuchos::RCP<Matrix>& A,
1844 const Teuchos::ArrayRCP<const bool>& dirichletCols,
1845 Scalar replaceWith) {
1846 TEUCHOS_ASSERT(static_cast<size_t>(dirichletCols.size()) == A->getColMap()->getLocalNumElements());
1847 for (size_t i = 0; i < A->getLocalNumRows(); i++) {
1848 Teuchos::ArrayView<const LocalOrdinal> indices;
1849 Teuchos::ArrayView<const Scalar> values;
1850 A->getLocalRowView(i, indices, values);
1851 // NOTE: This won't work with fancy node types.
1852 Scalar* valuesNC = const_cast<Scalar*>(values.getRawPtr());
1853 for (size_t j = 0; j < static_cast<size_t>(indices.size()); j++)
1854 if (dirichletCols[indices[j]])
1855 valuesNC[j] = replaceWith;
1856 }
1857}
1858
1859template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1861 ZeroDirichletCols(RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& A,
1862 const Kokkos::View<const bool*, typename Node::device_type>& dirichletCols,
1863 Scalar replaceWith) {
1864 using ATS = KokkosKernels::ArithTraits<Scalar>;
1865 using range_type = Kokkos::RangePolicy<LocalOrdinal, typename Node::execution_space>;
1866
1867 typename ATS::val_type impl_replaceWith = replaceWith;
1868
1869 auto localMatrix = A->getLocalMatrixDevice();
1870 LocalOrdinal numRows = A->getLocalNumRows();
1871
1872 Kokkos::parallel_for(
1873 "MueLu:Utils::ZeroDirichletCols", range_type(0, numRows),
1874 KOKKOS_LAMBDA(const LocalOrdinal row) {
1875 auto rowView = localMatrix.row(row);
1876 auto length = rowView.length;
1877 for (decltype(length) colID = 0; colID < length; colID++)
1878 if (dirichletCols(rowView.colidx(colID))) {
1879 rowView.value(colID) = impl_replaceWith;
1880 }
1881 });
1882}
1883
1884template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1886 FindDirichletRowsAndPropagateToCols(Teuchos::RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& A,
1887 Teuchos::RCP<Xpetra::Vector<int, LocalOrdinal, GlobalOrdinal, Node>>& isDirichletRow,
1888 Teuchos::RCP<Xpetra::Vector<int, LocalOrdinal, GlobalOrdinal, Node>>& isDirichletCol) {
1889 // Make sure A's RowMap == DomainMap
1890 if (!A->getRowMap()->isSameAs(*A->getDomainMap())) {
1891 throw std::runtime_error("UtilitiesBase::FindDirichletRowsAndPropagateToCols row and domain maps must match.");
1892 }
1893 RCP<const Xpetra::Import<LocalOrdinal, GlobalOrdinal, Node>> importer = A->getCrsGraph()->getImporter();
1894 bool has_import = !importer.is_null();
1895
1896 // Find the Dirichlet rows
1897 std::vector<LocalOrdinal> dirichletRows;
1898 FindDirichletRows(A, dirichletRows);
1899
1900#if 0
1901 printf("[%d] DirichletRow Ids = ",A->getRowMap()->getComm()->getRank());
1902 for(size_t i=0; i<(size_t) dirichletRows.size(); i++)
1903 printf("%d ",dirichletRows[i]);
1904 printf("\n");
1905 fflush(stdout);
1906#endif
1907 // Allocate all as non-Dirichlet
1908 isDirichletRow = Xpetra::VectorFactory<int, LocalOrdinal, GlobalOrdinal, Node>::Build(A->getRowMap(), true);
1909 isDirichletCol = Xpetra::VectorFactory<int, LocalOrdinal, GlobalOrdinal, Node>::Build(A->getColMap(), true);
1910
1911 {
1912 Teuchos::ArrayRCP<int> dr_rcp = isDirichletRow->getDataNonConst(0);
1913 Teuchos::ArrayView<int> dr = dr_rcp();
1914 Teuchos::ArrayRCP<int> dc_rcp = isDirichletCol->getDataNonConst(0);
1915 Teuchos::ArrayView<int> dc = dc_rcp();
1916 for (size_t i = 0; i < (size_t)dirichletRows.size(); i++) {
1917 dr[dirichletRows[i]] = 1;
1918 if (!has_import) dc[dirichletRows[i]] = 1;
1919 }
1920 }
1921
1922 if (has_import)
1923 isDirichletCol->doImport(*isDirichletRow, *importer, Xpetra::CombineMode::ADD);
1924}
1925
1926template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1927RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
1929 ReplaceNonZerosWithOnes(const RCP<Matrix>& original) {
1930 using ISC = typename KokkosKernels::ArithTraits<Scalar>::val_type;
1931 using range_type = Kokkos::RangePolicy<LocalOrdinal, typename Node::execution_space>;
1932 using local_matrix_type = typename CrsMatrix::local_matrix_device_type;
1933 using values_type = typename local_matrix_type::values_type;
1934
1935 const ISC ONE = KokkosKernels::ArithTraits<ISC>::one();
1936 const ISC ZERO = KokkosKernels::ArithTraits<ISC>::zero();
1937
1938 // Copy the values array of the old matrix to a new array, replacing all the non-zeros with one
1939 auto localMatrix = original->getLocalMatrixDevice();
1940 TEUCHOS_TEST_FOR_EXCEPTION(!original->hasCrsGraph(), Exceptions::RuntimeError, "ReplaceNonZerosWithOnes: Cannot get CrsGraph");
1941 values_type new_values("values", localMatrix.nnz());
1942
1943 Kokkos::parallel_for(
1944 "ReplaceNonZerosWithOnes", range_type(0, localMatrix.nnz()), KOKKOS_LAMBDA(const size_t i) {
1945 if (localMatrix.values(i) != ZERO)
1946 new_values(i) = ONE;
1947 else
1948 new_values(i) = ZERO;
1949 });
1950
1951 // Build the new matrix
1952 RCP<Matrix> NewMatrix = Xpetra::MatrixFactory<SC, LO, GO, NO>::Build(original->getCrsGraph(), new_values);
1953 TEUCHOS_TEST_FOR_EXCEPTION(NewMatrix.is_null(), Exceptions::RuntimeError, "ReplaceNonZerosWithOnes: MatrixFactory::Build() did not return matrix");
1954 NewMatrix->fillComplete(original->getDomainMap(), original->getRangeMap());
1955 return NewMatrix;
1956}
1957
1958template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1959RCP<const Xpetra::BlockedMap<LocalOrdinal, GlobalOrdinal, Node>>
1961 GeneratedBlockedTargetMap(const Xpetra::BlockedMap<LocalOrdinal, GlobalOrdinal, Node>& sourceBlockedMap,
1962 const Xpetra::Import<LocalOrdinal, GlobalOrdinal, Node>& Importer) {
1963 typedef Xpetra::Vector<int, LocalOrdinal, GlobalOrdinal, Node> IntVector;
1964 Xpetra::UnderlyingLib lib = sourceBlockedMap.lib();
1965
1966 // De-stride the map if we have to (might regret this later)
1967 RCP<const Map> fullMap = sourceBlockedMap.getMap();
1968 RCP<const Map> stridedMap = Teuchos::rcp_dynamic_cast<const Xpetra::StridedMap<LocalOrdinal, GlobalOrdinal, Node>>(fullMap);
1969 if (!stridedMap.is_null()) fullMap = stridedMap->getMap();
1970
1971 // Initial sanity checking for map compatibil
1972 const size_t numSubMaps = sourceBlockedMap.getNumMaps();
1973 if (!Importer.getSourceMap()->isCompatible(*fullMap))
1974 throw std::runtime_error("GenerateBlockedTargetMap(): Map compatibility error");
1975
1976 // Build an indicator vector
1977 RCP<IntVector> block_ids = Xpetra::VectorFactory<int, LocalOrdinal, GlobalOrdinal, Node>::Build(fullMap);
1978
1979 for (size_t i = 0; i < numSubMaps; i++) {
1980 RCP<const Map> map = sourceBlockedMap.getMap(i);
1981
1982 for (size_t j = 0; j < map->getLocalNumElements(); j++) {
1983 LocalOrdinal jj = fullMap->getLocalElement(map->getGlobalElement(j));
1984 block_ids->replaceLocalValue(jj, (int)i);
1985 }
1986 }
1987
1988 // Get the block ids for the new map
1989 RCP<const Map> targetMap = Importer.getTargetMap();
1990 RCP<IntVector> new_block_ids = Xpetra::VectorFactory<int, LocalOrdinal, GlobalOrdinal, Node>::Build(targetMap);
1991 new_block_ids->doImport(*block_ids, Importer, Xpetra::CombineMode::ADD);
1992 Teuchos::ArrayRCP<const int> dataRCP = new_block_ids->getData(0);
1993 Teuchos::ArrayView<const int> data = dataRCP();
1994
1995 // Get the GIDs for each subblock
1996 Teuchos::Array<Teuchos::Array<GlobalOrdinal>> elementsInSubMap(numSubMaps);
1997 for (size_t i = 0; i < targetMap->getLocalNumElements(); i++) {
1998 elementsInSubMap[data[i]].push_back(targetMap->getGlobalElement(i));
1999 }
2000
2001 // Generate the new submaps
2002 std::vector<RCP<const Map>> subMaps(numSubMaps);
2003 for (size_t i = 0; i < numSubMaps; i++) {
2004 subMaps[i] = Xpetra::MapFactory<LocalOrdinal, GlobalOrdinal, Node>::Build(lib, Teuchos::OrdinalTraits<GlobalOrdinal>::invalid(), elementsInSubMap[i](), targetMap->getIndexBase(), targetMap->getComm());
2005 }
2006
2007 // Build the BlockedMap
2008 return rcp(new BlockedMap(targetMap, subMaps));
2009}
2010
2011template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
2013 MapsAreNested(const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>& rowMap, const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>& colMap) {
2014 if ((rowMap.lib() == Xpetra::UseTpetra) && (colMap.lib() == Xpetra::UseTpetra)) {
2015 auto tpRowMap = toTpetra(rowMap);
2016 auto tpColMap = toTpetra(colMap);
2017 return tpColMap.isLocallyFitted(tpRowMap);
2018 }
2019
2020 ArrayView<const GlobalOrdinal> rowElements = rowMap.getLocalElementList();
2021 ArrayView<const GlobalOrdinal> colElements = colMap.getLocalElementList();
2022
2023 const size_t numElements = rowElements.size();
2024
2025 if (size_t(colElements.size()) < numElements)
2026 return false;
2027
2028 bool goodMap = true;
2029 for (size_t i = 0; i < numElements; i++)
2030 if (rowElements[i] != colElements[i]) {
2031 goodMap = false;
2032 break;
2033 }
2034
2035 return goodMap;
2036}
2037
2038template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
2039Teuchos::RCP<Xpetra::Vector<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>>
2041 ReverseCuthillMcKee(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op) {
2042 using local_matrix_type = typename Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_device_type;
2043 using local_graph_type = typename local_matrix_type::staticcrsgraph_type;
2044 using lno_nnz_view_t = typename local_graph_type::entries_type::non_const_type;
2045 using device = typename local_graph_type::device_type;
2046 using execution_space = typename local_matrix_type::execution_space;
2047 using ordinal_type = typename local_matrix_type::ordinal_type;
2048
2049 local_graph_type localGraph = Op.getLocalMatrixDevice().graph;
2050
2051 lno_nnz_view_t rcmOrder = KokkosGraph::Experimental::graph_rcm<device, typename local_graph_type::row_map_type, typename local_graph_type::entries_type, lno_nnz_view_t>(localGraph.row_map, localGraph.entries);
2052
2053 RCP<Xpetra::Vector<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>> retval =
2054 Xpetra::VectorFactory<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>::Build(Op.getRowMap());
2055
2056 // Copy out and reorder data
2057 auto view1D = Kokkos::subview(retval->getLocalViewDevice(Tpetra::Access::ReadWrite), Kokkos::ALL(), 0);
2058 Kokkos::parallel_for(
2059 "Utilities::ReverseCuthillMcKee",
2060 Kokkos::RangePolicy<ordinal_type, execution_space>(0, localGraph.numRows()),
2061 KOKKOS_LAMBDA(const ordinal_type rowIdx) {
2062 view1D(rcmOrder(rowIdx)) = rowIdx;
2063 });
2064 return retval;
2065}
2066
2067template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
2068Teuchos::RCP<Xpetra::Vector<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>>
2070 CuthillMcKee(const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op) {
2071 using local_matrix_type = typename Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_device_type;
2072 using local_graph_type = typename local_matrix_type::staticcrsgraph_type;
2073 using lno_nnz_view_t = typename local_graph_type::entries_type::non_const_type;
2074 using device = typename local_graph_type::device_type;
2075 using execution_space = typename local_matrix_type::execution_space;
2076 using ordinal_type = typename local_matrix_type::ordinal_type;
2077
2078 local_graph_type localGraph = Op.getLocalMatrixDevice().graph;
2079 LocalOrdinal numRows = localGraph.numRows();
2080
2081 lno_nnz_view_t rcmOrder = KokkosGraph::Experimental::graph_rcm<device, typename local_graph_type::row_map_type, typename local_graph_type::entries_type, lno_nnz_view_t>(localGraph.row_map, localGraph.entries);
2082
2083 RCP<Xpetra::Vector<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>> retval =
2084 Xpetra::VectorFactory<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>::Build(Op.getRowMap());
2085
2086 // Copy out data
2087 auto view1D = Kokkos::subview(retval->getLocalViewDevice(Tpetra::Access::ReadWrite), Kokkos::ALL(), 0);
2088 // Since KokkosKernels produced RCM, also reverse the order of the view to get CM
2089 Kokkos::parallel_for(
2090 "Utilities::ReverseCuthillMcKee",
2091 Kokkos::RangePolicy<ordinal_type, execution_space>(0, numRows),
2092 KOKKOS_LAMBDA(const ordinal_type rowIdx) {
2093 view1D(rcmOrder(numRows - 1 - rowIdx)) = rowIdx;
2094 });
2095 return retval;
2096}
2097
2098} // namespace MueLu
2099
2100#define MUELU_UTILITIESBASE_SHORT
2101#endif // MUELU_UTILITIESBASE_DEF_HPP
2102
2103// LocalWords: LocalOrdinal
#define MueLu_sumAll(rcpComm, in, out)
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Exception throws to report incompatible objects (like maps).
Exception throws to report errors in the internal logical of the program.
Teuchos::ScalarTraits< Scalar >::magnitudeType Magnitude
static Teuchos::RCP< Vector > GetInverse(Teuchos::RCP< const Vector > v, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps() *100, Scalar valReplacement=Teuchos::ScalarTraits< Scalar >::zero())
Return vector containing inverse of input vector.
static RCP< Vector > GetMatrixDiagonalInverse(const Matrix &A, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps() *100, Scalar valReplacement=Teuchos::ScalarTraits< Scalar >::zero(), const bool doLumped=false)
Extract Matrix Diagonal.
static Teuchos::ArrayRCP< Scalar > GetMatrixDiagonal_arcp(const Matrix &A)
Extract Matrix Diagonal.
static RCP< CrsMatrixWrap > GetThresholdedMatrix(const RCP< Matrix > &Ain, const Magnitude threshold, const bool keepDiagonal=true, const GlobalOrdinal expectedNNZperRow=-1)
Threshold a matrix.
static RCP< Matrix > Crs2Op(RCP< CrsMatrix > Op)
static RCP< Vector > GetMatrixDiagonal(const Matrix &A)
Extract Matrix Diagonal.
static RCP< Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > GetThresholdedGraph(const RCP< Matrix > &A, const Magnitude threshold, const GlobalOrdinal expectedNNZperRow=-1)
Threshold a graph.
Namespace for MueLu classes and methods.
Kokkos::View< bool *, memory_space > DetectDirichletRows_kokkos(const Xpetra::Matrix< SC, LO, GO, NO > &A, const typename Teuchos::ScalarTraits< SC >::magnitudeType &tol, const bool count_twos_as_dirichlet)
void ApplyRowSumCriterion(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const typename Teuchos::ScalarTraits< Scalar >::magnitudeType rowSumTol, Kokkos::View< bool *, memory_space > &dirichletRows)
Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > removeSmallEntries(Teuchos::RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const typename Teuchos::ScalarTraits< Scalar >::magnitudeType threshold, const bool keepDiagonal)
KOKKOS_FORCEINLINE_FUNCTION bool isDirichletRow(typename CrsMatrix::ordinal_type rowId, KokkosSparse::SparseRowViewConst< CrsMatrix > &row, const typename KokkosKernels::ArithTraits< typename CrsMatrix::value_type >::magnitudeType &tol, const bool count_twos_as_dirichlet)