26#ifndef DUMUX_MULTIDOMAIN_FACECENTERED_LOCAL_ASSEMBLER_HH
27#define DUMUX_MULTIDOMAIN_FACECENTERED_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 Element& element,
89 const SolutionVector&
curSol,
99 (element.partitionType() ==
Dune::GhostEntity))
107 template<
class JacobianMatrixRow,
class Gr
idVariablesTuple>
110 auto assembleCouplingBlocks = [&](
const auto& residual)
113 using namespace Dune::Hybrid;
114 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](
auto&& i)
123 ParentType::assembleJacobianAndResidual(jacRow[
domainId], res, *std::get<domainId>(gridVariables), noReassembler, assembleCouplingBlocks);
130 template<std::size_t otherId,
class JacRow,
class GridVariables,
131 typename std::enable_if_t<(otherId == id),
int> = 0>
139 template<std::size_t otherId,
class JacRow,
class GridVariables,
140 typename std::enable_if_t<(otherId != id),
int> = 0>
144 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
157 for (
auto&& scv : scvs(this->fvGeometry()))
159 const auto& curVolVars = elemVolVars[scv];
160 auto source = this->localResidual().computeSource(
problem(), element, this->fvGeometry(), elemVolVars, scv);
161 source *= -scv.volume()*curVolVars.extrusionFactor();
162 residual[scv.localDofIndex()] = std::move(source);
172 {
return this->
evalLocalSourceResidual(neighbor, implicit ? this->curElemVolVars() : this->prevElemVolVars()); }
180 const auto& element = this->element();
182 auto&& fvGeometry = this->fvGeometry();
183 auto&& curElemVolVars = this->curElemVolVars();
184 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
187 couplingManager_.bindCouplingContext(
domainId, element, this->assembler());
188 fvGeometry.bind(element);
192 curElemVolVars.bind(element, fvGeometry,
curSol);
193 elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars);
194 if (!this->assembler().isStationaryProblem())
195 this->prevElemVolVars().bindElement(element, fvGeometry, this->assembler().prevSol()[
domainId]);
199 auto& prevElemVolVars = this->prevElemVolVars();
200 const auto& prevSol = this->assembler().prevSol()[
domainId];
202 curElemVolVars.bindElement(element, fvGeometry,
curSol);
203 prevElemVolVars.bind(element, fvGeometry, prevSol);
204 elemFluxVarsCache.bind(element, fvGeometry, prevElemVolVars);
207 this->elemBcTypes().update(
problem(), this->element(), this->fvGeometry());
211 template<std::
size_t i = domainId>
213 {
return this->assembler().problem(
domainId); }
216 template<std::
size_t i = domainId>
218 {
return ParentType::curSol()[dId]; }
222 {
return couplingManager_; }
239template<std::
size_t id,
class TypeTag,
class Assembler, DiffMethod DM = DiffMethod::numeric,
bool implicit = true>
248template<std::
size_t id,
class TypeTag,
class Assembler>
251 id, TypeTag, 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);
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 const auto& curSolJ = this->curSol(domainJ);
319 auto updateCoupledVariables = [&] ()
323 if constexpr (enableGridFluxVarsCache)
325 if constexpr (enableGridVolVarsCache)
326 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
328 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, gridVariables.gridFluxVarsCache());
332 if constexpr (enableGridVolVarsCache)
333 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), elemFluxVarsCache);
335 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, elemFluxVarsCache);
339 for (
const auto& scv : scvs(fvGeometry))
341 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, scv, domainJ);
343 for (
const auto globalJ : stencil)
345 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *
this, scv, domainJ, globalJ);
347 const auto origPriVarsJ = curSolJ[globalJ];
348 auto priVarsJ = origPriVarsJ;
350 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
352 auto evalCouplingResidual = [&](Scalar priVar)
354 priVarsJ[pvIdx] = priVar;
355 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
356 updateCoupledVariables();
357 return this->couplingManager().evalCouplingResidual(domainI, *
this, scv, domainJ, globalJ);
363 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
364 static const int numDiffMethod = getParamFromGroup<int>(paramGroup,
"Assembly.NumericDifferenceMethod");
365 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
368 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
371 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
377 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
382 if (scv.boundary() && this->elemBcTypes().hasDirichlet())
384 const auto bcTypes = this->elemBcTypes()[fvGeometry.frontalScvfOnBoundary(scv).localIndex()];
385 if (bcTypes.hasDirichlet())
390 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
392 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
394 if (bcTypes.isCouplingDirichlet(eqIdx))
395 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][pvIdx];
396 else if (bcTypes.isDirichlet(eqIdx))
397 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
404 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
407 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
414 updateCoupledVariables();
426template<std::
size_t id,
class TypeTag,
class Assembler>
429 id, TypeTag, Assembler,
430 SubDomainFaceCenteredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>,
431 DiffMethod::numeric, false
437 using FVElementGeometry =
typename GridGeometry::LocalView;
438 using GridView =
typename GridGeometry::GridView;
439 using Element =
typename GridView::template Codim<0>::Entity;
440 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
441 static constexpr auto domainI = Dune::index_constant<id>();
450 template<
class ElemSol>
453 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, scv.dofIndex(), elemSol[scv.localDofIndex()], pvIdx);
459 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
467 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
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 (face-centered staggered methods)
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
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition: diffmethod.hh:37
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
An assembler for Jacobian and residual contribution per element (Face-centered methods)
Definition: fclocalassembler.hh:296
Definition: partialreassembler.hh:44
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 all face-centered staggered local assemblers.
Definition: subdomainfclocalassembler.hh:55
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition: subdomainfclocalassembler.hh:171
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: subdomainfclocalassembler.hh:221
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: subdomainfclocalassembler.hh:108
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: subdomainfclocalassembler.hh:177
const auto & curSol(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: subdomainfclocalassembler.hh:217
const Problem & problem(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: subdomainfclocalassembler.hh:212
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:150
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:132
SubDomainFaceCenteredLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: subdomainfclocalassembler.hh:87
static constexpr auto domainId
export the domain id of this sub-domain
Definition: subdomainfclocalassembler.hh:80
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: subdomainfclocalassembler.hh:84
The face-centered staggered scheme multidomain local assembler.
Definition: subdomainfclocalassembler.hh:240
Face-centered staggered scheme multi domain local assembler using numeric differentiation and implici...
Definition: subdomainfclocalassembler.hh:255
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:302
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: subdomainfclocalassembler.hh:277
void maybeUpdateCouplingContext(const SubControlVolume &scv, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: subdomainfclocalassembler.hh:283
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, const JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: subdomainfclocalassembler.hh:292
Face-centered staggered scheme multi domain local assembler using numeric differentiation and explici...
Definition: subdomainfclocalassembler.hh:433
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: subdomainfclocalassembler.hh:445
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:468
void maybeUpdateCouplingContext(const SubControlVolume &scv, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: subdomainfclocalassembler.hh:451
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, const JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: subdomainfclocalassembler.hh:460
Declares all properties used in Dumux.