15#ifndef DUMUX_EXPERIMENTAL_SUBDOMAIN_CC_LOCAL_ASSEMBLER_HH
16#define DUMUX_EXPERIMENTAL_SUBDOMAIN_CC_LOCAL_ASSEMBLER_HH
18#include <dune/common/indices.hh>
19#include <dune/common/reservedvector.hh>
20#include <dune/common/hybridutilities.hh>
21#include <dune/grid/common/gridenums.hh>
22#include <dune/istl/matrixindexset.hh>
49template<std::
size_t id,
class TypeTag,
class Assembler,
class Implementation, DiffMethod dm>
57 using SolutionVector =
typename Assembler::SolutionVector;
61 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
62 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
63 using ElementFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
67 using FVElementGeometry =
typename GridGeometry::LocalView;
68 using SubControlVolume =
typename GridGeometry::SubControlVolume;
69 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
71 using GridView =
typename GridGeometry::GridView;
72 using Element =
typename GridView::template Codim<0>::Entity;
74 using CouplingManager =
typename Assembler::CouplingManager;
75 using ElementResidualVector =
typename ParentType::LocalResidual::ElementResidualVector;
79 static constexpr auto domainId =
typename Dune::index_constant<id>();
83 using ParentType::ParentType;
87 const Element& element,
88 const SolutionVector&
curSol,
97 (element.partitionType() ==
Dune::GhostEntity),
98 assembler.isImplicit())
106 template<
class JacobianMatrixRow,
class SubRes
idualVector,
class Gr
idVariablesTuple,
class StageParams>
108 const StageParams& stageParams, SubResidualVector& temporal, SubResidualVector& spatial,
109 SubResidualVector& constrainedDofs)
111 auto assembleCouplingBlocks = [&](
const auto& residual)
114 using namespace Dune::Hybrid;
115 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](
auto&& i)
117 if (std::decay_t<
decltype(i)>{} != id)
124 ParentType::assembleJacobianAndResidual(
126 *std::get<domainId>(gridVariables),
127 stageParams, temporal, spatial,
130 assembleCouplingBlocks
138 template<std::
size_t otherId,
class JacRow,
class Gr
idVariables>
140 const LocalResidualValues& res, GridVariables& gridVariables)
142 if constexpr (
id != otherId)
143 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
152 const auto& element = this->element();
154 auto&& fvGeometry = this->fvGeometry();
155 auto&& curElemVolVars = this->curElemVolVars();
156 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
159 couplingManager_.bindCouplingContext(
domainId, element, this->assembler());
160 fvGeometry.bind(element);
161 curElemVolVars.bind(element, fvGeometry,
curSol);
162 elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars);
166 template<std::
size_t i = domainId>
168 {
return this->assembler().problem(
domainId); }
171 template<std::
size_t i = domainId>
173 {
return ParentType::curSol()[dId]; }
177 {
return couplingManager_; }
194template<std::
size_t id,
class TypeTag,
class Assembler, DiffMethod DM = DiffMethod::numeric>
204template<std::
size_t id,
class TypeTag,
class Assembler>
207 SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric>, DiffMethod::numeric>
217 using FVElementGeometry =
typename GridGeometry::LocalView;
218 using SubControlVolume =
typename GridGeometry::SubControlVolume;
219 using GridView =
typename GridGeometry::GridView;
220 using Element =
typename GridView::template Codim<0>::Entity;
223 enum { dim = GridView::dimension };
227 static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize;
228 static constexpr auto domainI = Dune::index_constant<id>();
237 template<
class ElemSol>
240 if (this->assembler().isImplicit())
241 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, scv.dofIndex(), elemSol[0], pvIdx);
247 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
250 if (this->assembler().isImplicit())
251 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *
this, origResiduals, A, gridVariables);
258 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
260 const LocalResidualValues& res, GridVariables& gridVariables)
262 if (!this->assembler().isImplicit())
270 const auto& element = this->element();
271 const auto& fvGeometry = this->fvGeometry();
272 const auto& gridGeometry = fvGeometry.gridGeometry();
273 auto&& curElemVolVars = this->curElemVolVars();
274 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
277 const auto globalI = gridGeometry.elementMapper().index(element);
278 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
279 const auto& curSolJ = this->curSol(domainJ);
282 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
285 auto updateCoupledVariables = [&] ()
289 if (enableGridFluxVarsCache)
291 if (enableGridVolVarsCache)
293 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
294 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), elemFluxVarsCache);
298 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, gridVariables.gridFluxVarsCache());
299 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, elemFluxVarsCache);
304 if (enableGridVolVarsCache)
305 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), elemFluxVarsCache);
307 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, elemFluxVarsCache);
311 for (
const auto& globalJ : stencil)
314 const auto origPriVarsJ = curSolJ[globalJ];
315 auto priVarsJ = origPriVarsJ;
317 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ)[0];
319 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
321 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)[0];
331 LocalResidualValues partialDeriv(0.0);
332 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
333 static const int numDiffMethod = getParamFromGroup<int>(paramGroup,
"Assembly.NumericDifferenceMethod");
334 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
336 epsCoupl(priVarsJ[pvIdx], pvIdx), numDiffMethod);
339 if constexpr (Problem::enableInternalDirichletConstraints())
341 const auto& scv = fvGeometry.scv(globalI);
342 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
343 if (internalDirichletConstraints.none())
345 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
346 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
350 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
352 if (internalDirichletConstraints[eqIdx])
353 A[globalI][globalJ][eqIdx][pvIdx] = 0.0;
355 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
361 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
362 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
366 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
369 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
375 updateCoupledVariables();
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:48
Definition: partialreassembler.hh:32
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition: experimental/assembly/cclocalassembler.hh:181
GG GridGeometry
export the grid geometry type
Definition: experimental/discretization/gridvariables.hh:38
Cell-centered scheme multidomain local assembler using numeric differentiation.
Definition: experimental/assembly/subdomaincclocalassembler.hh:208
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const LocalResidualValues &res, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: experimental/assembly/subdomaincclocalassembler.hh:259
void maybeUpdateCouplingContext(const SubControlVolume &scv, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: experimental/assembly/subdomaincclocalassembler.hh:238
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: experimental/assembly/subdomaincclocalassembler.hh:248
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
Definition: experimental/assembly/subdomaincclocalassembler.hh:232
A base class for all multidomain local assemblers.
Definition: experimental/assembly/subdomaincclocalassembler.hh:51
const auto & curSol(Dune::index_constant< i > dId=domainId) const
return reference to the subdomain solution
Definition: experimental/assembly/subdomaincclocalassembler.hh:172
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: experimental/assembly/subdomaincclocalassembler.hh:176
const Problem & problem(Dune::index_constant< i > dId=domainId) const
return reference to the subdomain problem
Definition: experimental/assembly/subdomaincclocalassembler.hh:167
SubDomainCCLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
the constructor
Definition: experimental/assembly/subdomaincclocalassembler.hh:86
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/subdomaincclocalassembler.hh:107
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacRow &jacRow, const LocalResidualValues &res, GridVariables &gridVariables)
Assemble the entries in a coupling block of the jacobian. There is no coupling block between a domain...
Definition: experimental/assembly/subdomaincclocalassembler.hh:139
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: experimental/assembly/subdomaincclocalassembler.hh:149
GetPropType< TypeTag, Properties::LocalResidual > LocalResidual
the local residual type of this domain
Definition: experimental/assembly/subdomaincclocalassembler.hh:81
static constexpr auto domainId
export the domain id of this sub-domain
Definition: experimental/assembly/subdomaincclocalassembler.hh:79
The cell-centered scheme multidomain local assembler.
Definition: experimental/assembly/subdomaincclocalassembler.hh:195
typename ScalarT< X >::type Scalar
export the underlying scalar type
Definition: 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:49
Defines all properties used in Dumux.
An enum class to define various differentiation methods available in order to compute the derivatives...
Element solution classes and factory functions.
An assembler for Jacobian and residual contribution per element (cell-centered 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
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.