version 3.8
multidomain/facet/box/darcyslaw.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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_BOX_FACET_COUPLING_DARCYS_LAW_HH
13#define DUMUX_DISCRETIZATION_BOX_FACET_COUPLING_DARCYS_LAW_HH
14
15#include <vector>
16#include <cmath>
17
18#include <dune/common/exceptions.hh>
19#include <dune/common/fvector.hh>
20#include <dune/common/float_cmp.hh>
21
24
28
29namespace Dumux {
30
37template<class Scalar, class GridGeometry>
39{
41
42 using FVElementGeometry = typename GridGeometry::LocalView;
43 using SubControlVolume = typename GridGeometry::SubControlVolume;
44 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
45 using Extrusion = Extrusion_t<GridGeometry>;
46 using GridView = typename GridGeometry::GridView;
47 using Element = typename GridView::template Codim<0>::Entity;
48 using CoordScalar = typename GridView::ctype;
49 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
50
51 static constexpr int dim = GridView::dimension;
52 static constexpr int dimWorld = GridView::dimensionworld;
53
54public:
55
56 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
57 static Scalar flux(const Problem& problem,
58 const Element& element,
59 const FVElementGeometry& fvGeometry,
60 const ElementVolumeVariables& elemVolVars,
61 const SubControlVolumeFace& scvf,
62 const int phaseIdx,
63 const ElementFluxVarsCache& elemFluxVarCache)
64 {
65 // if this scvf is not on an interior boundary, use the standard law
66 if (!scvf.interiorBoundary())
67 return DefaultDarcysLaw::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarCache);
68
69 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
70 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
71 DUNE_THROW(Dune::NotImplemented, "Xi != 1.0 cannot be used with the Box-Facet-Coupling scheme");
72
73 // get some references for convenience
74 const auto& fluxVarCache = elemFluxVarCache[scvf];
75 const auto& shapeValues = fluxVarCache.shapeValues();
76 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
77 const auto& insideVolVars = elemVolVars[insideScv];
78
79 // evaluate user-defined interior boundary types
80 const auto bcTypes = problem.interiorBoundaryTypes(element, scvf);
81
82 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
83
84 // on interior Neumann boundaries, evaluate the flux using the facet permeability
85 if (bcTypes.hasOnlyNeumann())
86 {
87 // interpolate pressure/density to scvf integration point
88 Scalar p = 0.0;
89 Scalar rho = 0.0;
90 for (const auto& scv : scvs(fvGeometry))
91 {
92 const auto& volVars = elemVolVars[scv];
93 p += volVars.pressure(phaseIdx)*shapeValues[scv.indexInElement()][0];
94 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
95 }
96
97 // compute tpfa flux from integration point to facet centerline
98 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
99
100 using std::sqrt;
101 // If this is a surface grid, use the square root of the facet extrusion factor
102 // as an approximate average distance from scvf ip to facet center
103 using std::sqrt;
104 const auto a = facetVolVars.extrusionFactor();
105 auto gradP = scvf.unitOuterNormal();
106 gradP *= dim == dimWorld ? 0.5*a : 0.5*sqrt(a);
107 gradP /= gradP.two_norm2();
108 gradP *= (facetVolVars.pressure(phaseIdx) - p);
109 if (enableGravity)
110 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
111
112 // apply facet permeability and return the flux
113 return -1.0*Extrusion::area(fvGeometry, scvf)
114 *insideVolVars.extrusionFactor()
115 *vtmv(scvf.unitOuterNormal(), facetVolVars.permeability(), gradP);
116 }
117
118 // on interior Dirichlet boundaries use the facet pressure and evaluate flux
119 else if (bcTypes.hasOnlyDirichlet())
120 {
121 // create vector with nodal pressures
122 std::vector<Scalar> pressures(element.subEntities(dim));
123 for (const auto& scv : scvs(fvGeometry))
124 pressures[scv.localDofIndex()] = elemVolVars[scv].pressure(phaseIdx);
125
126 // substitute with facet pressures for those scvs touching this facet
127 for (const auto& scvfJ : scvfs(fvGeometry))
128 if (scvfJ.interiorBoundary() && scvfJ.facetIndexInElement() == scvf.facetIndexInElement())
129 pressures[ fvGeometry.scv(scvfJ.insideScvIdx()).localDofIndex() ]
130 = problem.couplingManager().getLowDimVolVars(element, scvfJ).pressure(phaseIdx);
131
132 // evaluate gradP - rho*g at integration point
133 Scalar rho(0.0);
134 Dune::FieldVector<Scalar, dimWorld> gradP(0.0);
135 for (const auto& scv : scvs(fvGeometry))
136 {
137 rho += elemVolVars[scv].density(phaseIdx)*shapeValues[scv.indexInElement()][0];
138 gradP.axpy(pressures[scv.localDofIndex()], fluxVarCache.gradN(scv.indexInElement()));
139 }
140
141 if (enableGravity)
142 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
143
144 // apply matrix permeability and return the flux
145 return -1.0*Extrusion::area(fvGeometry, scvf)
146 *insideVolVars.extrusionFactor()
147 *vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), gradP);
148 }
149
150 // mixed boundary types are not supported
151 else
152 DUNE_THROW(Dune::NotImplemented, "Mixed boundary types are not supported");
153 }
154
155 // compute transmissibilities ti for analytical jacobians
156 template<class Problem, class ElementVolumeVariables, class FluxVarCache>
157 static std::vector<Scalar> calculateTransmissibilities(const Problem& problem,
158 const Element& element,
159 const FVElementGeometry& fvGeometry,
160 const ElementVolumeVariables& elemVolVars,
161 const SubControlVolumeFace& scvf,
162 const FluxVarCache& fluxVarCache)
163 {
164 DUNE_THROW(Dune::NotImplemented, "transmissibilty computation for BoxFacetCouplingDarcysLaw");
165 }
166};
167
168} // end namespace Dumux
169
170#endif
Darcy's law for the box scheme in the context of coupled models where coupling occurs across the face...
Definition: multidomain/facet/box/darcyslaw.hh:39
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVarsCache &elemFluxVarCache)
Definition: multidomain/facet/box/darcyslaw.hh:57
static std::vector< Scalar > calculateTransmissibilities(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const FluxVarCache &fluxVarCache)
Definition: multidomain/facet/box/darcyslaw.hh:157
Darcy's law for control-volume finite element schemes.
Definition: flux/cvfe/darcyslaw.hh:36
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVarsCache &elemFluxVarCache)
Returns the advective flux of a fluid phase across the given sub-control volume face.
Definition: flux/cvfe/darcyslaw.hh:59
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
Specialization of Darcy's Law for control-volume finite element schemes.
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:851
The available discretization methods in Dumux.
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.