version 3.8
experimental/assembly/subdomaincclocalassembler.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//
15#ifndef DUMUX_EXPERIMENTAL_SUBDOMAIN_CC_LOCAL_ASSEMBLER_HH
16#define DUMUX_EXPERIMENTAL_SUBDOMAIN_CC_LOCAL_ASSEMBLER_HH
17
18#include <dune/common/indices.hh>
19#include <dune/common/reservedvector.hh>
20#include <dune/common/hybridutilities.hh>
21#include <dune/grid/common/gridenums.hh> // for GhostEntity
22#include <dune/istl/matrixindexset.hh>
23
33
35
36namespace Dumux::Experimental {
37
49template<std::size_t id, class TypeTag, class Assembler, class Implementation, DiffMethod dm>
50class SubDomainCCLocalAssemblerBase : public Experimental::CCLocalAssembler<TypeTag, Assembler, dm, Implementation>
51{
53
57 using SolutionVector = typename Assembler::SolutionVector;
59
61 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
62 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
63 using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView;
64 using Scalar = typename GridVariables::Scalar;
65
66 using GridGeometry = typename GridVariables::GridGeometry;
67 using FVElementGeometry = typename GridGeometry::LocalView;
68 using SubControlVolume = typename GridGeometry::SubControlVolume;
69 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
70 using Extrusion = Extrusion_t<GridGeometry>;
71 using GridView = typename GridGeometry::GridView;
72 using Element = typename GridView::template Codim<0>::Entity;
73
74 using CouplingManager = typename Assembler::CouplingManager;
75 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
76
77public:
79 static constexpr auto domainId = typename Dune::index_constant<id>();
83 using ParentType::ParentType;
84
86 explicit SubDomainCCLocalAssemblerBase(const Assembler& assembler,
87 const Element& element,
88 const SolutionVector& curSol,
89 CouplingManager& couplingManager)
90 : ParentType(assembler,
91 element,
92 curSol,
93 localView(assembler.gridGeometry(domainId)),
94 localView(assembler.gridVariables(domainId).curGridVolVars()),
95 localView(assembler.gridVariables(domainId).gridFluxVarsCache()),
96 assembler.localResidual(domainId),
97 (element.partitionType() == Dune::GhostEntity),
98 assembler.isImplicit())
99 , couplingManager_(couplingManager)
100 {}
101
106 template<class JacobianMatrixRow, class SubResidualVector, class GridVariablesTuple, class StageParams>
107 void assembleJacobianAndResidual(JacobianMatrixRow& jacRow, SubResidualVector& res, GridVariablesTuple& gridVariables,
108 const StageParams& stageParams, SubResidualVector& temporal, SubResidualVector& spatial,
109 SubResidualVector& constrainedDofs)
110 {
111 auto assembleCouplingBlocks = [&](const auto& residual)
112 {
113 // for the coupling blocks
114 using namespace Dune::Hybrid;
115 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](auto&& i)
116 {
117 if (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(
125 jacRow[domainId], res,
126 *std::get<domainId>(gridVariables),
127 stageParams, temporal, spatial,
128 constrainedDofs,
129 noReassembler,
130 assembleCouplingBlocks
131 );
132 }
133
138 template<std::size_t otherId, class JacRow, class GridVariables>
139 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
140 const LocalResidualValues& res, GridVariables& gridVariables)
141 {
142 if constexpr (id != otherId)
143 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
144 }
145
150 {
151 // get some references for convenience
152 const auto& element = this->element();
153 const auto& curSol = this->curSol(domainId);
154 auto&& fvGeometry = this->fvGeometry();
155 auto&& curElemVolVars = this->curElemVolVars();
156 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
157
158 // bind the caches
159 couplingManager_.bindCouplingContext(domainId, element, this->assembler());
160 fvGeometry.bind(element);
161 curElemVolVars.bind(element, fvGeometry, curSol);
162 elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars);
163 }
164
166 template<std::size_t i = domainId>
167 const Problem& problem(Dune::index_constant<i> dId = domainId) const
168 { return this->assembler().problem(domainId); }
169
171 template<std::size_t i = domainId>
172 const auto& curSol(Dune::index_constant<i> dId = domainId) const
173 { return ParentType::curSol()[dId]; }
174
176 CouplingManager& couplingManager()
177 { return couplingManager_; }
178
179private:
180 CouplingManager& couplingManager_;
181};
182
194template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric>
196
204template<std::size_t id, class TypeTag, class Assembler>
205class SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric>
206: public SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler,
207 SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric>, DiffMethod::numeric>
208{
212
215
217 using FVElementGeometry = typename GridGeometry::LocalView;
218 using SubControlVolume = typename GridGeometry::SubControlVolume;
219 using GridView = typename GridGeometry::GridView;
220 using Element = typename GridView::template Codim<0>::Entity;
221
223 enum { dim = GridView::dimension };
224
227 static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize;
228 static constexpr auto domainI = Dune::index_constant<id>();
229
230public:
232 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
233
237 template<class ElemSol>
238 void maybeUpdateCouplingContext(const SubControlVolume& scv, ElemSol& elemSol, const int pvIdx)
239 {
240 if (this->assembler().isImplicit())
241 this->couplingManager().updateCouplingContext(domainI, *this, domainI, scv.dofIndex(), elemSol[0], pvIdx);
242 }
243
247 template<class JacobianMatrixDiagBlock, class GridVariables>
248 void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector& origResiduals, JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
249 {
250 if (this->assembler().isImplicit())
251 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *this, origResiduals, A, gridVariables);
252 }
253
258 template<std::size_t otherId, class JacobianBlock, class GridVariables>
259 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
260 const LocalResidualValues& res, GridVariables& gridVariables)
261 {
262 if (!this->assembler().isImplicit())
263 return;
264
266 // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. //
268
269 // get some aliases for convenience
270 const auto& element = this->element();
271 const auto& fvGeometry = this->fvGeometry();
272 const auto& gridGeometry = fvGeometry.gridGeometry();
273 auto&& curElemVolVars = this->curElemVolVars();
274 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
275
276 // get stencil information
277 const auto globalI = gridGeometry.elementMapper().index(element);
278 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
279 const auto& curSolJ = this->curSol(domainJ);
280
281 // make sure the flux variables cache does not contain any artifacts from prior differentiation
282 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
283
284 // convenience lambda for call to update self
285 auto updateCoupledVariables = [&] ()
286 {
287 // Update ourself after the context has been modified. Depending on the
288 // type of caching, other objects might have to be updated. All ifs can be optimized away.
289 if (enableGridFluxVarsCache)
290 {
291 if (enableGridVolVarsCache)
292 {
293 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
294 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache);
295 }
296 else
297 {
298 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, gridVariables.gridFluxVarsCache());
299 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache);
300 }
301 }
302 else
303 {
304 if (enableGridVolVarsCache)
305 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache);
306 else
307 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache);
308 }
309 };
310
311 for (const auto& globalJ : stencil)
312 {
313 // undeflected privars and privars to be deflected
314 const auto origPriVarsJ = curSolJ[globalJ];
315 auto priVarsJ = origPriVarsJ;
316
317 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ)[0];
318
319 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
320 {
321 auto evalCouplingResidual = [&](Scalar priVar)
322 {
323 // update the volume variables and the flux var cache
324 priVarsJ[pvIdx] = priVar;
325 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
326 updateCoupledVariables();
327 return this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ)[0];
328 };
329
330 // derive the residuals numerically
331 LocalResidualValues partialDeriv(0.0);
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);
335 NumericDifferentiation::partialDerivative(evalCouplingResidual, priVarsJ[pvIdx], partialDeriv, origResidual,
336 epsCoupl(priVarsJ[pvIdx], pvIdx), numDiffMethod);
337
338 // add the current partial derivatives to the global jacobian matrix
339 if constexpr (Problem::enableInternalDirichletConstraints())
340 {
341 const auto& scv = fvGeometry.scv(globalI);
342 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
343 if (internalDirichletConstraints.none())
344 {
345 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
346 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
347 }
348 else
349 {
350 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
351 {
352 if (internalDirichletConstraints[eqIdx])
353 A[globalI][globalJ][eqIdx][pvIdx] = 0.0;
354 else
355 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
356 }
357 }
358 }
359 else
360 {
361 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
362 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
363 }
364
365 // restore the current element solution
366 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
367
368 // restore the undeflected state of the coupling context
369 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
370 }
371
372 // restore original state of the flux vars cache and/or vol vars in case of global caching.
373 // This has to be done in order to guarantee that the undeflected residual computation done
374 // above is correct when jumping to the next coupled dof of domainJ.
375 updateCoupledVariables();
376 }
377 }
378};
379
380} // end namespace Dumux
381
382#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 (cell-centered methods)
Definition: experimental/assembly/cclocalassembler.hh:181
GG GridGeometry
export the grid geometry type
Definition: experimental/discretization/gridvariables.hh:38
Cell-centered scheme multidomain local assembler using numeric differentiation.
Definition: experimental/assembly/subdomaincclocalassembler.hh:208
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const LocalResidualValues &res, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: experimental/assembly/subdomaincclocalassembler.hh:259
void maybeUpdateCouplingContext(const SubControlVolume &scv, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: experimental/assembly/subdomaincclocalassembler.hh:238
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: experimental/assembly/subdomaincclocalassembler.hh:248
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
Definition: experimental/assembly/subdomaincclocalassembler.hh:232
A base class for all multidomain local assemblers.
Definition: experimental/assembly/subdomaincclocalassembler.hh:51
const auto & curSol(Dune::index_constant< i > dId=domainId) const
return reference to the subdomain solution
Definition: experimental/assembly/subdomaincclocalassembler.hh:172
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: experimental/assembly/subdomaincclocalassembler.hh:176
const Problem & problem(Dune::index_constant< i > dId=domainId) const
return reference to the subdomain problem
Definition: experimental/assembly/subdomaincclocalassembler.hh:167
SubDomainCCLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
the constructor
Definition: experimental/assembly/subdomaincclocalassembler.hh:86
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/subdomaincclocalassembler.hh:107
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacRow &jacRow, const LocalResidualValues &res, GridVariables &gridVariables)
Assemble the entries in a coupling block of the jacobian. There is no coupling block between a domain...
Definition: experimental/assembly/subdomaincclocalassembler.hh:139
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: experimental/assembly/subdomaincclocalassembler.hh:149
GetPropType< TypeTag, Properties::LocalResidual > LocalResidual
the local residual type of this domain
Definition: experimental/assembly/subdomaincclocalassembler.hh:81
static constexpr auto domainId
export the domain id of this sub-domain
Definition: experimental/assembly/subdomaincclocalassembler.hh:79
The cell-centered scheme multidomain local assembler.
Definition: experimental/assembly/subdomaincclocalassembler.hh:195
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...
Element solution classes and factory functions.
An assembler for Jacobian and residual contribution per element (cell-centered 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.