version 3.10.0
Loading...
Searching...
No Matches
projector.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
18#ifndef DUMUX_DISCRETIZATION_PROJECTOR_HH
19#define DUMUX_DISCRETIZATION_PROJECTOR_HH
20
21#include <algorithm>
22#include <string>
23#include <utility>
24#include <type_traits>
25
26#include <dune/common/timer.hh>
27#include <dune/common/fmatrix.hh>
28#include <dune/common/exceptions.hh>
29#include <dune/common/promotiontraits.hh>
30#include <dune/common/parametertree.hh>
31
32#include <dune/geometry/quadraturerules.hh>
33#include <dune/istl/matrixindexset.hh>
34#include <dune/istl/bcrsmatrix.hh>
35#include <dune/istl/bvector.hh>
36
43
44namespace Dumux {
45
52template<class ScalarType>
54{
55 using MatrixBlockType = Dune::FieldMatrix<ScalarType, 1, 1>;
56
57public:
59 using Scalar = ScalarType;
62
64 struct Params
65 {
66 std::size_t maxIterations{100};
68 int verbosity{0};
69 };
70
72 Projector() = delete;
73
79 Projector(Matrix&& massMatrix, Matrix&& projectionMatrix)
80 : massMat_(std::make_shared<Matrix>(std::move(massMatrix)))
81 , projMat_(std::make_shared<Matrix>(std::move(projectionMatrix)))
82 , numDofsTarget_(massMat_->N())
83 {
84 if (massMat_->N() != projMat_->N())
85 DUNE_THROW(Dune::InvalidStateException, "Matrix row size mismatch: " << massMat_->N() << " vs " << projMat_->N());
86
87 massMatrixSolver_.setMatrix(massMat_);
88 }
89
100 Projector(Matrix&& massMatrix,
101 Matrix&& projectionMatrix,
102 std::vector<std::size_t>&& indexMap,
103 std::size_t numDofsTarget)
104 : massMat_(std::make_shared<Matrix>(std::move(massMatrix)))
105 , projMat_(std::make_shared<Matrix>(std::move(projectionMatrix)))
106 , indexMapTarget_(std::move(indexMap))
107 , numDofsTarget_(numDofsTarget)
108 {
109 if (indexMapTarget_.size() != massMat_->N())
110 DUNE_THROW(Dune::InvalidStateException, "Target index map size mismatch: " << indexMapTarget_.size() << " vs " << massMat_->N());
111
112 if (massMat_->N() != projMat_->N())
113 DUNE_THROW(Dune::InvalidStateException, "Matrix row size mismatch: " << massMat_->N() << " vs " << projMat_->N());
114
115 if (!indexMapTarget_.empty())
116 if (*std::max_element(indexMapTarget_.begin(), indexMapTarget_.end()) > numDofsTarget_)
117 DUNE_THROW(Dune::InvalidStateException, "Index map exceeds provided number of dofs in target domain!");
118
119 massMatrixSolver_.setMatrix(massMat_);
120 }
121
128 template< class BlockType, std::enable_if_t<std::is_convertible<BlockType, ScalarType>::value, int> = 0 >
129 Dune::BlockVector<BlockType> project(const Dune::BlockVector<BlockType>& u, const Params& params = Params{}) const
130 {
131 // be picky about size of u
132 if ( u.size() != projMat_->M())
133 DUNE_THROW(Dune::InvalidStateException, "Vector size mismatch");
134
135 Dune::BlockVector<BlockType> up(massMat_->N());
136
137 auto rhs = up;
138 projMat_->mv(u, rhs);
139
140 Dune::ParameterTree solverParams;
141 solverParams["maxit"] = std::to_string(params.maxIterations);
142 solverParams["reduction"] = std::to_string(params.residualReduction);
143 solverParams["verbose"] = std::to_string(params.verbosity);
144 auto solver = massMatrixSolver_; // copy the solver to modify the parameters
145 solver.setParams(solverParams);
146 solver.solve(up, rhs);
147
148 // if target space occupies a larger region, fill missing entries with zero
149 if (!indexMapTarget_.empty())
150 {
151 Dune::BlockVector<BlockType> result(numDofsTarget_);
152
153 result = 0.0;
154 for (std::size_t i = 0; i < indexMapTarget_.size(); ++i)
155 result[indexMapTarget_[i]] = up[i];
156
157 return result;
158 }
159
160 return up;
161 }
162
169 template< class BlockType, std::enable_if_t<!std::is_convertible<BlockType, ScalarType>::value, int> = 0 >
170 Dune::BlockVector<BlockType> project(const Dune::BlockVector<BlockType>& u, const Params& params = Params{}) const
171 {
172 Dune::BlockVector<BlockType> result(numDofsTarget_);
173
174 for (int pvIdx = 0; pvIdx < BlockType::size(); ++pvIdx)
175 {
176 Dune::BlockVector<Dune::FieldVector<Scalar, 1>> tmp(u.size());
177 std::transform(u.begin(), u.end(), tmp.begin(), [pvIdx] (const auto& v) { return v[pvIdx]; });
178
179 const auto p = project(tmp, params);
180 for (std::size_t i = 0; i < p.size(); ++i)
181 result[i][pvIdx] = p[i];
182 }
183
184 return result;
185 }
186
190 static Params defaultParams()
191 { return {}; }
192
193private:
194 std::shared_ptr<Matrix> massMat_;
195 std::shared_ptr<Matrix> projMat_;
196
199 > massMatrixSolver_;
200
201 std::vector<std::size_t> indexMapTarget_;
202 std::size_t numDofsTarget_;
203};
204
208template<class FEBasisDomain, class FEBasisTarget>
210{
211 using FiniteElementDomain = typename FEBasisDomain::LocalView::Tree::FiniteElement;
212 using FiniteElementTarget = typename FEBasisTarget::LocalView::Tree::FiniteElement;
213 using ScalarDomain = typename FiniteElementDomain::Traits::LocalBasisType::Traits::RangeFieldType;
214 using ScalarTarget = typename FiniteElementTarget::Traits::LocalBasisType::Traits::RangeFieldType;
215 using Scalar = typename Dune::PromotionTraits<ScalarDomain, ScalarTarget>::PromotedType;
216public:
218};
219
220
221// Projector construction details
222namespace Detail {
223
231template<class Matrix>
232void setupReducedMatrices(const Matrix& massMatrix, const Matrix& projMatrix, const std::vector<bool>& dofIsVoid,
233 Matrix& reducedM, Matrix& reducedP, std::vector<std::size_t>& expansionMap)
234{
235 const std::size_t numNonVoidDofs = std::count_if(dofIsVoid.begin(), dofIsVoid.end(), [] (bool v) { return !v; });
236
237 // reduce matrices to only dofs that take part and create index map
238 std::vector<std::size_t> reductionMap(massMatrix.N());
239 expansionMap.resize(numNonVoidDofs);
240
241 std::size_t idxInReducedSpace = 0;
242 for (std::size_t dofIdx = 0; dofIdx < dofIsVoid.size(); ++dofIdx)
243 if (!dofIsVoid[dofIdx])
244 {
245 reductionMap[dofIdx] = idxInReducedSpace;
246 expansionMap[idxInReducedSpace] = dofIdx;
247 idxInReducedSpace++;
248 }
249
250 // create reduced M/P matrix
251 Dune::MatrixIndexSet patternMReduced, patternPReduced;
252 patternMReduced.resize(numNonVoidDofs, numNonVoidDofs);
253 patternPReduced.resize(numNonVoidDofs, projMatrix.M());
254 for (auto rowIt = massMatrix.begin(); rowIt != massMatrix.end(); ++rowIt)
255 if (!dofIsVoid[rowIt.index()])
256 {
257 const auto reducedRowIdx = reductionMap[rowIt.index()];
258 for (auto colIt = (*rowIt).begin(); colIt != (*rowIt).end(); ++colIt)
259 if (!dofIsVoid[colIt.index()])
260 patternMReduced.add(reducedRowIdx, reductionMap[colIt.index()]);
261 }
262
263 for (auto rowIt = projMatrix.begin(); rowIt != projMatrix.end(); ++rowIt)
264 if (!dofIsVoid[rowIt.index()])
265 {
266 const auto reducedRowIdx = reductionMap[rowIt.index()];
267 for (auto colIt = (*rowIt).begin(); colIt != (*rowIt).end(); ++colIt)
268 patternPReduced.add(reducedRowIdx, colIt.index());
269 }
270
271 patternMReduced.exportIdx(reducedM);
272 patternPReduced.exportIdx(reducedP);
273
274 // fill matrix entries
275 for (auto rowIt = massMatrix.begin(); rowIt != massMatrix.end(); ++rowIt)
276 if (!dofIsVoid[rowIt.index()])
277 {
278 const auto reducedRowIdx = reductionMap[rowIt.index()];
279 for (auto colIt = (*rowIt).begin(); colIt != (*rowIt).end(); ++colIt)
280 if (!dofIsVoid[colIt.index()])
281 reducedM[reducedRowIdx][reductionMap[colIt.index()]] = *colIt;
282 }
283
284 for (auto rowIt = projMatrix.begin(); rowIt != projMatrix.end(); ++rowIt)
285 if (!dofIsVoid[rowIt.index()])
286 {
287 const auto reducedRowIdx = reductionMap[rowIt.index()];
288 for (auto colIt = (*rowIt).begin(); colIt != (*rowIt).end(); ++colIt)
289 reducedP[reducedRowIdx][colIt.index()] = *colIt;
290 }
291}
292
309template<bool doBidirectional, class FEBasisDomain, class FEBasisTarget, class GlueType>
310auto createProjectionMatrices(const FEBasisDomain& feBasisDomain,
311 const FEBasisTarget& feBasisTarget,
312 const GlueType& glue,
313 bool treatDiagonalZeroes = true)
314{
315 // we assume that target dim <= domain dimension
316 static constexpr int domainDim = FEBasisDomain::GridView::dimension;
317 static constexpr int targetDim = FEBasisTarget::GridView::dimension;
318 static_assert(targetDim <= domainDim, "This expects target dim < domain dim, please swap arguments");
319
320 using ForwardProjector = typename ProjectorTraits<FEBasisDomain, FEBasisTarget>::Projector;
321 using BackwardProjector = typename ProjectorTraits<FEBasisTarget, FEBasisDomain>::Projector;
322
323 using ForwardProjectionMatrix = typename ForwardProjector::Matrix;
324 using BackwardProjectionMatrix = typename BackwardProjector::Matrix;
325
326 auto domainLocalView = feBasisDomain.localView();
327 auto targetLocalView = feBasisTarget.localView();
328
329 // determine mass matrix patterns (standard FE scheme pattern)
330 Dune::MatrixIndexSet backwardPatternM, forwardPatternM;
331 forwardPatternM = getFEJacobianPattern(feBasisTarget);
332 if (doBidirectional) backwardPatternM = getFEJacobianPattern(feBasisDomain);
333
334 // determine projection matrix patterns
335 Dune::MatrixIndexSet backwardPatternP, forwardPatternP;
336 forwardPatternP.resize(feBasisTarget.size(), feBasisDomain.size());
337 if (doBidirectional) backwardPatternP.resize(feBasisDomain.size(), feBasisTarget.size());
338
339 using std::max;
340 unsigned int maxBasisOrder = 0;
341 for (const auto& is : intersections(glue))
342 {
343 // since target dim <= domain dim there is maximum one!
344 targetLocalView.bind( is.targetEntity(0) );
345 const auto& targetLocalBasis = targetLocalView.tree().finiteElement().localBasis();
346
347 for (unsigned int nIdx = 0; nIdx < is.numDomainNeighbors(); ++nIdx)
348 {
349 domainLocalView.bind( is.domainEntity(nIdx) );
350 const auto& domainLocalBasis = domainLocalView.tree().finiteElement().localBasis();
351
352 // keep track of maximum basis order (used in integration)
353 maxBasisOrder = max(maxBasisOrder, max(domainLocalBasis.order(), targetLocalBasis.order()));
354
355 for (unsigned int i = 0; i < domainLocalBasis.size(); ++i)
356 for (unsigned int j = 0; j < targetLocalBasis.size(); ++j)
357 {
358 forwardPatternP.add(targetLocalView.index(j), domainLocalView.index(i));
359 if (doBidirectional) backwardPatternP.add(domainLocalView.index(i), targetLocalView.index(j));
360 }
361 }
362 }
363
364 // assemble matrices
365 ForwardProjectionMatrix forwardM, forwardP;
366 forwardPatternM.exportIdx(forwardM); forwardM = 0.0;
367 forwardPatternP.exportIdx(forwardP); forwardP = 0.0;
368
369 BackwardProjectionMatrix backwardM, backwardP;
370 if (doBidirectional)
371 {
372 backwardPatternM.exportIdx(backwardM); backwardM = 0.0;
373 backwardPatternP.exportIdx(backwardP); backwardP = 0.0;
374 }
375
376 for (const auto& is : intersections(glue))
377 {
378 const auto& targetElement = is.targetEntity(0);
379 const auto& targetElementGeometry = targetElement.geometry();
380
381 targetLocalView.bind( targetElement );
382 const auto& targetLocalBasis = targetLocalView.tree().finiteElement().localBasis();
383
384 // perform integration over intersection geometry
385 using IsGeometry = typename std::decay_t<decltype(is.geometry())>;
386 using ctype = typename IsGeometry::ctype;
387
388 const auto& isGeometry = is.geometry();
389 const int intOrder = maxBasisOrder + 1;
390 const auto& quad = Dune::QuadratureRules<ctype, IsGeometry::mydimension>::rule(isGeometry.type(), intOrder);
391 for (auto&& qp : quad)
392 {
393 const auto weight = qp.weight();
394 const auto ie = isGeometry.integrationElement(qp.position());
395 const auto globalPos = isGeometry.global(qp.position());
396
397 std::vector< Dune::FieldVector<ctype, 1> > targetShapeVals;
398 targetLocalBasis.evaluateFunction(targetElementGeometry.local(globalPos), targetShapeVals);
399
400 // mass matrix entries target domain
401 for (unsigned int i = 0; i < targetLocalBasis.size(); ++i)
402 {
403 const auto dofIdxI = targetLocalView.index(i);
404 forwardM[dofIdxI][dofIdxI] += ie*weight*targetShapeVals[i]*targetShapeVals[i];
405
406 for (unsigned int j = i+1; j < targetLocalBasis.size(); ++j)
407 {
408 const auto dofIdxJ = targetLocalView.index(j);
409 const auto value = ie*weight*targetShapeVals[i]*targetShapeVals[j];
410 forwardM[dofIdxI][dofIdxJ] += value;
411 forwardM[dofIdxJ][dofIdxI] += value;
412 }
413 }
414
415 // If targetDim < domainDim, there can be several "neighbors" if
416 // targetElement is aligned with a facet of domainElement. In
417 // this case make sure the basis functions are not added
418 // multiple times! (division by numNeighbors)
419 const auto numNeighbors = is.numDomainNeighbors();
420 for (unsigned int nIdx = 0; nIdx < numNeighbors; ++nIdx)
421 {
422 const auto& domainElement = is.domainEntity(nIdx);
423 domainLocalView.bind( domainElement );
424 const auto& domainLocalBasis = domainLocalView.tree().finiteElement().localBasis();
425
426 std::vector< Dune::FieldVector<ctype, 1> > domainShapeVals;
427 domainLocalBasis.evaluateFunction(domainElement.geometry().local(globalPos), domainShapeVals);
428
429 // add entries in matrices
430 for (unsigned int i = 0; i < domainLocalBasis.size(); ++i)
431 {
432 const auto dofIdxDomain = domainLocalView.index(i);
433 const auto domainShapeVal = domainShapeVals[i];
434 if (doBidirectional)
435 {
436 backwardM[dofIdxDomain][dofIdxDomain] += ie*weight*domainShapeVal*domainShapeVal;
437
438 for (unsigned int j = i+1; j < domainLocalBasis.size(); ++j)
439 {
440 const auto dofIdxDomainJ = domainLocalView.index(j);
441 const auto value = ie*weight*domainShapeVal*domainShapeVals[j];
442 backwardM[dofIdxDomain][dofIdxDomainJ] += value;
443 backwardM[dofIdxDomainJ][dofIdxDomain] += value;
444 }
445 }
446
447 for (unsigned int j = 0; j < targetLocalBasis.size(); ++j)
448 {
449 const auto dofIdxTarget = targetLocalView.index(j);
450 const auto entry = ie*weight*domainShapeVal*targetShapeVals[j];
451
452 forwardP[dofIdxTarget][dofIdxDomain] += entry/numNeighbors;
453 if (doBidirectional)
454 backwardP[dofIdxDomain][dofIdxTarget] += entry;
455 }
456 }
457 }
458 }
459 }
460
461 // maybe treat zeroes on the diagonal
462 if (treatDiagonalZeroes)
463 {
464 for (std::size_t dofIdxTarget = 0; dofIdxTarget < forwardM.N(); ++dofIdxTarget)
465 if (forwardM[dofIdxTarget][dofIdxTarget] == 0.0)
466 forwardM[dofIdxTarget][dofIdxTarget] = 1.0;
467
468 if (doBidirectional)
469 {
470 for (std::size_t dofIdxDomain = 0; dofIdxDomain < backwardM.N(); ++dofIdxDomain)
471 if (backwardM[dofIdxDomain][dofIdxDomain] == 0.0)
472 backwardM[dofIdxDomain][dofIdxDomain] = 1.0;
473 }
474 }
475
476 return std::make_pair( std::make_pair(std::move(forwardM), std::move(forwardP)),
477 std::make_pair(std::move(backwardM), std::move(backwardP)) );
478}
479
485template<bool doBidirectional, class FEBasisDomain, class FEBasisTarget, class GlueType>
486auto makeProjectorPair(const FEBasisDomain& feBasisDomain,
487 const FEBasisTarget& feBasisTarget,
488 const GlueType& glue)
489{
490 using ForwardProjector = typename ProjectorTraits<FEBasisDomain, FEBasisTarget>::Projector;
491 using BackwardProjector = typename ProjectorTraits<FEBasisTarget, FEBasisDomain>::Projector;
492
493 using ForwardProjectionMatrix = typename ForwardProjector::Matrix;
494 using BackwardProjectionMatrix = typename BackwardProjector::Matrix;
495
496 auto projectionMatrices = createProjectionMatrices<doBidirectional>(feBasisDomain, feBasisTarget, glue, false);
497 auto& forwardMatrices = projectionMatrices.first;
498 auto& backwardMatrices = projectionMatrices.second;
499
500 auto& forwardM = forwardMatrices.first;
501 auto& forwardP = forwardMatrices.second;
502
503 auto& backwardM = backwardMatrices.first;
504 auto& backwardP = backwardMatrices.second;
505
506 // determine the dofs that do not take part in intersections
507 std::vector<bool> isVoidTarget(forwardM.N(), false);
508 for (std::size_t dofIdxTarget = 0; dofIdxTarget < forwardM.N(); ++dofIdxTarget)
509 if (forwardM[dofIdxTarget][dofIdxTarget] == 0.0)
510 isVoidTarget[dofIdxTarget] = true;
511
512 std::vector<bool> isVoidDomain;
513 if (doBidirectional)
514 {
515 isVoidDomain.resize(backwardM.N(), false);
516 for (std::size_t dofIdxDomain = 0; dofIdxDomain < backwardM.N(); ++dofIdxDomain)
517 if (backwardM[dofIdxDomain][dofIdxDomain] == 0.0)
518 isVoidDomain[dofIdxDomain] = true;
519 }
520
521 const bool hasVoidTarget = std::any_of(isVoidTarget.begin(), isVoidTarget.end(), [] (bool v) { return v; });
522 const bool hasVoidDomain = std::any_of(isVoidDomain.begin(), isVoidDomain.end(), [] (bool v) { return v; });
523 if (!hasVoidDomain && !hasVoidTarget)
524 {
525 return std::make_pair(ForwardProjector(std::move(forwardM), std::move(forwardP)),
526 BackwardProjector(std::move(backwardM), std::move(backwardP)));
527 }
528 else if (!hasVoidDomain && hasVoidTarget)
529 {
530 std::vector<std::size_t> expansionMapTarget;
531 ForwardProjectionMatrix forwardMReduced, forwardPReduced;
532 setupReducedMatrices(forwardM, forwardP, isVoidTarget,
533 forwardMReduced, forwardPReduced, expansionMapTarget);
534
535 return std::make_pair( ForwardProjector(std::move(forwardMReduced),
536 std::move(forwardPReduced),
537 std::move(expansionMapTarget),
538 forwardM.N()),
539 BackwardProjector(std::move(backwardM), std::move(backwardP)) );
540 }
541 else if (hasVoidDomain && !hasVoidTarget)
542 {
543 if (doBidirectional)
544 {
545 std::vector<std::size_t> expansionMapDomain;
546 BackwardProjectionMatrix backwardMReduced, backwardPReduced;
547 setupReducedMatrices(backwardM, backwardP, isVoidDomain,
548 backwardMReduced, backwardPReduced, expansionMapDomain);
549
550 return std::make_pair( ForwardProjector(std::move(forwardM), std::move(forwardP)),
551 BackwardProjector(std::move(backwardMReduced),
552 std::move(backwardPReduced),
553 std::move(expansionMapDomain),
554 backwardM.N()) );
555 }
556 else
557 return std::make_pair( ForwardProjector(std::move(forwardM), std::move(forwardP)),
558 BackwardProjector(std::move(backwardM), std::move(backwardP)) );
559 }
560 else
561 {
562 std::vector<std::size_t> expansionMapTarget;
563 ForwardProjectionMatrix forwardMReduced, forwardPReduced;
564 setupReducedMatrices(forwardM, forwardP, isVoidTarget,
565 forwardMReduced, forwardPReduced, expansionMapTarget);
566
567 if (doBidirectional)
568 {
569 std::vector<std::size_t> expansionMapDomain;
570 BackwardProjectionMatrix backwardMReduced, backwardPReduced;
571 setupReducedMatrices(backwardM, backwardP, isVoidDomain,
572 backwardMReduced, backwardPReduced, expansionMapDomain);
573
574 return std::make_pair( ForwardProjector(std::move(forwardMReduced),
575 std::move(forwardPReduced),
576 std::move(expansionMapTarget),
577 forwardM.N()),
578 BackwardProjector(std::move(backwardMReduced),
579 std::move(backwardPReduced),
580 std::move(expansionMapDomain),
581 backwardM.N()) );
582 }
583 else
584 return std::make_pair( ForwardProjector(std::move(forwardMReduced),
585 std::move(forwardPReduced),
586 std::move(expansionMapTarget),
587 forwardM.N()),
588 BackwardProjector(std::move(backwardM), std::move(backwardP)) );
589 }
590}
591
592} // end namespace Detail
593
594
607template< class FEBasisDomain, class FEBasisTarget, class GlueType >
608auto makeProjectorPair(const FEBasisDomain& feBasisDomain,
609 const FEBasisTarget& feBasisTarget,
610 GlueType glue)
611{
612 // we assume that target dim <= domain dimension
613 static constexpr int domainDim = FEBasisDomain::GridView::dimension;
614 static constexpr int targetDim = FEBasisTarget::GridView::dimension;
615 static_assert(targetDim <= domainDim, "makeProjectorPair() expects targetDim < domainDim, please swap arguments");
616
617 return Detail::makeProjectorPair<true>(feBasisDomain, feBasisTarget, glue);
618}
619
630template< class FEBasisDomain, class FEBasisTarget, class GlueType >
631auto makeProjector(const FEBasisDomain& feBasisDomain,
632 const FEBasisTarget& feBasisTarget,
633 GlueType glue)
634{
635 // we assume that target dim <= domain dimension
636 static constexpr int domainDim = FEBasisDomain::GridView::dimension;
637 static constexpr int targetDim = FEBasisTarget::GridView::dimension;
638 static_assert(targetDim <= domainDim, "makeProjectorPair() expects targetDim < domainDim, please swap arguments");
639
640 return Detail::makeProjectorPair<false>(feBasisDomain, feBasisTarget, glue).first;
641}
642
654template< class FEBasisDomain, class FEBasisTarget, class GlueType >
655auto makeProjectionMatricesPair(const FEBasisDomain& feBasisDomain,
656 const FEBasisTarget& feBasisTarget,
657 GlueType glue)
658{
659 // we assume that target dim <= domain dimension
660 static constexpr int domainDim = FEBasisDomain::GridView::dimension;
661 static constexpr int targetDim = FEBasisTarget::GridView::dimension;
662 static_assert(targetDim <= domainDim, "makeProjectionMatrixPair() expects targetDim < domainDim, please swap arguments");
663
664 return Detail::createProjectionMatrices<true>(feBasisDomain, feBasisTarget, glue);
665}
666
675template< class FEBasisDomain, class FEBasisTarget, class GlueType >
676auto makeProjectionMatrices(const FEBasisDomain& feBasisDomain,
677 const FEBasisTarget& feBasisTarget,
678 GlueType glue)
679{
680 // we assume that target dim <= domain dimension
681 static constexpr int domainDim = FEBasisDomain::GridView::dimension;
682 static constexpr int targetDim = FEBasisTarget::GridView::dimension;
683 static_assert(targetDim <= domainDim, "makeProjectionMatrixPair() expects targetDim < domainDim, please swap arguments");
684
685 return Detail::createProjectionMatrices<false>(feBasisDomain, feBasisTarget, glue).first;
686}
687
688} // end namespace Dumux
689
690#endif
static Params defaultParams()
Returns the default parameters.
Definition projector.hh:190
Projector(Matrix &&massMatrix, Matrix &&projectionMatrix, std::vector< std::size_t > &&indexMap, std::size_t numDofsTarget)
Constructor for projection into a target space that occupies a larger geometric region than the domai...
Definition projector.hh:100
ScalarType Scalar
Export the scalar type.
Definition projector.hh:59
Projector(Matrix &&massMatrix, Matrix &&projectionMatrix)
Constructor. Receives the mass and projection matrix that define the linear system describing the L2-...
Definition projector.hh:79
Dune::BCRSMatrix< MatrixBlockType > Matrix
Export the type of the projection matrices.
Definition projector.hh:61
Projector()=delete
delete default constructor
Dune::BlockVector< BlockType > project(const Dune::BlockVector< BlockType > &u, const Params &params=Params{}) const
Project a solution u into up.
Definition projector.hh:129
Traits class stating the type of projector between to bases.
Definition projector.hh:210
Dumux::Projector< Scalar > Projector
Definition projector.hh:217
Provides helper aliases and functionality to obtain the types and instances of Dune::Functions functi...
Dune::MatrixIndexSet getFEJacobianPattern(const FEBasis &feBasis)
Helper function to generate Jacobian pattern for finite element scheme.
Definition jacobianpattern.hh:106
auto makeProjector(const FEBasisDomain &feBasisDomain, const FEBasisTarget &feBasisTarget, GlueType glue)
Creates a forward projector from the space feBasisDomain to the space with basis feBasisTarget.
Definition projector.hh:631
auto makeProjectorPair(const FEBasisDomain &feBasisDomain, const FEBasisTarget &feBasisTarget, GlueType glue)
Creates a pair of projectors between the space with basis feBasisDomain to the space with basis feBas...
Definition projector.hh:608
Detail::IstlIterativeLinearSolver< LSTraits, LATraits, Dune::CGSolver< typename LATraits::Vector >, Detail::IstlSolvers::IstlDefaultBlockLevelPreconditionerFactory< Dune::SeqSSOR > > SSORCGIstlSolver
An SSOR-preconditioned CG solver using dune-istl.
Definition istlsolvers.hh:688
Linear solvers from dune-istl.
Helper function to generate Jacobian pattern for different discretization methods.
Define traits for linear algebra backends.
Define traits for linear solvers.
Distance implementation details.
Definition cvfelocalresidual.hh:25
auto createProjectionMatrices(const FEBasisDomain &feBasisDomain, const FEBasisTarget &feBasisTarget, const GlueType &glue, bool treatDiagonalZeroes=true)
Creates the matrices underlying l2-projections.
Definition projector.hh:310
auto makeProjectorPair(const FEBasisDomain &feBasisDomain, const FEBasisTarget &feBasisTarget, const GlueType &glue)
Creates a projector class between two function space bases.
Definition projector.hh:486
void setupReducedMatrices(const Matrix &massMatrix, const Matrix &projMatrix, const std::vector< bool > &dofIsVoid, Matrix &reducedM, Matrix &reducedP, std::vector< std::size_t > &expansionMap)
Reduces a mass matrix and projection matrix such that they are composed of only those dofs that actua...
Definition projector.hh:232
Definition adapt.hh:17
auto makeProjectionMatricesPair(const FEBasisDomain &feBasisDomain, const FEBasisTarget &feBasisTarget, GlueType glue)
Creates the matrices underlying l2-projections.
Definition projector.hh:655
const Scalar PengRobinsonMixture< Scalar, StaticParameters >::u
Definition pengrobinsonmixture.hh:138
auto makeProjectionMatrices(const FEBasisDomain &feBasisDomain, const FEBasisTarget &feBasisTarget, GlueType glue)
Creates the matrices underlying l2-projections.
Definition projector.hh:676
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Definition linearalgebratraits.hh:51
Parameters that can be passed to project()
Definition projector.hh:65
Scalar residualReduction
Definition projector.hh:67
int verbosity
Definition projector.hh:68
std::size_t maxIterations
Definition projector.hh:66
Definition linearsolvertraits.hh:53