MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_TentativePFactory_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_TENTATIVEPFACTORY_DEF_HPP
47 #define MUELU_TENTATIVEPFACTORY_DEF_HPP
48 
49 #include <Xpetra_MapFactory.hpp>
50 #include <Xpetra_Map.hpp>
51 #include <Xpetra_CrsMatrix.hpp>
52 #include <Xpetra_Matrix.hpp>
53 #include <Xpetra_MultiVector.hpp>
54 #include <Xpetra_MultiVectorFactory.hpp>
55 #include <Xpetra_VectorFactory.hpp>
56 #include <Xpetra_Import.hpp>
57 #include <Xpetra_ImportFactory.hpp>
58 #include <Xpetra_CrsMatrixWrap.hpp>
59 #include <Xpetra_StridedMap.hpp>
60 #include <Xpetra_StridedMapFactory.hpp>
61 
63 
64 #include "MueLu_Aggregates.hpp"
65 #include "MueLu_AmalgamationFactory.hpp"
66 #include "MueLu_AmalgamationInfo.hpp"
67 #include "MueLu_CoarseMapFactory.hpp"
68 #include "MueLu_MasterList.hpp"
69 #include "MueLu_Monitor.hpp"
70 #include "MueLu_NullspaceFactory.hpp"
71 #include "MueLu_PerfUtils.hpp"
72 #include "MueLu_Utilities.hpp"
73 
74 namespace MueLu {
75 
76  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
78  RCP<ParameterList> validParamList = rcp(new ParameterList());
79 
80 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
81  SET_VALID_ENTRY("tentative: calculate qr");
82  SET_VALID_ENTRY("tentative: build coarse coordinates");
83  SET_VALID_ENTRY("tentative: constant column sums");
84 #undef SET_VALID_ENTRY
85  validParamList->set< std::string >("Nullspace name", "Nullspace", "Name for the input nullspace");
86 
87  validParamList->set< RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A");
88  validParamList->set< RCP<const FactoryBase> >("Aggregates", Teuchos::null, "Generating factory of the aggregates");
89  validParamList->set< RCP<const FactoryBase> >("Nullspace", Teuchos::null, "Generating factory of the nullspace");
90  validParamList->set< RCP<const FactoryBase> >("Scaled Nullspace", Teuchos::null, "Generating factory of the scaled nullspace");
91  validParamList->set< RCP<const FactoryBase> >("UnAmalgamationInfo", Teuchos::null, "Generating factory of UnAmalgamationInfo");
92  validParamList->set< RCP<const FactoryBase> >("CoarseMap", Teuchos::null, "Generating factory of the coarse map");
93  validParamList->set< RCP<const FactoryBase> >("Coordinates", Teuchos::null, "Generating factory of the coordinates");
94  validParamList->set< RCP<const FactoryBase> >("Node Comm", Teuchos::null, "Generating factory of the node level communicator");
95 
96  // Make sure we don't recursively validate options for the matrixmatrix kernels
97  ParameterList norecurse;
98  norecurse.disableRecursiveValidation();
99  validParamList->set<ParameterList> ("matrixmatrix: kernel params", norecurse, "MatrixMatrix kernel parameters");
100 
101  return validParamList;
102  }
103 
104  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
106 
107  const ParameterList& pL = GetParameterList();
108  // NOTE: This guy can only either be 'Nullspace' or 'Scaled Nullspace' or else the validator above will cause issues
109  std::string nspName = "Nullspace";
110  if(pL.isParameter("Nullspace name")) nspName = pL.get<std::string>("Nullspace name");
111 
112  Input(fineLevel, "A");
113  Input(fineLevel, "Aggregates");
114  Input(fineLevel, nspName);
115  Input(fineLevel, "UnAmalgamationInfo");
116  Input(fineLevel, "CoarseMap");
117  if( fineLevel.GetLevelID() == 0 &&
118  fineLevel.IsAvailable("Coordinates", NoFactory::get()) && // we have coordinates (provided by user app)
119  pL.get<bool>("tentative: build coarse coordinates") ) { // and we want coordinates on other levels
120  bTransferCoordinates_ = true; // then set the transfer coordinates flag to true
121  Input(fineLevel, "Coordinates");
122  } else if (bTransferCoordinates_) {
123  Input(fineLevel, "Coordinates");
124  }
125  }
126 
127  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
129  return BuildP(fineLevel, coarseLevel);
130  }
131 
132  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
134  FactoryMonitor m(*this, "Build", coarseLevel);
135  typedef typename Teuchos::ScalarTraits<Scalar>::coordinateType coordinate_type;
136  typedef Xpetra::MultiVector<coordinate_type,LO,GO,NO> RealValuedMultiVector;
137  typedef Xpetra::MultiVectorFactory<coordinate_type,LO,GO,NO> RealValuedMultiVectorFactory;
138 
139  const ParameterList& pL = GetParameterList();
140  std::string nspName = "Nullspace";
141  if(pL.isParameter("Nullspace name")) nspName = pL.get<std::string>("Nullspace name");
142 
143 
144  RCP<Matrix> A = Get< RCP<Matrix> > (fineLevel, "A");
145  RCP<Aggregates> aggregates = Get< RCP<Aggregates> > (fineLevel, "Aggregates");
146  RCP<AmalgamationInfo> amalgInfo = Get< RCP<AmalgamationInfo> > (fineLevel, "UnAmalgamationInfo");
147  RCP<MultiVector> fineNullspace = Get< RCP<MultiVector> > (fineLevel, nspName);
148  RCP<const Map> coarseMap = Get< RCP<const Map> > (fineLevel, "CoarseMap");
149  RCP<RealValuedMultiVector> fineCoords;
150  if(bTransferCoordinates_) {
151  fineCoords = Get< RCP<RealValuedMultiVector> >(fineLevel, "Coordinates");
152  }
153 
154  // FIXME: We should remove the NodeComm on levels past the threshold
155  if(fineLevel.IsAvailable("Node Comm")) {
156  RCP<const Teuchos::Comm<int> > nodeComm = Get<RCP<const Teuchos::Comm<int> > >(fineLevel,"Node Comm");
157  Set<RCP<const Teuchos::Comm<int> > >(coarseLevel, "Node Comm", nodeComm);
158  }
159 
160  TEUCHOS_TEST_FOR_EXCEPTION(A->getRowMap()->getNodeNumElements() != fineNullspace->getMap()->getNodeNumElements(),
161  Exceptions::RuntimeError,"MueLu::TentativePFactory::MakeTentative: Size mismatch between A and Nullspace");
162 
163  RCP<Matrix> Ptentative;
164  RCP<MultiVector> coarseNullspace;
165  RCP<RealValuedMultiVector> coarseCoords;
166 
167  if(bTransferCoordinates_) {
168  //*** Create the coarse coordinates ***
169  // First create the coarse map and coarse multivector
170  ArrayView<const GO> elementAList = coarseMap->getNodeElementList();
171  LO blkSize = 1;
172  if (rcp_dynamic_cast<const StridedMap>(coarseMap) != Teuchos::null) {
173  blkSize = rcp_dynamic_cast<const StridedMap>(coarseMap)->getFixedBlockSize();
174  }
175  GO indexBase = coarseMap->getIndexBase();
176  LO numCoarseNodes = Teuchos::as<LO>(elementAList.size() / blkSize);
177  Array<GO> nodeList(numCoarseNodes);
178  const int numDimensions = fineCoords->getNumVectors();
179 
180  for (LO i = 0; i < numCoarseNodes; i++) {
181  nodeList[i] = (elementAList[i*blkSize]-indexBase)/blkSize + indexBase;
182  }
183  RCP<const Map> coarseCoordsMap = MapFactory::Build(fineCoords->getMap()->lib(),
185  nodeList,
186  indexBase,
187  fineCoords->getMap()->getComm());
188  coarseCoords = RealValuedMultiVectorFactory::Build(coarseCoordsMap, numDimensions);
189 
190  // Create overlapped fine coordinates to reduce global communication
191  RCP<RealValuedMultiVector> ghostedCoords;
192  if (aggregates->AggregatesCrossProcessors()) {
193  RCP<const Map> aggMap = aggregates->GetMap();
194  RCP<const Import> importer = ImportFactory::Build(fineCoords->getMap(), aggMap);
195 
196  ghostedCoords = RealValuedMultiVectorFactory::Build(aggMap, numDimensions);
197  ghostedCoords->doImport(*fineCoords, *importer, Xpetra::INSERT);
198  } else {
199  ghostedCoords = fineCoords;
200  }
201 
202  // Get some info about aggregates
203  int myPID = coarseCoordsMap->getComm()->getRank();
204  LO numAggs = aggregates->GetNumAggregates();
205  ArrayRCP<LO> aggSizes = aggregates->ComputeAggregateSizes();
206  const ArrayRCP<const LO> vertex2AggID = aggregates->GetVertex2AggId()->getData(0);
207  const ArrayRCP<const LO> procWinner = aggregates->GetProcWinner()->getData(0);
208 
209  // Fill in coarse coordinates
210  for (int dim = 0; dim < numDimensions; ++dim) {
211  ArrayRCP<const coordinate_type> fineCoordsData = ghostedCoords->getData(dim);
212  ArrayRCP<coordinate_type> coarseCoordsData = coarseCoords->getDataNonConst(dim);
213 
214  for (LO lnode = 0; lnode < Teuchos::as<LO>(vertex2AggID.size()); lnode++) {
215  if (procWinner[lnode] == myPID &&
216  lnode < fineCoordsData.size() &&
217  vertex2AggID[lnode] < coarseCoordsData.size() &&
218  Teuchos::ScalarTraits<coordinate_type>::isnaninf(fineCoordsData[lnode]) == false) {
219  coarseCoordsData[vertex2AggID[lnode]] += fineCoordsData[lnode];
220  }
221  }
222  for (LO agg = 0; agg < numAggs; agg++) {
223  coarseCoordsData[agg] /= aggSizes[agg];
224  }
225  }
226  }
227 
228  if (!aggregates->AggregatesCrossProcessors())
229  BuildPuncoupled(A, aggregates, amalgInfo, fineNullspace, coarseMap, Ptentative, coarseNullspace,coarseLevel.GetLevelID());
230  else
231  BuildPcoupled (A, aggregates, amalgInfo, fineNullspace, coarseMap, Ptentative, coarseNullspace);
232 
233  // If available, use striding information of fine level matrix A for range
234  // map and coarseMap as domain map; otherwise use plain range map of
235  // Ptent = plain range map of A for range map and coarseMap as domain map.
236  // NOTE:
237  // The latter is not really safe, since there is no striding information
238  // for the range map. This is not really a problem, since striding
239  // information is always available on the intermedium levels and the
240  // coarsest levels.
241  if (A->IsView("stridedMaps") == true)
242  Ptentative->CreateView("stridedMaps", A->getRowMap("stridedMaps"), coarseMap);
243  else
244  Ptentative->CreateView("stridedMaps", Ptentative->getRangeMap(), coarseMap);
245 
246  if(bTransferCoordinates_) {
247  Set(coarseLevel, "Coordinates", coarseCoords);
248  }
249  Set(coarseLevel, "Nullspace", coarseNullspace);
250  Set(coarseLevel, "P", Ptentative);
251 
252  if (IsPrint(Statistics2)) {
253  RCP<ParameterList> params = rcp(new ParameterList());
254  params->set("printLoadBalancingInfo", true);
255  GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*Ptentative, "Ptent", params);
256  }
257  }
258 
259  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
262  RCP<const Map> coarseMap, RCP<Matrix>& Ptentative, RCP<MultiVector>& coarseNullspace, const int levelID) const {
263  RCP<const Map> rowMap = A->getRowMap();
264  RCP<const Map> colMap = A->getColMap();
265 
266  const size_t numRows = rowMap->getNodeNumElements();
267 
268  typedef Teuchos::ScalarTraits<SC> STS;
269  typedef typename STS::magnitudeType Magnitude;
270  const SC zero = STS::zero();
271  const SC one = STS::one();
272  const LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
273 
274  const GO numAggs = aggregates->GetNumAggregates();
275  const size_t NSDim = fineNullspace->getNumVectors();
276  ArrayRCP<LO> aggSizes = aggregates->ComputeAggregateSizes();
277 
278 
279  // Sanity checking
280  const ParameterList& pL = GetParameterList();
281  const bool &doQRStep = pL.get<bool>("tentative: calculate qr");
282  const bool &constantColSums = pL.get<bool>("tentative: constant column sums");
283 
284  TEUCHOS_TEST_FOR_EXCEPTION(doQRStep && constantColSums,Exceptions::RuntimeError,
285  "MueLu::TentativePFactory::MakeTentative: cannot use 'constant column sums' and 'calculate qr' at the same time");
286 
287  // Aggregates map is based on the amalgamated column map
288  // We can skip global-to-local conversion if LIDs in row map are
289  // same as LIDs in column map
290  bool goodMap = isGoodMap(*rowMap, *colMap);
291 
292  // Create a lookup table to determine the rows (fine DOFs) that belong to a given aggregate.
293  // aggStart is a pointer into aggToRowMapLO
294  // aggStart[i]..aggStart[i+1] are indices into aggToRowMapLO
295  // aggToRowMapLO[aggStart[i]]..aggToRowMapLO[aggStart[i+1]-1] are the DOFs in aggregate i
296  ArrayRCP<LO> aggStart;
297  ArrayRCP<LO> aggToRowMapLO;
298  ArrayRCP<GO> aggToRowMapGO;
299  if (goodMap) {
300  amalgInfo->UnamalgamateAggregatesLO(*aggregates, aggStart, aggToRowMapLO);
301  GetOStream(Runtime1) << "Column map is consistent with the row map, good." << std::endl;
302 
303  } else {
304  amalgInfo->UnamalgamateAggregates(*aggregates, aggStart, aggToRowMapGO);
305  GetOStream(Warnings0) << "Column map is not consistent with the row map\n"
306  << "using GO->LO conversion with performance penalty" << std::endl;
307  }
308 
309  coarseNullspace = MultiVectorFactory::Build(coarseMap, NSDim);
310 
311  // Pull out the nullspace vectors so that we can have random access.
312  ArrayRCP<ArrayRCP<const SC> > fineNS (NSDim);
313  ArrayRCP<ArrayRCP<SC> > coarseNS(NSDim);
314  for (size_t i = 0; i < NSDim; i++) {
315  fineNS[i] = fineNullspace->getData(i);
316  if (coarseMap->getNodeNumElements() > 0)
317  coarseNS[i] = coarseNullspace->getDataNonConst(i);
318  }
319 
320  size_t nnzEstimate = numRows * NSDim;
321 
322  // Time to construct the matrix and fill in the values
323  Ptentative = rcp(new CrsMatrixWrap(rowMap, coarseMap, 0));
324  RCP<CrsMatrix> PtentCrs = rcp_dynamic_cast<CrsMatrixWrap>(Ptentative)->getCrsMatrix();
325 
326  ArrayRCP<size_t> iaPtent;
327  ArrayRCP<LO> jaPtent;
328  ArrayRCP<SC> valPtent;
329 
330  PtentCrs->allocateAllValues(nnzEstimate, iaPtent, jaPtent, valPtent);
331 
332  ArrayView<size_t> ia = iaPtent();
333  ArrayView<LO> ja = jaPtent();
334  ArrayView<SC> val = valPtent();
335 
336  ia[0] = 0;
337  for (size_t i = 1; i <= numRows; i++)
338  ia[i] = ia[i-1] + NSDim;
339 
340  for (size_t j = 0; j < nnzEstimate; j++) {
341  ja [j] = INVALID;
342  val[j] = zero;
343  }
344 
345 
346  if (doQRStep) {
348  // Standard aggregate-wise QR //
350  for (GO agg = 0; agg < numAggs; agg++) {
351  LO aggSize = aggStart[agg+1] - aggStart[agg];
352 
353  Xpetra::global_size_t offset = agg*NSDim;
354 
355  // Extract the piece of the nullspace corresponding to the aggregate, and
356  // put it in the flat array, "localQR" (in column major format) for the
357  // QR routine.
358  Teuchos::SerialDenseMatrix<LO,SC> localQR(aggSize, NSDim);
359  if (goodMap) {
360  for (size_t j = 0; j < NSDim; j++)
361  for (LO k = 0; k < aggSize; k++)
362  localQR(k,j) = fineNS[j][aggToRowMapLO[aggStart[agg]+k]];
363  } else {
364  for (size_t j = 0; j < NSDim; j++)
365  for (LO k = 0; k < aggSize; k++)
366  localQR(k,j) = fineNS[j][rowMap->getLocalElement(aggToRowMapGO[aggStart[agg]+k])];
367  }
368 
369  // Test for zero columns
370  for (size_t j = 0; j < NSDim; j++) {
371  bool bIsZeroNSColumn = true;
372 
373  for (LO k = 0; k < aggSize; k++)
374  if (localQR(k,j) != zero)
375  bIsZeroNSColumn = false;
376 
377  TEUCHOS_TEST_FOR_EXCEPTION(bIsZeroNSColumn == true, Exceptions::RuntimeError,
378  "MueLu::TentativePFactory::MakeTentative: fine level NS part has a zero column");
379  }
380 
381  // Calculate QR decomposition (standard)
382  // NOTE: Q is stored in localQR and R is stored in coarseNS
383  if (aggSize >= Teuchos::as<LO>(NSDim)) {
384 
385  if (NSDim == 1) {
386  // Only one nullspace vector, calculate Q and R by hand
387  Magnitude norm = STS::magnitude(zero);
388  for (size_t k = 0; k < Teuchos::as<size_t>(aggSize); k++)
389  norm += STS::magnitude(localQR(k,0)*localQR(k,0));
391 
392  // R = norm
393  coarseNS[0][offset] = norm;
394 
395  // Q = localQR(:,0)/norm
396  for (LO i = 0; i < aggSize; i++)
397  localQR(i,0) /= norm;
398 
399  } else {
401  qrSolver.setMatrix(Teuchos::rcp(&localQR, false));
402  qrSolver.factor();
403 
404  // R = upper triangular part of localQR
405  for (size_t j = 0; j < NSDim; j++)
406  for (size_t k = 0; k <= j; k++)
407  coarseNS[j][offset+k] = localQR(k,j); //TODO is offset+k the correct local ID?!
408 
409  // Calculate Q, the tentative prolongator.
410  // The Lapack GEQRF call only works for myAggsize >= NSDim
411  qrSolver.formQ();
413  for (size_t j = 0; j < NSDim; j++)
414  for (size_t i = 0; i < Teuchos::as<size_t>(aggSize); i++)
415  localQR(i,j) = (*qFactor)(i,j);
416  }
417 
418  } else {
419  // Special handling for aggSize < NSDim (i.e. single node aggregates in structural mechanics)
420 
421  // The local QR decomposition is not possible in the "overconstrained"
422  // case (i.e. number of columns in localQR > number of rows), which
423  // corresponds to #DOFs in Aggregate < NSDim. For usual problems this
424  // is only possible for single node aggregates in structural mechanics.
425  // (Similar problems may arise in discontinuous Galerkin problems...)
426  // We bypass the QR decomposition and use an identity block in the
427  // tentative prolongator for the single node aggregate and transfer the
428  // corresponding fine level null space information 1-to-1 to the coarse
429  // level null space part.
430 
431  // NOTE: The resulting tentative prolongation operator has
432  // (aggSize*DofsPerNode-NSDim) zero columns leading to a singular
433  // coarse level operator A. To deal with that one has the following
434  // options:
435  // - Use the "RepairMainDiagonal" flag in the RAPFactory (default:
436  // false) to add some identity block to the diagonal of the zero rows
437  // in the coarse level operator A, such that standard level smoothers
438  // can be used again.
439  // - Use special (projection-based) level smoothers, which can deal
440  // with singular matrices (very application specific)
441  // - Adapt the code below to avoid zero columns. However, we do not
442  // support a variable number of DOFs per node in MueLu/Xpetra which
443  // makes the implementation really hard.
444 
445  // R = extended (by adding identity rows) localQR
446  for (size_t j = 0; j < NSDim; j++)
447  for (size_t k = 0; k < NSDim; k++)
448  if (k < as<size_t>(aggSize))
449  coarseNS[j][offset+k] = localQR(k,j);
450  else
451  coarseNS[j][offset+k] = (k == j ? one : zero);
452 
453  // Q = I (rectangular)
454  for (size_t i = 0; i < as<size_t>(aggSize); i++)
455  for (size_t j = 0; j < NSDim; j++)
456  localQR(i,j) = (j == i ? one : zero);
457  }
458 
459 
460  // Process each row in the local Q factor
461  // FIXME: What happens if maps are blocked?
462  for (LO j = 0; j < aggSize; j++) {
463  LO localRow = (goodMap ? aggToRowMapLO[aggStart[agg]+j] : rowMap->getLocalElement(aggToRowMapGO[aggStart[agg]+j]));
464 
465  size_t rowStart = ia[localRow];
466  for (size_t k = 0, lnnz = 0; k < NSDim; k++) {
467  // Skip zeros (there may be plenty of them, i.e., NSDim > 1 or boundary conditions)
468  if (localQR(j,k) != zero) {
469  ja [rowStart+lnnz] = offset + k;
470  val[rowStart+lnnz] = localQR(j,k);
471  lnnz++;
472  }
473  }
474  }
475  }
476 
477  } else {
478  GetOStream(Runtime1) << "TentativePFactory : bypassing local QR phase" << std::endl;
479  if (NSDim>1)
480  GetOStream(Warnings0) << "TentativePFactory : for nontrivial nullspace, this may degrade performance" << std::endl;
482  // "no-QR" option //
484  // Local Q factor is just the fine nullspace support over the current aggregate.
485  // Local R factor is the identity.
486  // TODO I have not implemented any special handling for aggregates that are too
487  // TODO small to locally support the nullspace, as is done in the standard QR
488  // TODO case above.
489  if (goodMap) {
490  for (GO agg = 0; agg < numAggs; agg++) {
491  const LO aggSize = aggStart[agg+1] - aggStart[agg];
492  Xpetra::global_size_t offset = agg*NSDim;
493 
494  // Process each row in the local Q factor
495  // FIXME: What happens if maps are blocked?
496  for (LO j = 0; j < aggSize; j++) {
497 
498  //TODO Here I do not check for a zero nullspace column on the aggregate.
499  // as is done in the standard QR case.
500 
501  const LO localRow = aggToRowMapLO[aggStart[agg]+j];
502 
503  const size_t rowStart = ia[localRow];
504 
505  for (size_t k = 0, lnnz = 0; k < NSDim; k++) {
506  // Skip zeros (there may be plenty of them, i.e., NSDim > 1 or boundary conditions)
507  SC qr_jk = fineNS[k][aggToRowMapLO[aggStart[agg]+j]];
508  if(constantColSums) qr_jk = qr_jk / (double)aggSizes[agg];
509  if (qr_jk != zero) {
510  ja [rowStart+lnnz] = offset + k;
511  val[rowStart+lnnz] = qr_jk;
512  lnnz++;
513  }
514  }
515  }
516  for (size_t j = 0; j < NSDim; j++)
517  coarseNS[j][offset+j] = one;
518  } //for (GO agg = 0; agg < numAggs; agg++)
519 
520  } else {
521  for (GO agg = 0; agg < numAggs; agg++) {
522  const LO aggSize = aggStart[agg+1] - aggStart[agg];
523  Xpetra::global_size_t offset = agg*NSDim;
524  for (LO j = 0; j < aggSize; j++) {
525 
526  const LO localRow = rowMap->getLocalElement(aggToRowMapGO[aggStart[agg]+j]);
527 
528  const size_t rowStart = ia[localRow];
529 
530  for (size_t k = 0, lnnz = 0; k < NSDim; k++) {
531  // Skip zeros (there may be plenty of them, i.e., NSDim > 1 or boundary conditions)
532  SC qr_jk = fineNS[k][rowMap->getLocalElement(aggToRowMapGO[aggStart[agg]+j])];
533  if(constantColSums) qr_jk = qr_jk / (double)aggSizes[agg];
534  if (qr_jk != zero) {
535  ja [rowStart+lnnz] = offset + k;
536  val[rowStart+lnnz] = qr_jk;
537  lnnz++;
538  }
539  }
540  }
541  for (size_t j = 0; j < NSDim; j++)
542  coarseNS[j][offset+j] = one;
543  } //for (GO agg = 0; agg < numAggs; agg++)
544 
545  } //if (goodmap) else ...
546 
547  } //if doQRStep ... else
548 
549  // Compress storage (remove all INVALID, which happen when we skip zeros)
550  // We do that in-place
551  size_t ia_tmp = 0, nnz = 0;
552  for (size_t i = 0; i < numRows; i++) {
553  for (size_t j = ia_tmp; j < ia[i+1]; j++)
554  if (ja[j] != INVALID) {
555  ja [nnz] = ja [j];
556  val[nnz] = val[j];
557  nnz++;
558  }
559  ia_tmp = ia[i+1];
560  ia[i+1] = nnz;
561  }
562  if (rowMap->lib() == Xpetra::UseTpetra) {
563  // - Cannot resize for Epetra, as it checks for same pointers
564  // - Need to resize for Tpetra, as it check ().size() == ia[numRows]
565  // NOTE: these invalidate ja and val views
566  jaPtent .resize(nnz);
567  valPtent.resize(nnz);
568  }
569 
570  GetOStream(Runtime1) << "TentativePFactory : aggregates do not cross process boundaries" << std::endl;
571 
572  PtentCrs->setAllValues(iaPtent, jaPtent, valPtent);
573 
574 
575  // Managing labels & constants for ESFC
576  RCP<ParameterList> FCparams;
577  if(pL.isSublist("matrixmatrix: kernel params"))
578  FCparams=rcp(new ParameterList(pL.sublist("matrixmatrix: kernel params")));
579  else
580  FCparams= rcp(new ParameterList);
581  // By default, we don't need global constants for TentativeP
582  FCparams->set("compute global constants",FCparams->get("compute global constants",false));
583  std::string levelIDs = toString(levelID);
584  FCparams->set("Timer Label",std::string("MueLu::TentativeP-")+levelIDs);
585  RCP<const Export> dummy_e;
586  RCP<const Import> dummy_i;
587 
588  PtentCrs->expertStaticFillComplete(coarseMap, A->getDomainMap(),dummy_i,dummy_e,FCparams);
589  }
590 
591  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
594  RCP<const Map> coarseMap, RCP<Matrix>& Ptentative, RCP<MultiVector>& coarseNullspace) const {
595  typedef Teuchos::ScalarTraits<SC> STS;
596  typedef typename STS::magnitudeType Magnitude;
597  const SC zero = STS::zero();
598  const SC one = STS::one();
599 
600  // number of aggregates
601  GO numAggs = aggregates->GetNumAggregates();
602 
603  // Create a lookup table to determine the rows (fine DOFs) that belong to a given aggregate.
604  // aggStart is a pointer into aggToRowMap
605  // aggStart[i]..aggStart[i+1] are indices into aggToRowMap
606  // aggToRowMap[aggStart[i]]..aggToRowMap[aggStart[i+1]-1] are the DOFs in aggregate i
607  ArrayRCP<LO> aggStart;
608  ArrayRCP< GO > aggToRowMap;
609  amalgInfo->UnamalgamateAggregates(*aggregates, aggStart, aggToRowMap);
610 
611  // find size of the largest aggregate
612  LO maxAggSize=0;
613  for (GO i=0; i<numAggs; ++i) {
614  LO sizeOfThisAgg = aggStart[i+1] - aggStart[i];
615  if (sizeOfThisAgg > maxAggSize) maxAggSize = sizeOfThisAgg;
616  }
617 
618  // dimension of fine level nullspace
619  const size_t NSDim = fineNullspace->getNumVectors();
620 
621  // index base for coarse Dof map (usually 0)
622  GO indexBase=A->getRowMap()->getIndexBase();
623 
624  const RCP<const Map> nonUniqueMap = amalgInfo->ComputeUnamalgamatedImportDofMap(*aggregates);
625  const RCP<const Map> uniqueMap = A->getDomainMap();
626  RCP<const Import> importer = ImportFactory::Build(uniqueMap, nonUniqueMap);
627  RCP<MultiVector> fineNullspaceWithOverlap = MultiVectorFactory::Build(nonUniqueMap,NSDim);
628  fineNullspaceWithOverlap->doImport(*fineNullspace,*importer,Xpetra::INSERT);
629 
630  // Pull out the nullspace vectors so that we can have random access.
631  ArrayRCP< ArrayRCP<const SC> > fineNS(NSDim);
632  for (size_t i=0; i<NSDim; ++i)
633  fineNS[i] = fineNullspaceWithOverlap->getData(i);
634 
635  //Allocate storage for the coarse nullspace.
636  coarseNullspace = MultiVectorFactory::Build(coarseMap, NSDim);
637 
638  ArrayRCP< ArrayRCP<SC> > coarseNS(NSDim);
639  for (size_t i=0; i<NSDim; ++i)
640  if (coarseMap->getNodeNumElements() > 0) coarseNS[i] = coarseNullspace->getDataNonConst(i);
641 
642  //This makes the rowmap of Ptent the same as that of A->
643  //This requires moving some parts of some local Q's to other processors
644  //because aggregates can span processors.
645  RCP<const Map > rowMapForPtent = A->getRowMap();
646  const Map& rowMapForPtentRef = *rowMapForPtent;
647 
648  // Set up storage for the rows of the local Qs that belong to other processors.
649  // FIXME This is inefficient and could be done within the main loop below with std::vector's.
650  RCP<const Map> colMap = A->getColMap();
651 
652  RCP<const Map > ghostQMap;
653  RCP<MultiVector> ghostQvalues;
654  Array<RCP<Xpetra::Vector<GO,LO,GO,Node> > > ghostQcolumns;
655  RCP<Xpetra::Vector<GO,LO,GO,Node> > ghostQrowNums;
656  ArrayRCP< ArrayRCP<SC> > ghostQvals;
657  ArrayRCP< ArrayRCP<GO> > ghostQcols;
658  ArrayRCP< GO > ghostQrows;
659 
660  Array<GO> ghostGIDs;
661  for (LO j=0; j<numAggs; ++j) {
662  for (LO k=aggStart[j]; k<aggStart[j+1]; ++k) {
663  if (rowMapForPtentRef.isNodeGlobalElement(aggToRowMap[k]) == false) {
664  ghostGIDs.push_back(aggToRowMap[k]);
665  }
666  }
667  }
668  ghostQMap = MapFactory::Build(A->getRowMap()->lib(),
670  ghostGIDs,
671  indexBase, A->getRowMap()->getComm()); //JG:Xpetra::global_size_t>?
672  //Vector to hold bits of Q that go to other processors.
673  ghostQvalues = MultiVectorFactory::Build(ghostQMap,NSDim);
674  //Note that Epetra does not support MultiVectors templated on Scalar != double.
675  //So to work around this, we allocate an array of Vectors. This shouldn't be too
676  //expensive, as the number of Vectors is NSDim.
677  ghostQcolumns.resize(NSDim);
678  for (size_t i=0; i<NSDim; ++i)
679  ghostQcolumns[i] = Xpetra::VectorFactory<GO,LO,GO,Node>::Build(ghostQMap);
680  ghostQrowNums = Xpetra::VectorFactory<GO,LO,GO,Node>::Build(ghostQMap);
681  if (ghostQvalues->getLocalLength() > 0) {
682  ghostQvals.resize(NSDim);
683  ghostQcols.resize(NSDim);
684  for (size_t i=0; i<NSDim; ++i) {
685  ghostQvals[i] = ghostQvalues->getDataNonConst(i);
686  ghostQcols[i] = ghostQcolumns[i]->getDataNonConst(0);
687  }
688  ghostQrows = ghostQrowNums->getDataNonConst(0);
689  }
690 
691  //importer to handle moving Q
692  importer = ImportFactory::Build(ghostQMap, A->getRowMap());
693 
694  // Dense QR solver
696 
697  //Allocate temporary storage for the tentative prolongator.
698  Array<GO> globalColPtr(maxAggSize*NSDim,0);
699  Array<LO> localColPtr(maxAggSize*NSDim,0);
700  Array<SC> valPtr(maxAggSize*NSDim,0.);
701 
702  //Create column map for Ptent, estimate local #nonzeros in Ptent, and create Ptent itself.
703  const Map& coarseMapRef = *coarseMap;
704 
705  // For the 3-arrays constructor
706  ArrayRCP<size_t> ptent_rowptr;
707  ArrayRCP<LO> ptent_colind;
708  ArrayRCP<Scalar> ptent_values;
709 
710  // Because ArrayRCPs are slow...
711  ArrayView<size_t> rowptr_v;
712  ArrayView<LO> colind_v;
713  ArrayView<Scalar> values_v;
714 
715  // For temporary usage
716  Array<size_t> rowptr_temp;
717  Array<LO> colind_temp;
718  Array<Scalar> values_temp;
719 
720  RCP<CrsMatrix> PtentCrs;
721 
722  RCP<CrsMatrixWrap> PtentCrsWrap = rcp(new CrsMatrixWrap(rowMapForPtent, NSDim));
723  PtentCrs = PtentCrsWrap->getCrsMatrix();
724  Ptentative = PtentCrsWrap;
725 
726  //*****************************************************************
727  //Loop over all aggregates and calculate local QR decompositions.
728  //*****************************************************************
729  GO qctr=0; //for indexing into Ptent data vectors
730  const Map& nonUniqueMapRef = *nonUniqueMap;
731 
732  size_t total_nnz_count=0;
733 
734  for (GO agg=0; agg<numAggs; ++agg)
735  {
736  LO myAggSize = aggStart[agg+1]-aggStart[agg];
737  // For each aggregate, extract the corresponding piece of the nullspace and put it in the flat array,
738  // "localQR" (in column major format) for the QR routine.
739  Teuchos::SerialDenseMatrix<LO,SC> localQR(myAggSize, NSDim);
740  for (size_t j=0; j<NSDim; ++j) {
741  bool bIsZeroNSColumn = true;
742  for (LO k=0; k<myAggSize; ++k)
743  {
744  // aggToRowMap[aggPtr[i]+k] is the kth DOF in the ith aggregate
745  // fineNS[j][n] is the nth entry in the jth NS vector
746  try{
747  SC nsVal = fineNS[j][ nonUniqueMapRef.getLocalElement(aggToRowMap[aggStart[agg]+k]) ]; // extract information from fine level NS
748  localQR(k,j) = nsVal;
749  if (nsVal != zero) bIsZeroNSColumn = false;
750  }
751  catch(...) {
752  GetOStream(Runtime1,-1) << "length of fine level nsp: " << fineNullspace->getGlobalLength() << std::endl;
753  GetOStream(Runtime1,-1) << "length of fine level nsp w overlap: " << fineNullspaceWithOverlap->getGlobalLength() << std::endl;
754  GetOStream(Runtime1,-1) << "(local?) aggId=" << agg << std::endl;
755  GetOStream(Runtime1,-1) << "aggSize=" << myAggSize << std::endl;
756  GetOStream(Runtime1,-1) << "agg DOF=" << k << std::endl;
757  GetOStream(Runtime1,-1) << "NS vector j=" << j << std::endl;
758  GetOStream(Runtime1,-1) << "j*myAggSize + k = " << j*myAggSize + k << std::endl;
759  GetOStream(Runtime1,-1) << "aggToRowMap["<<agg<<"][" << k << "] = " << aggToRowMap[aggStart[agg]+k] << std::endl;
760  GetOStream(Runtime1,-1) << "id aggToRowMap[agg][k]=" << aggToRowMap[aggStart[agg]+k] << " is global element in nonUniqueMap = " <<
761  nonUniqueMapRef.isNodeGlobalElement(aggToRowMap[aggStart[agg]+k]) << std::endl;
762  GetOStream(Runtime1,-1) << "colMap local id aggToRowMap[agg][k]=" << nonUniqueMapRef.getLocalElement(aggToRowMap[aggStart[agg]+k]) << std::endl;
763  GetOStream(Runtime1,-1) << "fineNS...=" << fineNS[j][ nonUniqueMapRef.getLocalElement(aggToRowMap[aggStart[agg]+k]) ] << std::endl;
764  GetOStream(Errors,-1) << "caught an error!" << std::endl;
765  }
766  } //for (LO k=0 ...
767  TEUCHOS_TEST_FOR_EXCEPTION(bIsZeroNSColumn == true, Exceptions::RuntimeError, "MueLu::TentativePFactory::MakeTentative: fine level NS part has a zero column. Error.");
768  } //for (LO j=0 ...
769 
770  Xpetra::global_size_t offset=agg*NSDim;
771 
772  if(myAggSize >= Teuchos::as<LocalOrdinal>(NSDim)) {
773  // calculate QR decomposition (standard)
774  // R is stored in localQR (size: myAggSize x NSDim)
775 
776  // Householder multiplier
777  SC tau = localQR(0,0);
778 
779  if (NSDim == 1) {
780  // Only one nullspace vector, so normalize by hand
781  Magnitude dtemp=0;
782  for (size_t k = 0; k < Teuchos::as<size_t>(myAggSize); ++k) {
783  Magnitude tmag = STS::magnitude(localQR(k,0));
784  dtemp += tmag*tmag;
785  }
787  tau = localQR(0,0);
788  localQR(0,0) = dtemp;
789  } else {
790  qrSolver.setMatrix( Teuchos::rcp(&localQR, false) );
791  qrSolver.factor();
792  }
793 
794  // Extract R, the coarse nullspace. This is stored in upper triangular part of localQR.
795  // Note: coarseNS[i][.] is the ith coarse nullspace vector, which may be counter to your intuition.
796  // This stores the (offset+k)th entry only if it is local according to the coarseMap.
797  for (size_t j=0; j<NSDim; ++j) {
798  for (size_t k=0; k<=j; ++k) {
799  try {
800  if (coarseMapRef.isNodeLocalElement(offset+k)) {
801  coarseNS[j][offset+k] = localQR(k, j); //TODO is offset+k the correct local ID?!
802  }
803  }
804  catch(...) {
805  GetOStream(Errors,-1) << "caught error in coarseNS insert, j="<<j<<", offset+k = "<<offset+k<<std::endl;
806  }
807  }
808  }
809 
810  // Calculate Q, the tentative prolongator.
811  // The Lapack GEQRF call only works for myAggsize >= NSDim
812 
813  if (NSDim == 1) {
814  // Only one nullspace vector, so calculate Q by hand
815  Magnitude dtemp = Teuchos::ScalarTraits<SC>::magnitude(localQR(0,0));
816  localQR(0,0) = tau;
817  dtemp = 1 / dtemp;
818  for (LocalOrdinal i=0; i<myAggSize; ++i) {
819  localQR(i,0) *= dtemp ;
820  }
821  } else {
822  qrSolver.formQ();
823  Teuchos::RCP<Teuchos::SerialDenseMatrix<LO,SC> > qFactor = qrSolver.getQ();
824  for (size_t j=0; j<NSDim; j++) {
825  for (size_t i = 0; i < Teuchos::as<size_t>(myAggSize); i++) {
826  localQR(i,j) = (*qFactor)(i,j);
827  }
828  }
829  }
830 
831  // end default case (myAggSize >= NSDim)
832  } else { // special handling for myAggSize < NSDim (i.e. 1pt nodes)
833  // See comments for the uncoupled case
834 
835  // R = extended (by adding identity rows) localQR
836  for (size_t j = 0; j < NSDim; j++)
837  for (size_t k = 0; k < NSDim; k++) {
838  TEUCHOS_TEST_FOR_EXCEPTION(!coarseMapRef.isNodeLocalElement(offset+k), Exceptions::RuntimeError,
839  "Caught error in coarseNS insert, j=" << j << ", offset+k = " << offset+k);
840 
841  if (k < as<size_t>(myAggSize))
842  coarseNS[j][offset+k] = localQR(k,j);
843  else
844  coarseNS[j][offset+k] = (k == j ? one : zero);
845  }
846 
847  // Q = I (rectangular)
848  for (size_t i = 0; i < as<size_t>(myAggSize); i++)
849  for (size_t j = 0; j < NSDim; j++)
850  localQR(i,j) = (j == i ? one : zero);
851  } // end else (special handling for 1pt aggregates)
852 
853  //Process each row in the local Q factor. If the row is local to the current processor
854  //according to the rowmap, insert it into Ptentative. Otherwise, save it in ghostQ
855  //to be communicated later to the owning processor.
856  //FIXME -- what happens if maps are blocked?
857  for (GO j=0; j<myAggSize; ++j) {
858  //This loop checks whether row associated with current DOF is local, according to rowMapForPtent.
859  //If it is, the row is inserted. If not, the row number, columns, and values are saved in
860  //MultiVectors that will be sent to other processors.
861  GO globalRow = aggToRowMap[aggStart[agg]+j];
862 
863  //TODO is the use of Xpetra::global_size_t below correct?
864  if (rowMapForPtentRef.isNodeGlobalElement(globalRow) == false ) {
865  ghostQrows[qctr] = globalRow;
866  for (size_t k=0; k<NSDim; ++k) {
867  ghostQcols[k][qctr] = coarseMapRef.getGlobalElement(agg*NSDim+k);
868  ghostQvals[k][qctr] = localQR(j,k);
869  }
870  ++qctr;
871  } else {
872  size_t nnz=0;
873  for (size_t k=0; k<NSDim; ++k) {
874  try{
875  if (localQR(j,k) != Teuchos::ScalarTraits<SC>::zero()) {
876  localColPtr[nnz] = agg * NSDim + k;
877  globalColPtr[nnz] = coarseMapRef.getGlobalElement(localColPtr[nnz]);
878  valPtr[nnz] = localQR(j,k);
879  ++total_nnz_count;
880  ++nnz;
881  }
882  }
883  catch(...) {
884  GetOStream(Errors,-1) << "caught error in colPtr/valPtr insert, current index="<<nnz<<std::endl;
885  }
886  } //for (size_t k=0; k<NSDim; ++k)
887 
888  try{
889  Ptentative->insertGlobalValues(globalRow,globalColPtr.view(0,nnz),valPtr.view(0,nnz));
890  }
891  catch(...) {
892  GetOStream(Errors,-1) << "pid " << A->getRowMap()->getComm()->getRank()
893  << "caught error during Ptent row insertion, global row "
894  << globalRow << std::endl;
895  }
896  }
897  } //for (GO j=0; j<myAggSize; ++j)
898 
899  } // for (LO agg=0; agg<numAggs; ++agg)
900 
901 
902  // ***********************************************************
903  // ************* end of aggregate-wise QR ********************
904  // ***********************************************************
905  GetOStream(Runtime1) << "TentativePFactory : aggregates may cross process boundaries" << std::endl;
906  // Import ghost parts of Q factors and insert into Ptentative.
907  // First import just the global row numbers.
908  RCP<Xpetra::Vector<GO,LO,GO,Node> > targetQrowNums = Xpetra::VectorFactory<GO,LO,GO,Node>::Build(rowMapForPtent);
909  targetQrowNums->putScalar(-1);
910  targetQrowNums->doImport(*ghostQrowNums,*importer,Xpetra::INSERT);
911  ArrayRCP< GO > targetQrows = targetQrowNums->getDataNonConst(0);
912 
913  // Now create map based on just the row numbers imported.
914  Array<GO> gidsToImport;
915  gidsToImport.reserve(targetQrows.size());
916  for (typename ArrayRCP<GO>::iterator r=targetQrows.begin(); r!=targetQrows.end(); ++r) {
917  if (*r > -1) {
918  gidsToImport.push_back(*r);
919  }
920  }
921  RCP<const Map > reducedMap = MapFactory::Build( A->getRowMap()->lib(),
923  gidsToImport, indexBase, A->getRowMap()->getComm() );
924 
925  // Import using the row numbers that this processor will receive.
926  importer = ImportFactory::Build(ghostQMap, reducedMap);
927 
928  Array<RCP<Xpetra::Vector<GO,LO,GO,Node> > > targetQcolumns(NSDim);
929  for (size_t i=0; i<NSDim; ++i) {
930  targetQcolumns[i] = Xpetra::VectorFactory<GO,LO,GO,Node>::Build(reducedMap);
931  targetQcolumns[i]->doImport(*(ghostQcolumns[i]),*importer,Xpetra::INSERT);
932  }
933  RCP<MultiVector> targetQvalues = MultiVectorFactory::Build(reducedMap,NSDim);
934  targetQvalues->doImport(*ghostQvalues,*importer,Xpetra::INSERT);
935 
936  ArrayRCP< ArrayRCP<SC> > targetQvals;
937  ArrayRCP<ArrayRCP<GO> > targetQcols;
938  if (targetQvalues->getLocalLength() > 0) {
939  targetQvals.resize(NSDim);
940  targetQcols.resize(NSDim);
941  for (size_t i=0; i<NSDim; ++i) {
942  targetQvals[i] = targetQvalues->getDataNonConst(i);
943  targetQcols[i] = targetQcolumns[i]->getDataNonConst(0);
944  }
945  }
946 
947  valPtr = Array<SC>(NSDim,0.);
948  globalColPtr = Array<GO>(NSDim,0);
949  for (typename Array<GO>::iterator r=gidsToImport.begin(); r!=gidsToImport.end(); ++r) {
950  if (targetQvalues->getLocalLength() > 0) {
951  for (size_t j=0; j<NSDim; ++j) {
952  valPtr[j] = targetQvals[j][reducedMap->getLocalElement(*r)];
953  globalColPtr[j] = targetQcols[j][reducedMap->getLocalElement(*r)];
954  }
955  Ptentative->insertGlobalValues(*r, globalColPtr.view(0,NSDim), valPtr.view(0,NSDim));
956  } //if (targetQvalues->getLocalLength() > 0)
957  }
958 
959  Ptentative->fillComplete(coarseMap, A->getDomainMap());
960  }
961 
962  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
963  bool TentativePFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::isGoodMap(const Map& rowMap, const Map& colMap) const {
964  ArrayView<const GO> rowElements = rowMap.getNodeElementList();
965  ArrayView<const GO> colElements = colMap.getNodeElementList();
966 
967  const size_t numElements = rowElements.size();
968 
969  if (size_t(colElements.size()) < numElements)
970  return false;
971 
972  bool goodMap = true;
973  for (size_t i = 0; i < numElements; i++)
974  if (rowElements[i] != colElements[i]) {
975  goodMap = false;
976  break;
977  }
978 
979  return goodMap;
980  }
981 
982 } //namespace MueLu
983 
984 // TODO ReUse: If only P or Nullspace is missing, TentativePFactory can be smart and skip part of the computation.
985 
986 #define MUELU_TENTATIVEPFACTORY_SHORT
987 #endif // MUELU_TENTATIVEPFACTORY_DEF_HPP
void BuildPuncoupled(RCP< Matrix > A, RCP< Aggregates > aggregates, RCP< AmalgamationInfo > amalgInfo, RCP< MultiVector > fineNullspace, RCP< const Map > coarseMap, RCP< Matrix > &Ptentative, RCP< MultiVector > &coarseNullspace, const int levelID) const
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
Important warning messages (one line)
void reserve(size_type n)
MueLu::DefaultLocalOrdinal LocalOrdinal
int setMatrix(const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &A)
void BuildPcoupled(RCP< Matrix > A, RCP< Aggregates > aggregates, RCP< AmalgamationInfo > amalgInfo, RCP< MultiVector > fineNullspace, RCP< const Map > coarseMap, RCP< Matrix > &Ptentative, RCP< MultiVector > &coarseNullspace) const
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
T & get(const std::string &name, T def_value)
void UnamalgamateAggregates(const Aggregates &aggregates, Teuchos::ArrayRCP< LocalOrdinal > &aggStart, Teuchos::ArrayRCP< GlobalOrdinal > &aggToRowMap) const
UnamalgamateAggregates.
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
ArrayView< T > view(size_type offset, size_type size)
static T squareroot(T x)
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
size_type size() const
void UnamalgamateAggregatesLO(const Aggregates &aggregates, Teuchos::ArrayRCP< LocalOrdinal > &aggStart, Teuchos::ArrayRCP< LO > &aggToRowMap) const
size_type size() const
Teuchos::ArrayRCP< LO > ComputeAggregateSizes(bool forceRecompute=false) const
Compute sizes of aggregates.
static const NoFactory * get()
Print even more statistics.
#define SET_VALID_ENTRY(name)
void resize(const size_type n, const T &val=T())
bool isParameter(const std::string &name) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
Teuchos::RCP< Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > ComputeUnamalgamatedImportDofMap(const Aggregates &aggregates) const
ComputeUnamalgamatedImportDofMap build overlapping dof row map from aggregates needed for overlapping...
bool isSublist(const std::string &name) const
LO GetNumAggregates() const
returns the number of aggregates of the current processor. Note: could/should be renamed to GetNumLoc...
void resize(size_type new_size, const value_type &x=value_type())
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
static bool isnaninf(const T &x)
bool isGoodMap(const Map &rowMap, const Map &colMap) const
iterator end()
static magnitudeType magnitude(T a)
void push_back(const value_type &x)
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
int GetLevelID() const
Return level number.
Definition: MueLu_Level.cpp:76
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getQ() const
Exception throws to report errors in the internal logical of the program.
iterator end() const
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
Description of what is happening (more verbose)
iterator begin()
iterator begin() const
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need&#39;s value has been saved.