3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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_DISCRETIZATION_BOX_FACET_COUPLING_DARCYS_LAW_HH
25#define DUMUX_DISCRETIZATION_BOX_FACET_COUPLING_DARCYS_LAW_HH
26
27#include <dune/common/exceptions.hh>
28#include <dune/common/fvector.hh>
29#include <dune/common/float_cmp.hh>
30#include <dune/geometry/multilineargeometry.hh>
31
32#include <dumux/common/math.hh>
37
40
41namespace Dumux {
42
49template<class Scalar, class GridGeometry>
51{
53
54 using FVElementGeometry = typename GridGeometry::LocalView;
55 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
56 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
57 using GridView = typename GridGeometry::GridView;
58 using Element = typename GridView::template Codim<0>::Entity;
59 using CoordScalar = typename GridView::ctype;
60 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
61
62 static constexpr int dim = GridView::dimension;
63 static constexpr int dimWorld = GridView::dimensionworld;
64
65public:
66
67 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
68 static Scalar flux(const Problem& problem,
69 const Element& element,
70 const FVElementGeometry& fvGeometry,
71 const ElementVolumeVariables& elemVolVars,
72 const SubControlVolumeFace& scvf,
73 const int phaseIdx,
74 const ElementFluxVarsCache& elemFluxVarCache)
75 {
76 // if this scvf is not on an interior boundary, use the standard law
77 if (!scvf.interiorBoundary())
78 return DefaultBoxDarcysLaw::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarCache);
79
80 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
81 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
82 DUNE_THROW(Dune::NotImplemented, "Xi != 1.0 cannot be used with the Box-Facet-Coupling scheme");
83
84 // get some references for convenience
85 const auto& fluxVarCache = elemFluxVarCache[scvf];
86 const auto& shapeValues = fluxVarCache.shapeValues();
87 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
88 const auto& insideVolVars = elemVolVars[insideScv];
89
90 // evaluate user-defined interior boundary types
91 const auto bcTypes = problem.interiorBoundaryTypes(element, scvf);
92
93 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
94
95 // on interior Neumann boundaries, evaluate the flux using the facet permeability
96 if (bcTypes.hasOnlyNeumann())
97 {
98 // interpolate pressure/density to scvf integration point
99 Scalar p = 0.0;
100 Scalar rho = 0.0;
101 for (const auto& scv : scvs(fvGeometry))
102 {
103 const auto& volVars = elemVolVars[scv];
104 p += volVars.pressure(phaseIdx)*shapeValues[scv.indexInElement()][0];
105 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
106 }
107
108 // compute tpfa flux from integration point to facet centerline
109 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
110
111 using std::sqrt;
112 // If this is a surface grid, use the square root of the facet extrusion factor
113 // as an approximate average distance from scvf ip to facet center
114 using std::sqrt;
115 const auto a = facetVolVars.extrusionFactor();
116 auto gradP = scvf.unitOuterNormal();
117 gradP *= dim == dimWorld ? 0.5*a : 0.5*sqrt(a);
118 gradP /= gradP.two_norm2();
119 gradP *= (facetVolVars.pressure(phaseIdx) - p);
120 if (enableGravity)
121 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
122
123 // apply facet permeability and return the flux
124 return -1.0*scvf.area()
125 *insideVolVars.extrusionFactor()
126 *vtmv(scvf.unitOuterNormal(), facetVolVars.permeability(), gradP);
127 }
128
129 // on interior Dirichlet boundaries use the facet pressure and evaluate flux
130 else if (bcTypes.hasOnlyDirichlet())
131 {
132 // create vector with nodal pressures
133 std::vector<Scalar> pressures(element.subEntities(dim));
134 for (const auto& scv : scvs(fvGeometry))
135 pressures[scv.localDofIndex()] = elemVolVars[scv].pressure(phaseIdx);
136
137 // substitute with facet pressures for those scvs touching this facet
138 for (const auto& scvfJ : scvfs(fvGeometry))
139 if (scvfJ.interiorBoundary() && scvfJ.facetIndexInElement() == scvf.facetIndexInElement())
140 pressures[ fvGeometry.scv(scvfJ.insideScvIdx()).localDofIndex() ]
141 = problem.couplingManager().getLowDimVolVars(element, scvfJ).pressure(phaseIdx);
142
143 // evaluate gradP - rho*g at integration point
144 Scalar rho(0.0);
145 Dune::FieldVector<Scalar, dimWorld> gradP(0.0);
146 for (const auto& scv : scvs(fvGeometry))
147 {
148 rho += elemVolVars[scv].density(phaseIdx)*shapeValues[scv.indexInElement()][0];
149 gradP.axpy(pressures[scv.localDofIndex()], fluxVarCache.gradN(scv.indexInElement()));
150 }
151
152 if (enableGravity)
153 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
154
155 // apply matrix permeability and return the flux
156 return -1.0*scvf.area()
157 *insideVolVars.extrusionFactor()
158 *vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), gradP);
159 }
160
161 // mixed boundary types are not supported
162 else
163 DUNE_THROW(Dune::NotImplemented, "Mixed boundary types are not supported");
164 }
165
166 // compute transmissibilities ti for analytical jacobians
167 template<class Problem, class ElementVolumeVariables, class FluxVarCache>
168 static std::vector<Scalar> calculateTransmissibilities(const Problem& problem,
169 const Element& element,
170 const FVElementGeometry& fvGeometry,
171 const ElementVolumeVariables& elemVolVars,
172 const SubControlVolumeFace& scvf,
173 const FluxVarCache& fluxVarCache)
174 {
175 DUNE_THROW(Dune::NotImplemented, "transmissibilty computation for BoxFacetCouplingDarcysLaw");
176 }
177};
178
179} // end namespace Dumux
180
181#endif
A function to compute a geometry's diameter, i.e. the longest distance between points of a geometry.
A class for collision detection of two geometries and computation of intersection corners.
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
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:840
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Darcy's law for box schemes.
Definition: flux/box/darcyslaw.hh:63
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: flux/box/darcyslaw.hh:79
Darcy's law for the box scheme scheme in the context of coupled models where coupling occurs across t...
Definition: multidomain/facet/box/darcyslaw.hh:51
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:68
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:168
Declares all properties used in Dumux.
Specialization of Darcy's Law for the box method.