232void setupReducedMatrices(
const Matrix& massMatrix, 
const Matrix& projMatrix, 
const std::vector<bool>& dofIsVoid,
 
  233                          Matrix& reducedM, Matrix& reducedP, std::vector<std::size_t>& expansionMap)
 
  235    const std::size_t numNonVoidDofs = std::count_if(dofIsVoid.begin(), dofIsVoid.end(), [] (
bool v) { return !v; });
 
  238    std::vector<std::size_t> reductionMap(massMatrix.N());
 
  239    expansionMap.resize(numNonVoidDofs);
 
  241    std::size_t idxInReducedSpace = 0;
 
  242    for (std::size_t dofIdx = 0; dofIdx < dofIsVoid.size(); ++dofIdx)
 
  243        if (!dofIsVoid[dofIdx])
 
  245            reductionMap[dofIdx] = idxInReducedSpace;
 
  246            expansionMap[idxInReducedSpace] = dofIdx;
 
  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()])
 
  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()]);
 
  263    for (
auto rowIt = projMatrix.begin(); rowIt != projMatrix.end(); ++rowIt)
 
  264        if (!dofIsVoid[rowIt.index()])
 
  266            const auto reducedRowIdx = reductionMap[rowIt.index()];
 
  267            for (
auto colIt = (*rowIt).begin(); colIt != (*rowIt).end(); ++colIt)
 
  268                patternPReduced.add(reducedRowIdx, colIt.index());
 
  271    patternMReduced.exportIdx(reducedM);
 
  272    patternPReduced.exportIdx(reducedP);
 
  275    for (
auto rowIt = massMatrix.begin(); rowIt != massMatrix.end(); ++rowIt)
 
  276        if (!dofIsVoid[rowIt.index()])
 
  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;
 
  284    for (
auto rowIt = projMatrix.begin(); rowIt != projMatrix.end(); ++rowIt)
 
  285        if (!dofIsVoid[rowIt.index()])
 
  287            const auto reducedRowIdx = reductionMap[rowIt.index()];
 
  288            for (
auto colIt = (*rowIt).begin(); colIt != (*rowIt).end(); ++colIt)
 
  289                reducedP[reducedRowIdx][colIt.index()] = *colIt;
 
 
  311                              const FEBasisTarget& feBasisTarget,
 
  312                              const GlueType& glue,
 
  313                              bool treatDiagonalZeroes = 
true)
 
  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");
 
  323    using ForwardProjectionMatrix = 
typename ForwardProjector::Matrix;
 
  324    using BackwardProjectionMatrix = 
typename BackwardProjector::Matrix;
 
  326    auto domainLocalView = feBasisDomain.localView();
 
  327    auto targetLocalView = feBasisTarget.localView();
 
  330    Dune::MatrixIndexSet backwardPatternM, forwardPatternM;
 
  335    Dune::MatrixIndexSet backwardPatternP, forwardPatternP;
 
  336    forwardPatternP.resize(feBasisTarget.size(), feBasisDomain.size());
 
  337    if (doBidirectional) backwardPatternP.resize(feBasisDomain.size(), feBasisTarget.size());
 
  340    unsigned int maxBasisOrder = 0;
 
  341    for (
const auto& is : intersections(glue))
 
  344        targetLocalView.bind( is.targetEntity(0) );
 
  345        const auto& targetLocalBasis = targetLocalView.tree().finiteElement().localBasis();
 
  347        for (
unsigned int nIdx = 0; nIdx < is.numDomainNeighbors(); ++nIdx)
 
  349            domainLocalView.bind( is.domainEntity(nIdx) );
 
  350            const auto& domainLocalBasis = domainLocalView.tree().finiteElement().localBasis();
 
  353            maxBasisOrder = max(maxBasisOrder, max(domainLocalBasis.order(), targetLocalBasis.order()));
 
  355            for (
unsigned int i = 0; i < domainLocalBasis.size(); ++i)
 
  356                for (
unsigned int j = 0; j < targetLocalBasis.size(); ++j)
 
  358                    forwardPatternP.add(targetLocalView.index(j), domainLocalView.index(i));
 
  359                    if (doBidirectional) backwardPatternP.add(domainLocalView.index(i), targetLocalView.index(j));
 
  365    ForwardProjectionMatrix forwardM, forwardP;
 
  366    forwardPatternM.exportIdx(forwardM); forwardM = 0.0;
 
  367    forwardPatternP.exportIdx(forwardP); forwardP = 0.0;
 
  369    BackwardProjectionMatrix backwardM, backwardP;
 
  372        backwardPatternM.exportIdx(backwardM); backwardM = 0.0;
 
  373        backwardPatternP.exportIdx(backwardP); backwardP = 0.0;
 
  376    for (
const auto& is : intersections(glue))
 
  378        const auto& targetElement = is.targetEntity(0);
 
  379        const auto& targetElementGeometry = targetElement.geometry();
 
  381        targetLocalView.bind( targetElement );
 
  382        const auto& targetLocalBasis = targetLocalView.tree().finiteElement().localBasis();
 
  385        using IsGeometry = 
