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