24#ifndef DUMUX_PNM_BOUNDARYFLUX_HH
25#define DUMUX_PNM_BOUNDARYFLUX_HH
31#include <unordered_map>
35#include <dune/common/exceptions.hh>
46template<
class Gr
idVariables,
class LocalRes
idual,
class SolutionVector>
49 using Problem = std::decay_t<decltype(std::declval<LocalResidual>().problem())>;
51 using GridView =
typename GridGeometry::GridView;
52 using Element =
typename GridView::template Codim<0>::Entity;
53 using FVElementGeometry =
typename GridGeometry::LocalView;
54 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
58 using NumEqVector =
typename SolutionVector::block_type;
59 static constexpr auto numEq = NumEqVector::dimension;
64 NumEqVector totalFlux;
65 std::unordered_map<std::size_t, NumEqVector> fluxPerPore;
68 friend std::ostream&
operator<< (std::ostream& stream,
const Result& result)
70 stream << result.totalFlux;
75 const auto& operator[] (
int eqIdx)
const
76 {
return totalFlux[eqIdx]; }
79 operator NumEqVector()
const
85 using Scalar =
typename GridVariables::Scalar;
88 const LocalResidual& localResidual,
89 const SolutionVector& sol)
90 : localResidual_(localResidual)
91 , gridVariables_(gridVariables)
93 , isStationary_(localResidual.isStationary())
95 const auto numDofs = localResidual_.problem().gridGeometry().numDofs();
96 isConsidered_.resize(numDofs,
false);
97 boundaryFluxes_.resize(numDofs);
106 template<
class Label>
107 Result
getFlux(
const std::vector<Label>& labels,
const bool verbose =
false)
const
110 auto restriction = [&] (
const SubControlVolume& scv)
112 const Label poreLabel = localResidual_.problem().gridGeometry().poreLabel(scv.dofIndex());
113 return std::any_of(labels.begin(), labels.end(),
114 [&](
const Label l){ return l == poreLabel; });
117 std::fill(boundaryFluxes_.begin(), boundaryFluxes_.end(), NumEqVector(0.0));
118 std::fill(isConsidered_.begin(), isConsidered_.end(),
false);
121 for (
const auto& element : elements(localResidual_.problem().gridGeometry().gridView()))
122 computeBoundaryFlux_(element, restriction, verbose);
125 result.totalFlux = std::accumulate(boundaryFluxes_.begin(), boundaryFluxes_.end(), NumEqVector(0.0));;
126 for (
int i = 0; i < isConsidered_.size(); ++i)
128 if (isConsidered_[i])
129 result.fluxPerPore[i] = boundaryFluxes_[i];
142 Result
getFlux(std::string_view minMax,
const int coord,
const bool verbose =
false)
const
144 if (!(minMax ==
"min" || minMax ==
"max"))
145 DUNE_THROW(Dune::InvalidStateException,
146 "second argument must be either 'min' or 'max' (string) !");
149 auto onMeasuringBoundary = [&] (
const Scalar pos)
151 return ( (minMax ==
"min" && pos < localResidual_.problem().gridGeometry().bBoxMin()[coord] + eps) ||
152 (minMax ==
"max" && pos > localResidual_.problem().gridGeometry().bBoxMax()[coord] - eps) );
156 auto restriction = [&] (
const SubControlVolume& scv)
158 bool considerAllDirections = coord < 0 ? true :
false;
161 bool considerScv = localResidual_.problem().gridGeometry().dofOnBoundary(scv.dofIndex()) && onMeasuringBoundary(scv.dofPosition()[coord]);
165 if (considerScv && !considerAllDirections)
167 const auto& pos = scv.dofPosition();
168 if (!(pos[coord] < localResidual_.problem().gridGeometry().bBoxMin()[coord] + eps || pos[coord] > localResidual_.problem().gridGeometry().bBoxMax()[coord] -eps ))
175 std::fill(boundaryFluxes_.begin(), boundaryFluxes_.end(), NumEqVector(0.0));
176 std::fill(isConsidered_.begin(), isConsidered_.end(),
false);
179 for (
const auto& element : elements(localResidual_.problem().gridGeometry().gridView()))
180 computeBoundaryFlux_(element, restriction, verbose);
183 result.totalFlux = std::accumulate(boundaryFluxes_.begin(), boundaryFluxes_.end(), NumEqVector(0.0));;
184 for (
int i = 0; i < isConsidered_.size(); ++i)
186 if (isConsidered_[i])
187 result.fluxPerPore[i] = boundaryFluxes_[i];
202 template<
class RestrictingFunction>
203 void computeBoundaryFlux_(
const Element& element,
204 RestrictingFunction considerScv,
205 const bool verbose =
false)
const
210 const auto fvGeometry =
localView(localResidual_.problem().gridGeometry()).bind(element);
211 const auto curElemVolVars =
localView(gridVariables_.curGridVolVars()).bind(element, fvGeometry, sol_);
213 auto prevElemVolVars =
localView(gridVariables_.prevGridVolVars());
215 prevElemVolVars.bindElement(element, fvGeometry, sol_);
217 ElementBoundaryTypes elemBcTypes;
218 elemBcTypes.update(localResidual_.problem(), element, fvGeometry);
219 const auto elemFluxVarsCache =
localView(gridVariables_.gridFluxVarsCache()).bindElement(element, fvGeometry, curElemVolVars);
220 auto residual = localResidual_.evalFluxAndSource(element, fvGeometry, curElemVolVars, elemFluxVarsCache, elemBcTypes);
223 residual += localResidual_.evalStorage(element, fvGeometry, prevElemVolVars, curElemVolVars);
225 for (
auto&& scv : scvs(fvGeometry))
228 if (considerScv(scv))
230 isConsidered_[scv.dofIndex()] =
true;
233 const auto& bcTypes = elemBcTypes.get(fvGeometry, scv);
236 for (
int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx)
242 if (!bcTypes.isDirichlet(eqIdx))
244 auto source = localResidual_.computeSource(localResidual_.problem(), element, fvGeometry, curElemVolVars, scv);
245 source *= scv.volume() * curElemVolVars[scv].extrusionFactor();
246 flux[eqIdx] = source[eqIdx];
249 flux[eqIdx] = residual[scv.indexInElement()][eqIdx];
255 boundaryFluxes_[scv.dofIndex()] -= flux;
258 std::cout <<
"SCV of element " << scv.elementIndex() <<
" at vertex " << scv.dofIndex() <<
" has flux: " << flux << std::endl;
263 const LocalResidual localResidual_;
264 const GridVariables& gridVariables_;
265 const SolutionVector& sol_;
267 mutable std::vector<bool> isConsidered_;
268 mutable std::vector<NumEqVector> boundaryFluxes_;
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
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:46
auto operator<<(std::ostream &os, const T &t) -> std::enable_if_t< decltype(Detail::hasName(t))::value, std::ostream & >
Return the class name of the tagged type calling t.name()
Definition: tag.hh:59
Definition: discretization/porenetwork/fvelementgeometry.hh:34
std::decay_t< decltype(std::declval< Problem >().gridGeometry())> GridGeometry
Definition: common/typetraits/problem.hh:45
typename Detail::template ProblemTraits< Problem, typename GridGeometry::DiscretizationMethod >::BoundaryTypes BoundaryTypes
Definition: common/typetraits/problem.hh:46
This class stores an array of BoundaryTypes objects.
Definition: cvfe/elementboundarytypes.hh:40
Class for the calculation of fluxes at the boundary of pore-network models.
Definition: boundaryflux.hh:48
BoundaryFlux(const GridVariables &gridVariables, const LocalResidual &localResidual, const SolutionVector &sol)
Definition: boundaryflux.hh:87
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:142
typename GridVariables::Scalar Scalar
Definition: boundaryflux.hh:85
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:107
Type traits for problem classes.
Boundary types gathered on an element.