typename std::decay_t<
decltype(is.geometry())>;
 
  386        using ctype = 
typename IsGeometry::ctype;
 
  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)
 
  393            const auto weight = qp.weight();
 
  394            const auto ie = isGeometry.integrationElement(qp.position());
 
  395            const auto globalPos = isGeometry.global(qp.position());
 
  397            std::vector< Dune::FieldVector<ctype, 1> > targetShapeVals;
 
  398            targetLocalBasis.evaluateFunction(targetElementGeometry.local(globalPos), targetShapeVals);
 
  401            for (
unsigned int i = 0; i < targetLocalBasis.size(); ++i)
 
  403                const auto dofIdxI = targetLocalView.index(i);
 
  404                forwardM[dofIdxI][dofIdxI] += ie*weight*targetShapeVals[i]*targetShapeVals[i];
 
  406                for (
unsigned int j = i+1; j < targetLocalBasis.size(); ++j)
 
  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;
 
  419            const auto numNeighbors = is.numDomainNeighbors();
 
  420            for (
unsigned int nIdx = 0; nIdx < numNeighbors; ++nIdx)
 
  422                const auto& domainElement = is.domainEntity(nIdx);
 
  423                domainLocalView.bind( domainElement );
 
  424                const auto& domainLocalBasis = domainLocalView.tree().finiteElement().localBasis();
 
  426                std::vector< Dune::FieldVector<ctype, 1> > domainShapeVals;
 
  427                domainLocalBasis.evaluateFunction(domainElement.geometry().local(globalPos), domainShapeVals);
 
  430                for (
unsigned int i = 0; i < domainLocalBasis.size(); ++i)
 
  432                    const auto dofIdxDomain = domainLocalView.index(i);
 
  433                    const auto domainShapeVal = domainShapeVals[i];
 
  436                        backwardM[dofIdxDomain][dofIdxDomain] += ie*weight*domainShapeVal*domainShapeVal;
 
  438                        for (
unsigned int j = i+1; j < domainLocalBasis.size(); ++j)
 
  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;
 
  447                    for (
unsigned int j = 0; j < targetLocalBasis.size(); ++j)
 
  449                        const auto dofIdxTarget = targetLocalView.index(j);
 
  450                        const auto entry = ie*weight*domainShapeVal*targetShapeVals[j];
 
  452                        forwardP[dofIdxTarget][dofIdxDomain] += entry/numNeighbors;
 
  454                            backwardP[dofIdxDomain][dofIdxTarget] += entry;
 
  462    if (treatDiagonalZeroes)
 
  464        for (std::size_t dofIdxTarget = 0; dofIdxTarget < forwardM.N(); ++dofIdxTarget)
 
  465            if (forwardM[dofIdxTarget][dofIdxTarget] == 0.0)
 
  466                forwardM[dofIdxTarget][dofIdxTarget] = 1.0;
 
  470            for (std::size_t dofIdxDomain = 0; dofIdxDomain < backwardM.N(); ++dofIdxDomain)
 
  471                if (backwardM[dofIdxDomain][dofIdxDomain] == 0.0)
 
  472                    backwardM[dofIdxDomain][dofIdxDomain] = 1.0;
 
  476    return std::make_pair( std::make_pair(std::move(forwardM), std::move(forwardP)),
 
  477                           std::make_pair(std::move(backwardM), std::move(backwardP)) );
 
 
  487                       const FEBasisTarget& feBasisTarget,
 
  488                       const GlueType& glue)
 
  493    using ForwardProjectionMatrix = 
