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>
48template<std::
size_t id,
class TypeTag,
class Assembler,
class Implementation, DiffMethod dm,
bool implicit>
56 using SolutionVector =
typename Assembler::SolutionVector;
60 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
61 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
62 using ElementFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
63 using Scalar =
typename GridVariables::Scalar;
65 using GridGeometry =
typename GridVariables::GridGeometry;
66 using FVElementGeometry =
typename GridGeometry::LocalView;
67 using SubControlVolume =
typename GridGeometry::SubControlVolume;
68 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
70 using GridView =
typename GridGeometry::GridView;
71 using Element =
typename GridView::template Codim<0>::Entity;
73 using CouplingManager =
typename Assembler::CouplingManager;
79 static constexpr auto domainId =
typename Dune::index_constant<id>();
81 using ParentType::ParentType;
87 const Assembler& assembler,
88 const Element& element,
89 const SolutionVector&
curSol,
100 (element.partitionType() ==
Dune::GhostEntity))
108 template<
class JacobianMatrixRow,
class SubRes
idualVector,
class Gr
idVariablesTuple>
111 auto assembleCouplingBlocks = [&](
const auto& residual)
114 using namespace Dune::Hybrid;
115 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](
auto&& i)
117 if constexpr (std::decay_t<
decltype(i)>{} != id)
124 ParentType::assembleJacobianAndResidual(jacRow[
domainId], res, *std::get<domainId>(gridVariables), noReassembler, assembleCouplingBlocks);
131 template<std::size_t otherId,
class JacRow,
class GridVariables,
132 typename std::enable_if_t<(otherId == id),
int> = 0>
140 template<std::size_t otherId,
class JacRow,
class GridVariables,
141 typename std::enable_if_t<(otherId != id),
int> = 0>
145 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
158 for (
const auto& scv : scvs(this->fvGeometry()))
160 const auto& curVolVars = elemVolVars[scv];
161 auto source = this->localResidual().computeSource(
problem(), element, this->fvGeometry(), elemVolVars, scv);
162 source *= -
Extrusion::volume(this->fvGeometry(), scv)*curVolVars.extrusionFactor();
163 residual[scv.localDofIndex()] = std::move(source);
173 {
return this->
evalLocalSourceResidual(neighbor, implicit ? this->curElemVolVars() : this->prevElemVolVars()); }
181 const auto& element = this->element();
183 auto&& fvGeometry = this->fvGeometry();
184 auto&& curElemVolVars = this->curElemVolVars();
185 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
188 couplingManager_.bindCouplingContext(
domainId, element, this->assembler());
189 fvGeometry.bind(element);
191 if constexpr (implicit)
193 curElemVolVars.bind(element, fvGeometry,
curSol);
194 elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars);
195 if (!this->assembler().isStationaryProblem())
196 this->prevElemVolVars().bindElement(element, fvGeometry, this->assembler().prevSol()[
domainId]);
200 auto& prevElemVolVars = this->prevElemVolVars();
201 const auto& prevSol = this->assembler().prevSol()[
domainId];
203 curElemVolVars.bindElement(element, fvGeometry,
curSol);
204 prevElemVolVars.bind(element, fvGeometry, prevSol);
205 elemFluxVarsCache.bind(element, fvGeometry, prevElemVolVars);
208 this->elemBcTypes().update(
problem(), this->element(), this->fvGeometry());
212 template<std::
size_t i = domainId>
214 {
return this->assembler().problem(
domainId); }
217 template<std::
size_t i = domainId>
219 {
return ParentType::curSol()[dId]; }
223 {
return couplingManager_; }
240template<std::
size_t id,
class TypeTag,
class Assembler, DiffMethod DM = DiffMethod::numeric,
bool implicit = true>
249template<std::
size_t id,
class TypeTag,
class Assembler>
261 using FVElementGeometry =
typename GridGeometry::LocalView;
262 using SubControlVolume =
typename GridGeometry::SubControlVolume;
263 using GridView =
typename GridGeometry::GridView;
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 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 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 (this->elemBcTypes().hasDirichlet())
383 const auto bcTypes = this->elemBcTypes().get(fvGeometry, scv);
384 if (bcTypes.isCouplingDirichlet(eqIdx))
385 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][eqIdx];
386 else if (bcTypes.isDirichlet(eqIdx))
387 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
391 if constexpr (Problem::enableInternalDirichletConstraints())
393 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
394 if (internalDirichletConstraints[eqIdx])
395 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();
423template<std::
size_t id,
class TypeTag,
class Assembler>
426 SubDomainCVFELocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>, DiffMethod::numeric, false >
442 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:299
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:48
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:49
CVFE scheme multi domain local assembler using numeric differentiation and explicit time discretizati...
Definition: multidomain/subdomaincvfelocalassembler.hh:427
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: multidomain/subdomaincvfelocalassembler.hh:434
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:443
CVFE scheme multi domain local assembler using numeric differentiation and implicit time discretizati...
Definition: multidomain/subdomaincvfelocalassembler.hh:253
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: multidomain/subdomaincvfelocalassembler.hh:292
void maybeUpdateCouplingContext(const SubControlVolume &scv, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: multidomain/subdomaincvfelocalassembler.hh:283
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
A base class for all CVFE subdomain local assemblers.
Definition: multidomain/subdomaincvfelocalassembler.hh:50
const auto & curSol(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: multidomain/subdomaincvfelocalassembler.hh:218
static constexpr auto domainId
export the domain id of this sub-domain
Definition: multidomain/subdomaincvfelocalassembler.hh:79
SubDomainCVFELocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: multidomain/subdomaincvfelocalassembler.hh:86
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:151
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: multidomain/subdomaincvfelocalassembler.hh:178
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: multidomain/subdomaincvfelocalassembler.hh:83
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition: multidomain/subdomaincvfelocalassembler.hh:172
const Problem & problem(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: multidomain/subdomaincvfelocalassembler.hh:213
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: multidomain/subdomaincvfelocalassembler.hh:222
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:133
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:109
The CVFE scheme multidomain local assembler.
Definition: multidomain/subdomaincvfelocalassembler.hh:241
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
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.