version 3.8
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-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_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
34
36
37namespace Dumux::Experimental {
38
50template<std::size_t id, class TypeTag, class Assembler, class Implementation, DiffMethod dm>
51class SubDomainCVFELocalAssemblerBase : public Experimental::CVFELocalAssembler<TypeTag, Assembler, dm, Implementation>
52{
54
58 using SolutionVector = typename Assembler::SolutionVector;
60
62 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
63 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
64 using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView;
65 using Scalar = typename GridVariables::Scalar;
66
67 using GridGeometry = typename GridVariables::GridGeometry;
68 using FVElementGeometry = typename GridGeometry::LocalView;
69 using SubControlVolume = typename GridGeometry::SubControlVolume;
70 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
71 using Extrusion = Extrusion_t<GridGeometry>;
72 using GridView = typename GridGeometry::GridView;
73 using Element = typename GridView::template Codim<0>::Entity;
74
75 using CouplingManager = typename Assembler::CouplingManager;
76
77 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
78
79public:
81 static constexpr auto domainId = typename Dune::index_constant<id>();
83 using ParentType::ParentType;
85 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
86
87 // the constructor
89 const Assembler& assembler,
90 const Element& element,
91 const SolutionVector& curSol,
92 CouplingManager& couplingManager
93 )
94 : ParentType(assembler,
95 element,
96 curSol,
97 localView(assembler.gridGeometry(domainId)),
98 localView(assembler.gridVariables(domainId).curGridVolVars()),
99 localView(assembler.gridVariables(domainId).gridFluxVarsCache()),
100 assembler.localResidual(domainId),
101 (element.partitionType() == Dune::GhostEntity),
102 assembler.isImplicit())
103 , couplingManager_(couplingManager)
104 {}
105
110 template<class JacobianMatrixRow, class SubResidualVector, class GridVariablesTuple, class StageParams>
111 void assembleJacobianAndResidual(JacobianMatrixRow& jacRow, SubResidualVector& res, GridVariablesTuple& gridVariables,
112 const StageParams& stageParams, SubResidualVector& temporal, SubResidualVector& spatial,
113 SubResidualVector& constrainedDofs)
114 {
115 auto assembleCouplingBlocks = [&](const auto& residual)
116 {
117 // assemble the coupling blocks
118 using namespace Dune::Hybrid;
119 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](auto&& i)
120 {
121 if constexpr (std::decay_t<decltype(i)>{} != id)
122 this->assembleJacobianCoupling(i, jacRow, residual, gridVariables);
123 });
124 };
125
126 // the coupled model does not support partial reassembly yet
127 const DefaultPartialReassembler* noReassembler = nullptr;
128 ParentType::assembleJacobianAndResidual(
129 jacRow[domainId], res,
130 *std::get<domainId>(gridVariables),
131 stageParams, temporal, spatial,
132 constrainedDofs,
133 noReassembler,
134 assembleCouplingBlocks
135 );
136 }
137
142 template<std::size_t otherId, class JacRow, class GridVariables>
143 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
144 const ElementResidualVector& res, GridVariables& gridVariables)
145 {
146 if constexpr (id != otherId)
147 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
148 }
149
154 {
155 // get some references for convenience
156 const auto& element = this->element();
157 const auto& curSol = this->curSol(domainId);
158 auto&& fvGeometry = this->fvGeometry();
159 auto&& curElemVolVars = this->curElemVolVars();
160 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
161
162 // bind the caches
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);
167 else
168 {
169 curElemVolVars.bind(element, fvGeometry, curSol);
170 elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars);
171 }
172
173 this->elemBcTypes().update(problem(), this->element(), this->fvGeometry());
174 }
175
177 template<std::size_t i = domainId>
178 const Problem& problem(Dune::index_constant<i> dId = domainId) const
179 { return this->assembler().problem(domainId); }
180
182 template<std::size_t i = domainId>
183 const auto& curSol(Dune::index_constant<i> dId = domainId) const
184 { return ParentType::curSol()[dId]; }
185
187 CouplingManager& couplingManager()
188 { return couplingManager_; }
189
190private:
191 CouplingManager& couplingManager_;
192};
193
205template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric>
207
215template<std::size_t id, class TypeTag, class Assembler>
216class SubDomainCVFELocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric>
217: public SubDomainCVFELocalAssemblerBase<id, TypeTag, Assembler,
218 SubDomainCVFELocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric>, DiffMethod::numeric>
219{
224
225 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
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;
232
233 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
234 static constexpr int dim = GridView::dimension;
235
238 static constexpr auto domainI = Dune::index_constant<id>();
239
240public:
244
248 template<class ElemSol>
249 void maybeUpdateCouplingContext(const SubControlVolume& scv, ElemSol& elemSol, const int pvIdx)
250 {
251 if (this->assembler().isImplicit())
252 this->couplingManager().updateCouplingContext(domainI, *this, domainI, scv.dofIndex(), elemSol[scv.localDofIndex()], pvIdx);
253 }
254
258 template<class JacobianMatrixDiagBlock, class GridVariables>
259 void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector& origResiduals, JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
260 {
261 if (this->assembler().isImplicit())
262 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *this, origResiduals, A, gridVariables);
263 }
264
269 template<std::size_t otherId, class JacobianBlock, class GridVariables>
270 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
271 const ElementResidualVector& res, GridVariables& gridVariables)
272 {
273 if (!this->assembler().isImplicit())
274 return;
275
277 // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. //
279
280 // get some aliases for convenience
281 const auto& element = this->element();
282 const auto& fvGeometry = this->fvGeometry();
283 auto&& curElemVolVars = this->curElemVolVars();
284 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
285
286 // convenience lambda for call to update self
287 auto updateCoupledVariables = [&] ()
288 {
289 // Update ourself after the context has been modified. Depending on the
290 // type of caching, other objects might have to be updated. All ifs can be optimized away.
291 if constexpr (enableGridFluxVarsCache)
292 {
293 if constexpr (enableGridVolVarsCache)
294 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
295 else
296 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, gridVariables.gridFluxVarsCache());
297 }
298 else
299 {
300 if constexpr (enableGridVolVarsCache)
301 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache);
302 else
303 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache);
304 }
305 };
306
307 // get element stencil information
308 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
309 const auto& curSolJ = this->curSol(domainJ);
310 for (const auto globalJ : stencil)
311 {
312 // undeflected privars and privars to be deflected
313 const auto origPriVarsJ = curSolJ[globalJ];
314 auto priVarsJ = origPriVarsJ;
315
316 // the undeflected coupling residual
317 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ);
318
319 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
320 {
321 auto evalCouplingResidual = [&](Scalar priVar)
322 {
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);
327 };
328
329 // derive the residuals numerically
330 ElementResidualVector partialDerivs(fvGeometry.numScv());
331
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
336 NumericDifferentiation::partialDerivative(evalCouplingResidual, origPriVarsJ[pvIdx], partialDerivs, origResidual,
337 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
338
339 // update the global stiffness matrix with the current partial derivatives
340 for (const auto& scv : scvs(fvGeometry))
341 {
342 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
343 {
344 // A[i][col][eqIdx][pvIdx] is the rate of change of
345 // the residual of equation 'eqIdx' at dof 'i'
346 // depending on the primary variable 'pvIdx' at dof
347 // 'col'.
348 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
349
350 // If the dof is coupled by a Dirichlet condition,
351 // set the derived value only once (i.e. overwrite existing values).
352 if (this->elemBcTypes().hasDirichlet())
353 {
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;
359 }
360
361 // enforce internal Dirichlet constraints
362 if constexpr (Problem::enableInternalDirichletConstraints())
363 {
364 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
365 if (internalDirichletConstraints[eqIdx])
366 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
367 }
368
369 }
370 }
371
372 // restore the current element solution
373 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
374
375 // restore the undeflected state of the coupling context
376 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
377 }
378
379 // Restore original state of the flux vars cache and/or vol vars.
380 // This has to be done in case they depend on variables of domainJ before
381 // we continue with the numeric derivative w.r.t the next globalJ. Otherwise,
382 // the next "origResidual" will be incorrect.
383 updateCoupledVariables();
384 }
385 }
386};
387
388} // end namespace Dumux
389
390#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 (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.