typename ForwardProjector::Matrix;
 
  494    using BackwardProjectionMatrix = 
typename BackwardProjector::Matrix;
 
  497    auto& forwardMatrices = projectionMatrices.first;
 
  498    auto& backwardMatrices = projectionMatrices.second;
 
  500    auto& forwardM = forwardMatrices.first;
 
  501    auto& forwardP = forwardMatrices.second;
 
  503    auto& backwardM = backwardMatrices.first;
 
  504    auto& backwardP = backwardMatrices.second;
 
  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;
 
  512    std::vector<bool> isVoidDomain;
 
  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;
 
  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)
 
  525        return std::make_pair(ForwardProjector(std::move(forwardM), std::move(forwardP)),
 
  526                              BackwardProjector(std::move(backwardM), std::move(backwardP)));
 
  528    else if (!hasVoidDomain && hasVoidTarget)
 
  530        std::vector<std::size_t> expansionMapTarget;
 
  531        ForwardProjectionMatrix forwardMReduced, forwardPReduced;
 
  533                             forwardMReduced, forwardPReduced, expansionMapTarget);
 
  535        return std::make_pair( ForwardProjector(std::move(forwardMReduced),
 
  536                                                std::move(forwardPReduced),
 
  537                                                std::move(expansionMapTarget),
 
  539                               BackwardProjector(std::move(backwardM), std::move(backwardP)) );
 
  541    else if (hasVoidDomain && !hasVoidTarget)
 
  545            std::vector<std::size_t> expansionMapDomain;
 
  546            BackwardProjectionMatrix backwardMReduced, backwardPReduced;
 
  548                                 backwardMReduced, backwardPReduced, expansionMapDomain);
 
  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),
 
  557            return std::make_pair( ForwardProjector(std::move(forwardM), std::move(forwardP)),
 
  558                                   BackwardProjector(std::move(backwardM), std::move(backwardP)) );
 
  562        std::vector<std::size_t> expansionMapTarget;
 
  563        ForwardProjectionMatrix forwardMReduced, forwardPReduced;
 
  565                             forwardMReduced, forwardPReduced, expansionMapTarget);
 
  569            std::vector<std::size_t> expansionMapDomain;
 
  570            BackwardProjectionMatrix backwardMReduced, backwardPReduced;
 
  572                                 backwardMReduced, backwardPReduced, expansionMapDomain);
 
  574            return std::make_pair( ForwardProjector(std::move(forwardMReduced),
 
  575                                                    std::move(forwardPReduced),
 
  576                                                    std::move(expansionMapTarget),
 
  578                                   BackwardProjector(std::move(backwardMReduced),
 
  579                                                     std::move(backwardPReduced),
 
  580                                                     std::move(expansionMapDomain),
 
  584            return std::make_pair( ForwardProjector(std::move(forwardMReduced),
 
  585                                                    std::move(forwardPReduced),
 
  586                                                    std::move(expansionMapTarget),
 
  588                                   BackwardProjector(std::move(backwardM), std::move(backwardP)) );