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