MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_GeometricInterpolationPFactory_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_DEF_HPP
47 #define MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEF_HPP
48 
49 #include "Xpetra_CrsGraph.hpp"
51 
52 #include "MueLu_MasterList.hpp"
53 #include "MueLu_Monitor.hpp"
54 #include "MueLu_Aggregates.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
71  validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null,
72  "Generating factory of the matrix A");
73  validParamList->set<RCP<const FactoryBase> >("Aggregates", Teuchos::null,
74  "Aggregates generated by StructuredAggregationFactory used to construct a piece-constant prolongator.");
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> >("coarseCoordinatesFineMap", Teuchos::null,
80  "map of the coarse coordinates' GIDs as indexed on the fine mesh.");
81  validParamList->set<RCP<const FactoryBase> >("coarseCoordinatesMap", Teuchos::null,
82  "map of the coarse coordinates' GIDs as indexed on the coarse mesh.");
83  validParamList->set<RCP<const FactoryBase> >("Nullspace", Teuchos::null,
84  "Fine level nullspace used to construct the coarse level nullspace.");
85  validParamList->set<RCP<const FactoryBase> >("numDimensions", Teuchos::null,
86  "Number of spacial dimensions in the problem.");
87  validParamList->set<RCP<const FactoryBase> >("lCoarseNodesPerDim", Teuchos::null,
88  "Number of nodes per spatial dimension on the coarse grid.");
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 
104  if( pL.get<bool>("interp: build coarse coordinates") ||
105  (pL.get<int>("interp: interpolation order") == 1) ) {
106  Input(fineLevel, "Coordinates");
107  Input(fineLevel, "coarseCoordinatesFineMap");
108  Input(fineLevel, "coarseCoordinatesMap");
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
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::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 = pL.get<int> ("interp: interpolation order");
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  RCP<const Map> coarseCoordsFineMap = Get< RCP<const Map> >(fineLevel, "coarseCoordinatesFineMap");
152  RCP<const Map> coarseCoordsMap = Get< RCP<const Map> >(fineLevel, "coarseCoordinatesMap");
153  fineCoordinates = Get< RCP<realvaluedmultivector_type> >(fineLevel, "Coordinates");
154  coarseCoordinates = Xpetra::MultiVectorFactory<real_type,LO,GO,Node>::Build(coarseCoordsFineMap,
155  fineCoordinates->getNumVectors());
156  RCP<const Import> coordsImporter = ImportFactory::Build(fineCoordinates->getMap(),
157  coarseCoordsFineMap);
158  coarseCoordinates->doImport(*fineCoordinates, *coordsImporter, Xpetra::INSERT);
159  coarseCoordinates->replaceMap(coarseCoordsMap);
160 
161  Set(coarseLevel, "Coordinates", coarseCoordinates);
162  }
163 
164  *out << "Fine and coarse coordinates have been loaded from the fine level and set on the coarse level." << std::endl;
165 
166  if(interpolationOrder == 0) {
167  SubFactoryMonitor sfm(*this, "BuildConstantP", coarseLevel);
168  // Compute the prolongator using piece-wise constant interpolation
169  BuildConstantP(P, prolongatorGraph, A);
170  } else if(interpolationOrder == 1) {
171  // Compute the prolongator using piece-wise linear interpolation
172  // First get all the required coordinates to compute the local part of P
173  RCP<realvaluedmultivector_type> ghostCoordinates
174  = Xpetra::MultiVectorFactory<real_type,LO,GO,NO>::Build(prolongatorGraph->getColMap(),
175  fineCoordinates->getNumVectors());
176  RCP<const Import> ghostImporter = ImportFactory::Build(coarseCoordinates->getMap(),
177  prolongatorGraph->getColMap());
178  ghostCoordinates->doImport(*coarseCoordinates, *ghostImporter, Xpetra::INSERT);
179 
180  SubFactoryMonitor sfm(*this, "BuildLinearP", coarseLevel);
181  BuildLinearP(A, prolongatorGraph, fineCoordinates, ghostCoordinates, numDimensions, P);
182  }
183 
184  *out << "The prolongator matrix has been built." << std::endl;
185 
186  {
187  SubFactoryMonitor sfm(*this, "BuildNullspace", coarseLevel);
188  // Build the coarse nullspace
189  RCP<MultiVector> fineNullspace = Get< RCP<MultiVector> > (fineLevel, "Nullspace");
190  RCP<MultiVector> coarseNullspace = MultiVectorFactory::Build(P->getDomainMap(),
191  fineNullspace->getNumVectors());
192  P->apply(*fineNullspace, *coarseNullspace, Teuchos::TRANS, Teuchos::ScalarTraits<SC>::one(),
194  Set(coarseLevel, "Nullspace", coarseNullspace);
195  }
196 
197  *out << "The coarse nullspace is constructed and set on the coarse level." << std::endl;
198 
199  Set(coarseLevel, "P", P);
200 
201  *out << "GeometricInterpolationPFactory::BuildP has completed." << std::endl;
202 
203  } // BuildP
204 
205  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
207  BuildConstantP(RCP<Matrix>& P, RCP<const CrsGraph>& prolongatorGraph, RCP<Matrix>& A) const {
208 
209  // Set debug outputs based on environment variable
211  if(const char* dbg = std::getenv("MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
212  out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
213  out->setShowAllFrontMatter(false).setShowProcRank(true);
214  } else {
215  out = Teuchos::getFancyOStream(rcp(new Teuchos::oblackholestream()));
216  }
217 
218  *out << "BuildConstantP" << std::endl;
219 
220  std::vector<size_t> strideInfo(1);
221  strideInfo[0] = A->GetFixedBlockSize();
222  RCP<const StridedMap> stridedDomainMap =
223  StridedMapFactory::Build(prolongatorGraph->getDomainMap(), strideInfo);
224 
225  *out << "Call prolongator constructor" << std::endl;
226 
227  // Create the prolongator matrix and its associated objects
228  RCP<ParameterList> dummyList = rcp(new ParameterList());
229  P = rcp(new CrsMatrixWrap(prolongatorGraph, dummyList));
230  RCP<CrsMatrix> PCrs = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
231  PCrs->setAllToScalar(1.0);
232  PCrs->fillComplete();
233 
234  // set StridingInformation of P
235  if (A->IsView("stridedMaps") == true) {
236  P->CreateView("stridedMaps", A->getRowMap("stridedMaps"), stridedDomainMap);
237  } else {
238  P->CreateView("stridedMaps", P->getRangeMap(), stridedDomainMap);
239  }
240 
241  } // BuildConstantP
242 
243  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
246  RCP<realvaluedmultivector_type>& fineCoordinates,
247  RCP<realvaluedmultivector_type>& ghostCoordinates,
248  const int numDimensions, RCP<Matrix>& P) const {
249 
250  // Set debug outputs based on environment variable
252  if(const char* dbg = std::getenv("MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
253  out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
254  out->setShowAllFrontMatter(false).setShowProcRank(true);
255  } else {
256  out = Teuchos::getFancyOStream(rcp(new Teuchos::oblackholestream()));
257  }
258 
259  *out << "Entering BuildLinearP" << std::endl;
260 
261  // Extract coordinates for interpolation stencil calculations
262  const LO numFineNodes = fineCoordinates->getLocalLength();
263  const LO numGhostNodes = ghostCoordinates->getLocalLength();
264  Array<ArrayRCP<const real_type> > fineCoords(3);
265  Array<ArrayRCP<const real_type> > ghostCoords(3);
266  const real_type realZero = Teuchos::as<real_type>(0.0);
267  ArrayRCP<real_type> fineZero(numFineNodes, realZero);
268  ArrayRCP<real_type> ghostZero(numGhostNodes, realZero);
269  for(int dim = 0; dim < 3; ++dim) {
270  if(dim < numDimensions) {
271  fineCoords[dim] = fineCoordinates->getData(dim);
272  ghostCoords[dim] = ghostCoordinates->getData(dim);
273  } else {
274  fineCoords[dim] = fineZero;
275  ghostCoords[dim] = ghostZero;
276  }
277  }
278 
279  *out << "Coordinates extracted from the multivectors!" << std::endl;
280 
281  // Compute 2^numDimensions using bit logic to avoid round-off errors
282  const int numInterpolationPoints = 1 << numDimensions;
283  const int dofsPerNode = A->GetFixedBlockSize();
284 
285  std::vector<size_t> strideInfo(1);
286  strideInfo[0] = dofsPerNode;
287  RCP<const StridedMap> stridedDomainMap =
288  StridedMapFactory::Build(prolongatorGraph->getDomainMap(), strideInfo);
289 
290  *out << "The maps of P have been computed" << std::endl;
291 
292  RCP<ParameterList> dummyList = rcp(new ParameterList());
293  P = rcp(new CrsMatrixWrap(prolongatorGraph, dummyList));
294  RCP<CrsMatrix> PCrs = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
295  PCrs->resumeFill(); // The Epetra matrix is considered filled at this point.
296 
297  LO interpolationNodeIdx = 0, rowIdx = 0;
298  ArrayView<const LO> colIndices;
299  Array<SC> values;
300  Array<Array<real_type> > coords(numInterpolationPoints + 1);
301  Array<real_type> stencil(numInterpolationPoints);
302  for(LO nodeIdx = 0; nodeIdx < numFineNodes; ++nodeIdx) {
303  if(PCrs->getNumEntriesInLocalRow(nodeIdx*dofsPerNode) == 1) {
304  values.resize(1);
305  values[0] = 1.0;
306  for(LO dof = 0; dof < dofsPerNode; ++dof) {
307  rowIdx = nodeIdx*dofsPerNode + dof;
308  prolongatorGraph->getLocalRowView(rowIdx, colIndices);
309  PCrs->replaceLocalValues(rowIdx, colIndices, values());
310  }
311  } else {
312  // Extract the coordinates associated with the current node
313  // and the neighboring coarse nodes
314  coords[0].resize(3);
315  for(int dim = 0; dim < 3; ++dim) {
316  coords[0][dim] = fineCoords[dim][nodeIdx];
317  }
318  prolongatorGraph->getLocalRowView(nodeIdx*dofsPerNode, colIndices);
319  for(int interpolationIdx=0; interpolationIdx < numInterpolationPoints; ++interpolationIdx) {
320  coords[interpolationIdx + 1].resize(3);
321  interpolationNodeIdx = colIndices[interpolationIdx] / dofsPerNode;
322  for(int dim = 0; dim < 3; ++dim) {
323  coords[interpolationIdx + 1][dim] = ghostCoords[dim][interpolationNodeIdx];
324  }
325  }
326  ComputeLinearInterpolationStencil(numDimensions, numInterpolationPoints, coords, stencil);
327  values.resize(numInterpolationPoints);
328  for(LO valueIdx = 0; valueIdx < numInterpolationPoints; ++valueIdx) {
329  values[valueIdx] = Teuchos::as<SC>(stencil[valueIdx]);
330  }
331 
332  // Set values in all the rows corresponding to nodeIdx
333  for(LO dof = 0; dof < dofsPerNode; ++dof) {
334  rowIdx = nodeIdx*dofsPerNode + dof;
335  prolongatorGraph->getLocalRowView(rowIdx, colIndices);
336  PCrs->replaceLocalValues(rowIdx, colIndices, values());
337  }
338  }
339  }
340 
341  *out << "The calculation of the interpolation stencils has completed." << std::endl;
342 
343  PCrs->fillComplete();
344 
345  *out << "All values in P have been set and expertStaticFillComplete has been performed." << std::endl;
346 
347  // set StridingInformation of P
348  if (A->IsView("stridedMaps") == true) {
349  P->CreateView("stridedMaps", A->getRowMap("stridedMaps"), stridedDomainMap);
350  } else {
351  P->CreateView("stridedMaps", P->getRangeMap(), stridedDomainMap);
352  }
353 
354  } // BuildLinearP
355 
356 
357  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
359  ComputeLinearInterpolationStencil(const int numDimensions, const int numInterpolationPoints,
360  const Array<Array<real_type> > coord,
361  Array<real_type>& stencil) const {
362 
363  // 7 8 Find xi, eta and zeta such that
364  // x---------x
365  // /| /| Rx = x_p - sum N_i(xi,eta,zeta)x_i = 0
366  // 5/ | 6/ | Ry = y_p - sum N_i(xi,eta,zeta)y_i = 0
367  // x---------x | Rz = z_p - sum N_i(xi,eta,zeta)z_i = 0
368  // | | *P | |
369  // | x------|--x We can do this with a Newton solver:
370  // | /3 | /4 We will start with initial guess (xi,eta,zeta) = (0,0,0)
371  // |/ |/ We compute the Jacobian and iterate until convergence...
372  // z y x---------x
373  // | / 1 2 Once we have (xi,eta,zeta), we can evaluate all N_i which
374  // |/ give us the weights for the interpolation stencil!
375  // o---x
376  //
377 
378  Teuchos::SerialDenseMatrix<LO,real_type> Jacobian(numDimensions, numDimensions);
379  Teuchos::SerialDenseVector<LO,real_type> residual(numDimensions);
380  Teuchos::SerialDenseVector<LO,real_type> solutionDirection(numDimensions);
381  Teuchos::SerialDenseVector<LO,real_type> paramCoords(numDimensions);
383  int iter = 0, max_iter = 5;
384  real_type functions[4][8], norm_ref = 1.0, norm2 = 1.0, tol = 1.0e-5;
385  paramCoords.size(numDimensions);
386 
387  while( (iter < max_iter) && (norm2 > tol*norm_ref) ) {
388  ++iter;
389  norm2 = 0.0;
390  solutionDirection.size(numDimensions);
391  residual.size(numDimensions);
392  Jacobian = 0.0;
393 
394  // Compute Jacobian and Residual
395  GetInterpolationFunctions(numDimensions, paramCoords, functions);
396  for(LO i = 0; i < numDimensions; ++i) {
397  residual(i) = coord[0][i]; // Add coordinates from point of interest
398  for(LO k = 0; k < numInterpolationPoints; ++k) {
399  residual(i) -= functions[0][k]*coord[k+1][i]; //Remove contribution from all coarse points
400  }
401  if(iter == 1) {
402  norm_ref += residual(i)*residual(i);
403  if(i == numDimensions - 1) {
404  norm_ref = std::sqrt(norm_ref);
405  }
406  }
407 
408  for(LO j = 0; j < numDimensions; ++j) {
409  for(LO k = 0; k < numInterpolationPoints; ++k) {
410  Jacobian(i,j) += functions[j+1][k]*coord[k+1][i];
411  }
412  }
413  }
414 
415  // Set Jacobian, Vectors and solve problem
416  problem.setMatrix(Teuchos::rcp(&Jacobian, false));
417  problem.setVectors(Teuchos::rcp(&solutionDirection, false), Teuchos::rcp(&residual, false));
418  if(problem.shouldEquilibrate()) {problem.factorWithEquilibration(true);}
419  problem.solve();
420 
421  for(LO i = 0; i < numDimensions; ++i) {
422  paramCoords(i) = paramCoords(i) + solutionDirection(i);
423  }
424 
425  // Recompute Residual norm
426  GetInterpolationFunctions(numDimensions, paramCoords, functions);
427  for(LO i = 0; i < numDimensions; ++i) {
428  real_type tmp = coord[0][i];
429  for(LO k = 0; k < numInterpolationPoints; ++k) {
430  tmp -= functions[0][k]*coord[k+1][i];
431  }
432  norm2 += tmp*tmp;
433  tmp = 0.0;
434  }
435  norm2 = std::sqrt(norm2);
436  }
437 
438  // Load the interpolation values onto the stencil.
439  for(LO i = 0; i < numInterpolationPoints; ++i) {
440  stencil[i] = functions[0][i];
441  }
442 
443  } // End ComputeLinearInterpolationStencil
444 
445  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
447  GetInterpolationFunctions(const LO numDimensions,
448  const Teuchos::SerialDenseVector<LO, real_type> parametricCoordinates,
449  real_type functions[4][8]) const {
450  real_type xi = 0.0, eta = 0.0, zeta = 0.0, denominator = 0.0;
451  if(numDimensions == 1) {
452  xi = parametricCoordinates[0];
453  denominator = 2.0;
454  } else if(numDimensions == 2) {
455  xi = parametricCoordinates[0];
456  eta = parametricCoordinates[1];
457  denominator = 4.0;
458  } else if(numDimensions == 3) {
459  xi = parametricCoordinates[0];
460  eta = parametricCoordinates[1];
461  zeta = parametricCoordinates[2];
462  denominator = 8.0;
463  }
464 
465  functions[0][0] = (1.0 - xi)*(1.0 - eta)*(1.0 - zeta) / denominator;
466  functions[0][1] = (1.0 + xi)*(1.0 - eta)*(1.0 - zeta) / denominator;
467  functions[0][2] = (1.0 - xi)*(1.0 + eta)*(1.0 - zeta) / denominator;
468  functions[0][3] = (1.0 + xi)*(1.0 + eta)*(1.0 - zeta) / denominator;
469  functions[0][4] = (1.0 - xi)*(1.0 - eta)*(1.0 + zeta) / denominator;
470  functions[0][5] = (1.0 + xi)*(1.0 - eta)*(1.0 + zeta) / denominator;
471  functions[0][6] = (1.0 - xi)*(1.0 + eta)*(1.0 + zeta) / denominator;
472  functions[0][7] = (1.0 + xi)*(1.0 + eta)*(1.0 + zeta) / denominator;
473 
474  functions[1][0] = -(1.0 - eta)*(1.0 - zeta) / denominator;
475  functions[1][1] = (1.0 - eta)*(1.0 - zeta) / denominator;
476  functions[1][2] = -(1.0 + eta)*(1.0 - zeta) / denominator;
477  functions[1][3] = (1.0 + eta)*(1.0 - zeta) / denominator;
478  functions[1][4] = -(1.0 - eta)*(1.0 + zeta) / denominator;
479  functions[1][5] = (1.0 - eta)*(1.0 + zeta) / denominator;
480  functions[1][6] = -(1.0 + eta)*(1.0 + zeta) / denominator;
481  functions[1][7] = (1.0 + eta)*(1.0 + zeta) / denominator;
482 
483  functions[2][0] = -(1.0 - xi)*(1.0 - zeta) / denominator;
484  functions[2][1] = -(1.0 + xi)*(1.0 - zeta) / denominator;
485  functions[2][2] = (1.0 - xi)*(1.0 - zeta) / denominator;
486  functions[2][3] = (1.0 + xi)*(1.0 - zeta) / denominator;
487  functions[2][4] = -(1.0 - xi)*(1.0 + zeta) / denominator;
488  functions[2][5] = -(1.0 + xi)*(1.0 + zeta) / denominator;
489  functions[2][6] = (1.0 - xi)*(1.0 + zeta) / denominator;
490  functions[2][7] = (1.0 + xi)*(1.0 + zeta) / denominator;
491 
492  functions[3][0] = -(1.0 - xi)*(1.0 - eta) / denominator;
493  functions[3][1] = -(1.0 + xi)*(1.0 - eta) / denominator;
494  functions[3][2] = -(1.0 - xi)*(1.0 + eta) / denominator;
495  functions[3][3] = -(1.0 + xi)*(1.0 + eta) / denominator;
496  functions[3][4] = (1.0 - xi)*(1.0 - eta) / denominator;
497  functions[3][5] = (1.0 + xi)*(1.0 - eta) / denominator;
498  functions[3][6] = (1.0 - xi)*(1.0 + eta) / denominator;
499  functions[3][7] = (1.0 + xi)*(1.0 + eta) / denominator;
500 
501  } // End GetInterpolationFunctions
502 
503 } // namespace MueLu
504 
505 #endif // MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEF_HPP
basic_FancyOStream & setShowProcRank(const bool showProcRank)
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
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
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)
#define SET_VALID_ENTRY(name)
void BuildLinearP(RCP< Matrix > &A, RCP< const CrsGraph > &prolongatorGraph, RCP< realvaluedmultivector_type > &fineCoordinates, RCP< realvaluedmultivector_type > &ghostCoordinates, const int numDimensions, RCP< Matrix > &P) const
basic_FancyOStream & setShowAllFrontMatter(const bool showAllFrontMatter)
void resize(size_type new_size, const value_type &x=value_type())
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
TransListIter iter
int size(OrdinalType length_in)
int setVectors(const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &X, const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &B)
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 BuildConstantP(RCP< Matrix > &P, RCP< const CrsGraph > &prolongatorGraph, RCP< Matrix > &A) const
typename Teuchos::ScalarTraits< SC >::coordinateType real_type
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
int setMatrix(const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &A)