26#ifndef DUMUX_MULTIDOMAIN_PQ1BUBBLE_SUBDOMAIN_LOCAL_ASSEMBLER_HH
27#define DUMUX_MULTIDOMAIN_PQ1BUBBLE_SUBDOMAIN_LOCAL_ASSEMBLER_HH
29#include <dune/common/indices.hh>
30#include <dune/common/hybridutilities.hh>
31#include <dune/grid/common/gridenums.hh>
53template<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;
70 using SubControlVolume =
typename GridGeometry::SubControlVolume;
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 Gr
idVariablesTuple>
112 auto assembleCouplingBlocks = [&](
const auto& residual)
115 using namespace Dune::Hybrid;
116 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](
auto&& i)
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));
159 for (
auto&& scv : scvs(this->fvGeometry()))
161 const auto& curVolVars = elemVolVars[scv];
162 auto source = this->localResidual().computeSource(
problem(), element, this->fvGeometry(), elemVolVars, scv);
163 source *= -scv.volume()*curVolVars.extrusionFactor();
164 residual[scv.localDofIndex()] = std::move(source);
174 {
return this->
evalLocalSourceResidual(neighbor, implicit ? this->curElemVolVars() : this->prevElemVolVars()); }
182 const auto& element = this->element();
184 auto&& fvGeometry = this->fvGeometry();
185 auto&& curElemVolVars = this->curElemVolVars();
186 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
189 couplingManager_.bindCouplingContext(
domainId, element, this->assembler());
190 fvGeometry.bind(element);
194 curElemVolVars.bind(element, fvGeometry,
curSol);
195 elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars);
196 if (!this->assembler().isStationaryProblem())
197 this->prevElemVolVars().bindElement(element, fvGeometry, this->assembler().prevSol()[
domainId]);
201 auto& prevElemVolVars = this->prevElemVolVars();
202 const auto& prevSol = this->assembler().prevSol()[
domainId];
204 curElemVolVars.bindElement(element, fvGeometry,
curSol);
205 prevElemVolVars.bind(element, fvGeometry, prevSol);
206 elemFluxVarsCache.bind(element, fvGeometry, prevElemVolVars);
209 this->elemBcTypes().update(
problem(), this->element(), this->fvGeometry());
213 template<std::
size_t i = domainId>
215 {
return this->assembler().problem(
domainId); }
218 template<std::
size_t i = domainId>
220 {
return ParentType::curSol()[dId]; }
224 {
return couplingManager_; }
241template<std::
size_t id,
class TypeTag,
class Assembler, DiffMethod DM = DiffMethod::numeric,
bool implicit = true>
250template<std::
size_t id,
class TypeTag,
class Assembler>
262 using FVElementGeometry =
typename GridGeometry::LocalView;
263 using GridView =
typename GridGeometry::GridView;
264 using Element =
typename GridView::template Codim<0>::Entity;
265 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
268 enum { dim = GridView::dimension };
272 static constexpr auto domainI = Dune::index_constant<id>();
282 template<
class ElemSol>
285 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, scv.dofIndex(), elemSol[scv.localDofIndex()], pvIdx);
291 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
294 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *
this, origResiduals, A, gridVariables);
303 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
312 const auto& element = this->element();
313 const auto& fvGeometry = this->fvGeometry();
314 auto&& curElemVolVars = this->curElemVolVars();
315 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
318 auto updateCoupledVariables = [&] ()
322 if constexpr (enableGridFluxVarsCache)
324 if constexpr (enableGridVolVarsCache)
325 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
327 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, gridVariables.gridFluxVarsCache());
331 if constexpr (enableGridVolVarsCache)
332 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), elemFluxVarsCache);
334 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, elemFluxVarsCache);
338 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
339 const auto& curSolJ = this->curSol(domainJ);
340 for (
const auto globalJ : stencil)
343 const auto origPriVarsJ = curSolJ[globalJ];
344 auto priVarsJ = origPriVarsJ;
347 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ);
349 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
351 auto evalCouplingResidual = [&](Scalar priVar)
353 priVarsJ[pvIdx] = priVar;
354 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
355 updateCoupledVariables();
356 return this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ);
362 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
363 static const int numDiffMethod = getParamFromGroup<int>(paramGroup,
"Assembly.NumericDifferenceMethod");
364 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
367 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
370 for (
const auto& scv : scvs(fvGeometry))
372 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
378 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
382 if (this->elemBcTypes().hasDirichlet())
384 const auto bcTypes = this->elemBcTypes().get(fvGeometry, scv);
385 if (bcTypes.isCouplingDirichlet(eqIdx))
386 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][eqIdx];
387 else if (bcTypes.isDirichlet(eqIdx))
388 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
392 if constexpr (Problem::enableInternalDirichletConstraints())
394 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
395 if (internalDirichletConstraints[eqIdx])
396 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
402 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
405 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
412 updateCoupledVariables();
An enum class to define various differentiation methods available in order to compute the derivatives...
An adapter class for local assemblers using numeric differentiation.
An assembler for Jacobian and residual contribution per element (pq1bubble method)
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A class for numeric differentiation.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
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
Definition: deprecated.hh:149
Definition: partialreassembler.hh:44
An assembler for Jacobian and residual contribution per element (PQ1Bubble methods)
Definition: pq1bubblelocalassembler.hh:304
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
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:60
A base class for PQ1Bubble local assemblers.
Definition: subdomainpq1bubblelocalassembler.hh:55
const Problem & problem(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: subdomainpq1bubblelocalassembler.hh:214
static constexpr auto domainId
export the domain id of this sub-domain
Definition: subdomainpq1bubblelocalassembler.hh:80
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: subdomainpq1bubblelocalassembler.hh:84
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: subdomainpq1bubblelocalassembler.hh:223
ElementResidualVector evalLocalSourceResidual(const Element &element, const ElementVolumeVariables &elemVolVars) const
Evaluates the local source term for an element and given element volume variables.
Definition: subdomainpq1bubblelocalassembler.hh:152
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: subdomainpq1bubblelocalassembler.hh:134
void assembleJacobianAndResidual(JacobianMatrixRow &jacRow, SubSolutionVector &res, GridVariablesTuple &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: subdomainpq1bubblelocalassembler.hh:110
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition: subdomainpq1bubblelocalassembler.hh:173
SubDomainPQ1BubbleLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: subdomainpq1bubblelocalassembler.hh:87
const auto & curSol(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: subdomainpq1bubblelocalassembler.hh:219
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: subdomainpq1bubblelocalassembler.hh:179
The PQ1Bubble scheme multidomain local assembler.
Definition: subdomainpq1bubblelocalassembler.hh:242
Control-volume fe staggered scheme multi domain local assembler using numeric differentiation and imp...
Definition: subdomainpq1bubblelocalassembler.hh:254
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, const JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: subdomainpq1bubblelocalassembler.hh:292
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: subdomainpq1bubblelocalassembler.hh:304
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: subdomainpq1bubblelocalassembler.hh:277
void maybeUpdateCouplingContext(const SubControlVolume &scv, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: subdomainpq1bubblelocalassembler.hh:283
Declares all properties used in Dumux.