version 3.10-dev
multidomain/facet/cellcentered/mpfa/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_CCMPFA_FACETCOUPLING_MANAGER_HH
13#define DUMUX_CCMPFA_FACETCOUPLING_MANAGER_HH
14
15#include <algorithm>
16#include <cassert>
17
23
28
29namespace Dumux {
30
42template<class MDTraits, class CouplingMapper, std::size_t bulkDomainId, std::size_t lowDimDomainId>
43class FacetCouplingManager<MDTraits, CouplingMapper, bulkDomainId, lowDimDomainId, DiscretizationMethods::CCMpfa>
44: public FacetCouplingManager<MDTraits, CouplingMapper, bulkDomainId, lowDimDomainId, DiscretizationMethods::CCTpfa>
45{
47
48 // domain id instances
49 using BulkIdType = typename MDTraits::template SubDomain<bulkDomainId>::Index;
50 using LowDimIdType = typename MDTraits::template SubDomain<lowDimDomainId>::Index;
51 static constexpr auto bulkId = BulkIdType();
52 static constexpr auto lowDimId = LowDimIdType();
53
54 // the sub-domain type tags
55 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
56
57 // further types specific to the sub-problems
58 template<std::size_t id> using Scalar = GetPropType<SubDomainTypeTag<id>, Properties::Scalar>;
59 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
60 template<std::size_t id> using NumEqVector = Dumux::NumEqVector<GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>>;
61 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
62 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
63 template<std::size_t id> using SubControlVolume = typename GridGeometry<id>::SubControlVolume;
64 template<std::size_t id> using SubControlVolumeFace = typename GridGeometry<id>::SubControlVolumeFace;
65 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
66 template<std::size_t id> using GridIndexType = typename IndexTraits< GridView<id> >::GridIndex;
67 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
68 template<std::size_t id> using LocalResidual = typename MDTraits::template SubDomain<id>::LocalResidual;
69
70 template<std::size_t id> using GridVariables = GetPropType<SubDomainTypeTag<id>, Properties::GridVariables>;
71 template<std::size_t id> using GridVolumeVariables = typename GridVariables<id>::GridVolumeVariables;
72 template<std::size_t id> using ElementVolumeVariables = typename GridVolumeVariables<id>::LocalView;
73
74 // grid ids
75 static constexpr int bulkDim = GridView<bulkDomainId>::dimension;
76 static constexpr int lowDimDim = GridView<lowDimDomainId>::dimension;
77 static constexpr auto bulkGridId = CouplingMapper::template gridId<bulkDim>();
78 static constexpr auto lowDimGridId = CouplingMapper::template gridId<lowDimDim>();
79
80 static constexpr bool lowDimUsesBox = GridGeometry<lowDimId>::discMethod == DiscretizationMethods::box;
81
82public:
83
85 using SolutionVector = typename MDTraits::SolutionVector;
86
95 void init(std::shared_ptr< Problem<bulkId> > bulkProblem,
96 std::shared_ptr< Problem<lowDimId> > lowDimProblem,
97 std::shared_ptr< CouplingMapper > couplingMapper,
98 const SolutionVector& curSol)
99 {
100 // Initialize the parent class
101 ParentType::init(bulkProblem, lowDimProblem, couplingMapper, curSol);
102
103 // determine all bulk scvfs that coincide with low dim elements
104 bulkScvfIsOnFacetElement_.assign(bulkProblem->gridGeometry().numScvf(), false);
105 const auto& bulkMap = couplingMapper->couplingMap(bulkGridId, lowDimGridId);
106 for (const auto& entry : bulkMap)
107 for (const auto& couplingEntry : entry.second.elementToScvfMap)
108 for (const auto& scvfIdx : couplingEntry.second)
109 bulkScvfIsOnFacetElement_[scvfIdx] = true;
110
111 // store pointer to mapper
112 couplingMapperPtr_ = couplingMapper;
113 }
114
118 bool isOnInteriorBoundary(const Element<bulkId>& element,
119 const SubControlVolumeFace<bulkId>& scvf) const
120 { return bulkScvfIsOnFacetElement_[scvf.index()]; }
121
122 using ParentType::evalCouplingResidual;
129 template< class LowDimLocalAssembler >
130 typename LocalResidual<lowDimId>::ElementResidualVector
132 const LowDimLocalAssembler& lowDimLocalAssembler,
133 BulkIdType,
134 GridIndexType<bulkId> dofIdxGlobalJ)
135 { return evalCouplingResidual(lowDimId, lowDimLocalAssembler, bulkId); }
136
142 template< class LowDimLocalAssembler >
143 typename LocalResidual<lowDimId>::ElementResidualVector
144 evalCouplingResidual(LowDimIdType, const LowDimLocalAssembler& lowDimLocalAssembler, BulkIdType)
145 {
146 // make sure this is called for the element for which the context was set
147 assert(this->lowDimCouplingContext().isSet);
148 assert(this->problem(lowDimId).gridGeometry().elementMapper().index(lowDimLocalAssembler.element()) == this->lowDimCouplingContext().elementIdx);
149
150 // fill element residual vector with the sources
151 typename LowDimLocalAssembler::LocalResidual::ElementResidualVector res(lowDimLocalAssembler.fvGeometry().numScv());
152 res = 0.0;
153 for (const auto& scv : scvs(lowDimLocalAssembler.fvGeometry()))
154 res[scv.localDofIndex()] -= evalSourcesFromBulk(lowDimLocalAssembler.element(),
155 lowDimLocalAssembler.fvGeometry(),
156 lowDimLocalAssembler.curElemVolVars(),
157 scv);
158 return res;
159 }
160
164 NumEqVector<lowDimId> evalSourcesFromBulk(const Element<lowDimId>& element,
165 const FVElementGeometry<lowDimId>& fvGeometry,
166 const ElementVolumeVariables<lowDimId>& elemVolVars,
167 const SubControlVolume<lowDimId>& scv) const
168 {
169 // make sure the this is called for the element of the context
170 assert(this->problem(lowDimId).gridGeometry().elementMapper().index(element) == this->lowDimCouplingContext().elementIdx);
171
172 NumEqVector<lowDimId> sources(0.0);
173
174 const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId);
175 auto it = map.find(this->lowDimCouplingContext().elementIdx);
176 if (it == map.end())
177 return sources;
178
179 assert(this->lowDimCouplingContext().isSet);
180 for (const auto& embedment : it->second.embedments)
181 {
182 // list of scvfs in the bulk domain whose fluxes enter this scv
183 // if low dim domain uses a cc scheme, this is all scvfs lying on this element
184 // if it uses box, it is the one scvf coinciding with the given scv
185 const auto& coincidingScvfs = embedment.second;
186 const auto& scvfList = lowDimUsesBox ? std::vector<GridIndexType<lowDimId>>{ coincidingScvfs[scv.localDofIndex()] }
187 : coincidingScvfs;
188
189 sources += this->evalBulkFluxes(this->problem(bulkId).gridGeometry().element(embedment.first),
190 *this->lowDimCouplingContext().bulkFvGeometry,
191 *this->lowDimCouplingContext().bulkElemVolVars,
192 *this->lowDimCouplingContext().bulkElemFluxVarsCache,
193 *this->lowDimCouplingContext().bulkLocalResidual,
194 scvfList);
195 }
196
197 return sources;
198 }
199
204 template<class JacobianPattern>
205 void extendJacobianPattern(LowDimIdType, JacobianPattern& pattern) const
206 {
207 const auto& lowDimFVGridGeometry = this->problem(lowDimId).gridGeometry();
208 for (const auto& element : elements(lowDimFVGridGeometry.gridView()))
209 {
210
211 const auto eIdx = lowDimFVGridGeometry.elementMapper().index(element);
212 const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId);
213 auto it = map.find(eIdx);
214
215 // if element is coupled, take one of the neighbors and add coupling stencil to pattern
216 if (it != map.end())
217 {
218 // coupling stencil of the first neighbor
219 const auto bulkElemIdx = it->second.embedments[0].first;
220 const auto& bulkMapEntry = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId).at(bulkElemIdx);
221 const auto& couplingStencil = bulkMapEntry.couplingStencil;
222
223 for (auto globalJ : couplingStencil)
224 {
225 if (lowDimUsesBox)
226 {
227 for (int i = 0; i < element.subEntities(lowDimDim); ++i)
228 pattern.add(lowDimFVGridGeometry.vertexMapper().subIndex(element, i, lowDimDim), globalJ);
229 }
230 else
231 pattern.add(eIdx, globalJ);
232 }
233 }
234 }
235 }
236
238 template<class JacobianPattern>
239 void extendJacobianPattern(BulkIdType, JacobianPattern& pattern) const
240 {}
241
248 template<class LowDimLocalAssembler, class JacobianMatrixDiagBlock, class GridVariables>
250 const LowDimLocalAssembler& lowDimLocalAssembler,
251 const typename LowDimLocalAssembler::LocalResidual::ElementResidualVector&,
252 JacobianMatrixDiagBlock& A,
253 GridVariables& gridVariables)
254 {
255 // Since coupling only occurs via the fluxes, there are no
256 // additional derivatives for explicit time integration schemes
257 if (!LowDimLocalAssembler::isImplicit())
258 return;
259
260 // lambda to update the coupling context for a given lowDim element/dofIdx
261 auto updateContext = [&] (auto elemIdx, auto dofIdx, auto priVars, auto pvIdx)
262 {
263 // deflect the solution
264 auto& ldSol = this->curSol(lowDimId);
265 ldSol[dofIdx][pvIdx] = priVars[pvIdx];
266
267 // update the corresponding vol vars in the bulk context
268 assert(this->bulkCouplingContext().isSet);
269 const auto& bulkMap = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
270 const auto& couplingElementStencil = bulkMap.find(this->bulkCouplingContext().elementIdx)->second.couplingElementStencil;
271
272 auto it = std::find(couplingElementStencil.begin(), couplingElementStencil.end(), elemIdx);
273 assert(it != couplingElementStencil.end());
274 const auto idxInContext = std::distance(couplingElementStencil.begin(), it);
275
276 auto& volVars = this->bulkCouplingContext().lowDimVolVars[idxInContext];
277 const auto& fvGeom = this->bulkCouplingContext().lowDimFvGeometries[idxInContext];
278 const auto& element = this->problem(lowDimId).gridGeometry().element(elemIdx);
279
280 // if low dim domain uses the box scheme, we have to create interpolated vol vars
281 if (lowDimUsesBox)
282 {
283 const auto elemGeom = element.geometry();
284 FacetCoupling::makeInterpolatedVolVars(volVars, this->problem(lowDimId), ldSol, fvGeom, element, elemGeom, elemGeom.center());
285 }
286 // if low dim domain uses a cc scheme we can directly update the vol vars
287 else
288 volVars.update( elementSolution(element, ldSol, this->problem(lowDimId).gridGeometry()),
289 this->problem(lowDimId),
290 element,
291 fvGeom.scv(elemIdx) );
292
293 // update the element flux variables cache (tij depend on low dim values in context)
294 const auto contextElem = this->problem(bulkId).gridGeometry().element(this->bulkCouplingContext().elementIdx);
295 this->lowDimCouplingContext().bulkElemFluxVarsCache->update(contextElem,
296 *this->lowDimCouplingContext().bulkFvGeometry,
297 *this->lowDimCouplingContext().bulkElemVolVars);
298 };
299
300 const auto eIdx = this->problem(lowDimId).gridGeometry().elementMapper().index(lowDimLocalAssembler.element());
301
302 // bug tracking
303 assert(this->lowDimCouplingContext().isSet);
304 assert(this->lowDimCouplingContext().elementIdx == eIdx);
305
306 // if the element is coupled, evaluate additional source derivatives
307 const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId);
308 auto it = map.find(eIdx);
309 if (it != map.end())
310 evalLowDimSourceDerivatives_(updateContext, lowDimLocalAssembler, A);
311 }
312
314 template<class LocalAssemblerI, class JacobianMatrixDiagBlock, class GridVariables>
316 const LocalAssemblerI& localAssemblerI,
317 const typename LocalAssemblerI::LocalResidual::ElementResidualVector& origResiduals,
318 JacobianMatrixDiagBlock& A,
319 GridVariables& gridVariables)
320 {}
321
322private:
324 template<class UpdateContext, class LowDimLocalAssembler, class JacobianMatrixDiagBlock>
325 void evalLowDimSourceDerivatives_(const UpdateContext& updateContext,
326 const LowDimLocalAssembler& lowDimLocalAssembler,
327 JacobianMatrixDiagBlock& A)
328 {
329 const auto& lowDimFVGridGeometry = this->problem(lowDimId).gridGeometry();
330 const auto eIdx = lowDimFVGridGeometry.elementMapper().index(lowDimLocalAssembler.element());
331
332 // coupling stencil of the first neighbor
333 const auto bulkElemIdx = this->bulkCouplingContext().elementIdx;
334 const auto& bulkMapEntry = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId).at(bulkElemIdx);
335 const auto& couplingStencil = bulkMapEntry.couplingStencil;
336 const auto& couplingElementStencil = bulkMapEntry.couplingElementStencil;
337
338 // compute the undeflected residual (reuse coupling residual function)
339 const auto origResidual = evalCouplingResidual(lowDimId, lowDimLocalAssembler, bulkId);
340
341 // container of dofs within this element
342 std::vector< std::decay_t<decltype(couplingStencil[0])> > elemDofs;
343 elemDofs.reserve(lowDimLocalAssembler.fvGeometry().numScv());
344 for (const auto& scv : scvs(lowDimLocalAssembler.fvGeometry()))
345 elemDofs.push_back(scv.dofIndex());
346
347 // compute derivate for all additional dofs in the stencil
348 for (const auto couplingElemIdx : couplingElementStencil)
349 {
350 // skip the same element
351 if (couplingElemIdx == eIdx)
352 continue;
353
354 // container of dofs within the other element
355 std::vector< std::decay_t<decltype(couplingStencil[0])> > elemDofsJ;
356 if (lowDimUsesBox)
357 {
358 const auto& elemJ = lowDimFVGridGeometry.element(couplingElemIdx);
359 for (int i = 0; i < elemJ.subEntities(lowDimDim); ++i)
360 elemDofsJ.push_back(lowDimFVGridGeometry.vertexMapper().subIndex(elemJ, i, lowDimDim));
361 }
362 else
363 elemDofsJ.push_back(couplingElemIdx);
364
365 for (auto dofIndex : elemDofsJ)
366 {
367 auto partialDerivs = origResidual;
368 const auto origPriVars = this->curSol(lowDimId)[dofIndex];
369
370 // calculate derivatives w.r.t to the privars at the dof at hand
371 static constexpr auto numEq = std::decay_t<decltype(origPriVars)>::dimension;
372 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
373 {
374 // reset partial derivatives
375 partialDerivs = 0.0;
376
377 auto evalResiduals = [&](Scalar<lowDimId> priVar)
378 {
379 auto priVars = origPriVars;
380 priVars[pvIdx] = priVar;
381
382 // Update context to deflected solution and reevaluate residual
383 updateContext(couplingElemIdx, dofIndex, priVars, pvIdx);
384 return this->evalCouplingResidual(lowDimId, lowDimLocalAssembler, bulkId);
385 };
386
387 static const int numDiffMethod = getParamFromGroup<int>(this->problem(lowDimId).paramGroup(), "Assembly.NumericDifferenceMethod");
388 static const NumericEpsilon< Scalar<lowDimId>, numEq > eps{this->problem(lowDimId).paramGroup()};
389 NumericDifferentiation::partialDerivative(evalResiduals, origPriVars[pvIdx], partialDerivs,
390 origResidual, eps(origPriVars[pvIdx], pvIdx), numDiffMethod);
391
392 // update the global stiffness matrix with the current partial derivatives
393 // A[i][col][eqIdx][pvIdx] is the rate of change of the residual of equation
394 // 'eqIdx' at dof 'i' depending on the primary variable 'pvIdx' at dof 'col'.
395 for (const auto& scv : scvs(lowDimLocalAssembler.fvGeometry()))
396 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
397 A[scv.dofIndex()][dofIndex][eqIdx][pvIdx] += partialDerivs[scv.indexInElement()][eqIdx];
398
399 // restore the original coupling context
400 updateContext(couplingElemIdx, dofIndex, origPriVars, pvIdx);
401 }
402 }
403 }
404 }
405
407 std::vector<bool> bulkScvfIsOnFacetElement_;
408
410 std::shared_ptr< CouplingMapper > couplingMapperPtr_;
411};
412
413} // end namespace Dumux
414
415#endif
NumEqVector< lowDimId > evalSourcesFromBulk(const Element< lowDimId > &element, const FVElementGeometry< lowDimId > &fvGeometry, const ElementVolumeVariables< lowDimId > &elemVolVars, const SubControlVolume< lowDimId > &scv) const
Computes the sources in a lower-dimensional sub-control volume stemming from the bulk domain.
Definition: multidomain/facet/cellcentered/mpfa/couplingmanager.hh:164
bool isOnInteriorBoundary(const Element< bulkId > &element, const SubControlVolumeFace< bulkId > &scvf) const
returns true if a bulk scvf coincides with a facet element.
Definition: multidomain/facet/cellcentered/mpfa/couplingmanager.hh:118
LocalResidual< lowDimId >::ElementResidualVector evalCouplingResidual(LowDimIdType, const LowDimLocalAssembler &lowDimLocalAssembler, BulkIdType)
Evaluates the coupling element residual of a lower-dimensional domain element with respect to a dof i...
Definition: multidomain/facet/cellcentered/mpfa/couplingmanager.hh:144
void extendJacobianPattern(LowDimIdType, JacobianPattern &pattern) const
Extend the jacobian pattern of the diagonal block of the lowdim domain by the elements that are in th...
Definition: multidomain/facet/cellcentered/mpfa/couplingmanager.hh:205
typename MDTraits::SolutionVector SolutionVector
the type of the solution vector
Definition: multidomain/facet/cellcentered/mpfa/couplingmanager.hh:85
LocalResidual< lowDimId >::ElementResidualVector evalCouplingResidual(LowDimIdType, const LowDimLocalAssembler &lowDimLocalAssembler, BulkIdType, GridIndexType< bulkId > dofIdxGlobalJ)
Evaluates the coupling element residual of a lower-dimensional domain element with respect to a dof i...
Definition: multidomain/facet/cellcentered/mpfa/couplingmanager.hh:131
void evalAdditionalDomainDerivatives(BulkIdType, const LocalAssemblerI &localAssemblerI, const typename LocalAssemblerI::LocalResidual::ElementResidualVector &origResiduals, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
The bulk domain has no additional derivatives.
Definition: multidomain/facet/cellcentered/mpfa/couplingmanager.hh:315
void init(std::shared_ptr< Problem< bulkId > > bulkProblem, std::shared_ptr< Problem< lowDimId > > lowDimProblem, std::shared_ptr< CouplingMapper > couplingMapper, const SolutionVector &curSol)
Initialize the coupling manager.
Definition: multidomain/facet/cellcentered/mpfa/couplingmanager.hh:95
void evalAdditionalDomainDerivatives(LowDimIdType, const LowDimLocalAssembler &lowDimLocalAssembler, const typename LowDimLocalAssembler::LocalResidual::ElementResidualVector &, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
evaluate additional derivatives of the element residual of the low-dim domain with respect to dofs in...
Definition: multidomain/facet/cellcentered/mpfa/couplingmanager.hh:249
void extendJacobianPattern(BulkIdType, JacobianPattern &pattern) const
The bulk domain has no extended jacobian pattern.
Definition: multidomain/facet/cellcentered/mpfa/couplingmanager.hh:239
Manages the coupling between bulk elements and lower dimensional elements where the coupling occurs a...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:42
Implementation for the coupling manager between two domains of dimension d and (d-1) for models consi...
Definition: multidomain/facet/couplingmanager.hh:83
static void partialDerivative(const Function &function, Scalar x0, FunctionEvalType &derivative, const FunctionEvalType &fx0, const int numericDifferenceMethod=1)
Computes the derivative of a function with respect to a function parameter.
Definition: numericdifferentiation.hh:50
Defines all properties used in Dumux.
Element solution classes and factory functions.
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:34
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
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:282
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.
The interface of the coupling manager for multi domain problems.
constexpr Box box
Definition: method.hh:147
Definition: adapt.hh:17
A helper to deduce a vector with the same size as numbers of equations.
A class for numeric differentiation.
An adapter class for local assemblers using numeric differentiation.
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:26