24#ifndef DUMUX_FLUX_OVER_SURFACE_STAGGERED_HH
25#define DUMUX_FLUX_OVER_SURFACE_STAGGERED_HH
32#include <dune/common/exceptions.hh>
33#include <dune/common/fvector.hh>
34#include <dune/geometry/type.hh>
35#include <dune/geometry/multilineargeometry.hh>
48template<
class Gr
idVariables,
class SolutionVector,
class ModelTraits,
class LocalRes
idual>
51 using Scalar =
typename GridVariables::Scalar;
52 using GridGeometry =
typename GridVariables::GridGeometry;
53 using FVElementGeometry =
typename GridGeometry::LocalView;
54 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
56 using GridView =
typename GridGeometry::GridView;
57 using VolumeVariables =
typename GridVariables::VolumeVariables;
58 using Element =
typename GridView::template Codim<0>::Entity;
60 using CellCenterPrimaryVariables = std::decay_t<decltype(std::declval<SolutionVector>()[GridGeometry::cellCenterIdx()][0])>;
64 dim = GridView::dimension,
65 dimWorld = GridView::dimensionworld
68 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
70 static constexpr auto surfaceDim = dimWorld - 1;
71 using SurfaceGeometryType = Dune::MultiLinearGeometry< Scalar, surfaceDim, dimWorld >;
76 template<
int mydim,
int coordDim,
class Scalar=
double>
79 using Geo = SurfaceGeometryType;
84 SurfaceData(Geo&& geo)
87 geometries_.push_back(std::move(geo));
90 SurfaceData(std::vector<Geo>&& geo)
92 values_.resize(geo.size());
93 geometries_ = std::forward<decltype(geo)>(geo);
96 void addValue(
int surfaceIdx,
const CellCenterPrimaryVariables& value)
98 values_[surfaceIdx] += value;
101 void addSubSurface(Geo&& geo)
103 values_.emplace_back(0.0);
104 geometries_.push_back(std::move(geo));
107 auto& subSurfaces()
const
112 CellCenterPrimaryVariables& value(
int surfaceIdx)
const
114 return values_[surfaceIdx];
122 void printSurfaceBoundaries(
int surfaceIdx)
const
124 const auto& geometry = geometries_[surfaceIdx];
125 for(
int i = 0; i < geometry.corners(); ++i)
126 std::cout <<
"(" << geometry.corner(i) <<
") ";
131 std::fill(values_.begin(), values_.end(), CellCenterPrimaryVariables(0.0));
136 std::vector<Geo> geometries_;
137 std::vector<CellCenterPrimaryVariables> values_;
150 : gridVariables_(gridVariables),
153 static_assert(std::is_same<Sol, SolutionVector>::value,
"Make sure that sol has the same type as SolutionVector."
154 "Use StaggeredVtkOutputModule<GridVariables, decltype(sol)> when calling the constructor.");
156 verbose_ = getParamFromGroup<bool>(problem_().paramGroup(),
"FluxOverSurface.Verbose",
false);
167 surfaces_[name] = SurfaceData<surfaceDim,dim>(std::forward<
decltype(surfaces)>(surfaces));
178 void addSurface(
const std::string& name,
const GlobalPosition& p0,
const GlobalPosition& p1)
180 surfaces_[name].addSubSurface(
makeSurface(std::vector<GlobalPosition>{p0, p1}));
194 const GlobalPosition& p0,
195 const GlobalPosition& p1,
196 const GlobalPosition& p2,
197 const GlobalPosition& p3)
199 surfaces_[name].addSubSurface(
makeSurface(std::vector<GlobalPosition>{p0, p1, p2, p3}));
207 static SurfaceGeometryType
makeSurface(
const std::vector<Dune::FieldVector<Scalar, 2>>& corners)
217 static SurfaceGeometryType
makeSurface(
const std::vector<Dune::FieldVector<Scalar, 3>>& corners)
227 auto fluxType = [
this](
const auto& element,
228 const auto& fvGeometry,
229 const auto& elemVolVars,
230 const auto& elemFaceVars,
232 const auto& elemFluxVarsCache)
234 LocalResidual localResidual(&problem_());
237 return localResidual.computeBoundaryFluxForCellCenter(problem_(), element, fvGeometry, scvf, elemVolVars, elemFaceVars, {}, elemFluxVarsCache);
239 return localResidual.computeFluxForCellCenter(problem_(), element, fvGeometry, elemVolVars, elemFaceVars, scvf, elemFluxVarsCache);
250 auto fluxType = [](
const auto& element,
251 const auto& fvGeometry,
252 const auto& elemVolVars,
253 const auto& elemFaceVars,
255 const auto& elemFluxVarsCache)
257 CellCenterPrimaryVariables result(0.0);
258 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
259 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
260 const Scalar extrusionFactor =
harmonicMean(insideVolVars.extrusionFactor(), outsideVolVars.extrusionFactor());
261 result[0] = elemFaceVars[scvf].velocitySelf() * Extrusion::area(fvGeometry, scvf) * extrusionFactor * scvf.directionSign();
280 template<
class FluxType>
284 for(
auto&&
surface : surfaces_)
288 std::vector<bool> dofVisited(problem_().gridGeometry().numFaceDofs(),
false);
290 auto fvGeometry =
localView(problem_().gridGeometry());
291 auto elemVolVars =
localView(gridVariables_.curGridVolVars());
292 auto elemFaceVars =
localView(gridVariables_.curGridFaceVars());
293 auto elemFluxVarsCache =
localView(gridVariables_.gridFluxVarsCache());
295 for(
auto&& element : elements(problem_().gridGeometry().gridView()))
297 fvGeometry.bind(element);
298 elemVolVars.bind(element, fvGeometry, sol_);
299 elemFaceVars.bind(element, fvGeometry, sol_);
300 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
302 for(
auto && scvf : scvfs(fvGeometry))
304 const auto dofIdx = scvf.dofIndex();
306 if(dofVisited[dofIdx])
309 dofVisited[dofIdx] =
true;
313 for(
auto&&
surface : surfaces_)
315 const auto& subSurfaces =
surface.second.subSurfaces();
317 for(
int surfaceIdx = 0; surfaceIdx < subSurfaces.size(); ++surfaceIdx)
321 const auto result = fluxType(element, fvGeometry, elemVolVars, elemFaceVars, scvf, elemFluxVarsCache);
323 surface.second.addValue(surfaceIdx, result);
327 std::cout <<
"Flux at face " << scvf.center() <<
" (" <<
surface.first <<
"): " << result;
328 std::cout <<
" (directionIndex: " << scvf.directionIndex() <<
"; surface boundaries: ";
329 surface.second.printSurfaceBoundaries(surfaceIdx);
330 std::cout <<
", surfaceIdx " << surfaceIdx <<
")" << std::endl;
344 auto&
values(
const std::string& name)
const
346 return surfaces_.at(name).values();
356 const auto& surfaceResults =
values(name);
357 return std::accumulate(surfaceResults.begin(), surfaceResults.end(), CellCenterPrimaryVariables(0.0));
362 const auto& problem_()
const {
return gridVariables_.curGridVolVars().problem(); }
364 std::map<std::string, SurfaceData<surfaceDim ,dim> > surfaces_;
365 const GridVariables& gridVariables_;
366 const SolutionVector& sol_;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Helper classes to compute the integration elements.
Detect if a point intersects a geometry.
Create Dune geometries from user-specified points.
auto makeDuneQuadrilaterial(const std::vector< Dune::FieldVector< CoordScalar, 3 > > &points)
Creates a dune quadrilateral geometry given 4 corner points.
Definition: makegeometry.hh:154
bool intersectsPointGeometry(const Dune::FieldVector< ctype, dimworld > &point, const Geometry &g)
Find out whether a point is inside a three-dimensional geometry.
Definition: intersectspointgeometry.hh:40
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:69
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
constexpr Line line
Definition: couplingmanager1d3d_line.hh:43
constexpr Surface surface
Definition: couplingmanager1d3d_surface.hh:49
Class used to calculate fluxes over surfaces. This only works for the staggered grid discretization.
Definition: fluxoversurface.hh:50
std::vector< SurfaceGeometryType > SurfaceList
Definition: fluxoversurface.hh:142
auto & values(const std::string &name) const
Return the fluxes of the individual sub surface of a given name.
Definition: fluxoversurface.hh:344
void addSurface(const std::string &name, const GlobalPosition &p0, const GlobalPosition &p1)
Add a surface under a given name, specifying the surface's corner points. This is a specialization fo...
Definition: fluxoversurface.hh:178
void calculateVolumeFluxes()
Calculate the volume fluxes over all surfaces.
Definition: fluxoversurface.hh:248
void addSurface(const std::string &name, SurfaceList &&surfaces)
Add a collection of sub surfaces under a given name.
Definition: fluxoversurface.hh:165
auto netFlux(const std::string &name) const
Return the cumulative net fluxes of a surface of a given name.
Definition: fluxoversurface.hh:354
static SurfaceGeometryType makeSurface(const std::vector< Dune::FieldVector< Scalar, 2 > > &corners)
Creates a geometrical surface object for (2D).
Definition: fluxoversurface.hh:207
void addSurface(const std::string &name, const GlobalPosition &p0, const GlobalPosition &p1, const GlobalPosition &p2, const GlobalPosition &p3)
Add a surface under a given name, specifying the surface's corner points. This is a specialization fo...
Definition: fluxoversurface.hh:193
void calculateMassOrMoleFluxes()
Calculate the mass or mole fluxes over all surfaces.
Definition: fluxoversurface.hh:225
FluxOverSurface(const GridVariables &gridVariables, const Sol &sol)
The constructor.
Definition: fluxoversurface.hh:148
static SurfaceGeometryType makeSurface(const std::vector< Dune::FieldVector< Scalar, 3 > > &corners)
Creates a geometrical surface object for (3D).
Definition: fluxoversurface.hh:217
void calculateFluxes(const FluxType &fluxType)
Calculate the fluxes over all surfaces for a given flux type.
Definition: fluxoversurface.hh:281