version 3.8
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// 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_FACETCOUPLING_MANAGER_HH
13#define DUMUX_FACETCOUPLING_MANAGER_HH
14
17
21
22namespace Dumux {
23namespace FacetCoupling{
24
39template<class VolumeVariables, class Problem, class SolutionVector, class FVGeometry>
40void makeInterpolatedVolVars(VolumeVariables& volVars,
41 const Problem& problem,
42 const SolutionVector& sol,
43 const FVGeometry& fvGeometry,
44 const typename FVGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
45 const typename FVGeometry::GridGeometry::GridView::template Codim<0>::Entity::Geometry& elemGeom,
46 const typename FVGeometry::GridGeometry::GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate& pos)
47{
48 // interpolate solution and set it for each entry in element solution
49 auto elemSol = elementSolution(element, sol, fvGeometry.gridGeometry());
50 const auto centerSol = evalSolution(element, elemGeom, fvGeometry.gridGeometry(), elemSol, pos);
51 for (unsigned int i = 0; i < fvGeometry.numScv(); ++i)
52 elemSol[i] = centerSol;
53
54 // Update volume variables with the interpolated solution. Note that this standard
55 // implementation only works for element-wise constant parameters as we simply use
56 // the first element scv for the vol var update. For heterogeneities within the element
57 // or more complex models (e.g. 2p with interface solver) a corresponding overload
58 // of this function has to be provided!
59 volVars.update(elemSol, problem, element, *scvs(fvGeometry).begin());
60}
61} // end namespace FacetCoupling
62
78template< class MDTraits,
79 class CouplingMapper,
80 std::size_t bulkDomainId = 0,
81 std::size_t lowDimDomainId = 1,
82 class DiscretizationMethod = typename GetPropType<typename MDTraits::template SubDomain<bulkDomainId>::TypeTag, Properties::GridGeometry>::DiscretizationMethod >
84
97template< class MDTraits,
98 class CouplingMapper,
99 std::size_t bulkDomainId = 0,
100 std::size_t facetDomainId = 1,
101 std::size_t edgeDomainId = 2 >
103: public FacetCouplingManager<MDTraits, CouplingMapper, bulkDomainId, facetDomainId>
104, public FacetCouplingManager<MDTraits, CouplingMapper, facetDomainId, edgeDomainId>
105{
108
109 // convenience aliases and instances of the domain ids
110 using BulkIdType = typename MDTraits::template SubDomain<bulkDomainId>::Index;
111 using FacetIdType = typename MDTraits::template SubDomain<facetDomainId>::Index;
112 using EdgeIdType = typename MDTraits::template SubDomain<edgeDomainId>::Index;
113 static constexpr auto bulkId = BulkIdType();
114 static constexpr auto facetId = FacetIdType();
115 static constexpr auto edgeId = EdgeIdType();
116
117 // the sub-domain type tags
118 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
119
120 // further types specific to the sub-problems
121 template<std::size_t id> using PrimaryVariables = GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>;
122 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
123
124 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
125 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
126 template<std::size_t id> using SubControlVolumeFace = typename GridGeometry<id>::SubControlVolumeFace;
127 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
128 template<std::size_t id> using GridIndexType = typename IndexTraits<GridView<id>>::GridIndex;
129 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
130
131 template<std::size_t id> using GridVariables = GetPropType<SubDomainTypeTag<id>, Properties::GridVariables>;
132 template<std::size_t id> using ElementVolumeVariables = typename GridVariables<id>::GridVolumeVariables::LocalView;
133 template<std::size_t id> using ElementFluxVariablesCache = typename GridVariables<id>::GridFluxVariablesCache::LocalView;
134
135 // helper function to check if a domain uses mpfa
136 template<std::size_t id>
137 static constexpr bool usesMpfa(Dune::index_constant<id> domainId)
138 { return GridGeometry<domainId>::discMethod == DiscretizationMethods::ccmpfa; }
139
140public:
142 template<std::size_t i, std::size_t j>
143 using CouplingStencilType = typename std::conditional< (j == edgeDomainId),
144 typename FacetEdgeManager::template CouplingStencilType<i, j>,
145 typename BulkFacetManager::template CouplingStencilType<i, j> >::type;
146
148 using SolutionVector = typename MDTraits::SolutionVector;
149
159 void init(std::shared_ptr< Problem<bulkId> > bulkProblem,
160 std::shared_ptr< Problem<facetId> > facetProblem,
161 std::shared_ptr< Problem<edgeId> > edgeProblem,
162 std::shared_ptr< CouplingMapper > couplingMapper,
163 const SolutionVector& curSol)
164 {
165 BulkFacetManager::init(bulkProblem, facetProblem, couplingMapper, curSol);
166 FacetEdgeManager::init(facetProblem, edgeProblem, couplingMapper, curSol);
167 }
168
170 using BulkFacetManager::couplingStencil;
171 using FacetEdgeManager::couplingStencil;
172
173 using BulkFacetManager::isCoupled;
174 using FacetEdgeManager::isCoupled;
175
176 using BulkFacetManager::isOnInteriorBoundary;
177 using FacetEdgeManager::isOnInteriorBoundary;
178
179 using BulkFacetManager::getLowDimVolVars;
180 using FacetEdgeManager::getLowDimVolVars;
181
182 using BulkFacetManager::getLowDimElement;
183 using FacetEdgeManager::getLowDimElement;
184
185 using BulkFacetManager::getLowDimElementIndex;
186 using FacetEdgeManager::getLowDimElementIndex;
187
188 using BulkFacetManager::evalSourcesFromBulk;
189 using FacetEdgeManager::evalSourcesFromBulk;
190
191 using BulkFacetManager::evalCouplingResidual;
192 using FacetEdgeManager::evalCouplingResidual;
193
194 using BulkFacetManager::bindCouplingContext;
195 using FacetEdgeManager::bindCouplingContext;
196
197 using BulkFacetManager::updateCouplingContext;
198 using FacetEdgeManager::updateCouplingContext;
199
200 using BulkFacetManager::updateCoupledVariables;
201 using FacetEdgeManager::updateCoupledVariables;
202
203 using BulkFacetManager::extendJacobianPattern;
204 using FacetEdgeManager::extendJacobianPattern;
205
206 using BulkFacetManager::evalAdditionalDomainDerivatives;
207 using FacetEdgeManager::evalAdditionalDomainDerivatives;
208
209 // extension of the jacobian pattern for the facet domain only occurs
210 // within the bulk-facet coupling & for mpfa being used in the bulk domain.
211 template<class JacobianPattern>
212 void extendJacobianPattern(FacetIdType, JacobianPattern& pattern) const
213 { BulkFacetManager::extendJacobianPattern(facetId, pattern); }
214
215 template<class FacetLocalAssembler, class JacobianMatrixDiagBlock, class GridVariables>
217 const FacetLocalAssembler& facetLocalAssembler,
218 const typename FacetLocalAssembler::LocalResidual::ElementResidualVector& origResiduals,
219 JacobianMatrixDiagBlock& A,
220 GridVariables& gridVariables)
221 { BulkFacetManager::evalAdditionalDomainDerivatives(facetId, facetLocalAssembler, origResiduals, A, gridVariables); }
222
227 const Element<bulkId>& element,
228 EdgeIdType domainJ) const
229 { return FacetEdgeManager::getEmptyStencil(edgeId); }
230
235 const Element<edgeId>& element,
236 BulkIdType domainJ) const
237 { return BulkFacetManager::getEmptyStencil(bulkId); }
238
244 {
245 BulkFacetManager::updateSolution(sol);
246 FacetEdgeManager::updateSolution(sol);
247 }
248
255 template<std::size_t i,
256 std::size_t j,
257 class LocalAssembler,
258 std::enable_if_t<((i==bulkId && j==edgeId) || ((i==edgeId && j==bulkId))), int> = 0>
259 typename LocalAssembler::LocalResidual::ElementResidualVector
260 evalCouplingResidual(Dune::index_constant<i> domainI,
261 const LocalAssembler& localAssembler,
262 Dune::index_constant<j> domainJ,
263 GridIndexType<j> dofIdxGlobalJ)
264 {
265 typename LocalAssembler::LocalResidual::ElementResidualVector res(1);
266 res = 0.0;
267 return res;
268 }
269
274 template< class Assembler >
275 void bindCouplingContext(FacetIdType, const Element<facetId>& element, const Assembler& assembler)
276 {
277 BulkFacetManager::bindCouplingContext(facetId, element, assembler);
278 FacetEdgeManager::bindCouplingContext(facetId, element, assembler);
279 }
280
285 template< class FacetLocalAssembler >
286 void updateCouplingContext(FacetIdType domainI,
287 const FacetLocalAssembler& facetLocalAssembler,
288 FacetIdType domainJ,
289 GridIndexType<facetId> dofIdxGlobalJ,
290 const PrimaryVariables<facetId>& priVarsJ,
291 unsigned int pvIdxJ)
292 {
293 BulkFacetManager::updateCouplingContext(domainI, facetLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
294 FacetEdgeManager::updateCouplingContext(domainI, facetLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
295 }
296
302 template<std::size_t i,
303 std::size_t j,
304 class LocalAssembler,
305 std::enable_if_t<((i==bulkId && j==edgeId) || (i==edgeId && j==bulkId)), int> = 0>
306 void updateCouplingContext(Dune::index_constant<i> domainI,
307 const LocalAssembler& localAssembler,
308 Dune::index_constant<j> domainJ,
309 GridIndexType<j> dofIdxGlobalJ,
310 const PrimaryVariables<j>& priVarsJ,
311 unsigned int pvIdxJ)
312 { /*do nothing here*/ }
313
319 template< class FacetLocalAssembler, class UpdatableElementVolVars, class UpdatableFluxVarCache>
320 void updateCoupledVariables(FacetIdType domainI,
321 const FacetLocalAssembler& facetLocalAssembler,
322 UpdatableElementVolVars& elemVolVars,
323 UpdatableFluxVarCache& elemFluxVarsCache)
324 {
325 BulkFacetManager::updateCoupledVariables(domainI, facetLocalAssembler, elemVolVars, elemFluxVarsCache);
326 FacetEdgeManager::updateCoupledVariables(domainI, facetLocalAssembler, elemVolVars, elemFluxVarsCache);
327 }
328
330 template<std::size_t id, std::enable_if_t<(id == bulkId || id == facetId), int> = 0>
331 const Problem<id>& problem() const { return BulkFacetManager::template problem<id>(); }
332
334 template<std::size_t id, std::enable_if_t<(id == bulkId || id == facetId), int> = 0>
335 Problem<id>& problem() { return BulkFacetManager::template problem<id>(); }
336
338 template<std::size_t id, std::enable_if_t<(id == edgeId), int> = 0>
339 const Problem<id>& problem() const { return FacetEdgeManager::template problem<id>(); }
340
342 template<std::size_t id, std::enable_if_t<(id == edgeId), int> = 0>
343 Problem<id>& problem() { return FacetEdgeManager::template problem<id>(); }
344};
345
346} // end namespace Dumux
347
348// Here, we have to include all available implementations
352
353#endif
Implementation for the coupling manager between two domains of dimension d and (d-1) for models consi...
Definition: multidomain/facet/couplingmanager.hh:83
Class that handles the coupling between three sub-domains in models where the coupling between the tw...
Definition: multidomain/facet/couplingmanager.hh:105
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:306
const Problem< id > & problem() const
Return a const reference to bulk or facet problem.
Definition: multidomain/facet/couplingmanager.hh:331
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:320
typename MDTraits::SolutionVector SolutionVector
the type of the solution vector
Definition: multidomain/facet/couplingmanager.hh:148
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:243
void evalAdditionalDomainDerivatives(FacetIdType, const FacetLocalAssembler &facetLocalAssembler, const typename FacetLocalAssembler::LocalResidual::ElementResidualVector &origResiduals, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Definition: multidomain/facet/couplingmanager.hh:216
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:286
Problem< id > & problem()
Return a reference to bulk or facet problem.
Definition: multidomain/facet/couplingmanager.hh:335
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:145
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:226
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:234
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:275
void extendJacobianPattern(FacetIdType, JacobianPattern &pattern) const
Definition: multidomain/facet/couplingmanager.hh:212
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:159
LocalAssembler::LocalResidual::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:260
Defines all properties used in Dumux.
Element solution classes and factory functions.
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 CVFEElementSolution< 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:152
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition: cellcentered/elementsolution.hh:101
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:40
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Defines the index types used for grid and local indices.
The available discretization methods in Dumux.
constexpr CCMpfa ccmpfa
Definition: method.hh:146
Definition: adapt.hh:17
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:26