14#ifndef DUMUX_MULTIDOMAIN_FACECENTERED_LOCAL_ASSEMBLER_HH
15#define DUMUX_MULTIDOMAIN_FACECENTERED_LOCAL_ASSEMBLER_HH
17#include <dune/common/indices.hh>
18#include <dune/common/hybridutilities.hh>
19#include <dune/grid/common/gridenums.hh>
41template<std::
size_t id,
class TypeTag,
class Assembler,
class Implementation, DiffMethod dm,
bool implicit>
47 using SolutionVector =
typename Assembler::SolutionVector;
50 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
51 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
52 using ElementFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
53 using Scalar =
typename GridVariables::Scalar;
55 using GridGeometry =
typename GridVariables::GridGeometry;
56 using FVElementGeometry =
typename GridGeometry::LocalView;
57 using SubControlVolume =
typename GridGeometry::SubControlVolume;
58 using GridView =
typename GridGeometry::GridView;
59 using Element =
typename GridView::template Codim<0>::Entity;
61 using CouplingManager =
typename Assembler::CouplingManager;
67 static constexpr auto domainId =
typename Dune::index_constant<id>();
69 using ParentType::ParentType;
75 const Element& element,
76 const SolutionVector&
curSol,
86 (element.partitionType() ==
Dune::GhostEntity))
94 template<
class JacobianMatrixRow,
class SubRes
idualVector,
class Gr
idVariablesTuple>
97 auto assembleCouplingBlocks = [&](
const auto& residual)
100 using namespace Dune::Hybrid;
101 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](
auto&& i)
110 ParentType::assembleJacobianAndResidual(jacRow[
domainId], res, *std::get<domainId>(gridVariables), noReassembler, assembleCouplingBlocks);
117 template<std::size_t otherId,
class JacRow,
class GridVariables,
118 typename std::enable_if_t<(otherId == id),
int> = 0>
126 template<std::size_t otherId,
class JacRow,
class GridVariables,
127 typename std::enable_if_t<(otherId != id),
int> = 0>
131 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
144 for (
auto&& scv : scvs(this->fvGeometry()))
146 const auto& curVolVars = elemVolVars[scv];
147 auto source = this->localResidual().computeSource(
problem(), element, this->fvGeometry(), elemVolVars, scv);
148 source *= -scv.volume()*curVolVars.extrusionFactor();
149 residual[scv.localDofIndex()] = std::move(source);
159 {
return this->
evalLocalSourceResidual(neighbor, implicit ? this->curElemVolVars() : this->prevElemVolVars()); }
167 const auto& element = this->element();
169 auto&& fvGeometry = this->fvGeometry();
170 auto&& curElemVolVars = this->curElemVolVars();
171 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
174 couplingManager_.bindCouplingContext(
domainId, element, this->assembler());
175 fvGeometry.bind(element);
179 curElemVolVars.bind(element, fvGeometry,
curSol);
180 elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars);
181 if (!this->assembler().isStationaryProblem())
182 this->prevElemVolVars().bindElement(element, fvGeometry, this->assembler().prevSol()[
domainId]);
186 auto& prevElemVolVars = this->prevElemVolVars();
187 const auto& prevSol = this->assembler().prevSol()[
domainId];
189 curElemVolVars.bindElement(element, fvGeometry,
curSol);
190 prevElemVolVars.bind(element, fvGeometry, prevSol);
191 elemFluxVarsCache.bind(element, fvGeometry, prevElemVolVars);
194 this->elemBcTypes().update(
problem(), this->element(), this->fvGeometry());
198 template<std::
size_t i = domainId>
200 {
return this->assembler().problem(
domainId); }
203 template<std::
size_t i = domainId>
205 {
return ParentType::curSol()[dId]; }
209 {
return couplingManager_; }
226template<std::
size_t id,
class TypeTag,
class Assembler, DiffMethod DM = DiffMethod::numeric,
bool implicit = true>
235template<std::
size_t id,
class TypeTag,
class Assembler>
238 id, TypeTag, Assembler,
249 using FVElementGeometry =
typename GridGeometry::LocalView;
250 using GridView =
typename GridGeometry::GridView;
251 using Element =
typename GridView::template Codim<0>::Entity;
252 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
255 enum { dim = GridView::dimension };
259 static constexpr auto domainI = Dune::index_constant<id>();
269 template<
class ElemSol>
272 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, scv.dofIndex(), elemSol[scv.localDofIndex()], pvIdx);
278 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
281 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *
this, origResiduals, A, gridVariables);
288 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
297 const auto& element = this->element();
298 const auto& fvGeometry = this->fvGeometry();
299 auto&& curElemVolVars = this->curElemVolVars();
300 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
303 const auto& curSolJ = this->curSol(domainJ);
306 auto updateCoupledVariables = [&] ()
310 if constexpr (enableGridFluxVarsCache)
312 if constexpr (enableGridVolVarsCache)
313 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
315 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, gridVariables.gridFluxVarsCache());
319 if constexpr (enableGridVolVarsCache)
320 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), elemFluxVarsCache);
322 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, elemFluxVarsCache);
326 for (
const auto& scv : scvs(fvGeometry))
328 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, scv, domainJ);
330 for (
const auto globalJ : stencil)
332 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *
this, scv, domainJ, globalJ);
334 const auto origPriVarsJ = curSolJ[globalJ];
335 auto priVarsJ = origPriVarsJ;
337 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
339 auto evalCouplingResidual = [&](Scalar priVar)
341 priVarsJ[pvIdx] = priVar;
342 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
343 updateCoupledVariables();
344 return this->couplingManager().evalCouplingResidual(domainI, *
this, scv, domainJ, globalJ);
350 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
351 static const int numDiffMethod = getParamFromGroup<int>(paramGroup,
"Assembly.NumericDifferenceMethod");
352 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
355 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
358 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
364 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
369 if (scv.boundary() && this->elemBcTypes().hasDirichlet())
371 const auto bcTypes = this->elemBcTypes()[fvGeometry.frontalScvfOnBoundary(scv).localIndex()];
372 if (bcTypes.hasDirichlet())
377 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
379 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
381 if (bcTypes.isCouplingDirichlet(eqIdx))
382 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][pvIdx];
383 else if (bcTypes.isDirichlet(eqIdx))
384 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
391 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
394 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
401 updateCoupledVariables();
413template<std::
size_t id,
class TypeTag,
class Assembler>
416 id, TypeTag, Assembler,
417 SubDomainFaceCenteredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>,
418 DiffMethod::numeric, false
424 using FVElementGeometry =
typename GridGeometry::LocalView;
425 using GridView =
typename GridGeometry::GridView;
426 using Element =
typename GridView::template Codim<0>::Entity;
427 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
428 static constexpr auto domainI = Dune::index_constant<id>();
437 template<
class ElemSol>
440 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, scv.dofIndex(), elemSol[scv.localDofIndex()], pvIdx);
446 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
454 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:37
Definition: partialreassembler.hh:32
An assembler for Jacobian and residual contribution per element (Face-centered methods)
Definition: fclocalassembler.hh:284
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
Face-centered staggered scheme multi domain local assembler using numeric differentiation and implici...
Definition: subdomainfclocalassembler.hh:242
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const ElementResidualVector &res, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainfclocalassembler.hh:289
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: subdomainfclocalassembler.hh:264
void maybeUpdateCouplingContext(const SubControlVolume &scv, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: subdomainfclocalassembler.hh:270
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, const JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: subdomainfclocalassembler.hh:279
Face-centered staggered scheme multi domain local assembler using numeric differentiation and explici...
Definition: subdomainfclocalassembler.hh:420
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: subdomainfclocalassembler.hh:432
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const ElementResidualVector &res, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainfclocalassembler.hh:455
void maybeUpdateCouplingContext(const SubControlVolume &scv, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: subdomainfclocalassembler.hh:438
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, const JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: subdomainfclocalassembler.hh:447
A base class for all face-centered staggered local assemblers.
Definition: subdomainfclocalassembler.hh:43
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition: subdomainfclocalassembler.hh:158
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: subdomainfclocalassembler.hh:208
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: subdomainfclocalassembler.hh:164
const auto & curSol(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: subdomainfclocalassembler.hh:204
const Problem & problem(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: subdomainfclocalassembler.hh:199
void assembleJacobianAndResidual(JacobianMatrixRow &jacRow, SubResidualVector &res, GridVariablesTuple &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: subdomainfclocalassembler.hh:95
ElementResidualVector evalLocalSourceResidual(const Element &element, const ElementVolumeVariables &elemVolVars) const
Evaluates the local source term for an element and given element volume variables.
Definition: subdomainfclocalassembler.hh:137
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacRow &jacRow, const ElementResidualVector &res, GridVariables &gridVariables)
Assemble the entries in a coupling block of the jacobian. There is no coupling block between a domain...
Definition: subdomainfclocalassembler.hh:119
SubDomainFaceCenteredLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: subdomainfclocalassembler.hh:74
static constexpr auto domainId
export the domain id of this sub-domain
Definition: subdomainfclocalassembler.hh:67
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: subdomainfclocalassembler.hh:71
The face-centered staggered scheme multidomain local assembler.
Definition: subdomainfclocalassembler.hh:227
Defines all properties used in Dumux.
An enum class to define various differentiation methods available in order to compute the derivatives...
An assembler for Jacobian and residual contribution per element (face-centered staggered methods)
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition: diffmethod.hh:25
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Definition: common/pdesolver.hh:24
A class for numeric differentiation.
An adapter class for local assemblers using numeric differentiation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.