37    using Problem = std::decay_t<decltype(std::declval<LocalResidual>().problem())>;
 
   39    using GridView = 
typename GridGeometry::GridView;
 
   40    using Element = 
typename GridView::template Codim<0>::Entity;
 
   41    using FVElementGeometry = 
typename GridGeometry::LocalView;
 
   42    using SubControlVolume = 
typename FVElementGeometry::SubControlVolume;
 
   46    using NumEqVector = 
typename SolutionVector::block_type;
 
   47    static constexpr auto numEq = NumEqVector::dimension;
 
   52        NumEqVector totalFlux;
 
   53        std::unordered_map<std::size_t, NumEqVector> fluxPerPore;
 
   56        friend std::ostream& operator<< (std::ostream& stream, 
const Result& result)
 
   58            stream << result.totalFlux;
 
   63        const auto& operator[] (
int eqIdx)
 const 
   64        { 
return totalFlux[eqIdx]; }
 
   67        operator NumEqVector()
 const 
   73    using Scalar = 
typename GridVariables::Scalar;
 
   76                 const LocalResidual& localResidual,
 
   77                 const SolutionVector& sol)
 
   78    : localResidual_(localResidual)
 
   79    , gridVariables_(gridVariables)
 
   81    , isStationary_(localResidual.isStationary())
 
   83        const auto numDofs = localResidual_.problem().gridGeometry().numDofs();
 
   84        isConsidered_.resize(numDofs, 
false);
 
   85        boundaryFluxes_.resize(numDofs);
 
 
   95    Result 
getFlux(
const std::vector<Label>& labels, 
const bool verbose = 
false)
 const 
   98        auto restriction = [&] (
const SubControlVolume& scv)
 
  100            const Label poreLabel = localResidual_.problem().gridGeometry().poreLabel(scv.dofIndex());
 
  101            return std::any_of(labels.begin(), labels.end(),
 
  102                               [&](
const Label l){ return l == poreLabel; });
 
  105        std::fill(boundaryFluxes_.begin(), boundaryFluxes_.end(), NumEqVector(0.0));
 
  106        std::fill(isConsidered_.begin(), isConsidered_.end(), 
false);
 
  109        for (
const auto& element : elements(localResidual_.problem().gridGeometry().gridView()))
 
  110            computeBoundaryFlux_(element, restriction, verbose);
 
  113        result.totalFlux = std::accumulate(boundaryFluxes_.begin(), boundaryFluxes_.end(), NumEqVector(0.0));;
 
  114        for (
int i = 0; i < isConsidered_.size(); ++i)
 
  116            if (isConsidered_[i])
 
  117                result.fluxPerPore[i] = boundaryFluxes_[i];
 
 
  130    Result 
getFlux(std::string_view minMax, 
const int coord, 
const bool verbose = 
false)
 const 
  132        if (!(minMax == 
"min" || minMax == 
"max"))
 
  133            DUNE_THROW(Dune::InvalidStateException,
 
  134                    "second argument must be either 'min' or 'max' (string) !");
 
  137        auto onMeasuringBoundary = [&] (
const Scalar pos)
 
  139            return ( (minMax == 
"min" && pos < localResidual_.problem().gridGeometry().bBoxMin()[coord] + eps) ||
 
  140                     (minMax == 
"max" && pos > localResidual_.problem().gridGeometry().bBoxMax()[coord] - eps) );
 
  144        auto restriction = [&] (
const SubControlVolume& scv)
 
  146            bool considerAllDirections = coord < 0 ? true : 
