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