18#ifndef DUMUX_DISCRETIZATION_PROJECTOR_HH
19#define DUMUX_DISCRETIZATION_PROJECTOR_HH
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>
32#include <dune/geometry/quadraturerules.hh>
33#include <dune/istl/matrixindexset.hh>
34#include <dune/istl/bcrsmatrix.hh>
35#include <dune/istl/bvector.hh>
52template<
class ScalarType>
55 using MatrixBlockType = Dune::FieldMatrix<ScalarType, 1, 1>;
80 : massMat_(std::make_shared<
Matrix>(std::move(massMatrix)))
81 , projMat_(std::make_shared<
Matrix>(std::move(projectionMatrix)))
82 , numDofsTarget_(massMat_->N())
84 if (massMat_->N() != projMat_->N())
85 DUNE_THROW(Dune::InvalidStateException,
"Matrix row size mismatch: " << massMat_->N() <<
" vs " << projMat_->N());
87 massMatrixSolver_.setMatrix(massMat_);
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)
109 if (indexMapTarget_.size() != massMat_->N())
110 DUNE_THROW(Dune::InvalidStateException,
"Target index map size mismatch: " << indexMapTarget_.size() <<
" vs " << massMat_->N());
112 if (massMat_->N() != projMat_->N())
113 DUNE_THROW(Dune::InvalidStateException,
"Matrix row size mismatch: " << massMat_->N() <<
" vs " << projMat_->N());
115 if (*std::max_element(indexMapTarget_.begin(), indexMapTarget_.end()) > numDofsTarget_)
116 DUNE_THROW(Dune::InvalidStateException,
"Index map exceeds provided number of dofs in target domain!");
118 massMatrixSolver_.setMatrix(massMat_);
127 template< class BlockType, std::enable_if_t<std::is_convertible<BlockType, ScalarType>::value,
int> = 0 >
128 Dune::BlockVector<BlockType>
project(
const Dune::BlockVector<BlockType>& u,
const Params& params =
Params{})
const
131 if ( u.size() != projMat_->M())
132 DUNE_THROW(Dune::InvalidStateException,
"Vector size mismatch");
134 Dune::BlockVector<BlockType> up(massMat_->N());
137 projMat_->mv(u, rhs);
139 Dune::ParameterTree solverParams;
140 solverParams[
"maxit"] = std::to_string(params.maxIterations);
141 solverParams[
"reduction"] = std::to_string(params.residualReduction);
142 solverParams[
"verbose"] = std::to_string(params.verbosity);
143 auto solver = massMatrixSolver_;
144 solver.setParams(solverParams);
145 solver.solve(up, rhs);
148 if (!indexMapTarget_.empty())
150 Dune::BlockVector<BlockType> result(numDofsTarget_);
153 for (std::size_t i = 0; i < indexMapTarget_.size(); ++i)
154 result[indexMapTarget_[i]] = up[i];
168 template< class BlockType, std::enable_if_t<!std::is_convertible<BlockType, ScalarType>::value,
int> = 0 >
169 Dune::BlockVector<BlockType>
project(
const Dune::BlockVector<BlockType>& u,
const Params& params =
Params{})
const
171 Dune::BlockVector<BlockType> result(numDofsTarget_);
173 for (
int pvIdx = 0; pvIdx < BlockType::size(); ++pvIdx)
175 Dune::BlockVector<Dune::FieldVector<Scalar, 1>> tmp(u.size());
176 std::transform(u.begin(), u.end(), tmp.begin(), [pvIdx] (
const auto& v) { return v[pvIdx]; });
178 const auto p =
project(tmp, params);
179 for (std::size_t i = 0; i < p.size(); ++i)
180 result[i][pvIdx] = p[i];
193 std::shared_ptr<Matrix> massMat_;
194 std::shared_ptr<Matrix> projMat_;
200 std::vector<std::size_t> indexMapTarget_;
201 std::size_t numDofsTarget_;
207template<
class FEBasisDomain,
class FEBasisTarget>
210 using FiniteElementDomain =
typename FEBasisDomain::LocalView::Tree::FiniteElement;
211 using FiniteElementTarget =
typename FEBasisTarget::LocalView::Tree::FiniteElement;
212 using ScalarDomain =
typename FiniteElementDomain::Traits::LocalBasisType::Traits::RangeFieldType;
213 using ScalarTarget =
typename FiniteElementTarget::Traits::LocalBasisType::Traits::RangeFieldType;
214 using Scalar =
typename Dune::PromotionTraits<ScalarDomain, ScalarTarget>::PromotedType;
230template<
class Matrix>
231void setupReducedMatrices(
const Matrix& massMatrix,
const Matrix& projMatrix,
const std::vector<bool>& dofIsVoid,
232 Matrix& reducedM, Matrix& reducedP, std::vector<std::size_t>& expansionMap)
234 const std::size_t numNonVoidDofs = std::count_if(dofIsVoid.begin(), dofIsVoid.end(), [] (
bool v) { return !v; });
237 std::vector<std::size_t> reductionMap(massMatrix.N());
238 expansionMap.resize(numNonVoidDofs);
240 std::size_t idxInReducedSpace = 0;
241 for (std::size_t dofIdx = 0; dofIdx < dofIsVoid.size(); ++dofIdx)
242 if (!dofIsVoid[dofIdx])
244 reductionMap[dofIdx] = idxInReducedSpace;
245 expansionMap[idxInReducedSpace] = dofIdx;
250 Dune::MatrixIndexSet patternMReduced, patternPReduced;
251 patternMReduced.resize(numNonVoidDofs, numNonVoidDofs);
252 patternPReduced.resize(numNonVoidDofs, projMatrix.M());
253 for (
auto rowIt = massMatrix.begin(); rowIt != massMatrix.end(); ++rowIt)
254 if (!dofIsVoid[rowIt.index()])
256 const auto reducedRowIdx = reductionMap[rowIt.index()];
257 for (
auto colIt = (*rowIt).begin(); colIt != (*rowIt).end(); ++colIt)
258 if (!dofIsVoid[colIt.index()])
259 patternMReduced.add(reducedRowIdx, reductionMap[colIt.index()]);
262 for (
auto rowIt = projMatrix.begin(); rowIt != projMatrix.end(); ++rowIt)
263 if (!dofIsVoid[rowIt.index()])
265 const auto reducedRowIdx = reductionMap[rowIt.index()];
266 for (
auto colIt = (*rowIt).begin(); colIt != (*rowIt).end(); ++colIt)
267 patternPReduced.add(reducedRowIdx, colIt.index());
270 patternMReduced.exportIdx(reducedM);
271 patternPReduced.exportIdx(reducedP);
274 for (
auto rowIt = massMatrix.begin(); rowIt != massMatrix.end(); ++rowIt)
275 if (!dofIsVoid[rowIt.index()])
277 const auto reducedRowIdx = reductionMap[rowIt.index()];
278 for (
auto colIt = (*rowIt).begin(); colIt != (*rowIt).end(); ++colIt)
279 if (!dofIsVoid[colIt.index()])
280 reducedM[reducedRowIdx][reductionMap[colIt.index()]] = *colIt;
283 for (
auto rowIt = projMatrix.begin(); rowIt != projMatrix.end(); ++rowIt)
284 if (!dofIsVoid[rowIt.index()])
286 const auto reducedRowIdx = reductionMap[rowIt.index()];
287 for (
auto colIt = (*rowIt).begin(); colIt != (*rowIt).end(); ++colIt)
288 reducedP[reducedRowIdx][colIt.index()] = *colIt;
308template<
bool doB
idirectional,
class FEBasisDomain,
class FEBasisTarget,
class GlueType>
310 const FEBasisTarget& feBasisTarget,
311 const GlueType& glue,
312 bool treatDiagonalZeroes =
true)
315 static constexpr int domainDim = FEBasisDomain::GridView::dimension;
316 static constexpr int targetDim = FEBasisTarget::GridView::dimension;
317 static_assert(targetDim <= domainDim,
"This expects target dim < domain dim, please swap arguments");
322 using ForwardProjectionMatrix =
typename ForwardProjector::Matrix;
323 using BackwardProjectionMatrix =
typename BackwardProjector::Matrix;
325 auto domainLocalView = feBasisDomain.localView();
326 auto targetLocalView = feBasisTarget.localView();
329 Dune::MatrixIndexSet backwardPatternM, forwardPatternM;
334 Dune::MatrixIndexSet backwardPatternP, forwardPatternP;
335 forwardPatternP.resize(feBasisTarget.size(), feBasisDomain.size());
336 if (doBidirectional) backwardPatternP.resize(feBasisDomain.size(), feBasisTarget.size());
339 unsigned int maxBasisOrder = 0;
340 for (
const auto& is : intersections(glue))
343 targetLocalView.bind( is.targetEntity(0) );
344 const auto& targetLocalBasis = targetLocalView.tree().finiteElement().localBasis();
346 for (
unsigned int nIdx = 0; nIdx < is.numDomainNeighbors(); ++nIdx)
348 domainLocalView.bind( is.domainEntity(nIdx) );
349 const auto& domainLocalBasis = domainLocalView.tree().finiteElement().localBasis();
352 maxBasisOrder = max(maxBasisOrder, max(domainLocalBasis.order(), targetLocalBasis.order()));
354 for (
unsigned int i = 0; i < domainLocalBasis.size(); ++i)
355 for (
unsigned int j = 0; j < targetLocalBasis.size(); ++j)
357 forwardPatternP.add(targetLocalView.index(j), domainLocalView.index(i));
358 if (doBidirectional) backwardPatternP.add(domainLocalView.index(i), targetLocalView.index(j));
364 ForwardProjectionMatrix forwardM, forwardP;
365 forwardPatternM.exportIdx(forwardM); forwardM = 0.0;
366 forwardPatternP.exportIdx(forwardP); forwardP = 0.0;
368 BackwardProjectionMatrix backwardM, backwardP;
371 backwardPatternM.exportIdx(backwardM); backwardM = 0.0;
372 backwardPatternP.exportIdx(backwardP); backwardP = 0.0;
375 for (
const auto& is : intersections(glue))
377 const auto& targetElement = is.targetEntity(0);
378 const auto& targetElementGeometry = targetElement.geometry();
380 targetLocalView.bind( targetElement );
381 const auto& targetLocalBasis = targetLocalView.tree().finiteElement().localBasis();
384 using IsGeometry =
typename std::decay_t<
decltype(is.geometry())>;
385 using ctype =
typename IsGeometry::ctype;
387 const auto& isGeometry = is.geometry();
388 const int intOrder = maxBasisOrder + 1;
389 const auto& quad = Dune::QuadratureRules<ctype, IsGeometry::mydimension>::rule(isGeometry.type(), intOrder);
390 for (
auto&& qp : quad)
392 const auto weight = qp.weight();
393 const auto ie = isGeometry.integrationElement(qp.position());
394 const auto globalPos = isGeometry.global(qp.position());
396 std::vector< Dune::FieldVector<ctype, 1> > targetShapeVals;
397 targetLocalBasis.evaluateFunction(targetElementGeometry.local(globalPos), targetShapeVals);
400 for (
unsigned int i = 0; i < targetLocalBasis.size(); ++i)
402 const auto dofIdxI = targetLocalView.index(i);
403 forwardM[dofIdxI][dofIdxI] += ie*weight*targetShapeVals[i]*targetShapeVals[i];
405 for (
unsigned int j = i+1; j < targetLocalBasis.size(); ++j)
407 const auto dofIdxJ = targetLocalView.index(j);
408 const auto value = ie*weight*targetShapeVals[i]*targetShapeVals[j];
409 forwardM[dofIdxI][dofIdxJ] += value;
410 forwardM[dofIdxJ][dofIdxI] += value;
418 const auto numNeighbors = is.numDomainNeighbors();
419 for (
unsigned int nIdx = 0; nIdx < numNeighbors; ++nIdx)
421 const auto& domainElement = is.domainEntity(nIdx);
422 domainLocalView.bind( domainElement );
423 const auto& domainLocalBasis = domainLocalView.tree().finiteElement().localBasis();
425 std::vector< Dune::FieldVector<ctype, 1> > domainShapeVals;
426 domainLocalBasis.evaluateFunction(domainElement.geometry().local(globalPos), domainShapeVals);
429 for (
unsigned int i = 0; i < domainLocalBasis.size(); ++i)
431 const auto dofIdxDomain = domainLocalView.index(i);
432 const auto domainShapeVal = domainShapeVals[i];
435 backwardM[dofIdxDomain][dofIdxDomain] += ie*weight*domainShapeVal*domainShapeVal;
437 for (
unsigned int j = i+1; j < domainLocalBasis.size(); ++j)
439 const auto dofIdxDomainJ = domainLocalView.index(j);
440 const auto value = ie*weight*domainShapeVal*domainShapeVals[j];
441 backwardM[dofIdxDomain][dofIdxDomainJ] += value;
442 backwardM[dofIdxDomainJ][dofIdxDomain] += value;
446 for (
unsigned int j = 0; j < targetLocalBasis.size(); ++j)
448 const auto dofIdxTarget = targetLocalView.index(j);
449 const auto entry = ie*weight*domainShapeVal*targetShapeVals[j];
451 forwardP[dofIdxTarget][dofIdxDomain] += entry/numNeighbors;
453 backwardP[dofIdxDomain][dofIdxTarget] += entry;
461 if (treatDiagonalZeroes)
463 for (std::size_t dofIdxTarget = 0; dofIdxTarget < forwardM.N(); ++dofIdxTarget)
464 if (forwardM[dofIdxTarget][dofIdxTarget] == 0.0)
465 forwardM[dofIdxTarget][dofIdxTarget] = 1.0;
469 for (std::size_t dofIdxDomain = 0; dofIdxDomain < backwardM.N(); ++dofIdxDomain)
470 if (backwardM[dofIdxDomain][dofIdxDomain] == 0.0)
471 backwardM[dofIdxDomain][dofIdxDomain] = 1.0;
475 return std::make_pair( std::make_pair(std::move(forwardM), std::move(forwardP)),
476 std::make_pair(std::move(backwardM), std::move(backwardP)) );
484template<
bool doB
idirectional,
class FEBasisDomain,
class FEBasisTarget,
class GlueType>
486 const FEBasisTarget& feBasisTarget,
487 const GlueType& glue)
492 using ForwardProjectionMatrix =
typename ForwardProjector::Matrix;
493 using BackwardProjectionMatrix =
typename BackwardProjector::Matrix;
495 auto projectionMatrices = createProjectionMatrices<doBidirectional>(feBasisDomain, feBasisTarget, glue,
false);
496 auto& forwardMatrices = projectionMatrices.first;
497 auto& backwardMatrices = projectionMatrices.second;
499 auto& forwardM = forwardMatrices.first;
500 auto& forwardP = forwardMatrices.second;
502 auto& backwardM = backwardMatrices.first;
503 auto& backwardP = backwardMatrices.second;
506 std::vector<bool> isVoidTarget(forwardM.N(),
false);
507 for (std::size_t dofIdxTarget = 0; dofIdxTarget < forwardM.N(); ++dofIdxTarget)
508 if (forwardM[dofIdxTarget][dofIdxTarget] == 0.0)
509 isVoidTarget[dofIdxTarget] =
true;
511 std::vector<bool> isVoidDomain;
514 isVoidDomain.resize(backwardM.N(),
false);
515 for (std::size_t dofIdxDomain = 0; dofIdxDomain < backwardM.N(); ++dofIdxDomain)
516 if (backwardM[dofIdxDomain][dofIdxDomain] == 0.0)
517 isVoidDomain[dofIdxDomain] =
true;
520 const bool hasVoidTarget = std::any_of(isVoidTarget.begin(), isVoidTarget.end(), [] (
bool v) { return v; });
521 const bool hasVoidDomain = std::any_of(isVoidDomain.begin(), isVoidDomain.end(), [] (
bool v) { return v; });
522 if (!hasVoidDomain && !hasVoidTarget)
524 return std::make_pair(ForwardProjector(std::move(forwardM), std::move(forwardP)),
525 BackwardProjector(std::move(backwardM), std::move(backwardP)));
527 else if (!hasVoidDomain && hasVoidTarget)
529 std::vector<std::size_t> expansionMapTarget;
530 ForwardProjectionMatrix forwardMReduced, forwardPReduced;
532 forwardMReduced, forwardPReduced, expansionMapTarget);
534 return std::make_pair( ForwardProjector(std::move(forwardMReduced),
535 std::move(forwardPReduced),
536 std::move(expansionMapTarget),
538 BackwardProjector(std::move(backwardM), std::move(backwardP)) );
540 else if (hasVoidDomain && !hasVoidTarget)
544 std::vector<std::size_t> expansionMapDomain;
545 BackwardProjectionMatrix backwardMReduced, backwardPReduced;
547 backwardMReduced, backwardPReduced, expansionMapDomain);
549 return std::make_pair( ForwardProjector(std::move(forwardM), std::move(forwardP)),
550 BackwardProjector(std::move(backwardMReduced),
551 std::move(backwardPReduced),
552 std::move(expansionMapDomain),
556 return std::make_pair( ForwardProjector(std::move(forwardM), std::move(forwardP)),
557 BackwardProjector(std::move(backwardM), std::move(backwardP)) );
561 std::vector<std::size_t> expansionMapTarget;
562 ForwardProjectionMatrix forwardMReduced, forwardPReduced;
564 forwardMReduced, forwardPReduced, expansionMapTarget);
568 std::vector<std::size_t> expansionMapDomain;
569 BackwardProjectionMatrix backwardMReduced, backwardPReduced;
571 backwardMReduced, backwardPReduced, expansionMapDomain);
573 return std::make_pair( ForwardProjector(std::move(forwardMReduced),
574 std::move(forwardPReduced),
575 std::move(expansionMapTarget),
577 BackwardProjector(std::move(backwardMReduced),
578 std::move(backwardPReduced),
579 std::move(expansionMapDomain),
583 return std::make_pair( ForwardProjector(std::move(forwardMReduced),
584 std::move(forwardPReduced),
585 std::move(expansionMapTarget),
587 BackwardProjector(std::move(backwardM), std::move(backwardP)) );
606template<
class FEBasisDomain,
class FEBasisTarget,
class GlueType >
608 const FEBasisTarget& feBasisTarget,
612 static constexpr int domainDim = FEBasisDomain::GridView::dimension;
613 static constexpr int targetDim = FEBasisTarget::GridView::dimension;
614 static_assert(targetDim <= domainDim,
"makeProjectorPair() expects targetDim < domainDim, please swap arguments");
616 return Detail::makeProjectorPair<true>(feBasisDomain, feBasisTarget, glue);
629template<
class FEBasisDomain,
class FEBasisTarget,
class GlueType >
631 const FEBasisTarget& feBasisTarget,
635 static constexpr int domainDim = FEBasisDomain::GridView::dimension;
636 static constexpr int targetDim = FEBasisTarget::GridView::dimension;
637 static_assert(targetDim <= domainDim,
"makeProjectorPair() expects targetDim < domainDim, please swap arguments");
639 return Detail::makeProjectorPair<false>(feBasisDomain, feBasisTarget, glue).first;
653template<
class FEBasisDomain,
class FEBasisTarget,
class GlueType >
655 const FEBasisTarget& feBasisTarget,
659 static constexpr int domainDim = FEBasisDomain::GridView::dimension;
660 static constexpr int targetDim = FEBasisTarget::GridView::dimension;
661 static_assert(targetDim <= domainDim,
"makeProjectionMatrixPair() expects targetDim < domainDim, please swap arguments");
663 return Detail::createProjectionMatrices<true>(feBasisDomain, feBasisTarget, glue);
674template<
class FEBasisDomain,
class FEBasisTarget,
class GlueType >
676 const FEBasisTarget& feBasisTarget,
680 static constexpr int domainDim = FEBasisDomain::GridView::dimension;
681 static constexpr int targetDim = FEBasisTarget::GridView::dimension;
682 static_assert(targetDim <= domainDim,
"makeProjectionMatrixPair() expects targetDim < domainDim, please swap arguments");
684 return Detail::createProjectionMatrices<false>(feBasisDomain, feBasisTarget, glue).first;
Standard dune-istl iterative linear solvers.
Definition: istlsolvers.hh:159
Does an L2-projection from one discrete function space into another. The convenience functions makePr...
Definition: projector.hh:54
static Params defaultParams()
Returns the default parameters.
Definition: projector.hh:189
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
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:128
Traits class stating the type of projector between to bases.
Definition: projector.hh:209
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:630
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:607
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.
auto createProjectionMatrices(const FEBasisDomain &feBasisDomain, const FEBasisTarget &feBasisTarget, const GlueType &glue, bool treatDiagonalZeroes=true)
Creates the matrices underlying l2-projections.
Definition: projector.hh:309
auto makeProjectorPair(const FEBasisDomain &feBasisDomain, const FEBasisTarget &feBasisTarget, const GlueType &glue)
Creates a projector class between two function space bases.
Definition: projector.hh:485
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:231
auto makeProjectionMatricesPair(const FEBasisDomain &feBasisDomain, const FEBasisTarget &feBasisTarget, GlueType glue)
Creates the matrices underlying l2-projections.
Definition: projector.hh:654
auto makeProjectionMatrices(const FEBasisDomain &feBasisDomain, const FEBasisTarget &feBasisTarget, GlueType glue)
Creates the matrices underlying l2-projections.
Definition: projector.hh:675
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