false;
 
  149            bool considerScv = localResidual_.problem().gridGeometry().dofOnBoundary(scv.dofIndex()) && onMeasuringBoundary(scv.dofPosition()[coord]);
 
  153            if (considerScv && !considerAllDirections)
 
  155                const auto& pos = scv.dofPosition();
 
  156                if (!(pos[coord] < localResidual_.problem().gridGeometry().bBoxMin()[coord] + eps || pos[coord] > localResidual_.problem().gridGeometry().bBoxMax()[coord] -eps ))
 
  163        std::fill(boundaryFluxes_.begin(), boundaryFluxes_.end(), NumEqVector(0.0));
 
  164        std::fill(isConsidered_.begin(), isConsidered_.end(), 
false);
 
  167        for (
const auto& element : elements(localResidual_.problem().gridGeometry().gridView()))
 
  168            computeBoundaryFlux_(element, restriction, verbose);
 
  171        result.totalFlux = std::accumulate(boundaryFluxes_.begin(), boundaryFluxes_.end(), NumEqVector(0.0));;
 
  172        for (
int i = 0; i < isConsidered_.size(); ++i)
 
  174            if (isConsidered_[i])
 
  175                result.fluxPerPore[i] = boundaryFluxes_[i];
 
 
  190    template<
class RestrictingFunction>
 
  191    void computeBoundaryFlux_(
const Element& element,
 
  192                              RestrictingFunction considerScv,
 
  193                              const bool verbose = 
false)
 const 
  198        const auto fvGeometry = 
localView(localResidual_.problem().gridGeometry()).bind(element);
 
  199        const auto curElemVolVars = 
localView(gridVariables_.curGridVolVars()).bind(element, fvGeometry, sol_);
 
  201        auto prevElemVolVars = 
localView(gridVariables_.prevGridVolVars());
 
  203            prevElemVolVars.bindElement(element, fvGeometry, sol_);
 
  205        ElementBoundaryTypes elemBcTypes;
 
  206        elemBcTypes.
update(localResidual_.problem(), element, fvGeometry);
 
  207        const auto elemFluxVarsCache = 
localView(gridVariables_.gridFluxVarsCache()).bindElement(element, fvGeometry, curElemVolVars);
 
  208        auto residual = localResidual_.evalFluxAndSource(element, fvGeometry, curElemVolVars, elemFluxVarsCache, elemBcTypes);
 
  211            residual += localResidual_.evalStorage(element, fvGeometry, prevElemVolVars, curElemVolVars);
 
  213        for (
auto&& scv : scvs(fvGeometry))
 
  216            if (considerScv(scv))
 
  218                isConsidered_[scv.dofIndex()] = 
true;
 
  221                const auto& bcTypes = elemBcTypes.get(fvGeometry, scv);
 
  224                for (
int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx)
 
  230                    if (!bcTypes.isDirichlet(eqIdx))
 
  232                        auto source = localResidual_.computeSource(localResidual_.problem(), element, fvGeometry, curElemVolVars, scv);
 
  233                        source *= scv.volume() * curElemVolVars[scv].extrusionFactor();
 
  234                        flux[eqIdx] = source[eqIdx];
 
  237                        flux[eqIdx] = residual[scv.indexInElement()][eqIdx];
 
  243                boundaryFluxes_[scv.dofIndex()] -= flux;
 
  246                    std::cout << 
"SCV of element " << scv.elementIndex()  << 
" at vertex " << scv.dofIndex() << 
" has flux: " << flux << std::endl;
 
  251    const LocalResidual localResidual_; 
 
  252    const GridVariables& gridVariables_;
 
  253    const SolutionVector& sol_;
 
  255    mutable std::vector<bool> isConsidered_;
 
  256    mutable std::vector<NumEqVector> boundaryFluxes_;
 
 
void update(const Problem &problem, const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry)
Update the boundary types for all vertices of an element.
Definition cvfe/elementboundarytypes.hh:41
Result getFlux(std::string_view minMax, const int coord, const bool verbose=false) const
Returns the cumulative flux in ,  or  of several pores at a given location on the boundary.
Definition boundaryflux.hh:130
Result getFlux(const std::vector< Label > &labels, const bool verbose=false) const
Returns the cumulative flux in ,  or  of several pores for a given list of pore labels to consider.
Definition boundaryflux.hh:95
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition numeqvector.hh:34
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition localview.hh:26