Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup FacetCoupling | ||
10 | * \copydoc Dumux::BoxFacetCouplingFouriersLaw | ||
11 | */ | ||
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> | ||
23 | #include <dumux/common/parameters.hh> | ||
24 | #include <dumux/common/properties.hh> | ||
25 | |||
26 | #include <dumux/flux/referencesystemformulation.hh> | ||
27 | #include <dumux/flux/box/fourierslaw.hh> | ||
28 | #include <dumux/discretization/method.hh> | ||
29 | #include <dumux/discretization/extrusion.hh> | ||
30 | |||
31 | namespace Dumux { | ||
32 | |||
33 | /*! | ||
34 | * \ingroup FacetCoupling | ||
35 | * \brief Fourier's law for the box scheme in the context of coupled models | ||
36 | * where coupling occurs across the facets of the bulk domain elements | ||
37 | * with a lower-dimensional domain living on these facets. | ||
38 | */ | ||
39 | template<class TypeTag> | ||
40 | class BoxFacetCouplingFouriersLaw | ||
41 | : public FouriersLawImplementation<TypeTag, DiscretizationMethods::Box> | ||
42 | { | ||
43 | using ParentType = FouriersLawImplementation<TypeTag, DiscretizationMethods::Box>; | ||
44 | |||
45 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
46 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
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 | |||
59 | public: | ||
60 | |||
61 | /*! | ||
62 | * \brief Compute the conductive heat flux across a sub-control volume face. | ||
63 | */ | ||
64 | template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache> | ||
65 | 4612608 | 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 |
2/2✓ Branch 0 taken 4239872 times.
✓ Branch 1 taken 372736 times.
|
4612608 | if (!scvf.interiorBoundary()) |
74 | 4239872 | return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarCache); | |
75 | |||
76 |
5/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 372734 times.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
|
372736 | static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0); |
77 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 372736 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 372736 times.
|
372736 | 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 | 372736 | const auto& fluxVarCache = elemFluxVarCache[scvf]; | |
82 | 372736 | const auto& shapeValues = fluxVarCache.shapeValues(); | |
83 | 745472 | const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx()); | |
84 | 372736 | const auto& insideVolVars = elemVolVars[insideScv]; | |
85 | |||
86 | // on interior Neumann boundaries, evaluate the flux using the facet thermal conductivity | ||
87 | 372736 | const auto bcTypes = problem.interiorBoundaryTypes(element, scvf); | |
88 |
2/4✓ Branch 0 taken 372736 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 372736 times.
✗ Branch 3 not taken.
|
745472 | if (bcTypes.hasOnlyNeumann()) |
89 | { | ||
90 | // compute tpfa flux from integration point to facet centerline | ||
91 | 745472 | const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf); | |
92 | |||
93 | // interpolate temperature to scvf integration point | ||
94 | 372736 | Scalar T = 0.0; | |
95 |
4/4✓ Branch 0 taken 1490944 times.
✓ Branch 1 taken 372736 times.
✓ Branch 2 taken 1490944 times.
✓ Branch 3 taken 372736 times.
|
3727360 | for (const auto& scv : scvs(fvGeometry)) |
96 | 5963776 | 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 | 372736 | const auto a = facetVolVars.extrusionFactor(); | |
103 | 372736 | auto gradT = scvf.unitOuterNormal(); | |
104 | 372736 | gradT *= dim == dimWorld ? 0.5*a : 0.5*sqrt(a); | |
105 | 372736 | gradT /= gradT.two_norm2(); | |
106 | 745472 | gradT *= (facetVolVars.temperature() - T); | |
107 | |||
108 | 372736 | return -1.0*Extrusion::area(fvGeometry, scvf) | |
109 | 372736 | *insideVolVars.extrusionFactor() | |
110 | 1118208 | *vtmv(scvf.unitOuterNormal(), | |
111 | facetVolVars.effectiveThermalConductivity(), | ||
112 | 372736 | 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 | ||
164 |