MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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"
51 
52 #include "MueLu_MasterList.hpp"
53 #include "MueLu_Monitor.hpp"
54 #include "MueLu_IndexManager_kokkos.hpp"
55 
56 // Including this one last ensure that the short names of the above headers are defined properly
58 
59 namespace MueLu {
60 
61  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
63  RCP<ParameterList> validParamList = rcp(new ParameterList());
64 
65 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
66  SET_VALID_ENTRY("interp: interpolation order");
67  SET_VALID_ENTRY("interp: build coarse coordinates");
68 #undef SET_VALID_ENTRY
69 
70  // general variables needed in GeometricInterpolationPFactory_kokkos
71  validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null,
72  "Generating factory of the matrix A");
73  validParamList->set<RCP<const FactoryBase> >("prolongatorGraph", Teuchos::null,
74  "Graph generated by StructuredAggregationFactory used to construct a piece-linear prolongator.");
75  validParamList->set<RCP<const FactoryBase> >("Coordinates", Teuchos::null,
76  "Fine level coordinates used to construct piece-wise linear prolongator and coarse level coordinates.");
77  validParamList->set<RCP<const FactoryBase> >("Nullspace", Teuchos::null,
78  "Fine level nullspace used to construct the coarse level nullspace.");
79  validParamList->set<RCP<const FactoryBase> >("numDimensions", Teuchos::null,
80  "Number of spacial dimensions in the problem.");
81  validParamList->set<RCP<const FactoryBase> >("lCoarseNodesPerDim", Teuchos::null,
82  "Number of nodes per spatial dimension on the coarse grid.");
83  validParamList->set<RCP<const FactoryBase> >("indexManager", Teuchos::null,
84  "The index manager associated with the local mesh.");
85 
86  return validParamList;
87  }
88 
89  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
91  DeclareInput(Level& fineLevel, Level& coarseLevel) const {
92  const ParameterList& pL = GetParameterList();
93 
94  Input(fineLevel, "A");
95  Input(fineLevel, "Nullspace");
96  Input(fineLevel, "numDimensions");
97  Input(fineLevel, "prolongatorGraph");
98  Input(fineLevel, "lCoarseNodesPerDim");
99 
100  if( pL.get<bool>("interp: build coarse coordinates") ||
101  (pL.get<int>("interp: interpolation order") == 1) ) {
102  Input(fineLevel, "Coordinates");
103  Input(fineLevel, "indexManager");
104  }
105 
106  }
107 
108  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
110  Build(Level& fineLevel, Level &coarseLevel) const {
111  return BuildP(fineLevel, coarseLevel);
112  }
113 
114  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
116  BuildP(Level &fineLevel, Level &coarseLevel) const {
117  FactoryMonitor m(*this, "BuildP", coarseLevel);
118 
119  // Set debug outputs based on environment variable
121  if(const char* dbg = std::getenv("MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
122  out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
123  out->setShowAllFrontMatter(false).setShowProcRank(true);
124  } else {
125  out = Teuchos::getFancyOStream(rcp(new Teuchos::oblackholestream()));
126  }
127 
128  *out << "Starting GeometricInterpolationPFactory_kokkos::BuildP." << std::endl;
129 
130  // Get inputs from the parameter list
131  const ParameterList& pL = GetParameterList();
132  const bool buildCoarseCoordinates = pL.get<bool>("interp: build coarse coordinates");
133  const int interpolationOrder = pL.get<int> ("interp: interpolation order");
134  const int numDimensions = Get<int>(fineLevel, "numDimensions");
135 
136  // Declared main input/outputs to be retrieved and placed on the fine resp. coarse level
137  RCP<Matrix> A = Get<RCP<Matrix> >(fineLevel, "A");
138  RCP<const CrsGraph> prolongatorGraph = Get<RCP<CrsGraph> >(fineLevel, "prolongatorGraph");
139  RCP<realvaluedmultivector_type> fineCoordinates, coarseCoordinates;
140  RCP<Matrix> P;
141 
142  // Check if we need to build coarse coordinates as they are used if we construct
143  // a linear interpolation prolongator
144  if(buildCoarseCoordinates || (interpolationOrder == 1)) {
145  SubFactoryMonitor sfm(*this, "BuildCoordinates", coarseLevel);
146 
147  // Extract data from fine level
148  RCP<IndexManager_kokkos> geoData = Get<RCP<IndexManager_kokkos> >(fineLevel, "indexManager");
149  fineCoordinates = Get< RCP<realvaluedmultivector_type> >(fineLevel, "Coordinates");
150 
151  // Build coarse coordinates map/multivector
152  RCP<const Map> coarseCoordsMap = MapFactory::Build(fineCoordinates->getMap()->lib(),
154  geoData->getNumCoarseNodes(),
155  fineCoordinates->getMap()->getIndexBase(),
156  fineCoordinates->getMap()->getComm());
158  Build(coarseCoordsMap, fineCoordinates->getNumVectors());
159 
160  // Construct and launch functor to fill coarse coordinates values
161  coarseCoordinatesBuilderFunctor myCoarseCoordinatesBuilder(geoData,
162  fineCoordinates-> template getLocalView<device_type>(),
163  coarseCoordinates->template getLocalView<device_type>());
164  Kokkos::parallel_for("GeometricInterpolation: build coarse coordinates",
165  Kokkos::RangePolicy<execution_space>(0, geoData->getNumCoarseNodes()),
166  myCoarseCoordinatesBuilder);
167 
168  Set(coarseLevel, "Coordinates", coarseCoordinates);
169  }
170 
171  *out << "Fine and coarse coordinates have been loaded from the fine level and set on the coarse level." << std::endl;
172 
173  if(interpolationOrder == 0) {
174  SubFactoryMonitor sfm(*this, "BuildConstantP", coarseLevel);
175  // Compute the prolongator using piece-wise constant interpolation
176  BuildConstantP(P, prolongatorGraph, A);
177  } else if(interpolationOrder == 1) {
178  // Compute the prolongator using piece-wise linear interpolation
179  // First get all the required coordinates to compute the local part of P
180  RCP<realvaluedmultivector_type> ghostCoordinates
181  = Xpetra::MultiVectorFactory<real_type,LO,GO,NO>::Build(prolongatorGraph->getColMap(),
182  fineCoordinates->getNumVectors());
183  RCP<const Import> ghostImporter = ImportFactory::Build(coarseCoordinates->getMap(),
184  prolongatorGraph->getColMap());
185  ghostCoordinates->doImport(*coarseCoordinates, *ghostImporter, Xpetra::INSERT);
186 
187  SubFactoryMonitor sfm(*this, "BuildLinearP", coarseLevel);
188  BuildLinearP(A, prolongatorGraph, fineCoordinates, ghostCoordinates, numDimensions, P);
189  }
190 
191  *out << "The prolongator matrix has been built." << std::endl;
192 
193  {
194  SubFactoryMonitor sfm(*this, "BuildNullspace", coarseLevel);
195  // Build the coarse nullspace
196  RCP<MultiVector> fineNullspace = Get< RCP<MultiVector> > (fineLevel, "Nullspace");
197  RCP<MultiVector> coarseNullspace = MultiVectorFactory::Build(P->getDomainMap(),
198  fineNullspace->getNumVectors());
199  P->apply(*fineNullspace, *coarseNullspace, Teuchos::TRANS, Teuchos::ScalarTraits<SC>::one(),
201  Set(coarseLevel, "Nullspace", coarseNullspace);
202  }
203 
204  *out << "The coarse nullspace is constructed and set on the coarse level." << std::endl;
205 
206  Array<LO> lNodesPerDir = Get<Array<LO> >(fineLevel, "lCoarseNodesPerDim");
207  Set(coarseLevel, "numDimensions", numDimensions);
208  Set(coarseLevel, "lNodesPerDim", lNodesPerDir);
209  Set(coarseLevel, "P", P);
210 
211  *out << "GeometricInterpolationPFactory_kokkos::BuildP has completed." << std::endl;
212 
213  } // BuildP
214 
215  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
217  BuildConstantP(RCP<Matrix>& P, RCP<const CrsGraph>& prolongatorGraph, RCP<Matrix>& A) const {
218 
219  // Set debug outputs based on environment variable
221  if(const char* dbg = std::getenv("MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
222  out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
223  out->setShowAllFrontMatter(false).setShowProcRank(true);
224  } else {
225  out = Teuchos::getFancyOStream(rcp(new Teuchos::oblackholestream()));
226  }
227 
228  *out << "BuildConstantP" << std::endl;
229 
230  std::vector<size_t> strideInfo(1);
231  strideInfo[0] = A->GetFixedBlockSize();
232  RCP<const StridedMap> stridedDomainMap =
233  StridedMapFactory::Build(prolongatorGraph->getDomainMap(), strideInfo);
234 
235  *out << "Call prolongator constructor" << std::endl;
236 
237  // Create the prolongator matrix and its associated objects
238  RCP<ParameterList> dummyList = rcp(new ParameterList());
239  P = rcp(new CrsMatrixWrap(prolongatorGraph, dummyList));
240  RCP<CrsMatrix> PCrs = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
241  PCrs->setAllToScalar(1.0);
242  PCrs->fillComplete();
243 
244  // set StridingInformation of P
245  if (A->IsView("stridedMaps") == true) {
246  P->CreateView("stridedMaps", A->getRowMap("stridedMaps"), stridedDomainMap);
247  } else {
248  P->CreateView("stridedMaps", P->getRangeMap(), stridedDomainMap);
249  }
250 
251  } // BuildConstantP
252 
253  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
256  RCP<realvaluedmultivector_type>& fineCoordinates,
257  RCP<realvaluedmultivector_type>& ghostCoordinates,
258  const int numDimensions, RCP<Matrix>& P) const {
259 
260  // Set debug outputs based on environment variable
262  if(const char* dbg = std::getenv("MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
263  out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
264  out->setShowAllFrontMatter(false).setShowProcRank(true);
265  } else {
266  out = Teuchos::getFancyOStream(rcp(new Teuchos::oblackholestream()));
267  }
268 
269  *out << "Entering BuildLinearP" << std::endl;
270 
271  // Extract coordinates for interpolation stencil calculations
272  const LO numFineNodes = fineCoordinates->getLocalLength();
273  const LO numGhostNodes = ghostCoordinates->getLocalLength();
274  Array<ArrayRCP<const real_type> > fineCoords(3);
275  Array<ArrayRCP<const real_type> > ghostCoords(3);
276  const real_type realZero = Teuchos::as<real_type>(0.0);
277  ArrayRCP<real_type> fineZero(numFineNodes, realZero);
278  ArrayRCP<real_type> ghostZero(numGhostNodes, realZero);
279  for(int dim = 0; dim < 3; ++dim) {
280  if(dim < numDimensions) {
281  fineCoords[dim] = fineCoordinates->getData(dim);
282  ghostCoords[dim] = ghostCoordinates->getData(dim);
283  } else {
284  fineCoords[dim] = fineZero;
285  ghostCoords[dim] = ghostZero;
286  }
287  }
288 
289  *out << "Coordinates extracted from the multivectors!" << std::endl;
290 
291  // Compute 2^numDimensions using bit logic to avoid round-off errors
292  const int numInterpolationPoints = 1 << numDimensions;
293  const int dofsPerNode = A->GetFixedBlockSize();
294 
295  std::vector<size_t> strideInfo(1);
296  strideInfo[0] = dofsPerNode;
297  RCP<const StridedMap> stridedDomainMap =
298  StridedMapFactory::Build(prolongatorGraph->getDomainMap(), strideInfo);
299 
300  *out << "The maps of P have been computed" << std::endl;
301 
302  RCP<ParameterList> dummyList = rcp(new ParameterList());
303  P = rcp(new CrsMatrixWrap(prolongatorGraph, dummyList));
304  RCP<CrsMatrix> PCrs = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
305  PCrs->resumeFill(); // The Epetra matrix is considered filled at this point.
306 
307  LO interpolationNodeIdx = 0, rowIdx = 0;
308  ArrayView<const LO> colIndices;
309  Array<SC> values;
310  Array<Array<real_type> > coords(numInterpolationPoints + 1);
311  Array<real_type> stencil(numInterpolationPoints);
312  for(LO nodeIdx = 0; nodeIdx < numFineNodes; ++nodeIdx) {
313  if(PCrs->getNumEntriesInLocalRow(nodeIdx*dofsPerNode) == 1) {
314  values.resize(1);
315  values[0] = 1.0;
316  for(LO dof = 0; dof < dofsPerNode; ++dof) {
317  rowIdx = nodeIdx*dofsPerNode + dof;
318  prolongatorGraph->getLocalRowView(rowIdx, colIndices);
319  PCrs->replaceLocalValues(rowIdx, colIndices, values());
320  }
321  } else {
322  // Extract the coordinates associated with the current node
323  // and the neighboring coarse nodes
324  coords[0].resize(3);
325  for(int dim = 0; dim < 3; ++dim) {
326  coords[0][dim] = fineCoords[dim][nodeIdx];
327  }
328  prolongatorGraph->getLocalRowView(nodeIdx*dofsPerNode, colIndices);
329  for(int interpolationIdx=0; interpolationIdx < numInterpolationPoints; ++interpolationIdx) {
330  coords[interpolationIdx + 1].resize(3);
331  interpolationNodeIdx = colIndices[interpolationIdx] / dofsPerNode;
332  for(int dim = 0; dim < 3; ++dim) {
333  coords[interpolationIdx + 1][dim] = ghostCoords[dim][interpolationNodeIdx];
334  }
335  }
336  ComputeLinearInterpolationStencil(numDimensions, numInterpolationPoints, coords, stencil);
337  values.resize(numInterpolationPoints);
338  for(LO valueIdx = 0; valueIdx < numInterpolationPoints; ++valueIdx) {
339  values[valueIdx] = Teuchos::as<SC>(stencil[valueIdx]);
340  }
341 
342  // Set values in all the rows corresponding to nodeIdx
343  for(LO dof = 0; dof < dofsPerNode; ++dof) {
344  rowIdx = nodeIdx*dofsPerNode + dof;
345  prolongatorGraph->getLocalRowView(rowIdx, colIndices);
346  PCrs->replaceLocalValues(rowIdx, colIndices, values());
347  }
348  }
349  }
350 
351  *out << "The calculation of the interpolation stencils has completed." << std::endl;
352 
353  PCrs->fillComplete();
354 
355  *out << "All values in P have been set and expertStaticFillComplete has been performed." << std::endl;
356 
357  // set StridingInformation of P
358  if (A->IsView("stridedMaps") == true) {
359  P->CreateView("stridedMaps", A->getRowMap("stridedMaps"), stridedDomainMap);
360  } else {
361  P->CreateView("stridedMaps", P->getRangeMap(), stridedDomainMap);
362  }
363 
364  } // BuildLinearP
365 
366 
367  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
369  ComputeLinearInterpolationStencil(const int numDimensions, const int numInterpolationPoints,
370  const Array<Array<real_type> > coord,
371  Array<real_type>& stencil) const {
372 
373  // 7 8 Find xi, eta and zeta such that
374  // x---------x
375  // /| /| Rx = x_p - sum N_i(xi,eta,zeta)x_i = 0
376  // 5/ | 6/ | Ry = y_p - sum N_i(xi,eta,zeta)y_i = 0
377  // x---------x | Rz = z_p - sum N_i(xi,eta,zeta)z_i = 0
378  // | | *P | |
379  // | x------|--x We can do this with a Newton solver:
380  // | /3 | /4 We will start with initial guess (xi,eta,zeta) = (0,0,0)
381  // |/ |/ We compute the Jacobian and iterate until convergence...
382  // z y x---------x
383  // | / 1 2 Once we have (xi,eta,zeta), we can evaluate all N_i which
384  // |/ give us the weights for the interpolation stencil!
385  // o---x
386  //
387 
388  Teuchos::SerialDenseMatrix<LO,real_type> Jacobian(numDimensions, numDimensions);
389  Teuchos::SerialDenseVector<LO,real_type> residual(numDimensions);
390  Teuchos::SerialDenseVector<LO,real_type> solutionDirection(numDimensions);
391  Teuchos::SerialDenseVector<LO,real_type> paramCoords(numDimensions);
393  int iter = 0, max_iter = 5;
394  real_type functions[4][8], norm_ref = 1.0, norm2 = 1.0, tol = 1.0e-5;
395  paramCoords.size(numDimensions);
396 
397  while( (iter < max_iter) && (norm2 > tol*norm_ref) ) {
398  ++iter;
399  norm2 = 0.0;
400  solutionDirection.size(numDimensions);
401  residual.size(numDimensions);
402  Jacobian = 0.0;
403 
404  // Compute Jacobian and Residual
405  GetInterpolationFunctions(numDimensions, paramCoords, functions);
406  for(LO i = 0; i < numDimensions; ++i) {
407  residual(i) = coord[0][i]; // Add coordinates from point of interest
408  for(LO k = 0; k < numInterpolationPoints; ++k) {
409  residual(i) -= functions[0][k]*coord[k+1][i]; //Remove contribution from all coarse points
410  }
411  if(iter == 1) {
412  norm_ref += residual(i)*residual(i);
413  if(i == numDimensions - 1) {
414  norm_ref = std::sqrt(norm_ref);
415  }
416  }
417 
418  for(LO j = 0; j < numDimensions; ++j) {
419  for(LO k = 0; k < numInterpolationPoints; ++k) {
420  Jacobian(i,j) += functions[j+1][k]*coord[k+1][i];
421  }
422  }
423  }
424 
425  // Set Jacobian, Vectors and solve problem
426  problem.setMatrix(Teuchos::rcp(&Jacobian, false));
427  problem.setVectors(Teuchos::rcp(&solutionDirection, false), Teuchos::rcp(&residual, false));
428  if(problem.shouldEquilibrate()) {problem.factorWithEquilibration(true);}
429  problem.solve();
430 
431  for(LO i = 0; i < numDimensions; ++i) {
432  paramCoords(i) = paramCoords(i) + solutionDirection(i);
433  }
434 
435  // Recompute Residual norm
436  GetInterpolationFunctions(numDimensions, paramCoords, functions);
437  for(LO i = 0; i < numDimensions; ++i) {
438  real_type tmp = coord[0][i];
439  for(LO k = 0; k < numInterpolationPoints; ++k) {
440  tmp -= functions[0][k]*coord[k+1][i];
441  }
442  norm2 += tmp*tmp;
443  tmp = 0.0;
444  }
445  norm2 = std::sqrt(norm2);
446  }
447 
448  // Load the interpolation values onto the stencil.
449  for(LO i = 0; i < numInterpolationPoints; ++i) {
450  stencil[i] = functions[0][i];
451  }
452 
453  } // End ComputeLinearInterpolationStencil
454 
455  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
457  GetInterpolationFunctions(const LO numDimensions,
458  const Teuchos::SerialDenseVector<LO, real_type> parametricCoordinates,
459  real_type functions[4][8]) const {
460  real_type xi = 0.0, eta = 0.0, zeta = 0.0, denominator = 0.0;
461  if(numDimensions == 1) {
462  xi = parametricCoordinates[0];
463  denominator = 2.0;
464  } else if(numDimensions == 2) {
465  xi = parametricCoordinates[0];
466  eta = parametricCoordinates[1];
467  denominator = 4.0;
468  } else if(numDimensions == 3) {
469  xi = parametricCoordinates[0];
470  eta = parametricCoordinates[1];
471  zeta = parametricCoordinates[2];
472  denominator = 8.0;
473  }
474 
475  functions[0][0] = (1.0 - xi)*(1.0 - eta)*(1.0 - zeta) / denominator;
476  functions[0][1] = (1.0 + xi)*(1.0 - eta)*(1.0 - zeta) / denominator;
477  functions[0][2] = (1.0 - xi)*(1.0 + eta)*(1.0 - zeta) / denominator;
478  functions[0][3] = (1.0 + xi)*(1.0 + eta)*(1.0 - zeta) / denominator;
479  functions[0][4] = (1.0 - xi)*(1.0 - eta)*(1.0 + zeta) / denominator;
480  functions[0][5] = (1.0 + xi)*(1.0 - eta)*(1.0 + zeta) / denominator;
481  functions[0][6] = (1.0 - xi)*(1.0 + eta)*(1.0 + zeta) / denominator;
482  functions[0][7] = (1.0 + xi)*(1.0 + eta)*(1.0 + zeta) / denominator;
483 
484  functions[1][0] = -(1.0 - eta)*(1.0 - zeta) / denominator;
485  functions[1][1] = (1.0 - eta)*(1.0 - zeta) / denominator;
486  functions[1][2] = -(1.0 + eta)*(1.0 - zeta) / denominator;
487  functions[1][3] = (1.0 + eta)*(1.0 - zeta) / denominator;
488  functions[1][4] = -(1.0 - eta)*(1.0 + zeta) / denominator;
489  functions[1][5] = (1.0 - eta)*(1.0 + zeta) / denominator;
490  functions[1][6] = -(1.0 + eta)*(1.0 + zeta) / denominator;
491  functions[1][7] = (1.0 + eta)*(1.0 + zeta) / denominator;
492 
493  functions[2][0] = -(1.0 - xi)*(1.0 - zeta) / denominator;
494  functions[2][1] = -(1.0 + xi)*(1.0 - zeta) / denominator;
495  functions[2][2] = (1.0 - xi)*(1.0 - zeta) / denominator;
496  functions[2][3] = (1.0 + xi)*(1.0 - zeta) / denominator;
497  functions[2][4] = -(1.0 - xi)*(1.0 + zeta) / denominator;
498  functions[2][5] = -(1.0 + xi)*(1.0 + zeta) / denominator;
499  functions[2][6] = (1.0 - xi)*(1.0 + zeta) / denominator;
500  functions[2][7] = (1.0 + xi)*(1.0 + zeta) / denominator;
501 
502  functions[3][0] = -(1.0 - xi)*(1.0 - eta) / denominator;
503  functions[3][1] = -(1.0 + xi)*(1.0 - eta) / denominator;
504  functions[3][2] = -(1.0 - xi)*(1.0 + eta) / denominator;
505  functions[3][3] = -(1.0 + xi)*(1.0 + eta) / denominator;
506  functions[3][4] = (1.0 - xi)*(1.0 - eta) / denominator;
507  functions[3][5] = (1.0 + xi)*(1.0 - eta) / denominator;
508  functions[3][6] = (1.0 - xi)*(1.0 + eta) / denominator;
509  functions[3][7] = (1.0 + xi)*(1.0 + eta) / denominator;
510 
511  } // End GetInterpolationFunctions
512 
513  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
516  coord_view_type fineCoordView,
517  coord_view_type coarseCoordView)
518  : geoData_(*geoData), fineCoordView_(fineCoordView), coarseCoordView_(coarseCoordView) {
519 
520  } // coarseCoordinatesBuilderFunctor()
521 
522  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
523  KOKKOS_INLINE_FUNCTION
526 
527  LO fineNodeIdx;
528  LO nodeCoarseTuple[3] = {0, 0, 0};
529  LO nodeFineTuple[3] = {0, 0, 0};
530  auto coarseningRate = geoData_.getCoarseningRates();
531  auto fineNodesPerDir = geoData_.getLocalFineNodesPerDir();
532  auto coarseNodesPerDir = geoData_.getCoarseNodesPerDir();
533  geoData_.getCoarseLID2CoarseTuple(coarseNodeIdx, nodeCoarseTuple);
534  for(int dim = 0; dim < 3; ++dim) {
535  if(nodeCoarseTuple[dim] == coarseNodesPerDir(dim) - 1) {
536  nodeFineTuple[dim] = fineNodesPerDir(dim) - 1;
537  } else {
538  nodeFineTuple[dim] = nodeCoarseTuple[dim]*coarseningRate(dim);
539  }
540  }
541 
542  fineNodeIdx = nodeFineTuple[2]*fineNodesPerDir(1)*fineNodesPerDir(0)
543  + nodeFineTuple[1]*fineNodesPerDir(0) + nodeFineTuple[0];
544 
545  for(int dim = 0; dim < fineCoordView_.extent_int(1); ++dim) {
546  coarseCoordView_(coarseNodeIdx, dim) = fineCoordView_(fineNodeIdx, dim);
547  }
548  }
549 
550 } // namespace MueLu
551 
552 #endif // MUELU_GEOMETRICINTERPOLATIONPFACTORY_KOKKOS_DEF_HPP
void parallel_for(const ExecPolicy &policy, const FunctorType &functor, const std::string &str="", typename Impl::enable_if< Kokkos::Impl::is_execution_policy< ExecPolicy >::value >::type *=0)
basic_FancyOStream & setShowProcRank(const bool showProcRank)
coarseCoordinatesBuilderFunctor(RCP< IndexManager_kokkos > geoData, coord_view_type fineCoordView, coord_view_type coarseCoordView)
typename Kokkos::View< impl_scalar_type **, Kokkos::LayoutLeft, execution_space > coord_view_type
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Timer to be used in factories. Similar to Monitor but with additional timers.
LocalOrdinal LO
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.
typename Teuchos::ScalarTraits< SC >::coordinateType real_type
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
void ComputeLinearInterpolationStencil(const int numDimensions, const int numInterpolationPoints, const Array< Array< real_type > > coord, Array< real_type > &stencil) const
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
void factorWithEquilibration(bool flag)
void BuildConstantP(RCP< Matrix > &P, RCP< const CrsGraph > &prolongatorGraph, RCP< Matrix > &A) const
#define SET_VALID_ENTRY(name)
basic_FancyOStream & setShowAllFrontMatter(const bool showAllFrontMatter)
void resize(size_type new_size, const value_type &x=value_type())
TransListIter iter
int size(OrdinalType length_in)
int setVectors(const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &X, const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &B)
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
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.
static Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Build(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
int setMatrix(const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &A)