MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_AggregateQualityEstimateFactory_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_AGGREGATEQUALITYESTIMATEFACTORY_DEF_HPP
11#define MUELU_AGGREGATEQUALITYESTIMATEFACTORY_DEF_HPP
12#include <iomanip>
14
15#include "MueLu_Level.hpp"
16
17#include <Teuchos_SerialDenseVector.hpp>
18#include <Teuchos_LAPACK.hpp>
19
21#include <Xpetra_IO.hpp>
22#include "MueLu_MasterList.hpp"
23#include "MueLu_Monitor.hpp"
24#include "MueLu_Utilities.hpp"
25
26#include <vector>
27
28namespace MueLu {
29
30template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
32
33template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
35
36template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
38 Input(fineLevel, "A");
39 Input(fineLevel, "Aggregates");
40 Input(fineLevel, "CoarseMap");
41}
42
43template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
45 RCP<ParameterList> validParamList = rcp(new ParameterList());
46
47#define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
48 SET_VALID_ENTRY("aggregate qualities: good aggregate threshold");
49 SET_VALID_ENTRY("aggregate qualities: file output");
50 SET_VALID_ENTRY("aggregate qualities: file base");
51 SET_VALID_ENTRY("aggregate qualities: check symmetry");
52 SET_VALID_ENTRY("aggregate qualities: algorithm");
53 SET_VALID_ENTRY("aggregate qualities: zero threshold");
54 SET_VALID_ENTRY("aggregate qualities: percentiles");
55 SET_VALID_ENTRY("aggregate qualities: mode");
56
57#undef SET_VALID_ENTRY
58
59 validParamList->set<RCP<const FactoryBase>>("A", Teuchos::null, "Generating factory of the matrix A");
60 validParamList->set<RCP<const FactoryBase>>("Aggregates", Teuchos::null, "Generating factory of the aggregates");
61 validParamList->set<RCP<const FactoryBase>>("CoarseMap", Teuchos::null, "Generating factory of the coarse map");
62
63 return validParamList;
64}
65
66template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
68 FactoryMonitor m(*this, "Build", fineLevel);
69
70 RCP<Matrix> A = Get<RCP<Matrix>>(fineLevel, "A");
71 RCP<Aggregates> aggregates = Get<RCP<Aggregates>>(fineLevel, "Aggregates");
72
73 RCP<const Map> map = Get<RCP<const Map>>(fineLevel, "CoarseMap");
74
75 assert(!aggregates->AggregatesCrossProcessors());
76 ParameterList pL = GetParameterList();
77 std::string mode = pL.get<std::string>("aggregate qualities: mode");
78 GetOStream(Statistics1) << "AggregateQuality: mode " << mode << std::endl;
79
80 RCP<Xpetra::MultiVector<magnitudeType, LO, GO, NO>> aggregate_qualities;
81 if (mode == "eigenvalue" || mode == "both") {
82 aggregate_qualities = Xpetra::MultiVectorFactory<magnitudeType, LO, GO, NO>::Build(map, 1);
83 ComputeAggregateQualities(A, aggregates, aggregate_qualities);
84 OutputAggQualities(fineLevel, aggregate_qualities);
85 }
86 if (mode == "size" || mode == "both") {
87 RCP<LocalOrdinalVector> aggregate_sizes = Xpetra::VectorFactory<LO, LO, GO, NO>::Build(map);
88 ComputeAggregateSizes(A, aggregates, aggregate_sizes);
89 Set(fineLevel, "AggregateSizes", aggregate_sizes);
90 OutputAggSizes(fineLevel, aggregate_sizes);
91 }
92 Set(coarseLevel, "AggregateQualities", aggregate_qualities);
93}
94
95template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
96void AggregateQualityEstimateFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::ConvertAggregatesData(RCP<const Aggregates> aggs, ArrayRCP<LO>& aggSortedVertices, ArrayRCP<LO>& aggsToIndices, ArrayRCP<LO>& aggSizes) {
97 // Reorder local aggregate information into a format amenable to computing
98 // per-aggregate quantities. Specifically, we compute a format
99 // similar to compressed sparse row format for sparse matrices in which
100 // we store all the local vertices in a single array in blocks corresponding
101 // to aggregates. (This array is aggSortedVertices.) We then store a second
102 // array (aggsToIndices) whose k-th element stores the index of the first
103 // vertex in aggregate k in the array aggSortedVertices.
104
105 const LO LO_ZERO = Teuchos::OrdinalTraits<LO>::zero();
106 const LO LO_ONE = Teuchos::OrdinalTraits<LO>::one();
107
108 LO numAggs = aggs->GetNumAggregates();
109 aggSizes = aggs->ComputeAggregateSizesArrayRCP();
110
111 aggsToIndices = ArrayRCP<LO>(numAggs + LO_ONE, LO_ZERO);
112
113 for (LO i = 0; i < numAggs; ++i) {
114 aggsToIndices[i + LO_ONE] = aggsToIndices[i] + aggSizes[i];
115 }
116
117 const RCP<LOMultiVector> vertex2AggId = aggs->GetVertex2AggId();
118 const ArrayRCP<const LO> vertex2AggIdData = vertex2AggId->getData(0);
119
120 LO numNodes = vertex2AggId->getLocalLength();
121 aggSortedVertices = ArrayRCP<LO>(numNodes, -LO_ONE);
122 std::vector<LO> vertexInsertionIndexByAgg(numNodes, LO_ZERO);
123
124 for (LO i = 0; i < numNodes; ++i) {
125 LO aggId = vertex2AggIdData[i];
126 if (aggId < 0 || aggId >= numAggs) continue;
127
128 aggSortedVertices[aggsToIndices[aggId] + vertexInsertionIndexByAgg[aggId]] = i;
129 vertexInsertionIndexByAgg[aggId]++;
130 }
131}
132
133template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
134void AggregateQualityEstimateFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::ComputeAggregateQualities(RCP<const Matrix> A, RCP<const Aggregates> aggs, RCP<Xpetra::MultiVector<magnitudeType, LO, GO, Node>> agg_qualities) const {
135 const SC SCALAR_ONE = Teuchos::ScalarTraits<SC>::one();
136 const SC SCALAR_TWO = SCALAR_ONE + SCALAR_ONE;
137
138 const LO LO_ZERO = Teuchos::OrdinalTraits<LO>::zero();
139 const LO LO_ONE = Teuchos::OrdinalTraits<LO>::one();
140
141 using MT = magnitudeType;
142 const MT MT_ZERO = Teuchos::ScalarTraits<MT>::zero();
143 const MT MT_ONE = Teuchos::ScalarTraits<MT>::one();
144 ParameterList pL = GetParameterList();
145
146 RCP<const Matrix> AT = A;
147
148 // Algorithm check
149 std::string algostr = pL.get<std::string>("aggregate qualities: algorithm");
150 MT zeroThreshold = Teuchos::as<MT>(pL.get<double>("aggregate qualities: zero threshold"));
151 enum AggAlgo { ALG_FORWARD = 0,
152 ALG_REVERSE };
153 AggAlgo algo;
154 if (algostr == "forward") {
155 algo = ALG_FORWARD;
156 GetOStream(Statistics1) << "AggregateQuality: Using 'forward' algorithm" << std::endl;
157 } else if (algostr == "reverse") {
158 algo = ALG_REVERSE;
159 GetOStream(Statistics1) << "AggregateQuality: Using 'reverse' algorithm" << std::endl;
160 } else {
161 TEUCHOS_TEST_FOR_EXCEPTION(1, Exceptions::RuntimeError, "\"algorithm\" must be one of (forward|reverse)");
162 }
163
164 bool check_symmetry = pL.get<bool>("aggregate qualities: check symmetry");
165 if (check_symmetry) {
166 RCP<MultiVector> x = MultiVectorFactory::Build(A->getMap(), 1, false);
167 x->Xpetra_randomize();
168
169 RCP<MultiVector> tmp = MultiVectorFactory::Build(A->getMap(), 1, false);
170
171 A->apply(*x, *tmp, Teuchos::NO_TRANS); // tmp now stores A*x
172 A->apply(*x, *tmp, Teuchos::TRANS, -SCALAR_ONE, SCALAR_ONE); // tmp now stores A*x - A^T*x
173
174 Array<magnitudeType> tmp_norm(1);
175 tmp->norm2(tmp_norm());
176
177 Array<magnitudeType> x_norm(1);
178 tmp->norm2(x_norm());
179
180 if (tmp_norm[0] > 1e-10 * x_norm[0]) {
181 std::string transpose_string = "transpose";
182 RCP<ParameterList> whatever;
183 AT = Utilities::Transpose(*rcp_const_cast<Matrix>(A), true, transpose_string, whatever);
184
185 assert(A->getMap()->isSameAs(*(AT->getMap())));
186 }
187 }
188
189 // Reorder local aggregate information into a format amenable to computing
190 // per-aggregate quantities. Specifically, we compute a format
191 // similar to compressed sparse row format for sparse matrices in which
192 // we store all the local vertices in a single array in blocks corresponding
193 // to aggregates. (This array is aggSortedVertices.) We then store a second
194 // array (aggsToIndices) whose k-th element stores the index of the first
195 // vertex in aggregate k in the array aggSortedVertices.
196
197 ArrayRCP<LO> aggSortedVertices, aggsToIndices, aggSizes;
198 ConvertAggregatesData(aggs, aggSortedVertices, aggsToIndices, aggSizes);
199
200 LO numAggs = aggs->GetNumAggregates();
201
202 // Compute the per-aggregate quality estimate
203
204 typedef Teuchos::SerialDenseMatrix<LO, MT> DenseMatrix;
205 typedef Teuchos::SerialDenseVector<LO, MT> DenseVector;
206
207 ArrayView<const LO> rowIndices;
208 ArrayView<const SC> rowValues;
209 ArrayView<const SC> colValues;
210 Teuchos::LAPACK<LO, MT> myLapack;
211
212 // Iterate over each aggregate to compute the quality estimate
213 for (LO aggId = LO_ZERO; aggId < numAggs; ++aggId) {
214 LO aggSize = aggSizes[aggId];
215 DenseMatrix A_aggPart(aggSize, aggSize, true);
216 DenseVector offDiagonalAbsoluteSums(aggSize, true);
217
218 // Iterate over each node in the aggregate
219 for (LO idx = LO_ZERO; idx < aggSize; ++idx) {
220 LO nodeId = aggSortedVertices[idx + aggsToIndices[aggId]];
221 A->getLocalRowView(nodeId, rowIndices, rowValues);
222 AT->getLocalRowView(nodeId, rowIndices, colValues);
223
224 // Iterate over each element in the row corresponding to the current node
225 for (LO elem = LO_ZERO; elem < rowIndices.size(); ++elem) {
226 LO nodeId2 = rowIndices[elem];
227 SC val = (rowValues[elem] + colValues[elem]) / SCALAR_TWO;
228
229 LO idxInAgg = -LO_ONE; // -1 if element is not in aggregate
230
231 // Check whether the element belongs in the aggregate. If it does
232 // find, its index. Otherwise, add it's value to the off diagonal
233 // sums
234 for (LO idx2 = LO_ZERO; idx2 < aggSize; ++idx2) {
235 if (aggSortedVertices[idx2 + aggsToIndices[aggId]] == nodeId2) {
236 // Element does belong to aggregate
237 idxInAgg = idx2;
238 break;
239 }
240 }
241
242 if (idxInAgg == -LO_ONE) { // Element does not belong to aggregate
243
244 offDiagonalAbsoluteSums[idx] += Teuchos::ScalarTraits<SC>::magnitude(val);
245
246 } else { // Element does belong to aggregate
247
248 A_aggPart(idx, idxInAgg) = Teuchos::ScalarTraits<SC>::real(val);
249 }
250 }
251 }
252
253 // Construct a diagonal matrix consisting of the diagonal
254 // of A_aggPart
255 DenseMatrix A_aggPartDiagonal(aggSize, aggSize, true);
256 MT diag_sum = MT_ZERO;
257 for (int i = 0; i < aggSize; ++i) {
258 A_aggPartDiagonal(i, i) = Teuchos::ScalarTraits<SC>::real(A_aggPart(i, i));
259 diag_sum += Teuchos::ScalarTraits<SC>::real(A_aggPart(i, i));
260 }
261
262 DenseMatrix ones(aggSize, aggSize, false);
263 ones.putScalar(MT_ONE);
264
265 // Compute matrix on top of generalized Rayleigh quotient
266 // topMatrix = A_aggPartDiagonal - A_aggPartDiagonal*ones*A_aggPartDiagonal/diag_sum;
267 DenseMatrix tmp(aggSize, aggSize, false);
268 DenseMatrix topMatrix(A_aggPartDiagonal);
269
270 tmp.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, MT_ONE, ones, A_aggPartDiagonal, MT_ZERO);
271 topMatrix.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, -MT_ONE / diag_sum, A_aggPartDiagonal, tmp, MT_ONE);
272
273 // Compute matrix on bottom of generalized Rayleigh quotient
274 DenseMatrix bottomMatrix(A_aggPart);
275 MT matrixNorm = A_aggPart.normInf();
276
277 // Forward mode: Include a small perturbation to the bottom matrix to make it nonsingular
278 const MT boost = (algo == ALG_FORWARD) ? (-1e4 * Teuchos::ScalarTraits<MT>::eps() * matrixNorm) : MT_ZERO;
279
280 for (int i = 0; i < aggSize; ++i) {
281 bottomMatrix(i, i) -= offDiagonalAbsoluteSums(i) + boost;
282 }
283
284 // Compute generalized eigenvalues
285 LO sdim, info;
286 DenseVector alpha_real(aggSize, false);
287 DenseVector alpha_imag(aggSize, false);
288 DenseVector beta(aggSize, false);
289
290 DenseVector workArray(14 * (aggSize + 1), false);
291
292 LO(*ptr2func)
293 (MT*, MT*, MT*);
294 ptr2func = NULL;
295 LO* bwork = NULL;
296 MT* vl = NULL;
297 MT* vr = NULL;
298
299 const char compute_flag = 'N';
300 if (algo == ALG_FORWARD) {
301 // Forward: Solve the generalized eigenvalue problem as is
302 myLapack.GGES(compute_flag, compute_flag, compute_flag, ptr2func, aggSize,
303 topMatrix.values(), aggSize, bottomMatrix.values(), aggSize, &sdim,
304 alpha_real.values(), alpha_imag.values(), beta.values(), vl, aggSize,
305 vr, aggSize, workArray.values(), workArray.length(), bwork,
306 &info);
307 TEUCHOS_ASSERT(info == LO_ZERO);
308
309 MT maxEigenVal = MT_ZERO;
310 for (int i = LO_ZERO; i < aggSize; ++i) {
311 // NOTE: In theory, the eigenvalues should be nearly real
312 // TEUCHOS_ASSERT(fabs(alpha_imag[i]) <= 1e-8*fabs(alpha_real[i])); // Eigenvalues should be nearly real
313 maxEigenVal = std::max(maxEigenVal, alpha_real[i] / beta[i]);
314 }
315
316 (agg_qualities->getDataNonConst(0))[aggId] = (MT_ONE + MT_ONE) * maxEigenVal;
317 } else {
318 // Reverse: Swap the top and bottom matrices for the generalized eigenvalue problem
319 // This is trickier, since we need to grab the smallest non-zero eigenvalue and invert it.
320 myLapack.GGES(compute_flag, compute_flag, compute_flag, ptr2func, aggSize,
321 bottomMatrix.values(), aggSize, topMatrix.values(), aggSize, &sdim,
322 alpha_real.values(), alpha_imag.values(), beta.values(), vl, aggSize,
323 vr, aggSize, workArray.values(), workArray.length(), bwork,
324 &info);
325
326 TEUCHOS_ASSERT(info == LO_ZERO);
327
328 MT minEigenVal = MT_ZERO;
329
330 for (int i = LO_ZERO; i < aggSize; ++i) {
331 MT ev = alpha_real[i] / beta[i];
332 if (ev > zeroThreshold) {
333 if (minEigenVal == MT_ZERO)
334 minEigenVal = ev;
335 else
336 minEigenVal = std::min(minEigenVal, ev);
337 }
338 }
339 if (minEigenVal == MT_ZERO)
340 (agg_qualities->getDataNonConst(0))[aggId] = Teuchos::ScalarTraits<MT>::rmax();
341 else
342 (agg_qualities->getDataNonConst(0))[aggId] = (MT_ONE + MT_ONE) / minEigenVal;
343 }
344 } // end aggId loop
345}
346
347template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
348void AggregateQualityEstimateFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::OutputAggQualities(const Level& level, RCP<const Xpetra::MultiVector<magnitudeType, LO, GO, Node>> agg_qualities) const {
349 ParameterList pL = GetParameterList();
350
351 magnitudeType good_agg_thresh = Teuchos::as<magnitudeType>(pL.get<double>("aggregate qualities: good aggregate threshold"));
352 using MT = magnitudeType;
353
354 ArrayRCP<const MT> data = agg_qualities->getData(0);
355
356 LO num_bad_aggs = 0;
357 MT worst_agg = 0.0;
358
359 MT mean_bad_agg = 0.0;
360 MT mean_good_agg = 0.0;
361
362 for (size_t i = 0; i < agg_qualities->getLocalLength(); ++i) {
363 if (data[i] > good_agg_thresh) {
364 num_bad_aggs++;
365 mean_bad_agg += data[i];
366 } else {
367 mean_good_agg += data[i];
368 }
369 worst_agg = std::max(worst_agg, data[i]);
370 }
371
372 if (num_bad_aggs > 0) mean_bad_agg /= num_bad_aggs;
373 mean_good_agg /= agg_qualities->getLocalLength() - num_bad_aggs;
374
375 if (num_bad_aggs == 0) {
376 GetOStream(Statistics1) << "All aggregates passed the quality measure. Worst aggregate had quality " << worst_agg << ". Mean aggregate quality " << mean_good_agg << "." << std::endl;
377 } else {
378 GetOStream(Statistics1) << num_bad_aggs << " out of " << agg_qualities->getLocalLength() << " did not pass the quality measure. Worst aggregate had quality " << worst_agg << ". "
379 << "Mean bad aggregate quality " << mean_bad_agg << ". Mean good aggregate quality " << mean_good_agg << "." << std::endl;
380 }
381
382 if (pL.get<bool>("aggregate qualities: file output")) {
383 std::string filename = pL.get<std::string>("aggregate qualities: file base") + "." + std::to_string(level.GetLevelID());
384 Xpetra::IO<magnitudeType, LO, GO, Node>::Write(filename, *agg_qualities);
385 }
386
387 {
388 const auto n = size_t(agg_qualities->getLocalLength());
389
390 std::vector<MT> tmp;
391 tmp.reserve(n);
392
393 for (size_t i = 0; i < n; ++i) {
394 tmp.push_back(data[i]);
395 }
396
397 std::sort(tmp.begin(), tmp.end());
398
399 Teuchos::ArrayView<const double> percents = pL.get<Teuchos::Array<double>>("aggregate qualities: percentiles")();
400
401 GetOStream(Statistics1) << "AGG QUALITY HEADER : | LEVEL | TOTAL |";
402 for (auto percent : percents) {
403 GetOStream(Statistics1) << std::fixed << std::setprecision(4) << 100.0 * percent << "% |";
404 }
405 GetOStream(Statistics1) << std::endl;
406
407 GetOStream(Statistics1) << "AGG QUALITY PERCENTILES: | " << level.GetLevelID() << " | " << n << "|";
408 for (auto percent : percents) {
409 size_t i = size_t(n * percent);
410 i = i < n ? i : n - 1u;
411 i = i > 0u ? i : 0u;
412 GetOStream(Statistics1) << std::fixed << std::setprecision(4) << tmp[i] << " |";
413 }
414 GetOStream(Statistics1) << std::endl;
415 }
416}
417
418template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
419void AggregateQualityEstimateFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::ComputeAggregateSizes(RCP<const Matrix> A, RCP<const Aggregates> aggs, RCP<LocalOrdinalVector> agg_sizes) const {
420 ArrayRCP<LO> aggSortedVertices, aggsToIndices, aggSizes;
421 ConvertAggregatesData(aggs, aggSortedVertices, aggsToIndices, aggSizes);
422
423 // Iterate over each node in the aggregate
424 auto data = agg_sizes->getDataNonConst(0);
425 for (LO i = 0; i < (LO)aggSizes.size(); i++)
426 data[i] = aggSizes[i];
427}
428
429template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
430void AggregateQualityEstimateFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::OutputAggSizes(const Level& level, RCP<const LocalOrdinalVector> agg_sizes) const {
431 ParameterList pL = GetParameterList();
432 using MT = magnitudeType;
433
434 ArrayRCP<const LO> data = agg_sizes->getData(0);
435
436 if (pL.get<bool>("aggregate qualities: file output")) {
437 std::string filename = pL.get<std::string>("aggregate qualities: file base") + ".sizes." + std::to_string(level.GetLevelID());
438 Xpetra::IO<SC, LO, GO, Node>::WriteLOMV(filename, *agg_sizes);
439 }
440
441 {
442 size_t n = (size_t)agg_sizes->getLocalLength();
443
444 std::vector<MT> tmp;
445 tmp.reserve(n);
446
447 for (size_t i = 0; i < n; ++i) {
448 tmp.push_back(Teuchos::as<MT>(data[i]));
449 }
450
451 std::sort(tmp.begin(), tmp.end());
452
453 Teuchos::ArrayView<const double> percents = pL.get<Teuchos::Array<double>>("aggregate qualities: percentiles")();
454
455 GetOStream(Statistics1) << "AGG SIZE HEADER : | LEVEL | TOTAL |";
456 for (auto percent : percents) {
457 GetOStream(Statistics1) << std::fixed << std::setprecision(4) << 100.0 * percent << "% |";
458 }
459 GetOStream(Statistics1) << std::endl;
460
461 GetOStream(Statistics1) << "AGG SIZE PERCENTILES: | " << level.GetLevelID() << " | " << n << "|";
462 for (auto percent : percents) {
463 size_t i = size_t(n * percent);
464 i = i < n ? i : n - 1u;
465 i = i > 0u ? i : 0u;
466 GetOStream(Statistics1) << std::fixed << std::setprecision(4) << tmp[i] << " |";
467 }
468 GetOStream(Statistics1) << std::endl;
469 }
470}
471
472} // namespace MueLu
473
474#endif // MUELU_AGGREGATEQUALITYESTIMATE_DEF_HPP
#define SET_VALID_ENTRY(name)
void OutputAggQualities(const Level &level, RCP< const Xpetra::MultiVector< magnitudeType, LO, GO, Node > > agg_qualities) const
void OutputAggSizes(const Level &level, RCP< const LocalOrdinalVector > agg_sizes) const
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void ComputeAggregateQualities(RCP< const Matrix > A, RCP< const Aggregates > aggs, RCP< Xpetra::MultiVector< magnitudeType, LO, GO, Node > > agg_qualities) const
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Specifies the data that this class needs, and the factories that generate that data.
static void ConvertAggregatesData(RCP< const Aggregates > aggs, ArrayRCP< LO > &aggSortedVertices, ArrayRCP< LO > &aggsToIndices, ArrayRCP< LO > &aggSizes)
Build aggregate quality esimates with this factory.
Teuchos::ScalarTraits< Scalar >::magnitudeType magnitudeType
void Build(Level &fineLevel, Level &coarseLevel) const
Build aggregate quality esimates with this factory.
void ComputeAggregateSizes(RCP< const Matrix > A, RCP< const Aggregates > aggs, RCP< LocalOrdinalVector > agg_sizes) const
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
void Input(Level &level, const std::string &varName) const
T Get(Level &level, const std::string &varName) const
void Set(Level &level, const std::string &varName, const T &data) const
Class that holds all level-specific information.
int GetLevelID() const
Return level number.
virtual const Teuchos::ParameterList & GetParameterList() const
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Transpose(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, bool optimizeTranspose=false, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
Namespace for MueLu classes and methods.
@ Statistics1
Print more statistics.