3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
multidomain/facet/box/fourierslaw.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_FOURIERS_LAW_HH
25#define DUMUX_DISCRETIZATION_BOX_FACET_COUPLING_FOURIERS_LAW_HH
26
27#include <cmath>
28#include <vector>
29
30#include <dune/common/exceptions.hh>
31#include <dune/common/fvector.hh>
32#include <dune/common/float_cmp.hh>
33
34#include <dumux/common/math.hh>
37
41
42namespace Dumux {
43
50template<class TypeTag>
52: public FouriersLawImplementation<TypeTag, DiscretizationMethod::box>
53{
55
58 using FVElementGeometry = typename GridGeometry::LocalView;
59 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
60 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
61 using GridView = typename GridGeometry::GridView;
62 using Element = typename GridView::template Codim<0>::Entity;
63 using CoordScalar = typename GridView::ctype;
64 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
65
66 static constexpr int dim = GridView::dimension;
67 static constexpr int dimWorld = GridView::dimensionworld;
68
69public:
70
74 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
75 static Scalar flux(const Problem& problem,
76 const Element& element,
77 const FVElementGeometry& fvGeometry,
78 const ElementVolumeVariables& elemVolVars,
79 const SubControlVolumeFace& scvf,
80 const ElementFluxVarsCache& elemFluxVarCache)
81 {
82 // if this scvf is not on an interior boundary, use the standard law
83 if (!scvf.interiorBoundary())
84 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarCache);
85
86 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
87 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
88 DUNE_THROW(Dune::NotImplemented, "Xi != 1.0 cannot be used with the Box-Facet-Coupling scheme");
89
90 // get some references for convenience
91 const auto& fluxVarCache = elemFluxVarCache[scvf];
92 const auto& shapeValues = fluxVarCache.shapeValues();
93 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
94 const auto& insideVolVars = elemVolVars[insideScv];
95
96 // on interior Neumann boundaries, evaluate the flux using the facet thermal conductivity
97 const auto bcTypes = problem.interiorBoundaryTypes(element, scvf);
98 if (bcTypes.hasOnlyNeumann())
99 {
100 // compute tpfa flux from integration point to facet centerline
101 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
102
103 // interpolate temperature to scvf integration point
104 Scalar T = 0.0;
105 for (const auto& scv : scvs(fvGeometry))
106 T += elemVolVars[scv].temperature()*shapeValues[scv.indexInElement()][0];
107
108 using std::sqrt;
109 // If this is a surface grid, use the square root of the facet extrusion factor
110 // as an approximate average distance from scvf ip to facet center
111 using std::sqrt;
112 const auto a = facetVolVars.extrusionFactor();
113 auto gradT = scvf.unitOuterNormal();
114 gradT *= dim == dimWorld ? 0.5*a : 0.5*sqrt(a);
115 gradT /= gradT.two_norm2();
116 gradT *= (facetVolVars.temperature() - T);
117
118 return -1.0*scvf.area()
119 *insideVolVars.extrusionFactor()
120 *vtmv(scvf.unitOuterNormal(),
121 facetVolVars.effectiveThermalConductivity(),
122 gradT);
123 }
124
125 // on interior Dirichlet boundaries use the facet temperature and evaluate flux
126 else if (bcTypes.hasOnlyDirichlet())
127 {
128 // create vector with nodal temperatures
129 std::vector<Scalar> temperatures(element.subEntities(dim));
130 for (const auto& scv : scvs(fvGeometry))
131 temperatures[scv.localDofIndex()] = elemVolVars[scv].temperature();
132
133 // substitute with facet temperatures for those scvs touching this facet
134 for (const auto& scvfJ : scvfs(fvGeometry))
135 if (scvfJ.interiorBoundary() && scvfJ.facetIndexInElement() == scvf.facetIndexInElement())
136 temperatures[ fvGeometry.scv(scvfJ.insideScvIdx()).localDofIndex() ]
137 = problem.couplingManager().getLowDimVolVars(element, scvfJ).temperature();
138
139 // evaluate gradT at integration point
140 Dune::FieldVector<Scalar, dimWorld> gradT(0.0);
141 for (const auto& scv : scvs(fvGeometry))
142 gradT.axpy(temperatures[scv.localDofIndex()], fluxVarCache.gradN(scv.indexInElement()));
143
144 // apply matrix thermal conductivity and return the flux
145 return -1.0*scvf.area()
146 *insideVolVars.extrusionFactor()
147 *vtmv(scvf.unitOuterNormal(),
148 insideVolVars.effectiveThermalConductivity(),
149 gradT);
150 }
151
152 // mixed boundary types are not supported
153 else
154 DUNE_THROW(Dune::NotImplemented, "Mixed boundary types are not supported");
155 }
156
157 // compute transmissibilities ti for analytical jacobians
158 template<class Problem, class ElementVolumeVariables, class FluxVarCache>
159 static std::vector<Scalar> calculateTransmissibilities(const Problem& problem,
160 const Element& element,
161 const FVElementGeometry& fvGeometry,
162 const ElementVolumeVariables& elemVolVars,
163 const SubControlVolumeFace& scvf,
164 const FluxVarCache& fluxVarCache,
165 unsigned int phaseIdx)
166 {
167 DUNE_THROW(Dune::NotImplemented, "transmissibilty computation for BoxFacetCouplingFouriersLaw");
168 }
169};
170
171} // end namespace Dumux
172
173#endif
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
forward declaration of the method-specific implementation
Definition: flux/fourierslaw.hh:37
Specialization of Fourier's Law for the box method.
Definition: flux/box/fourierslaw.hh:45
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache)
Definition: flux/box/fourierslaw.hh:56
Fourier's law for the box scheme in the context of coupled models where coupling occurs across the fa...
Definition: multidomain/facet/box/fourierslaw.hh:53
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVarsCache &elemFluxVarCache)
Compute the conductive heat flux across a sub-control volume face.
Definition: multidomain/facet/box/fourierslaw.hh:75
static std::vector< Scalar > calculateTransmissibilities(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const FluxVarCache &fluxVarCache, unsigned int phaseIdx)
Definition: multidomain/facet/box/fourierslaw.hh:159
Declares all properties used in Dumux.
This file contains the data which is required to calculate energy fluxes due to molecular diffusion w...