MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_AggregateQualityEstimateFactory_def.hpp
Go to the documentation of this file.
1// @HEADER
2//
3// ***********************************************************************
4//
5// MueLu: A package for multigrid based preconditioning
6// Copyright 2012 Sandia Corporation
7//
8// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9// the U.S. Government retains certain rights in this software.
10//
11// Redistribution and use in source and binary forms, with or without
12// modification, are permitted provided that the following conditions are
13// met:
14//
15// 1. Redistributions of source code must retain the above copyright
16// notice, this list of conditions and the following disclaimer.
17//
18// 2. Redistributions in binary form must reproduce the above copyright
19// notice, this list of conditions and the following disclaimer in the
20// documentation and/or other materials provided with the distribution.
21//
22// 3. Neither the name of the Corporation nor the names of the
23// contributors may be used to endorse or promote products derived from
24// this software without specific prior written permission.
25//
26// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37//
38// Questions? Contact
39// Jonathan Hu (jhu@sandia.gov)
40// Andrey Prokopenko (aprokop@sandia.gov)
41// Ray Tuminaro (rstumin@sandia.gov)
42//
43// ***********************************************************************
44//
45// @HEADER
46#ifndef MUELU_AGGREGATEQUALITYESTIMATEFACTORY_DEF_HPP
47#define MUELU_AGGREGATEQUALITYESTIMATEFACTORY_DEF_HPP
48#include <iomanip>
50
51#include "MueLu_Level.hpp"
52
53#include <Teuchos_SerialDenseVector.hpp>
54#include <Teuchos_LAPACK.hpp>
55
57#include <Xpetra_IO.hpp>
58#include "MueLu_MasterList.hpp"
59#include "MueLu_Monitor.hpp"
60#include "MueLu_Utilities.hpp"
61
62#include <vector>
63
64namespace MueLu {
65
66 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
69
70 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
72
73 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
75
76 Input(currentLevel, "A");
77 Input(currentLevel, "Aggregates");
78 Input(currentLevel, "CoarseMap");
79
80 }
81
82 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
84 RCP<ParameterList> validParamList = rcp(new ParameterList());
85
86#define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
87 SET_VALID_ENTRY("aggregate qualities: good aggregate threshold");
88 SET_VALID_ENTRY("aggregate qualities: file output");
89 SET_VALID_ENTRY("aggregate qualities: file base");
90 SET_VALID_ENTRY("aggregate qualities: check symmetry");
91 SET_VALID_ENTRY("aggregate qualities: algorithm");
92 SET_VALID_ENTRY("aggregate qualities: zero threshold");
93 SET_VALID_ENTRY("aggregate qualities: percentiles");
94 SET_VALID_ENTRY("aggregate qualities: mode");
95
96#undef SET_VALID_ENTRY
97
98 validParamList->set< RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A");
99 validParamList->set< RCP<const FactoryBase> >("Aggregates", Teuchos::null, "Generating factory of the aggregates");
100 validParamList->set< RCP<const FactoryBase> >("CoarseMap", Teuchos::null, "Generating factory of the coarse map");
101
102 return validParamList;
103 }
104
105
106 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
108
109 FactoryMonitor m(*this, "Build", currentLevel);
110
111 RCP<Matrix> A = Get<RCP<Matrix>>(currentLevel, "A");
112 RCP<Aggregates> aggregates = Get<RCP<Aggregates>>(currentLevel, "Aggregates");
113
114 RCP<const Map> map = Get< RCP<const Map> >(currentLevel, "CoarseMap");
115
116
117 assert(!aggregates->AggregatesCrossProcessors());
118 ParameterList pL = GetParameterList();
119 std::string mode = pL.get<std::string>("aggregate qualities: mode");
120 GetOStream(Statistics1) << "AggregateQuality: mode "<<mode << std::endl;
121
122 RCP<Xpetra::MultiVector<magnitudeType,LO,GO,NO>> aggregate_qualities;
123 if(mode == "eigenvalue" || mode == "both") {
124 aggregate_qualities = Xpetra::MultiVectorFactory<magnitudeType,LO,GO,NO>::Build(map, 1);
125 ComputeAggregateQualities(A, aggregates, aggregate_qualities);
126 OutputAggQualities(currentLevel, aggregate_qualities);
127 }
128 if(mode == "size" || mode =="both") {
129 RCP<LocalOrdinalVector> aggregate_sizes = Xpetra::VectorFactory<LO,LO,GO,NO>::Build(map);
130 ComputeAggregateSizes(A,aggregates,aggregate_sizes);
131 Set(currentLevel, "AggregateSizes",aggregate_sizes);
132 OutputAggSizes(currentLevel, aggregate_sizes);
133 }
134 Set(currentLevel, "AggregateQualities", aggregate_qualities);
135
136
137 }
138
139 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
140 void AggregateQualityEstimateFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::ConvertAggregatesData(RCP<const Aggregates> aggs, ArrayRCP<LO>& aggSortedVertices, ArrayRCP<LO>& aggsToIndices, ArrayRCP<LO>& aggSizes) {
141
142 // Reorder local aggregate information into a format amenable to computing
143 // per-aggregate quantities. Specifically, we compute a format
144 // similar to compressed sparse row format for sparse matrices in which
145 // we store all the local vertices in a single array in blocks corresponding
146 // to aggregates. (This array is aggSortedVertices.) We then store a second
147 // array (aggsToIndices) whose k-th element stores the index of the first
148 // vertex in aggregate k in the array aggSortedVertices.
149
150 const LO LO_ZERO = Teuchos::OrdinalTraits<LO>::zero();
151 const LO LO_ONE = Teuchos::OrdinalTraits<LO>::one();
152
153 LO numAggs = aggs->GetNumAggregates();
154 aggSizes = aggs->ComputeAggregateSizesArrayRCP();
155
156 aggsToIndices = ArrayRCP<LO>(numAggs+LO_ONE,LO_ZERO);
157
158 for (LO i=0;i<numAggs;++i) {
159 aggsToIndices[i+LO_ONE] = aggsToIndices[i] + aggSizes[i];
160 }
161
162 const RCP<LOMultiVector> vertex2AggId = aggs->GetVertex2AggId();
163 const ArrayRCP<const LO> vertex2AggIdData = vertex2AggId->getData(0);
164
165 LO numNodes = vertex2AggId->getLocalLength();
166 aggSortedVertices = ArrayRCP<LO>(numNodes,-LO_ONE);
167 std::vector<LO> vertexInsertionIndexByAgg(numNodes,LO_ZERO);
168
169 for (LO i=0;i<numNodes;++i) {
170
171 LO aggId = vertex2AggIdData[i];
172 if (aggId<0 || aggId>=numAggs) continue;
173
174 aggSortedVertices[aggsToIndices[aggId]+vertexInsertionIndexByAgg[aggId]] = i;
175 vertexInsertionIndexByAgg[aggId]++;
176
177 }
178
179
180 }
181
182 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
183 void AggregateQualityEstimateFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::ComputeAggregateQualities(RCP<const Matrix> A, RCP<const Aggregates> aggs, RCP<Xpetra::MultiVector<magnitudeType,LO,GO,Node>> agg_qualities) const {
184
185 const SC SCALAR_ONE = Teuchos::ScalarTraits<SC>::one();
186 const SC SCALAR_TWO = SCALAR_ONE + SCALAR_ONE;
187
188 const LO LO_ZERO = Teuchos::OrdinalTraits<LO>::zero();
189 const LO LO_ONE = Teuchos::OrdinalTraits<LO>::one();
190
191 using MT = magnitudeType;
192 const MT MT_ZERO = Teuchos::ScalarTraits<MT>::zero();
193 const MT MT_ONE = Teuchos::ScalarTraits<MT>::one();
194 ParameterList pL = GetParameterList();
195
196 RCP<const Matrix> AT = A;
197
198 // Algorithm check
199 std::string algostr = pL.get<std::string>("aggregate qualities: algorithm");
200 MT zeroThreshold = Teuchos::as<MT>(pL.get<double>("aggregate qualities: zero threshold"));
201 enum AggAlgo {ALG_FORWARD=0, ALG_REVERSE};
202 AggAlgo algo;
203 if(algostr == "forward") {algo = ALG_FORWARD; GetOStream(Statistics1) << "AggregateQuality: Using 'forward' algorithm" << std::endl;}
204 else if(algostr == "reverse") {algo = ALG_REVERSE; GetOStream(Statistics1) << "AggregateQuality: Using 'reverse' algorithm" << std::endl;}
205 else {
206 TEUCHOS_TEST_FOR_EXCEPTION(1, Exceptions::RuntimeError, "\"algorithm\" must be one of (forward|reverse)");
207 }
208
209 bool check_symmetry = pL.get<bool>("aggregate qualities: check symmetry");
210 if (check_symmetry) {
211
212 RCP<MultiVector> x = MultiVectorFactory::Build(A->getMap(), 1, false);
213 x->Xpetra_randomize();
214
215 RCP<MultiVector> tmp = MultiVectorFactory::Build(A->getMap(), 1, false);
216
217 A->apply(*x, *tmp, Teuchos::NO_TRANS); // tmp now stores A*x
218 A->apply(*x, *tmp, Teuchos::TRANS, -SCALAR_ONE, SCALAR_ONE); // tmp now stores A*x - A^T*x
219
220 Array<magnitudeType> tmp_norm(1);
221 tmp->norm2(tmp_norm());
222
223 Array<magnitudeType> x_norm(1);
224 tmp->norm2(x_norm());
225
226 if (tmp_norm[0] > 1e-10*x_norm[0]) {
227 std::string transpose_string = "transpose";
228 RCP<ParameterList> whatever;
229 AT = Utilities::Transpose(*rcp_const_cast<Matrix>(A), true, transpose_string, whatever);
230
231 assert(A->getMap()->isSameAs( *(AT->getMap()) ));
232 }
233
234 }
235
236 // Reorder local aggregate information into a format amenable to computing
237 // per-aggregate quantities. Specifically, we compute a format
238 // similar to compressed sparse row format for sparse matrices in which
239 // we store all the local vertices in a single array in blocks corresponding
240 // to aggregates. (This array is aggSortedVertices.) We then store a second
241 // array (aggsToIndices) whose k-th element stores the index of the first
242 // vertex in aggregate k in the array aggSortedVertices.
243
244 ArrayRCP<LO> aggSortedVertices, aggsToIndices, aggSizes;
245 ConvertAggregatesData(aggs, aggSortedVertices, aggsToIndices, aggSizes);
246
247 LO numAggs = aggs->GetNumAggregates();
248
249 // Compute the per-aggregate quality estimate
250
251 typedef Teuchos::SerialDenseMatrix<LO,MT> DenseMatrix;
252 typedef Teuchos::SerialDenseVector<LO,MT> DenseVector;
253
254 ArrayView<const LO> rowIndices;
255 ArrayView<const SC> rowValues;
256 ArrayView<const SC> colValues;
257 Teuchos::LAPACK<LO,MT> myLapack;
258
259 // Iterate over each aggregate to compute the quality estimate
260 for (LO aggId=LO_ZERO; aggId<numAggs; ++aggId) {
261
262 LO aggSize = aggSizes[aggId];
263 DenseMatrix A_aggPart(aggSize, aggSize, true);
264 DenseVector offDiagonalAbsoluteSums(aggSize, true);
265
266 // Iterate over each node in the aggregate
267 for (LO idx=LO_ZERO; idx<aggSize; ++idx) {
268
269 LO nodeId = aggSortedVertices[idx+aggsToIndices[aggId]];
270 A->getLocalRowView(nodeId, rowIndices, rowValues);
271 AT->getLocalRowView(nodeId, rowIndices, colValues);
272
273 // Iterate over each element in the row corresponding to the current node
274 for (LO elem=LO_ZERO; elem<rowIndices.size();++elem) {
275
276 LO nodeId2 = rowIndices[elem];
277 SC val = (rowValues[elem]+colValues[elem])/SCALAR_TWO;
278
279 LO idxInAgg = -LO_ONE; // -1 if element is not in aggregate
280
281 // Check whether the element belongs in the aggregate. If it does
282 // find, its index. Otherwise, add it's value to the off diagonal
283 // sums
284 for (LO idx2=LO_ZERO; idx2<aggSize; ++idx2) {
285
286 if (aggSortedVertices[idx2+aggsToIndices[aggId]] == nodeId2) {
287
288 // Element does belong to aggregate
289 idxInAgg = idx2;
290 break;
291
292 }
293
294 }
295
296 if (idxInAgg == -LO_ONE) { // Element does not belong to aggregate
297
298 offDiagonalAbsoluteSums[idx] += Teuchos::ScalarTraits<SC>::magnitude(val);
299
300 } else { // Element does belong to aggregate
301
302 A_aggPart(idx,idxInAgg) = Teuchos::ScalarTraits<SC>::real(val);
303
304 }
305
306 }
307
308 }
309
310 // Construct a diagonal matrix consisting of the diagonal
311 // of A_aggPart
312 DenseMatrix A_aggPartDiagonal(aggSize, aggSize, true);
313 MT diag_sum = MT_ZERO;
314 for (int i=0;i<aggSize;++i) {
315 A_aggPartDiagonal(i,i) = Teuchos::ScalarTraits<SC>::real(A_aggPart(i,i));
316 diag_sum += Teuchos::ScalarTraits<SC>::real(A_aggPart(i,i));
317 }
318
319 DenseMatrix ones(aggSize, aggSize, false);
320 ones.putScalar(MT_ONE);
321
322 // Compute matrix on top of generalized Rayleigh quotient
323 // topMatrix = A_aggPartDiagonal - A_aggPartDiagonal*ones*A_aggPartDiagonal/diag_sum;
324 DenseMatrix tmp(aggSize, aggSize, false);
325 DenseMatrix topMatrix(A_aggPartDiagonal);
326
327 tmp.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, MT_ONE, ones, A_aggPartDiagonal, MT_ZERO);
328 topMatrix.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, -MT_ONE/diag_sum, A_aggPartDiagonal, tmp, MT_ONE);
329
330 // Compute matrix on bottom of generalized Rayleigh quotient
331 DenseMatrix bottomMatrix(A_aggPart);
332 MT matrixNorm = A_aggPart.normInf();
333
334 // Forward mode: Include a small perturbation to the bottom matrix to make it nonsingular
335 const MT boost = (algo == ALG_FORWARD) ? (-1e4*Teuchos::ScalarTraits<MT>::eps()*matrixNorm) : MT_ZERO;
336
337 for (int i=0;i<aggSize;++i){
338 bottomMatrix(i,i) -= offDiagonalAbsoluteSums(i) + boost;
339 }
340
341 // Compute generalized eigenvalues
342 LO sdim, info;
343 DenseVector alpha_real(aggSize, false);
344 DenseVector alpha_imag(aggSize, false);
345 DenseVector beta(aggSize, false);
346
347 DenseVector workArray(14*(aggSize+1), false);
348
349 LO (*ptr2func)(MT*, MT*, MT*);
350 ptr2func = NULL;
351 LO* bwork = NULL;
352 MT* vl = NULL;
353 MT* vr = NULL;
354
355 const char compute_flag ='N';
356 if(algo == ALG_FORWARD) {
357 // Forward: Solve the generalized eigenvalue problem as is
358 myLapack.GGES(compute_flag,compute_flag,compute_flag,ptr2func,aggSize,
359 topMatrix.values(),aggSize,bottomMatrix.values(),aggSize,&sdim,
360 alpha_real.values(),alpha_imag.values(),beta.values(),vl,aggSize,
361 vr,aggSize,workArray.values(),workArray.length(),bwork,
362 &info);
363 TEUCHOS_ASSERT(info == LO_ZERO);
364
365 MT maxEigenVal = MT_ZERO;
366 for (int i=LO_ZERO;i<aggSize;++i) {
367 // NOTE: In theory, the eigenvalues should be nearly real
368 //TEUCHOS_ASSERT(fabs(alpha_imag[i]) <= 1e-8*fabs(alpha_real[i])); // Eigenvalues should be nearly real
369 maxEigenVal = std::max(maxEigenVal, alpha_real[i]/beta[i]);
370 }
371
372 (agg_qualities->getDataNonConst(0))[aggId] = (MT_ONE+MT_ONE)*maxEigenVal;
373 }
374 else {
375 // Reverse: Swap the top and bottom matrices for the generalized eigenvalue problem
376 // This is trickier, since we need to grab the smallest non-zero eigenvalue and invert it.
377 myLapack.GGES(compute_flag,compute_flag,compute_flag,ptr2func,aggSize,
378 bottomMatrix.values(),aggSize,topMatrix.values(),aggSize,&sdim,
379 alpha_real.values(),alpha_imag.values(),beta.values(),vl,aggSize,
380 vr,aggSize,workArray.values(),workArray.length(),bwork,
381 &info);
382
383 TEUCHOS_ASSERT(info == LO_ZERO);
384
385 MT minEigenVal = MT_ZERO;
386
387 for (int i=LO_ZERO;i<aggSize;++i) {
388 MT ev = alpha_real[i] / beta[i];
389 if(ev > zeroThreshold) {
390 if (minEigenVal == MT_ZERO) minEigenVal = ev;
391 else minEigenVal = std::min(minEigenVal,ev);
392 }
393 }
394 if(minEigenVal == MT_ZERO) (agg_qualities->getDataNonConst(0))[aggId] = Teuchos::ScalarTraits<MT>::rmax();
395 else (agg_qualities->getDataNonConst(0))[aggId] = (MT_ONE+MT_ONE) / minEigenVal;
396 }
397 }//end aggId loop
398 }
399
400 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
401 void AggregateQualityEstimateFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::OutputAggQualities(const Level& level, RCP<const Xpetra::MultiVector<magnitudeType,LO,GO,Node>> agg_qualities) const {
402
403 ParameterList pL = GetParameterList();
404
405 magnitudeType good_agg_thresh = Teuchos::as<magnitudeType>(pL.get<double>("aggregate qualities: good aggregate threshold"));
406 using MT = magnitudeType;
407
408 ArrayRCP<const MT> data = agg_qualities->getData(0);
409
410 LO num_bad_aggs = 0;
411 MT worst_agg = 0.0;
412
413 MT mean_bad_agg = 0.0;
414 MT mean_good_agg = 0.0;
415
416
417 for (size_t i=0;i<agg_qualities->getLocalLength();++i) {
418
419 if (data[i] > good_agg_thresh) {
420 num_bad_aggs++;
421 mean_bad_agg += data[i];
422 }
423 else {
424 mean_good_agg += data[i];
425 }
426 worst_agg = std::max(worst_agg, data[i]);
427 }
428
429
430 if (num_bad_aggs > 0) mean_bad_agg /= num_bad_aggs;
431 mean_good_agg /= agg_qualities->getLocalLength() - num_bad_aggs;
432
433 if (num_bad_aggs == 0) {
434 GetOStream(Statistics1) << "All aggregates passed the quality measure. Worst aggregate had quality " << worst_agg << ". Mean aggregate quality " << mean_good_agg << "." << std::endl;
435 } else {
436 GetOStream(Statistics1) << num_bad_aggs << " out of " << agg_qualities->getLocalLength() << " did not pass the quality measure. Worst aggregate had quality " << worst_agg << ". "
437 << "Mean bad aggregate quality " << mean_bad_agg << ". Mean good aggregate quality " << mean_good_agg << "." << std::endl;
438 }
439
440 if (pL.get<bool>("aggregate qualities: file output")) {
441 std::string filename = pL.get<std::string>("aggregate qualities: file base")+"."+std::to_string(level.GetLevelID());
442 Xpetra::IO<magnitudeType,LO,GO,Node>::Write(filename, *agg_qualities);
443 }
444
445 {
446 const auto n = size_t(agg_qualities->getLocalLength());
447
448 std::vector<MT> tmp;
449 tmp.reserve(n);
450
451 for (size_t i=0; i<n; ++i) {
452 tmp.push_back(data[i]);
453 }
454
455 std::sort(tmp.begin(), tmp.end());
456
457 Teuchos::ArrayView<const double> percents = pL.get<Teuchos::Array<double> >("aggregate qualities: percentiles")();
458
459 GetOStream(Statistics1) << "AGG QUALITY HEADER : | LEVEL | TOTAL |";
460 for (auto percent : percents) {
461 GetOStream(Statistics1) << std::fixed << std::setprecision(4) <<100.0*percent << "% |";
462 }
463 GetOStream(Statistics1) << std::endl;
464
465 GetOStream(Statistics1) << "AGG QUALITY PERCENTILES: | " << level.GetLevelID() << " | " << n << "|";
466 for (auto percent : percents) {
467 size_t i = size_t(n*percent);
468 i = i < n ? i : n-1u;
469 i = i > 0u ? i : 0u;
470 GetOStream(Statistics1) << std::fixed <<std::setprecision(4) << tmp[i] << " |";
471 }
472 GetOStream(Statistics1) << std::endl;
473
474 }
475 }
476
477
478
479template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
480 void AggregateQualityEstimateFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::ComputeAggregateSizes(RCP<const Matrix> A, RCP<const Aggregates> aggs, RCP<LocalOrdinalVector> agg_sizes) const {
481
482 ArrayRCP<LO> aggSortedVertices, aggsToIndices, aggSizes;
483 ConvertAggregatesData(aggs, aggSortedVertices, aggsToIndices, aggSizes);
484
485 // Iterate over each node in the aggregate
486 auto data = agg_sizes->getDataNonConst(0);
487 for (LO i=0; i<(LO)aggSizes.size(); i++)
488 data[i] = aggSizes[i];
489}
490
491
492
493template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
494 void AggregateQualityEstimateFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::OutputAggSizes(const Level& level, RCP<const LocalOrdinalVector> agg_sizes) const {
495
496 ParameterList pL = GetParameterList();
497 using MT = magnitudeType;
498
499 ArrayRCP<const LO> data = agg_sizes->getData(0);
500
501
502 if (pL.get<bool>("aggregate qualities: file output")) {
503 std::string filename = pL.get<std::string>("aggregate qualities: file base")+".sizes."+std::to_string(level.GetLevelID());
504 Xpetra::IO<LO,LO,GO,Node>::Write(filename, *agg_sizes);
505 }
506
507 {
508 size_t n = (size_t)agg_sizes->getLocalLength();
509
510 std::vector<MT> tmp;
511 tmp.reserve(n);
512
513 for (size_t i=0; i<n; ++i) {
514 tmp.push_back(Teuchos::as<MT>(data[i]));
515 }
516
517 std::sort(tmp.begin(), tmp.end());
518
519 Teuchos::ArrayView<const double> percents = pL.get<Teuchos::Array<double> >("aggregate qualities: percentiles")();
520
521 GetOStream(Statistics1) << "AGG SIZE HEADER : | LEVEL | TOTAL |";
522 for (auto percent : percents) {
523 GetOStream(Statistics1) << std::fixed << std::setprecision(4) <<100.0*percent << "% |";
524 }
525 GetOStream(Statistics1) << std::endl;
526
527 GetOStream(Statistics1) << "AGG SIZE PERCENTILES: | " << level.GetLevelID() << " | " << n << "|";
528 for (auto percent : percents) {
529 size_t i = size_t(n*percent);
530 i = i < n ? i : n-1u;
531 i = i > 0u ? i : 0u;
532 GetOStream(Statistics1) << std::fixed <<std::setprecision(4) << tmp[i] << " |";
533 }
534 GetOStream(Statistics1) << std::endl;
535
536 }
537 }
538
539
540
541} // namespace MueLu
542
543#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
void Build(Level &currentLevel) const
Build aggregate quality esimates with this factory.
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 &currentLevel) 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 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.