GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: dumux/dumux/experimental/assembly/cvfelocalassembler.hh
Date: 2025-06-14 19:21:29
Exec Total Coverage
Lines: 108 161 67.1%
Functions: 24 37 64.9%
Branches: 73 195 37.4%

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-FileCopyrightText: 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 CVFEDiscretization
12 * \brief An assembler for Jacobian and residual contribution per element (CVFE methods)
13 */
14 #ifndef DUMUX_EXPERIMENTAL_CVFE_LOCAL_ASSEMBLER_HH
15 #define DUMUX_EXPERIMENTAL_CVFE_LOCAL_ASSEMBLER_HH
16
17 #include <dune/common/exceptions.hh>
18 #include <dune/common/hybridutilities.hh>
19 #include <dune/common/reservedvector.hh>
20 #include <dune/grid/common/gridenums.hh>
21 #include <dune/istl/matrixindexset.hh>
22 #include <dune/istl/bvector.hh>
23
24 #include <dumux/common/typetraits/typetraits.hh>
25 #include <dumux/common/properties.hh>
26 #include <dumux/common/parameters.hh>
27 #include <dumux/common/numericdifferentiation.hh>
28 #include <dumux/common/typetraits/localdofs_.hh>
29
30 #include <dumux/assembly/numericepsilon.hh>
31 #include <dumux/assembly/diffmethod.hh>
32 #include <dumux/experimental/assembly/fvlocalassemblerbase.hh>
33 #include <dumux/assembly/partialreassembler.hh>
34 #include <dumux/assembly/entitycolor.hh>
35
36 #include <dumux/discretization/cvfe/elementsolution.hh>
37 #include <dumux/discretization/cvfe/localdof.hh>
38
39 #include <dumux/assembly/cvfevolvarsdeflectionpolicy_.hh>
40
41 namespace Dumux::Experimental {
42
43 /*!
44 * \ingroup Experimental
45 * \ingroup Assembly
46 * \ingroup CVFEDiscretization
47 * \brief A base class for all local CVFE assemblers
48 * \tparam TypeTag The TypeTag
49 * \tparam Assembler The assembler type
50 * \tparam Implementation The actual implementation
51 */
52 template<class TypeTag, class Assembler, class Implementation>
53 51968 class CVFELocalAssemblerBase : public Experimental::FVLocalAssemblerBase<TypeTag, Assembler, Implementation>
54 {
55 using ParentType = Experimental::FVLocalAssemblerBase<TypeTag, Assembler, Implementation>;
56 using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>;
57 using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
58 using Problem = GetPropType<TypeTag, Properties::Problem>;
59 using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
60 using SolutionVector = typename Assembler::SolutionVector;
61
62 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
63 static constexpr int dim = GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension;
64
65 public:
66
67 38528 using ParentType::ParentType;
68
69 using LocalResidual = typename ParentType::LocalResidual;
70 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
71
72
73 33792 void bindLocalViews()
74 {
75 33792 ParentType::bindLocalViews();
76 33792 this->elemBcTypes().update(this->asImp_().problem(), this->element(), this->fvGeometry());
77 33792 }
78
79 /*!
80 * \brief Computes the derivatives with respect to the given element and adds them
81 * to the global matrix. The element residual is written into the right hand side.
82 */
83 template <class ResidualVector, class StageParams, class PartialReassembler = DefaultPartialReassembler, class CouplingFunction = Noop>
84 29568 void assembleJacobianAndResidual(JacobianMatrix& jac, ResidualVector& res, GridVariables& gridVariables,
85 const StageParams& stageParams, ResidualVector& temporal, ResidualVector& spatial,
86 ResidualVector& constrainedDofs,
87 const PartialReassembler* partialReassembler = nullptr,
88 const CouplingFunction& maybeAssembleCouplingBlocks = noop)
89 {
90 29568 this->asImp_().bindLocalViews();
91 29568 const auto eIdxGlobal = this->asImp_().problem().gridGeometry().elementMapper().index(this->element());
92
93 29568 this->localResidual().spatialWeight(1.0);
94 29568 this->localResidual().temporalWeight(1.0);
95
96 29568 const auto sWeight = stageParams.spatialWeight(stageParams.size()-1);
97 29568 const auto tWeight = stageParams.temporalWeight(stageParams.size()-1);
98
99 29568 const auto flux = this->evalLocalFluxAndSourceResidual(this->curElemVolVars());
100 29568 const auto storage = this->localResidual().evalStorage(this->fvGeometry(), this->curElemVolVars());
101 29568 ElementResidualVector origResidual(flux.size());
102 29568 origResidual = 0.0;
103
2/2
✓ Branch 0 taken 132096 times.
✓ Branch 1 taken 29568 times.
161664 for (const auto& scv : scvs(this->fvGeometry()))
104 {
105 132096 spatial[scv.dofIndex()] += flux[scv.localDofIndex()];
106 132096 temporal[scv.dofIndex()] += storage[scv.localDofIndex()];
107 528384 origResidual[scv.localDofIndex()] += flux[scv.localDofIndex()]*sWeight + storage[scv.localDofIndex()]*tWeight;
108 264192 res[scv.dofIndex()] += origResidual[scv.localDofIndex()];
109 }
110
111
2/2
✓ Branch 0 taken 3456 times.
✓ Branch 1 taken 26112 times.
29568 this->localResidual().spatialWeight(sWeight);
112
2/2
✓ Branch 0 taken 3456 times.
✓ Branch 1 taken 26112 times.
29568 this->localResidual().temporalWeight(tWeight);
113
114
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 26112 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
26112 if (partialReassembler && partialReassembler->elementColor(eIdxGlobal) == EntityColor::green)
115 {
116 // assemble the coupling blocks for coupled models (does nothing if not coupled)
117 26112 maybeAssembleCouplingBlocks(origResidual);
118 }
119
1/2
✓ Branch 0 taken 29568 times.
✗ Branch 1 not taken.
29568 else if (!this->elementIsGhost())
120 {
121 29568 this->asImp_().assembleJacobian(jac, gridVariables, origResidual, partialReassembler); // forward to the internal implementation
122
123 // assemble the coupling blocks for coupled models (does nothing if not coupled)
124 29568 maybeAssembleCouplingBlocks(origResidual);
125 }
126 else
127 {
128 // Treatment of ghost elements
129 assert(this->elementIsGhost());
130
131 // handle dofs per codimension
132 const auto& gridGeometry = this->asImp_().problem().gridGeometry();
133 Dune::Hybrid::forEach(std::make_integer_sequence<int, dim>{}, [&](auto d)
134 {
135 constexpr int codim = dim - d;
136 const auto& localCoeffs = gridGeometry.feCache().get(this->element().type()).localCoefficients();
137 for (int idx = 0; idx < localCoeffs.size(); ++idx)
138 {
139 const auto& localKey = localCoeffs.localKey(idx);
140
141 // skip if we are not handling this codim right now
142 if (localKey.codim() != codim)
143 continue;
144
145 // do not change the non-ghost entities
146 auto entity = this->element().template subEntity<codim>(localKey.subEntity());
147 if (entity.partitionType() == Dune::InteriorEntity || entity.partitionType() == Dune::BorderEntity)
148 continue;
149
150 // WARNING: this only works if the mapping from codim+subEntity to
151 // global dofIndex is unique (on dof per entity of this codim).
152 // For more general mappings, we should use a proper local-global mapping here.
153 // For example through dune-functions.
154 const auto dofIndex = gridGeometry.dofMapper().index(entity);
155
156 // this might be a vector-valued dof
157 using BlockType = typename JacobianMatrix::block_type;
158 BlockType &J = jac[dofIndex][dofIndex];
159 for (int j = 0; j < BlockType::rows; ++j)
160 J[j][j] = 1.0;
161
162 // set residual for the ghost dof
163 res[dofIndex] = 0;
164 constrainedDofs[dofIndex] = 1;
165 }
166 });
167 }
168
169 134176 auto applyDirichlet = [&] (const auto& scvI,
170 const auto& dirichletValues,
171 const auto eqIdx,
172 const auto pvIdx)
173 {
174 52304 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
175 52304 constrainedDofs[scvI.dofIndex()][eqIdx] = 1;
176
177 52304 auto& row = jac[scvI.dofIndex()];
178
3/4
✗ Branch 0 not taken.
✓ Branch 1 taken 858064 times.
✓ Branch 2 taken 805760 times.
✓ Branch 3 taken 52304 times.
858064 for (auto col = row.begin(); col != row.end(); ++col)
179 1611520 row[col.index()][eqIdx] = 0.0;
180
181 52304 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
182
183 // if a periodic dof has Dirichlet values also apply the same Dirichlet values to the other dof
184 104608 if (this->asImp_().problem().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex()))
185 {
186 const auto periodicDof = this->asImp_().problem().gridGeometry().periodicallyMappedDof(scvI.dofIndex());
187 res[periodicDof][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
188 constrainedDofs[periodicDof][eqIdx] = 1;
189 const auto end = jac[periodicDof].end();
190 for (auto it = jac[periodicDof].begin(); it != end; ++it)
191 (*it) = periodicDof != it.index() ? 0.0 : 1.0;
192 }
193 };
194
195 59136 this->asImp_().enforceDirichletConstraints(applyDirichlet);
196 29568 }
197
198 /*!
199 * \brief Computes the derivatives with respect to the given element and adds them
200 * to the global matrix.
201 */
202 void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables)
203 {
204 this->asImp_().bindLocalViews();
205 this->asImp_().assembleJacobian(jac, gridVariables); // forward to the internal implementation
206
207 auto applyDirichlet = [&] (const auto& scvI,
208 const auto& dirichletValues,
209 const auto eqIdx,
210 const auto pvIdx)
211 {
212 auto& row = jac[scvI.dofIndex()];
213 for (auto col = row.begin(); col != row.end(); ++col)
214 row[col.index()][eqIdx] = 0.0;
215
216 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
217 };
218
219 this->asImp_().enforceDirichletConstraints(applyDirichlet);
220 }
221
222 /*!
223 * \brief Assemble the residual only
224 */
225 template <class ResidualVector>
226 void assembleResidual(ResidualVector& res)
227 {
228 this->asImp_().bindLocalViews();
229 const auto residual = this->evalLocalResidual();
230
231 for (const auto& localDof : localDofs(this->fvGeometry()))
232 res[localDof.dofIndex()] += residual[localDof.index()];
233
234 auto applyDirichlet = [&] (const auto& scvI,
235 const auto& dirichletValues,
236 const auto eqIdx,
237 const auto pvIdx)
238 {
239 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
240 };
241
242 this->asImp_().enforceDirichletConstraints(applyDirichlet);
243 }
244
245 /*!
246 * \brief Assemble the residual only
247 */
248 template<class ResidualVector>
249 8960 void assembleCurrentResidual(ResidualVector& spatialRes, ResidualVector& temporalRes)
250 {
251 8960 this->asImp_().bindLocalViews();
252 8960 const auto flux = this->evalLocalFluxAndSourceResidual(this->curElemVolVars());
253 8960 const auto storage = this->localResidual().evalStorage(this->fvGeometry(), this->curElemVolVars());
254
2/2
✓ Branch 0 taken 40960 times.
✓ Branch 1 taken 8960 times.
49920 for (const auto& scv : scvs(this->fvGeometry()))
255 {
256 40960 spatialRes[scv.dofIndex()] += flux[scv.localDofIndex()];
257 81920 temporalRes[scv.dofIndex()] += storage[scv.localDofIndex()];
258 }
259 8960 }
260
261 //! Enforce Dirichlet constraints
262 template<typename ApplyFunction>
263 29568 void enforceDirichletConstraints(const ApplyFunction& applyDirichlet)
264 {
265 // enforce Dirichlet boundary conditions
266 29568 this->asImp_().evalDirichletBoundaries(applyDirichlet);
267 // take care of internal Dirichlet constraints (if enabled)
268 29568 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
269 }
270
271 /*!
272 * \brief Evaluates Dirichlet boundaries
273 */
274 template< typename ApplyDirichletFunctionType >
275 29568 void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
276 {
277 // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0
278 // and set the residual to (privar - dirichletvalue)
279
2/2
✓ Branch 0 taken 4112 times.
✓ Branch 1 taken 25456 times.
29568 if (this->elemBcTypes().hasDirichlet())
280 {
281
2/2
✓ Branch 0 taken 28544 times.
✓ Branch 1 taken 4112 times.
32656 for (const auto& scvI : scvs(this->fvGeometry()))
282 {
283
2/2
✓ Branch 1 taken 18160 times.
✓ Branch 2 taken 10384 times.
28544 const auto bcTypes = this->elemBcTypes().get(this->fvGeometry(), scvI);
284
2/2
✓ Branch 0 taken 18160 times.
✓ Branch 1 taken 10384 times.
28544 if (bcTypes.hasDirichlet())
285 {
286
1/2
✓ Branch 1 taken 2176 times.
✗ Branch 2 not taken.
18160 const auto dirichletValues = this->asImp_().problem().dirichlet(this->element(), scvI);
287
288 // set the Dirichlet conditions in residual and jacobian
289
2/2
✓ Branch 0 taken 52304 times.
✓ Branch 1 taken 18160 times.
70464 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
290 {
291
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 52304 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
52304 if (bcTypes.isDirichlet(eqIdx))
292 {
293
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 52304 times.
52304 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
294
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 52304 times.
52304 assert(0 <= pvIdx && pvIdx < numEq);
295
1/2
✓ Branch 1 taken 4352 times.
✗ Branch 2 not taken.
52304 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
296 }
297 }
298 }
299 }
300 }
301 29568 }
302
303 /*!
304 * \brief Update the coupling context for coupled models.
305 * \note This does nothing per default (not a coupled model).
306 */
307 template<class... Args>
308 void maybeUpdateCouplingContext(Args&&...) {}
309
310 /*!
311 * \brief Update the additional domain derivatives for coupled models.
312 * \note This does nothing per default (not a coupled model).
313 */
314 template<class... Args>
315 void maybeEvalAdditionalDomainDerivatives(Args&&...) {}
316 };
317
318 /*!
319 * \ingroup Experimental
320 * \ingroup Assembly
321 * \ingroup CVFEDiscretization
322 * \brief An assembler for Jacobian and residual contribution per element (CVFE methods)
323 * \tparam TypeTag The TypeTag
324 * \tparam diffMethod The differentiation method to residual compute derivatives
325 * \tparam Implementation via CRTP
326 */
327 template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, class Implementation = void>
328 class CVFELocalAssembler;
329
330 /*!
331 * \ingroup Experimental
332 * \ingroup Assembly
333 * \ingroup CVFEDiscretization
334 * \brief Control volume finite element local assembler using numeric differentiation
335 */
336 template<class TypeTag, class Assembler, class Implementation>
337
2/8
✓ Branch 0 taken 10368 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 3072 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
51968 class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, Implementation>
338 : public CVFELocalAssemblerBase<TypeTag, Assembler,
339 NonVoidOr<CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, Implementation>, Implementation>>
340 {
341 using ThisType = CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, Implementation>;
342 using ParentType = CVFELocalAssemblerBase<TypeTag, Assembler, NonVoidOr<ThisType, Implementation>>;
343 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
344 using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
345 using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
346 using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>;
347 using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>;
348
349 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
350 static constexpr int dim = GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension;
351
352 static constexpr bool enableGridFluxVarsCache
353 = GridVariables::GridFluxVariablesCache::cachingEnabled;
354 static constexpr bool solutionDependentFluxVarsCache
355 = GridVariables::GridFluxVariablesCache::FluxVariablesCache::isSolDependent;
356
357 public:
358
359 using LocalResidual = typename ParentType::LocalResidual;
360 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
361
3/5
✓ Branch 2 taken 26112 times.
✗ Branch 3 not taken.
✓ Branch 6 taken 7680 times.
✗ Branch 7 not taken.
✓ Branch 1 taken 4736 times.
38528 using ParentType::ParentType;
362
363 template <class PartialReassembler = DefaultPartialReassembler>
364 29568 void assembleJacobian(JacobianMatrix& A, GridVariables& gridVariables,
365 const ElementResidualVector& origResiduals,
366 const PartialReassembler* partialReassembler = nullptr)
367 {
368
1/2
✓ Branch 0 taken 29568 times.
✗ Branch 1 not taken.
29568 if (this->isImplicit())
369 29568 assembleJacobianImplicit_(A, gridVariables, origResiduals, partialReassembler);
370 else
371 assembleJacobianExplicit_(A, gridVariables, origResiduals, partialReassembler);
372 29568 }
373
374 private:
375 /*!
376 * \brief Computes the derivatives with respect to the given element and adds them
377 * to the global matrix.
378 *
379 * \return The element residual at the current solution.
380 */
381 template <class PartialReassembler = DefaultPartialReassembler>
382 29568 void assembleJacobianImplicit_(JacobianMatrix& A, GridVariables& gridVariables,
383 const ElementResidualVector& origResiduals,
384 const PartialReassembler* partialReassembler = nullptr)
385 {
386 // get some aliases for convenience
387
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 29560 times.
29568 const auto& element = this->element();
388
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 29560 times.
29568 const auto& fvGeometry = this->fvGeometry();
389
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 29560 times.
29568 const auto& curSol = this->asImp_().curSol();
390
391
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 29560 times.
29568 auto&& curElemVolVars = this->curElemVolVars();
392
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 29560 times.
29568 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
393
394 //////////////////////////////////////////////////////////////////////////////////////////////////
395 // //
396 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
397 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
398 // actual element. In the actual element we evaluate the derivative of the entire residual. //
399 // //
400 //////////////////////////////////////////////////////////////////////////////////////////////////
401
402 // if all volvars in the stencil have to be updated or if it's enough to only update the
403 // volVars for the scv whose associated dof has been deflected
404
5/6
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 29560 times.
✓ Branch 3 taken 5 times.
✓ Branch 4 taken 3 times.
✓ Branch 6 taken 5 times.
✗ Branch 7 not taken.
29568 static const bool updateAllVolVars = getParamFromGroup<bool>(
405
1/2
✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
10 this->asImp_().problem().paramGroup(), "Assembly.BoxVolVarsDependOnAllElementDofs", false
406 );
407
408 // create the element solution
409 29568 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
410
411 // create the vector storing the partial derivatives
412 29568 ElementResidualVector partialDerivs(Dumux::Detail::LocalDofs::numLocalDofs(fvGeometry));
413
414 29568 auto deflectionPolicy = Detail::CVFE::makeVariablesDeflectionPolicy(
415 29568 gridVariables.curGridVolVars(),
416 curElemVolVars,
417 fvGeometry,
418 updateAllVolVars
419 );
420
421 293760 auto assembleDerivative = [&, this](const auto& localDof)
422 {
423 // dof index and corresponding actual pri vars
424 132096 const auto dofIdx = localDof.dofIndex();
425 132096 const auto localIdx = localDof.index();
426 132096 deflectionPolicy.store(localDof);
427
428 // calculate derivatives w.r.t to the privars at the dof at hand
429
2/2
✓ Branch 0 taken 291840 times.
✓ Branch 1 taken 132096 times.
423936 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
430 {
431 291840 partialDerivs = 0.0;
432
433 583680 auto evalResiduals = [&](Scalar priVar)
434 {
435 // update the volume variables and compute element residual
436 291840 elemSol[localIdx][pvIdx] = priVar;
437 291840 deflectionPolicy.update(elemSol, localDof, this->asImp_().problem());
438 if constexpr (solutionDependentFluxVarsCache)
439 {
440 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
441 if constexpr (enableGridFluxVarsCache)
442 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
443 }
444 291840 this->asImp_().maybeUpdateCouplingContext(localDof, elemSol, pvIdx);
445 291840 return this->evalLocalResidual();
446 };
447
448 // derive the residuals numerically
449
5/6
✓ Branch 0 taken 6 times.
✓ Branch 1 taken 291834 times.
✓ Branch 3 taken 5 times.
✓ Branch 4 taken 1 times.
✓ Branch 6 taken 5 times.
✗ Branch 7 not taken.
291840 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
450
5/6
✓ Branch 0 taken 7 times.
✓ Branch 1 taken 291833 times.
✓ Branch 3 taken 5 times.
✓ Branch 4 taken 2 times.
✓ Branch 6 taken 5 times.
✗ Branch 7 not taken.
291840 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
451
1/2
✓ Branch 1 taken 291840 times.
✗ Branch 2 not taken.
291840 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[localIdx][pvIdx], partialDerivs, origResiduals,
452 291840 eps_(elemSol[localIdx][pvIdx], pvIdx), numDiffMethod);
453
454 // update the global stiffness matrix with the current partial derivatives
455
2/2
✓ Branch 0 taken 1499136 times.
✓ Branch 1 taken 291840 times.
1790976 for (const auto& localDofJ : localDofs(fvGeometry))
456 {
457 // don't add derivatives for green dofs
458 835584 if (!partialReassembler
459
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 835584 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
835584 || partialReassembler->dofColor(localDofJ.dofIndex()) != EntityColor::green)
460 {
461
2/2
✓ Branch 0 taken 3661824 times.
✓ Branch 1 taken 1499136 times.
5160960 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
462 {
463 // A[i][col][eqIdx][pvIdx] is the rate of change of
464 // the residual of equation 'eqIdx' at dof 'i'
465 // depending on the primary variable 'pvIdx' at dof
466 // 'col'.
467
1/2
✓ Branch 1 taken 3661824 times.
✗ Branch 2 not taken.
3661824 A[localDofJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[localDofJ.index()][eqIdx];
468 }
469 }
470 }
471
472 // restore the original state of the scv's volume variables
473 291840 deflectionPolicy.restore(localDof);
474
475 // restore the original element solution
476
1/2
✓ Branch 1 taken 82944 times.
✗ Branch 2 not taken.
291840 elemSol[localIdx][pvIdx] = curSol[localDof.dofIndex()][pvIdx];
477
1/2
✓ Branch 1 taken 82944 times.
✗ Branch 2 not taken.
291840 this->asImp_().maybeUpdateCouplingContext(localDof, elemSol, pvIdx);
478 }
479 };
480
481 // calculation of the derivatives
482
2/2
✓ Branch 0 taken 132096 times.
✓ Branch 1 taken 29568 times.
161664 for (const auto& localDof : localDofs(fvGeometry))
483 132096 assembleDerivative(localDof);
484
485 // restore original state of the flux vars cache in case of global caching.
486 // In the case of local caching this is obsolete because the elemFluxVarsCache used here goes out of scope after this.
487 if constexpr (enableGridFluxVarsCache)
488 10368 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
489
490 // evaluate additional derivatives that might arise from the coupling (no-op if not coupled)
491 29568 this->asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables);
492 29568 }
493
494 /*!
495 * \brief Computes the derivatives with respect to the given element and adds them
496 * to the global matrix.
497 * \return The element residual at the current solution.
498 */
499 template <class PartialReassembler = DefaultPartialReassembler>
500 void assembleJacobianExplicit_(JacobianMatrix& A, GridVariables& gridVariables,
501 const ElementResidualVector& origResiduals,
502 const PartialReassembler* partialReassembler = nullptr)
503 {
504 if (partialReassembler)
505 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
506
507 // get some aliases for convenience
508 const auto& element = this->element();
509 const auto& fvGeometry = this->fvGeometry();
510 const auto& curSol = this->asImp_().curSol();
511 auto&& curElemVolVars = this->curElemVolVars();
512
513 // create the element solution
514 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
515
516 // create the vector storing the partial derivatives
517 ElementResidualVector partialDerivs(Dumux::Detail::LocalDofs::numLocalDofs(fvGeometry));
518
519 // calculation of the derivatives
520 for (const auto& scv : scvs(fvGeometry))
521 {
522 // dof index and corresponding actual pri vars
523 const auto dofIdx = scv.dofIndex();
524 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
525 const VolumeVariables origVolVars(curVolVars);
526
527 // calculate derivatives w.r.t to the privars at the dof at hand
528 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
529 {
530 partialDerivs = 0.0;
531
532 auto evalStorage = [&](Scalar priVar)
533 {
534 elemSol[scv.localDofIndex()][pvIdx] = priVar;
535 curVolVars.update(elemSol, this->asImp_().problem(), element, scv);
536 return this->evalStorage();
537 };
538
539 // derive the residuals numerically
540 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
541 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
542 NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals,
543 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
544
545 // update the global stiffness matrix with the current partial derivatives
546 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
547 {
548 // A[i][col][eqIdx][pvIdx] is the rate of change of
549 // the residual of equation 'eqIdx' at dof 'i'
550 // depending on the primary variable 'pvIdx' at dof
551 // 'col'.
552 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
553 }
554
555 // restore the original state of the scv's volume variables
556 curVolVars = origVolVars;
557
558 // restore the original element solution
559 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
560 }
561 }
562 }
563 };
564
565 } // end namespace Dumux
566
567 #endif
568