MueLu Version of the Day
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"
50#include "Xpetra_CrsMatrixUtils.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
59namespace 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: build coarse coordinates");
67#undef SET_VALID_ENTRY
68
69 // general variables needed in GeometricInterpolationPFactory
70 validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null,
71 "Generating factory of the matrix A");
72 validParamList->set<RCP<const FactoryBase> >("Aggregates", Teuchos::null,
73 "Aggregates generated by StructuredAggregationFactory used to construct a piece-constant prolongator.");
74 validParamList->set<RCP<const FactoryBase> >("prolongatorGraph", Teuchos::null,
75 "Graph generated by StructuredAggregationFactory used to construct a piece-linear prolongator.");
76 validParamList->set<RCP<const FactoryBase> >("Coordinates", Teuchos::null,
77 "Fine level coordinates used to construct piece-wise linear prolongator and coarse level coordinates.");
78 validParamList->set<RCP<const FactoryBase> >("coarseCoordinatesFineMap", Teuchos::null,
79 "map of the coarse coordinates' GIDs as indexed on the fine mesh.");
80 validParamList->set<RCP<const FactoryBase> >("coarseCoordinatesMap", Teuchos::null,
81 "map of the coarse coordinates' GIDs as indexed on the coarse mesh.");
82 validParamList->set<RCP<const FactoryBase> >("Nullspace", Teuchos::null,
83 "Fine level nullspace used to construct the coarse level nullspace.");
84 validParamList->set<RCP<const FactoryBase> >("numDimensions", Teuchos::null,
85 "Number of spacial dimensions in the problem.");
86 validParamList->set<RCP<const FactoryBase> >("lCoarseNodesPerDim", Teuchos::null,
87 "Number of nodes per spatial dimension on the coarse grid.");
88 validParamList->set<RCP<const FactoryBase> >("structuredInterpolationOrder", Teuchos::null,
89 "Interpolation order for constructing the prolongator.");
90 validParamList->set<bool> ("keep coarse coords", false, "Flag to keep coordinates for special coarse grid solve");
91 validParamList->set<bool> ("interp: remove small entries", true, "Remove small interpolation coeficient from prolongator to reduce fill-in on coarse level");
92
93 return validParamList;
94 }
95
96 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
98 DeclareInput(Level& fineLevel, Level& /* coarseLevel */) const {
99 const ParameterList& pL = GetParameterList();
100
101 Input(fineLevel, "A");
102 Input(fineLevel, "Nullspace");
103 Input(fineLevel, "numDimensions");
104 Input(fineLevel, "prolongatorGraph");
105 Input(fineLevel, "lCoarseNodesPerDim");
106 Input(fineLevel, "structuredInterpolationOrder");
107
108 if( pL.get<bool>("interp: build coarse coordinates") ||
109 Get<int>(fineLevel, "structuredInterpolationOrder") == 1) {
110 Input(fineLevel, "Coordinates");
111 Input(fineLevel, "coarseCoordinatesFineMap");
112 Input(fineLevel, "coarseCoordinatesMap");
113 }
114
115 }
116
117 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
119 Build(Level& fineLevel, Level &coarseLevel) const {
120 return BuildP(fineLevel, coarseLevel);
121 }
122
123 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
125 BuildP(Level &fineLevel, Level &coarseLevel) const {
126 FactoryMonitor m(*this, "BuildP", coarseLevel);
127
128 // Set debug outputs based on environment variable
129 RCP<Teuchos::FancyOStream> out;
130 if(const char* dbg = std::getenv("MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
131 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
132 out->setShowAllFrontMatter(false).setShowProcRank(true);
133 } else {
134 out = Teuchos::getFancyOStream(rcp(new Teuchos::oblackholestream()));
135 }
136
137 *out << "Starting GeometricInterpolationPFactory::BuildP." << std::endl;
138
139 // Get inputs from the parameter list
140 const ParameterList& pL = GetParameterList();
141 const bool removeSmallEntries = pL.get<bool>("interp: remove small entries");
142 const bool buildCoarseCoordinates = pL.get<bool>("interp: build coarse coordinates");
143 const int interpolationOrder = Get<int>(fineLevel, "structuredInterpolationOrder");
144 const int numDimensions = Get<int>(fineLevel, "numDimensions");
145
146 // Declared main input/outputs to be retrieved and placed on the fine resp. coarse level
147 RCP<Matrix> A = Get<RCP<Matrix> >(fineLevel, "A");
148 RCP<const CrsGraph> prolongatorGraph = Get<RCP<CrsGraph> >(fineLevel, "prolongatorGraph");
149 RCP<realvaluedmultivector_type> fineCoordinates, coarseCoordinates;
150 RCP<Matrix> P;
151
152 // Check if we need to build coarse coordinates as they are used if we construct
153 // a linear interpolation prolongator
154 if(buildCoarseCoordinates || (interpolationOrder == 1)) {
155 SubFactoryMonitor sfm(*this, "BuildCoordinates", coarseLevel);
156 RCP<const Map> coarseCoordsFineMap = Get< RCP<const Map> >(fineLevel, "coarseCoordinatesFineMap");
157 RCP<const Map> coarseCoordsMap = Get< RCP<const Map> >(fineLevel, "coarseCoordinatesMap");
158 fineCoordinates = Get< RCP<realvaluedmultivector_type> >(fineLevel, "Coordinates");
159 coarseCoordinates = Xpetra::MultiVectorFactory<real_type,LO,GO,Node>::Build(coarseCoordsFineMap,
160 fineCoordinates->getNumVectors());
161 RCP<const Import> coordsImporter = ImportFactory::Build(fineCoordinates->getMap(),
162 coarseCoordsFineMap);
163 coarseCoordinates->doImport(*fineCoordinates, *coordsImporter, Xpetra::INSERT);
164 coarseCoordinates->replaceMap(coarseCoordsMap);
165
166 Set(coarseLevel, "Coordinates", coarseCoordinates);
167
168 if(pL.get<bool>("keep coarse coords")) {
169 coarseLevel.Set<RCP<realvaluedmultivector_type> >("Coordinates2", coarseCoordinates, NoFactory::get());
170 }
171 }
172
173 *out << "Fine and coarse coordinates have been loaded from the fine level and set on the coarse level." << std::endl;
174
175 if(interpolationOrder == 0) {
176 SubFactoryMonitor sfm(*this, "BuildConstantP", coarseLevel);
177 // Compute the prolongator using piece-wise constant interpolation
178 BuildConstantP(P, prolongatorGraph, A);
179 } else if(interpolationOrder == 1) {
180 // Compute the prolongator using piece-wise linear interpolation
181 // First get all the required coordinates to compute the local part of P
182 RCP<const Map> prolongatorColMap = prolongatorGraph->getColMap();
183
184 const size_t dofsPerNode = static_cast<size_t>(A->GetFixedBlockSize());
185 const size_t numColIndices = prolongatorColMap->getNodeNumElements();
186 TEUCHOS_TEST_FOR_EXCEPTION((numColIndices % dofsPerNode) != 0,
188 "Something went wrong, the number of columns in the prolongator is not a multiple of dofsPerNode!");
189 const size_t numGhostCoords = numColIndices / dofsPerNode;
190 const GO indexBase = prolongatorColMap->getIndexBase();
191 const GO coordIndexBase = fineCoordinates->getMap()->getIndexBase();
192
193 ArrayView<const GO> prolongatorColIndices = prolongatorColMap->getNodeElementList();
194 Array<GO> ghostCoordIndices(numGhostCoords);
195 for(size_t ghostCoordIdx = 0; ghostCoordIdx < numGhostCoords; ++ghostCoordIdx) {
196 ghostCoordIndices[ghostCoordIdx]
197 = (prolongatorColIndices[ghostCoordIdx*dofsPerNode] - indexBase) / dofsPerNode
198 + coordIndexBase;
199 }
200 RCP<Map> ghostCoordMap = MapFactory::Build(fineCoordinates->getMap()->lib(),
201 prolongatorColMap->getGlobalNumElements() / dofsPerNode,
202 ghostCoordIndices(),
203 coordIndexBase,
204 fineCoordinates->getMap()->getComm());
205
206 RCP<realvaluedmultivector_type> ghostCoordinates
207 = Xpetra::MultiVectorFactory<real_type,LO,GO,NO>::Build(ghostCoordMap,
208 fineCoordinates->getNumVectors());
209 RCP<const Import> ghostImporter = ImportFactory::Build(coarseCoordinates->getMap(),
210 ghostCoordMap);
211 ghostCoordinates->doImport(*coarseCoordinates, *ghostImporter, Xpetra::INSERT);
212
213 BuildLinearP(coarseLevel, A, prolongatorGraph, fineCoordinates,
214 ghostCoordinates, numDimensions, removeSmallEntries, P);
215 }
216
217 *out << "The prolongator matrix has been built." << std::endl;
218
219 {
220 SubFactoryMonitor sfm(*this, "BuildNullspace", coarseLevel);
221 // Build the coarse nullspace
222 RCP<MultiVector> fineNullspace = Get< RCP<MultiVector> > (fineLevel, "Nullspace");
223 RCP<MultiVector> coarseNullspace = MultiVectorFactory::Build(P->getDomainMap(),
224 fineNullspace->getNumVectors());
225 P->apply(*fineNullspace, *coarseNullspace, Teuchos::TRANS, Teuchos::ScalarTraits<SC>::one(),
226 Teuchos::ScalarTraits<SC>::zero());
227 Set(coarseLevel, "Nullspace", coarseNullspace);
228 }
229
230 *out << "The coarse nullspace is constructed and set on the coarse level." << std::endl;
231
232 Set(coarseLevel, "P", P);
233
234 *out << "GeometricInterpolationPFactory::BuildP has completed." << std::endl;
235
236 } // BuildP
237
238 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
240 BuildConstantP(RCP<Matrix>& P, RCP<const CrsGraph>& prolongatorGraph, RCP<Matrix>& A) const {
241
242 // Set debug outputs based on environment variable
243 RCP<Teuchos::FancyOStream> out;
244 if(const char* dbg = std::getenv("MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
245 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
246 out->setShowAllFrontMatter(false).setShowProcRank(true);
247 } else {
248 out = Teuchos::getFancyOStream(rcp(new Teuchos::oblackholestream()));
249 }
250
251 *out << "BuildConstantP" << std::endl;
252
253 std::vector<size_t> strideInfo(1);
254 strideInfo[0] = A->GetFixedBlockSize();
255 RCP<const StridedMap> stridedDomainMap =
256 StridedMapFactory::Build(prolongatorGraph->getDomainMap(), strideInfo);
257
258 *out << "Call prolongator constructor" << std::endl;
259
260 // Create the prolongator matrix and its associated objects
261 RCP<ParameterList> dummyList = rcp(new ParameterList());
262 P = rcp(new CrsMatrixWrap(prolongatorGraph, dummyList));
263 RCP<CrsMatrix> PCrs = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
264 PCrs->setAllToScalar(1.0);
265 PCrs->fillComplete();
266
267 // set StridingInformation of P
268 if (A->IsView("stridedMaps") == true) {
269 P->CreateView("stridedMaps", A->getRowMap("stridedMaps"), stridedDomainMap);
270 } else {
271 P->CreateView("stridedMaps", P->getRangeMap(), stridedDomainMap);
272 }
273
274 } // BuildConstantP
275
276 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
278 BuildLinearP(Level& coarseLevel,
279 RCP<Matrix>& A, RCP<const CrsGraph>& prolongatorGraph,
280 RCP<realvaluedmultivector_type>& fineCoordinates,
281 RCP<realvaluedmultivector_type>& ghostCoordinates,
282 const int numDimensions, const bool removeSmallEntries,
283 RCP<Matrix>& P) const {
284
285 // Set debug outputs based on environment variable
286 RCP<Teuchos::FancyOStream> out;
287 if(const char* dbg = std::getenv("MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEBUG")) {
288 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
289 out->setShowAllFrontMatter(false).setShowProcRank(true);
290 } else {
291 out = Teuchos::getFancyOStream(rcp(new Teuchos::oblackholestream()));
292 }
293
294 // Compute 2^numDimensions using bit logic to avoid round-off errors
295 const int numInterpolationPoints = 1 << numDimensions;
296 const int dofsPerNode = A->GetFixedBlockSize();
297
298 RCP<ParameterList> dummyList = rcp(new ParameterList());
299 P = rcp(new CrsMatrixWrap(prolongatorGraph, dummyList));
300 RCP<CrsMatrix> PCrs = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
301 PCrs->resumeFill(); // The Epetra matrix is considered filled at this point.
302
303 {
304 *out << "Entering BuildLinearP" << std::endl;
305 SubFactoryMonitor sfm(*this, "BuildLinearP", coarseLevel);
306
307 // Extract coordinates for interpolation stencil calculations
308 const LO numFineNodes = fineCoordinates->getLocalLength();
309 const LO numGhostNodes = ghostCoordinates->getLocalLength();
310 Array<ArrayRCP<const real_type> > fineCoords(3);
311 Array<ArrayRCP<const real_type> > ghostCoords(3);
312 const real_type realZero = Teuchos::as<real_type>(0.0);
313 ArrayRCP<real_type> fineZero(numFineNodes, realZero);
314 ArrayRCP<real_type> ghostZero(numGhostNodes, realZero);
315 for(int dim = 0; dim < 3; ++dim) {
316 if(dim < numDimensions) {
317 fineCoords[dim] = fineCoordinates->getData(dim);
318 ghostCoords[dim] = ghostCoordinates->getData(dim);
319 } else {
320 fineCoords[dim] = fineZero;
321 ghostCoords[dim] = ghostZero;
322 }
323 }
324
325 *out << "Coordinates extracted from the multivectors!" << std::endl;
326
327 { // Construct the linear interpolation prolongator
328 LO interpolationNodeIdx = 0, rowIdx = 0;
329 ArrayView<const LO> colIndices;
330 Array<SC> values;
331 Array<Array<real_type> > coords(numInterpolationPoints + 1);
332 Array<real_type> stencil(numInterpolationPoints);
333 for(LO nodeIdx = 0; nodeIdx < numFineNodes; ++nodeIdx) {
334 if(PCrs->getNumEntriesInLocalRow(nodeIdx*dofsPerNode) == 1) {
335 values.resize(1);
336 values[0] = 1.0;
337 for(LO dof = 0; dof < dofsPerNode; ++dof) {
338 rowIdx = nodeIdx*dofsPerNode + dof;
339 prolongatorGraph->getLocalRowView(rowIdx, colIndices);
340 PCrs->replaceLocalValues(rowIdx, colIndices, values());
341 }
342 } else {
343 // Extract the coordinates associated with the current node
344 // and the neighboring coarse nodes
345 coords[0].resize(3);
346 for(int dim = 0; dim < 3; ++dim) {
347 coords[0][dim] = fineCoords[dim][nodeIdx];
348 }
349 prolongatorGraph->getLocalRowView(nodeIdx*dofsPerNode, colIndices);
350 for(int interpolationIdx=0; interpolationIdx < numInterpolationPoints; ++interpolationIdx) {
351 coords[interpolationIdx + 1].resize(3);
352 interpolationNodeIdx = colIndices[interpolationIdx] / dofsPerNode;
353 for(int dim = 0; dim < 3; ++dim) {
354 coords[interpolationIdx + 1][dim] = ghostCoords[dim][interpolationNodeIdx];
355 }
356 }
357 RCP<Teuchos::TimeMonitor> tm = rcp(new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer("Compute Linear Interpolation")));
358 ComputeLinearInterpolationStencil(numDimensions, numInterpolationPoints, coords, stencil);
359 tm = Teuchos::null;
360 values.resize(numInterpolationPoints);
361 for(LO valueIdx = 0; valueIdx < numInterpolationPoints; ++valueIdx) {
362 values[valueIdx] = Teuchos::as<SC>(stencil[valueIdx]);
363 }
364
365 // Set values in all the rows corresponding to nodeIdx
366 for(LO dof = 0; dof < dofsPerNode; ++dof) {
367 rowIdx = nodeIdx*dofsPerNode + dof;
368 prolongatorGraph->getLocalRowView(rowIdx, colIndices);
369 PCrs->replaceLocalValues(rowIdx, colIndices, values());
370 }
371 } // Check for Coarse vs. Fine point
372 } // Loop over fine nodes
373 } // Construct the linear interpolation prolongator
374
375 *out << "The calculation of the interpolation stencils has completed." << std::endl;
376
377 PCrs->fillComplete();
378 }
379
380 *out << "All values in P have been set and fillComplete has been performed." << std::endl;
381
382 // Note lbv Jan 29 2019: this should be handle at aggregation level
383 // if the user really does not want potential d2 neighbors on coarse grid
384 // that way we would avoid a new graph construction...
385
386 // Check if we want to remove small entries from P
387 // to reduce stencil growth on next level.
388 if(removeSmallEntries) {
389 *out << "Entering remove small entries" << std::endl;
390 SubFactoryMonitor sfm(*this, "remove small entries", coarseLevel);
391
392 ArrayRCP<const size_t> rowptrOrig;
393 ArrayRCP<const LO> colindOrig;
394 ArrayRCP<const Scalar> valuesOrig;
395 PCrs->getAllValues(rowptrOrig, colindOrig, valuesOrig);
396
397 const size_t numRows = static_cast<size_t>(rowptrOrig.size() - 1);
398 ArrayRCP<size_t> rowPtr(numRows + 1);
399 ArrayRCP<size_t> nnzOnRows(numRows);
400 rowPtr[0] = 0;
401 size_t countRemovedEntries = 0;
402 for(size_t rowIdx = 0; rowIdx < numRows; ++rowIdx) {
403 for(size_t entryIdx = rowptrOrig[rowIdx]; entryIdx < rowptrOrig[rowIdx + 1]; ++entryIdx) {
404 if(Teuchos::ScalarTraits<Scalar>::magnitude(valuesOrig[entryIdx]) < 1e-6) {++countRemovedEntries;}
405 }
406 rowPtr[rowIdx + 1] = rowptrOrig[rowIdx + 1] - countRemovedEntries;
407 nnzOnRows[rowIdx] = rowPtr[rowIdx + 1] - rowPtr[rowIdx];
408 }
409 GetOStream(Statistics1) << "interp: number of small entries removed= " << countRemovedEntries << " / " << rowptrOrig[numRows] << std::endl;
410
411 size_t countKeptEntries = 0;
412 ArrayRCP<LO> colInd(rowPtr[numRows]);
413 ArrayRCP<SC> values(rowPtr[numRows]);
414 for(size_t entryIdx = 0; entryIdx < rowptrOrig[numRows]; ++entryIdx) {
415 if(Teuchos::ScalarTraits<Scalar>::magnitude(valuesOrig[entryIdx]) > 1e-6) {
416 colInd[countKeptEntries] = colindOrig[entryIdx];
417 values[countKeptEntries] = valuesOrig[entryIdx];
418 ++countKeptEntries;
419 }
420 }
421
422 P = rcp(new CrsMatrixWrap(prolongatorGraph->getRowMap(),
423 prolongatorGraph->getColMap(),
424 nnzOnRows));
425 RCP<CrsMatrix> PCrsSqueezed = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
426 PCrsSqueezed->resumeFill(); // The Epetra matrix is considered filled at this point.
427 PCrsSqueezed->setAllValues(rowPtr, colInd, values);
428 PCrsSqueezed->expertStaticFillComplete(prolongatorGraph->getDomainMap(),
429 prolongatorGraph->getRangeMap());
430 }
431
432 std::vector<size_t> strideInfo(1);
433 strideInfo[0] = dofsPerNode;
434 RCP<const StridedMap> stridedDomainMap =
435 StridedMapFactory::Build(prolongatorGraph->getDomainMap(), strideInfo);
436
437 *out << "The strided maps of P have been computed" << std::endl;
438
439 // set StridingInformation of P
440 if (A->IsView("stridedMaps") == true) {
441 P->CreateView("stridedMaps", A->getRowMap("stridedMaps"), stridedDomainMap);
442 } else {
443 P->CreateView("stridedMaps", P->getRangeMap(), stridedDomainMap);
444 }
445
446 } // BuildLinearP
447
448
449 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
451 ComputeLinearInterpolationStencil(const int numDimensions, const int numInterpolationPoints,
452 const Array<Array<real_type> > coord,
453 Array<real_type>& stencil) const {
454
455 // 7 8 Find xi, eta and zeta such that
456 // x---------x
457 // /| /| Rx = x_p - sum N_i(xi,eta,zeta)x_i = 0
458 // 5/ | 6/ | Ry = y_p - sum N_i(xi,eta,zeta)y_i = 0
459 // x---------x | Rz = z_p - sum N_i(xi,eta,zeta)z_i = 0
460 // | | *P | |
461 // | x------|--x We can do this with a Newton solver:
462 // | /3 | /4 We will start with initial guess (xi,eta,zeta) = (0,0,0)
463 // |/ |/ We compute the Jacobian and iterate until convergence...
464 // z y x---------x
465 // | / 1 2 Once we have (xi,eta,zeta), we can evaluate all N_i which
466 // |/ give us the weights for the interpolation stencil!
467 // o---x
468 //
469
470 Teuchos::SerialDenseMatrix<LO,real_type> Jacobian(numDimensions, numDimensions);
471 Teuchos::SerialDenseVector<LO,real_type> residual(numDimensions);
472 Teuchos::SerialDenseVector<LO,real_type> solutionDirection(numDimensions);
473 Teuchos::SerialDenseVector<LO,real_type> paramCoords(numDimensions);
474 Teuchos::SerialDenseSolver<LO,real_type> problem;
475 int iter = 0, max_iter = 5;
476 real_type functions[4][8], norm_ref = 1.0, norm2 = 1.0, tol = 1.0e-5;
477 paramCoords.size(numDimensions);
478
479 while( (iter < max_iter) && (norm2 > tol*norm_ref) ) {
480 ++iter;
481 norm2 = 0.0;
482 solutionDirection.size(numDimensions);
483 residual.size(numDimensions);
484 Jacobian = 0.0;
485
486 // Compute Jacobian and Residual
487 GetInterpolationFunctions(numDimensions, paramCoords, functions);
488 for(LO i = 0; i < numDimensions; ++i) {
489 residual(i) = coord[0][i]; // Add coordinates from point of interest
490 for(LO k = 0; k < numInterpolationPoints; ++k) {
491 residual(i) -= functions[0][k]*coord[k+1][i]; //Remove contribution from all coarse points
492 }
493 if(iter == 1) {
494 norm_ref += residual(i)*residual(i);
495 if(i == numDimensions - 1) {
496 norm_ref = std::sqrt(norm_ref);
497 }
498 }
499
500 for(LO j = 0; j < numDimensions; ++j) {
501 for(LO k = 0; k < numInterpolationPoints; ++k) {
502 Jacobian(i,j) += functions[j+1][k]*coord[k+1][i];
503 }
504 }
505 }
506
507 // Set Jacobian, Vectors and solve problem
508 problem.setMatrix(Teuchos::rcp(&Jacobian, false));
509 problem.setVectors(Teuchos::rcp(&solutionDirection, false), Teuchos::rcp(&residual, false));
510 if(problem.shouldEquilibrate()) {problem.factorWithEquilibration(true);}
511 problem.solve();
512
513 for(LO i = 0; i < numDimensions; ++i) {
514 paramCoords(i) = paramCoords(i) + solutionDirection(i);
515 }
516
517 // Recompute Residual norm
518 GetInterpolationFunctions(numDimensions, paramCoords, functions);
519 for(LO i = 0; i < numDimensions; ++i) {
520 real_type tmp = coord[0][i];
521 for(LO k = 0; k < numInterpolationPoints; ++k) {
522 tmp -= functions[0][k]*coord[k+1][i];
523 }
524 norm2 += tmp*tmp;
525 tmp = 0.0;
526 }
527 norm2 = std::sqrt(norm2);
528 }
529
530 // Load the interpolation values onto the stencil.
531 for(LO i = 0; i < numInterpolationPoints; ++i) {
532 stencil[i] = functions[0][i];
533 }
534
535 } // End ComputeLinearInterpolationStencil
536
537 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
539 GetInterpolationFunctions(const LO numDimensions,
540 const Teuchos::SerialDenseVector<LO, real_type> parametricCoordinates,
541 real_type functions[4][8]) const {
542 real_type xi = 0.0, eta = 0.0, zeta = 0.0, denominator = 0.0;
543 if(numDimensions == 1) {
544 xi = parametricCoordinates[0];
545 denominator = 2.0;
546 } else if(numDimensions == 2) {
547 xi = parametricCoordinates[0];
548 eta = parametricCoordinates[1];
549 denominator = 4.0;
550 } else if(numDimensions == 3) {
551 xi = parametricCoordinates[0];
552 eta = parametricCoordinates[1];
553 zeta = parametricCoordinates[2];
554 denominator = 8.0;
555 }
556
557 functions[0][0] = (1.0 - xi)*(1.0 - eta)*(1.0 - zeta) / denominator;
558 functions[0][1] = (1.0 + xi)*(1.0 - eta)*(1.0 - zeta) / denominator;
559 functions[0][2] = (1.0 - xi)*(1.0 + eta)*(1.0 - zeta) / denominator;
560 functions[0][3] = (1.0 + xi)*(1.0 + eta)*(1.0 - zeta) / denominator;
561 functions[0][4] = (1.0 - xi)*(1.0 - eta)*(1.0 + zeta) / denominator;
562 functions[0][5] = (1.0 + xi)*(1.0 - eta)*(1.0 + zeta) / denominator;
563 functions[0][6] = (1.0 - xi)*(1.0 + eta)*(1.0 + zeta) / denominator;
564 functions[0][7] = (1.0 + xi)*(1.0 + eta)*(1.0 + zeta) / denominator;
565
566 functions[1][0] = -(1.0 - eta)*(1.0 - zeta) / denominator;
567 functions[1][1] = (1.0 - eta)*(1.0 - zeta) / denominator;
568 functions[1][2] = -(1.0 + eta)*(1.0 - zeta) / denominator;
569 functions[1][3] = (1.0 + eta)*(1.0 - zeta) / denominator;
570 functions[1][4] = -(1.0 - eta)*(1.0 + zeta) / denominator;
571 functions[1][5] = (1.0 - eta)*(1.0 + zeta) / denominator;
572 functions[1][6] = -(1.0 + eta)*(1.0 + zeta) / denominator;
573 functions[1][7] = (1.0 + eta)*(1.0 + zeta) / denominator;
574
575 functions[2][0] = -(1.0 - xi)*(1.0 - zeta) / denominator;
576 functions[2][1] = -(1.0 + xi)*(1.0 - zeta) / denominator;
577 functions[2][2] = (1.0 - xi)*(1.0 - zeta) / denominator;
578 functions[2][3] = (1.0 + xi)*(1.0 - zeta) / denominator;
579 functions[2][4] = -(1.0 - xi)*(1.0 + zeta) / denominator;
580 functions[2][5] = -(1.0 + xi)*(1.0 + zeta) / denominator;
581 functions[2][6] = (1.0 - xi)*(1.0 + zeta) / denominator;
582 functions[2][7] = (1.0 + xi)*(1.0 + zeta) / denominator;
583
584 functions[3][0] = -(1.0 - xi)*(1.0 - eta) / denominator;
585 functions[3][1] = -(1.0 + xi)*(1.0 - eta) / denominator;
586 functions[3][2] = -(1.0 - xi)*(1.0 + eta) / denominator;
587 functions[3][3] = -(1.0 + xi)*(1.0 + eta) / denominator;
588 functions[3][4] = (1.0 - xi)*(1.0 - eta) / denominator;
589 functions[3][5] = (1.0 + xi)*(1.0 - eta) / denominator;
590 functions[3][6] = (1.0 - xi)*(1.0 + eta) / denominator;
591 functions[3][7] = (1.0 + xi)*(1.0 + eta) / denominator;
592
593 } // End GetInterpolationFunctions
594
595} // namespace MueLu
596
597#endif // MUELU_GEOMETRICINTERPOLATIONPFACTORY_DEF_HPP
#define SET_VALID_ENTRY(name)
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
typename Teuchos::ScalarTraits< SC >::coordinateType real_type
void BuildLinearP(Level &coarseLevel, RCP< Matrix > &A, RCP< const CrsGraph > &prolongatorGraph, RCP< realvaluedmultivector_type > &fineCoordinates, RCP< realvaluedmultivector_type > &ghostCoordinates, const int numDimensions, const bool keepD2, RCP< Matrix > &P) const
void BuildConstantP(RCP< Matrix > &P, RCP< const CrsGraph > &prolongatorGraph, RCP< Matrix > &A) const
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
void ComputeLinearInterpolationStencil(const int numDimensions, const int numInterpolationPoints, const Array< Array< real_type > > coord, Array< real_type > &stencil) const
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void GetInterpolationFunctions(const LO numDimensions, const Teuchos::SerialDenseVector< LO, real_type > parametricCoordinates, real_type functions[4][8]) const
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
void Set(const std::string &ename, const T &entry, const FactoryBase *factory=NoFactory::get())
static const NoFactory * get()
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
Namespace for MueLu classes and methods.
@ Statistics1
Print more statistics.