129 RCP<Teuchos::FancyOStream> out;
130 if(
const char* dbg = std::getenv(
"MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
131 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
132 out->setShowAllFrontMatter(
false).setShowProcRank(
true);
134 out = Teuchos::getFancyOStream(rcp(
new Teuchos::oblackholestream()));
137 *out <<
"Starting GeometricInterpolationPFactory::BuildP." << std::endl;
141 const bool removeSmallEntries = pL.get<
bool>(
"interp: remove small entries");
142 const bool buildCoarseCoordinates = pL.get<
bool>(
"interp: build coarse coordinates");
143 const int interpolationOrder =
Get<int>(fineLevel,
"structuredInterpolationOrder");
144 const int numDimensions =
Get<int>(fineLevel,
"numDimensions");
148 RCP<const CrsGraph> prolongatorGraph =
Get<RCP<CrsGraph> >(fineLevel,
"prolongatorGraph");
149 RCP<realvaluedmultivector_type> fineCoordinates, coarseCoordinates;
154 if(buildCoarseCoordinates || (interpolationOrder == 1)) {
156 RCP<const Map> coarseCoordsFineMap =
Get< RCP<const Map> >(fineLevel,
"coarseCoordinatesFineMap");
160 coarseCoordinates = Xpetra::MultiVectorFactory<real_type,LO,GO,Node>::Build(coarseCoordsFineMap,
161 fineCoordinates->getNumVectors());
162 RCP<const Import> coordsImporter = ImportFactory::Build(fineCoordinates->getMap(),
163 coarseCoordsFineMap);
164 coarseCoordinates->doImport(*fineCoordinates, *coordsImporter, Xpetra::INSERT);
165 coarseCoordinates->replaceMap(coarseCoordsMap);
167 Set(coarseLevel,
"Coordinates", coarseCoordinates);
169 if(pL.get<
bool>(
"keep coarse coords")) {
170 coarseLevel.
Set<RCP<realvaluedmultivector_type> >(
"Coordinates2", coarseCoordinates,
NoFactory::get());
174 *out <<
"Fine and coarse coordinates have been loaded from the fine level and set on the coarse level." << std::endl;
177 if(interpolationOrder == 0) {
181 }
else if(interpolationOrder == 1) {
184 RCP<const Map> prolongatorColMap = prolongatorGraph->getColMap();
186 const size_t dofsPerNode =
static_cast<size_t>(A->GetFixedBlockSize());
187 const size_t numColIndices = prolongatorColMap->getLocalNumElements();
188 TEUCHOS_TEST_FOR_EXCEPTION((numColIndices % dofsPerNode) != 0,
190 "Something went wrong, the number of columns in the prolongator is not a multiple of dofsPerNode!");
191 const size_t numGhostCoords = numColIndices / dofsPerNode;
192 const GO indexBase = prolongatorColMap->getIndexBase();
193 const GO coordIndexBase = fineCoordinates->getMap()->getIndexBase();
195 ArrayView<const GO> prolongatorColIndices = prolongatorColMap->getLocalElementList();
196 Array<GO> ghostCoordIndices(numGhostCoords);
197 for(
size_t ghostCoordIdx = 0; ghostCoordIdx < numGhostCoords; ++ghostCoordIdx) {
198 ghostCoordIndices[ghostCoordIdx]
199 = (prolongatorColIndices[ghostCoordIdx*dofsPerNode] - indexBase) / dofsPerNode
202 RCP<Map> ghostCoordMap = MapFactory::Build(fineCoordinates->getMap()->lib(),
203 prolongatorColMap->getGlobalNumElements() / dofsPerNode,
206 fineCoordinates->getMap()->getComm());
208 RCP<realvaluedmultivector_type> ghostCoordinates
209 = Xpetra::MultiVectorFactory<real_type,LO,GO,NO>::Build(ghostCoordMap,
210 fineCoordinates->getNumVectors());
211 RCP<const Import> ghostImporter = ImportFactory::Build(coarseCoordinates->getMap(),
213 ghostCoordinates->doImport(*coarseCoordinates, *ghostImporter, Xpetra::INSERT);
215 BuildLinearP(coarseLevel, A, prolongatorGraph, fineCoordinates,
216 ghostCoordinates, numDimensions, removeSmallEntries, P);
219 *out <<
"The prolongator matrix has been built." << std::endl;
225 RCP<MultiVector> coarseNullspace = MultiVectorFactory::Build(P->getDomainMap(),
226 fineNullspace->getNumVectors());
228 using helpers=Xpetra::Helpers<Scalar,LocalOrdinal,GlobalOrdinal,Node>;
229 if(helpers::isTpetraBlockCrs(A)) {
232 Ptrans->apply(*fineNullspace, *coarseNullspace, Teuchos::NO_TRANS, Teuchos::ScalarTraits<SC>::one(),
233 Teuchos::ScalarTraits<SC>::zero());
236 P->apply(*fineNullspace, *coarseNullspace, Teuchos::TRANS, Teuchos::ScalarTraits<SC>::one(),
237 Teuchos::ScalarTraits<SC>::zero());
240 Set(coarseLevel,
"Nullspace", coarseNullspace);
243 *out <<
"The coarse nullspace is constructed and set on the coarse level." << std::endl;
245 Set(coarseLevel,
"P", P);
247 *out <<
"GeometricInterpolationPFactory::BuildP has completed." << std::endl;
253 BuildConstantP(RCP<Matrix>& P, RCP<const CrsGraph>& prolongatorGraph, RCP<Matrix>& A)
const {
256 RCP<Teuchos::FancyOStream> out;
257 if(
const char* dbg = std::getenv(
"MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
258 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
259 out->setShowAllFrontMatter(
false).setShowProcRank(
true);
261 out = Teuchos::getFancyOStream(rcp(
new Teuchos::oblackholestream()));
264 *out <<
"BuildConstantP" << std::endl;
266 std::vector<size_t> strideInfo(1);
267 strideInfo[0] = A->GetFixedBlockSize();
268 RCP<const StridedMap> stridedDomainMap =
269 StridedMapFactory::Build(prolongatorGraph->getDomainMap(), strideInfo);
271 *out <<
"Call prolongator constructor" << std::endl;
273 using helpers=Xpetra::Helpers<Scalar,LocalOrdinal,GlobalOrdinal,Node>;
274 if(helpers::isTpetraBlockCrs(A)) {
275 SC one = Teuchos::ScalarTraits<SC>::one();
276 SC zero = Teuchos::ScalarTraits<SC>::zero();
277 LO NSDim = A->GetStorageBlockSize();
280 RCP<const Map> BlockMap = prolongatorGraph->getDomainMap();
281 Teuchos::ArrayView<const GO> block_dofs = BlockMap->getLocalElementList();
282 Teuchos::Array<GO> point_dofs(block_dofs.size()*NSDim);
283 for(LO i=0, ct=0; i<block_dofs.size(); i++) {
284 for(LO j=0; j<NSDim; j++) {
285 point_dofs[ct] = block_dofs[i]*NSDim + j;
290 RCP<const Map> PointMap = MapFactory::Build(BlockMap->lib(),
291 BlockMap->getGlobalNumElements() *NSDim,
293 BlockMap->getIndexBase(),
294 BlockMap->getComm());
295 strideInfo[0] = A->GetFixedBlockSize();
296 RCP<const StridedMap> stridedPointMap = StridedMapFactory::Build(PointMap, strideInfo);
298 RCP<Xpetra::CrsMatrix<SC,LO,GO,NO> > P_xpetra = Xpetra::CrsMatrixFactory<SC,LO,GO,NO>::BuildBlock(prolongatorGraph, PointMap, A->getRangeMap(),NSDim);
299 RCP<Xpetra::TpetraBlockCrsMatrix<SC,LO,GO,NO> > P_tpetra = rcp_dynamic_cast<Xpetra::TpetraBlockCrsMatrix<SC,LO,GO,NO> >(P_xpetra);
300 if(P_tpetra.is_null())
throw std::runtime_error(
"BuildConstantP Matrix factory did not return a Tpetra::BlockCrsMatrix");
301 RCP<CrsMatrixWrap> P_wrap = rcp(
new CrsMatrixWrap(P_xpetra));
304 Teuchos::Array<LO> temp(1);
305 Teuchos::ArrayView<const LO> indices;
306 Teuchos::Array<Scalar> block(NSDim*NSDim, zero);
307 for(LO i=0; i<NSDim; i++)
308 block[i*NSDim+i] = one;
309 for(LO i=0; i<(int)prolongatorGraph->getLocalNumRows(); i++) {
310 prolongatorGraph->getLocalRowView(i,indices);
311 for(LO j=0; j<(LO)indices.size();j++) {
312 temp[0] = indices[j];
313 P_tpetra->replaceLocalValues(i,temp(),block());
318 if (A->IsView(
"stridedMaps") ==
true) {
319 P->CreateView(
"stridedMaps", A->getRowMap(
"stridedMaps"), stridedPointMap);
322 P->CreateView(
"stridedMaps", P->getRangeMap(), PointMap);
327 RCP<ParameterList> dummyList = rcp(
new ParameterList());
328 P = rcp(
new CrsMatrixWrap(prolongatorGraph, dummyList));
329 RCP<CrsMatrix> PCrs = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
330 PCrs->setAllToScalar(1.0);
331 PCrs->fillComplete();
334 if (A->IsView(
"stridedMaps") ==
true)
335 P->CreateView(
"stridedMaps", A->getRowMap(
"stridedMaps"), stridedDomainMap);
337 P->CreateView(
"stridedMaps", P->getRangeMap(), stridedDomainMap);
345 RCP<Matrix>& A, RCP<const CrsGraph>& prolongatorGraph,
346 RCP<realvaluedmultivector_type>& fineCoordinates,
347 RCP<realvaluedmultivector_type>& ghostCoordinates,
348 const int numDimensions,
const bool removeSmallEntries,
349 RCP<Matrix>& P)
const {
352 RCP<Teuchos::FancyOStream> out;
353 if(
const char* dbg = std::getenv(
"MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
354 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
355 out->setShowAllFrontMatter(
false).setShowProcRank(
true);
357 out = Teuchos::getFancyOStream(rcp(
new Teuchos::oblackholestream()));
361 const int numInterpolationPoints = 1 << numDimensions;
362 const int dofsPerNode = A->GetFixedBlockSize()/ A->GetStorageBlockSize();;
364 RCP<ParameterList> dummyList = rcp(
new ParameterList());
365 P = rcp(
new CrsMatrixWrap(prolongatorGraph, dummyList));
366 RCP<CrsMatrix> PCrs = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
370 *out <<
"Entering BuildLinearP" << std::endl;
374 const LO numFineNodes = fineCoordinates->getLocalLength();
375 const LO numGhostNodes = ghostCoordinates->getLocalLength();
376 Array<ArrayRCP<const real_type> > fineCoords(3);
377 Array<ArrayRCP<const real_type> > ghostCoords(3);
378 const real_type realZero = Teuchos::as<real_type>(0.0);
379 ArrayRCP<real_type> fineZero(numFineNodes, realZero);
380 ArrayRCP<real_type> ghostZero(numGhostNodes, realZero);
381 for(
int dim = 0; dim < 3; ++dim) {
382 if(dim < numDimensions) {
383 fineCoords[dim] = fineCoordinates->getData(dim);
384 ghostCoords[dim] = ghostCoordinates->getData(dim);
386 fineCoords[dim] = fineZero;
387 ghostCoords[dim] = ghostZero;
391 *out <<
"Coordinates extracted from the multivectors!" << std::endl;
394 LO interpolationNodeIdx = 0, rowIdx = 0;
395 ArrayView<const LO> colIndices;
397 Array<Array<real_type> > coords(numInterpolationPoints + 1);
398 Array<real_type> stencil(numInterpolationPoints);
399 for(LO nodeIdx = 0; nodeIdx < numFineNodes; ++nodeIdx) {
400 if(PCrs->getNumEntriesInLocalRow(nodeIdx*dofsPerNode) == 1) {
403 for(LO dof = 0; dof < dofsPerNode; ++dof) {
404 rowIdx = nodeIdx*dofsPerNode + dof;
405 prolongatorGraph->getLocalRowView(rowIdx, colIndices);
406 PCrs->replaceLocalValues(rowIdx, colIndices, values());
412 for(
int dim = 0; dim < 3; ++dim) {
413 coords[0][dim] = fineCoords[dim][nodeIdx];
415 prolongatorGraph->getLocalRowView(nodeIdx*dofsPerNode, colIndices);
416 for(
int interpolationIdx=0; interpolationIdx < numInterpolationPoints; ++interpolationIdx) {
417 coords[interpolationIdx + 1].resize(3);
418 interpolationNodeIdx = colIndices[interpolationIdx] / dofsPerNode;
419 for(
int dim = 0; dim < 3; ++dim) {
420 coords[interpolationIdx + 1][dim] = ghostCoords[dim][interpolationNodeIdx];
423 RCP<Teuchos::TimeMonitor> tm = rcp(
new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer(
"Compute Linear Interpolation")));
426 values.resize(numInterpolationPoints);
427 for(LO valueIdx = 0; valueIdx < numInterpolationPoints; ++valueIdx) {
428 values[valueIdx] = Teuchos::as<SC>(stencil[valueIdx]);
432 for(LO dof = 0; dof < dofsPerNode; ++dof) {
433 rowIdx = nodeIdx*dofsPerNode + dof;
434 prolongatorGraph->getLocalRowView(rowIdx, colIndices);
435 PCrs->replaceLocalValues(rowIdx, colIndices, values());
441 *out <<
"The calculation of the interpolation stencils has completed." << std::endl;
443 PCrs->fillComplete();
446 *out <<
"All values in P have been set and fillComplete has been performed." << std::endl;
454 if(removeSmallEntries) {
455 *out <<
"Entering remove small entries" << std::endl;
458 ArrayRCP<const size_t> rowptrOrig;
459 ArrayRCP<const LO> colindOrig;
460 ArrayRCP<const Scalar> valuesOrig;
461 PCrs->getAllValues(rowptrOrig, colindOrig, valuesOrig);
463 const size_t numRows =
static_cast<size_t>(rowptrOrig.size() - 1);
464 ArrayRCP<size_t> rowPtr(numRows + 1);
465 ArrayRCP<size_t> nnzOnRows(numRows);
467 size_t countRemovedEntries = 0;
468 for(
size_t rowIdx = 0; rowIdx < numRows; ++rowIdx) {
469 for(
size_t entryIdx = rowptrOrig[rowIdx]; entryIdx < rowptrOrig[rowIdx + 1]; ++entryIdx) {
470 if(Teuchos::ScalarTraits<Scalar>::magnitude(valuesOrig[entryIdx]) < 1e-6) {++countRemovedEntries;}
472 rowPtr[rowIdx + 1] = rowptrOrig[rowIdx + 1] - countRemovedEntries;
473 nnzOnRows[rowIdx] = rowPtr[rowIdx + 1] - rowPtr[rowIdx];
475 GetOStream(
Statistics1) <<
"interp: number of small entries removed= " << countRemovedEntries <<
" / " << rowptrOrig[numRows] << std::endl;
477 size_t countKeptEntries = 0;
478 ArrayRCP<LO> colInd(rowPtr[numRows]);
479 ArrayRCP<SC> values(rowPtr[numRows]);
480 for(
size_t entryIdx = 0; entryIdx < rowptrOrig[numRows]; ++entryIdx) {
481 if(Teuchos::ScalarTraits<Scalar>::magnitude(valuesOrig[entryIdx]) > 1e-6) {
482 colInd[countKeptEntries] = colindOrig[entryIdx];
483 values[countKeptEntries] = valuesOrig[entryIdx];
488 P = rcp(
new CrsMatrixWrap(prolongatorGraph->getRowMap(),
489 prolongatorGraph->getColMap(),
491 RCP<CrsMatrix> PCrsSqueezed = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
492 PCrsSqueezed->resumeFill();
493 PCrsSqueezed->setAllValues(rowPtr, colInd, values);
494 PCrsSqueezed->expertStaticFillComplete(prolongatorGraph->getDomainMap(),
495 prolongatorGraph->getRangeMap());
498 std::vector<size_t> strideInfo(1);
499 strideInfo[0] = dofsPerNode;
500 RCP<const StridedMap> stridedDomainMap =
501 StridedMapFactory::Build(prolongatorGraph->getDomainMap(), strideInfo);
503 *out <<
"The strided maps of P have been computed" << std::endl;
506 if (A->IsView(
"stridedMaps") ==
true) {
507 P->CreateView(
"stridedMaps", A->getRowMap(
"stridedMaps"), stridedDomainMap);
509 P->CreateView(
"stridedMaps", P->getRangeMap(), stridedDomainMap);
518 const Array<Array<real_type> > coord,
519 Array<real_type>& stencil)
const {
536 Teuchos::SerialDenseMatrix<LO,real_type> Jacobian(numDimensions, numDimensions);
537 Teuchos::SerialDenseVector<LO,real_type> residual(numDimensions);
538 Teuchos::SerialDenseVector<LO,real_type> solutionDirection(numDimensions);
539 Teuchos::SerialDenseVector<LO,real_type> paramCoords(numDimensions);
540 Teuchos::SerialDenseSolver<LO,real_type> problem;
541 int iter = 0, max_iter = 5;
542 real_type functions[4][8], norm_ref = 1.0, norm2 = 1.0, tol = 1.0e-5;
543 paramCoords.size(numDimensions);
545 while( (iter < max_iter) && (norm2 > tol*norm_ref) ) {
548 solutionDirection.size(numDimensions);
549 residual.size(numDimensions);
554 for(LO i = 0; i < numDimensions; ++i) {
555 residual(i) = coord[0][i];
556 for(LO k = 0; k < numInterpolationPoints; ++k) {
557 residual(i) -= functions[0][k]*coord[k+1][i];
560 norm_ref += residual(i)*residual(i);
561 if(i == numDimensions - 1) {
562 norm_ref = std::sqrt(norm_ref);
566 for(LO j = 0; j < numDimensions; ++j) {
567 for(LO k = 0; k < numInterpolationPoints; ++k) {
568 Jacobian(i,j) += functions[j+1][k]*coord[k+1][i];
574 problem.setMatrix(Teuchos::rcp(&Jacobian,
false));
575 problem.setVectors(Teuchos::rcp(&solutionDirection,
false), Teuchos::rcp(&residual,
false));
576 if(problem.shouldEquilibrate()) {problem.factorWithEquilibration(
true);}
579 for(LO i = 0; i < numDimensions; ++i) {
580 paramCoords(i) = paramCoords(i) + solutionDirection(i);
585 for(LO i = 0; i < numDimensions; ++i) {
587 for(LO k = 0; k < numInterpolationPoints; ++k) {
588 tmp -= functions[0][k]*coord[k+1][i];
593 norm2 = std::sqrt(norm2);
597 for(LO i = 0; i < numInterpolationPoints; ++i) {
598 stencil[i] = functions[0][i];
606 const Teuchos::SerialDenseVector<LO, real_type> parametricCoordinates,
608 real_type xi = 0.0, eta = 0.0, zeta = 0.0, denominator = 0.0;
609 if(numDimensions == 1) {
610 xi = parametricCoordinates[0];
612 }
else if(numDimensions == 2) {
613 xi = parametricCoordinates[0];
614 eta = parametricCoordinates[1];
616 }
else if(numDimensions == 3) {
617 xi = parametricCoordinates[0];
618 eta = parametricCoordinates[1];
619 zeta = parametricCoordinates[2];
623 functions[0][0] = (1.0 - xi)*(1.0 - eta)*(1.0 - zeta) / denominator;
624 functions[0][1] = (1.0 + xi)*(1.0 - eta)*(1.0 - zeta) / denominator;
625 functions[0][2] = (1.0 - xi)*(1.0 + eta)*(1.0 - zeta) / denominator;
626 functions[0][3] = (1.0 + xi)*(1.0 + eta)*(1.0 - zeta) / denominator;
627 functions[0][4] = (1.0 - xi)*(1.0 - eta)*(1.0 + zeta) / denominator;
628 functions[0][5] = (1.0 + xi)*(1.0 - eta)*(1.0 + zeta) / denominator;
629 functions[0][6] = (1.0 - xi)*(1.0 + eta)*(1.0 + zeta) / denominator;
630 functions[0][7] = (1.0 + xi)*(1.0 + eta)*(1.0 + zeta) / denominator;
632 functions[1][0] = -(1.0 - eta)*(1.0 - zeta) / denominator;
633 functions[1][1] = (1.0 - eta)*(1.0 - zeta) / denominator;
634 functions[1][2] = -(1.0 + eta)*(1.0 - zeta) / denominator;
635 functions[1][3] = (1.0 + eta)*(1.0 - zeta) / denominator;
636 functions[1][4] = -(1.0 - eta)*(1.0 + zeta) / denominator;
637 functions[1][5] = (1.0 - eta)*(1.0 + zeta) / denominator;
638 functions[1][6] = -(1.0 + eta)*(1.0 + zeta) / denominator;
639 functions[1][7] = (1.0 + eta)*(1.0 + zeta) / denominator;
641 functions[2][0] = -(1.0 - xi)*(1.0 - zeta) / denominator;
642 functions[2][1] = -(1.0 + xi)*(1.0 - zeta) / denominator;
643 functions[2][2] = (1.0 - xi)*(1.0 - zeta) / denominator;
644 functions[2][3] = (1.0 + xi)*(1.0 - zeta) / denominator;
645 functions[2][4] = -(1.0 - xi)*(1.0 + zeta) / denominator;
646 functions[2][5] = -(1.0 + xi)*(1.0 + zeta) / denominator;
647 functions[2][6] = (1.0 - xi)*(1.0 + zeta) / denominator;
648 functions[2][7] = (1.0 + xi)*(1.0 + zeta) / denominator;
650 functions[3][0] = -(1.0 - xi)*(1.0 - eta) / denominator;
651 functions[3][1] = -(1.0 + xi)*(1.0 - eta) / denominator;
652 functions[3][2] = -(1.0 - xi)*(1.0 + eta) / denominator;
653 functions[3][3] = -(1.0 + xi)*(1.0 + eta) / denominator;
654 functions[3][4] = (1.0 - xi)*(1.0 - eta) / denominator;
655 functions[3][5] = (1.0 + xi)*(1.0 - eta) / denominator;
656 functions[3][6] = (1.0 - xi)*(1.0 + eta) / denominator;
657 functions[3][7] = (1.0 + xi)*(1.0 + eta) / denominator;