14#ifndef DUMUX_MULTIDOMAIN_SUBDOMAIN_CVFE_LOCAL_ASSEMBLER_HH
15#define DUMUX_MULTIDOMAIN_SUBDOMAIN_CVFE_LOCAL_ASSEMBLER_HH
19#include <dune/common/reservedvector.hh>
20#include <dune/common/indices.hh>
21#include <dune/common/hybridutilities.hh>
22#include <dune/grid/common/gridenums.hh>
23#include <dune/istl/matrixindexset.hh>
25#include <dumux/common/typetraits/localdofs_.hh>
26#include <dumux/common/typetraits/boundary_.hh>
51template<std::
size_t id,
class TypeTag,
class Assembler,
class Implementation, DiffMethod dm,
bool implicit>
59 using SolutionVector =
typename Assembler::SolutionVector;
63 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
64 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
65 using ElementFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
66 using Scalar =
typename GridVariables::Scalar;
68 using GridGeometry =
typename GridVariables::GridGeometry;
69 using FVElementGeometry =
typename GridGeometry::LocalView;
71 using GridView =
typename GridGeometry::GridView;
72 using Element =
typename GridView::template Codim<0>::Entity;
74 using CouplingManager =
typename Assembler::CouplingManager;
80 static constexpr auto domainId =
typename Dune::index_constant<id>();
82 using ParentType::ParentType;
88 const Assembler& assembler,
89 const Element& element,
90 const SolutionVector&
curSol,
101 (element.partitionType() ==
Dune::GhostEntity))
109 template<
class JacobianMatrixRow,
class SubRes
idualVector,
class Gr
idVariablesTuple>
112 auto assembleCouplingBlocks = [&](
const auto& residual)
115 using namespace Dune::Hybrid;
116 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](
auto&& i)
118 if constexpr (std::decay_t<
decltype(i)>{} != id)
125 ParentType::assembleJacobianAndResidual(jacRow[
domainId], res, *std::get<domainId>(gridVariables), noReassembler, assembleCouplingBlocks);
132 template<std::size_t otherId,
class JacRow,
class GridVariables,
133 typename std::enable_if_t<(otherId == id),
int> = 0>
141 template<std::size_t otherId,
class JacRow,
class GridVariables,
142 typename std::enable_if_t<(otherId != id),
int> = 0>
146 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
154 static_assert(!Detail::LocalDofs::hasNonCVLocalDofsInterface<FVElementGeometry>(),
"Separate source calculation not implemented for hybrid schemes.");
161 for (
const auto& scv :
scvs(this->fvGeometry()))
163 const auto& curVolVars = elemVolVars[scv];
164 auto source = this->localResidual().computeSource(
problem(), element, this->fvGeometry(), elemVolVars, scv);
165 source *= -
Extrusion::volume(this->fvGeometry(), scv)*curVolVars.extrusionFactor();
166 residual[scv.localDofIndex()] = std::move(source);
176 {
return this->
evalLocalSourceResidual(neighbor, implicit ? this->curElemVolVars() : this->prevElemVolVars()); }
184 const auto& element = this->element();
186 auto&& fvGeometry = this->fvGeometry();
187 auto&& curElemVolVars = this->curElemVolVars();
188 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
191 couplingManager_.bindCouplingContext(
domainId, element, this->assembler());
192 fvGeometry.bind(element);
194 if constexpr (implicit)
196 curElemVolVars.bind(element, fvGeometry,
curSol);
197 elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars);
198 if (!this->assembler().isStationaryProblem())
199 this->prevElemVolVars().bindElement(element, fvGeometry, this->assembler().prevSol()[
domainId]);
203 auto& prevElemVolVars = this->prevElemVolVars();
204 const auto& prevSol = this->assembler().prevSol()[
domainId];
206 curElemVolVars.bindElement(element, fvGeometry,
curSol);
207 prevElemVolVars.bind(element, fvGeometry, prevSol);
208 elemFluxVarsCache.bind(element, fvGeometry, prevElemVolVars);
211 this->elemBcTypes().update(
problem(), this->element(), this->fvGeometry());
215 template<std::
size_t i = domainId>
217 {
return this->assembler().problem(
domainId); }
220 template<std::
size_t i = domainId>
222 {
return ParentType::curSol()[dId]; }
226 {
return couplingManager_; }
243template<std::
size_t id,
class TypeTag,
class Assembler, DiffMethod DM = DiffMethod::numeric,
bool implicit = true>
252template<std::
size_t id,
class TypeTag,
class Assembler>
262 using GridView =
typename GridGeometry::GridView;
263 using FVElementGeometry =
typename GridGeometry::LocalView;
264 using Element =
typename GridView::template Codim<0>::Entity;
268 static constexpr int dim = GridView::dimension;
272 static constexpr auto domainI = Dune::index_constant<id>();
282 template<
class ScvOrLocalDof,
class ElemSol>
285 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, scvOrLocalDof.dofIndex(), elemSol[Detail::LocalDofs::index(scvOrLocalDof)], pvIdx);
291 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
294 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *
this, origResiduals, A, gridVariables);
301 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
310 const auto& element = this->element();
311 const auto& fvGeometry = this->fvGeometry();
312 auto&& curElemVolVars = this->curElemVolVars();
313 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
316 auto updateCoupledVariables = [&] ()
320 if constexpr (enableGridFluxVarsCache)
322 if constexpr (enableGridVolVarsCache)
323 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
325 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, gridVariables.gridFluxVarsCache());
329 if constexpr (enableGridVolVarsCache)
330 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), elemFluxVarsCache);
332 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, elemFluxVarsCache);
337 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
338 const auto& curSolJ = this->curSol(domainJ);
339 for (
const auto globalJ : stencil)
342 const auto origPriVarsJ = curSolJ[globalJ];
343 auto priVarsJ = origPriVarsJ;
346 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ);
348 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
350 auto evalCouplingResidual = [&](Scalar priVar)
352 priVarsJ[pvIdx] = priVar;
353 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
354 updateCoupledVariables();
355 return this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ);
361 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
362 static const int numDiffMethod = getParamFromGroup<int>(paramGroup,
"Assembly.NumericDifferenceMethod");
363 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
366 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
369 for (
const auto& scv :
scvs(fvGeometry))
371 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
377 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
381 if constexpr (!Detail::hasProblemBoundaryTypesForIntersectionFunction<Problem, FVElementGeometry, typename GridGeometry::GridView::Intersection>())
385 if (this->elemBcTypes().hasDirichlet())
387 const auto bcTypes = this->elemBcTypes().get(fvGeometry, scv);
388 if (bcTypes.isCouplingDirichlet(eqIdx))
389 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][eqIdx];
390 else if (bcTypes.isDirichlet(eqIdx))
391 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
396 if constexpr (Problem::enableInternalDirichletConstraints())
398 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
399 if (internalDirichletConstraints[eqIdx])
400 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
406 if constexpr (Detail::LocalDofs::hasNonCVLocalDofsInterface<FVElementGeometry>())
408 for (
const auto& localDof : nonCVLocalDofs(fvGeometry))
410 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
416 A[localDof.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[localDof.index()][eqIdx];
422 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
425 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
432 updateCoupledVariables();
443template<std::
size_t id,
class TypeTag,
class Assembler>
446 SubDomainCVFELocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>, DiffMethod::numeric, false >
462 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
An assembler for Jacobian and residual contribution per element (CVFE methods)
An assembler for Jacobian and residual contribution per element (CVFE methods)
Definition: assembly/cvfelocalassembler.hh:317
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:37
Definition: partialreassembler.hh:32
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
CVFE scheme multi domain local assembler using numeric differentiation and explicit time discretizati...
Definition: multidomain/subdomaincvfelocalassembler.hh:447
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: multidomain/subdomaincvfelocalassembler.hh:454
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const ElementResidualVector &res, GridVariables &gridVariables)
Computes the coupling derivatives with respect to the given element and adds them to the global matri...
Definition: multidomain/subdomaincvfelocalassembler.hh:463
CVFE scheme multi domain local assembler using numeric differentiation and implicit time discretizati...
Definition: multidomain/subdomaincvfelocalassembler.hh:256
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: multidomain/subdomaincvfelocalassembler.hh:292
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: multidomain/subdomaincvfelocalassembler.hh:277
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: multidomain/subdomaincvfelocalassembler.hh:302
void maybeUpdateCouplingContext(const ScvOrLocalDof &scvOrLocalDof, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: multidomain/subdomaincvfelocalassembler.hh:283
A base class for all CVFE subdomain local assemblers.
Definition: multidomain/subdomaincvfelocalassembler.hh:53
const auto & curSol(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: multidomain/subdomaincvfelocalassembler.hh:221
static constexpr auto domainId
export the domain id of this sub-domain
Definition: multidomain/subdomaincvfelocalassembler.hh:80
SubDomainCVFELocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: multidomain/subdomaincvfelocalassembler.hh:87
ElementResidualVector evalLocalSourceResidual(const Element &element, const ElementVolumeVariables &elemVolVars) const
Evaluates the local source term for an element and given element volume variables.
Definition: multidomain/subdomaincvfelocalassembler.hh:152
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: multidomain/subdomaincvfelocalassembler.hh:181
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: multidomain/subdomaincvfelocalassembler.hh:84
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition: multidomain/subdomaincvfelocalassembler.hh:175
const Problem & problem(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: multidomain/subdomaincvfelocalassembler.hh:216
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: multidomain/subdomaincvfelocalassembler.hh:225
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: multidomain/subdomaincvfelocalassembler.hh:134
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: multidomain/subdomaincvfelocalassembler.hh:110
The CVFE scheme multidomain local assembler.
Definition: multidomain/subdomaincvfelocalassembler.hh:244
Defines all properties used in Dumux.
An enum class to define various differentiation methods available in order to compute the derivatives...
Helper classes to compute the integration elements.
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition: diffmethod.hh:25
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
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Class representing dofs on elements for control-volume finite element schemes.
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition: localdof.hh:79
Definition: common/pdesolver.hh:24
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.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A arithmetic block vector type based on DUNE's reserved vector.