46#ifndef MUELU_GEOMETRICINTERPOLATIONPFACTORY_KOKKOS_DEF_HPP
47#define MUELU_GEOMETRICINTERPOLATIONPFACTORY_KOKKOS_DEF_HPP
49#include "Xpetra_CrsGraph.hpp"
50#include "Xpetra_CrsMatrixUtils.hpp"
54#include "MueLu_IndexManager_kokkos.hpp"
56#ifdef HAVE_MUELU_TPETRA
57#include "Xpetra_TpetraCrsMatrix.hpp"
66 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
68 RCP<ParameterList> validParamList = rcp(
new ParameterList());
70#define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
75 validParamList->set<RCP<const FactoryBase> >(
"A", Teuchos::null,
76 "Generating factory of the matrix A");
77 validParamList->set<RCP<const FactoryBase> >(
"prolongatorGraph", Teuchos::null,
78 "Graph generated by StructuredAggregationFactory used to construct a piece-linear prolongator.");
79 validParamList->set<RCP<const FactoryBase> >(
"Coordinates", Teuchos::null,
80 "Fine level coordinates used to construct piece-wise linear prolongator and coarse level coordinates.");
81 validParamList->set<RCP<const FactoryBase> >(
"Nullspace", Teuchos::null,
82 "Fine level nullspace used to construct the coarse level nullspace.");
83 validParamList->set<RCP<const FactoryBase> >(
"numDimensions", Teuchos::null,
84 "Number of spacial dimensions in the problem.");
85 validParamList->set<RCP<const FactoryBase> >(
"lCoarseNodesPerDim", Teuchos::null,
86 "Number of nodes per spatial dimension on the coarse grid.");
87 validParamList->set<RCP<const FactoryBase> >(
"indexManager", Teuchos::null,
88 "The index manager associated with the local mesh.");
89 validParamList->set<RCP<const FactoryBase> >(
"structuredInterpolationOrder", Teuchos::null,
90 "Interpolation order for constructing the prolongator.");
92 return validParamList;
95 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
98 const ParameterList& pL = GetParameterList();
100 Input(fineLevel,
"A");
101 Input(fineLevel,
"Nullspace");
102 Input(fineLevel,
"numDimensions");
103 Input(fineLevel,
"prolongatorGraph");
104 Input(fineLevel,
"lCoarseNodesPerDim");
105 Input(fineLevel,
"structuredInterpolationOrder");
107 if( pL.get<
bool>(
"interp: build coarse coordinates") ||
108 Get<int>(fineLevel,
"structuredInterpolationOrder") == 1) {
109 Input(fineLevel,
"Coordinates");
110 Input(fineLevel,
"indexManager");
115 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
118 return BuildP(fineLevel, coarseLevel);
121 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
127 RCP<Teuchos::FancyOStream> out;
128 if(
const char* dbg = std::getenv(
"MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
129 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
130 out->setShowAllFrontMatter(
false).setShowProcRank(
true);
132 out = Teuchos::getFancyOStream(rcp(
new Teuchos::oblackholestream()));
135 *out <<
"Starting GeometricInterpolationPFactory_kokkos::BuildP." << std::endl;
138 const ParameterList& pL = GetParameterList();
139 const bool buildCoarseCoordinates = pL.get<
bool>(
"interp: build coarse coordinates");
140 const int interpolationOrder = Get<int>(fineLevel,
"structuredInterpolationOrder");
141 const int numDimensions = Get<int>(fineLevel,
"numDimensions");
144 RCP<Matrix> A = Get<RCP<Matrix> >(fineLevel,
"A");
145 RCP<const CrsGraph> prolongatorGraph = Get<RCP<CrsGraph> >(fineLevel,
"prolongatorGraph");
146 RCP<realvaluedmultivector_type> fineCoordinates, coarseCoordinates;
151 if(buildCoarseCoordinates || (interpolationOrder == 1)) {
155 RCP<IndexManager_kokkos> geoData = Get<RCP<IndexManager_kokkos> >(fineLevel,
"indexManager");
156 fineCoordinates = Get< RCP<realvaluedmultivector_type> >(fineLevel,
"Coordinates");
159 RCP<const Map> coarseCoordsMap = MapFactory::Build(fineCoordinates->getMap()->lib(),
160 Teuchos::OrdinalTraits<GO>::invalid(),
161 geoData->getNumCoarseNodes(),
162 fineCoordinates->getMap()->getIndexBase(),
163 fineCoordinates->getMap()->getComm());
164 coarseCoordinates = Xpetra::MultiVectorFactory<real_type,LO,GO,Node>::
165 Build(coarseCoordsMap, fineCoordinates->getNumVectors());
170 fineCoordinates-> getDeviceLocalView(Xpetra::Access::ReadWrite),
171 coarseCoordinates->getDeviceLocalView(Xpetra::Access::OverwriteAll));
172 Kokkos::parallel_for(
"GeometricInterpolation: build coarse coordinates",
173 Kokkos::RangePolicy<execution_space>(0, geoData->getNumCoarseNodes()),
174 myCoarseCoordinatesBuilder);
176 Set(coarseLevel,
"Coordinates", coarseCoordinates);
179 *out <<
"Fine and coarse coordinates have been loaded from the fine level and set on the coarse level." << std::endl;
181 if(interpolationOrder == 0) {
184 BuildConstantP(P, prolongatorGraph, A);
185 }
else if(interpolationOrder == 1) {
188 RCP<realvaluedmultivector_type> ghostCoordinates
189 = Xpetra::MultiVectorFactory<real_type,LO,GO,NO>::Build(prolongatorGraph->getColMap(),
190 fineCoordinates->getNumVectors());
191 RCP<const Import> ghostImporter = ImportFactory::Build(coarseCoordinates->getMap(),
192 prolongatorGraph->getColMap());
193 ghostCoordinates->doImport(*coarseCoordinates, *ghostImporter, Xpetra::INSERT);
196 BuildLinearP(A, prolongatorGraph, fineCoordinates, ghostCoordinates, numDimensions, P);
199 *out <<
"The prolongator matrix has been built." << std::endl;
204 RCP<MultiVector> fineNullspace = Get< RCP<MultiVector> > (fineLevel,
"Nullspace");
205 RCP<MultiVector> coarseNullspace = MultiVectorFactory::Build(P->getDomainMap(),
206 fineNullspace->getNumVectors());
207 P->apply(*fineNullspace, *coarseNullspace, Teuchos::TRANS, Teuchos::ScalarTraits<SC>::one(),
208 Teuchos::ScalarTraits<SC>::zero());
209 Set(coarseLevel,
"Nullspace", coarseNullspace);
212 *out <<
"The coarse nullspace is constructed and set on the coarse level." << std::endl;
214 Array<LO> lNodesPerDir = Get<Array<LO> >(fineLevel,
"lCoarseNodesPerDim");
215 Set(coarseLevel,
"numDimensions", numDimensions);
216 Set(coarseLevel,
"lNodesPerDim", lNodesPerDir);
217 Set(coarseLevel,
"P", P);
219 *out <<
"GeometricInterpolationPFactory_kokkos::BuildP has completed." << std::endl;
223 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
225 BuildConstantP(RCP<Matrix>& P, RCP<const CrsGraph>& prolongatorGraph, RCP<Matrix>& A)
const {
228 RCP<Teuchos::FancyOStream> out;
229 if(
const char* dbg = std::getenv(
"MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
230 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
231 out->setShowAllFrontMatter(
false).setShowProcRank(
true);
233 out = Teuchos::getFancyOStream(rcp(
new Teuchos::oblackholestream()));
236 *out <<
"BuildConstantP" << std::endl;
238 std::vector<size_t> strideInfo(1);
239 strideInfo[0] = A->GetFixedBlockSize();
240 RCP<const StridedMap> stridedDomainMap =
241 StridedMapFactory::Build(prolongatorGraph->getDomainMap(), strideInfo);
243 *out <<
"Call prolongator constructor" << std::endl;
244 using helpers=Xpetra::Helpers<Scalar,LocalOrdinal,GlobalOrdinal,Node>;
245 if(helpers::isTpetraBlockCrs(A)) {
246#ifdef HAVE_MUELU_TPETRA
247 LO NSDim = A->GetStorageBlockSize();
251 RCP<const Map> BlockMap = prolongatorGraph->getDomainMap();
252 Teuchos::ArrayView<const GO> block_dofs = BlockMap->getLocalElementList();
253 Teuchos::Array<GO> point_dofs(block_dofs.size()*NSDim);
254 for(LO i=0, ct=0; i<block_dofs.size(); i++) {
255 for(LO j=0; j<NSDim; j++) {
256 point_dofs[ct] = block_dofs[i]*NSDim + j;
261 RCP<const Map> PointMap = MapFactory::Build(BlockMap->lib(),
262 BlockMap->getGlobalNumElements() *NSDim,
264 BlockMap->getIndexBase(),
265 BlockMap->getComm());
266 strideInfo[0] = A->GetFixedBlockSize();
267 RCP<const StridedMap> stridedPointMap = StridedMapFactory::Build(PointMap, strideInfo);
269 RCP<Xpetra::CrsMatrix<SC,LO,GO,NO> > P_xpetra = Xpetra::CrsMatrixFactory<SC,LO,GO,NO>::BuildBlock(prolongatorGraph, PointMap, A->getRangeMap(),NSDim);
270 RCP<Xpetra::TpetraBlockCrsMatrix<SC,LO,GO,NO> > P_tpetra = rcp_dynamic_cast<Xpetra::TpetraBlockCrsMatrix<SC,LO,GO,NO> >(P_xpetra);
271 if(P_tpetra.is_null())
throw std::runtime_error(
"BuildConstantP: Matrix factory did not return a Tpetra::BlockCrsMatrix");
272 RCP<CrsMatrixWrap> P_wrap = rcp(
new CrsMatrixWrap(P_xpetra));
274 const LO stride = strideInfo[0]*strideInfo[0];
275 const LO in_stride = strideInfo[0];
276 typename CrsMatrix::local_graph_type localGraph = prolongatorGraph->getLocalGraphDevice();
277 auto rowptr = localGraph.row_map;
278 auto indices = localGraph.entries;
279 auto values = P_tpetra->getTpetra_BlockCrsMatrix()->getValuesDeviceNonConst();
281 using ISC =
typename Tpetra::BlockCrsMatrix<SC,LO,GO,NO>::impl_scalar_type;
282 ISC one = Teuchos::ScalarTraits<ISC>::one();
284 const Kokkos::TeamPolicy<execution_space> policy(prolongatorGraph->getLocalNumRows(), 1);
286 Kokkos::parallel_for(
"MueLu:GeoInterpFact::BuildConstantP::fill", policy,
287 KOKKOS_LAMBDA(
const typename Kokkos::TeamPolicy<execution_space>::member_type &thread) {
288 auto row = thread.league_rank();
289 for(LO j = (LO)rowptr[row]; j<(LO) rowptr[row+1]; j++) {
290 LO block_offset = j*stride;
291 for(LO k=0; k<in_stride; k++)
292 values[block_offset + k*(in_stride+1) ] = one;
297 if (A->IsView(
"stridedMaps") ==
true) {
298 P->CreateView(
"stridedMaps", A->getRowMap(
"stridedMaps"), stridedPointMap);
301 P->CreateView(
"stridedMaps", P->getRangeMap(), PointMap);
305 throw std::runtime_error(
"GeometricInteroplationFactory::BuildConstantP(): BlockCrs requires Tpetra");
311 RCP<ParameterList> dummyList = rcp(
new ParameterList());
312 P = rcp(
new CrsMatrixWrap(prolongatorGraph, dummyList));
313 RCP<CrsMatrix> PCrs = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
314 PCrs->setAllToScalar(1.0);
315 PCrs->fillComplete();
318 if (A->IsView(
"stridedMaps") ==
true) {
319 P->CreateView(
"stridedMaps", A->getRowMap(
"stridedMaps"), stridedDomainMap);
321 P->CreateView(
"stridedMaps", P->getRangeMap(), stridedDomainMap);
327 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
329 BuildLinearP(RCP<Matrix>& A, RCP<const CrsGraph>& prolongatorGraph,
330 RCP<realvaluedmultivector_type>& fineCoordinates,
331 RCP<realvaluedmultivector_type>& ghostCoordinates,
332 const int numDimensions, RCP<Matrix>& P)
const {
335 RCP<Teuchos::FancyOStream> out;
336 if(
const char* dbg = std::getenv(
"MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
337 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
338 out->setShowAllFrontMatter(
false).setShowProcRank(
true);
340 out = Teuchos::getFancyOStream(rcp(
new Teuchos::oblackholestream()));
343 *out <<
"Entering BuildLinearP" << std::endl;
346 const LO numFineNodes = fineCoordinates->getLocalLength();
347 const LO numGhostNodes = ghostCoordinates->getLocalLength();
348 Array<ArrayRCP<const real_type> > fineCoords(3);
349 Array<ArrayRCP<const real_type> > ghostCoords(3);
350 const real_type realZero = Teuchos::as<real_type>(0.0);
351 ArrayRCP<real_type> fineZero(numFineNodes, realZero);
352 ArrayRCP<real_type> ghostZero(numGhostNodes, realZero);
353 for(
int dim = 0; dim < 3; ++dim) {
354 if(dim < numDimensions) {
355 fineCoords[dim] = fineCoordinates->getData(dim);
356 ghostCoords[dim] = ghostCoordinates->getData(dim);
358 fineCoords[dim] = fineZero;
359 ghostCoords[dim] = ghostZero;
363 *out <<
"Coordinates extracted from the multivectors!" << std::endl;
366 const int numInterpolationPoints = 1 << numDimensions;
367 const int dofsPerNode = A->GetFixedBlockSize();
369 std::vector<size_t> strideInfo(1);
370 strideInfo[0] = dofsPerNode;
371 RCP<const StridedMap> stridedDomainMap =
372 StridedMapFactory::Build(prolongatorGraph->getDomainMap(), strideInfo);
374 *out <<
"The maps of P have been computed" << std::endl;
376 RCP<ParameterList> dummyList = rcp(
new ParameterList());
377 P = rcp(
new CrsMatrixWrap(prolongatorGraph, dummyList));
378 RCP<CrsMatrix> PCrs = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
381 LO interpolationNodeIdx = 0, rowIdx = 0;
382 ArrayView<const LO> colIndices;
384 Array<Array<real_type> > coords(numInterpolationPoints + 1);
385 Array<real_type> stencil(numInterpolationPoints);
386 for(LO nodeIdx = 0; nodeIdx < numFineNodes; ++nodeIdx) {
387 if(PCrs->getNumEntriesInLocalRow(nodeIdx*dofsPerNode) == 1) {
390 for(LO dof = 0; dof < dofsPerNode; ++dof) {
391 rowIdx = nodeIdx*dofsPerNode + dof;
392 prolongatorGraph->getLocalRowView(rowIdx, colIndices);
393 PCrs->replaceLocalValues(rowIdx, colIndices, values());
399 for(
int dim = 0; dim < 3; ++dim) {
400 coords[0][dim] = fineCoords[dim][nodeIdx];
402 prolongatorGraph->getLocalRowView(nodeIdx*dofsPerNode, colIndices);
403 for(
int interpolationIdx=0; interpolationIdx < numInterpolationPoints; ++interpolationIdx) {
404 coords[interpolationIdx + 1].resize(3);
405 interpolationNodeIdx = colIndices[interpolationIdx] / dofsPerNode;
406 for(
int dim = 0; dim < 3; ++dim) {
407 coords[interpolationIdx + 1][dim] = ghostCoords[dim][interpolationNodeIdx];
410 ComputeLinearInterpolationStencil(numDimensions, numInterpolationPoints, coords, stencil);
411 values.resize(numInterpolationPoints);
412 for(LO valueIdx = 0; valueIdx < numInterpolationPoints; ++valueIdx) {
413 values[valueIdx] = Teuchos::as<SC>(stencil[valueIdx]);
417 for(LO dof = 0; dof < dofsPerNode; ++dof) {
418 rowIdx = nodeIdx*dofsPerNode + dof;
419 prolongatorGraph->getLocalRowView(rowIdx, colIndices);
420 PCrs->replaceLocalValues(rowIdx, colIndices, values());
425 *out <<
"The calculation of the interpolation stencils has completed." << std::endl;
427 PCrs->fillComplete();
429 *out <<
"All values in P have been set and expertStaticFillComplete has been performed." << std::endl;
432 if (A->IsView(
"stridedMaps") ==
true) {
433 P->CreateView(
"stridedMaps", A->getRowMap(
"stridedMaps"), stridedDomainMap);
435 P->CreateView(
"stridedMaps", P->getRangeMap(), stridedDomainMap);
441 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
444 const Array<Array<real_type> > coord,
445 Array<real_type>& stencil)
const {
462 Teuchos::SerialDenseMatrix<LO,real_type> Jacobian(numDimensions, numDimensions);
463 Teuchos::SerialDenseVector<LO,real_type> residual(numDimensions);
464 Teuchos::SerialDenseVector<LO,real_type> solutionDirection(numDimensions);
465 Teuchos::SerialDenseVector<LO,real_type> paramCoords(numDimensions);
466 Teuchos::SerialDenseSolver<LO,real_type> problem;
467 int iter = 0, max_iter = 5;
468 real_type functions[4][8], norm_ref = 1.0, norm2 = 1.0, tol = 1.0e-5;
469 paramCoords.size(numDimensions);
471 while( (iter < max_iter) && (norm2 > tol*norm_ref) ) {
474 solutionDirection.size(numDimensions);
475 residual.size(numDimensions);
479 GetInterpolationFunctions(numDimensions, paramCoords, functions);
480 for(LO i = 0; i < numDimensions; ++i) {
481 residual(i) = coord[0][i];
482 for(LO k = 0; k < numInterpolationPoints; ++k) {
483 residual(i) -= functions[0][k]*coord[k+1][i];
486 norm_ref += residual(i)*residual(i);
487 if(i == numDimensions - 1) {
488 norm_ref = std::sqrt(norm_ref);
492 for(LO j = 0; j < numDimensions; ++j) {
493 for(LO k = 0; k < numInterpolationPoints; ++k) {
494 Jacobian(i,j) += functions[j+1][k]*coord[k+1][i];
500 problem.setMatrix(Teuchos::rcp(&Jacobian,
false));
501 problem.setVectors(Teuchos::rcp(&solutionDirection,
false), Teuchos::rcp(&residual,
false));
502 if(problem.shouldEquilibrate()) {problem.factorWithEquilibration(
true);}
505 for(LO i = 0; i < numDimensions; ++i) {
506 paramCoords(i) = paramCoords(i) + solutionDirection(i);
510 GetInterpolationFunctions(numDimensions, paramCoords, functions);
511 for(LO i = 0; i < numDimensions; ++i) {
513 for(LO k = 0; k < numInterpolationPoints; ++k) {
514 tmp -= functions[0][k]*coord[k+1][i];
519 norm2 = std::sqrt(norm2);
523 for(LO i = 0; i < numInterpolationPoints; ++i) {
524 stencil[i] = functions[0][i];
529 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
532 const Teuchos::SerialDenseVector<LO, real_type> parametricCoordinates,
534 real_type xi = 0.0, eta = 0.0, zeta = 0.0, denominator = 0.0;
535 if(numDimensions == 1) {
536 xi = parametricCoordinates[0];
538 }
else if(numDimensions == 2) {
539 xi = parametricCoordinates[0];
540 eta = parametricCoordinates[1];
542 }
else if(numDimensions == 3) {
543 xi = parametricCoordinates[0];
544 eta = parametricCoordinates[1];
545 zeta = parametricCoordinates[2];
549 functions[0][0] = (1.0 - xi)*(1.0 - eta)*(1.0 - zeta) / denominator;
550 functions[0][1] = (1.0 + xi)*(1.0 - eta)*(1.0 - zeta) / denominator;
551 functions[0][2] = (1.0 - xi)*(1.0 + eta)*(1.0 - zeta) / denominator;
552 functions[0][3] = (1.0 + xi)*(1.0 + eta)*(1.0 - zeta) / denominator;
553 functions[0][4] = (1.0 - xi)*(1.0 - eta)*(1.0 + zeta) / denominator;
554 functions[0][5] = (1.0 + xi)*(1.0 - eta)*(1.0 + zeta) / denominator;
555 functions[0][6] = (1.0 - xi)*(1.0 + eta)*(1.0 + zeta) / denominator;
556 functions[0][7] = (1.0 + xi)*(1.0 + eta)*(1.0 + zeta) / denominator;
558 functions[1][0] = -(1.0 - eta)*(1.0 - zeta) / denominator;
559 functions[1][1] = (1.0 - eta)*(1.0 - zeta) / denominator;
560 functions[1][2] = -(1.0 + eta)*(1.0 - zeta) / denominator;
561 functions[1][3] = (1.0 + eta)*(1.0 - zeta) / denominator;
562 functions[1][4] = -(1.0 - eta)*(1.0 + zeta) / denominator;
563 functions[1][5] = (1.0 - eta)*(1.0 + zeta) / denominator;
564 functions[1][6] = -(1.0 + eta)*(1.0 + zeta) / denominator;
565 functions[1][7] = (1.0 + eta)*(1.0 + zeta) / denominator;
567 functions[2][0] = -(1.0 - xi)*(1.0 - zeta) / denominator;
568 functions[2][1] = -(1.0 + xi)*(1.0 - zeta) / denominator;
569 functions[2][2] = (1.0 - xi)*(1.0 - zeta) / denominator;
570 functions[2][3] = (1.0 + xi)*(1.0 - zeta) / denominator;
571 functions[2][4] = -(1.0 - xi)*(1.0 + zeta) / denominator;
572 functions[2][5] = -(1.0 + xi)*(1.0 + zeta) / denominator;
573 functions[2][6] = (1.0 - xi)*(1.0 + zeta) / denominator;
574 functions[2][7] = (1.0 + xi)*(1.0 + zeta) / denominator;
576 functions[3][0] = -(1.0 - xi)*(1.0 - eta) / denominator;
577 functions[3][1] = -(1.0 + xi)*(1.0 - eta) / denominator;
578 functions[3][2] = -(1.0 - xi)*(1.0 + eta) / denominator;
579 functions[3][3] = -(1.0 + xi)*(1.0 + eta) / denominator;
580 functions[3][4] = (1.0 - xi)*(1.0 - eta) / denominator;
581 functions[3][5] = (1.0 + xi)*(1.0 - eta) / denominator;
582 functions[3][6] = (1.0 - xi)*(1.0 + eta) / denominator;
583 functions[3][7] = (1.0 + xi)*(1.0 + eta) / denominator;
587 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
592 : geoData_(*geoData), fineCoordView_(fineCoordView), coarseCoordView_(coarseCoordView) {
596 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
597 KOKKOS_INLINE_FUNCTION
602 LO nodeCoarseTuple[3] = {0, 0, 0};
603 LO nodeFineTuple[3] = {0, 0, 0};
604 auto coarseningRate = geoData_.getCoarseningRates();
605 auto fineNodesPerDir = geoData_.getLocalFineNodesPerDir();
606 auto coarseNodesPerDir = geoData_.getCoarseNodesPerDir();
607 geoData_.getCoarseLID2CoarseTuple(coarseNodeIdx, nodeCoarseTuple);
608 for(
int dim = 0; dim < 3; ++dim) {
609 if(nodeCoarseTuple[dim] == coarseNodesPerDir(dim) - 1) {
610 nodeFineTuple[dim] = fineNodesPerDir(dim) - 1;
612 nodeFineTuple[dim] = nodeCoarseTuple[dim]*coarseningRate(dim);
616 fineNodeIdx = nodeFineTuple[2]*fineNodesPerDir(1)*fineNodesPerDir(0)
617 + nodeFineTuple[1]*fineNodesPerDir(0) + nodeFineTuple[0];
619 for(
int dim = 0; dim < fineCoordView_.extent_int(1); ++dim) {
620 coarseCoordView_(coarseNodeIdx, dim) = fineCoordView_(fineNodeIdx, dim);
#define SET_VALID_ENTRY(name)
Timer to be used in factories. Similar to Monitor but with additional timers.
typename Teuchos::ScalarTraits< SC >::coordinateType real_type
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void ComputeLinearInterpolationStencil(const int numDimensions, const int numInterpolationPoints, const Array< Array< real_type > > coord, Array< real_type > &stencil) const
void GetInterpolationFunctions(const LO numDimensions, const Teuchos::SerialDenseVector< LO, real_type > parametricCoordinates, real_type functions[4][8]) const
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
typename Kokkos::View< impl_scalar_type **, Kokkos::LayoutLeft, device_type > coord_view_type
void BuildConstantP(RCP< Matrix > &P, RCP< const CrsGraph > &prolongatorGraph, RCP< Matrix > &A) const
void BuildLinearP(RCP< Matrix > &A, RCP< const CrsGraph > &prolongatorGraph, RCP< realvaluedmultivector_type > &fineCoordinates, RCP< realvaluedmultivector_type > &ghostCoordinates, const int numDimensions, RCP< Matrix > &P) const
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Class that holds all level-specific information.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
Namespace for MueLu classes and methods.
KOKKOS_INLINE_FUNCTION void operator()(const LO nodeIdx) const
coarseCoordinatesBuilderFunctor(RCP< IndexManager_kokkos > geoData, coord_view_type fineCoordView, coord_view_type coarseCoordView)