GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: /builds/dumux-repositories/dumux/dumux/multidomain/subdomainfclocalassembler.hh
Date: 2024-09-21 20:52:54
Exec Total Coverage
Lines: 78 87 89.7%
Functions: 76 249 30.5%
Branches: 89 172 51.7%

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 Assembly
10 * \ingroup StaggeredDiscretization
11 * \ingroup MultiDomain
12 * \brief An assembler for Jacobian and residual contribution per element (face-centered staggered methods) for multidomain problems
13 */
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
21 #include <dumux/common/properties.hh>
22 #include <dumux/common/parameters.hh>
23 #include <dumux/common/numericdifferentiation.hh>
24 #include <dumux/assembly/numericepsilon.hh>
25 #include <dumux/assembly/diffmethod.hh>
26 #include <dumux/assembly/fclocalassembler.hh>
27
28 namespace Dumux {
29
30 /*!
31 * \ingroup Assembly
32 * \ingroup StaggeredDiscretization
33 * \ingroup MultiDomain
34 * \brief A base class for all face-centered staggered local assemblers
35 * \tparam id the id of the sub domain
36 * \tparam TypeTag the TypeTag
37 * \tparam Assembler the assembler type
38 * \tparam Implementation the actual implementation type
39 * \tparam implicit Specifies whether the time discretization is implicit or not not (i.e. explicit)
40 */
41 template<std::size_t id, class TypeTag, class Assembler, class Implementation, DiffMethod dm, bool implicit>
42 5198220 class SubDomainFaceCenteredLocalAssemblerBase : public FaceCenteredLocalAssembler<TypeTag, Assembler, dm, implicit, Implementation>
43 {
44 using ParentType = FaceCenteredLocalAssembler<TypeTag, Assembler, dm, implicit, Implementation>;
45
46 using Problem = GetPropType<TypeTag, Properties::Problem>;
47 using SolutionVector = typename Assembler::SolutionVector;
48
49 using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
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
65 public:
66 //! export the domain id of this sub-domain
67 static constexpr auto domainId = typename Dune::index_constant<id>();
68 //! pull up constructor of parent class
69 using ParentType::ParentType;
70 //! export element residual vector type
71 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
72
73 // the constructor
74 2599110 explicit SubDomainFaceCenteredLocalAssemblerBase(const Assembler& assembler,
75 const Element& element,
76 const SolutionVector& curSol,
77 CouplingManager& couplingManager)
78 : ParentType(assembler,
79 element,
80 curSol,
81
1/2
✓ Branch 1 taken 2568074 times.
✗ Branch 2 not taken.
2568074 localView(assembler.gridGeometry(domainId)),
82
1/2
✓ Branch 1 taken 2568074 times.
✗ Branch 2 not taken.
2599110 localView(assembler.gridVariables(domainId).curGridVolVars()),
83
1/2
✓ Branch 1 taken 2568074 times.
✗ Branch 2 not taken.
2599110 localView(assembler.gridVariables(domainId).prevGridVolVars()),
84
3/6
✓ Branch 1 taken 2568074 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 2568074 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 2568074 times.
✗ Branch 8 not taken.
7735258 localView(assembler.gridVariables(domainId).gridFluxVarsCache()),
85 assembler.localResidual(domainId),
86 5198220 (element.partitionType() == Dune::GhostEntity))
87
8/16
✓ Branch 2 taken 2568074 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 2568074 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 2568074 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 2599110 times.
✗ Branch 12 not taken.
✓ Branch 14 taken 2568074 times.
✗ Branch 15 not taken.
✓ Branch 16 taken 31036 times.
✓ Branch 17 taken 2568074 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
✓ Branch 20 taken 2568074 times.
✗ Branch 21 not taken.
2599110 , couplingManager_(couplingManager)
88 2599110 {}
89
90 /*!
91 * \brief Computes the derivatives with respect to the given element and adds them
92 * to the global matrix. The element residual is written into the right hand side.
93 */
94 template<class JacobianMatrixRow, class SubResidualVector, class GridVariablesTuple>
95 void assembleJacobianAndResidual(JacobianMatrixRow& jacRow, SubResidualVector& res, GridVariablesTuple& gridVariables)
96 {
97 3788220 auto assembleCouplingBlocks = [&](const auto& residual)
98 {
99 // assemble the coupling blocks
100 using namespace Dune::Hybrid;
101 1894110 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](auto&& i)
102 {
103 if (i != id)
104
0/2
✗ Branch 2 not taken.
✗ Branch 3 not taken.
1894110 this->assembleJacobianCoupling(i, jacRow, residual, gridVariables);
105 });
106 };
107
108 // the coupled model does not support partial reassembly yet
109 1894110 const DefaultPartialReassembler* noReassembler = nullptr;
110
4/8
✓ Branch 1 taken 1894110 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1894110 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1894110 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1894110 times.
✗ Branch 11 not taken.
7576440 ParentType::assembleJacobianAndResidual(jacRow[domainId], res, *std::get<domainId>(gridVariables), noReassembler, assembleCouplingBlocks);
111 }
112
113 /*!
114 * \brief Assemble the entries in a coupling block of the jacobian.
115 * There is no coupling block between a domain and itself.
116 */
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
123 /*!
124 * \brief Assemble the entries in a coupling block of the jacobian.
125 */
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
0/35
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✗ Branch 23 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.
6087776 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
132 }
133
134 /*!
135 * \brief Evaluates the local source term for an element and given element volume variables
136 */
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
155 /*!
156 * \brief Evaluates the local source term depending on time discretization scheme
157 */
158 ElementResidualVector evalLocalSourceResidual(const Element& neighbor) const
159 { return this->evalLocalSourceResidual(neighbor, implicit ? this->curElemVolVars() : this->prevElemVolVars()); }
160
161 /*!
162 * \brief Prepares all local views necessary for local assembly.
163 */
164 2599110 void bindLocalViews()
165 {
166 // get some references for convenience
167 2599110 const auto& element = this->element();
168 5198220 const auto& curSol = this->curSol(domainId);
169 5198220 auto&& fvGeometry = this->fvGeometry();
170 2599110 auto&& curElemVolVars = this->curElemVolVars();
171 2599110 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
172
173 // bind the caches
174 2599110 couplingManager_.bindCouplingContext(domainId, element, this->assembler());
175 2599110 fvGeometry.bind(element);
176
177 if (implicit)
178 {
179 2599110 curElemVolVars.bind(element, fvGeometry, curSol);
180 2599110 elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars);
181
4/4
✓ Branch 0 taken 27500 times.
✓ Branch 1 taken 3536 times.
✓ Branch 2 taken 27500 times.
✓ Branch 3 taken 3536 times.
2630146 if (!this->assembler().isStationaryProblem())
182 110000 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 5198220 this->elemBcTypes().update(problem(), this->element(), this->fvGeometry());
195 2599110 }
196
197 //! return reference to the underlying problem
198 template<std::size_t i = domainId>
199 const Problem& problem(Dune::index_constant<i> dId = domainId) const
200
6/8
✗ Branch 1 not taken.
✓ Branch 2 taken 171 times.
✓ Branch 3 taken 37 times.
✓ Branch 4 taken 7 times.
✓ Branch 5 taken 11 times.
✓ Branch 6 taken 37 times.
✓ Branch 7 taken 7 times.
✗ Branch 8 not taken.
30000142 { return this->assembler().problem(domainId); }
201
202 //! return reference to the underlying problem
203 template<std::size_t i = domainId>
204 const auto& curSol(Dune::index_constant<i> dId = domainId) const
205
2/4
✗ Branch 0 not taken.
✓ Branch 1 taken 160 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 160 times.
8986440 { return ParentType::curSol()[dId]; }
206
207 //! return reference to the coupling manager
208 CouplingManager& couplingManager()
209 { return couplingManager_; }
210
211 private:
212 CouplingManager& couplingManager_; //!< the coupling manager
213 };
214
215 /*!
216 * \ingroup Assembly
217 * \ingroup StaggeredDiscretization
218 * \ingroup MultiDomain
219 * \brief The face-centered staggered scheme multidomain local assembler
220 * \tparam id the id of the sub domain
221 * \tparam TypeTag the TypeTag
222 * \tparam Assembler the assembler type
223 * \tparam DM the numeric differentiation method
224 * \tparam implicit whether the assembler is explicit or implicit in time
225 */
226 template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric, bool implicit = true>
227 class SubDomainFaceCenteredLocalAssembler;
228
229 /*!
230 * \ingroup Assembly
231 * \ingroup StaggeredDiscretization
232 * \ingroup MultiDomain
233 * \brief Face-centered staggered scheme multi domain local assembler using numeric differentiation and implicit time discretization
234 */
235 template<std::size_t id, class TypeTag, class Assembler>
236 2599110 class SubDomainFaceCenteredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>
237 : public SubDomainFaceCenteredLocalAssemblerBase<
238 id, TypeTag, Assembler,
239 SubDomainFaceCenteredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, true>,
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>;
245 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
246 using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>;
247
248 using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
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
254 enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() };
255 enum { dim = GridView::dimension };
256
257 static constexpr bool enableGridFluxVarsCache = GetPropType<TypeTag, Properties::GridVariables>::GridFluxVariablesCache::cachingEnabled;
258 static constexpr bool enableGridVolVarsCache = GetPropType<TypeTag, Properties::GridVariables>::GridVolumeVariables::cachingEnabled;
259 static constexpr auto domainI = Dune::index_constant<id>();
260
261 public:
262 2599110 using ParentType::ParentType;
263 //! export element residual vector type
264 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
265
266 /*!
267 * \brief Update the coupling context for coupled models.
268 */
269 template<class ElemSol>
270 59056 void maybeUpdateCouplingContext(const SubControlVolume& scv, ElemSol& elemSol, const int pvIdx)
271 {
272 361893008 this->couplingManager().updateCouplingContext(domainI, *this, domainI, scv.dofIndex(), elemSol[scv.localDofIndex()], pvIdx);
273 59056 }
274
275 /*!
276 * \brief Update the additional domain derivatives for coupled models.
277 */
278 template<class JacobianMatrixDiagBlock, class GridVariables>
279 void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector& origResiduals, const JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
280 {
281 1894110 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *this, origResiduals, A, gridVariables);
282 }
283
284 /*!
285 * \brief Computes the derivatives with respect to the given element and adds them
286 * to the global matrix.
287 */
288 template<std::size_t otherId, class JacobianBlock, class GridVariables>
289 3010608 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
290 const ElementResidualVector& res, GridVariables& gridVariables)
291 {
292 ////////////////////////////////////////////////////////////////////////////////////////////////////////
293 // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. //
294 ////////////////////////////////////////////////////////////////////////////////////////////////////////
295
296 // get some aliases for convenience
297 3010608 const auto& element = this->element();
298 3010608 const auto& fvGeometry = this->fvGeometry();
299 3010608 auto&& curElemVolVars = this->curElemVolVars();
300 3010608 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
301
302 // the solution vector of the other domain
303 3010608 const auto& curSolJ = this->curSol(domainJ);
304
305 // convenience lambda for call to update self
306 27988654 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 28456020 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 998488 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache);
323 }
324 };
325
326
5/5
✓ Branch 0 taken 124144 times.
✓ Branch 1 taken 8980332 times.
✓ Branch 2 taken 2235240 times.
✓ Branch 3 taken 8949296 times.
✓ Branch 4 taken 2235240 times.
15092412 for (const auto& scv : scvs(fvGeometry))
327 {
328
1/2
✓ Branch 1 taken 1488664 times.
✗ Branch 2 not taken.
12050768 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, scv, domainJ);
329
330
4/4
✓ Branch 0 taken 9583226 times.
✓ Branch 1 taken 9073440 times.
✓ Branch 2 taken 9583226 times.
✓ Branch 3 taken 9073440 times.
47565228 for (const auto globalJ : stencil)
331 {
332
1/2
✓ Branch 1 taken 9583226 times.
✗ Branch 2 not taken.
11075380 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, scv, domainJ, globalJ); // TODO is this necessary?
333 // undeflected privars and privars to be deflected
334 11075380 const auto origPriVarsJ = curSolJ[globalJ];
335 11075380 auto priVarsJ = origPriVarsJ;
336
337
2/2
✓ Branch 0 taken 13902666 times.
✓ Branch 1 taken 9583226 times.
26470200 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
338 {
339 71005484 auto evalCouplingResidual = [&](Scalar priVar)
340 {
341 19974532 priVarsJ[pvIdx] = priVar;
342 13902666 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
343 13902666 updateCoupledVariables();
344 13902666 return this->couplingManager().evalCouplingResidual(domainI, *this, scv, domainJ, globalJ);
345 };
346
347 // derive the residuals numerically
348
4/4
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 28607 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 319 times.
30760552 ElementResidualVector partialDerivs(element.subEntities(1));
349
350
4/4
✓ Branch 0 taken 71 times.
✓ Branch 1 taken 13902595 times.
✓ Branch 2 taken 71 times.
✓ Branch 3 taken 13902595 times.
30789640 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
351
4/6
✓ Branch 0 taken 71 times.
✓ Branch 1 taken 13902595 times.
✓ Branch 3 taken 71 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 71 times.
✗ Branch 7 not taken.
15394820 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
352
4/6
✓ Branch 0 taken 71 times.
✓ Branch 1 taken 13902595 times.
✓ Branch 3 taken 71 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 71 times.
✗ Branch 7 not taken.
15394820 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
353
354
2/4
✓ Branch 1 taken 13902666 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 7830800 times.
✗ Branch 5 not taken.
23225620 NumericDifferentiation::partialDerivative(evalCouplingResidual, origPriVarsJ[pvIdx], partialDerivs, origResidual,
355 23225620 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
356
357 // update the global stiffness matrix with the current partial derivatives
358
2/2
✓ Branch 0 taken 13902666 times.
✓ Branch 1 taken 13902666 times.
30789640 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
3/6
✓ Branch 1 taken 13902666 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 13902666 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 13902666 times.
✗ Branch 8 not taken.
46184460 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
365 }
366
367 // handle Dirichlet boundary conditions
368 // TODO internal constraints
369
6/6
✓ Branch 0 taken 350500 times.
✓ Branch 1 taken 13552166 times.
✓ Branch 2 taken 314711 times.
✓ Branch 3 taken 35789 times.
✓ Branch 4 taken 314711 times.
✓ Branch 5 taken 35789 times.
15394820 if (scv.boundary() && this->elemBcTypes().hasDirichlet())
370 {
371 641656 const auto bcTypes = this->elemBcTypes()[fvGeometry.frontalScvfOnBoundary(scv).localIndex()];
372
4/4
✓ Branch 0 taken 311336 times.
✓ Branch 1 taken 3375 times.
✓ Branch 2 taken 311336 times.
✓ Branch 3 taken 3375 times.
641656 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
2/2
✓ Branch 0 taken 311336 times.
✓ Branch 1 taken 311336 times.
634680 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
378 {
379
2/2
✓ Branch 0 taken 311336 times.
✓ Branch 1 taken 646372 times.
969716 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
380 {
381
2/4
✗ Branch 0 not taken.
✓ Branch 1 taken 646372 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 646372 times.
1304752 if (bcTypes.isCouplingDirichlet(eqIdx))
382 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][pvIdx];
383
2/2
✓ Branch 0 taken 582072 times.
✓ Branch 1 taken 64300 times.
652376 else if (bcTypes.isDirichlet(eqIdx))
384
2/4
✓ Branch 1 taken 646172 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 646172 times.
✗ Branch 5 not taken.
1304352 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
385 }
386 }
387 }
388 }
389
390 // restore the current element solution
391
2/4
✓ Branch 1 taken 7830800 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 7830800 times.
✗ Branch 5 not taken.
23225620 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
392
393 // restore the undeflected state of the coupling context
394
1/2
✓ Branch 1 taken 13902666 times.
✗ Branch 2 not taken.
15394820 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 11075380 updateCoupledVariables();
402 }
403 }
404 3010608 }
405 };
406
407 /*!
408 * \ingroup Assembly
409 * \ingroup StaggeredDiscretization
410 * \ingroup MultiDomain
411 * \brief Face-centered staggered scheme multi domain local assembler using numeric differentiation and explicit time discretization
412 */
413 template<std::size_t id, class TypeTag, class Assembler>
414 class SubDomainFaceCenteredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>
415 : public SubDomainFaceCenteredLocalAssemblerBase<
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>;
423 using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
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>();
429 public:
430 using ParentType::ParentType;
431 //! export element residual vector type
432 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
433
434 /*!
435 * \brief Update the coupling context for coupled models.
436 */
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
443 /*!
444 * \brief Update the additional domain derivatives for coupled models.
445 */
446 template<class JacobianMatrixDiagBlock, class GridVariables>
447 void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector& origResiduals, const JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
448 {}
449
450 /*!
451 * \brief Computes the derivatives with respect to the given element and adds them
452 * to the global matrix.
453 */
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
463