3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
multidomain/facet/couplingmanager.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_FACETCOUPLING_MANAGER_HH
25#define DUMUX_FACETCOUPLING_MANAGER_HH
26
29
33
34namespace Dumux {
35namespace FacetCoupling{
36
51template<class VolumeVariables, class Problem, class SolutionVector, class FVGeometry>
52void makeInterpolatedVolVars(VolumeVariables& volVars,
53 const Problem& problem,
54 const SolutionVector& sol,
55 const FVGeometry& fvGeometry,
56 const typename FVGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
57 const typename FVGeometry::GridGeometry::GridView::template Codim<0>::Entity::Geometry& elemGeom,
58 const typename FVGeometry::GridGeometry::GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate& pos)
59{
60 // interpolate solution and set it for each entry in element solution
61 auto elemSol = elementSolution(element, sol, fvGeometry.gridGeometry());
62 const auto centerSol = evalSolution(element, elemGeom, fvGeometry.gridGeometry(), elemSol, pos);
63 for (unsigned int i = 0; i < fvGeometry.numScv(); ++i)
64 elemSol[i] = centerSol;
65
66 // Update volume variables with the interpolated solution. Note that this standard
67 // implementation only works for element-wise constant parameters as we simply use
68 // the first element scv for the vol var update. For heterogeneities within the element
69 // or more complex models (e.g. 2p with interface solver) a corresponding overload
70 // of this function has to be provided!
71 volVars.update(elemSol, problem, element, *scvs(fvGeometry).begin());
72}
73} // end namespace FacetCoupling
74
90template< class MDTraits,
91 class CouplingMapper,
92 std::size_t bulkDomainId = 0,
93 std::size_t lowDimDomainId = 1,
94 DiscretizationMethod bulkDM = GetPropType<typename MDTraits::template SubDomain<bulkDomainId>::TypeTag, Properties::GridGeometry>::discMethod >
96
109template< class MDTraits,
110 class CouplingMapper,
111 std::size_t bulkDomainId = 0,
112 std::size_t facetDomainId = 1,
113 std::size_t edgeDomainId = 2 >
115: public FacetCouplingManager<MDTraits, CouplingMapper, bulkDomainId, facetDomainId>
116, public FacetCouplingManager<MDTraits, CouplingMapper, facetDomainId, edgeDomainId>
117{
120
121 // convenience aliases and instances of the domain ids
122 using BulkIdType = typename MDTraits::template SubDomain<bulkDomainId>::Index;
123 using FacetIdType = typename MDTraits::template SubDomain<facetDomainId>::Index;
124 using EdgeIdType = typename MDTraits::template SubDomain<edgeDomainId>::Index;
125 static constexpr auto bulkId = BulkIdType();
126 static constexpr auto facetId = FacetIdType();
127 static constexpr auto edgeId = EdgeIdType();
128
129 // the sub-domain type tags
130 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
131
132 // further types specific to the sub-problems
133 template<std::size_t id> using LocalResidual = GetPropType<SubDomainTypeTag<id>, Properties::LocalResidual>;
134 template<std::size_t id> using PrimaryVariables = GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>;
135 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
136
137 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
138 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
139 template<std::size_t id> using SubControlVolumeFace = typename GridGeometry<id>::SubControlVolumeFace;
140 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
141 template<std::size_t id> using GridIndexType = typename IndexTraits<GridView<id>>::GridIndex;
142 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
143
144 template<std::size_t id> using GridVariables = GetPropType<SubDomainTypeTag<id>, Properties::GridVariables>;
145 template<std::size_t id> using ElementVolumeVariables = typename GridVariables<id>::GridVolumeVariables::LocalView;
146 template<std::size_t id> using ElementFluxVariablesCache = typename GridVariables<id>::GridFluxVariablesCache::LocalView;
147
148 // helper function to check if a domain uses mpfa
149 template<std::size_t id>
150 static constexpr bool usesMpfa(Dune::index_constant<id> domainId)
151 { return GridGeometry<domainId>::discMethod == DiscretizationMethod::ccmpfa; }
152
153public:
155 template<std::size_t i, std::size_t j>
156 using CouplingStencilType = typename std::conditional< (j == edgeDomainId),
157 typename FacetEdgeManager::template CouplingStencilType<i, j>,
158 typename BulkFacetManager::template CouplingStencilType<i, j> >::type;
159
161 using SolutionVector = typename MDTraits::SolutionVector;
162
172 void init(std::shared_ptr< Problem<bulkId> > bulkProblem,
173 std::shared_ptr< Problem<facetId> > facetProblem,
174 std::shared_ptr< Problem<edgeId> > edgeProblem,
175 std::shared_ptr< CouplingMapper > couplingMapper,
176 const SolutionVector& curSol)
177 {
178 BulkFacetManager::init(bulkProblem, facetProblem, couplingMapper, curSol);
179 FacetEdgeManager::init(facetProblem, edgeProblem, couplingMapper, curSol);
180 }
181
183 using BulkFacetManager::couplingStencil;
184 using FacetEdgeManager::couplingStencil;
185
186 using BulkFacetManager::isCoupled;
187 using FacetEdgeManager::isCoupled;
188
189 using BulkFacetManager::isOnInteriorBoundary;
190 using FacetEdgeManager::isOnInteriorBoundary;
191
192 using BulkFacetManager::getLowDimVolVars;
193 using FacetEdgeManager::getLowDimVolVars;
194
195 using BulkFacetManager::getLowDimElement;
196 using FacetEdgeManager::getLowDimElement;
197
198 using BulkFacetManager::getLowDimElementIndex;
199 using FacetEdgeManager::getLowDimElementIndex;
200
201 using BulkFacetManager::evalSourcesFromBulk;
202 using FacetEdgeManager::evalSourcesFromBulk;
203
204 using BulkFacetManager::evalCouplingResidual;
205 using FacetEdgeManager::evalCouplingResidual;
206
207 using BulkFacetManager::bindCouplingContext;
208 using FacetEdgeManager::bindCouplingContext;
209
210 using BulkFacetManager::updateCouplingContext;
211 using FacetEdgeManager::updateCouplingContext;
212
213 using BulkFacetManager::updateCoupledVariables;
214 using FacetEdgeManager::updateCoupledVariables;
215
216 using BulkFacetManager::extendJacobianPattern;
217 using FacetEdgeManager::extendJacobianPattern;
218
219 using BulkFacetManager::evalAdditionalDomainDerivatives;
220 using FacetEdgeManager::evalAdditionalDomainDerivatives;
221
222 // extension of the jacobian pattern for the facet domain only occurs
223 // within the bulk-facet coupling & for mpfa being used in the bulk domain.
224 template<class JacobianPattern>
225 void extendJacobianPattern(FacetIdType, JacobianPattern& pattern) const
226 { BulkFacetManager::extendJacobianPattern(facetId, pattern); }
227
228 template<class FacetLocalAssembler, class JacobianMatrixDiagBlock, class GridVariables>
230 const FacetLocalAssembler& facetLocalAssembler,
231 const typename FacetLocalAssembler::LocalResidual::ElementResidualVector& origResiduals,
232 JacobianMatrixDiagBlock& A,
233 GridVariables& gridVariables)
234 { BulkFacetManager::evalAdditionalDomainDerivatives(facetId, facetLocalAssembler, origResiduals, A, gridVariables); }
235
240 const Element<bulkId>& element,
241 EdgeIdType domainJ) const
242 { return FacetEdgeManager::getEmptyStencil(edgeId); }
243
248 const Element<edgeId>& element,
249 BulkIdType domainJ) const
250 { return BulkFacetManager::getEmptyStencil(bulkId); }
251
257 {
258 BulkFacetManager::updateSolution(sol);
259 FacetEdgeManager::updateSolution(sol);
260 }
261
268 template<std::size_t i,
269 std::size_t j,
270 class LocalAssembler,
271 std::enable_if_t<((i==bulkId && j==edgeId) || ((i==edgeId && j==bulkId))), int> = 0>
272 typename LocalResidual<i>::ElementResidualVector
273 evalCouplingResidual(Dune::index_constant<i> domainI,
274 const LocalAssembler& localAssembler,
275 Dune::index_constant<j> domainJ,
276 GridIndexType<j> dofIdxGlobalJ)
277 {
278 typename LocalResidual<i>::ElementResidualVector res(1);
279 res = 0.0;
280 return res;
281 }
282
287 template< class Assembler >
288 void bindCouplingContext(FacetIdType, const Element<facetId>& element, const Assembler& assembler)
289 {
290 BulkFacetManager::bindCouplingContext(facetId, element, assembler);
291 FacetEdgeManager::bindCouplingContext(facetId, element, assembler);
292 }
293
298 template< class FacetLocalAssembler >
299 void updateCouplingContext(FacetIdType domainI,
300 const FacetLocalAssembler& facetLocalAssembler,
301 FacetIdType domainJ,
302 GridIndexType<facetId> dofIdxGlobalJ,
303 const PrimaryVariables<facetId>& priVarsJ,
304 unsigned int pvIdxJ)
305 {
306 BulkFacetManager::updateCouplingContext(domainI, facetLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
307 FacetEdgeManager::updateCouplingContext(domainI, facetLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
308 }
309
315 template<std::size_t i,
316 std::size_t j,
317 class LocalAssembler,
318 std::enable_if_t<((i==bulkId && j==edgeId) || (i==edgeId && j==bulkId)), int> = 0>
319 void updateCouplingContext(Dune::index_constant<i> domainI,
320 const LocalAssembler& localAssembler,
321 Dune::index_constant<j> domainJ,
322 GridIndexType<j> dofIdxGlobalJ,
323 const PrimaryVariables<j>& priVarsJ,
324 unsigned int pvIdxJ)
325 { /*do nothing here*/ }
326
332 template< class FacetLocalAssembler, class UpdatableElementVolVars, class UpdatableFluxVarCache>
333 void updateCoupledVariables(FacetIdType domainI,
334 const FacetLocalAssembler& facetLocalAssembler,
335 UpdatableElementVolVars& elemVolVars,
336 UpdatableFluxVarCache& elemFluxVarsCache)
337 {
338 BulkFacetManager::updateCoupledVariables(domainI, facetLocalAssembler, elemVolVars, elemFluxVarsCache);
339 FacetEdgeManager::updateCoupledVariables(domainI, facetLocalAssembler, elemVolVars, elemFluxVarsCache);
340 }
341
343 template<std::size_t id, std::enable_if_t<(id == bulkId || id == facetId), int> = 0>
344 const Problem<id>& problem() const { return BulkFacetManager::template problem<id>(); }
345
347 template<std::size_t id, std::enable_if_t<(id == bulkId || id == facetId), int> = 0>
348 Problem<id>& problem() { return BulkFacetManager::template problem<id>(); }
349
351 template<std::size_t id, std::enable_if_t<(id == edgeId), int> = 0>
352 const Problem<id>& problem() const { return FacetEdgeManager::template problem<id>(); }
353
355 template<std::size_t id, std::enable_if_t<(id == edgeId), int> = 0>
356 Problem<id>& problem() { return FacetEdgeManager::template problem<id>(); }
357};
358
359} // end namespace Dumux
360
361// Here, we have to include all available implementations
365
366#endif
Defines the index types used for grid and local indices.
Element solution classes and factory functions.
The available discretization methods in Dumux.
free functions for the evaluation of primary variables inside elements.
PrimaryVariables evalSolution(const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const BoxElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
Interpolates a given box element solution at a given global position. Uses the finite element cache o...
Definition: evalsolution.hh:94
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethod::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:115
void makeInterpolatedVolVars(VolumeVariables &volVars, const Problem &problem, const SolutionVector &sol, const FVGeometry &fvGeometry, const typename FVGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const typename FVGeometry::GridGeometry::GridView::template Codim< 0 >::Entity::Geometry &elemGeom, const typename FVGeometry::GridGeometry::GridView::template Codim< 0 >::Entity::Geometry::GlobalCoordinate &pos)
Free function that allows the creation of a volume variables object interpolated to a given position ...
Definition: multidomain/facet/couplingmanager.hh:52
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
Struture to define the index types used for grid and local indices.
Definition: indextraits.hh:38
A vector of primary variables.
Definition: common/properties.hh:60
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:70
Definition: common/properties.hh:89
Definition: common/properties.hh:113
The grid variables object managing variable data on the grid (volvars/fluxvars cache)
Definition: common/properties.hh:134
Implementation for the coupling manager between two domains of dimension d and (d-1) for models consi...
Definition: multidomain/facet/couplingmanager.hh:95
Class that handles the coupling between three sub-domains in models where the coupling between the tw...
Definition: multidomain/facet/couplingmanager.hh:117
void updateCouplingContext(Dune::index_constant< i > domainI, const LocalAssembler &localAssembler, Dune::index_constant< j > domainJ, GridIndexType< j > dofIdxGlobalJ, const PrimaryVariables< j > &priVarsJ, unsigned int pvIdxJ)
Interface for updating the coupling context between the bulk and the edge domain. We do nothing here ...
Definition: multidomain/facet/couplingmanager.hh:319
const Problem< id > & problem() const
Return a const reference to bulk or facet problem.
Definition: multidomain/facet/couplingmanager.hh:344
void updateCoupledVariables(FacetIdType domainI, const FacetLocalAssembler &facetLocalAssembler, UpdatableElementVolVars &elemVolVars, UpdatableFluxVarCache &elemFluxVarsCache)
Interface for updating the local views of the facet domain after updateCouplingContext the coupling c...
Definition: multidomain/facet/couplingmanager.hh:333
typename MDTraits::SolutionVector SolutionVector
the type of the solution vector
Definition: multidomain/facet/couplingmanager.hh:161
void updateSolution(const SolutionVector &sol)
updates the current solution. We have to overload this here to avoid ambiguity and update the solutio...
Definition: multidomain/facet/couplingmanager.hh:256
void evalAdditionalDomainDerivatives(FacetIdType, const FacetLocalAssembler &facetLocalAssembler, const typename FacetLocalAssembler::LocalResidual::ElementResidualVector &origResiduals, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Definition: multidomain/facet/couplingmanager.hh:229
void updateCouplingContext(FacetIdType domainI, const FacetLocalAssembler &facetLocalAssembler, FacetIdType domainJ, GridIndexType< facetId > dofIdxGlobalJ, const PrimaryVariables< facetId > &priVarsJ, unsigned int pvIdxJ)
Interface for updating the coupling context of the facet domain. In this case we have to update both ...
Definition: multidomain/facet/couplingmanager.hh:299
Problem< id > & problem()
Return a reference to bulk or facet problem.
Definition: multidomain/facet/couplingmanager.hh:348
typename std::conditional<(j==edgeDomainId), typename FacetEdgeManager::template CouplingStencilType< i, j >, typename BulkFacetManager::template CouplingStencilType< i, j > >::type CouplingStencilType
types used for coupling stencils
Definition: multidomain/facet/couplingmanager.hh:158
const CouplingStencilType< bulkId, edgeId > & couplingStencil(BulkIdType domainI, const Element< bulkId > &element, EdgeIdType domainJ) const
The coupling stencil of the bulk with the edge domain (empty stencil).
Definition: multidomain/facet/couplingmanager.hh:239
LocalResidual< i >::ElementResidualVector evalCouplingResidual(Dune::index_constant< i > domainI, const LocalAssembler &localAssembler, Dune::index_constant< j > domainJ, GridIndexType< j > dofIdxGlobalJ)
Interface for evaluating the coupling residual between the bulk and the edge domain....
Definition: multidomain/facet/couplingmanager.hh:273
const CouplingStencilType< edgeId, bulkId > & couplingStencil(EdgeIdType domainI, const Element< edgeId > &element, BulkIdType domainJ) const
The coupling stencil of the edge with the bulk domain (empty stencil).
Definition: multidomain/facet/couplingmanager.hh:247
void bindCouplingContext(FacetIdType, const Element< facetId > &element, const Assembler &assembler)
Interface for binding the coupling context for the facet domain. In this case we have to bind both th...
Definition: multidomain/facet/couplingmanager.hh:288
void extendJacobianPattern(FacetIdType, JacobianPattern &pattern) const
Definition: multidomain/facet/couplingmanager.hh:225
void init(std::shared_ptr< Problem< bulkId > > bulkProblem, std::shared_ptr< Problem< facetId > > facetProblem, std::shared_ptr< Problem< edgeId > > edgeProblem, std::shared_ptr< CouplingMapper > couplingMapper, const SolutionVector &curSol)
Initialize the coupling manager.
Definition: multidomain/facet/couplingmanager.hh:172
Declares all properties used in Dumux.