46 #ifndef MUELU_TENTATIVEPFACTORY_DEF_HPP
47 #define MUELU_TENTATIVEPFACTORY_DEF_HPP
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>
64 #include "MueLu_Aggregates.hpp"
65 #include "MueLu_AmalgamationFactory.hpp"
66 #include "MueLu_AmalgamationInfo.hpp"
67 #include "MueLu_CoarseMapFactory.hpp"
70 #include "MueLu_NullspaceFactory.hpp"
71 #include "MueLu_PerfUtils.hpp"
72 #include "MueLu_Utilities.hpp"
76 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
78 RCP<ParameterList> validParamList = rcp(
new ParameterList());
80 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
83 #undef SET_VALID_ENTRY
84 validParamList->set< std::string >(
"Nullspace name",
"Nullspace",
"Name for the input nullspace");
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");
96 ParameterList norecurse;
97 norecurse.disableRecursiveValidation();
98 validParamList->set<ParameterList> (
"matrixmatrix: kernel params", norecurse,
"MatrixMatrix kernel parameters");
100 return validParamList;
103 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
106 const ParameterList& pL = GetParameterList();
108 std::string nspName =
"Nullspace";
109 if(pL.isParameter(
"Nullspace name")) nspName = pL.get<std::string>(
"Nullspace name");
111 Input(fineLevel,
"A");
112 Input(fineLevel,
"Aggregates");
113 Input(fineLevel, nspName);
114 Input(fineLevel,
"UnAmalgamationInfo");
115 Input(fineLevel,
"CoarseMap");
118 pL.get<
bool>(
"tentative: build coarse coordinates") ) {
119 bTransferCoordinates_ =
true;
120 Input(fineLevel,
"Coordinates");
121 }
else if (bTransferCoordinates_) {
122 Input(fineLevel,
"Coordinates");
126 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
128 return BuildP(fineLevel, coarseLevel);
131 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
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;
138 const ParameterList& pL = GetParameterList();
139 std::string nspName =
"Nullspace";
140 if(pL.isParameter(
"Nullspace name")) nspName = pL.get<std::string>(
"Nullspace name");
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");
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);
159 TEUCHOS_TEST_FOR_EXCEPTION(A->getRowMap()->getNodeNumElements() != fineNullspace->getMap()->getNodeNumElements(),
162 RCP<Matrix> Ptentative;
163 RCP<MultiVector> coarseNullspace;
164 RCP<RealValuedMultiVector> coarseCoords;
166 if(bTransferCoordinates_) {
169 ArrayView<const GO> elementAList = coarseMap->getNodeElementList();
171 if (rcp_dynamic_cast<const StridedMap>(coarseMap) != Teuchos::null) {
172 blkSize = rcp_dynamic_cast<const StridedMap>(coarseMap)->getFixedBlockSize();
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();
179 for (LO i = 0; i < numCoarseNodes; i++) {
180 nodeList[i] = (elementAList[i*blkSize]-indexBase)/blkSize + indexBase;
182 RCP<const Map> coarseCoordsMap = MapFactory::Build(fineCoords->getMap()->lib(),
183 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
186 fineCoords->getMap()->getComm());
187 coarseCoords = RealValuedMultiVectorFactory::Build(coarseCoordsMap, numDimensions);
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);
195 ghostedCoords = RealValuedMultiVectorFactory::Build(aggMap, numDimensions);
196 ghostedCoords->doImport(*fineCoords, *importer, Xpetra::INSERT);
198 ghostedCoords = fineCoords;
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);
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);
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];
221 for (LO agg = 0; agg < numAggs; agg++) {
222 coarseCoordsData[agg] /= aggSizes[agg];
227 if (!aggregates->AggregatesCrossProcessors())
228 BuildPuncoupled(A, aggregates, amalgInfo, fineNullspace, coarseMap, Ptentative, coarseNullspace,coarseLevel.
GetLevelID());
230 BuildPcoupled (A, aggregates, amalgInfo, fineNullspace, coarseMap, Ptentative, coarseNullspace);
240 if (A->IsView(
"stridedMaps") ==
true)
241 Ptentative->CreateView(
"stridedMaps", A->getRowMap(
"stridedMaps"), coarseMap);
243 Ptentative->CreateView(
"stridedMaps", Ptentative->getRangeMap(), coarseMap);
245 if(bTransferCoordinates_) {
246 Set(coarseLevel,
"Coordinates", coarseCoords);
248 Set(coarseLevel,
"Nullspace", coarseNullspace);
249 Set(coarseLevel,
"P", Ptentative);
252 RCP<ParameterList> params = rcp(
new ParameterList());
253 params->set(
"printLoadBalancingInfo",
true);
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();
265 const size_t numRows = rowMap->getNodeNumElements();
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();
273 const GO numAggs = aggregates->GetNumAggregates();
274 const size_t NSDim = fineNullspace->getNumVectors();
279 bool goodMap = isGoodMap(*rowMap, *colMap);
285 ArrayRCP<LO> aggStart;
286 ArrayRCP<LO> aggToRowMapLO;
287 ArrayRCP<GO> aggToRowMapGO;
289 amalgInfo->UnamalgamateAggregatesLO(*aggregates, aggStart, aggToRowMapLO);
290 GetOStream(
Runtime1) <<
"Column map is consistent with the row map, good." << std::endl;
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;
298 coarseNullspace = MultiVectorFactory::Build(coarseMap, NSDim);
300 const ParameterList& pL = GetParameterList();
301 const bool &doQRStep = pL.get<
bool>(
"tentative: calculate qr");
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);
312 size_t nnzEstimate = numRows * NSDim;
315 Ptentative = rcp(
new CrsMatrixWrap(rowMap, coarseMap, 0, Xpetra::StaticProfile));
316 RCP<CrsMatrix> PtentCrs = rcp_dynamic_cast<CrsMatrixWrap>(Ptentative)->getCrsMatrix();
318 ArrayRCP<size_t> iaPtent;
319 ArrayRCP<LO> jaPtent;
320 ArrayRCP<SC> valPtent;
322 PtentCrs->allocateAllValues(nnzEstimate, iaPtent, jaPtent, valPtent);
324 ArrayView<size_t> ia = iaPtent();
325 ArrayView<LO> ja = jaPtent();
326 ArrayView<SC> val = valPtent();
329 for (
size_t i = 1; i <= numRows; i++)
330 ia[i] = ia[i-1] + NSDim;
332 for (
size_t j = 0; j < nnzEstimate; j++) {
342 for (GO agg = 0; agg < numAggs; agg++) {
343 LO aggSize = aggStart[agg+1] - aggStart[agg];
345 Xpetra::global_size_t offset = agg*NSDim;
350 Teuchos::SerialDenseMatrix<LO,SC> localQR(aggSize, NSDim);
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]];
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])];
362 for (
size_t j = 0; j < NSDim; j++) {
363 bool bIsZeroNSColumn =
true;
365 for (LO k = 0; k < aggSize; k++)
366 if (localQR(k,j) != zero)
367 bIsZeroNSColumn =
false;
370 "MueLu::TentativePFactory::MakeTentative: fine level NS part has a zero column");
375 if (aggSize >= Teuchos::as<LO>(NSDim)) {
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);
385 coarseNS[0][offset] = norm;
388 for (LO i = 0; i < aggSize; i++)
389 localQR(i,0) /= norm;
392 Teuchos::SerialQRDenseSolver<LO,SC> qrSolver;
393 qrSolver.setMatrix(Teuchos::rcp(&localQR,
false));
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);
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);
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);
443 coarseNS[j][offset+k] = (k == j ? one : zero);
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);
454 for (LO j = 0; j < aggSize; j++) {
455 LO localRow = (goodMap ? aggToRowMapLO[aggStart[agg]+j] : rowMap->getLocalElement(aggToRowMapGO[aggStart[agg]+j]));
457 size_t rowStart = ia[localRow];
458 for (
size_t k = 0, lnnz = 0; k < NSDim; k++) {
460 if (localQR(j,k) != zero) {
461 ja [rowStart+lnnz] = offset + k;
462 val[rowStart+lnnz] = localQR(j,k);
470 GetOStream(
Runtime1) <<
"TentativePFactory : bypassing local QR phase" << std::endl;
472 GetOStream(
Warnings0) <<
"TentativePFactor : for nontrivial nullspace, this may degrade performance" << std::endl;
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;
488 for (LO j = 0; j < aggSize; j++) {
493 const LO localRow = aggToRowMapLO[aggStart[agg]+j];
495 const size_t rowStart = ia[localRow];
497 for (
size_t k = 0, lnnz = 0; k < NSDim; k++) {
499 const SC qr_jk = fineNS[k][aggToRowMapLO[aggStart[agg]+j]];
501 ja [rowStart+lnnz] = offset + k;
502 val[rowStart+lnnz] = qr_jk;
507 for (
size_t j = 0; j < NSDim; j++)
508 coarseNS[j][offset+j] = one;
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++) {
517 const LO localRow = rowMap->getLocalElement(aggToRowMapGO[aggStart[agg]+j]);
519 const size_t rowStart = ia[localRow];
521 for (
size_t k = 0, lnnz = 0; k < NSDim; k++) {
523 const SC qr_jk = fineNS[k][rowMap->getLocalElement(aggToRowMapGO[aggStart[agg]+j])];
525 ja [rowStart+lnnz] = offset + k;
526 val[rowStart+lnnz] = qr_jk;
531 for (
size_t j = 0; j < NSDim; j++)
532 coarseNS[j][offset+j] = one;
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) {
552 if (rowMap->lib() == Xpetra::UseTpetra) {
556 jaPtent .resize(nnz);
557 valPtent.resize(nnz);
560 GetOStream(
Runtime1) <<
"TentativePFactory : aggregates do not cross process boundaries" << std::endl;
562 PtentCrs->setAllValues(iaPtent, jaPtent, valPtent);
566 RCP<ParameterList> FCparams;
567 if(pL.isSublist(
"matrixmatrix: kernel params"))
568 FCparams=rcp(
new ParameterList(pL.sublist(
"matrixmatrix: kernel params")));
570 FCparams= rcp(
new ParameterList);
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;
578 PtentCrs->expertStaticFillComplete(coarseMap, A->getDomainMap(),dummy_i,dummy_e,FCparams);
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();
591 GO numAggs = aggregates->GetNumAggregates();
597 ArrayRCP<LO> aggStart;
598 ArrayRCP< GO > aggToRowMap;
599 amalgInfo->UnamalgamateAggregates(*aggregates, aggStart, aggToRowMap);
603 for (GO i=0; i<numAggs; ++i) {
604 LO sizeOfThisAgg = aggStart[i+1] - aggStart[i];
605 if (sizeOfThisAgg > maxAggSize) maxAggSize = sizeOfThisAgg;
609 const size_t NSDim = fineNullspace->getNumVectors();
612 GO indexBase=A->getRowMap()->getIndexBase();
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);
621 ArrayRCP< ArrayRCP<const SC> > fineNS(NSDim);
622 for (
size_t i=0; i<NSDim; ++i)
623 fineNS[i] = fineNullspaceWithOverlap->getData(i);
626 coarseNullspace = MultiVectorFactory::Build(coarseMap, NSDim);
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);
635 RCP<const Map > rowMapForPtent = A->getRowMap();
636 const Map& rowMapForPtentRef = *rowMapForPtent;
640 RCP<const Map> colMap = A->getColMap();
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;
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]);
658 ghostQMap = MapFactory::Build(A->getRowMap()->lib(),
659 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
661 indexBase, A->getRowMap()->getComm());
663 ghostQvalues = MultiVectorFactory::Build(ghostQMap,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);
678 ghostQrows = ghostQrowNums->getDataNonConst(0);
682 importer = ImportFactory::Build(ghostQMap, A->getRowMap());
685 Teuchos::SerialQRDenseSolver<LO,SC> qrSolver;
688 Array<GO> globalColPtr(maxAggSize*NSDim,0);
689 Array<LO> localColPtr(maxAggSize*NSDim,0);
690 Array<SC> valPtr(maxAggSize*NSDim,0.);
693 const Map& coarseMapRef = *coarseMap;
696 ArrayRCP<size_t> ptent_rowptr;
697 ArrayRCP<LO> ptent_colind;
698 ArrayRCP<Scalar> ptent_values;
701 ArrayView<size_t> rowptr_v;
702 ArrayView<LO> colind_v;
703 ArrayView<Scalar> values_v;
706 Array<size_t> rowptr_temp;
707 Array<LO> colind_temp;
708 Array<Scalar> values_temp;
710 RCP<CrsMatrix> PtentCrs;
712 RCP<CrsMatrixWrap> PtentCrsWrap = rcp(
new CrsMatrixWrap(rowMapForPtent, NSDim, Xpetra::StaticProfile));
713 PtentCrs = PtentCrsWrap->getCrsMatrix();
714 Ptentative = PtentCrsWrap;
720 const Map& nonUniqueMapRef = *nonUniqueMap;
722 size_t total_nnz_count=0;
724 for (GO agg=0; agg<numAggs; ++agg)
726 LO myAggSize = aggStart[agg+1]-aggStart[agg];
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)
737 SC nsVal = fineNS[j][ nonUniqueMapRef.getLocalElement(aggToRowMap[aggStart[agg]+k]) ];
738 localQR(k,j) = nsVal;
739 if (nsVal != zero) bIsZeroNSColumn =
false;
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;
757 TEUCHOS_TEST_FOR_EXCEPTION(bIsZeroNSColumn ==
true,
Exceptions::RuntimeError,
"MueLu::TentativePFactory::MakeTentative: fine level NS part has a zero column. Error.");
760 Xpetra::global_size_t offset=agg*NSDim;
762 if(myAggSize >= Teuchos::as<LocalOrdinal>(NSDim)) {
767 SC tau = localQR(0,0);
772 for (
size_t k = 0; k < Teuchos::as<size_t>(myAggSize); ++k) {
773 Magnitude tmag = STS::magnitude(localQR(k,0));
776 dtemp = Teuchos::ScalarTraits<Magnitude>::squareroot(dtemp);
778 localQR(0,0) = dtemp;
780 qrSolver.setMatrix( Teuchos::rcp(&localQR,
false) );
787 for (
size_t j=0; j<NSDim; ++j) {
788 for (
size_t k=0; k<=j; ++k) {
790 if (coarseMapRef.isNodeLocalElement(offset+k)) {
791 coarseNS[j][offset+k] = localQR(k, j);
795 GetOStream(
Errors,-1) <<
"caught error in coarseNS insert, j="<<j<<
", offset+k = "<<offset+k<<std::endl;
805 Magnitude dtemp = Teuchos::ScalarTraits<SC>::magnitude(localQR(0,0));
809 localQR(i,0) *= dtemp ;
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);
826 for (
size_t j = 0; j < NSDim; j++)
827 for (
size_t k = 0; k < NSDim; k++) {
829 "Caught error in coarseNS insert, j=" << j <<
", offset+k = " << offset+k);
831 if (k < as<size_t>(myAggSize))
832 coarseNS[j][offset+k] = localQR(k,j);
834 coarseNS[j][offset+k] = (k == j ? one : zero);
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);
847 for (GO j=0; j<myAggSize; ++j) {
851 GO globalRow = aggToRowMap[aggStart[agg]+j];
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);
863 for (
size_t k=0; k<NSDim; ++k) {
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);
874 GetOStream(
Errors,-1) <<
"caught error in colPtr/valPtr insert, current index="<<nnz<<std::endl;
879 Ptentative->insertGlobalValues(globalRow,globalColPtr.view(0,nnz),valPtr.view(0,nnz));
882 GetOStream(
Errors,-1) <<
"pid " << A->getRowMap()->getComm()->getRank()
883 <<
"caught error during Ptent row insertion, global row "
884 << globalRow << std::endl;
895 GetOStream(
Runtime1) <<
"TentativePFactory : aggregates may cross process boundaries" << std::endl;
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);
904 Array<GO> gidsToImport;
905 gidsToImport.reserve(targetQrows.size());
906 for (
typename ArrayRCP<GO>::iterator r=targetQrows.begin(); r!=targetQrows.end(); ++r) {
908 gidsToImport.push_back(*r);
911 RCP<const Map > reducedMap = MapFactory::Build( A->getRowMap()->lib(),
912 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
913 gidsToImport, indexBase, A->getRowMap()->getComm() );
916 importer = ImportFactory::Build(ghostQMap, reducedMap);
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);
923 RCP<MultiVector> targetQvalues = MultiVectorFactory::Build(reducedMap,NSDim);
924 targetQvalues->doImport(*ghostQvalues,*importer,Xpetra::INSERT);
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);
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)];
945 Ptentative->insertGlobalValues(*r, globalColPtr.view(0,NSDim), valPtr.view(0,NSDim));
949 Ptentative->fillComplete(coarseMap, A->getDomainMap());
952 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
954 ArrayView<const GO> rowElements = rowMap.getNodeElementList();
955 ArrayView<const GO> colElements = colMap.getNodeElementList();
957 const size_t numElements = rowElements.size();
960 for (
size_t i = 0; i < numElements; i++)
961 if (rowElements[i] != colElements[i]) {
973 #define MUELU_TENTATIVEPFACTORY_SHORT
#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.
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()
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.
@ Runtime1
Description of what is happening (more verbose)
std::string toString(const T &what)
Little helper function to convert non-string types to strings.