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)) 82 #undef SET_VALID_ENTRY 84 validParamList->set< RCP<const FactoryBase> >(
"A", Teuchos::null,
"Generating factory of the matrix A");
85 validParamList->set< RCP<const FactoryBase> >(
"Aggregates", Teuchos::null,
"Generating factory of the aggregates");
86 validParamList->set< RCP<const FactoryBase> >(
"Nullspace", Teuchos::null,
"Generating factory of the nullspace");
87 validParamList->set< RCP<const FactoryBase> >(
"UnAmalgamationInfo", Teuchos::null,
"Generating factory of UnAmalgamationInfo");
88 validParamList->set< RCP<const FactoryBase> >(
"CoarseMap", Teuchos::null,
"Generating factory of the coarse map");
91 ParameterList norecurse;
92 norecurse.disableRecursiveValidation();
93 validParamList->set<ParameterList> (
"matrixmatrix: kernel params", norecurse,
"MatrixMatrix kernel parameters");
95 return validParamList;
98 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
100 Input(fineLevel,
"A");
101 Input(fineLevel,
"Aggregates");
102 Input(fineLevel,
"Nullspace");
103 Input(fineLevel,
"UnAmalgamationInfo");
104 Input(fineLevel,
"CoarseMap");
107 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
109 return BuildP(fineLevel, coarseLevel);
112 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
116 RCP<Matrix> A = Get< RCP<Matrix> > (fineLevel,
"A");
117 RCP<Aggregates> aggregates = Get< RCP<Aggregates> > (fineLevel,
"Aggregates");
118 RCP<AmalgamationInfo> amalgInfo = Get< RCP<AmalgamationInfo> >(fineLevel,
"UnAmalgamationInfo");
119 RCP<MultiVector> fineNullspace = Get< RCP<MultiVector> > (fineLevel,
"Nullspace");
120 RCP<const Map> coarseMap = Get< RCP<const Map> > (fineLevel,
"CoarseMap");
122 RCP<Matrix> Ptentative;
123 RCP<MultiVector> coarseNullspace;
124 if (!aggregates->AggregatesCrossProcessors())
125 BuildPuncoupled(A, aggregates, amalgInfo, fineNullspace, coarseMap, Ptentative, coarseNullspace,coarseLevel.
GetLevelID());
127 BuildPcoupled (A, aggregates, amalgInfo, fineNullspace, coarseMap, Ptentative, coarseNullspace);
137 if (A->IsView(
"stridedMaps") ==
true)
138 Ptentative->CreateView(
"stridedMaps", A->getRowMap(
"stridedMaps"), coarseMap);
140 Ptentative->CreateView(
"stridedMaps", Ptentative->getRangeMap(), coarseMap);
142 Set(coarseLevel,
"Nullspace", coarseNullspace);
143 Set(coarseLevel,
"P", Ptentative);
146 RCP<ParameterList> params = rcp(
new ParameterList());
147 params->set(
"printLoadBalancingInfo",
true);
152 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
154 BuildPuncoupled(RCP<Matrix> A, RCP<Aggregates> aggregates, RCP<AmalgamationInfo> amalgInfo, RCP<MultiVector> fineNullspace,
155 RCP<const Map> coarseMap, RCP<Matrix>& Ptentative, RCP<MultiVector>& coarseNullspace,
const int levelID)
const {
156 RCP<const Map> rowMap = A->getRowMap();
157 RCP<const Map> colMap = A->getColMap();
159 const size_t numRows = rowMap->getNodeNumElements();
161 typedef Teuchos::ScalarTraits<SC> STS;
162 typedef typename STS::magnitudeType Magnitude;
163 const SC zero = STS::zero();
164 const SC one = STS::one();
165 const LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
167 const GO numAggs = aggregates->GetNumAggregates();
168 const size_t NSDim = fineNullspace->getNumVectors();
173 bool goodMap = isGoodMap(*rowMap, *colMap);
179 ArrayRCP<LO> aggStart;
180 ArrayRCP<LO> aggToRowMapLO;
181 ArrayRCP<GO> aggToRowMapGO;
183 amalgInfo->UnamalgamateAggregatesLO(*aggregates, aggStart, aggToRowMapLO);
184 GetOStream(
Runtime1) <<
"Column map is consistent with the row map, good." << std::endl;
187 amalgInfo->UnamalgamateAggregates(*aggregates, aggStart, aggToRowMapGO);
188 GetOStream(
Warnings0) <<
"Column map is not consistent with the row map\n" 189 <<
"using GO->LO conversion with performance penalty" << std::endl;
192 coarseNullspace = MultiVectorFactory::Build(coarseMap, NSDim);
194 const ParameterList& pL = GetParameterList();
195 const bool &doQRStep = pL.get<
bool>(
"tentative: calculate qr");
198 ArrayRCP<ArrayRCP<const SC> > fineNS (NSDim);
199 ArrayRCP<ArrayRCP<SC> > coarseNS(NSDim);
200 for (
size_t i = 0; i < NSDim; i++) {
201 fineNS[i] = fineNullspace->getData(i);
202 if (coarseMap->getNodeNumElements() > 0)
203 coarseNS[i] = coarseNullspace->getDataNonConst(i);
206 size_t nnzEstimate = numRows * NSDim;
209 Ptentative = rcp(
new CrsMatrixWrap(rowMap, coarseMap, 0, Xpetra::StaticProfile));
210 RCP<CrsMatrix> PtentCrs = rcp_dynamic_cast<CrsMatrixWrap>(Ptentative)->getCrsMatrix();
212 ArrayRCP<size_t> iaPtent;
213 ArrayRCP<LO> jaPtent;
214 ArrayRCP<SC> valPtent;
216 PtentCrs->allocateAllValues(nnzEstimate, iaPtent, jaPtent, valPtent);
218 ArrayView<size_t> ia = iaPtent();
219 ArrayView<LO> ja = jaPtent();
220 ArrayView<SC> val = valPtent();
223 for (
size_t i = 1; i <= numRows; i++)
224 ia[i] = ia[i-1] + NSDim;
226 for (
size_t j = 0; j < nnzEstimate; j++) {
236 for (GO agg = 0; agg < numAggs; agg++) {
237 LO aggSize = aggStart[agg+1] - aggStart[agg];
239 Xpetra::global_size_t offset = agg*NSDim;
244 Teuchos::SerialDenseMatrix<LO,SC> localQR(aggSize, NSDim);
246 for (
size_t j = 0; j < NSDim; j++)
247 for (LO k = 0; k < aggSize; k++)
248 localQR(k,j) = fineNS[j][aggToRowMapLO[aggStart[agg]+k]];
250 for (
size_t j = 0; j < NSDim; j++)
251 for (LO k = 0; k < aggSize; k++)
252 localQR(k,j) = fineNS[j][rowMap->getLocalElement(aggToRowMapGO[aggStart[agg]+k])];
256 for (
size_t j = 0; j < NSDim; j++) {
257 bool bIsZeroNSColumn =
true;
259 for (LO k = 0; k < aggSize; k++)
260 if (localQR(k,j) != zero)
261 bIsZeroNSColumn =
false;
264 "MueLu::TentativePFactory::MakeTentative: fine level NS part has a zero column");
269 if (aggSize >= Teuchos::as<LO>(NSDim)) {
273 Magnitude norm = STS::magnitude(zero);
274 for (
size_t k = 0; k < Teuchos::as<size_t>(aggSize); k++)
275 norm += STS::magnitude(localQR(k,0)*localQR(k,0));
276 norm = Teuchos::ScalarTraits<Magnitude>::squareroot(norm);
279 coarseNS[0][offset] = norm;
282 for (LO i = 0; i < aggSize; i++)
283 localQR(i,0) /= norm;
286 Teuchos::SerialQRDenseSolver<LO,SC> qrSolver;
287 qrSolver.setMatrix(Teuchos::rcp(&localQR,
false));
291 for (
size_t j = 0; j < NSDim; j++)
292 for (
size_t k = 0; k <= j; k++)
293 coarseNS[j][offset+k] = localQR(k,j);
298 Teuchos::RCP<Teuchos::SerialDenseMatrix<LO,SC> > qFactor = qrSolver.getQ();
299 for (
size_t j = 0; j < NSDim; j++)
300 for (
size_t i = 0; i < Teuchos::as<size_t>(aggSize); i++)
301 localQR(i,j) = (*qFactor)(i,j);
332 for (
size_t j = 0; j < NSDim; j++)
333 for (
size_t k = 0; k < NSDim; k++)
334 if (k < as<size_t>(aggSize))
335 coarseNS[j][offset+k] = localQR(k,j);
337 coarseNS[j][offset+k] = (k == j ? one : zero);
340 for (
size_t i = 0; i < as<size_t>(aggSize); i++)
341 for (
size_t j = 0; j < NSDim; j++)
342 localQR(i,j) = (j == i ? one : zero);
348 for (LO j = 0; j < aggSize; j++) {
349 LO localRow = (goodMap ? aggToRowMapLO[aggStart[agg]+j] : rowMap->getLocalElement(aggToRowMapGO[aggStart[agg]+j]));
351 size_t rowStart = ia[localRow];
352 for (
size_t k = 0, lnnz = 0; k < NSDim; k++) {
354 if (localQR(j,k) != zero) {
355 ja [rowStart+lnnz] = offset + k;
356 val[rowStart+lnnz] = localQR(j,k);
364 GetOStream(
Runtime1) <<
"TentativePFactory : bypassing local QR phase" << std::endl;
366 GetOStream(
Warnings0) <<
"TentativePFactor : for nontrivial nullspace, this may degrade performance" << std::endl;
376 for (GO agg = 0; agg < numAggs; agg++) {
377 const LO aggSize = aggStart[agg+1] - aggStart[agg];
378 Xpetra::global_size_t offset = agg*NSDim;
382 for (LO j = 0; j < aggSize; j++) {
387 const LO localRow = aggToRowMapLO[aggStart[agg]+j];
389 const size_t rowStart = ia[localRow];
391 for (
size_t k = 0, lnnz = 0; k < NSDim; k++) {
393 const SC qr_jk = fineNS[k][aggToRowMapLO[aggStart[agg]+j]];
395 ja [rowStart+lnnz] = offset + k;
396 val[rowStart+lnnz] = qr_jk;
401 for (
size_t j = 0; j < NSDim; j++)
402 coarseNS[j][offset+j] = one;
406 for (GO agg = 0; agg < numAggs; agg++) {
407 const LO aggSize = aggStart[agg+1] - aggStart[agg];
408 Xpetra::global_size_t offset = agg*NSDim;
409 for (LO j = 0; j < aggSize; j++) {
411 const LO localRow = rowMap->getLocalElement(aggToRowMapGO[aggStart[agg]+j]);
413 const size_t rowStart = ia[localRow];
415 for (
size_t k = 0, lnnz = 0; k < NSDim; k++) {
417 const SC qr_jk = fineNS[k][rowMap->getLocalElement(aggToRowMapGO[aggStart[agg]+j])];
419 ja [rowStart+lnnz] = offset + k;
420 val[rowStart+lnnz] = qr_jk;
425 for (
size_t j = 0; j < NSDim; j++)
426 coarseNS[j][offset+j] = one;
435 size_t ia_tmp = 0, nnz = 0;
436 for (
size_t i = 0; i < numRows; i++) {
437 for (
size_t j = ia_tmp; j < ia[i+1]; j++)
438 if (ja[j] != INVALID) {
446 if (rowMap->lib() == Xpetra::UseTpetra) {
450 jaPtent .resize(nnz);
451 valPtent.resize(nnz);
454 GetOStream(
Runtime1) <<
"TentativePFactory : aggregates do not cross process boundaries" << std::endl;
456 PtentCrs->setAllValues(iaPtent, jaPtent, valPtent);
460 RCP<ParameterList> FCparams;
461 if(pL.isSublist(
"matrixmatrix: kernel params"))
462 FCparams=rcp(
new ParameterList(pL.sublist(
"matrixmatrix: kernel params")));
464 FCparams= rcp(
new ParameterList);
466 FCparams->set(
"compute global constants",FCparams->get(
"compute global constants",
false));
467 std::string levelIDs =
toString(levelID);
468 FCparams->set(
"Timer Label",std::string(
"MueLu::TentativeP-")+levelIDs);
469 RCP<const Export> dummy_e;
470 RCP<const Import> dummy_i;
472 PtentCrs->expertStaticFillComplete(coarseMap, A->getDomainMap(),dummy_i,dummy_e,FCparams);
475 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
477 BuildPcoupled(RCP<Matrix> A, RCP<Aggregates> aggregates, RCP<AmalgamationInfo> amalgInfo, RCP<MultiVector> fineNullspace,
478 RCP<const Map> coarseMap, RCP<Matrix>& Ptentative, RCP<MultiVector>& coarseNullspace)
const {
479 typedef Teuchos::ScalarTraits<SC> STS;
480 typedef typename STS::magnitudeType Magnitude;
481 const SC zero = STS::zero();
482 const SC one = STS::one();
485 GO numAggs = aggregates->GetNumAggregates();
491 ArrayRCP<LO> aggStart;
492 ArrayRCP< GO > aggToRowMap;
493 amalgInfo->UnamalgamateAggregates(*aggregates, aggStart, aggToRowMap);
497 for (GO i=0; i<numAggs; ++i) {
498 LO sizeOfThisAgg = aggStart[i+1] - aggStart[i];
499 if (sizeOfThisAgg > maxAggSize) maxAggSize = sizeOfThisAgg;
503 const size_t NSDim = fineNullspace->getNumVectors();
506 GO indexBase=A->getRowMap()->getIndexBase();
508 const RCP<const Map> nonUniqueMap = amalgInfo->ComputeUnamalgamatedImportDofMap(*aggregates);
509 const RCP<const Map> uniqueMap = A->getDomainMap();
510 RCP<const Import> importer = ImportFactory::Build(uniqueMap, nonUniqueMap);
511 RCP<MultiVector> fineNullspaceWithOverlap = MultiVectorFactory::Build(nonUniqueMap,NSDim);
512 fineNullspaceWithOverlap->doImport(*fineNullspace,*importer,Xpetra::INSERT);
515 ArrayRCP< ArrayRCP<const SC> > fineNS(NSDim);
516 for (
size_t i=0; i<NSDim; ++i)
517 fineNS[i] = fineNullspaceWithOverlap->getData(i);
520 coarseNullspace = MultiVectorFactory::Build(coarseMap, NSDim);
522 ArrayRCP< ArrayRCP<SC> > coarseNS(NSDim);
523 for (
size_t i=0; i<NSDim; ++i)
524 if (coarseMap->getNodeNumElements() > 0) coarseNS[i] = coarseNullspace->getDataNonConst(i);
529 RCP<const Map > rowMapForPtent = A->getRowMap();
530 const Map& rowMapForPtentRef = *rowMapForPtent;
534 RCP<const Map> colMap = A->getColMap();
536 RCP<const Map > ghostQMap;
537 RCP<MultiVector> ghostQvalues;
538 Array<RCP<Xpetra::Vector<GO,LO,GO,Node> > > ghostQcolumns;
539 RCP<Xpetra::Vector<GO,LO,GO,Node> > ghostQrowNums;
540 ArrayRCP< ArrayRCP<SC> > ghostQvals;
541 ArrayRCP< ArrayRCP<GO> > ghostQcols;
542 ArrayRCP< GO > ghostQrows;
545 for (LO j=0; j<numAggs; ++j) {
546 for (LO k=aggStart[j]; k<aggStart[j+1]; ++k) {
547 if (rowMapForPtentRef.isNodeGlobalElement(aggToRowMap[k]) ==
false) {
548 ghostGIDs.push_back(aggToRowMap[k]);
552 ghostQMap = MapFactory::Build(A->getRowMap()->lib(),
553 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
555 indexBase, A->getRowMap()->getComm());
557 ghostQvalues = MultiVectorFactory::Build(ghostQMap,NSDim);
561 ghostQcolumns.resize(NSDim);
562 for (
size_t i=0; i<NSDim; ++i)
563 ghostQcolumns[i] = Xpetra::VectorFactory<GO,LO,GO,Node>::Build(ghostQMap);
564 ghostQrowNums = Xpetra::VectorFactory<GO,LO,GO,Node>::Build(ghostQMap);
565 if (ghostQvalues->getLocalLength() > 0) {
566 ghostQvals.resize(NSDim);
567 ghostQcols.resize(NSDim);
568 for (
size_t i=0; i<NSDim; ++i) {
569 ghostQvals[i] = ghostQvalues->getDataNonConst(i);
570 ghostQcols[i] = ghostQcolumns[i]->getDataNonConst(0);
572 ghostQrows = ghostQrowNums->getDataNonConst(0);
576 importer = ImportFactory::Build(ghostQMap, A->getRowMap());
579 Teuchos::SerialQRDenseSolver<LO,SC> qrSolver;
582 Array<GO> globalColPtr(maxAggSize*NSDim,0);
583 Array<LO> localColPtr(maxAggSize*NSDim,0);
584 Array<SC> valPtr(maxAggSize*NSDim,0.);
587 const Map& coarseMapRef = *coarseMap;
590 ArrayRCP<size_t> ptent_rowptr;
591 ArrayRCP<LO> ptent_colind;
592 ArrayRCP<Scalar> ptent_values;
595 ArrayView<size_t> rowptr_v;
596 ArrayView<LO> colind_v;
597 ArrayView<Scalar> values_v;
600 Array<size_t> rowptr_temp;
601 Array<LO> colind_temp;
602 Array<Scalar> values_temp;
604 RCP<CrsMatrix> PtentCrs;
606 RCP<CrsMatrixWrap> PtentCrsWrap = rcp(
new CrsMatrixWrap(rowMapForPtent, NSDim, Xpetra::StaticProfile));
607 PtentCrs = PtentCrsWrap->getCrsMatrix();
608 Ptentative = PtentCrsWrap;
614 const Map& nonUniqueMapRef = *nonUniqueMap;
616 size_t total_nnz_count=0;
618 for (GO agg=0; agg<numAggs; ++agg)
620 LO myAggSize = aggStart[agg+1]-aggStart[agg];
623 Teuchos::SerialDenseMatrix<LO,SC> localQR(myAggSize, NSDim);
624 for (
size_t j=0; j<NSDim; ++j) {
625 bool bIsZeroNSColumn =
true;
626 for (LO k=0; k<myAggSize; ++k)
631 SC nsVal = fineNS[j][ nonUniqueMapRef.getLocalElement(aggToRowMap[aggStart[agg]+k]) ];
632 localQR(k,j) = nsVal;
633 if (nsVal != zero) bIsZeroNSColumn =
false;
636 GetOStream(
Runtime1,-1) <<
"length of fine level nsp: " << fineNullspace->getGlobalLength() << std::endl;
637 GetOStream(
Runtime1,-1) <<
"length of fine level nsp w overlap: " << fineNullspaceWithOverlap->getGlobalLength() << std::endl;
638 GetOStream(
Runtime1,-1) <<
"(local?) aggId=" << agg << std::endl;
639 GetOStream(
Runtime1,-1) <<
"aggSize=" << myAggSize << std::endl;
640 GetOStream(
Runtime1,-1) <<
"agg DOF=" << k << std::endl;
641 GetOStream(
Runtime1,-1) <<
"NS vector j=" << j << std::endl;
642 GetOStream(
Runtime1,-1) <<
"j*myAggSize + k = " << j*myAggSize + k << std::endl;
643 GetOStream(
Runtime1,-1) <<
"aggToRowMap["<<agg<<
"][" << k <<
"] = " << aggToRowMap[aggStart[agg]+k] << std::endl;
644 GetOStream(
Runtime1,-1) <<
"id aggToRowMap[agg][k]=" << aggToRowMap[aggStart[agg]+k] <<
" is global element in nonUniqueMap = " <<
645 nonUniqueMapRef.isNodeGlobalElement(aggToRowMap[aggStart[agg]+k]) << std::endl;
646 GetOStream(
Runtime1,-1) <<
"colMap local id aggToRowMap[agg][k]=" << nonUniqueMapRef.getLocalElement(aggToRowMap[aggStart[agg]+k]) << std::endl;
647 GetOStream(
Runtime1,-1) <<
"fineNS...=" << fineNS[j][ nonUniqueMapRef.getLocalElement(aggToRowMap[aggStart[agg]+k]) ] << std::endl;
648 GetOStream(
Errors,-1) <<
"caught an error!" << std::endl;
651 TEUCHOS_TEST_FOR_EXCEPTION(bIsZeroNSColumn ==
true,
Exceptions::RuntimeError,
"MueLu::TentativePFactory::MakeTentative: fine level NS part has a zero column. Error.");
654 Xpetra::global_size_t offset=agg*NSDim;
656 if(myAggSize >= Teuchos::as<LocalOrdinal>(NSDim)) {
661 SC tau = localQR(0,0);
666 for (
size_t k = 0; k < Teuchos::as<size_t>(myAggSize); ++k) {
667 Magnitude tmag = STS::magnitude(localQR(k,0));
670 dtemp = Teuchos::ScalarTraits<Magnitude>::squareroot(dtemp);
672 localQR(0,0) = dtemp;
674 qrSolver.setMatrix( Teuchos::rcp(&localQR,
false) );
681 for (
size_t j=0; j<NSDim; ++j) {
682 for (
size_t k=0; k<=j; ++k) {
684 if (coarseMapRef.isNodeLocalElement(offset+k)) {
685 coarseNS[j][offset+k] = localQR(k, j);
689 GetOStream(
Errors,-1) <<
"caught error in coarseNS insert, j="<<j<<
", offset+k = "<<offset+k<<std::endl;
699 Magnitude dtemp = Teuchos::ScalarTraits<SC>::magnitude(localQR(0,0));
702 for (LocalOrdinal i=0; i<myAggSize; ++i) {
703 localQR(i,0) *= dtemp ;
707 Teuchos::RCP<Teuchos::SerialDenseMatrix<LO,SC> > qFactor = qrSolver.getQ();
708 for (
size_t j=0; j<NSDim; j++) {
709 for (
size_t i = 0; i < Teuchos::as<size_t>(myAggSize); i++) {
710 localQR(i,j) = (*qFactor)(i,j);
720 for (
size_t j = 0; j < NSDim; j++)
721 for (
size_t k = 0; k < NSDim; k++) {
723 "Caught error in coarseNS insert, j=" << j <<
", offset+k = " << offset+k);
725 if (k < as<size_t>(myAggSize))
726 coarseNS[j][offset+k] = localQR(k,j);
728 coarseNS[j][offset+k] = (k == j ? one : zero);
732 for (
size_t i = 0; i < as<size_t>(myAggSize); i++)
733 for (
size_t j = 0; j < NSDim; j++)
734 localQR(i,j) = (j == i ? one : zero);
741 for (GO j=0; j<myAggSize; ++j) {
745 GO globalRow = aggToRowMap[aggStart[agg]+j];
748 if (rowMapForPtentRef.isNodeGlobalElement(globalRow) == false ) {
749 ghostQrows[qctr] = globalRow;
750 for (
size_t k=0; k<NSDim; ++k) {
751 ghostQcols[k][qctr] = coarseMapRef.getGlobalElement(agg*NSDim+k);
752 ghostQvals[k][qctr] = localQR(j,k);
757 for (
size_t k=0; k<NSDim; ++k) {
759 if (localQR(j,k) != Teuchos::ScalarTraits<SC>::zero()) {
760 localColPtr[nnz] = agg * NSDim + k;
761 globalColPtr[nnz] = coarseMapRef.getGlobalElement(localColPtr[nnz]);
762 valPtr[nnz] = localQR(j,k);
768 GetOStream(
Errors,-1) <<
"caught error in colPtr/valPtr insert, current index="<<nnz<<std::endl;
773 Ptentative->insertGlobalValues(globalRow,globalColPtr.view(0,nnz),valPtr.view(0,nnz));
776 GetOStream(
Errors,-1) <<
"pid " << A->getRowMap()->getComm()->getRank()
777 <<
"caught error during Ptent row insertion, global row " 778 << globalRow << std::endl;
789 GetOStream(
Runtime1) <<
"TentativePFactory : aggregates may cross process boundaries" << std::endl;
792 RCP<Xpetra::Vector<GO,LO,GO,Node> > targetQrowNums = Xpetra::VectorFactory<GO,LO,GO,Node>::Build(rowMapForPtent);
793 targetQrowNums->putScalar(-1);
794 targetQrowNums->doImport(*ghostQrowNums,*importer,Xpetra::INSERT);
795 ArrayRCP< GO > targetQrows = targetQrowNums->getDataNonConst(0);
798 Array<GO> gidsToImport;
799 gidsToImport.reserve(targetQrows.size());
800 for (
typename ArrayRCP<GO>::iterator r=targetQrows.begin(); r!=targetQrows.end(); ++r) {
802 gidsToImport.push_back(*r);
805 RCP<const Map > reducedMap = MapFactory::Build( A->getRowMap()->lib(),
806 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
807 gidsToImport, indexBase, A->getRowMap()->getComm() );
810 importer = ImportFactory::Build(ghostQMap, reducedMap);
812 Array<RCP<Xpetra::Vector<GO,LO,GO,Node> > > targetQcolumns(NSDim);
813 for (
size_t i=0; i<NSDim; ++i) {
814 targetQcolumns[i] = Xpetra::VectorFactory<GO,LO,GO,Node>::Build(reducedMap);
815 targetQcolumns[i]->doImport(*(ghostQcolumns[i]),*importer,Xpetra::INSERT);
817 RCP<MultiVector> targetQvalues = MultiVectorFactory::Build(reducedMap,NSDim);
818 targetQvalues->doImport(*ghostQvalues,*importer,Xpetra::INSERT);
820 ArrayRCP< ArrayRCP<SC> > targetQvals;
821 ArrayRCP<ArrayRCP<GO> > targetQcols;
822 if (targetQvalues->getLocalLength() > 0) {
823 targetQvals.resize(NSDim);
824 targetQcols.resize(NSDim);
825 for (
size_t i=0; i<NSDim; ++i) {
826 targetQvals[i] = targetQvalues->getDataNonConst(i);
827 targetQcols[i] = targetQcolumns[i]->getDataNonConst(0);
831 valPtr = Array<SC>(NSDim,0.);
832 globalColPtr = Array<GO>(NSDim,0);
833 for (
typename Array<GO>::iterator r=gidsToImport.begin(); r!=gidsToImport.end(); ++r) {
834 if (targetQvalues->getLocalLength() > 0) {
835 for (
size_t j=0; j<NSDim; ++j) {
836 valPtr[j] = targetQvals[j][reducedMap->getLocalElement(*r)];
837 globalColPtr[j] = targetQcols[j][reducedMap->getLocalElement(*r)];
839 Ptentative->insertGlobalValues(*r, globalColPtr.view(0,NSDim), valPtr.view(0,NSDim));
843 Ptentative->fillComplete(coarseMap, A->getDomainMap());
846 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
848 ArrayView<const GO> rowElements = rowMap.getNodeElementList();
849 ArrayView<const GO> colElements = colMap.getNodeElementList();
851 const size_t numElements = rowElements.size();
854 for (
size_t i = 0; i < numElements; i++)
855 if (rowElements[i] != colElements[i]) {
867 #define MUELU_TENTATIVEPFACTORY_SHORT 868 #endif // MUELU_TENTATIVEPFACTORY_DEF_HPP Important warning messages (one line)
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
Timer to be used in factories. Similar to Monitor but with additional timers.
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Namespace for MueLu classes and methods.
int GetLevelID() const
Return level number.
#define SET_VALID_ENTRY(name)
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
Class that holds all level-specific information.
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
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
Exception throws to report errors in the internal logical of the program.
Description of what is happening (more verbose)
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
bool isGoodMap(const Map &rowMap, const Map &colMap) const