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