version 3.11-dev
experimental/assembly/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-FileCopyrightText: 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_CVFE_LOCAL_ASSEMBLER_HH
16#define DUMUX_EXPERIMENTAL_SUBDOMAIN_CVFE_LOCAL_ASSEMBLER_HH
17
18#include <type_traits>
19
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> // for GhostEntity
24#include <dune/istl/matrixindexset.hh>
25
26#include <dumux/common/typetraits/localdofs_.hh>
35
37
38namespace Dumux::Experimental {
39
51template<std::size_t id, class TypeTag, class Assembler, class Implementation, DiffMethod dm>
52class SubDomainCVFELocalAssemblerBase : public Experimental::CVFELocalAssembler<TypeTag, Assembler, dm, Implementation>
53{
55
59 using SolutionVector = typename Assembler::SolutionVector;
61
63 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
64 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
65 using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView;
66 using Scalar = typename GridVariables::Scalar;
67
68 using GridGeometry = typename GridVariables::GridGeometry;
69 using FVElementGeometry = typename GridGeometry::LocalView;
70 using SubControlVolume = typename GridGeometry::SubControlVolume;
71 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
72 using Extrusion = Extrusion_t<GridGeometry>;
73 using GridView = typename GridGeometry::GridView;
74 using Element = typename GridView::template Codim<0>::Entity;
75
76 using CouplingManager = typename Assembler::CouplingManager;
77
78 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
79
80public:
82 static constexpr auto domainId = typename Dune::index_constant<id>();
84 using ParentType::ParentType;
86 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
87
88 // the constructor
90 const Assembler& assembler,
91 const Element& element,
92 const SolutionVector& curSol,
93 CouplingManager& couplingManager
94 )
95 : ParentType(assembler,
96 element,
97 curSol,
98 localView(assembler.gridGeometry(domainId)),
99 localView(assembler.gridVariables(domainId).curGridVolVars()),
100 localView(assembler.gridVariables(domainId).gridFluxVarsCache()),
101 assembler.localResidual(domainId),
102 (element.partitionType() == Dune::GhostEntity),
103 assembler.isImplicit())
104 , couplingManager_(couplingManager)
105 {}
106
111 template<class JacobianMatrixRow, class SubResidualVector, class GridVariablesTuple, class StageParams>
112 void assembleJacobianAndResidual(JacobianMatrixRow& jacRow, SubResidualVector& res, GridVariablesTuple& gridVariables,
113 const StageParams& stageParams, SubResidualVector& temporal, SubResidualVector& spatial,
114 SubResidualVector& constrainedDofs)
115 {
116 auto assembleCouplingBlocks = [&](const auto& residual)
117 {
118 // assemble the coupling blocks
119 using namespace Dune::Hybrid;
120 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](auto&& i)
121 {
122 if constexpr (std::decay_t<decltype(i)>{} != id)
123 this->assembleJacobianCoupling(i, jacRow, residual, gridVariables);
124 });
125 };
126
127 // the coupled model does not support partial reassembly yet
128 const DefaultPartialReassembler* noReassembler = nullptr;
129 ParentType::assembleJacobianAndResidual(
130 jacRow[domainId], res,
131 *std::get<domainId>(gridVariables),
132 stageParams, temporal, spatial,
133 constrainedDofs,
134 noReassembler,
135 assembleCouplingBlocks
136 );
137 }
138
143 template<std::size_t otherId, class JacRow, class GridVariables>
144 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
145 const ElementResidualVector& res, GridVariables& gridVariables)
146 {
147 if constexpr (id != otherId)
148 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
149 }
150
155 {
156 // get some references for convenience
157 const auto& element = this->element();
158 const auto& curSol = this->curSol(domainId);
159 auto&& fvGeometry = this->fvGeometry();
160 auto&& curElemVolVars = this->curElemVolVars();
161 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
162
163 // bind the caches
164 couplingManager_.bindCouplingContext(domainId, element, this->assembler());
165 fvGeometry.bind(element);
166 if (std::abs(this->localResidual().spatialWeight()) < 1e-6)
167 curElemVolVars.bindElement(element, fvGeometry, curSol);
168 else
169 {
170 curElemVolVars.bind(element, fvGeometry, curSol);
171 elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars);
172 }
173
174 this->elemBcTypes().update(problem(), this->element(), this->fvGeometry());
175 }
176
178 template<std::size_t i = domainId>
179 const Problem& problem(Dune::index_constant<i> dId = domainId) const
180 { return this->assembler().problem(domainId); }
181
183 template<std::size_t i = domainId>
184 const auto& curSol(Dune::index_constant<i> dId = domainId) const
185 { return ParentType::curSol()[dId]; }
186
188 CouplingManager& couplingManager()
189 { return couplingManager_; }
190
191private:
192 CouplingManager& couplingManager_;
193};
194
206template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric>
208
216template<std::size_t id, class TypeTag, class Assembler>
217class SubDomainCVFELocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric>
218: public SubDomainCVFELocalAssemblerBase<id, TypeTag, Assembler,
219 SubDomainCVFELocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric>, DiffMethod::numeric>
220{
225
226 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
228 using FVElementGeometry = typename GridGeometry::LocalView;
229 using SubControlVolume = typename GridGeometry::SubControlVolume;
230 using GridView = typename GridGeometry::GridView;
231 using Element = typename GridView::template Codim<0>::Entity;
233
234 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
235 static constexpr int dim = GridView::dimension;
236
239 static constexpr auto domainI = Dune::index_constant<id>();
240
241public:
245
249 template<class ScvOrLocalDof, class ElemSol>
250 void maybeUpdateCouplingContext(const ScvOrLocalDof& scvOrLocalDof, ElemSol& elemSol, const int pvIdx)
251 {
252 if (this->assembler().isImplicit())
253 this->couplingManager().updateCouplingContext(domainI, *this, domainI, scvOrLocalDof.dofIndex(), elemSol[Dumux::Detail::LocalDofs::index(scvOrLocalDof)], pvIdx);
254 }
255
259 template<class JacobianMatrixDiagBlock, class GridVariables>
260 void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector& origResiduals, JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
261 {
262 if (this->assembler().isImplicit())
263 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *this, origResiduals, A, gridVariables);
264 }
265
270 template<std::size_t otherId, class JacobianBlock, class GridVariables>
271 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
272 const ElementResidualVector& res, GridVariables& gridVariables)
273 {
274 if (!this->assembler().isImplicit())
275 return;
276
278 // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. //
280
281 // get some aliases for convenience
282 const auto& element = this->element();
283 const auto& fvGeometry = this->fvGeometry();
284 auto&& curElemVolVars = this->curElemVolVars();
285 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
286
287 // convenience lambda for call to update self
288 auto updateCoupledVariables = [&] ()
289 {
290 // Update ourself after the context has been modified. Depending on the
291 // type of caching, other objects might have to be updated. All ifs can be optimized away.
292 if constexpr (enableGridFluxVarsCache)
293 {
294 if constexpr (enableGridVolVarsCache)
295 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
296 else
297 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, gridVariables.gridFluxVarsCache());
298 }
299 else
300 {
301 if constexpr (enableGridVolVarsCache)
302 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache);
303 else
304 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache);
305 }
306 };
307
308 // get element stencil information
309 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
310 const auto& curSolJ = this->curSol(domainJ);
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 // the undeflected coupling residual
318 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ);
319
320 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
321 {
322 auto evalCouplingResidual = [&](Scalar priVar)
323 {
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);
328 };
329
330 // derive the residuals numerically
331 ElementResidualVector partialDerivs(fvGeometry.numScv());
332
333 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
334 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
335 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
336
337 NumericDifferentiation::partialDerivative(evalCouplingResidual, origPriVarsJ[pvIdx], partialDerivs, origResidual,
338 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
339
340 // update the global stiffness matrix with the current partial derivatives
341 for (const auto& scv : scvs(fvGeometry))
342 {
343 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
344 {
345 // A[i][col][eqIdx][pvIdx] is the rate of change of
346 // the residual of equation 'eqIdx' at dof 'i'
347 // depending on the primary variable 'pvIdx' at dof
348 // 'col'.
349 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
350
351 // If the dof is coupled by a Dirichlet condition,
352 // set the derived value only once (i.e. overwrite existing values).
353 if (this->elemBcTypes().hasDirichlet())
354 {
355 const auto bcTypes = this->elemBcTypes().get(fvGeometry, scv);
356 if (bcTypes.isCouplingDirichlet(eqIdx))
357 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][eqIdx];
358 else if (bcTypes.isDirichlet(eqIdx))
359 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
360 }
361
362 // enforce internal Dirichlet constraints
363 if constexpr (Problem::enableInternalDirichletConstraints())
364 {
365 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
366 if (internalDirichletConstraints[eqIdx])
367 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
368 }
369
370 }
371 }
372
373 // restore the current element solution
374 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
375
376 // restore the undeflected state of the coupling context
377 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
378 }
379
380 // Restore original state of the flux vars cache and/or vol vars.
381 // This has to be done in case they depend on variables of domainJ before
382 // we continue with the numeric derivative w.r.t the next globalJ. Otherwise,
383 // the next "origResidual" will be incorrect.
384 updateCoupledVariables();
385 }
386 }
387};
388
389} // end namespace Dumux
390
391#endif
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:37
Definition: partialreassembler.hh:32
An assembler for Jacobian and residual contribution per element (CVFE methods)
Definition: experimental/assembly/cvfelocalassembler.hh:328
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:220
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:271
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:260
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:244
void maybeUpdateCouplingContext(const ScvOrLocalDof &scvOrLocalDof, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:250
A base class for all CVFE subdomain local assemblers.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:53
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:154
SubDomainCVFELocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:89
static constexpr auto domainId
export the domain id of this sub-domain
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:82
const auto & curSol(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:184
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:144
const Problem & problem(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:179
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:112
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:86
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:188
The CVFE scheme multidomain local assembler.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:207
typename ScalarT< X >::type Scalar
export the underlying scalar type
Definition: experimental/common/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:50
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
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition: localdof.hh:79
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.