24#ifndef DUMUX_CCMPFA_FACETCOUPLING_MANAGER_HH
25#define DUMUX_CCMPFA_FACETCOUPLING_MANAGER_HH
54template<
class MDTraits,
class CouplingMapper, std::
size_t bulkDomainId, std::
size_t lowDimDomainId>
55class FacetCouplingManager<MDTraits, CouplingMapper, bulkDomainId, lowDimDomainId, DiscretizationMethods::CCMpfa>
56:
public FacetCouplingManager<MDTraits, CouplingMapper, bulkDomainId, lowDimDomainId, DiscretizationMethods::CCTpfa>
61 using BulkIdType =
typename MDTraits::template SubDomain<bulkDomainId>::Index;
62 using LowDimIdType =
typename MDTraits::template SubDomain<lowDimDomainId>::Index;
63 static constexpr auto bulkId = BulkIdType();
64 static constexpr auto lowDimId = LowDimIdType();
67 template<std::
size_t id>
using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
74 template<std::
size_t id>
using FVElementGeometry =
typename GridGeometry<id>::LocalView;
75 template<std::
size_t id>
using SubControlVolume =
typename GridGeometry<id>::SubControlVolume;
76 template<std::
size_t id>
using SubControlVolumeFace =
typename GridGeometry<id>::SubControlVolumeFace;
77 template<std::
size_t id>
using GridView =
typename GridGeometry<id>::GridView;
79 template<std::
size_t id>
using Element =
typename GridView<id>::template Codim<0>::Entity;
83 template<std::
size_t id>
using GridVolumeVariables =
typename GridVariables<id>::GridVolumeVariables;
84 template<std::
size_t id>
using ElementVolumeVariables =
typename GridVolumeVariables<id>::LocalView;
87 static constexpr int bulkDim = GridView<bulkDomainId>::dimension;
88 static constexpr int lowDimDim = GridView<lowDimDomainId>::dimension;
89 static constexpr auto bulkGridId = CouplingMapper::template gridId<bulkDim>();
90 static constexpr auto lowDimGridId = CouplingMapper::template gridId<lowDimDim>();
107 void init(std::shared_ptr< Problem<bulkId> > bulkProblem,
108 std::shared_ptr< Problem<lowDimId> > lowDimProblem,
109 std::shared_ptr< CouplingMapper > couplingMapper,
113 ParentType::init(bulkProblem, lowDimProblem, couplingMapper, curSol);
116 bulkScvfIsOnFacetElement_.assign(bulkProblem->gridGeometry().numScvf(),
false);
117 const auto& bulkMap = couplingMapper->couplingMap(bulkGridId, lowDimGridId);
118 for (
const auto& entry : bulkMap)
119 for (
const auto& couplingEntry : entry.second.elementToScvfMap)
120 for (
const auto& scvfIdx : couplingEntry.second)
121 bulkScvfIsOnFacetElement_[scvfIdx] =
true;
124 couplingMapperPtr_ = couplingMapper;
131 const SubControlVolumeFace<bulkId>& scvf)
const
132 {
return bulkScvfIsOnFacetElement_[scvf.index()]; }
134 using ParentType::evalCouplingResidual;
141 template<
class LowDimLocalAssembler >
142 typename LocalResidual<lowDimId>::ElementResidualVector
144 const LowDimLocalAssembler& lowDimLocalAssembler,
146 GridIndexType<bulkId> dofIdxGlobalJ)
147 {
return evalCouplingResidual(lowDimId, lowDimLocalAssembler, bulkId); }
154 template<
class LowDimLocalAssembler >
155 typename LocalResidual<lowDimId>::ElementResidualVector
159 assert(this->lowDimCouplingContext().isSet);
160 assert(this->problem(lowDimId).gridGeometry().elementMapper().index(lowDimLocalAssembler.element()) == this->lowDimCouplingContext().elementIdx);
163 typename LowDimLocalAssembler::LocalResidual::ElementResidualVector res(lowDimLocalAssembler.fvGeometry().numScv());
165 for (
const auto& scv : scvs(lowDimLocalAssembler.fvGeometry()))
166 res[scv.localDofIndex()] -= evalSourcesFromBulk(lowDimLocalAssembler.element(),
167 lowDimLocalAssembler.fvGeometry(),
168 lowDimLocalAssembler.curElemVolVars(),
177 const FVElementGeometry<lowDimId>& fvGeometry,
178 const ElementVolumeVariables<lowDimId>& elemVolVars,
179 const SubControlVolume<lowDimId>& scv)
182 assert(this->problem(lowDimId).gridGeometry().elementMapper().index(element) == this->lowDimCouplingContext().elementIdx);
184 NumEqVector<lowDimId> sources(0.0);
186 const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId);
187 auto it = map.find(this->lowDimCouplingContext().elementIdx);
191 assert(this->lowDimCouplingContext().isSet);
192 for (
const auto& embedment : it->second.embedments)
197 const auto& coincidingScvfs = embedment.second;
198 const auto& scvfList = lowDimUsesBox ? std::vector<GridIndexType<lowDimId>>{ coincidingScvfs[scv.localDofIndex()] }
201 sources += this->evalBulkFluxes(this->problem(bulkId).gridGeometry().element(embedment.first),
202 *this->lowDimCouplingContext().bulkFvGeometry,
203 *this->lowDimCouplingContext().bulkElemVolVars,
204 *this->lowDimCouplingContext().bulkElemFluxVarsCache,
205 *this->lowDimCouplingContext().bulkLocalResidual,
216 template<
class JacobianPattern>
219 const auto& lowDimFVGridGeometry = this->problem(lowDimId).gridGeometry();
220 for (
const auto& element : elements(lowDimFVGridGeometry.gridView()))
223 const auto eIdx = lowDimFVGridGeometry.elementMapper().index(element);
224 const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId);
225 auto it = map.find(eIdx);
231 const auto bulkElemIdx = it->second.embedments[0].first;
232 const auto& bulkMapEntry = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId).at(bulkElemIdx);
233 const auto& couplingStencil = bulkMapEntry.couplingStencil;
235 for (
auto globalJ : couplingStencil)
239 for (
int i = 0; i < element.subEntities(lowDimDim); ++i)
240 pattern.add(lowDimFVGridGeometry.vertexMapper().subIndex(element, i, lowDimDim), globalJ);
243 pattern.add(eIdx, globalJ);
250 template<
class JacobianPattern>
260 template<
class LowDimLocalAssembler,
class JacobianMatrixDiagBlock,
class Gr
idVariables>
262 const LowDimLocalAssembler& lowDimLocalAssembler,
263 const typename LowDimLocalAssembler::LocalResidual::ElementResidualVector&,
264 JacobianMatrixDiagBlock& A,
265 GridVariables& gridVariables)
269 if (!LowDimLocalAssembler::isImplicit())
273 auto updateContext = [&] (
auto elemIdx,
auto dofIdx,
auto priVars,
auto pvIdx)
276 auto& ldSol = this->curSol(lowDimId);
277 ldSol[dofIdx][pvIdx] = priVars[pvIdx];
280 assert(this->bulkCouplingContext().isSet);
281 const auto& bulkMap = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
282 const auto& couplingElementStencil = bulkMap.find(this->bulkCouplingContext().elementIdx)->second.couplingElementStencil;
284 auto it = std::find(couplingElementStencil.begin(), couplingElementStencil.end(), elemIdx);
285 assert(it != couplingElementStencil.end());
286 const auto idxInContext =
std::distance(couplingElementStencil.begin(), it);
288 auto& volVars = this->bulkCouplingContext().lowDimVolVars[idxInContext];
289 const auto& fvGeom = this->bulkCouplingContext().lowDimFvGeometries[idxInContext];
290 const auto& element = this->problem(lowDimId).gridGeometry().element(elemIdx);
295 const auto elemGeom = element.geometry();
300 volVars.update(
elementSolution(element, ldSol, this->problem(lowDimId).gridGeometry()),
301 this->problem(lowDimId),
303 fvGeom.scv(elemIdx) );
306 const auto contextElem = this->problem(bulkId).gridGeometry().element(this->bulkCouplingContext().elementIdx);
307 this->lowDimCouplingContext().bulkElemFluxVarsCache->update(contextElem,
308 *this->lowDimCouplingContext().bulkFvGeometry,
309 *this->lowDimCouplingContext().bulkElemVolVars);
312 const auto eIdx = this->problem(lowDimId).gridGeometry().elementMapper().index(lowDimLocalAssembler.element());
315 assert(this->lowDimCouplingContext().isSet);
316 assert(this->lowDimCouplingContext().elementIdx == eIdx);
319 const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId);
320 auto it = map.find(eIdx);
322 evalLowDimSourceDerivatives_(updateContext, lowDimLocalAssembler, A);
326 template<
class LocalAssemblerI,
class JacobianMatrixDiagBlock,
class Gr
idVariables>
328 const LocalAssemblerI& localAssemblerI,
329 const typename LocalAssemblerI::LocalResidual::ElementResidualVector& origResiduals,
330 JacobianMatrixDiagBlock& A,
331 GridVariables& gridVariables)
336 template<
class UpdateContext,
class LowDimLocalAssembler,
class JacobianMatrixDiagBlock>
337 void evalLowDimSourceDerivatives_(
const UpdateContext& updateContext,
338 const LowDimLocalAssembler& lowDimLocalAssembler,
339 JacobianMatrixDiagBlock& A)
341 const auto& lowDimFVGridGeometry = this->problem(lowDimId).gridGeometry();
342 const auto eIdx = lowDimFVGridGeometry.elementMapper().index(lowDimLocalAssembler.element());
345 const auto bulkElemIdx = this->bulkCouplingContext().elementIdx;
346 const auto& bulkMapEntry = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId).at(bulkElemIdx);
347 const auto& couplingStencil = bulkMapEntry.couplingStencil;
348 const auto& couplingElementStencil = bulkMapEntry.couplingElementStencil;
351 const auto origResidual = evalCouplingResidual(lowDimId, lowDimLocalAssembler, bulkId);
354 std::vector< std::decay_t<
decltype(couplingStencil[0])> > elemDofs;
355 elemDofs.reserve(lowDimLocalAssembler.fvGeometry().numScv());
356 for (
const auto& scv : scvs(lowDimLocalAssembler.fvGeometry()))
357 elemDofs.push_back(scv.dofIndex());
360 for (
const auto couplingElemIdx : couplingElementStencil)
363 if (couplingElemIdx == eIdx)
367 std::vector< std::decay_t<
decltype(couplingStencil[0])> > elemDofsJ;
370 const auto& elemJ = lowDimFVGridGeometry.element(couplingElemIdx);
371 for (
int i = 0; i < elemJ.subEntities(lowDimDim); ++i)
372 elemDofsJ.push_back(lowDimFVGridGeometry.vertexMapper().subIndex(elemJ, i, lowDimDim));
375 elemDofsJ.push_back(couplingElemIdx);
377 for (
auto dofIndex : elemDofsJ)
379 auto partialDerivs = origResidual;
380 const auto origPriVars = this->curSol(lowDimId)[dofIndex];
383 static constexpr auto numEq = std::decay_t<
decltype(origPriVars)>::dimension;
384 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
389 auto evalResiduals = [&](Scalar<lowDimId> priVar)
391 auto priVars = origPriVars;
392 priVars[pvIdx] = priVar;
395 updateContext(couplingElemIdx, dofIndex, priVars, pvIdx);
396 return this->evalCouplingResidual(lowDimId, lowDimLocalAssembler, bulkId);
399 static const int numDiffMethod = getParamFromGroup<int>(this->problem(lowDimId).paramGroup(),
"Assembly.NumericDifferenceMethod");
400 static const NumericEpsilon< Scalar<lowDimId>, numEq > eps{this->problem(lowDimId).paramGroup()};
402 origResidual, eps(origPriVars[pvIdx], pvIdx), numDiffMethod);
407 for (
const auto& scv : scvs(lowDimLocalAssembler.fvGeometry()))
408 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
409 A[scv.dofIndex()][dofIndex][eqIdx][pvIdx] += partialDerivs[scv.indexInElement()][eqIdx];
412 updateContext(couplingElemIdx, dofIndex, origPriVars, pvIdx);
419 std::vector<bool> bulkScvfIsOnFacetElement_;
422 std::shared_ptr< CouplingMapper > couplingMapperPtr_;
A helper to deduce a vector with the same size as numbers of equations.
A class for numeric differentiation.
Defines the index types used for grid and local indices.
An adapter class for local assemblers using numeric differentiation.
The available discretization methods in Dumux.
Element solution classes and factory functions.
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:294
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
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:46
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
constexpr Box box
Definition: method.hh:136
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:38
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:61
Property to specify the type of scalar values.
Definition: common/properties.hh:43
A vector of primary variables.
Definition: common/properties.hh:49
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:55
Definition: common/properties.hh:72
Definition: common/properties.hh:100
The grid variables object managing variable data on the grid (volvars/fluxvars cache)
Definition: common/properties.hh:121
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:130
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:156
NumEqVector< lowDimId > evalSourcesFromBulk(const Element< lowDimId > &element, const FVElementGeometry< lowDimId > &fvGeometry, const ElementVolumeVariables< lowDimId > &elemVolVars, const SubControlVolume< lowDimId > &scv)
Computes the sources in a lower-dimensional sub-control volume stemming from the bulk domain.
Definition: multidomain/facet/cellcentered/mpfa/couplingmanager.hh:176
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:217
typename MDTraits::SolutionVector SolutionVector
the type of the solution vector
Definition: multidomain/facet/cellcentered/mpfa/couplingmanager.hh:97
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:143
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:327
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:107
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:261
void extendJacobianPattern(BulkIdType, JacobianPattern &pattern) const
The bulk domain has no extended jacobian pattern.
Definition: multidomain/facet/cellcentered/mpfa/couplingmanager.hh:251
Manages the coupling between bulk elements and lower dimensional elements where the coupling occurs a...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:54
Implementation for the coupling manager between two domains of dimension d and (d-1) for models consi...
Definition: multidomain/facet/couplingmanager.hh:95
The interface of the coupling manager for multi domain problems.
Declares all properties used in Dumux.