15#ifndef DUMUX_EXPERIMENTAL_SUBDOMAIN_CVFE_LOCAL_ASSEMBLER_HH
16#define DUMUX_EXPERIMENTAL_SUBDOMAIN_CVFE_LOCAL_ASSEMBLER_HH
20#include <dune/common/reservedvector.hh>
21#include <dune/common/indices.hh>
22#include <dune/common/hybridutilities.hh>
23#include <dune/grid/common/gridenums.hh>
24#include <dune/istl/matrixindexset.hh>
26#include <dumux/common/typetraits/localdofs_.hh>
51template<std::
size_t id,
class TypeTag,
class Assembler,
class Implementation, DiffMethod dm>
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;
69 using FVElementGeometry =
typename GridGeometry::LocalView;
70 using SubControlVolume =
typename GridGeometry::SubControlVolume;
71 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
73 using GridView =
typename GridGeometry::GridView;
74 using Element =
typename GridView::template Codim<0>::Entity;
76 using CouplingManager =
typename Assembler::CouplingManager;
82 static constexpr auto domainId =
typename Dune::index_constant<id>();
84 using ParentType::ParentType;
90 const Assembler& assembler,
91 const Element& element,
92 const SolutionVector&
curSol,
102 (element.partitionType() ==
Dune::GhostEntity),
103 assembler.isImplicit())
111 template<
class JacobianMatrixRow,
class SubRes
idualVector,
class Gr
idVariablesTuple,
class StageParams>
113 const StageParams& stageParams, SubResidualVector& temporal, SubResidualVector& spatial,
114 SubResidualVector& constrainedDofs)
116 auto assembleCouplingBlocks = [&](
const auto& residual)
119 using namespace Dune::Hybrid;
120 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](
auto&& i)
122 if constexpr (std::decay_t<
decltype(i)>{} != id)
129 ParentType::assembleJacobianAndResidual(
131 *std::get<domainId>(gridVariables),
132 stageParams, temporal, spatial,
135 assembleCouplingBlocks
143 template<std::
size_t otherId,
class JacRow,
class Gr
idVariables>
147 if constexpr (
id != otherId)
148 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
157 const auto& element = this->element();
159 auto&& fvGeometry = this->fvGeometry();
160 auto&& curElemVolVars = this->curElemVolVars();
161 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
164 couplingManager_.bindCouplingContext(
domainId, element, this->assembler());
165 fvGeometry.bind(element);
166 if (std::abs(this->localResidual().spatialWeight()) < 1e-6)
167 curElemVolVars.bindElement(element, fvGeometry,
curSol);
170 curElemVolVars.bind(element, fvGeometry,
curSol);
171 elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars);
174 this->elemBcTypes().update(
problem(), this->element(), this->fvGeometry());
178 template<std::
size_t i = domainId>
180 {
return this->assembler().problem(
domainId); }
183 template<std::
size_t i = domainId>
185 {
return ParentType::curSol()[dId]; }
189 {
return couplingManager_; }
206template<std::
size_t id,
class TypeTag,
class Assembler, DiffMethod DM = DiffMethod::numeric>
216template<std::
size_t id,
class TypeTag,
class Assembler>
219 SubDomainCVFELocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric>, DiffMethod::numeric>
228 using FVElementGeometry =
typename GridGeometry::LocalView;
229 using SubControlVolume =
typename GridGeometry::SubControlVolume;
230 using GridView =
typename GridGeometry::GridView;
231 using Element =
typename GridView::template Codim<0>::Entity;
235 static constexpr int dim = GridView::dimension;
239 static constexpr auto domainI = Dune::index_constant<id>();
249 template<
class ScvOrLocalDof,
class ElemSol>
252 if (this->assembler().isImplicit())
253 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, scvOrLocalDof.dofIndex(), elemSol[Dumux::Detail::LocalDofs::index(scvOrLocalDof)], pvIdx);
259 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
262 if (this->assembler().isImplicit())
263 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *
this, origResiduals, A, gridVariables);
270 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
274 if (!this->assembler().isImplicit())
282 const auto& element = this->element();
283 const auto& fvGeometry = this->fvGeometry();
284 auto&& curElemVolVars = this->curElemVolVars();
285 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
288 auto updateCoupledVariables = [&] ()
292 if constexpr (enableGridFluxVarsCache)
294 if constexpr (enableGridVolVarsCache)
295 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
297 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, gridVariables.gridFluxVarsCache());
301 if constexpr (enableGridVolVarsCache)
302 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), elemFluxVarsCache);
304 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, elemFluxVarsCache);
309 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
310 const auto& curSolJ = this->curSol(domainJ);
311 for (
const auto globalJ : stencil)
314 const auto origPriVarsJ = curSolJ[globalJ];
315 auto priVarsJ = origPriVarsJ;
318 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ);
320 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
322 auto evalCouplingResidual = [&](Scalar priVar)
324 priVarsJ[pvIdx] = priVar;
325 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
326 updateCoupledVariables();
327 return this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ);
333 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
334 static const int numDiffMethod = getParamFromGroup<int>(paramGroup,
"Assembly.NumericDifferenceMethod");
335 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
338 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
341 for (
const auto& scv :
scvs(fvGeometry))
343 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
349 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
353 if (this->elemBcTypes().hasDirichlet())
355 const auto bcTypes = this->elemBcTypes().get(fvGeometry, scv);
356 if (bcTypes.isCouplingDirichlet(eqIdx))
357 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][eqIdx];
358 else if (bcTypes.isDirichlet(eqIdx))
359 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
363 if constexpr (Problem::enableInternalDirichletConstraints())
365 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
366 if (internalDirichletConstraints[eqIdx])
367 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
374 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
377 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
384 updateCoupledVariables();
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 (CVFE methods)
Definition: experimental/assembly/cvfelocalassembler.hh:328
GG GridGeometry
export the grid geometry type
Definition: experimental/discretization/gridvariables.hh:38
CVFE scheme multi domain local assembler using numeric differentiation.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:220
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: experimental/assembly/subdomaincvfelocalassembler.hh:271
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:260
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:244
void maybeUpdateCouplingContext(const ScvOrLocalDof &scvOrLocalDof, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:250
A base class for all CVFE subdomain local assemblers.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:53
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:154
SubDomainCVFELocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:89
static constexpr auto domainId
export the domain id of this sub-domain
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:82
const auto & curSol(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:184
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: experimental/assembly/subdomaincvfelocalassembler.hh:144
const Problem & problem(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:179
void assembleJacobianAndResidual(JacobianMatrixRow &jacRow, SubResidualVector &res, GridVariablesTuple &gridVariables, const StageParams &stageParams, SubResidualVector &temporal, SubResidualVector &spatial, SubResidualVector &constrainedDofs)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:112
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:86
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:188
The CVFE scheme multidomain local assembler.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:207
typename ScalarT< X >::type Scalar
export the underlying scalar type
Definition: experimental/common/variables.hh:57
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.
An enum class to define various differentiation methods available in order to compute the derivatives...
An assembler for Jacobian and residual contribution per element (CVFE methods)
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
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Definition: experimental/assembly/cclocalassembler.hh:36
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.