MueLu Version of the Day
Loading...
Searching...
No Matches
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_CrsGraphFactory.hpp>
53#include <Xpetra_Matrix.hpp>
54#include <Xpetra_MatrixMatrix.hpp>
55#include <Xpetra_MultiVector.hpp>
56#include <Xpetra_MultiVectorFactory.hpp>
57#include <Xpetra_VectorFactory.hpp>
58#include <Xpetra_Import.hpp>
59#include <Xpetra_ImportFactory.hpp>
60#include <Xpetra_CrsMatrixWrap.hpp>
61#include <Xpetra_StridedMap.hpp>
62#include <Xpetra_StridedMapFactory.hpp>
63#include <Xpetra_IO.hpp>
64
65#include "Xpetra_TpetraBlockCrsMatrix.hpp"
66
68
69#include "MueLu_Aggregates.hpp"
70#include "MueLu_AmalgamationInfo.hpp"
71#include "MueLu_MasterList.hpp"
72#include "MueLu_Monitor.hpp"
73#include "MueLu_PerfUtils.hpp"
74#include "MueLu_Utilities.hpp"
75
76
77
78
79namespace MueLu {
80
81 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
83 RCP<ParameterList> validParamList = rcp(new ParameterList());
84
85#define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
86 SET_VALID_ENTRY("tentative: calculate qr");
87 SET_VALID_ENTRY("tentative: build coarse coordinates");
88 SET_VALID_ENTRY("tentative: constant column sums");
89#undef SET_VALID_ENTRY
90 validParamList->set< std::string >("Nullspace name", "Nullspace", "Name for the input nullspace");
91
92 validParamList->set< RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A");
93 validParamList->set< RCP<const FactoryBase> >("Aggregates", Teuchos::null, "Generating factory of the aggregates");
94 validParamList->set< RCP<const FactoryBase> >("Nullspace", Teuchos::null, "Generating factory of the nullspace");
95 validParamList->set< RCP<const FactoryBase> >("Scaled Nullspace", Teuchos::null, "Generating factory of the scaled nullspace");
96 validParamList->set< RCP<const FactoryBase> >("UnAmalgamationInfo", Teuchos::null, "Generating factory of UnAmalgamationInfo");
97 validParamList->set< RCP<const FactoryBase> >("CoarseMap", Teuchos::null, "Generating factory of the coarse map");
98 validParamList->set< RCP<const FactoryBase> >("Coordinates", Teuchos::null, "Generating factory of the coordinates");
99 validParamList->set< RCP<const FactoryBase> >("Node Comm", Teuchos::null, "Generating factory of the node level communicator");
100
101 // Make sure we don't recursively validate options for the matrixmatrix kernels
102 ParameterList norecurse;
103 norecurse.disableRecursiveValidation();
104 validParamList->set<ParameterList> ("matrixmatrix: kernel params", norecurse, "MatrixMatrix kernel parameters");
105
106 return validParamList;
107 }
108
109 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
111
112 const ParameterList& pL = GetParameterList();
113 // NOTE: This guy can only either be 'Nullspace' or 'Scaled Nullspace' or else the validator above will cause issues
114 std::string nspName = "Nullspace";
115 if(pL.isParameter("Nullspace name")) nspName = pL.get<std::string>("Nullspace name");
116
117 Input(fineLevel, "A");
118 Input(fineLevel, "Aggregates");
119 Input(fineLevel, nspName);
120 Input(fineLevel, "UnAmalgamationInfo");
121 Input(fineLevel, "CoarseMap");
122 if( fineLevel.GetLevelID() == 0 &&
123 fineLevel.IsAvailable("Coordinates", NoFactory::get()) && // we have coordinates (provided by user app)
124 pL.get<bool>("tentative: build coarse coordinates") ) { // and we want coordinates on other levels
125 bTransferCoordinates_ = true; // then set the transfer coordinates flag to true
126 Input(fineLevel, "Coordinates");
127 } else if (bTransferCoordinates_) {
128 Input(fineLevel, "Coordinates");
129 }
130 }
131
132 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
134 return BuildP(fineLevel, coarseLevel);
135 }
136
137 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
139 FactoryMonitor m(*this, "Build", coarseLevel);
140 typedef typename Teuchos::ScalarTraits<Scalar>::coordinateType coordinate_type;
141 typedef Xpetra::MultiVector<coordinate_type,LO,GO,NO> RealValuedMultiVector;
142 typedef Xpetra::MultiVectorFactory<coordinate_type,LO,GO,NO> RealValuedMultiVectorFactory;
143
144 const ParameterList& pL = GetParameterList();
145 std::string nspName = "Nullspace";
146 if(pL.isParameter("Nullspace name")) nspName = pL.get<std::string>("Nullspace name");
147
148
149 RCP<Matrix> Ptentative;
150 RCP<Matrix> A = Get< RCP<Matrix> > (fineLevel, "A");
151 RCP<Aggregates> aggregates = Get< RCP<Aggregates> > (fineLevel, "Aggregates");
152 // No coarse DoFs so we need to bail by setting Ptentattive to null and returning
153 // This level will ultimately be removed in MueLu_Hierarchy_defs.h via a resize()
154 if ( aggregates->GetNumGlobalAggregatesComputeIfNeeded() == 0) {
155 Ptentative = Teuchos::null;
156 Set(coarseLevel, "P", Ptentative);
157 return;
158 }
159
160 RCP<AmalgamationInfo> amalgInfo = Get< RCP<AmalgamationInfo> > (fineLevel, "UnAmalgamationInfo");
161 RCP<MultiVector> fineNullspace = Get< RCP<MultiVector> > (fineLevel, nspName);
162 RCP<const Map> coarseMap = Get< RCP<const Map> > (fineLevel, "CoarseMap");
163 RCP<RealValuedMultiVector> fineCoords;
165 fineCoords = Get< RCP<RealValuedMultiVector> >(fineLevel, "Coordinates");
166 }
167
168 // FIXME: We should remove the NodeComm on levels past the threshold
169 if(fineLevel.IsAvailable("Node Comm")) {
170 RCP<const Teuchos::Comm<int> > nodeComm = Get<RCP<const Teuchos::Comm<int> > >(fineLevel,"Node Comm");
171 Set<RCP<const Teuchos::Comm<int> > >(coarseLevel, "Node Comm", nodeComm);
172 }
173
174 // NOTE: We check DomainMap here rather than RowMap because those are different for BlockCrs matrices
175 TEUCHOS_TEST_FOR_EXCEPTION( A->getDomainMap()->getLocalNumElements() != fineNullspace->getMap()->getLocalNumElements(),
176 Exceptions::RuntimeError,"MueLu::TentativePFactory::MakeTentative: Size mismatch between A and Nullspace");
177
178 RCP<MultiVector> coarseNullspace;
179 RCP<RealValuedMultiVector> coarseCoords;
180
182 //*** Create the coarse coordinates ***
183 // First create the coarse map and coarse multivector
184 ArrayView<const GO> elementAList = coarseMap->getLocalElementList();
185 LO blkSize = 1;
186 if (rcp_dynamic_cast<const StridedMap>(coarseMap) != Teuchos::null) {
187 blkSize = rcp_dynamic_cast<const StridedMap>(coarseMap)->getFixedBlockSize();
188 }
189 GO indexBase = coarseMap->getIndexBase();
190 LO numCoarseNodes = Teuchos::as<LO>(elementAList.size() / blkSize);
191 Array<GO> nodeList(numCoarseNodes);
192 const int numDimensions = fineCoords->getNumVectors();
193
194 for (LO i = 0; i < numCoarseNodes; i++) {
195 nodeList[i] = (elementAList[i*blkSize]-indexBase)/blkSize + indexBase;
196 }
197 RCP<const Map> coarseCoordsMap = MapFactory::Build(fineCoords->getMap()->lib(),
198 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
199 nodeList,
200 indexBase,
201 fineCoords->getMap()->getComm());
202 coarseCoords = RealValuedMultiVectorFactory::Build(coarseCoordsMap, numDimensions);
203
204 // Create overlapped fine coordinates to reduce global communication
205 RCP<RealValuedMultiVector> ghostedCoords;
206 if (aggregates->AggregatesCrossProcessors()) {
207 RCP<const Map> aggMap = aggregates->GetMap();
208 RCP<const Import> importer = ImportFactory::Build(fineCoords->getMap(), aggMap);
209
210 ghostedCoords = RealValuedMultiVectorFactory::Build(aggMap, numDimensions);
211 ghostedCoords->doImport(*fineCoords, *importer, Xpetra::INSERT);
212 } else {
213 ghostedCoords = fineCoords;
214 }
215
216 // Get some info about aggregates
217 int myPID = coarseCoordsMap->getComm()->getRank();
218 LO numAggs = aggregates->GetNumAggregates();
219 ArrayRCP<LO> aggSizes = aggregates->ComputeAggregateSizesArrayRCP();
220 const ArrayRCP<const LO> vertex2AggID = aggregates->GetVertex2AggId()->getData(0);
221 const ArrayRCP<const LO> procWinner = aggregates->GetProcWinner()->getData(0);
222
223 // Fill in coarse coordinates
224 for (int dim = 0; dim < numDimensions; ++dim) {
225 ArrayRCP<const coordinate_type> fineCoordsData = ghostedCoords->getData(dim);
226 ArrayRCP<coordinate_type> coarseCoordsData = coarseCoords->getDataNonConst(dim);
227
228 for (LO lnode = 0; lnode < Teuchos::as<LO>(vertex2AggID.size()); lnode++) {
229 if (procWinner[lnode] == myPID &&
230 lnode < fineCoordsData.size() &&
231 vertex2AggID[lnode] < coarseCoordsData.size() &&
232 Teuchos::ScalarTraits<coordinate_type>::isnaninf(fineCoordsData[lnode]) == false) {
233 coarseCoordsData[vertex2AggID[lnode]] += fineCoordsData[lnode];
234 }
235 }
236 for (LO agg = 0; agg < numAggs; agg++) {
237 coarseCoordsData[agg] /= aggSizes[agg];
238 }
239 }
240 }
241
242 if (!aggregates->AggregatesCrossProcessors()) {
243 if(Xpetra::Helpers<SC,LO,GO,NO>::isTpetraBlockCrs(A)) {
244 BuildPuncoupledBlockCrs(A, aggregates, amalgInfo, fineNullspace, coarseMap, Ptentative, coarseNullspace,coarseLevel.GetLevelID());
245 }
246 else {
247 BuildPuncoupled(A, aggregates, amalgInfo, fineNullspace, coarseMap, Ptentative, coarseNullspace,coarseLevel.GetLevelID());
248 }
249 }
250 else
251 BuildPcoupled(A, aggregates, amalgInfo, fineNullspace, coarseMap, Ptentative, coarseNullspace);
252
253 // If available, use striding information of fine level matrix A for range
254 // map and coarseMap as domain map; otherwise use plain range map of
255 // Ptent = plain range map of A for range map and coarseMap as domain map.
256 // NOTE:
257 // The latter is not really safe, since there is no striding information
258 // for the range map. This is not really a problem, since striding
259 // information is always available on the intermedium levels and the
260 // coarsest levels.
261 if (A->IsView("stridedMaps") == true)
262 Ptentative->CreateView("stridedMaps", A->getRowMap("stridedMaps"), coarseMap);
263 else
264 Ptentative->CreateView("stridedMaps", Ptentative->getRangeMap(), coarseMap);
265
267 Set(coarseLevel, "Coordinates", coarseCoords);
268 }
269 Set(coarseLevel, "Nullspace", coarseNullspace);
270 Set(coarseLevel, "P", Ptentative);
271
272 if (IsPrint(Statistics2)) {
273 RCP<ParameterList> params = rcp(new ParameterList());
274 params->set("printLoadBalancingInfo", true);
275 GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*Ptentative, "Ptent", params);
276 }
277 }
278
279 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
281 BuildPuncoupledBlockCrs(RCP<Matrix> A, RCP<Aggregates> aggregates, RCP<AmalgamationInfo> amalgInfo, RCP<MultiVector> fineNullspace,
282 RCP<const Map> coarsePointMap, RCP<Matrix>& Ptentative, RCP<MultiVector>& coarseNullspace, const int levelID) const {
283 /* This routine generates a BlockCrs P for a BlockCrs A. There are a few assumptions here, which meet the use cases we care about, but could
284 be generalized later, if we ever need to do so:
285 1) Null space dimension === block size of matrix: So no elasticity right now
286 2) QR is not supported: Under assumption #1, this shouldn't cause problems.
287 3) Maps are "good": Aka the first chunk of the ColMap is the RowMap.
288
289 These assumptions keep our code way simpler and still support the use cases we actually care about.
290 */
291
292 RCP<const Map> rowMap = A->getRowMap();
293 RCP<const Map> rangeMap = A->getRangeMap();
294 RCP<const Map> colMap = A->getColMap();
295 // const size_t numFinePointRows = rangeMap->getLocalNumElements();
296 const size_t numFineBlockRows = rowMap->getLocalNumElements();
297
298 typedef Teuchos::ScalarTraits<SC> STS;
299 // typedef typename STS::magnitudeType Magnitude;
300 const SC zero = STS::zero();
301 const SC one = STS::one();
302 const LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
303
304 const GO numAggs = aggregates->GetNumAggregates();
305 const size_t NSDim = fineNullspace->getNumVectors();
306 ArrayRCP<LO> aggSizes = aggregates->ComputeAggregateSizesArrayRCP();
307
308 // Need to generate the coarse block map
309 // NOTE: We assume NSDim == block size here
310 // NOTE: We also assume that coarseMap has contiguous GIDs
311 //const size_t numCoarsePointRows = coarsePointMap->getLocalNumElements();
312 const size_t numCoarseBlockRows = coarsePointMap->getLocalNumElements() / NSDim;
313 RCP<const Map> coarseBlockMap = MapFactory::Build(coarsePointMap->lib(),
314 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
315 numCoarseBlockRows,
316 coarsePointMap->getIndexBase(),
317 coarsePointMap->getComm());
318 // Sanity checking
319 const ParameterList& pL = GetParameterList();
320 const bool &doQRStep = pL.get<bool>("tentative: calculate qr");
321 const bool &constantColSums = pL.get<bool>("tentative: constant column sums");
322
323 TEUCHOS_TEST_FOR_EXCEPTION(doQRStep && constantColSums,Exceptions::RuntimeError,
324 "MueLu::TentativePFactory::MakeTentative: cannot use 'constant column sums' and 'calculate qr' at the same time");
325
326 // The aggregates use the amalgamated column map, which in this case is what we want
327
328 // Aggregates map is based on the amalgamated column map
329 // We can skip global-to-local conversion if LIDs in row map are
330 // same as LIDs in column map
331 bool goodMap = MueLu::Utilities<SC,LO,GO,NO>::MapsAreNested(*rowMap, *colMap);
332
333 // Create a lookup table to determine the rows (fine DOFs) that belong to a given aggregate.
334 // aggStart is a pointer into aggToRowMapLO
335 // aggStart[i]..aggStart[i+1] are indices into aggToRowMapLO
336 // aggToRowMapLO[aggStart[i]]..aggToRowMapLO[aggStart[i+1]-1] are the DOFs in aggregate i
337 ArrayRCP<LO> aggStart;
338 ArrayRCP<LO> aggToRowMapLO;
339 ArrayRCP<GO> aggToRowMapGO;
340 if (goodMap) {
341 amalgInfo->UnamalgamateAggregatesLO(*aggregates, aggStart, aggToRowMapLO);
342 GetOStream(Runtime1) << "Column map is consistent with the row map, good." << std::endl;
343 } else {
344 throw std::runtime_error("TentativePFactory::PuncoupledBlockCrs: Inconsistent maps not currently supported");
345 }
346
347 coarseNullspace = MultiVectorFactory::Build(coarsePointMap, NSDim);
348
349 // Pull out the nullspace vectors so that we can have random access.
350 ArrayRCP<ArrayRCP<const SC> > fineNS (NSDim);
351 ArrayRCP<ArrayRCP<SC> > coarseNS(NSDim);
352 for (size_t i = 0; i < NSDim; i++) {
353 fineNS[i] = fineNullspace->getData(i);
354 if (coarsePointMap->getLocalNumElements() > 0)
355 coarseNS[i] = coarseNullspace->getDataNonConst(i);
356 }
357
358 // BlockCrs requires that we build the (block) graph first, so let's do that...
359 // NOTE: Because we're assuming that the NSDim == BlockSize, we only have one
360 // block non-zero per row in the matrix;
361 RCP<CrsGraph> BlockGraph = CrsGraphFactory::Build(rowMap,coarseBlockMap,0);
362 ArrayRCP<size_t> iaPtent;
363 ArrayRCP<LO> jaPtent;
364 BlockGraph->allocateAllIndices(numFineBlockRows, iaPtent, jaPtent);
365 ArrayView<size_t> ia = iaPtent();
366 ArrayView<LO> ja = jaPtent();
367
368 for (size_t i = 0; i < numFineBlockRows; i++) {
369 ia[i] = i;
370 ja[i] = INVALID;
371 }
372 ia[numCoarseBlockRows] = numCoarseBlockRows;
373
374
375 for (GO agg = 0; agg < numAggs; agg++) {
376 LO aggSize = aggStart[agg+1] - aggStart[agg];
377 Xpetra::global_size_t offset = agg;
378
379 for (LO j = 0; j < aggSize; j++) {
380 // FIXME: Allow for bad maps
381 const LO localRow = aggToRowMapLO[aggStart[agg]+j];
382 const size_t rowStart = ia[localRow];
383 ja[rowStart] = offset;
384 }
385 }
386
387 // Compress storage (remove all INVALID, which happen when we skip zeros)
388 // We do that in-place
389 size_t ia_tmp = 0, nnz = 0;
390 for (size_t i = 0; i < numFineBlockRows; i++) {
391 for (size_t j = ia_tmp; j < ia[i+1]; j++)
392 if (ja[j] != INVALID) {
393 ja [nnz] = ja [j];
394 nnz++;
395 }
396 ia_tmp = ia[i+1];
397 ia[i+1] = nnz;
398 }
399
400 if (rowMap->lib() == Xpetra::UseTpetra) {
401 // - Cannot resize for Epetra, as it checks for same pointers
402 // - Need to resize for Tpetra, as it check ().size() == ia[numRows]
403 // NOTE: these invalidate ja and val views
404 jaPtent .resize(nnz);
405 }
406
407 GetOStream(Runtime1) << "TentativePFactory : generating block graph" << std::endl;
408 BlockGraph->setAllIndices(iaPtent, jaPtent);
409
410 // Managing labels & constants for ESFC
411 {
412 RCP<ParameterList> FCparams;
413 if(pL.isSublist("matrixmatrix: kernel params"))
414 FCparams=rcp(new ParameterList(pL.sublist("matrixmatrix: kernel params")));
415 else
416 FCparams= rcp(new ParameterList);
417 // By default, we don't need global constants for TentativeP, but we do want it for the graph
418 // if we're printing statistics, so let's leave it on for now.
419 FCparams->set("compute global constants",FCparams->get("compute global constants",true));
420 std::string levelIDs = toString(levelID);
421 FCparams->set("Timer Label",std::string("MueLu::TentativeP-")+levelIDs);
422 RCP<const Export> dummy_e;
423 RCP<const Import> dummy_i;
424 BlockGraph->expertStaticFillComplete(coarseBlockMap,rowMap,dummy_i,dummy_e,FCparams);
425 }
426
427 // Now let's make a BlockCrs Matrix
428 // NOTE: Assumes block size== NSDim
429 RCP<Xpetra::CrsMatrix<SC,LO,GO,NO> > P_xpetra = Xpetra::CrsMatrixFactory<SC,LO,GO,NO>::BuildBlock(BlockGraph, coarsePointMap, rangeMap,NSDim);
430 RCP<Xpetra::TpetraBlockCrsMatrix<SC,LO,GO,NO> > P_tpetra = rcp_dynamic_cast<Xpetra::TpetraBlockCrsMatrix<SC,LO,GO,NO> >(P_xpetra);
431 if(P_tpetra.is_null()) throw std::runtime_error("BuildPUncoupled: Matrix factory did not return a Tpetra::BlockCrsMatrix");
432 RCP<CrsMatrixWrap> P_wrap = rcp(new CrsMatrixWrap(P_xpetra));
433
435 // "no-QR" option //
437 // Local Q factor is just the fine nullspace support over the current aggregate.
438 // Local R factor is the identity.
439 // NOTE: We're not going to do a QR here as we're assuming that blocksize == NSDim
440 // NOTE: "goodMap" case only
441 Teuchos::Array<Scalar> block(NSDim*NSDim, zero);
442 Teuchos::Array<LO> bcol(1);
443
444 GetOStream(Runtime1) << "TentativePFactory : bypassing local QR phase" << std::endl;
445 for (LO agg = 0; agg < numAggs; agg++) {
446 bcol[0] = agg;
447 const LO aggSize = aggStart[agg+1] - aggStart[agg];
448 Xpetra::global_size_t offset = agg*NSDim;
449
450 // Process each row in the local Q factor
451 // NOTE: Blocks are in row-major order
452 for (LO j = 0; j < aggSize; j++) {
453 const LO localBlockRow = aggToRowMapLO[aggStart[agg]+j];
454
455 for (size_t r = 0; r < NSDim; r++) {
456 LO localPointRow = localBlockRow*NSDim + r;
457 for (size_t c = 0; c < NSDim; c++)
458 block[r*NSDim+c] = fineNS[c][localPointRow];
459 }
460 // NOTE: Assumes columns==aggs and are ordered sequentially
461 P_tpetra->replaceLocalValues(localBlockRow,bcol(),block());
462
463 }//end aggSize
464
465 for (size_t j = 0; j < NSDim; j++)
466 coarseNS[j][offset+j] = one;
467
468 } //for (GO agg = 0; agg < numAggs; agg++)
469
470 Ptentative = P_wrap;
471 }
472
473 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
475 BuildPcoupled(RCP<Matrix> A, RCP<Aggregates> aggregates, RCP<AmalgamationInfo> amalgInfo, RCP<MultiVector> fineNullspace,
476 RCP<const Map> coarseMap, RCP<Matrix>& Ptentative, RCP<MultiVector>& coarseNullspace) const {
477 typedef Teuchos::ScalarTraits<SC> STS;
478 typedef typename STS::magnitudeType Magnitude;
479 const SC zero = STS::zero();
480 const SC one = STS::one();
481
482 // number of aggregates
483 GO numAggs = aggregates->GetNumAggregates();
484
485 // Create a lookup table to determine the rows (fine DOFs) that belong to a given aggregate.
486 // aggStart is a pointer into aggToRowMap
487 // aggStart[i]..aggStart[i+1] are indices into aggToRowMap
488 // aggToRowMap[aggStart[i]]..aggToRowMap[aggStart[i+1]-1] are the DOFs in aggregate i
489 ArrayRCP<LO> aggStart;
490 ArrayRCP< GO > aggToRowMap;
491 amalgInfo->UnamalgamateAggregates(*aggregates, aggStart, aggToRowMap);
492
493 // find size of the largest aggregate
494 LO maxAggSize=0;
495 for (GO i=0; i<numAggs; ++i) {
496 LO sizeOfThisAgg = aggStart[i+1] - aggStart[i];
497 if (sizeOfThisAgg > maxAggSize) maxAggSize = sizeOfThisAgg;
498 }
499
500 // dimension of fine level nullspace
501 const size_t NSDim = fineNullspace->getNumVectors();
502
503 // index base for coarse Dof map (usually 0)
504 GO indexBase=A->getRowMap()->getIndexBase();
505
506 const RCP<const Map> nonUniqueMap = amalgInfo->ComputeUnamalgamatedImportDofMap(*aggregates);
507 const RCP<const Map> uniqueMap = A->getDomainMap();
508 RCP<const Import> importer = ImportFactory::Build(uniqueMap, nonUniqueMap);
509 RCP<MultiVector> fineNullspaceWithOverlap = MultiVectorFactory::Build(nonUniqueMap,NSDim);
510 fineNullspaceWithOverlap->doImport(*fineNullspace,*importer,Xpetra::INSERT);
511
512 // Pull out the nullspace vectors so that we can have random access.
513 ArrayRCP< ArrayRCP<const SC> > fineNS(NSDim);
514 for (size_t i=0; i<NSDim; ++i)
515 fineNS[i] = fineNullspaceWithOverlap->getData(i);
516
517 //Allocate storage for the coarse nullspace.
518 coarseNullspace = MultiVectorFactory::Build(coarseMap, NSDim);
519
520 ArrayRCP< ArrayRCP<SC> > coarseNS(NSDim);
521 for (size_t i=0; i<NSDim; ++i)
522 if (coarseMap->getLocalNumElements() > 0) coarseNS[i] = coarseNullspace->getDataNonConst(i);
523
524 //This makes the rowmap of Ptent the same as that of A->
525 //This requires moving some parts of some local Q's to other processors
526 //because aggregates can span processors.
527 RCP<const Map > rowMapForPtent = A->getRowMap();
528 const Map& rowMapForPtentRef = *rowMapForPtent;
529
530 // Set up storage for the rows of the local Qs that belong to other processors.
531 // FIXME This is inefficient and could be done within the main loop below with std::vector's.
532 RCP<const Map> colMap = A->getColMap();
533
534 RCP<const Map > ghostQMap;
535 RCP<MultiVector> ghostQvalues;
536 Array<RCP<Xpetra::Vector<GO,LO,GO,Node> > > ghostQcolumns;
537 RCP<Xpetra::Vector<GO,LO,GO,Node> > ghostQrowNums;
538 ArrayRCP< ArrayRCP<SC> > ghostQvals;
539 ArrayRCP< ArrayRCP<GO> > ghostQcols;
540 ArrayRCP< GO > ghostQrows;
541
542 Array<GO> ghostGIDs;
543 for (LO j=0; j<numAggs; ++j) {
544 for (LO k=aggStart[j]; k<aggStart[j+1]; ++k) {
545 if (rowMapForPtentRef.isNodeGlobalElement(aggToRowMap[k]) == false) {
546 ghostGIDs.push_back(aggToRowMap[k]);
547 }
548 }
549 }
550 ghostQMap = MapFactory::Build(A->getRowMap()->lib(),
551 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
552 ghostGIDs,
553 indexBase, A->getRowMap()->getComm()); //JG:Xpetra::global_size_t>?
554 //Vector to hold bits of Q that go to other processors.
555 ghostQvalues = MultiVectorFactory::Build(ghostQMap,NSDim);
556 //Note that Epetra does not support MultiVectors templated on Scalar != double.
557 //So to work around this, we allocate an array of Vectors. This shouldn't be too
558 //expensive, as the number of Vectors is NSDim.
559 ghostQcolumns.resize(NSDim);
560 for (size_t i=0; i<NSDim; ++i)
561 ghostQcolumns[i] = Xpetra::VectorFactory<GO,LO,GO,Node>::Build(ghostQMap);
562 ghostQrowNums = Xpetra::VectorFactory<GO,LO,GO,Node>::Build(ghostQMap);
563 if (ghostQvalues->getLocalLength() > 0) {
564 ghostQvals.resize(NSDim);
565 ghostQcols.resize(NSDim);
566 for (size_t i=0; i<NSDim; ++i) {
567 ghostQvals[i] = ghostQvalues->getDataNonConst(i);
568 ghostQcols[i] = ghostQcolumns[i]->getDataNonConst(0);
569 }
570 ghostQrows = ghostQrowNums->getDataNonConst(0);
571 }
572
573 //importer to handle moving Q
574 importer = ImportFactory::Build(ghostQMap, A->getRowMap());
575
576 // Dense QR solver
577 Teuchos::SerialQRDenseSolver<LO,SC> qrSolver;
578
579 //Allocate temporary storage for the tentative prolongator.
580 Array<GO> globalColPtr(maxAggSize*NSDim,0);
581 Array<LO> localColPtr(maxAggSize*NSDim,0);
582 Array<SC> valPtr(maxAggSize*NSDim,0.);
583
584 //Create column map for Ptent, estimate local #nonzeros in Ptent, and create Ptent itself.
585 const Map& coarseMapRef = *coarseMap;
586
587 // For the 3-arrays constructor
588 ArrayRCP<size_t> ptent_rowptr;
589 ArrayRCP<LO> ptent_colind;
590 ArrayRCP<Scalar> ptent_values;
591
592 // Because ArrayRCPs are slow...
593 ArrayView<size_t> rowptr_v;
594 ArrayView<LO> colind_v;
595 ArrayView<Scalar> values_v;
596
597 // For temporary usage
598 Array<size_t> rowptr_temp;
599 Array<LO> colind_temp;
600 Array<Scalar> values_temp;
601
602 RCP<CrsMatrix> PtentCrs;
603
604 RCP<CrsMatrixWrap> PtentCrsWrap = rcp(new CrsMatrixWrap(rowMapForPtent, NSDim));
605 PtentCrs = PtentCrsWrap->getCrsMatrix();
606 Ptentative = PtentCrsWrap;
607
608 //*****************************************************************
609 //Loop over all aggregates and calculate local QR decompositions.
610 //*****************************************************************
611 GO qctr=0; //for indexing into Ptent data vectors
612 const Map& nonUniqueMapRef = *nonUniqueMap;
613
614 size_t total_nnz_count=0;
615
616 for (GO agg=0; agg<numAggs; ++agg)
617 {
618 LO myAggSize = aggStart[agg+1]-aggStart[agg];
619 // For each aggregate, extract the corresponding piece of the nullspace and put it in the flat array,
620 // "localQR" (in column major format) for the QR routine.
621 Teuchos::SerialDenseMatrix<LO,SC> localQR(myAggSize, NSDim);
622 for (size_t j=0; j<NSDim; ++j) {
623 bool bIsZeroNSColumn = true;
624 for (LO k=0; k<myAggSize; ++k)
625 {
626 // aggToRowMap[aggPtr[i]+k] is the kth DOF in the ith aggregate
627 // fineNS[j][n] is the nth entry in the jth NS vector
628 try{
629 SC nsVal = fineNS[j][ nonUniqueMapRef.getLocalElement(aggToRowMap[aggStart[agg]+k]) ]; // extract information from fine level NS
630 localQR(k,j) = nsVal;
631 if (nsVal != zero) bIsZeroNSColumn = false;
632 }
633 catch(...) {
634 GetOStream(Runtime1,-1) << "length of fine level nsp: " << fineNullspace->getGlobalLength() << std::endl;
635 GetOStream(Runtime1,-1) << "length of fine level nsp w overlap: " << fineNullspaceWithOverlap->getGlobalLength() << std::endl;
636 GetOStream(Runtime1,-1) << "(local?) aggId=" << agg << std::endl;
637 GetOStream(Runtime1,-1) << "aggSize=" << myAggSize << std::endl;
638 GetOStream(Runtime1,-1) << "agg DOF=" << k << std::endl;
639 GetOStream(Runtime1,-1) << "NS vector j=" << j << std::endl;
640 GetOStream(Runtime1,-1) << "j*myAggSize + k = " << j*myAggSize + k << std::endl;
641 GetOStream(Runtime1,-1) << "aggToRowMap["<<agg<<"][" << k << "] = " << aggToRowMap[aggStart[agg]+k] << std::endl;
642 GetOStream(Runtime1,-1) << "id aggToRowMap[agg][k]=" << aggToRowMap[aggStart[agg]+k] << " is global element in nonUniqueMap = " <<
643 nonUniqueMapRef.isNodeGlobalElement(aggToRowMap[aggStart[agg]+k]) << std::endl;
644 GetOStream(Runtime1,-1) << "colMap local id aggToRowMap[agg][k]=" << nonUniqueMapRef.getLocalElement(aggToRowMap[aggStart[agg]+k]) << std::endl;
645 GetOStream(Runtime1,-1) << "fineNS...=" << fineNS[j][ nonUniqueMapRef.getLocalElement(aggToRowMap[aggStart[agg]+k]) ] << std::endl;
646 GetOStream(Errors,-1) << "caught an error!" << std::endl;
647 }
648 } //for (LO k=0 ...
649 TEUCHOS_TEST_FOR_EXCEPTION(bIsZeroNSColumn == true, Exceptions::RuntimeError, "MueLu::TentativePFactory::MakeTentative: fine level NS part has a zero column. Error.");
650 } //for (LO j=0 ...
651
652 Xpetra::global_size_t offset=agg*NSDim;
653
654 if(myAggSize >= Teuchos::as<LocalOrdinal>(NSDim)) {
655 // calculate QR decomposition (standard)
656 // R is stored in localQR (size: myAggSize x NSDim)
657
658 // Householder multiplier
659 SC tau = localQR(0,0);
660
661 if (NSDim == 1) {
662 // Only one nullspace vector, so normalize by hand
663 Magnitude dtemp=0;
664 for (size_t k = 0; k < Teuchos::as<size_t>(myAggSize); ++k) {
665 Magnitude tmag = STS::magnitude(localQR(k,0));
666 dtemp += tmag*tmag;
667 }
668 dtemp = Teuchos::ScalarTraits<Magnitude>::squareroot(dtemp);
669 tau = localQR(0,0);
670 localQR(0,0) = dtemp;
671 } else {
672 qrSolver.setMatrix( Teuchos::rcp(&localQR, false) );
673 qrSolver.factor();
674 }
675
676 // Extract R, the coarse nullspace. This is stored in upper triangular part of localQR.
677 // Note: coarseNS[i][.] is the ith coarse nullspace vector, which may be counter to your intuition.
678 // This stores the (offset+k)th entry only if it is local according to the coarseMap.
679 for (size_t j=0; j<NSDim; ++j) {
680 for (size_t k=0; k<=j; ++k) {
681 try {
682 if (coarseMapRef.isNodeLocalElement(offset+k)) {
683 coarseNS[j][offset+k] = localQR(k, j); //TODO is offset+k the correct local ID?!
684 }
685 }
686 catch(...) {
687 GetOStream(Errors,-1) << "caught error in coarseNS insert, j="<<j<<", offset+k = "<<offset+k<<std::endl;
688 }
689 }
690 }
691
692 // Calculate Q, the tentative prolongator.
693 // The Lapack GEQRF call only works for myAggsize >= NSDim
694
695 if (NSDim == 1) {
696 // Only one nullspace vector, so calculate Q by hand
697 Magnitude dtemp = Teuchos::ScalarTraits<SC>::magnitude(localQR(0,0));
698 localQR(0,0) = tau;
699 dtemp = 1 / dtemp;
700 for (LocalOrdinal i=0; i<myAggSize; ++i) {
701 localQR(i,0) *= dtemp ;
702 }
703 } else {
704 qrSolver.formQ();
705 Teuchos::RCP<Teuchos::SerialDenseMatrix<LO,SC> > qFactor = qrSolver.getQ();
706 for (size_t j=0; j<NSDim; j++) {
707 for (size_t i = 0; i < Teuchos::as<size_t>(myAggSize); i++) {
708 localQR(i,j) = (*qFactor)(i,j);
709 }
710 }
711 }
712
713 // end default case (myAggSize >= NSDim)
714 } else { // special handling for myAggSize < NSDim (i.e. 1pt nodes)
715 // See comments for the uncoupled case
716
717 // R = extended (by adding identity rows) localQR
718 for (size_t j = 0; j < NSDim; j++)
719 for (size_t k = 0; k < NSDim; k++) {
720 TEUCHOS_TEST_FOR_EXCEPTION(!coarseMapRef.isNodeLocalElement(offset+k), Exceptions::RuntimeError,
721 "Caught error in coarseNS insert, j=" << j << ", offset+k = " << offset+k);
722
723 if (k < as<size_t>(myAggSize))
724 coarseNS[j][offset+k] = localQR(k,j);
725 else
726 coarseNS[j][offset+k] = (k == j ? one : zero);
727 }
728
729 // Q = I (rectangular)
730 for (size_t i = 0; i < as<size_t>(myAggSize); i++)
731 for (size_t j = 0; j < NSDim; j++)
732 localQR(i,j) = (j == i ? one : zero);
733 } // end else (special handling for 1pt aggregates)
734
735 //Process each row in the local Q factor. If the row is local to the current processor
736 //according to the rowmap, insert it into Ptentative. Otherwise, save it in ghostQ
737 //to be communicated later to the owning processor.
738 //FIXME -- what happens if maps are blocked?
739 for (GO j=0; j<myAggSize; ++j) {
740 //This loop checks whether row associated with current DOF is local, according to rowMapForPtent.
741 //If it is, the row is inserted. If not, the row number, columns, and values are saved in
742 //MultiVectors that will be sent to other processors.
743 GO globalRow = aggToRowMap[aggStart[agg]+j];
744
745 //TODO is the use of Xpetra::global_size_t below correct?
746 if (rowMapForPtentRef.isNodeGlobalElement(globalRow) == false ) {
747 ghostQrows[qctr] = globalRow;
748 for (size_t k=0; k<NSDim; ++k) {
749 ghostQcols[k][qctr] = coarseMapRef.getGlobalElement(agg*NSDim+k);
750 ghostQvals[k][qctr] = localQR(j,k);
751 }
752 ++qctr;
753 } else {
754 size_t nnz=0;
755 for (size_t k=0; k<NSDim; ++k) {
756 try{
757 if (localQR(j,k) != Teuchos::ScalarTraits<SC>::zero()) {
758 localColPtr[nnz] = agg * NSDim + k;
759 globalColPtr[nnz] = coarseMapRef.getGlobalElement(localColPtr[nnz]);
760 valPtr[nnz] = localQR(j,k);
761 ++total_nnz_count;
762 ++nnz;
763 }
764 }
765 catch(...) {
766 GetOStream(Errors,-1) << "caught error in colPtr/valPtr insert, current index="<<nnz<<std::endl;
767 }
768 } //for (size_t k=0; k<NSDim; ++k)
769
770 try{
771 Ptentative->insertGlobalValues(globalRow,globalColPtr.view(0,nnz),valPtr.view(0,nnz));
772 }
773 catch(...) {
774 GetOStream(Errors,-1) << "pid " << A->getRowMap()->getComm()->getRank()
775 << "caught error during Ptent row insertion, global row "
776 << globalRow << std::endl;
777 }
778 }
779 } //for (GO j=0; j<myAggSize; ++j)
780
781 } // for (LO agg=0; agg<numAggs; ++agg)
782
783
784 // ***********************************************************
785 // ************* end of aggregate-wise QR ********************
786 // ***********************************************************
787 GetOStream(Runtime1) << "TentativePFactory : aggregates may cross process boundaries" << std::endl;
788 // Import ghost parts of Q factors and insert into Ptentative.
789 // First import just the global row numbers.
790 RCP<Xpetra::Vector<GO,LO,GO,Node> > targetQrowNums = Xpetra::VectorFactory<GO,LO,GO,Node>::Build(rowMapForPtent);
791 targetQrowNums->putScalar(-1);
792 targetQrowNums->doImport(*ghostQrowNums,*importer,Xpetra::INSERT);
793 ArrayRCP< GO > targetQrows = targetQrowNums->getDataNonConst(0);
794
795 // Now create map based on just the row numbers imported.
796 Array<GO> gidsToImport;
797 gidsToImport.reserve(targetQrows.size());
798 for (typename ArrayRCP<GO>::iterator r=targetQrows.begin(); r!=targetQrows.end(); ++r) {
799 if (*r > -1) {
800 gidsToImport.push_back(*r);
801 }
802 }
803 RCP<const Map > reducedMap = MapFactory::Build( A->getRowMap()->lib(),
804 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
805 gidsToImport, indexBase, A->getRowMap()->getComm() );
806
807 // Import using the row numbers that this processor will receive.
808 importer = ImportFactory::Build(ghostQMap, reducedMap);
809
810 Array<RCP<Xpetra::Vector<GO,LO,GO,Node> > > targetQcolumns(NSDim);
811 for (size_t i=0; i<NSDim; ++i) {
812 targetQcolumns[i] = Xpetra::VectorFactory<GO,LO,GO,Node>::Build(reducedMap);
813 targetQcolumns[i]->doImport(*(ghostQcolumns[i]),*importer,Xpetra::INSERT);
814 }
815 RCP<MultiVector> targetQvalues = MultiVectorFactory::Build(reducedMap,NSDim);
816 targetQvalues->doImport(*ghostQvalues,*importer,Xpetra::INSERT);
817
818 ArrayRCP< ArrayRCP<SC> > targetQvals;
819 ArrayRCP<ArrayRCP<GO> > targetQcols;
820 if (targetQvalues->getLocalLength() > 0) {
821 targetQvals.resize(NSDim);
822 targetQcols.resize(NSDim);
823 for (size_t i=0; i<NSDim; ++i) {
824 targetQvals[i] = targetQvalues->getDataNonConst(i);
825 targetQcols[i] = targetQcolumns[i]->getDataNonConst(0);
826 }
827 }
828
829 valPtr = Array<SC>(NSDim,0.);
830 globalColPtr = Array<GO>(NSDim,0);
831 for (typename Array<GO>::iterator r=gidsToImport.begin(); r!=gidsToImport.end(); ++r) {
832 if (targetQvalues->getLocalLength() > 0) {
833 for (size_t j=0; j<NSDim; ++j) {
834 valPtr[j] = targetQvals[j][reducedMap->getLocalElement(*r)];
835 globalColPtr[j] = targetQcols[j][reducedMap->getLocalElement(*r)];
836 }
837 Ptentative->insertGlobalValues(*r, globalColPtr.view(0,NSDim), valPtr.view(0,NSDim));
838 } //if (targetQvalues->getLocalLength() > 0)
839 }
840
841 Ptentative->fillComplete(coarseMap, A->getDomainMap());
842 }
843
844
845
846 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
848 BuildPuncoupled(RCP<Matrix> A, RCP<Aggregates> aggregates, RCP<AmalgamationInfo> amalgInfo, RCP<MultiVector> fineNullspace,
849 RCP<const Map> coarseMap, RCP<Matrix>& Ptentative, RCP<MultiVector>& coarseNullspace, const int levelID) const {
850 RCP<const Map> rowMap = A->getRowMap();
851 RCP<const Map> colMap = A->getColMap();
852 const size_t numRows = rowMap->getLocalNumElements();
853
854 typedef Teuchos::ScalarTraits<SC> STS;
855 typedef typename STS::magnitudeType Magnitude;
856 const SC zero = STS::zero();
857 const SC one = STS::one();
858 const LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
859
860 const GO numAggs = aggregates->GetNumAggregates();
861 const size_t NSDim = fineNullspace->getNumVectors();
862 ArrayRCP<LO> aggSizes = aggregates->ComputeAggregateSizesArrayRCP();
863
864
865 // Sanity checking
866 const ParameterList& pL = GetParameterList();
867 const bool &doQRStep = pL.get<bool>("tentative: calculate qr");
868 const bool &constantColSums = pL.get<bool>("tentative: constant column sums");
869
870 TEUCHOS_TEST_FOR_EXCEPTION(doQRStep && constantColSums,Exceptions::RuntimeError,
871 "MueLu::TentativePFactory::MakeTentative: cannot use 'constant column sums' and 'calculate qr' at the same time");
872
873 // Aggregates map is based on the amalgamated column map
874 // We can skip global-to-local conversion if LIDs in row map are
875 // same as LIDs in column map
876 bool goodMap = MueLu::Utilities<SC,LO,GO,NO>::MapsAreNested(*rowMap, *colMap);
877
878 // Create a lookup table to determine the rows (fine DOFs) that belong to a given aggregate.
879 // aggStart is a pointer into aggToRowMapLO
880 // aggStart[i]..aggStart[i+1] are indices into aggToRowMapLO
881 // aggToRowMapLO[aggStart[i]]..aggToRowMapLO[aggStart[i+1]-1] are the DOFs in aggregate i
882 ArrayRCP<LO> aggStart;
883 ArrayRCP<LO> aggToRowMapLO;
884 ArrayRCP<GO> aggToRowMapGO;
885 if (goodMap) {
886 amalgInfo->UnamalgamateAggregatesLO(*aggregates, aggStart, aggToRowMapLO);
887 GetOStream(Runtime1) << "Column map is consistent with the row map, good." << std::endl;
888
889 } else {
890 amalgInfo->UnamalgamateAggregates(*aggregates, aggStart, aggToRowMapGO);
891 GetOStream(Warnings0) << "Column map is not consistent with the row map\n"
892 << "using GO->LO conversion with performance penalty" << std::endl;
893 }
894 coarseNullspace = MultiVectorFactory::Build(coarseMap, NSDim);
895
896 // Pull out the nullspace vectors so that we can have random access.
897 ArrayRCP<ArrayRCP<const SC> > fineNS (NSDim);
898 ArrayRCP<ArrayRCP<SC> > coarseNS(NSDim);
899 for (size_t i = 0; i < NSDim; i++) {
900 fineNS[i] = fineNullspace->getData(i);
901 if (coarseMap->getLocalNumElements() > 0)
902 coarseNS[i] = coarseNullspace->getDataNonConst(i);
903 }
904
905 size_t nnzEstimate = numRows * NSDim;
906
907 // Time to construct the matrix and fill in the values
908 Ptentative = rcp(new CrsMatrixWrap(rowMap, coarseMap, 0));
909 RCP<CrsMatrix> PtentCrs = rcp_dynamic_cast<CrsMatrixWrap>(Ptentative)->getCrsMatrix();
910
911 ArrayRCP<size_t> iaPtent;
912 ArrayRCP<LO> jaPtent;
913 ArrayRCP<SC> valPtent;
914
915 PtentCrs->allocateAllValues(nnzEstimate, iaPtent, jaPtent, valPtent);
916
917 ArrayView<size_t> ia = iaPtent();
918 ArrayView<LO> ja = jaPtent();
919 ArrayView<SC> val = valPtent();
920
921 ia[0] = 0;
922 for (size_t i = 1; i <= numRows; i++)
923 ia[i] = ia[i-1] + NSDim;
924
925 for (size_t j = 0; j < nnzEstimate; j++) {
926 ja [j] = INVALID;
927 val[j] = zero;
928 }
929
930
931 if (doQRStep) {
933 // Standard aggregate-wise QR //
935 for (GO agg = 0; agg < numAggs; agg++) {
936 LO aggSize = aggStart[agg+1] - aggStart[agg];
937
938 Xpetra::global_size_t offset = agg*NSDim;
939
940 // Extract the piece of the nullspace corresponding to the aggregate, and
941 // put it in the flat array, "localQR" (in column major format) for the
942 // QR routine.
943 Teuchos::SerialDenseMatrix<LO,SC> localQR(aggSize, NSDim);
944 if (goodMap) {
945 for (size_t j = 0; j < NSDim; j++)
946 for (LO k = 0; k < aggSize; k++)
947 localQR(k,j) = fineNS[j][aggToRowMapLO[aggStart[agg]+k]];
948 } else {
949 for (size_t j = 0; j < NSDim; j++)
950 for (LO k = 0; k < aggSize; k++)
951 localQR(k,j) = fineNS[j][rowMap->getLocalElement(aggToRowMapGO[aggStart[agg]+k])];
952 }
953
954 // Test for zero columns
955 for (size_t j = 0; j < NSDim; j++) {
956 bool bIsZeroNSColumn = true;
957
958 for (LO k = 0; k < aggSize; k++)
959 if (localQR(k,j) != zero)
960 bIsZeroNSColumn = false;
961
962 TEUCHOS_TEST_FOR_EXCEPTION(bIsZeroNSColumn == true, Exceptions::RuntimeError,
963 "MueLu::TentativePFactory::MakeTentative: fine level NS part has a zero column in NS column " << j);
964 }
965
966 // Calculate QR decomposition (standard)
967 // NOTE: Q is stored in localQR and R is stored in coarseNS
968 if (aggSize >= Teuchos::as<LO>(NSDim)) {
969
970 if (NSDim == 1) {
971 // Only one nullspace vector, calculate Q and R by hand
972 Magnitude norm = STS::magnitude(zero);
973 for (size_t k = 0; k < Teuchos::as<size_t>(aggSize); k++)
974 norm += STS::magnitude(localQR(k,0)*localQR(k,0));
975 norm = Teuchos::ScalarTraits<Magnitude>::squareroot(norm);
976
977 // R = norm
978 coarseNS[0][offset] = norm;
979
980 // Q = localQR(:,0)/norm
981 for (LO i = 0; i < aggSize; i++)
982 localQR(i,0) /= norm;
983
984 } else {
985 Teuchos::SerialQRDenseSolver<LO,SC> qrSolver;
986 qrSolver.setMatrix(Teuchos::rcp(&localQR, false));
987 qrSolver.factor();
988
989 // R = upper triangular part of localQR
990 for (size_t j = 0; j < NSDim; j++)
991 for (size_t k = 0; k <= j; k++)
992 coarseNS[j][offset+k] = localQR(k,j); //TODO is offset+k the correct local ID?!
993
994 // Calculate Q, the tentative prolongator.
995 // The Lapack GEQRF call only works for myAggsize >= NSDim
996 qrSolver.formQ();
997 Teuchos::RCP<Teuchos::SerialDenseMatrix<LO,SC> > qFactor = qrSolver.getQ();
998 for (size_t j = 0; j < NSDim; j++)
999 for (size_t i = 0; i < Teuchos::as<size_t>(aggSize); i++)
1000 localQR(i,j) = (*qFactor)(i,j);
1001 }
1002
1003 } else {
1004 // Special handling for aggSize < NSDim (i.e. single node aggregates in structural mechanics)
1005
1006 // The local QR decomposition is not possible in the "overconstrained"
1007 // case (i.e. number of columns in localQR > number of rows), which
1008 // corresponds to #DOFs in Aggregate < NSDim. For usual problems this
1009 // is only possible for single node aggregates in structural mechanics.
1010 // (Similar problems may arise in discontinuous Galerkin problems...)
1011 // We bypass the QR decomposition and use an identity block in the
1012 // tentative prolongator for the single node aggregate and transfer the
1013 // corresponding fine level null space information 1-to-1 to the coarse
1014 // level null space part.
1015
1016 // NOTE: The resulting tentative prolongation operator has
1017 // (aggSize*DofsPerNode-NSDim) zero columns leading to a singular
1018 // coarse level operator A. To deal with that one has the following
1019 // options:
1020 // - Use the "RepairMainDiagonal" flag in the RAPFactory (default:
1021 // false) to add some identity block to the diagonal of the zero rows
1022 // in the coarse level operator A, such that standard level smoothers
1023 // can be used again.
1024 // - Use special (projection-based) level smoothers, which can deal
1025 // with singular matrices (very application specific)
1026 // - Adapt the code below to avoid zero columns. However, we do not
1027 // support a variable number of DOFs per node in MueLu/Xpetra which
1028 // makes the implementation really hard.
1029
1030 // R = extended (by adding identity rows) localQR
1031 for (size_t j = 0; j < NSDim; j++)
1032 for (size_t k = 0; k < NSDim; k++)
1033 if (k < as<size_t>(aggSize))
1034 coarseNS[j][offset+k] = localQR(k,j);
1035 else
1036 coarseNS[j][offset+k] = (k == j ? one : zero);
1037
1038 // Q = I (rectangular)
1039 for (size_t i = 0; i < as<size_t>(aggSize); i++)
1040 for (size_t j = 0; j < NSDim; j++)
1041 localQR(i,j) = (j == i ? one : zero);
1042 }
1043
1044
1045 // Process each row in the local Q factor
1046 // FIXME: What happens if maps are blocked?
1047 for (LO j = 0; j < aggSize; j++) {
1048 LO localRow = (goodMap ? aggToRowMapLO[aggStart[agg]+j] : rowMap->getLocalElement(aggToRowMapGO[aggStart[agg]+j]));
1049
1050 size_t rowStart = ia[localRow];
1051 for (size_t k = 0, lnnz = 0; k < NSDim; k++) {
1052 // Skip zeros (there may be plenty of them, i.e., NSDim > 1 or boundary conditions)
1053 if (localQR(j,k) != zero) {
1054 ja [rowStart+lnnz] = offset + k;
1055 val[rowStart+lnnz] = localQR(j,k);
1056 lnnz++;
1057 }
1058 }
1059 }
1060 }
1061
1062 } else {
1063 GetOStream(Runtime1) << "TentativePFactory : bypassing local QR phase" << std::endl;
1064 if (NSDim>1)
1065 GetOStream(Warnings0) << "TentativePFactory : for nontrivial nullspace, this may degrade performance" << std::endl;
1067 // "no-QR" option //
1069 // Local Q factor is just the fine nullspace support over the current aggregate.
1070 // Local R factor is the identity.
1071 // TODO I have not implemented any special handling for aggregates that are too
1072 // TODO small to locally support the nullspace, as is done in the standard QR
1073 // TODO case above.
1074 if (goodMap) {
1075 for (GO agg = 0; agg < numAggs; agg++) {
1076 const LO aggSize = aggStart[agg+1] - aggStart[agg];
1077 Xpetra::global_size_t offset = agg*NSDim;
1078
1079 // Process each row in the local Q factor
1080 // FIXME: What happens if maps are blocked?
1081 for (LO j = 0; j < aggSize; j++) {
1082
1083 //TODO Here I do not check for a zero nullspace column on the aggregate.
1084 // as is done in the standard QR case.
1085
1086 const LO localRow = aggToRowMapLO[aggStart[agg]+j];
1087
1088 const size_t rowStart = ia[localRow];
1089
1090 for (size_t k = 0, lnnz = 0; k < NSDim; k++) {
1091 // Skip zeros (there may be plenty of them, i.e., NSDim > 1 or boundary conditions)
1092 SC qr_jk = fineNS[k][aggToRowMapLO[aggStart[agg]+j]];
1093 if(constantColSums) qr_jk = qr_jk / (Magnitude)aggSizes[agg];
1094 if (qr_jk != zero) {
1095 ja [rowStart+lnnz] = offset + k;
1096 val[rowStart+lnnz] = qr_jk;
1097 lnnz++;
1098 }
1099 }
1100 }
1101 for (size_t j = 0; j < NSDim; j++)
1102 coarseNS[j][offset+j] = one;
1103 } //for (GO agg = 0; agg < numAggs; agg++)
1104
1105 } else {
1106 for (GO agg = 0; agg < numAggs; agg++) {
1107 const LO aggSize = aggStart[agg+1] - aggStart[agg];
1108 Xpetra::global_size_t offset = agg*NSDim;
1109 for (LO j = 0; j < aggSize; j++) {
1110
1111 const LO localRow = rowMap->getLocalElement(aggToRowMapGO[aggStart[agg]+j]);
1112
1113 const size_t rowStart = ia[localRow];
1114
1115 for (size_t k = 0, lnnz = 0; k < NSDim; ++k) {
1116 // Skip zeros (there may be plenty of them, i.e., NSDim > 1 or boundary conditions)
1117 SC qr_jk = fineNS[k][rowMap->getLocalElement(aggToRowMapGO[aggStart[agg]+j])];
1118 if(constantColSums) qr_jk = qr_jk / (Magnitude)aggSizes[agg];
1119 if (qr_jk != zero) {
1120 ja [rowStart+lnnz] = offset + k;
1121 val[rowStart+lnnz] = qr_jk;
1122 lnnz++;
1123 }
1124 }
1125 }
1126 for (size_t j = 0; j < NSDim; j++)
1127 coarseNS[j][offset+j] = one;
1128 } //for (GO agg = 0; agg < numAggs; agg++)
1129
1130 } //if (goodmap) else ...
1131
1132 } //if doQRStep ... else
1133
1134 // Compress storage (remove all INVALID, which happen when we skip zeros)
1135 // We do that in-place
1136 size_t ia_tmp = 0, nnz = 0;
1137 for (size_t i = 0; i < numRows; i++) {
1138 for (size_t j = ia_tmp; j < ia[i+1]; j++)
1139 if (ja[j] != INVALID) {
1140 ja [nnz] = ja [j];
1141 val[nnz] = val[j];
1142 nnz++;
1143 }
1144 ia_tmp = ia[i+1];
1145 ia[i+1] = nnz;
1146 }
1147 if (rowMap->lib() == Xpetra::UseTpetra) {
1148 // - Cannot resize for Epetra, as it checks for same pointers
1149 // - Need to resize for Tpetra, as it check ().size() == ia[numRows]
1150 // NOTE: these invalidate ja and val views
1151 jaPtent .resize(nnz);
1152 valPtent.resize(nnz);
1153 }
1154
1155 GetOStream(Runtime1) << "TentativePFactory : aggregates do not cross process boundaries" << std::endl;
1156
1157 PtentCrs->setAllValues(iaPtent, jaPtent, valPtent);
1158
1159
1160 // Managing labels & constants for ESFC
1161 RCP<ParameterList> FCparams;
1162 if(pL.isSublist("matrixmatrix: kernel params"))
1163 FCparams=rcp(new ParameterList(pL.sublist("matrixmatrix: kernel params")));
1164 else
1165 FCparams= rcp(new ParameterList);
1166 // By default, we don't need global constants for TentativeP
1167 FCparams->set("compute global constants",FCparams->get("compute global constants",false));
1168 std::string levelIDs = toString(levelID);
1169 FCparams->set("Timer Label",std::string("MueLu::TentativeP-")+levelIDs);
1170 RCP<const Export> dummy_e;
1171 RCP<const Import> dummy_i;
1172
1173 PtentCrs->expertStaticFillComplete(coarseMap, A->getDomainMap(),dummy_i,dummy_e,FCparams);
1174 }
1175
1176
1177
1178} //namespace MueLu
1179
1180// TODO ReUse: If only P or Nullspace is missing, TentativePFactory can be smart and skip part of the computation.
1181
1182#define MUELU_TENTATIVEPFACTORY_SHORT
1183#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.
void Input(Level &level, const std::string &varName) const
T Get(Level &level, const std::string &varName) const
void Set(Level &level, const std::string &varName, const T &data) const
Class that holds all level-specific information.
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.
static const NoFactory * get()
virtual const Teuchos::ParameterList & GetParameterList() const
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 BuildPuncoupledBlockCrs(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 DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
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
static bool MapsAreNested(const Xpetra::Map< DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &rowMap, const Xpetra::Map< DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &colMap)
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
bool IsPrint(MsgType type, int thisProcRankOnly=-1) const
Find out whether we need to print out information for a specific message type.
Namespace for MueLu classes and methods.
@ Warnings0
Important warning messages (one line).
@ Statistics2
Print even more statistics.
@ Runtime1
Description of what is happening (more verbose).
std::string toString(const T &what)
Little helper function to convert non-string types to strings.