30#ifndef DUMUX_DISCRETIZATION_PROJECTOR_HH
31#define DUMUX_DISCRETIZATION_PROJECTOR_HH
38#include <dune/common/timer.hh>
39#include <dune/common/fmatrix.hh>
40#include <dune/common/exceptions.hh>
41#include <dune/common/promotiontraits.hh>
43#include <dune/geometry/quadraturerules.hh>
44#include <dune/istl/matrixindexset.hh>
45#include <dune/istl/bcrsmatrix.hh>
46#include <dune/istl/bvector.hh>
61template<
class ScalarType>
64 using MatrixBlockType = Dune::FieldMatrix<ScalarType, 1, 1>;
88 : massMat_(std::move(massMatrix))
89 , projMat_(std::move(projectionMatrix))
90 , numDofsTarget_(massMat_.N())
92 if (massMat_.N() != projMat_.N())
93 DUNE_THROW(Dune::InvalidStateException,
"Matrix row size mismatch: " << massMat_.N() <<
" vs " << projMat_.N());
107 Matrix&& projectionMatrix,
108 std::vector<std::size_t>&& indexMap,
109 std::size_t numDofsTarget)
110 : massMat_(std::move(massMatrix))
111 , projMat_(std::move(projectionMatrix))
112 , indexMapTarget_(std::move(indexMap))
113 , numDofsTarget_(numDofsTarget)
115 if (indexMapTarget_.size() != massMat_.N())
116 DUNE_THROW(Dune::InvalidStateException,
"Target index map size mismatch: " << indexMapTarget_.size() <<
" vs " << massMat_.N());
118 if (massMat_.N() != projMat_.N())
119 DUNE_THROW(Dune::InvalidStateException,
"Matrix row size mismatch: " << massMat_.N() <<
" vs " << projMat_.N());
121 if (*std::max_element(indexMapTarget_.begin(), indexMapTarget_.end()) > numDofsTarget_)
122 DUNE_THROW(Dune::InvalidStateException,
"Index map exceeds provided number of dofs in target domain!");
131 template< class BlockType, std::enable_if_t<std::is_convertible<BlockType, ScalarType>::value,
int> = 0 >
132 Dune::BlockVector<BlockType>
project(
const Dune::BlockVector<BlockType>& u,
const Params& params =
Params{})
const
135 if ( u.size() != projMat_.M())
136 DUNE_THROW(Dune::InvalidStateException,
"Vector size mismatch");
138 Dune::BlockVector<BlockType> up(massMat_.N());
146 solver.
solve(massMat_, up, rhs);
149 if (!indexMapTarget_.empty())
151 Dune::BlockVector<BlockType> result(numDofsTarget_);
154 for (std::size_t i = 0; i < indexMapTarget_.size(); ++i)
155 result[indexMapTarget_[i]] = up[i];
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
172 Dune::BlockVector<BlockType> result(numDofsTarget_);
174 for (
int pvIdx = 0; pvIdx < BlockType::size(); ++pvIdx)
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]; });
179 const auto p =
project(tmp, params);
180 for (std::size_t i = 0; i < p.size(); ++i)
181 result[i][pvIdx] = p[i];
197 std::vector<std::size_t> indexMapTarget_;
198 std::size_t numDofsTarget_;
204template<
class FEBasisDomain,
class FEBasisTarget>
207 using FiniteElementDomain =
typename FEBasisDomain::LocalView::Tree::FiniteElement;
208 using FiniteElementTarget =
typename FEBasisTarget::LocalView::Tree::FiniteElement;
209 using ScalarDomain =
typename FiniteElementDomain::Traits::LocalBasisType::Traits::RangeFieldType;
210 using ScalarTarget =
typename FiniteElementTarget::Traits::LocalBasisType::Traits::RangeFieldType;
211 using Scalar =
typename Dune::PromotionTraits<ScalarDomain, ScalarTarget>::PromotedType;
227template<
class Matrix>
228void setupReducedMatrices(
const Matrix& massMatrix,
const Matrix& projMatrix,
const std::vector<bool>& dofIsVoid,
229 Matrix& reducedM, Matrix& reducedP, std::vector<std::size_t>& expansionMap)
231 const std::size_t numNonVoidDofs = std::count_if(dofIsVoid.begin(), dofIsVoid.end(), [] (
bool v) { return !v; });
234 std::vector<std::size_t> reductionMap(massMatrix.N());
235 expansionMap.resize(numNonVoidDofs);
237 std::size_t idxInReducedSpace = 0;
238 for (std::size_t dofIdx = 0; dofIdx < dofIsVoid.size(); ++dofIdx)
239 if (!dofIsVoid[dofIdx])
241 reductionMap[dofIdx] = idxInReducedSpace;
242 expansionMap[idxInReducedSpace] = dofIdx;
247 Dune::MatrixIndexSet patternMReduced, patternPReduced;
248 patternMReduced.resize(numNonVoidDofs, numNonVoidDofs);
249 patternPReduced.resize(numNonVoidDofs, projMatrix.M());
250 for (
auto rowIt = massMatrix.begin(); rowIt != massMatrix.end(); ++rowIt)
251 if (!dofIsVoid[rowIt.index()])
253 const auto reducedRowIdx = reductionMap[rowIt.index()];
254 for (
auto colIt = (*rowIt).begin(); colIt != (*rowIt).end(); ++colIt)
255 if (!dofIsVoid[colIt.index()])
256 patternMReduced.add(reducedRowIdx, reductionMap[colIt.index()]);
259 for (
auto rowIt = projMatrix.begin(); rowIt != projMatrix.end(); ++rowIt)
260 if (!dofIsVoid[rowIt.index()])
262 const auto reducedRowIdx = reductionMap[rowIt.index()];
263 for (
auto colIt = (*rowIt).begin(); colIt != (*rowIt).end(); ++colIt)
264 patternPReduced.add(reducedRowIdx, colIt.index());
267 patternMReduced.exportIdx(reducedM);
268 patternPReduced.exportIdx(reducedP);
271 for (
auto rowIt = massMatrix.begin(); rowIt != massMatrix.end(); ++rowIt)
272 if (!dofIsVoid[rowIt.index()])
274 const auto reducedRowIdx = reductionMap[rowIt.index()];
275 for (
auto colIt = (*rowIt).begin(); colIt != (*rowIt).end(); ++colIt)
276 if (!dofIsVoid[colIt.index()])
277 reducedM[reducedRowIdx][reductionMap[colIt.index()]] = *colIt;
280 for (
auto rowIt = projMatrix.begin(); rowIt != projMatrix.end(); ++rowIt)
281 if (!dofIsVoid[rowIt.index()])
283 const auto reducedRowIdx = reductionMap[rowIt.index()];
284 for (
auto colIt = (*rowIt).begin(); colIt != (*rowIt).end(); ++colIt)
285 reducedP[reducedRowIdx][colIt.index()] = *colIt;
305template<
bool doB
idirectional,
class FEBasisDomain,
class FEBasisTarget,
class GlueType>
306auto createProjectionMatrices(
const FEBasisDomain& feBasisDomain,
307 const FEBasisTarget& feBasisTarget,
308 const GlueType& glue,
309 bool treatDiagonalZeroes =
true)
312 static constexpr int domainDim = FEBasisDomain::GridView::dimension;
313 static constexpr int targetDim = FEBasisTarget::GridView::dimension;
314 static_assert(targetDim <= domainDim,
"This expects target dim < domain dim, please swap arguments");
319 using ForwardProjectionMatrix =
typename ForwardProjector::Matrix;
320 using BackwardProjectionMatrix =
typename BackwardProjector::Matrix;
322 auto domainLocalView = feBasisDomain.localView();
323 auto targetLocalView = feBasisTarget.localView();
326 Dune::MatrixIndexSet backwardPatternM, forwardPatternM;
331 Dune::MatrixIndexSet backwardPatternP, forwardPatternP;
332 forwardPatternP.resize(feBasisTarget.size(), feBasisDomain.size());
333 if (doBidirectional) backwardPatternP.resize(feBasisDomain.size(), feBasisTarget.size());
336 unsigned int maxBasisOrder = 0;
337 for (
const auto& is : intersections(glue))
340 targetLocalView.bind( is.targetEntity(0) );
341 const auto& targetLocalBasis = targetLocalView.tree().finiteElement().localBasis();
343 for (
unsigned int nIdx = 0; nIdx < is.numDomainNeighbors(); ++nIdx)
345 domainLocalView.bind( is.domainEntity(nIdx) );
346 const auto& domainLocalBasis = domainLocalView.tree().finiteElement().localBasis();
349 maxBasisOrder = max(maxBasisOrder, max(domainLocalBasis.order(), targetLocalBasis.order()));
351 for (
unsigned int i = 0; i < domainLocalBasis.size(); ++i)
352 for (
unsigned int j = 0; j < targetLocalBasis.size(); ++j)
354 forwardPatternP.add(targetLocalView.index(j), domainLocalView.index(i));
355 if (doBidirectional) backwardPatternP.add(domainLocalView.index(i), targetLocalView.index(j));
361 ForwardProjectionMatrix forwardM, forwardP;
362 forwardPatternM.exportIdx(forwardM); forwardM = 0.0;
363 forwardPatternP.exportIdx(forwardP); forwardP = 0.0;
365 BackwardProjectionMatrix backwardM, backwardP;
368 backwardPatternM.exportIdx(backwardM); backwardM = 0.0;
369 backwardPatternP.exportIdx(backwardP); backwardP = 0.0;
372 for (
const auto& is : intersections(glue))
374 const auto& targetElement = is.targetEntity(0);
375 const auto& targetElementGeometry = targetElement.geometry();
377 targetLocalView.bind( targetElement );
378 const auto& targetLocalBasis = targetLocalView.tree().finiteElement().localBasis();
381 using IsGeometry =
typename std::decay_t<
decltype(is.geometry())>;
382 using ctype =
typename IsGeometry::ctype;
384 const auto& isGeometry = is.geometry();
385 const int intOrder = maxBasisOrder + 1;
386 const auto& quad = Dune::QuadratureRules<ctype, IsGeometry::mydimension>::rule(isGeometry.type(), intOrder);
387 for (
auto&& qp : quad)
389 const auto weight = qp.weight();
390 const auto ie = isGeometry.integrationElement(qp.position());
391 const auto globalPos = isGeometry.global(qp.position());
393 std::vector< Dune::FieldVector<ctype, 1> > targetShapeVals;
394 targetLocalBasis.evaluateFunction(targetElementGeometry.local(globalPos), targetShapeVals);
397 for (
unsigned int i = 0; i < targetLocalBasis.size(); ++i)
399 const auto dofIdxI = targetLocalView.index(i);
400 forwardM[dofIdxI][dofIdxI] += ie*weight*targetShapeVals[i]*targetShapeVals[i];
402 for (
unsigned int j = i+1; j < targetLocalBasis.size(); ++j)
404 const auto dofIdxJ = targetLocalView.index(j);
405 const auto value = ie*weight*targetShapeVals[i]*targetShapeVals[j];
406 forwardM[dofIdxI][dofIdxJ] += value;
407 forwardM[dofIdxJ][dofIdxI] += value;
415 const auto numNeighbors = is.numDomainNeighbors();
416 for (
unsigned int nIdx = 0; nIdx < numNeighbors; ++nIdx)
418 const auto& domainElement = is.domainEntity(nIdx);
419 domainLocalView.bind( domainElement );
420 const auto& domainLocalBasis = domainLocalView.tree().finiteElement().localBasis();
422 std::vector< Dune::FieldVector<ctype, 1> > domainShapeVals;
423 domainLocalBasis.evaluateFunction(domainElement.geometry().local(globalPos), domainShapeVals);
426 for (
unsigned int i = 0; i < domainLocalBasis.size(); ++i)
428 const auto dofIdxDomain = domainLocalView.index(i);
429 const auto domainShapeVal = domainShapeVals[i];
432 backwardM[dofIdxDomain][dofIdxDomain] += ie*weight*domainShapeVal*domainShapeVal;
434 for (
unsigned int j = i+1; j < domainLocalBasis.size(); ++j)
436 const auto dofIdxDomainJ = domainLocalView.index(j);
437 const auto value = ie*weight*domainShapeVal*domainShapeVals[j];
438 backwardM[dofIdxDomain][dofIdxDomainJ] += value;
439 backwardM[dofIdxDomainJ][dofIdxDomain] += value;
443 for (
unsigned int j = 0; j < targetLocalBasis.size(); ++j)
445 const auto dofIdxTarget = targetLocalView.index(j);
446 const auto entry = ie*weight*domainShapeVal*targetShapeVals[j];
448 forwardP[dofIdxTarget][dofIdxDomain] += entry/numNeighbors;
450 backwardP[dofIdxDomain][dofIdxTarget] += entry;
458 if (treatDiagonalZeroes)
460 for (std::size_t dofIdxTarget = 0; dofIdxTarget < forwardM.N(); ++dofIdxTarget)
461 if (forwardM[dofIdxTarget][dofIdxTarget] == 0.0)
462 forwardM[dofIdxTarget][dofIdxTarget] = 1.0;
466 for (std::size_t dofIdxDomain = 0; dofIdxDomain < backwardM.N(); ++dofIdxDomain)
467 if (backwardM[dofIdxDomain][dofIdxDomain] == 0.0)
468 backwardM[dofIdxDomain][dofIdxDomain] = 1.0;
472 return std::make_pair( std::make_pair(std::move(forwardM), std::move(forwardP)),
473 std::make_pair(std::move(backwardM), std::move(backwardP)) );
481template<
bool doB
idirectional,
class FEBasisDomain,
class FEBasisTarget,
class GlueType>
483 const FEBasisTarget& feBasisTarget,
484 const GlueType& glue)
489 using ForwardProjectionMatrix =
typename ForwardProjector::Matrix;
490 using BackwardProjectionMatrix =
typename BackwardProjector::Matrix;
492 auto projectionMatrices = createProjectionMatrices<doBidirectional>(feBasisDomain, feBasisTarget, glue,
false);
493 auto& forwardMatrices = projectionMatrices.first;
494 auto& backwardMatrices = projectionMatrices.second;
496 auto& forwardM = forwardMatrices.first;
497 auto& forwardP = forwardMatrices.second;
499 auto& backwardM = backwardMatrices.first;
500 auto& backwardP = backwardMatrices.second;
503 std::vector<bool> isVoidTarget(forwardM.N(),
false);
504 for (std::size_t dofIdxTarget = 0; dofIdxTarget < forwardM.N(); ++dofIdxTarget)
505 if (forwardM[dofIdxTarget][dofIdxTarget] == 0.0)
506 isVoidTarget[dofIdxTarget] =
true;
508 std::vector<bool> isVoidDomain;
511 isVoidDomain.resize(backwardM.N(),
false);
512 for (std::size_t dofIdxDomain = 0; dofIdxDomain < backwardM.N(); ++dofIdxDomain)
513 if (backwardM[dofIdxDomain][dofIdxDomain] == 0.0)
514 isVoidDomain[dofIdxDomain] =
true;
517 const bool hasVoidTarget = std::any_of(isVoidTarget.begin(), isVoidTarget.end(), [] (
bool v) { return v; });
518 const bool hasVoidDomain = std::any_of(isVoidDomain.begin(), isVoidDomain.end(), [] (
bool v) { return v; });
519 if (!hasVoidDomain && !hasVoidTarget)
521 return std::make_pair(ForwardProjector(std::move(forwardM), std::move(forwardP)),
522 BackwardProjector(std::move(backwardM), std::move(backwardP)));
524 else if (!hasVoidDomain && hasVoidTarget)
526 std::vector<std::size_t> expansionMapTarget;
527 ForwardProjectionMatrix forwardMReduced, forwardPReduced;
528 setupReducedMatrices(forwardM, forwardP, isVoidTarget,
529 forwardMReduced, forwardPReduced, expansionMapTarget);
531 return std::make_pair( ForwardProjector(std::move(forwardMReduced),
532 std::move(forwardPReduced),
533 std::move(expansionMapTarget),
535 BackwardProjector(std::move(backwardM), std::move(backwardP)) );
537 else if (hasVoidDomain && !hasVoidTarget)
541 std::vector<std::size_t> expansionMapDomain;
542 BackwardProjectionMatrix backwardMReduced, backwardPReduced;
543 setupReducedMatrices(backwardM, backwardP, isVoidDomain,
544 backwardMReduced, backwardPReduced, expansionMapDomain);
546 return std::make_pair( ForwardProjector(std::move(forwardM), std::move(forwardP)),
547 BackwardProjector(std::move(backwardMReduced),
548 std::move(backwardPReduced),
549 std::move(expansionMapDomain),
553 return std::make_pair( ForwardProjector(std::move(forwardM), std::move(forwardP)),
554 BackwardProjector(std::move(backwardM), std::move(backwardP)) );
558 std::vector<std::size_t> expansionMapTarget;
559 ForwardProjectionMatrix forwardMReduced, forwardPReduced;
560 setupReducedMatrices(forwardM, forwardP, isVoidTarget,
561 forwardMReduced, forwardPReduced, expansionMapTarget);
565 std::vector<std::size_t> expansionMapDomain;
566 BackwardProjectionMatrix backwardMReduced, backwardPReduced;
567 setupReducedMatrices(backwardM, backwardP, isVoidDomain,
568 backwardMReduced, backwardPReduced, expansionMapDomain);
570 return std::make_pair( ForwardProjector(std::move(forwardMReduced),
571 std::move(forwardPReduced),
572 std::move(expansionMapTarget),
574 BackwardProjector(std::move(backwardMReduced),
575 std::move(backwardPReduced),
576 std::move(expansionMapDomain),
580 return std::make_pair( ForwardProjector(std::move(forwardMReduced),
581 std::move(forwardPReduced),
582 std::move(expansionMapTarget),
584 BackwardProjector(std::move(backwardM), std::move(backwardP)) );
603template<
class FEBasisDomain,
class FEBasisTarget,
class GlueType >
605 const FEBasisTarget& feBasisTarget,
609 static constexpr int domainDim = FEBasisDomain::GridView::dimension;
610 static constexpr int targetDim = FEBasisTarget::GridView::dimension;
611 static_assert(targetDim <= domainDim,
"makeProjectorPair() expects targetDim < domainDim, please swap arguments");
613 return Impl::makeProjectorPair<true>(feBasisDomain, feBasisTarget, glue);
626template<
class FEBasisDomain,
class FEBasisTarget,
class GlueType >
628 const FEBasisTarget& feBasisTarget,
632 static constexpr int domainDim = FEBasisDomain::GridView::dimension;
633 static constexpr int targetDim = FEBasisTarget::GridView::dimension;
634 static_assert(targetDim <= domainDim,
"makeProjectorPair() expects targetDim < domainDim, please swap arguments");
636 return Impl::makeProjectorPair<false>(feBasisDomain, feBasisTarget, glue).first;
650template<
class FEBasisDomain,
class FEBasisTarget,
class GlueType >
652 const FEBasisTarget& feBasisTarget,
656 static constexpr int domainDim = FEBasisDomain::GridView::dimension;
657 static constexpr int targetDim = FEBasisTarget::GridView::dimension;
658 static_assert(targetDim <= domainDim,
"makeProjectionMatrixPair() expects targetDim < domainDim, please swap arguments");
660 return Impl::createProjectionMatrices<true>(feBasisDomain, feBasisTarget, glue);
671template<
class FEBasisDomain,
class FEBasisTarget,
class GlueType >
673 const FEBasisTarget& feBasisTarget,
677 static constexpr int domainDim = FEBasisDomain::GridView::dimension;
678 static constexpr int targetDim = FEBasisTarget::GridView::dimension;
679 static_assert(targetDim <= domainDim,
"makeProjectionMatrixPair() expects targetDim < domainDim, please swap arguments");
681 return Impl::createProjectionMatrices<false>(feBasisDomain, feBasisTarget, glue).first;
Helper function to generate Jacobian pattern for different discretization methods.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Provides helper aliases and functionality to obtain the types and instances of Dune::Functions functi...
Dumux sequential linear solver backends.
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:627
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:604
Dune::MatrixIndexSet getFEJacobianPattern(const FEBasis &feBasis)
Helper function to generate Jacobian pattern for finite element scheme.
Definition: jacobianpattern.hh:167
auto makeProjectionMatricesPair(const FEBasisDomain &feBasisDomain, const FEBasisTarget &feBasisTarget, GlueType glue)
Creates the matrices underlying l2-projections.
Definition: projector.hh:651
auto makeProjectionMatrices(const FEBasisDomain &feBasisDomain, const FEBasisTarget &feBasisTarget, GlueType glue)
Creates the matrices underlying l2-projections.
Definition: projector.hh:672
Does an L2-projection from one discrete function space into another. The convenience functions makePr...
Definition: projector.hh:63
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:106
ScalarType Scalar
Export the scalar type.
Definition: projector.hh:68
Projector(Matrix &&massMatrix, Matrix &&projectionMatrix)
Constructor. Receives the mass and projection matrix that define the linear system describing the L2-...
Definition: projector.hh:87
Dune::BCRSMatrix< MatrixBlockType > Matrix
Export the type of the projection matrices.
Definition: projector.hh:70
Projector()=delete
delete default constructor
Dune::BlockVector< BlockType > project(const Dune::BlockVector< BlockType > &u, const Params ¶ms=Params{}) const
Project a solution u into up.
Definition: projector.hh:132
Parameters that can be passed to project()
Definition: projector.hh:74
Scalar residualReduction
Definition: projector.hh:76
std::size_t maxIterations
Definition: projector.hh:75
Traits class stating the type of projector between to bases.
Definition: projector.hh:206
Dumux::Projector< Scalar > Projector
Definition: projector.hh:213
Sequential SSOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:474
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:479
void setResidualReduction(double r)
set the linear solver residual reduction
Definition: solver.hh:135
void setMaxIter(int i)
set the maximum number of linear solver iterations
Definition: solver.hh:127