GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: /builds/dumux-repositories/dumux/dumux/experimental/assembly/subdomaincclocalassembler.hh
Date: 2024-05-04 19:09:25
Exec Total Coverage
Lines: 71 80 88.8%
Functions: 5 12 41.7%
Branches: 46 104 44.2%

Line Branch Exec Source
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 //
7 /*!
8 * \file
9 * \ingroup Experimental
10 * \ingroup Assembly
11 * \ingroup CCDiscretization
12 * \ingroup MultiDomain
13 * \brief A multidomain local assembler for Jacobian and residual contribution per element (cell-centered methods)
14 */
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
24 #include <dumux/common/reservedblockvector.hh>
25 #include <dumux/common/properties.hh>
26 #include <dumux/common/parameters.hh>
27 #include <dumux/common/numericdifferentiation.hh>
28 #include <dumux/common/numeqvector.hh>
29 #include <dumux/discretization/elementsolution.hh>
30 #include <dumux/discretization/extrusion.hh>
31 #include <dumux/assembly/numericepsilon.hh>
32 #include <dumux/assembly/diffmethod.hh>
33
34 #include <dumux/experimental/assembly/cclocalassembler.hh>
35
36 namespace Dumux::Experimental {
37
38 /*!
39 * \ingroup Experimental
40 * \ingroup Assembly
41 * \ingroup CCDiscretization
42 * \ingroup MultiDomain
43 * \brief A base class for all multidomain local assemblers
44 * \tparam id the id of the sub domain
45 * \tparam TypeTag the TypeTag
46 * \tparam Assembler the assembler type
47 * \tparam Implementation the actual assembler implementation
48 */
49 template<std::size_t id, class TypeTag, class Assembler, class Implementation, DiffMethod dm>
50 9472 class SubDomainCCLocalAssemblerBase : public Experimental::CCLocalAssembler<TypeTag, Assembler, dm, Implementation>
51 {
52 using ParentType = Experimental::CCLocalAssembler<TypeTag, Assembler, dm, Implementation>;
53
54 using Problem = GetPropType<TypeTag, Properties::Problem>;
55 using LocalResidualValues = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>;
56 using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>;
57 using SolutionVector = typename Assembler::SolutionVector;
58 using ElementBoundaryTypes = GetPropType<TypeTag, Properties::ElementBoundaryTypes>;
59
60 using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
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
77 public:
78 //! export the domain id of this sub-domain
79 static constexpr auto domainId = typename Dune::index_constant<id>();
80 //! the local residual type of this domain
81 using LocalResidual = GetPropType<TypeTag, Properties::LocalResidual>;
82 //! pull up constructor of parent class
83 using ParentType::ParentType;
84
85 //! the constructor
86 4736 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
1/2
✓ Branch 1 taken 4736 times.
✗ Branch 2 not taken.
4736 localView(assembler.gridVariables(domainId).curGridVolVars()),
95
2/4
✓ Branch 1 taken 4736 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 4736 times.
✗ Branch 5 not taken.
9472 localView(assembler.gridVariables(domainId).gridFluxVarsCache()),
96 assembler.localResidual(domainId),
97 9472 (element.partitionType() == Dune::GhostEntity),
98 assembler.isImplicit())
99
7/14
✓ Branch 3 taken 4736 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 4736 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 4736 times.
✗ Branch 10 not taken.
✓ Branch 12 taken 4736 times.
✗ Branch 13 not taken.
✓ Branch 15 taken 4736 times.
✗ Branch 16 not taken.
✓ Branch 18 taken 4736 times.
✗ Branch 19 not taken.
✓ Branch 21 taken 4736 times.
✗ Branch 22 not taken.
4736 , couplingManager_(couplingManager)
100 4736 {}
101
102 /*!
103 * \brief Computes the derivatives with respect to the given element and adds them
104 * to the global matrix. The element residual is written into the right hand side.
105 */
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 6912 auto assembleCouplingBlocks = [&](const auto& residual)
112 {
113 // for the coupling blocks
114 using namespace Dune::Hybrid;
115 10368 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](auto&& i)
116 {
117 if (std::decay_t<decltype(i)>{} != id)
118 3456 this->assembleJacobianCoupling(i, jacRow, residual, gridVariables);
119 });
120 };
121
122 // the coupled model does not support partial reassembly yet
123 3456 const DefaultPartialReassembler* noReassembler = nullptr;
124 6912 ParentType::assembleJacobianAndResidual(
125
1/2
✓ Branch 1 taken 3456 times.
✗ Branch 2 not taken.
3456 jacRow[domainId], res,
126
3/6
✓ Branch 1 taken 3456 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 3456 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 3456 times.
✗ Branch 8 not taken.
10368 *std::get<domainId>(gridVariables),
127 stageParams, temporal, spatial,
128 constrainedDofs,
129 noReassembler,
130 assembleCouplingBlocks
131 );
132 }
133
134 /*!
135 * \brief Assemble the entries in a coupling block of the jacobian.
136 * There is no coupling block between a domain and itself.
137 */
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
0/24
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
✗ Branch 33 not taken.
✗ Branch 34 not taken.
✗ Branch 36 not taken.
✗ Branch 37 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
13824 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
144 }
145
146 /*!
147 * \brief Prepares all local views necessary for local assembly.
148 */
149 4736 void bindLocalViews()
150 {
151 // get some references for convenience
152 4736 const auto& element = this->element();
153 9472 const auto& curSol = this->curSol(domainId);
154 9472 auto&& fvGeometry = this->fvGeometry();
155 4736 auto&& curElemVolVars = this->curElemVolVars();
156 4736 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
157
158 // bind the caches
159 4736 couplingManager_.bindCouplingContext(domainId, element, this->assembler());
160 4736 fvGeometry.bind(element);
161 4736 curElemVolVars.bind(element, fvGeometry, curSol);
162 4736 elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars);
163 4736 }
164
165 //! return reference to the subdomain problem
166 template<std::size_t i = domainId>
167 const Problem& problem(Dune::index_constant<i> dId = domainId) const
168
2/8
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
53572 { return this->assembler().problem(domainId); }
169
170 //! return reference to the subdomain solution
171 template<std::size_t i = domainId>
172 const auto& curSol(Dune::index_constant<i> dId = domainId) const
173 16384 { return ParentType::curSol()[dId]; }
174
175 //! return reference to the coupling manager
176 CouplingManager& couplingManager()
177 { return couplingManager_; }
178
179 private:
180 CouplingManager& couplingManager_; //!< the coupling manager
181 };
182
183 /*!
184 * \ingroup Experimental
185 * \ingroup Assembly
186 * \ingroup CCDiscretization
187 * \ingroup MultiDomain
188 * \brief The cell-centered scheme multidomain local assembler
189 * \tparam id the id of the sub domain
190 * \tparam TypeTag the TypeTag
191 * \tparam Assembler the assembler type
192 * \tparam DM the numeric differentiation method
193 */
194 template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric>
195 class SubDomainCCLocalAssembler;
196
197 /*!
198 * \ingroup Experimental
199 * \ingroup Assembly
200 * \ingroup CCDiscretization
201 * \ingroup MultiDomain
202 * \brief Cell-centered scheme multidomain local assembler using numeric differentiation
203 */
204 template<std::size_t id, class TypeTag, class Assembler>
205 4736 class SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric>
206 : public SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler,
207 SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric>, DiffMethod::numeric>
208 {
209 using ThisType = SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric>;
210 using ParentType = SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler, ThisType, DiffMethod::numeric>;
211 using Problem = GetPropType<TypeTag, Properties::Problem>;
212
213 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
214 using LocalResidualValues = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>;
215
216 using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
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
222 enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() };
223 enum { dim = GridView::dimension };
224
225 static constexpr bool enableGridFluxVarsCache = GetPropType<TypeTag, Properties::GridVariables>::GridFluxVariablesCache::cachingEnabled;
226 static constexpr bool enableGridVolVarsCache = GetPropType<TypeTag, Properties::GridVariables>::GridVolumeVariables::cachingEnabled;
227 static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize;
228 static constexpr auto domainI = Dune::index_constant<id>();
229
230 public:
231 4736 using ParentType::ParentType;
232 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
233
234 /*!
235 * \brief Update the coupling context for coupled models.
236 */
237 template<class ElemSol>
238 13824 void maybeUpdateCouplingContext(const SubControlVolume& scv, ElemSol& elemSol, const int pvIdx)
239 {
240
1/2
✓ Branch 1 taken 13824 times.
✗ Branch 2 not taken.
13824 if (this->assembler().isImplicit())
241 41472 this->couplingManager().updateCouplingContext(domainI, *this, domainI, scv.dofIndex(), elemSol[0], pvIdx);
242 13824 }
243
244 /*!
245 * \brief Update the additional domain derivatives for coupled models.
246 */
247 template<class JacobianMatrixDiagBlock, class GridVariables>
248 void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector& origResiduals, JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
249 {
250 3456 if (this->assembler().isImplicit())
251 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *this, origResiduals, A, gridVariables);
252 }
253
254 /*!
255 * \brief Computes the derivatives with respect to the given element and adds them
256 * to the global matrix.
257 */
258 template<std::size_t otherId, class JacobianBlock, class GridVariables>
259 3456 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
260 const LocalResidualValues& res, GridVariables& gridVariables)
261 {
262
1/2
✗ Branch 1 not taken.
✓ Branch 2 taken 3456 times.
3456 if (!this->assembler().isImplicit())
263 return;
264
265 ////////////////////////////////////////////////////////////////////////////////////////////////////////
266 // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. //
267 ////////////////////////////////////////////////////////////////////////////////////////////////////////
268
269 // get some aliases for convenience
270 3456 const auto& element = this->element();
271 3456 const auto& fvGeometry = this->fvGeometry();
272 3456 const auto& gridGeometry = fvGeometry.gridGeometry();
273 3456 auto&& curElemVolVars = this->curElemVolVars();
274 3456 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
275
276 // get stencil information
277 6912 const auto globalI = gridGeometry.elementMapper().index(element);
278 3456 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
279 3456 const auto& curSolJ = this->curSol(domainJ);
280
281 // make sure the flux variables cache does not contain any artifacts from prior differentiation
282 3456 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
283
284 // convenience lambda for call to update self
285 362880 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 359424 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache);
295 }
296 else
297 {
298 359424 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 718848 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache);
308 }
309 };
310
311
4/4
✓ Branch 0 taken 89856 times.
✓ Branch 1 taken 3456 times.
✓ Branch 2 taken 89856 times.
✓ Branch 3 taken 3456 times.
190080 for (const auto& globalJ : stencil)
312 {
313 // undeflected privars and privars to be deflected
314
1/2
✓ Branch 1 taken 89856 times.
✗ Branch 2 not taken.
89856 const auto origPriVarsJ = curSolJ[globalJ];
315 89856 auto priVarsJ = origPriVarsJ;
316
317
1/2
✓ Branch 1 taken 89856 times.
✗ Branch 2 not taken.
89856 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ)[0];
318
319
2/2
✓ Branch 0 taken 269568 times.
✓ Branch 1 taken 89856 times.
359424 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
320 {
321 1347840 auto evalCouplingResidual = [&](Scalar priVar)
322 {
323 // update the volume variables and the flux var cache
324 269568 priVarsJ[pvIdx] = priVar;
325 269568 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
326 269568 updateCoupledVariables();
327 269568 return this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ)[0];
328 };
329
330 // derive the residuals numerically
331 269568 LocalResidualValues partialDeriv(0.0);
332
4/4
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 269566 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 269566 times.
539136 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
333
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 269566 times.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
269568 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
334
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 269566 times.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
269568 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
335
2/4
✓ Branch 1 taken 269568 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 269568 times.
✗ Branch 5 not taken.
539136 NumericDifferentiation::partialDerivative(evalCouplingResidual, priVarsJ[pvIdx], partialDeriv, origResidual,
336 539136 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
2/2
✓ Branch 0 taken 539136 times.
✓ Branch 1 taken 269568 times.
808704 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
362
3/6
✓ Branch 1 taken 539136 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 539136 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 539136 times.
✗ Branch 8 not taken.
1617408 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
363 }
364
365 // restore the current element solution
366 539136 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
367
368 // restore the undeflected state of the coupling context
369 539136 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
1/2
✓ Branch 1 taken 89856 times.
✗ Branch 2 not taken.
89856 updateCoupledVariables();
376 }
377 }
378 };
379
380 } // end namespace Dumux
381
382 #endif
383