3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
boundaryflux.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_PNM_BOUNDARYFLUX_HH
25#define DUMUX_PNM_BOUNDARYFLUX_HH
26
27#include <algorithm>
28#include <numeric>
29#include <vector>
30#include <type_traits>
31#include <unordered_map>
32#include <string_view>
33#include <iostream>
34
35#include <dune/common/exceptions.hh>
38
39namespace Dumux::PoreNetwork {
40
46template<class GridVariables, class LocalResidual, class SolutionVector>
48{
49 using Problem = std::decay_t<decltype(std::declval<LocalResidual>().problem())>;
50 using GridGeometry = typename ProblemTraits<Problem>::GridGeometry;
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;
55 using BoundaryTypes = typename ProblemTraits<Problem>::BoundaryTypes;
57
58 using NumEqVector = typename SolutionVector::block_type;
59 static constexpr auto numEq = NumEqVector::dimension;
60
62 struct Result
63 {
64 NumEqVector totalFlux;
65 std::unordered_map<std::size_t, NumEqVector> fluxPerPore;
66
68 friend std::ostream& operator<< (std::ostream& stream, const Result& result)
69 {
70 stream << result.totalFlux;
71 return stream;
72 }
73
75 const auto& operator[] (int eqIdx) const
76 { return totalFlux[eqIdx]; }
77
79 operator NumEqVector() const
80 { return totalFlux; }
81 };
82
83public:
84 // export the Scalar type
85 using Scalar = typename GridVariables::Scalar;
86
87 BoundaryFlux(const GridVariables& gridVariables,
88 const LocalResidual& localResidual,
89 const SolutionVector& sol)
90 : localResidual_(localResidual)
91 , gridVariables_(gridVariables)
92 , sol_(sol)
93 , isStationary_(localResidual.isStationary())
94 {
95 const auto numDofs = localResidual_.problem().gridGeometry().numDofs();
96 isConsidered_.resize(numDofs, false);
97 boundaryFluxes_.resize(numDofs);
98 }
99
106 template<class Label>
107 Result getFlux(const std::vector<Label>& labels, const bool verbose = false) const
108 {
109 // helper lambda to decide which scvs to consider for flux calculation
110 auto restriction = [&] (const SubControlVolume& scv)
111 {
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; });
115 };
116
117 std::fill(boundaryFluxes_.begin(), boundaryFluxes_.end(), NumEqVector(0.0));
118 std::fill(isConsidered_.begin(), isConsidered_.end(), false);
119
120 // sum up the fluxes
121 for (const auto& element : elements(localResidual_.problem().gridGeometry().gridView()))
122 getFlux(element, restriction, verbose);
123
124 Result result;
125 result.totalFlux = std::accumulate(boundaryFluxes_.begin(), boundaryFluxes_.end(), NumEqVector(0.0));;
126 for (int i = 0; i < isConsidered_.size(); ++i)
127 {
128 if (isConsidered_[i])
129 result.fluxPerPore[i] = boundaryFluxes_[i];
130 }
131
132 return result;
133 }
134
142 Result getFlux(std::string_view minMax, const int coord, const bool verbose = false) const
143 {
144 if (!(minMax == "min" || minMax == "max"))
145 DUNE_THROW(Dune::InvalidStateException,
146 "second argument must be either 'min' or 'max' (string) !");
147
148 const Scalar eps = 1e-6; //TODO
149 auto onMeasuringBoundary = [&] (const Scalar pos)
150 {
151 return ( (minMax == "min" && pos < localResidual_.problem().gridGeometry().bBoxMin()[coord] + eps) ||
152 (minMax == "max" && pos > localResidual_.problem().gridGeometry().bBoxMax()[coord] - eps) );
153 };
154
155 // helper lambda to decide which scvs to consider for flux calculation
156 auto restriction = [&] (const SubControlVolume& scv)
157 {
158 bool considerAllDirections = coord < 0 ? true : false;
159
160 //only consider SCVs on the boundary
161 bool considerScv = localResidual_.problem().gridGeometry().dofOnBoundary(scv.dofIndex()) && onMeasuringBoundary(scv.dofPosition()[coord]);
162
163 //check whether a vertex lies on a boundary and also check whether this boundary shall be
164 // considered for the flux calculation
165 if (considerScv && !considerAllDirections)
166 {
167 const auto& pos = scv.dofPosition();
168 if (!(pos[coord] < localResidual_.problem().gridGeometry().bBoxMin()[coord] + eps || pos[coord] > localResidual_.problem().gridGeometry().bBoxMax()[coord] -eps ))
169 considerScv = false;
170 }
171
172 return considerScv;
173 };
174
175 std::fill(boundaryFluxes_.begin(), boundaryFluxes_.end(), NumEqVector(0.0));
176 std::fill(isConsidered_.begin(), isConsidered_.end(), false);
177
178 // sum up the fluxes
179 for (const auto& element : elements(localResidual_.problem().gridGeometry().gridView()))
180 getFlux(element, restriction, verbose);
181
182 Result result;
183 result.totalFlux = std::accumulate(boundaryFluxes_.begin(), boundaryFluxes_.end(), NumEqVector(0.0));;
184 for (int i = 0; i < isConsidered_.size(); ++i)
185 {
186 if (isConsidered_[i])
187 result.fluxPerPore[i] = boundaryFluxes_[i];
188 }
189
190 return result;
191
192 }
193
201 template<class RestrictingFunction>
202 NumEqVector getFlux(const Element& element,
203 RestrictingFunction considerScv,
204 const bool verbose = false) const
205 {
206 NumEqVector flux(0.0);
207
208 // by default, all coordinate directions are considered for the definition of a boundary
209
210 // make sure FVElementGeometry and volume variables are bound to the element
211 auto fvGeometry = localView(localResidual_.problem().gridGeometry());
212 fvGeometry.bind(element);
213
214 auto curElemVolVars = localView(gridVariables_.curGridVolVars());
215 curElemVolVars.bind(element, fvGeometry, sol_);
216
217 auto prevElemVolVars = localView(gridVariables_.prevGridVolVars());
218 if (!isStationary_)
219 prevElemVolVars.bindElement(element, fvGeometry, sol_);
220
221 auto elemFluxVarsCache = localView(gridVariables_.gridFluxVarsCache());
222 elemFluxVarsCache.bindElement(element, fvGeometry, curElemVolVars);
223
224 ElementBoundaryTypes elemBcTypes;
225 elemBcTypes.update(localResidual_.problem(), element, fvGeometry);
226
227 auto residual = localResidual_.evalFluxAndSource(element, fvGeometry, curElemVolVars, elemFluxVarsCache, elemBcTypes);
228
229 if (!isStationary_)
230 residual += localResidual_.evalStorage(element, fvGeometry, prevElemVolVars, curElemVolVars);
231
232 for (auto&& scv : scvs(fvGeometry))
233 {
234 // compute the boundary flux using the local residual of the element's scv on the boundary
235 if (considerScv(scv))
236 {
237 isConsidered_[scv.dofIndex()] = true;
238
239 // The flux must be substracted:
240 // On an inlet boundary, the flux part of the local residual will be positive, since all fluxes will leave the SCV towards to interior domain.
241 // For the domain itself, however, the sign has to be negative, since mass is entering the system.
242 flux -= residual[scv.indexInElement()];
243
244 boundaryFluxes_[scv.dofIndex()] -= residual[scv.indexInElement()];
245
246 if (verbose)
247 std::cout << "SCV of element " << scv.elementIndex() << " at vertex " << scv.dofIndex() << " has flux: " << residual[scv.indexInElement()] << std::endl;
248 }
249 }
250 return flux;
251 }
252
253private:
254 const LocalResidual localResidual_; // store a copy of the local residual
255 const GridVariables& gridVariables_;
256 const SolutionVector& sol_;
257 bool isStationary_;
258 mutable std::vector<bool> isConsidered_;
259 mutable std::vector<NumEqVector> boundaryFluxes_;
260};
261
262} // end Dumux::PoreNetwork
263
264#endif
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
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:33
std::decay_t< decltype(std::declval< Problem >().gridGeometry())> GridGeometry
Definition: common/typetraits/problem.hh:45
typename Detail::template ProblemTraits< Problem, GridGeometry::discMethod >::BoundaryTypes BoundaryTypes
Definition: common/typetraits/problem.hh:46
This class stores an array of BoundaryTypes objects.
Definition: box/elementboundarytypes.hh:39
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry)
Update the boundary types for all vertices of an element.
Definition: box/elementboundarytypes.hh:52
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
NumEqVector getFlux(const Element &element, RestrictingFunction considerScv, const bool verbose=false) const
Returns the cumulative flux in , or of several pore throats at a given location on the boundary.
Definition: boundaryflux.hh:202
Result getFlux(std::string_view minMax, const int coord, const bool verbose=false) const
Returns the cumulative flux in , or of several pore throats 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 pore throats for a given list of pore labels to cons...
Definition: boundaryflux.hh:107
Type traits for problem classes.
Boundary types gathered on an element.