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>
50template<std::
size_t id,
class TypeTag,
class Assembler,
class Implementation, DiffMethod dm>
58 using SolutionVector =
typename Assembler::SolutionVector;
62 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
63 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
64 using ElementFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
68 using FVElementGeometry =
typename GridGeometry::LocalView;
69 using SubControlVolume =
typename GridGeometry::SubControlVolume;
70 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
72 using GridView =
typename GridGeometry::GridView;
73 using Element =
typename GridView::template Codim<0>::Entity;
75 using CouplingManager =
typename Assembler::CouplingManager;
81 static constexpr auto domainId =
typename Dune::index_constant<id>();
83 using ParentType::ParentType;
89 const Assembler& assembler,
90 const Element& element,
91 const SolutionVector&
curSol,
101 (element.partitionType() ==
Dune::GhostEntity),
102 assembler.isImplicit())
110 template<
class JacobianMatrixRow,
class SubRes
idualVector,
class Gr
idVariablesTuple,
class StageParams>
112 const StageParams& stageParams, SubResidualVector& temporal, SubResidualVector& spatial,
113 SubResidualVector& constrainedDofs)
115 auto assembleCouplingBlocks = [&](
const auto& residual)
118 using namespace Dune::Hybrid;
119 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](
auto&& i)
121 if constexpr (std::decay_t<
decltype(i)>{} != id)
128 ParentType::assembleJacobianAndResidual(
130 *std::get<domainId>(gridVariables),
131 stageParams, temporal, spatial,
134 assembleCouplingBlocks
142 template<std::
size_t otherId,
class JacRow,
class Gr
idVariables>
146 if constexpr (
id != otherId)
147 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
156 const auto& element = this->element();
158 auto&& fvGeometry = this->fvGeometry();
159 auto&& curElemVolVars = this->curElemVolVars();
160 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
163 couplingManager_.bindCouplingContext(
domainId, element, this->assembler());
164 fvGeometry.bind(element);
165 if (std::abs(this->localResidual().spatialWeight()) < 1e-6)
166 curElemVolVars.bindElement(element, fvGeometry,
curSol);
169 curElemVolVars.bind(element, fvGeometry,
curSol);
170 elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars);
173 this->elemBcTypes().update(
problem(), this->element(), this->fvGeometry());
177 template<std::
size_t i = domainId>
179 {
return this->assembler().problem(
domainId); }
182 template<std::
size_t i = domainId>
184 {
return ParentType::curSol()[dId]; }
188 {
return couplingManager_; }
205template<std::
size_t id,
class TypeTag,
class Assembler, DiffMethod DM = DiffMethod::numeric>
215template<std::
size_t id,
class TypeTag,
class Assembler>
218 SubDomainCVFELocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric>, DiffMethod::numeric>
227 using FVElementGeometry =
typename GridGeometry::LocalView;
228 using SubControlVolume =
typename GridGeometry::SubControlVolume;
229 using GridView =
typename GridGeometry::GridView;
230 using Element =
typename GridView::template Codim<0>::Entity;
234 static constexpr int dim = GridView::dimension;
238 static constexpr auto domainI = Dune::index_constant<id>();
248 template<
class ElemSol>
251 if (this->assembler().isImplicit())
252 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, scv.dofIndex(), elemSol[scv.localDofIndex()], pvIdx);
258 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
261 if (this->assembler().isImplicit())
262 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *
this, origResiduals, A, gridVariables);
269 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
273 if (!this->assembler().isImplicit())
281 const auto& element = this->element();
282 const auto& fvGeometry = this->fvGeometry();
283 auto&& curElemVolVars = this->curElemVolVars();
284 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
287 auto updateCoupledVariables = [&] ()
291 if constexpr (enableGridFluxVarsCache)
293 if constexpr (enableGridVolVarsCache)
294 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
296 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, gridVariables.gridFluxVarsCache());
300 if constexpr (enableGridVolVarsCache)
301 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), elemFluxVarsCache);
303 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, elemFluxVarsCache);
308 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
309 const auto& curSolJ = this->curSol(domainJ);
310 for (
const auto globalJ : stencil)
313 const auto origPriVarsJ = curSolJ[globalJ];
314 auto priVarsJ = origPriVarsJ;
317 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ);
319 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
321 auto evalCouplingResidual = [&](Scalar priVar)
323 priVarsJ[pvIdx] = priVar;
324 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
325 updateCoupledVariables();
326 return this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ);
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);
337 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
340 for (
const auto& scv : scvs(fvGeometry))
342 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
348 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
352 if (this->elemBcTypes().hasDirichlet())
354 const auto bcTypes = this->elemBcTypes().get(fvGeometry, scv);
355 if (bcTypes.isCouplingDirichlet(eqIdx))
356 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][eqIdx];
357 else if (bcTypes.isDirichlet(eqIdx))
358 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
362 if constexpr (Problem::enableInternalDirichletConstraints())
364 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
365 if (internalDirichletConstraints[eqIdx])
366 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
373 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
376 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
383 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 (CVFE methods)
Definition: experimental/assembly/cvfelocalassembler.hh:326
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:219
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:270
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:259
void maybeUpdateCouplingContext(const SubControlVolume &scv, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:249
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:243
A base class for all CVFE subdomain local assemblers.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:52
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:153
SubDomainCVFELocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:88
static constexpr auto domainId
export the domain id of this sub-domain
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:81
const auto & curSol(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:183
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:143
const Problem & problem(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:178
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:111
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:85
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:187
The CVFE scheme multidomain local assembler.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:206
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...
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
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.