version 3.8
subdomainfclocalassembler.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_FACECENTERED_LOCAL_ASSEMBLER_HH
15#define DUMUX_MULTIDOMAIN_FACECENTERED_LOCAL_ASSEMBLER_HH
16
17#include <dune/common/indices.hh>
18#include <dune/common/hybridutilities.hh>
19#include <dune/grid/common/gridenums.hh>
20
27
28namespace Dumux {
29
41template<std::size_t id, class TypeTag, class Assembler, class Implementation, DiffMethod dm, bool implicit>
42class SubDomainFaceCenteredLocalAssemblerBase : public FaceCenteredLocalAssembler<TypeTag, Assembler, dm, implicit, Implementation>
43{
45
47 using SolutionVector = typename Assembler::SolutionVector;
48
50 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
51 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
52 using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView;
53 using Scalar = typename GridVariables::Scalar;
54
55 using GridGeometry = typename GridVariables::GridGeometry;
56 using FVElementGeometry = typename GridGeometry::LocalView;
57 using SubControlVolume = typename GridGeometry::SubControlVolume;
58 using GridView = typename GridGeometry::GridView;
59 using Element = typename GridView::template Codim<0>::Entity;
60
61 using CouplingManager = typename Assembler::CouplingManager;
62
63 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
64
65public:
67 static constexpr auto domainId = typename Dune::index_constant<id>();
69 using ParentType::ParentType;
71 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
72
73 // the constructor
74 explicit SubDomainFaceCenteredLocalAssemblerBase(const Assembler& assembler,
75 const Element& element,
76 const SolutionVector& curSol,
77 CouplingManager& couplingManager)
78 : ParentType(assembler,
79 element,
80 curSol,
81 localView(assembler.gridGeometry(domainId)),
82 localView(assembler.gridVariables(domainId).curGridVolVars()),
83 localView(assembler.gridVariables(domainId).prevGridVolVars()),
84 localView(assembler.gridVariables(domainId).gridFluxVarsCache()),
85 assembler.localResidual(domainId),
86 (element.partitionType() == Dune::GhostEntity))
87 , couplingManager_(couplingManager)
88 {}
89
94 template<class JacobianMatrixRow, class SubResidualVector, class GridVariablesTuple>
95 void assembleJacobianAndResidual(JacobianMatrixRow& jacRow, SubResidualVector& res, GridVariablesTuple& gridVariables)
96 {
97 auto assembleCouplingBlocks = [&](const auto& residual)
98 {
99 // assemble the coupling blocks
100 using namespace Dune::Hybrid;
101 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](auto&& i)
102 {
103 if (i != id)
104 this->assembleJacobianCoupling(i, jacRow, residual, gridVariables);
105 });
106 };
107
108 // the coupled model does not support partial reassembly yet
109 const DefaultPartialReassembler* noReassembler = nullptr;
110 ParentType::assembleJacobianAndResidual(jacRow[domainId], res, *std::get<domainId>(gridVariables), noReassembler, assembleCouplingBlocks);
111 }
112
117 template<std::size_t otherId, class JacRow, class GridVariables,
118 typename std::enable_if_t<(otherId == id), int> = 0>
119 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
120 const ElementResidualVector& res, GridVariables& gridVariables)
121 {}
122
126 template<std::size_t otherId, class JacRow, class GridVariables,
127 typename std::enable_if_t<(otherId != id), int> = 0>
128 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
129 const ElementResidualVector& res, GridVariables& gridVariables)
130 {
131 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
132 }
133
137 ElementResidualVector evalLocalSourceResidual(const Element& element, const ElementVolumeVariables& elemVolVars) const
138 {
139 // initialize the residual vector for all scvs in this element
140 ElementResidualVector residual(this->fvGeometry().numScv());
141
142 // evaluate the volume terms (storage + source terms)
143 // forward to the local residual specialized for the discretization methods
144 for (auto&& scv : scvs(this->fvGeometry()))
145 {
146 const auto& curVolVars = elemVolVars[scv];
147 auto source = this->localResidual().computeSource(problem(), element, this->fvGeometry(), elemVolVars, scv);
148 source *= -scv.volume()*curVolVars.extrusionFactor();
149 residual[scv.localDofIndex()] = std::move(source);
150 }
151
152 return residual;
153 }
154
158 ElementResidualVector evalLocalSourceResidual(const Element& neighbor) const
159 { return this->evalLocalSourceResidual(neighbor, implicit ? this->curElemVolVars() : this->prevElemVolVars()); }
160
165 {
166 // get some references for convenience
167 const auto& element = this->element();
168 const auto& curSol = this->curSol(domainId);
169 auto&& fvGeometry = this->fvGeometry();
170 auto&& curElemVolVars = this->curElemVolVars();
171 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
172
173 // bind the caches
174 couplingManager_.bindCouplingContext(domainId, element, this->assembler());
175 fvGeometry.bind(element);
176
177 if (implicit)
178 {
179 curElemVolVars.bind(element, fvGeometry, curSol);
180 elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars);
181 if (!this->assembler().isStationaryProblem())
182 this->prevElemVolVars().bindElement(element, fvGeometry, this->assembler().prevSol()[domainId]);
183 }
184 else
185 {
186 auto& prevElemVolVars = this->prevElemVolVars();
187 const auto& prevSol = this->assembler().prevSol()[domainId];
188
189 curElemVolVars.bindElement(element, fvGeometry, curSol);
190 prevElemVolVars.bind(element, fvGeometry, prevSol);
191 elemFluxVarsCache.bind(element, fvGeometry, prevElemVolVars);
192 }
193
194 this->elemBcTypes().update(problem(), this->element(), this->fvGeometry());
195 }
196
198 template<std::size_t i = domainId>
199 const Problem& problem(Dune::index_constant<i> dId = domainId) const
200 { return this->assembler().problem(domainId); }
201
203 template<std::size_t i = domainId>
204 const auto& curSol(Dune::index_constant<i> dId = domainId) const
205 { return ParentType::curSol()[dId]; }
206
208 CouplingManager& couplingManager()
209 { return couplingManager_; }
210
211private:
212 CouplingManager& couplingManager_;
213};
214
226template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric, bool implicit = true>
228
235template<std::size_t id, class TypeTag, class Assembler>
236class SubDomainFaceCenteredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>
238 id, TypeTag, Assembler,
240 DiffMethod::numeric, /*implicit=*/true
241>
242{
243 using ThisType = SubDomainFaceCenteredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>;
244 using ParentType = SubDomainFaceCenteredLocalAssemblerBase<id, TypeTag, Assembler, ThisType, DiffMethod::numeric, /*implicit=*/true>;
247
249 using FVElementGeometry = typename GridGeometry::LocalView;
250 using GridView = typename GridGeometry::GridView;
251 using Element = typename GridView::template Codim<0>::Entity;
252 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
253
255 enum { dim = GridView::dimension };
256
259 static constexpr auto domainI = Dune::index_constant<id>();
260
261public:
265
269 template<class ElemSol>
270 void maybeUpdateCouplingContext(const SubControlVolume& scv, ElemSol& elemSol, const int pvIdx)
271 {
272 this->couplingManager().updateCouplingContext(domainI, *this, domainI, scv.dofIndex(), elemSol[scv.localDofIndex()], pvIdx);
273 }
274
278 template<class JacobianMatrixDiagBlock, class GridVariables>
279 void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector& origResiduals, const JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
280 {
281 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *this, origResiduals, A, gridVariables);
282 }
283
288 template<std::size_t otherId, class JacobianBlock, class GridVariables>
289 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
290 const ElementResidualVector& res, GridVariables& gridVariables)
291 {
293 // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. //
295
296 // get some aliases for convenience
297 const auto& element = this->element();
298 const auto& fvGeometry = this->fvGeometry();
299 auto&& curElemVolVars = this->curElemVolVars();
300 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
301
302 // the solution vector of the other domain
303 const auto& curSolJ = this->curSol(domainJ);
304
305 // convenience lambda for call to update self
306 auto updateCoupledVariables = [&] ()
307 {
308 // Update ourself after the context has been modified. Depending on the
309 // type of caching, other objects might have to be updated. All ifs can be optimized away.
310 if constexpr (enableGridFluxVarsCache)
311 {
312 if constexpr (enableGridVolVarsCache)
313 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
314 else
315 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, gridVariables.gridFluxVarsCache());
316 }
317 else
318 {
319 if constexpr (enableGridVolVarsCache)
320 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache);
321 else
322 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache);
323 }
324 };
325
326 for (const auto& scv : scvs(fvGeometry))
327 {
328 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, scv, domainJ);
329
330 for (const auto globalJ : stencil)
331 {
332 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, scv, domainJ, globalJ); // TODO is this necessary?
333 // undeflected privars and privars to be deflected
334 const auto origPriVarsJ = curSolJ[globalJ];
335 auto priVarsJ = origPriVarsJ;
336
337 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
338 {
339 auto evalCouplingResidual = [&](Scalar priVar)
340 {
341 priVarsJ[pvIdx] = priVar;
342 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
343 updateCoupledVariables();
344 return this->couplingManager().evalCouplingResidual(domainI, *this, scv, domainJ, globalJ);
345 };
346
347 // derive the residuals numerically
348 ElementResidualVector partialDerivs(element.subEntities(1));
349
350 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
351 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
352 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
353
354 NumericDifferentiation::partialDerivative(evalCouplingResidual, origPriVarsJ[pvIdx], partialDerivs, origResidual,
355 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
356
357 // update the global stiffness matrix with the current partial derivatives
358 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
359 {
360 // A[i][col][eqIdx][pvIdx] is the rate of change of
361 // the residual of equation 'eqIdx' at dof 'i'
362 // depending on the primary variable 'pvIdx' at dof
363 // 'col'.
364 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
365 }
366
367 // handle Dirichlet boundary conditions
368 // TODO internal constraints
369 if (scv.boundary() && this->elemBcTypes().hasDirichlet())
370 {
371 const auto bcTypes = this->elemBcTypes()[fvGeometry.frontalScvfOnBoundary(scv).localIndex()];
372 if (bcTypes.hasDirichlet())
373 {
374 // If the dof is coupled by a Dirichlet condition,
375 // set the derived value only once (i.e. overwrite existing values).
376 // For other dofs, add the contribution of the partial derivative.
377 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
378 {
379 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
380 {
381 if (bcTypes.isCouplingDirichlet(eqIdx))
382 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][pvIdx];
383 else if (bcTypes.isDirichlet(eqIdx))
384 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
385 }
386 }
387 }
388 }
389
390 // restore the current element solution
391 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
392
393 // restore the undeflected state of the coupling context
394 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
395 }
396
397 // Restore original state of the flux vars cache and/or vol vars.
398 // This has to be done in case they depend on variables of domainJ before
399 // we continue with the numeric derivative w.r.t the next globalJ. Otherwise,
400 // the next "origResidual" will be incorrect.
401 updateCoupledVariables();
402 }
403 }
404 }
405};
406
413template<std::size_t id, class TypeTag, class Assembler>
414class SubDomainFaceCenteredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>
416 id, TypeTag, Assembler,
417 SubDomainFaceCenteredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>,
418 DiffMethod::numeric, /*implicit=*/false
419>
420{
421 using ThisType = SubDomainFaceCenteredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>;
422 using ParentType = SubDomainFaceCenteredLocalAssemblerBase<id, TypeTag, Assembler, ThisType, DiffMethod::numeric, /*implicit=*/false>;
424 using FVElementGeometry = typename GridGeometry::LocalView;
425 using GridView = typename GridGeometry::GridView;
426 using Element = typename GridView::template Codim<0>::Entity;
427 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
428 static constexpr auto domainI = Dune::index_constant<id>();
429public:
433
437 template<class ElemSol>
438 void maybeUpdateCouplingContext(const SubControlVolume& scv, ElemSol& elemSol, const int pvIdx)
439 {
440 this->couplingManager().updateCouplingContext(domainI, *this, domainI, scv.dofIndex(), elemSol[scv.localDofIndex()], pvIdx);
441 }
442
446 template<class JacobianMatrixDiagBlock, class GridVariables>
447 void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector& origResiduals, const JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
448 {}
449
454 template<std::size_t otherId, class JacobianBlock, class GridVariables>
455 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
456 const ElementResidualVector& res, GridVariables& gridVariables)
457 {}
458};
459
460} // end namespace Dumux
461
462#endif
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 (Face-centered methods)
Definition: fclocalassembler.hh:284
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
Face-centered staggered scheme multi domain local assembler using numeric differentiation and implici...
Definition: subdomainfclocalassembler.hh:242
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:289
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: subdomainfclocalassembler.hh:264
void maybeUpdateCouplingContext(const SubControlVolume &scv, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: subdomainfclocalassembler.hh:270
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, const JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: subdomainfclocalassembler.hh:279
Face-centered staggered scheme multi domain local assembler using numeric differentiation and explici...
Definition: subdomainfclocalassembler.hh:420
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: subdomainfclocalassembler.hh:432
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:455
void maybeUpdateCouplingContext(const SubControlVolume &scv, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: subdomainfclocalassembler.hh:438
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, const JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: subdomainfclocalassembler.hh:447
A base class for all face-centered staggered local assemblers.
Definition: subdomainfclocalassembler.hh:43
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition: subdomainfclocalassembler.hh:158
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: subdomainfclocalassembler.hh:208
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: subdomainfclocalassembler.hh:164
const auto & curSol(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: subdomainfclocalassembler.hh:204
const Problem & problem(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: subdomainfclocalassembler.hh:199
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: subdomainfclocalassembler.hh:95
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:137
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:119
SubDomainFaceCenteredLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: subdomainfclocalassembler.hh:74
static constexpr auto domainId
export the domain id of this sub-domain
Definition: subdomainfclocalassembler.hh:67
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: subdomainfclocalassembler.hh:71
The face-centered staggered scheme multidomain local assembler.
Definition: subdomainfclocalassembler.hh:227
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 (face-centered staggered methods)
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition: diffmethod.hh:25
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: adapt.hh:17
Definition: common/pdesolver.hh:24
A class for numeric differentiation.
An adapter class for local assemblers using numeric differentiation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.