version 3.10-dev
multidomain/subdomaincvfelocalassembler.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_MULTIDOMAIN_SUBDOMAIN_CVFE_LOCAL_ASSEMBLER_HH
15#define DUMUX_MULTIDOMAIN_SUBDOMAIN_CVFE_LOCAL_ASSEMBLER_HH
16
17#include <type_traits>
18
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> // for GhostEntity
23#include <dune/istl/matrixindexset.hh>
24
34
35namespace Dumux {
36
48template<std::size_t id, class TypeTag, class Assembler, class Implementation, DiffMethod dm, bool implicit>
49class SubDomainCVFELocalAssemblerBase : public CVFELocalAssembler<TypeTag, Assembler, dm, implicit, Implementation>
50{
52
56 using SolutionVector = typename Assembler::SolutionVector;
58
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;
64
65 using GridGeometry = typename GridVariables::GridGeometry;
66 using FVElementGeometry = typename GridGeometry::LocalView;
67 using SubControlVolume = typename GridGeometry::SubControlVolume;
68 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
69 using Extrusion = Extrusion_t<GridGeometry>;
70 using GridView = typename GridGeometry::GridView;
71 using Element = typename GridView::template Codim<0>::Entity;
72
73 using CouplingManager = typename Assembler::CouplingManager;
74
75 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
76
77public:
79 static constexpr auto domainId = typename Dune::index_constant<id>();
81 using ParentType::ParentType;
83 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
84
85 // the constructor
87 const Assembler& assembler,
88 const Element& element,
89 const SolutionVector& curSol,
90 CouplingManager& couplingManager
91 )
92 : ParentType(assembler,
93 element,
94 curSol,
95 localView(assembler.gridGeometry(domainId)),
96 localView(assembler.gridVariables(domainId).curGridVolVars()),
97 localView(assembler.gridVariables(domainId).prevGridVolVars()),
98 localView(assembler.gridVariables(domainId).gridFluxVarsCache()),
99 assembler.localResidual(domainId),
100 (element.partitionType() == Dune::GhostEntity))
101 , couplingManager_(couplingManager)
102 {}
103
108 template<class JacobianMatrixRow, class SubResidualVector, class GridVariablesTuple>
109 void assembleJacobianAndResidual(JacobianMatrixRow& jacRow, SubResidualVector& res, GridVariablesTuple& gridVariables)
110 {
111 auto assembleCouplingBlocks = [&](const auto& residual)
112 {
113 // assemble the coupling blocks
114 using namespace Dune::Hybrid;
115 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](auto&& i)
116 {
117 if constexpr (std::decay_t<decltype(i)>{} != id)
118 this->assembleJacobianCoupling(i, jacRow, residual, gridVariables);
119 });
120 };
121
122 // the coupled model does not support partial reassembly yet
123 const DefaultPartialReassembler* noReassembler = nullptr;
124 ParentType::assembleJacobianAndResidual(jacRow[domainId], res, *std::get<domainId>(gridVariables), noReassembler, assembleCouplingBlocks);
125 }
126
131 template<std::size_t otherId, class JacRow, class GridVariables,
132 typename std::enable_if_t<(otherId == id), int> = 0>
133 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
134 const ElementResidualVector& res, GridVariables& gridVariables)
135 {}
136
140 template<std::size_t otherId, class JacRow, class GridVariables,
141 typename std::enable_if_t<(otherId != id), int> = 0>
142 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
143 const ElementResidualVector& res, GridVariables& gridVariables)
144 {
145 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
146 }
147
151 ElementResidualVector evalLocalSourceResidual(const Element& element, const ElementVolumeVariables& elemVolVars) const
152 {
153 // initialize the residual vector for all scvs in this element
154 ElementResidualVector residual(this->fvGeometry().numScv());
155
156 // evaluate the volume terms (storage + source terms)
157 // forward to the local residual specialized for the discretization methods
158 for (const auto& scv : scvs(this->fvGeometry()))
159 {
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);
164 }
165
166 return residual;
167 }
168
172 ElementResidualVector evalLocalSourceResidual(const Element& neighbor) const
173 { return this->evalLocalSourceResidual(neighbor, implicit ? this->curElemVolVars() : this->prevElemVolVars()); }
174
179 {
180 // get some references for convenience
181 const auto& element = this->element();
182 const auto& curSol = this->curSol(domainId);
183 auto&& fvGeometry = this->fvGeometry();
184 auto&& curElemVolVars = this->curElemVolVars();
185 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
186
187 // bind the caches
188 couplingManager_.bindCouplingContext(domainId, element, this->assembler());
189 fvGeometry.bind(element);
190
191 if constexpr (implicit)
192 {
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]);
197 }
198 else
199 {
200 auto& prevElemVolVars = this->prevElemVolVars();
201 const auto& prevSol = this->assembler().prevSol()[domainId];
202
203 curElemVolVars.bindElement(element, fvGeometry, curSol);
204 prevElemVolVars.bind(element, fvGeometry, prevSol);
205 elemFluxVarsCache.bind(element, fvGeometry, prevElemVolVars);
206 }
207
208 this->elemBcTypes().update(problem(), this->element(), this->fvGeometry());
209 }
210
212 template<std::size_t i = domainId>
213 const Problem& problem(Dune::index_constant<i> dId = domainId) const
214 { return this->assembler().problem(domainId); }
215
217 template<std::size_t i = domainId>
218 const auto& curSol(Dune::index_constant<i> dId = domainId) const
219 { return ParentType::curSol()[dId]; }
220
222 CouplingManager& couplingManager()
223 { return couplingManager_; }
224
225private:
226 CouplingManager& couplingManager_;
227};
228
240template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric, bool implicit = true>
242
249template<std::size_t id, class TypeTag, class Assembler>
250class SubDomainCVFELocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>
251: public SubDomainCVFELocalAssemblerBase<id, TypeTag, Assembler,
253{
254 using ThisType = SubDomainCVFELocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>;
255 using ParentType = SubDomainCVFELocalAssemblerBase<id, TypeTag, Assembler, ThisType, DiffMethod::numeric, /*implicit=*/true>;
258
259 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
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;
266
267 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
268 static constexpr int dim = GridView::dimension;
269
272 static constexpr auto domainI = Dune::index_constant<id>();
273
274public:
278
282 template<class ElemSol>
283 void maybeUpdateCouplingContext(const SubControlVolume& scv, ElemSol& elemSol, const int pvIdx)
284 {
285 this->couplingManager().updateCouplingContext(domainI, *this, domainI, scv.dofIndex(), elemSol[scv.localDofIndex()], pvIdx);
286 }
287
291 template<class JacobianMatrixDiagBlock, class GridVariables>
292 void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector& origResiduals, JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
293 {
294 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *this, origResiduals, A, gridVariables);
295 }
296
301 template<std::size_t otherId, class JacobianBlock, class GridVariables>
302 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
303 const ElementResidualVector& res, GridVariables& gridVariables)
304 {
306 // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. //
308
309 // get some aliases for convenience
310 const auto& element = this->element();
311 const auto& fvGeometry = this->fvGeometry();
312 auto&& curElemVolVars = this->curElemVolVars();
313 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
314
315 // convenience lambda for call to update self
316 auto updateCoupledVariables = [&] ()
317 {
318 // Update ourself after the context has been modified. Depending on the
319 // type of caching, other objects might have to be updated. All ifs can be optimized away.
320 if constexpr (enableGridFluxVarsCache)
321 {
322 if constexpr (enableGridVolVarsCache)
323 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
324 else
325 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, gridVariables.gridFluxVarsCache());
326 }
327 else
328 {
329 if constexpr (enableGridVolVarsCache)
330 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache);
331 else
332 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache);
333 }
334 };
335
336 // get element stencil information
337 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
338 const auto& curSolJ = this->curSol(domainJ);
339 for (const auto globalJ : stencil)
340 {
341 // undeflected privars and privars to be deflected
342 const auto origPriVarsJ = curSolJ[globalJ];
343 auto priVarsJ = origPriVarsJ;
344
345 // the undeflected coupling residual
346 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ);
347
348 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
349 {
350 auto evalCouplingResidual = [&](Scalar priVar)
351 {
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);
356 };
357
358 // derive the residuals numerically
359 ElementResidualVector partialDerivs(fvGeometry.numScv());
360
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);
364
365 NumericDifferentiation::partialDerivative(evalCouplingResidual, origPriVarsJ[pvIdx], partialDerivs, origResidual,
366 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
367
368 // update the global stiffness matrix with the current partial derivatives
369 for (const auto& scv : scvs(fvGeometry))
370 {
371 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
372 {
373 // A[i][col][eqIdx][pvIdx] is the rate of change of
374 // the residual of equation 'eqIdx' at dof 'i'
375 // depending on the primary variable 'pvIdx' at dof
376 // 'col'.
377 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
378
379 // If the dof is coupled by a Dirichlet condition,
380 // set the derived value only once (i.e. overwrite existing values).
381 if (this->elemBcTypes().hasDirichlet())
382 {
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;
388 }
389
390 // enforce internal Dirichlet constraints
391 if constexpr (Problem::enableInternalDirichletConstraints())
392 {
393 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
394 if (internalDirichletConstraints[eqIdx])
395 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
396 }
397
398 }
399 }
400
401 // restore the current element solution
402 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
403
404 // restore the undeflected state of the coupling context
405 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
406 }
407
408 // Restore original state of the flux vars cache and/or vol vars.
409 // This has to be done in case they depend on variables of domainJ before
410 // we continue with the numeric derivative w.r.t the next globalJ. Otherwise,
411 // the next "origResidual" will be incorrect.
412 updateCoupledVariables();
413 }
414 }
415};
416
423template<std::size_t id, class TypeTag, class Assembler>
424class SubDomainCVFELocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>
425: public SubDomainCVFELocalAssemblerBase<id, TypeTag, Assembler,
426 SubDomainCVFELocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>, DiffMethod::numeric, false >
427{
428 using ThisType = SubDomainCVFELocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>;
429 using ParentType = SubDomainCVFELocalAssemblerBase<id, TypeTag, Assembler, ThisType, DiffMethod::numeric, /*implicit=*/false>;
430
431public:
435
442 template<std::size_t otherId, class JacobianBlock, class GridVariables>
443 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
444 const ElementResidualVector& res, GridVariables& gridVariables)
445 {}
446};
447
448} // end namespace Dumux
449
450#endif
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:302
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:37
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:50
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
Definition: adapt.hh:17
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.