MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_GeometricInterpolationPFactory_kokkos_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_GEOMETRICINTERPOLATIONPFACTORY_KOKKOS_DEF_HPP
47#define MUELU_GEOMETRICINTERPOLATIONPFACTORY_KOKKOS_DEF_HPP
48
49#include "Xpetra_CrsGraph.hpp"
50#include "Xpetra_CrsMatrixUtils.hpp"
51
52#include "MueLu_MasterList.hpp"
53#include "MueLu_Monitor.hpp"
54#include "MueLu_IndexManager_kokkos.hpp"
55
56#include "Xpetra_TpetraCrsMatrix.hpp"
57
58
59// Including this one last ensure that the short names of the above headers are defined properly
61
62namespace MueLu {
63
64 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
66 RCP<ParameterList> validParamList = rcp(new ParameterList());
67
68#define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
69 SET_VALID_ENTRY("interp: build coarse coordinates");
70#undef SET_VALID_ENTRY
71
72 // general variables needed in GeometricInterpolationPFactory_kokkos
73 validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null,
74 "Generating factory of the matrix A");
75 validParamList->set<RCP<const FactoryBase> >("prolongatorGraph", Teuchos::null,
76 "Graph generated by StructuredAggregationFactory used to construct a piece-linear prolongator.");
77 validParamList->set<RCP<const FactoryBase> >("Coordinates", Teuchos::null,
78 "Fine level coordinates used to construct piece-wise linear prolongator and coarse level coordinates.");
79 validParamList->set<RCP<const FactoryBase> >("Nullspace", Teuchos::null,
80 "Fine level nullspace used to construct the coarse level nullspace.");
81 validParamList->set<RCP<const FactoryBase> >("numDimensions", Teuchos::null,
82 "Number of spacial dimensions in the problem.");
83 validParamList->set<RCP<const FactoryBase> >("lCoarseNodesPerDim", Teuchos::null,
84 "Number of nodes per spatial dimension on the coarse grid.");
85 validParamList->set<RCP<const FactoryBase> >("indexManager", Teuchos::null,
86 "The index manager associated with the local mesh.");
87 validParamList->set<RCP<const FactoryBase> >("structuredInterpolationOrder", Teuchos::null,
88 "Interpolation order for constructing the prolongator.");
89
90 return validParamList;
91 }
92
93 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
95 DeclareInput(Level& fineLevel, Level& coarseLevel) const {
96 const ParameterList& pL = GetParameterList();
97
98 Input(fineLevel, "A");
99 Input(fineLevel, "Nullspace");
100 Input(fineLevel, "numDimensions");
101 Input(fineLevel, "prolongatorGraph");
102 Input(fineLevel, "lCoarseNodesPerDim");
103 Input(fineLevel, "structuredInterpolationOrder");
104
105 if( pL.get<bool>("interp: build coarse coordinates") ||
106 Get<int>(fineLevel, "structuredInterpolationOrder") == 1) {
107 Input(fineLevel, "Coordinates");
108 Input(fineLevel, "indexManager");
109 }
110
111 }
112
113 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
115 Build(Level& fineLevel, Level &coarseLevel) const {
116 return BuildP(fineLevel, coarseLevel);
117 }
118
119 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
121 BuildP(Level &fineLevel, Level &coarseLevel) const {
122 FactoryMonitor m(*this, "BuildP", coarseLevel);
123
124 // Set debug outputs based on environment variable
125 RCP<Teuchos::FancyOStream> out;
126 if(const char* dbg = std::getenv("MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
127 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
128 out->setShowAllFrontMatter(false).setShowProcRank(true);
129 } else {
130 out = Teuchos::getFancyOStream(rcp(new Teuchos::oblackholestream()));
131 }
132
133 *out << "Starting GeometricInterpolationPFactory_kokkos::BuildP." << std::endl;
134
135 // Get inputs from the parameter list
136 const ParameterList& pL = GetParameterList();
137 const bool buildCoarseCoordinates = pL.get<bool>("interp: build coarse coordinates");
138 const int interpolationOrder = Get<int>(fineLevel, "structuredInterpolationOrder");
139 const int numDimensions = Get<int>(fineLevel, "numDimensions");
140
141 // Declared main input/outputs to be retrieved and placed on the fine resp. coarse level
142 RCP<Matrix> A = Get<RCP<Matrix> >(fineLevel, "A");
143 RCP<const CrsGraph> prolongatorGraph = Get<RCP<CrsGraph> >(fineLevel, "prolongatorGraph");
144 RCP<realvaluedmultivector_type> fineCoordinates, coarseCoordinates;
145 RCP<Matrix> P;
146
147 // Check if we need to build coarse coordinates as they are used if we construct
148 // a linear interpolation prolongator
149 if(buildCoarseCoordinates || (interpolationOrder == 1)) {
150 SubFactoryMonitor sfm(*this, "BuildCoordinates", coarseLevel);
151
152 // Extract data from fine level
153 RCP<IndexManager_kokkos> geoData = Get<RCP<IndexManager_kokkos> >(fineLevel, "indexManager");
154 fineCoordinates = Get< RCP<realvaluedmultivector_type> >(fineLevel, "Coordinates");
155
156 // Build coarse coordinates map/multivector
157 RCP<const Map> coarseCoordsMap = MapFactory::Build(fineCoordinates->getMap()->lib(),
158 Teuchos::OrdinalTraits<GO>::invalid(),
159 geoData->getNumCoarseNodes(),
160 fineCoordinates->getMap()->getIndexBase(),
161 fineCoordinates->getMap()->getComm());
162 coarseCoordinates = Xpetra::MultiVectorFactory<real_type,LO,GO,Node>::
163 Build(coarseCoordsMap, fineCoordinates->getNumVectors());
164
165 // Construct and launch functor to fill coarse coordinates values
166 // function should take a const view really
167 coarseCoordinatesBuilderFunctor myCoarseCoordinatesBuilder(geoData,
168 fineCoordinates-> getDeviceLocalView(Xpetra::Access::ReadWrite),
169 coarseCoordinates->getDeviceLocalView(Xpetra::Access::OverwriteAll));
170 Kokkos::parallel_for("GeometricInterpolation: build coarse coordinates",
171 Kokkos::RangePolicy<execution_space>(0, geoData->getNumCoarseNodes()),
172 myCoarseCoordinatesBuilder);
173
174 Set(coarseLevel, "Coordinates", coarseCoordinates);
175 }
176
177 *out << "Fine and coarse coordinates have been loaded from the fine level and set on the coarse level." << std::endl;
178
179 if(interpolationOrder == 0) {
180 SubFactoryMonitor sfm(*this, "BuildConstantP", coarseLevel);
181 // Compute the prolongator using piece-wise constant interpolation
182 BuildConstantP(P, prolongatorGraph, A);
183 } else if(interpolationOrder == 1) {
184 // Compute the prolongator using piece-wise linear interpolation
185 // First get all the required coordinates to compute the local part of P
186 RCP<realvaluedmultivector_type> ghostCoordinates
187 = Xpetra::MultiVectorFactory<real_type,LO,GO,NO>::Build(prolongatorGraph->getColMap(),
188 fineCoordinates->getNumVectors());
189 RCP<const Import> ghostImporter = ImportFactory::Build(coarseCoordinates->getMap(),
190 prolongatorGraph->getColMap());
191 ghostCoordinates->doImport(*coarseCoordinates, *ghostImporter, Xpetra::INSERT);
192
193 SubFactoryMonitor sfm(*this, "BuildLinearP", coarseLevel);
194 BuildLinearP(A, prolongatorGraph, fineCoordinates, ghostCoordinates, numDimensions, P);
195 }
196
197 *out << "The prolongator matrix has been built." << std::endl;
198
199 {
200 SubFactoryMonitor sfm(*this, "BuildNullspace", coarseLevel);
201 // Build the coarse nullspace
202 RCP<MultiVector> fineNullspace = Get< RCP<MultiVector> > (fineLevel, "Nullspace");
203 RCP<MultiVector> coarseNullspace = MultiVectorFactory::Build(P->getDomainMap(),
204 fineNullspace->getNumVectors());
205 P->apply(*fineNullspace, *coarseNullspace, Teuchos::TRANS, Teuchos::ScalarTraits<SC>::one(),
206 Teuchos::ScalarTraits<SC>::zero());
207 Set(coarseLevel, "Nullspace", coarseNullspace);
208 }
209
210 *out << "The coarse nullspace is constructed and set on the coarse level." << std::endl;
211
212 Array<LO> lNodesPerDir = Get<Array<LO> >(fineLevel, "lCoarseNodesPerDim");
213 Set(coarseLevel, "numDimensions", numDimensions);
214 Set(coarseLevel, "lNodesPerDim", lNodesPerDir);
215 Set(coarseLevel, "P", P);
216
217 *out << "GeometricInterpolationPFactory_kokkos::BuildP has completed." << std::endl;
218
219 } // BuildP
220
221 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
223 BuildConstantP(RCP<Matrix>& P, RCP<const CrsGraph>& prolongatorGraph, RCP<Matrix>& A) const {
224
225 // Set debug outputs based on environment variable
226 RCP<Teuchos::FancyOStream> out;
227 if(const char* dbg = std::getenv("MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
228 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
229 out->setShowAllFrontMatter(false).setShowProcRank(true);
230 } else {
231 out = Teuchos::getFancyOStream(rcp(new Teuchos::oblackholestream()));
232 }
233
234 *out << "BuildConstantP" << std::endl;
235
236 std::vector<size_t> strideInfo(1);
237 strideInfo[0] = A->GetFixedBlockSize();
238 RCP<const StridedMap> stridedDomainMap =
239 StridedMapFactory::Build(prolongatorGraph->getDomainMap(), strideInfo);
240
241 *out << "Call prolongator constructor" << std::endl;
242 using helpers=Xpetra::Helpers<Scalar,LocalOrdinal,GlobalOrdinal,Node>;
243 if(helpers::isTpetraBlockCrs(A)) {
244 LO NSDim = A->GetStorageBlockSize();
245
246 // Build the exploded Map
247 // FIXME: Should look at doing this on device
248 RCP<const Map> BlockMap = prolongatorGraph->getDomainMap();
249 Teuchos::ArrayView<const GO> block_dofs = BlockMap->getLocalElementList();
250 Teuchos::Array<GO> point_dofs(block_dofs.size()*NSDim);
251 for(LO i=0, ct=0; i<block_dofs.size(); i++) {
252 for(LO j=0; j<NSDim; j++) {
253 point_dofs[ct] = block_dofs[i]*NSDim + j;
254 ct++;
255 }
256 }
257
258 RCP<const Map> PointMap = MapFactory::Build(BlockMap->lib(),
259 BlockMap->getGlobalNumElements() *NSDim,
260 point_dofs(),
261 BlockMap->getIndexBase(),
262 BlockMap->getComm());
263 strideInfo[0] = A->GetFixedBlockSize();
264 RCP<const StridedMap> stridedPointMap = StridedMapFactory::Build(PointMap, strideInfo);
265
266 RCP<Xpetra::CrsMatrix<SC,LO,GO,NO> > P_xpetra = Xpetra::CrsMatrixFactory<SC,LO,GO,NO>::BuildBlock(prolongatorGraph, PointMap, A->getRangeMap(),NSDim);
267 RCP<Xpetra::TpetraBlockCrsMatrix<SC,LO,GO,NO> > P_tpetra = rcp_dynamic_cast<Xpetra::TpetraBlockCrsMatrix<SC,LO,GO,NO> >(P_xpetra);
268 if(P_tpetra.is_null()) throw std::runtime_error("BuildConstantP: Matrix factory did not return a Tpetra::BlockCrsMatrix");
269 RCP<CrsMatrixWrap> P_wrap = rcp(new CrsMatrixWrap(P_xpetra));
270
271 const LO stride = strideInfo[0]*strideInfo[0];
272 const LO in_stride = strideInfo[0];
273 typename CrsMatrix::local_graph_type localGraph = prolongatorGraph->getLocalGraphDevice();
274 auto rowptr = localGraph.row_map;
275 auto indices = localGraph.entries;
276 auto values = P_tpetra->getTpetra_BlockCrsMatrix()->getValuesDeviceNonConst();
277
278 using ISC = typename Tpetra::BlockCrsMatrix<SC,LO,GO,NO>::impl_scalar_type;
279 ISC one = Teuchos::ScalarTraits<ISC>::one();
280
281 const Kokkos::TeamPolicy<execution_space> policy(prolongatorGraph->getLocalNumRows(), 1);
282
283 Kokkos::parallel_for("MueLu:GeoInterpFact::BuildConstantP::fill", policy,
284 KOKKOS_LAMBDA(const typename Kokkos::TeamPolicy<execution_space>::member_type &thread) {
285 auto row = thread.league_rank();
286 for(LO j = (LO)rowptr[row]; j<(LO) rowptr[row+1]; j++) {
287 LO block_offset = j*stride;
288 for(LO k=0; k<in_stride; k++)
289 values[block_offset + k*(in_stride+1) ] = one;
290 }
291 });
292
293 P = P_wrap;
294 if (A->IsView("stridedMaps") == true) {
295 P->CreateView("stridedMaps", A->getRowMap("stridedMaps"), stridedPointMap);
296 }
297 else {
298 P->CreateView("stridedMaps", P->getRangeMap(), PointMap);
299 }
300 }
301 else {
302 // Create the prolongator matrix and its associated objects
303 RCP<ParameterList> dummyList = rcp(new ParameterList());
304 P = rcp(new CrsMatrixWrap(prolongatorGraph, dummyList));
305 RCP<CrsMatrix> PCrs = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
306 PCrs->setAllToScalar(1.0);
307 PCrs->fillComplete();
308
309 // set StridingInformation of P
310 if (A->IsView("stridedMaps") == true) {
311 P->CreateView("stridedMaps", A->getRowMap("stridedMaps"), stridedDomainMap);
312 } else {
313 P->CreateView("stridedMaps", P->getRangeMap(), stridedDomainMap);
314 }
315 }
316
317 } // BuildConstantP
318
319 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
321 BuildLinearP(RCP<Matrix>& A, RCP<const CrsGraph>& prolongatorGraph,
322 RCP<realvaluedmultivector_type>& fineCoordinates,
323 RCP<realvaluedmultivector_type>& ghostCoordinates,
324 const int numDimensions, RCP<Matrix>& P) const {
325
326 // Set debug outputs based on environment variable
327 RCP<Teuchos::FancyOStream> out;
328 if(const char* dbg = std::getenv("MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
329 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
330 out->setShowAllFrontMatter(false).setShowProcRank(true);
331 } else {
332 out = Teuchos::getFancyOStream(rcp(new Teuchos::oblackholestream()));
333 }
334
335 *out << "Entering BuildLinearP" << std::endl;
336
337 // Extract coordinates for interpolation stencil calculations
338 const LO numFineNodes = fineCoordinates->getLocalLength();
339 const LO numGhostNodes = ghostCoordinates->getLocalLength();
340 Array<ArrayRCP<const real_type> > fineCoords(3);
341 Array<ArrayRCP<const real_type> > ghostCoords(3);
342 const real_type realZero = Teuchos::as<real_type>(0.0);
343 ArrayRCP<real_type> fineZero(numFineNodes, realZero);
344 ArrayRCP<real_type> ghostZero(numGhostNodes, realZero);
345 for(int dim = 0; dim < 3; ++dim) {
346 if(dim < numDimensions) {
347 fineCoords[dim] = fineCoordinates->getData(dim);
348 ghostCoords[dim] = ghostCoordinates->getData(dim);
349 } else {
350 fineCoords[dim] = fineZero;
351 ghostCoords[dim] = ghostZero;
352 }
353 }
354
355 *out << "Coordinates extracted from the multivectors!" << std::endl;
356
357 // Compute 2^numDimensions using bit logic to avoid round-off errors
358 const int numInterpolationPoints = 1 << numDimensions;
359 const int dofsPerNode = A->GetFixedBlockSize();
360
361 std::vector<size_t> strideInfo(1);
362 strideInfo[0] = dofsPerNode;
363 RCP<const StridedMap> stridedDomainMap =
364 StridedMapFactory::Build(prolongatorGraph->getDomainMap(), strideInfo);
365
366 *out << "The maps of P have been computed" << std::endl;
367
368 RCP<ParameterList> dummyList = rcp(new ParameterList());
369 P = rcp(new CrsMatrixWrap(prolongatorGraph, dummyList));
370 RCP<CrsMatrix> PCrs = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
371 PCrs->resumeFill(); // The Epetra matrix is considered filled at this point.
372
373 LO interpolationNodeIdx = 0, rowIdx = 0;
374 ArrayView<const LO> colIndices;
375 Array<SC> values;
376 Array<Array<real_type> > coords(numInterpolationPoints + 1);
377 Array<real_type> stencil(numInterpolationPoints);
378 for(LO nodeIdx = 0; nodeIdx < numFineNodes; ++nodeIdx) {
379 if(PCrs->getNumEntriesInLocalRow(nodeIdx*dofsPerNode) == 1) {
380 values.resize(1);
381 values[0] = 1.0;
382 for(LO dof = 0; dof < dofsPerNode; ++dof) {
383 rowIdx = nodeIdx*dofsPerNode + dof;
384 prolongatorGraph->getLocalRowView(rowIdx, colIndices);
385 PCrs->replaceLocalValues(rowIdx, colIndices, values());
386 }
387 } else {
388 // Extract the coordinates associated with the current node
389 // and the neighboring coarse nodes
390 coords[0].resize(3);
391 for(int dim = 0; dim < 3; ++dim) {
392 coords[0][dim] = fineCoords[dim][nodeIdx];
393 }
394 prolongatorGraph->getLocalRowView(nodeIdx*dofsPerNode, colIndices);
395 for(int interpolationIdx=0; interpolationIdx < numInterpolationPoints; ++interpolationIdx) {
396 coords[interpolationIdx + 1].resize(3);
397 interpolationNodeIdx = colIndices[interpolationIdx] / dofsPerNode;
398 for(int dim = 0; dim < 3; ++dim) {
399 coords[interpolationIdx + 1][dim] = ghostCoords[dim][interpolationNodeIdx];
400 }
401 }
402 ComputeLinearInterpolationStencil(numDimensions, numInterpolationPoints, coords, stencil);
403 values.resize(numInterpolationPoints);
404 for(LO valueIdx = 0; valueIdx < numInterpolationPoints; ++valueIdx) {
405 values[valueIdx] = Teuchos::as<SC>(stencil[valueIdx]);
406 }
407
408 // Set values in all the rows corresponding to nodeIdx
409 for(LO dof = 0; dof < dofsPerNode; ++dof) {
410 rowIdx = nodeIdx*dofsPerNode + dof;
411 prolongatorGraph->getLocalRowView(rowIdx, colIndices);
412 PCrs->replaceLocalValues(rowIdx, colIndices, values());
413 }
414 }
415 }
416
417 *out << "The calculation of the interpolation stencils has completed." << std::endl;
418
419 PCrs->fillComplete();
420
421 *out << "All values in P have been set and expertStaticFillComplete has been performed." << std::endl;
422
423 // set StridingInformation of P
424 if (A->IsView("stridedMaps") == true) {
425 P->CreateView("stridedMaps", A->getRowMap("stridedMaps"), stridedDomainMap);
426 } else {
427 P->CreateView("stridedMaps", P->getRangeMap(), stridedDomainMap);
428 }
429
430 } // BuildLinearP
431
432
433 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
435 ComputeLinearInterpolationStencil(const int numDimensions, const int numInterpolationPoints,
436 const Array<Array<real_type> > coord,
437 Array<real_type>& stencil) const {
438
439 // 7 8 Find xi, eta and zeta such that
440 // x---------x
441 // /| /| Rx = x_p - sum N_i(xi,eta,zeta)x_i = 0
442 // 5/ | 6/ | Ry = y_p - sum N_i(xi,eta,zeta)y_i = 0
443 // x---------x | Rz = z_p - sum N_i(xi,eta,zeta)z_i = 0
444 // | | *P | |
445 // | x------|--x We can do this with a Newton solver:
446 // | /3 | /4 We will start with initial guess (xi,eta,zeta) = (0,0,0)
447 // |/ |/ We compute the Jacobian and iterate until convergence...
448 // z y x---------x
449 // | / 1 2 Once we have (xi,eta,zeta), we can evaluate all N_i which
450 // |/ give us the weights for the interpolation stencil!
451 // o---x
452 //
453
454 Teuchos::SerialDenseMatrix<LO,real_type> Jacobian(numDimensions, numDimensions);
455 Teuchos::SerialDenseVector<LO,real_type> residual(numDimensions);
456 Teuchos::SerialDenseVector<LO,real_type> solutionDirection(numDimensions);
457 Teuchos::SerialDenseVector<LO,real_type> paramCoords(numDimensions);
458 Teuchos::SerialDenseSolver<LO,real_type> problem;
459 int iter = 0, max_iter = 5;
460 real_type functions[4][8], norm_ref = 1.0, norm2 = 1.0, tol = 1.0e-5;
461 paramCoords.size(numDimensions);
462
463 while( (iter < max_iter) && (norm2 > tol*norm_ref) ) {
464 ++iter;
465 norm2 = 0.0;
466 solutionDirection.size(numDimensions);
467 residual.size(numDimensions);
468 Jacobian = 0.0;
469
470 // Compute Jacobian and Residual
471 GetInterpolationFunctions(numDimensions, paramCoords, functions);
472 for(LO i = 0; i < numDimensions; ++i) {
473 residual(i) = coord[0][i]; // Add coordinates from point of interest
474 for(LO k = 0; k < numInterpolationPoints; ++k) {
475 residual(i) -= functions[0][k]*coord[k+1][i]; //Remove contribution from all coarse points
476 }
477 if(iter == 1) {
478 norm_ref += residual(i)*residual(i);
479 if(i == numDimensions - 1) {
480 norm_ref = std::sqrt(norm_ref);
481 }
482 }
483
484 for(LO j = 0; j < numDimensions; ++j) {
485 for(LO k = 0; k < numInterpolationPoints; ++k) {
486 Jacobian(i,j) += functions[j+1][k]*coord[k+1][i];
487 }
488 }
489 }
490
491 // Set Jacobian, Vectors and solve problem
492 problem.setMatrix(Teuchos::rcp(&Jacobian, false));
493 problem.setVectors(Teuchos::rcp(&solutionDirection, false), Teuchos::rcp(&residual, false));
494 if(problem.shouldEquilibrate()) {problem.factorWithEquilibration(true);}
495 problem.solve();
496
497 for(LO i = 0; i < numDimensions; ++i) {
498 paramCoords(i) = paramCoords(i) + solutionDirection(i);
499 }
500
501 // Recompute Residual norm
502 GetInterpolationFunctions(numDimensions, paramCoords, functions);
503 for(LO i = 0; i < numDimensions; ++i) {
504 real_type tmp = coord[0][i];
505 for(LO k = 0; k < numInterpolationPoints; ++k) {
506 tmp -= functions[0][k]*coord[k+1][i];
507 }
508 norm2 += tmp*tmp;
509 tmp = 0.0;
510 }
511 norm2 = std::sqrt(norm2);
512 }
513
514 // Load the interpolation values onto the stencil.
515 for(LO i = 0; i < numInterpolationPoints; ++i) {
516 stencil[i] = functions[0][i];
517 }
518
519 } // End ComputeLinearInterpolationStencil
520
521 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
523 GetInterpolationFunctions(const LO numDimensions,
524 const Teuchos::SerialDenseVector<LO, real_type> parametricCoordinates,
525 real_type functions[4][8]) const {
526 real_type xi = 0.0, eta = 0.0, zeta = 0.0, denominator = 0.0;
527 if(numDimensions == 1) {
528 xi = parametricCoordinates[0];
529 denominator = 2.0;
530 } else if(numDimensions == 2) {
531 xi = parametricCoordinates[0];
532 eta = parametricCoordinates[1];
533 denominator = 4.0;
534 } else if(numDimensions == 3) {
535 xi = parametricCoordinates[0];
536 eta = parametricCoordinates[1];
537 zeta = parametricCoordinates[2];
538 denominator = 8.0;
539 }
540
541 functions[0][0] = (1.0 - xi)*(1.0 - eta)*(1.0 - zeta) / denominator;
542 functions[0][1] = (1.0 + xi)*(1.0 - eta)*(1.0 - zeta) / denominator;
543 functions[0][2] = (1.0 - xi)*(1.0 + eta)*(1.0 - zeta) / denominator;
544 functions[0][3] = (1.0 + xi)*(1.0 + eta)*(1.0 - zeta) / denominator;
545 functions[0][4] = (1.0 - xi)*(1.0 - eta)*(1.0 + zeta) / denominator;
546 functions[0][5] = (1.0 + xi)*(1.0 - eta)*(1.0 + zeta) / denominator;
547 functions[0][6] = (1.0 - xi)*(1.0 + eta)*(1.0 + zeta) / denominator;
548 functions[0][7] = (1.0 + xi)*(1.0 + eta)*(1.0 + zeta) / denominator;
549
550 functions[1][0] = -(1.0 - eta)*(1.0 - zeta) / denominator;
551 functions[1][1] = (1.0 - eta)*(1.0 - zeta) / denominator;
552 functions[1][2] = -(1.0 + eta)*(1.0 - zeta) / denominator;
553 functions[1][3] = (1.0 + eta)*(1.0 - zeta) / denominator;
554 functions[1][4] = -(1.0 - eta)*(1.0 + zeta) / denominator;
555 functions[1][5] = (1.0 - eta)*(1.0 + zeta) / denominator;
556 functions[1][6] = -(1.0 + eta)*(1.0 + zeta) / denominator;
557 functions[1][7] = (1.0 + eta)*(1.0 + zeta) / denominator;
558
559 functions[2][0] = -(1.0 - xi)*(1.0 - zeta) / denominator;
560 functions[2][1] = -(1.0 + xi)*(1.0 - zeta) / denominator;
561 functions[2][2] = (1.0 - xi)*(1.0 - zeta) / denominator;
562 functions[2][3] = (1.0 + xi)*(1.0 - zeta) / denominator;
563 functions[2][4] = -(1.0 - xi)*(1.0 + zeta) / denominator;
564 functions[2][5] = -(1.0 + xi)*(1.0 + zeta) / denominator;
565 functions[2][6] = (1.0 - xi)*(1.0 + zeta) / denominator;
566 functions[2][7] = (1.0 + xi)*(1.0 + zeta) / denominator;
567
568 functions[3][0] = -(1.0 - xi)*(1.0 - eta) / denominator;
569 functions[3][1] = -(1.0 + xi)*(1.0 - eta) / denominator;
570 functions[3][2] = -(1.0 - xi)*(1.0 + eta) / denominator;
571 functions[3][3] = -(1.0 + xi)*(1.0 + eta) / denominator;
572 functions[3][4] = (1.0 - xi)*(1.0 - eta) / denominator;
573 functions[3][5] = (1.0 + xi)*(1.0 - eta) / denominator;
574 functions[3][6] = (1.0 - xi)*(1.0 + eta) / denominator;
575 functions[3][7] = (1.0 + xi)*(1.0 + eta) / denominator;
576
577 } // End GetInterpolationFunctions
578
579 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
582 coord_view_type fineCoordView,
583 coord_view_type coarseCoordView)
584 : geoData_(*geoData), fineCoordView_(fineCoordView), coarseCoordView_(coarseCoordView) {
585
586 } // coarseCoordinatesBuilderFunctor()
587
588 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
589 KOKKOS_INLINE_FUNCTION
591 coarseCoordinatesBuilderFunctor::operator() (const LO coarseNodeIdx) const {
592
593 LO fineNodeIdx;
594 LO nodeCoarseTuple[3] = {0, 0, 0};
595 LO nodeFineTuple[3] = {0, 0, 0};
596 auto coarseningRate = geoData_.getCoarseningRates();
597 auto fineNodesPerDir = geoData_.getLocalFineNodesPerDir();
598 auto coarseNodesPerDir = geoData_.getCoarseNodesPerDir();
599 geoData_.getCoarseLID2CoarseTuple(coarseNodeIdx, nodeCoarseTuple);
600 for(int dim = 0; dim < 3; ++dim) {
601 if(nodeCoarseTuple[dim] == coarseNodesPerDir(dim) - 1) {
602 nodeFineTuple[dim] = fineNodesPerDir(dim) - 1;
603 } else {
604 nodeFineTuple[dim] = nodeCoarseTuple[dim]*coarseningRate(dim);
605 }
606 }
607
608 fineNodeIdx = nodeFineTuple[2]*fineNodesPerDir(1)*fineNodesPerDir(0)
609 + nodeFineTuple[1]*fineNodesPerDir(0) + nodeFineTuple[0];
610
611 for(int dim = 0; dim < fineCoordView_.extent_int(1); ++dim) {
612 coarseCoordView_(coarseNodeIdx, dim) = fineCoordView_(fineNodeIdx, dim);
613 }
614 }
615
616} // namespace MueLu
617
618#endif // MUELU_GEOMETRICINTERPOLATIONPFACTORY_KOKKOS_DEF_HPP
#define SET_VALID_ENTRY(name)
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
typename Kokkos::View< impl_scalar_type **, Kokkos::LayoutLeft, device_type > coord_view_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.
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.
virtual const Teuchos::ParameterList & GetParameterList() const
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
Namespace for MueLu classes and methods.
coarseCoordinatesBuilderFunctor(RCP< IndexManager_kokkos > geoData, coord_view_type fineCoordView, coord_view_type coarseCoordView)