GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: /builds/dumux-repositories/dumux/dumux/assembly/cclocalassembler.hh
Date: 2024-05-04 19:09:25
Exec Total Coverage
Lines: 128 137 93.4%
Functions: 358 448 79.9%
Branches: 133 224 59.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-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 CCDiscretization
11 * \brief An assembler for Jacobian and residual contribution per element (cell-centered methods)
12 */
13 #ifndef DUMUX_CC_LOCAL_ASSEMBLER_HH
14 #define DUMUX_CC_LOCAL_ASSEMBLER_HH
15
16 #include <dune/common/reservedvector.hh>
17 #include <dune/grid/common/gridenums.hh> // for GhostEntity
18 #include <dune/istl/matrixindexset.hh>
19
20 #include <dumux/common/reservedblockvector.hh>
21 #include <dumux/common/properties.hh>
22 #include <dumux/common/parameters.hh>
23 #include <dumux/common/numericdifferentiation.hh>
24 #include <dumux/common/numeqvector.hh>
25 #include <dumux/assembly/numericepsilon.hh>
26 #include <dumux/assembly/diffmethod.hh>
27 #include <dumux/assembly/fvlocalassemblerbase.hh>
28 #include <dumux/assembly/entitycolor.hh>
29 #include <dumux/assembly/partialreassembler.hh>
30 #include <dumux/discretization/fluxstencil.hh>
31 #include <dumux/discretization/cellcentered/elementsolution.hh>
32
33 namespace Dumux {
34
35 /*!
36 * \ingroup Assembly
37 * \ingroup CCDiscretization
38 * \brief A base class for all local cell-centered assemblers
39 * \tparam TypeTag The TypeTag
40 * \tparam Assembler The assembler type
41 * \tparam Implementation The actual implementation
42 * \tparam implicit Specifies whether the time discretization is implicit or not (i.e. explicit)
43 */
44 template<class TypeTag, class Assembler, class Implementation, bool implicit>
45 31534617 class CCLocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit>
46 {
47 using ParentType = FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit>;
48 using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView;
49 using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>;
50 using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
51 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
52 using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>;
53
54 public:
55
56 27578891 using ParentType::ParentType;
57
58 /*!
59 * \brief Computes the derivatives with respect to the given element and adds them
60 * to the global matrix. The element residual is written into the right hand side.
61 */
62 template <class ResidualVector, class PartialReassembler = DefaultPartialReassembler>
63 27917227 void assembleJacobianAndResidual(JacobianMatrix& jac, ResidualVector& res, GridVariables& gridVariables,
64 const PartialReassembler* partialReassembler)
65 {
66 27917227 this->asImp_().bindLocalViews();
67
3/6
✗ Branch 0 not taken.
✓ Branch 1 taken 9518 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 9518 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 9518 times.
83751681 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
68 if (partialReassembler
69
6/6
✓ Branch 0 taken 5405370 times.
✓ Branch 1 taken 19555825 times.
✓ Branch 2 taken 1927564 times.
✓ Branch 3 taken 3477806 times.
✓ Branch 4 taken 1927564 times.
✓ Branch 5 taken 3477806 times.
25194667 && partialReassembler->elementColor(globalI) == EntityColor::green)
70 {
71
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 174808 times.
2718904 res[globalI] = this->asImp_().evalLocalResidual()[0]; // forward to the internal implementation
72 }
73 else
74 {
75 25922934 res[globalI] = this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation
76 }
77 27917227 }
78
79 /*!
80 * \brief Computes the derivatives with respect to the given element and adds them
81 * to the global matrix.
82 */
83 void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables)
84 {
85
1/4
✓ Branch 1 taken 180764 times.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
180764 this->asImp_().bindLocalViews();
86
1/4
✓ Branch 1 taken 180764 times.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
180764 this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation
87 }
88
89 /*!
90 * \brief Assemble the residual only
91 */
92 template <class ResidualVector>
93 1043672 void assembleResidual(ResidualVector& res)
94 {
95 1043672 this->asImp_().bindLocalViews();
96 3131016 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
97
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 501640 times.
1043672 res[globalI] = this->asImp_().evalLocalResidual()[0]; // forward to the internal implementation
98
99 using Problem = GetPropType<TypeTag, Properties::Problem>;
100 if constexpr (Problem::enableInternalDirichletConstraints())
101 {
102 204 const auto applyDirichlet = [&] (const auto& scvI,
103 const auto& dirichletValues,
104 const auto eqIdx,
105 4 const auto pvIdx)
106 {
107 4 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
108 };
109
110 200 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
111 }
112 1043672 }
113 };
114
115 /*!
116 * \ingroup Assembly
117 * \ingroup CCDiscretization
118 * \brief An assembler for Jacobian and residual contribution per element (cell-centered methods)
119 * \tparam TypeTag The TypeTag
120 * \tparam diffMethod The differentiation method to residual compute derivatives
121 * \tparam implicit Specifies whether the time discretization is implicit or not not (i.e. explicit)
122 */
123 template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, bool implicit = true>
124 class CCLocalAssembler;
125
126 /*!
127 * \ingroup Assembly
128 * \ingroup CCDiscretization
129 * \brief Cell-centered scheme local assembler using numeric differentiation and implicit time discretization
130 */
131 template<class TypeTag, class Assembler>
132 18619569 class CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>
133 : public CCLocalAssemblerBase<TypeTag, Assembler,
134 CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true>, true >
135 {
136 using ThisType = CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true>;
137 using ParentType = CCLocalAssemblerBase<TypeTag, Assembler, ThisType, true>;
138 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
139 using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>;
140 using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity;
141 using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
142 using FVElementGeometry = typename GridGeometry::LocalView;
143 using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
144 using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>;
145 using Problem = typename GridVariables::GridVolumeVariables::Problem;
146
147 enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() };
148 enum { dim = GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension };
149
150 using FluxStencil = Dumux::FluxStencil<FVElementGeometry>;
151 static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize;
152 static constexpr bool enableGridFluxVarsCache = GetPropType<TypeTag, Properties::GridVariables>::GridFluxVariablesCache::cachingEnabled;
153
154 public:
155
156
2/8
✓ Branch 2 taken 16886209 times.
✗ Branch 3 not taken.
✓ Branch 6 taken 43604 times.
✗ Branch 7 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
16929813 using ParentType::ParentType;
157
158 /*!
159 * \brief Computes the derivatives with respect to the given element and adds them
160 * to the global matrix.
161 *
162 * \return The element residual at the current solution.
163 */
164 14959905 NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables)
165 {
166 //////////////////////////////////////////////////////////////////////////////////////////////////
167 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
168 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
169 // actual element. In the actual element we evaluate the derivative of the entire residual. //
170 //////////////////////////////////////////////////////////////////////////////////////////////////
171
172 // get some aliases for convenience
173 14959905 const auto& element = this->element();
174 14959905 const auto& fvGeometry = this->fvGeometry();
175 14959905 const auto& gridGeometry = this->assembler().gridGeometry();
176 14959905 auto&& curElemVolVars = this->curElemVolVars();
177 14959905 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
178
179 // get stencil information
180 29919810 const auto globalI = gridGeometry.elementMapper().index(element);
181 14959905 const auto& connectivityMap = gridGeometry.connectivityMap();
182 29919810 const auto numNeighbors = connectivityMap[globalI].size();
183
184 // container to store the neighboring elements
185 14959905 Dune::ReservedVector<Element, maxElementStencilSize> neighborElements;
186 14959905 neighborElements.resize(numNeighbors);
187
188 // assemble the undeflected residual
189 using Residuals = ReservedBlockVector<NumEqVector, maxElementStencilSize>;
190 29919810 Residuals origResiduals(numNeighbors + 1); origResiduals = 0.0;
191
2/3
✓ Branch 0 taken 31680 times.
✓ Branch 1 taken 7243173 times.
✗ Branch 2 not taken.
26918947 origResiduals[0] = this->evalLocalResidual()[0];
192
193 // lambda for convenient evaluation of the fluxes across scvfs in the neighbors
194 // if the neighbor is a ghost we don't want to add anything to their residual
195 // so we return 0 and omit computing the flux
196 485751317 auto evalNeighborFlux = [&] (const auto& neighbor, const auto& scvf)
197 {
198
3/5
✓ Branch 0 taken 5803176 times.
✓ Branch 1 taken 251156 times.
✓ Branch 2 taken 229357432 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
235411764 if (neighbor.partitionType() == Dune::GhostEntity)
199 346196 return NumEqVector(0.0);
200 else
201 359624524 return this->localResidual().evalFlux(this->problem(),
202 neighbor,
203 this->fvGeometry(),
204 this->curElemVolVars(),
205 1798122620 this->elemFluxVarsCache(), scvf);
206 };
207
208 // get the elements in which we need to evaluate the fluxes
209 // and calculate these in the undeflected state
210 14959905 unsigned int j = 1;
211
4/4
✓ Branch 0 taken 69118123 times.
✓ Branch 1 taken 14959905 times.
✓ Branch 2 taken 69118123 times.
✓ Branch 3 taken 14959905 times.
198075866 for (const auto& dataJ : connectivityMap[globalI])
212 {
213
1/2
✓ Branch 1 taken 20825638 times.
✗ Branch 2 not taken.
106665679 neighborElements[j-1] = gridGeometry.element(dataJ.globalJ);
214
4/4
✓ Branch 0 taken 113380381 times.
✓ Branch 1 taken 69118123 times.
✓ Branch 2 taken 70211802 times.
✓ Branch 3 taken 25949544 times.
390946552 for (const auto scvfIdx : dataJ.scvfsJ)
215
5/8
✓ Branch 1 taken 32105200 times.
✓ Branch 2 taken 8044478 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 32105200 times.
✓ Branch 5 taken 8044478 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 32105200 times.
✗ Branch 8 not taken.
410432362 origResiduals[j] += evalNeighborFlux(neighborElements[j-1], fvGeometry.scvf(scvfIdx));
216
217 69118123 ++j;
218 }
219
220 // reference to the element's scv (needed later) and corresponding vol vars
221
2/3
✓ Branch 0 taken 4672609 times.
✓ Branch 1 taken 2363892 times.
✗ Branch 2 not taken.
14959905 const auto& scv = fvGeometry.scv(globalI);
222
2/4
✓ Branch 1 taken 2363892 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 2363892 times.
✗ Branch 5 not taken.
29919810 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
223
224 // save a copy of the original privars and vol vars in order
225 // to restore the original solution after deflection
226 14959905 const auto& curSol = this->curSol();
227
1/2
✓ Branch 1 taken 4104085 times.
✗ Branch 2 not taken.
14959905 const auto origPriVars = curSol[globalI];
228 14959905 const auto origVolVars = curVolVars;
229
230 // element solution container to be deflected
231
1/2
✓ Branch 1 taken 4104085 times.
✗ Branch 2 not taken.
14959905 auto elemSol = elementSolution(element, curSol, gridGeometry);
232
233 // derivatives in the neighbors with respect to the current elements
234 // in index 0 we save the derivative of the element residual with respect to it's own dofs
235 14959905 Residuals partialDerivs(numNeighbors + 1);
236
237
2/2
✓ Branch 0 taken 34717902 times.
✓ Branch 1 taken 14959905 times.
49677807 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
238 {
239 34717902 partialDerivs = 0.0;
240
241 161852928 auto evalResiduals = [&](Scalar priVar)
242 {
243 34910862 Residuals partialDerivsTmp(numNeighbors + 1);
244 34910862 partialDerivsTmp = 0.0;
245 // update the volume variables and the flux var cache
246 57308408 elemSol[0][pvIdx] = priVar;
247 69821724 curVolVars.update(elemSol, this->problem(), element, scv);
248 107960087 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
249 if (enableGridFluxVarsCache)
250
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2341340 times.
11996812 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
251
252 // calculate the residual with the deflected primary variables
253
2/2
✓ Branch 0 taken 63360 times.
✓ Branch 1 taken 11895682 times.
58828946 partialDerivsTmp[0] = this->evalLocalResidual()[0];
254
255 // calculate the fluxes in the neighbors with the deflected primary variables
256
2/2
✓ Branch 0 taken 158859089 times.
✓ Branch 1 taken 34910862 times.
193769951 for (std::size_t k = 0; k < numNeighbors; ++k)
257
4/4
✓ Branch 0 taken 246495299 times.
✓ Branch 1 taken 158859089 times.
✓ Branch 2 taken 139011648 times.
✓ Branch 3 taken 51375438 times.
1179802392 for (auto scvfIdx : connectivityMap[globalI][k].scvfsJ)
258 903256818 partialDerivsTmp[k+1] += evalNeighborFlux(neighborElements[k], fvGeometry.scvf(scvfIdx));
259
260 34910862 return partialDerivsTmp;
261 };
262
263 // derive the residuals numerically
264
7/10
✓ Branch 0 taken 199 times.
✓ Branch 1 taken 34717703 times.
✓ Branch 3 taken 139 times.
✓ Branch 4 taken 60 times.
✓ Branch 6 taken 139 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 139 times.
✗ Branch 10 not taken.
✓ Branch 12 taken 139 times.
✗ Branch 13 not taken.
34717902 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
265
7/10
✓ Branch 0 taken 184 times.
✓ Branch 1 taken 34717718 times.
✓ Branch 3 taken 139 times.
✓ Branch 4 taken 45 times.
✓ Branch 6 taken 139 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 139 times.
✗ Branch 10 not taken.
✓ Branch 12 taken 139 times.
✗ Branch 13 not taken.
34717902 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
266
2/4
✓ Branch 1 taken 34717902 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 32526578 times.
✗ Branch 5 not taken.
67244480 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[0][pvIdx], partialDerivs, origResiduals,
267 67244480 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
268
269 // Correct derivative for ghost elements, i.e. set a 1 for the derivative w.r.t. the
270 // current primary variable and a 0 elsewhere. As we always solve for a delta of the
271 // solution with respect to the initial one, this results in a delta of 0 for ghosts.
272
2/2
✓ Branch 0 taken 171215 times.
✓ Branch 1 taken 34546687 times.
34717902 if (this->elementIsGhost())
273 {
274 171215 partialDerivs[0] = 0.0;
275 510685 partialDerivs[0][pvIdx] = 1.0;
276 }
277
278 // restore the original state of the scv's volume variables
279 34717902 curVolVars = origVolVars;
280
281 // restore the current element solution
282 67244480 elemSol[0][pvIdx] = origPriVars[pvIdx];
283
284 // add the current partial derivatives to the global jacobian matrix
285 // no special treatment is needed if globalJ is a ghost because then derivatives have been assembled to 0 above
286 if constexpr (Problem::enableInternalDirichletConstraints())
287 {
288 // check if own element has internal Dirichlet constraint
289 const auto internalDirichletConstraintsOwnElement = this->problem().hasInternalDirichletConstraint(this->element(), scv);
290 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
291
292 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
293 {
294 if (internalDirichletConstraintsOwnElement[eqIdx])
295 {
296 origResiduals[0][eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
297 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
298 }
299 else
300 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
301 }
302
303 // off-diagonal entries
304 j = 1;
305 for (const auto& dataJ : connectivityMap[globalI])
306 {
307 const auto& neighborElement = neighborElements[j-1];
308 const auto& neighborScv = fvGeometry.scv(dataJ.globalJ);
309 const auto internalDirichletConstraintsNeighbor = this->problem().hasInternalDirichletConstraint(neighborElement, neighborScv);
310
311 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
312 {
313 if (internalDirichletConstraintsNeighbor[eqIdx])
314 A[dataJ.globalJ][globalI][eqIdx][pvIdx] = 0.0;
315 else
316 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j][eqIdx];
317 }
318
319 ++j;
320 }
321 }
322 else // no internal Dirichlet constraints specified
323 {
324
2/2
✓ Branch 0 taken 88630900 times.
✓ Branch 1 taken 34717902 times.
123348802 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
325 {
326 // the diagonal entries
327
4/8
✓ Branch 1 taken 88630900 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 88630900 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 88630900 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 86439576 times.
✗ Branch 11 not taken.
352332276 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
328
329 // off-diagonal entries
330 88630900 j = 1;
331
4/4
✓ Branch 0 taken 386572681 times.
✓ Branch 1 taken 92580950 times.
✓ Branch 2 taken 386572681 times.
✓ Branch 3 taken 92580950 times.
1135569062 for (const auto& dataJ : connectivityMap[globalI])
332
4/8
✓ Branch 1 taken 390522731 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 390522731 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 390522731 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 384381357 times.
✗ Branch 11 not taken.
1555949550 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j++][eqIdx];
333 }
334 }
335 }
336
337 // restore original state of the flux vars cache in case of global caching.
338 // This has to be done in order to guarantee that everything is in an undeflected
339 // state before the assembly of another element is called. In the case of local caching
340 // this is obsolete because the elemFluxVarsCache used here goes out of scope after this.
341 // We only have to do this for the last primary variable, for all others the flux var cache
342 // is updated with the correct element volume variables before residual evaluations
343 if (enableGridFluxVarsCache)
344
2/4
✓ Branch 1 taken 579386 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 579386 times.
✗ Branch 5 not taken.
5862642 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
345
346 // return the original residual
347 34023895 return origResiduals[0];
348 }
349 };
350
351
352 /*!
353 * \ingroup Assembly
354 * \ingroup CCDiscretization
355 * \brief Cell-centered scheme local assembler using numeric differentiation and explicit time discretization
356 */
357 template<class TypeTag, class Assembler>
358 class CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>
359 : public CCLocalAssemblerBase<TypeTag, Assembler,
360 CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false>, false>
361 {
362 using ThisType = CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false>;
363 using ParentType = CCLocalAssemblerBase<TypeTag, Assembler, ThisType, false>;
364 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
365 using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>;
366 using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity;
367 using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
368 using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>;
369 using Problem = typename GridVariables::GridVolumeVariables::Problem;
370
371 enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() };
372
373 public:
374 using ParentType::ParentType;
375
376 /*!
377 * \brief Computes the derivatives with respect to the given element and adds them
378 * to the global matrix.
379 *
380 * \return The element residual at the current solution.
381 */
382 NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables)
383 {
384 if (this->assembler().isStationaryProblem())
385 DUNE_THROW(Dune::InvalidStateException, "Using explicit jacobian assembler with stationary local residual");
386
387 // assemble the undeflected residual
388 auto residual = this->evalLocalResidual()[0];
389 const auto storageResidual = this->evalLocalStorageResidual()[0];
390
391 //////////////////////////////////////////////////////////////////////////////////////////////////
392 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
393 // neighboring elements all derivatives are zero. For the assembled element only the storage //
394 // derivatives are non-zero. //
395 //////////////////////////////////////////////////////////////////////////////////////////////////
396
397 // get some aliases for convenience
398 const auto& element = this->element();
399 const auto& fvGeometry = this->fvGeometry();
400 const auto& gridGeometry = this->assembler().gridGeometry();
401 auto&& curElemVolVars = this->curElemVolVars();
402
403 // reference to the element's scv (needed later) and corresponding vol vars
404 const auto globalI = gridGeometry.elementMapper().index(element);
405 const auto& scv = fvGeometry.scv(globalI);
406 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
407
408 // save a copy of the original privars and vol vars in order
409 // to restore the original solution after deflection
410 const auto& curSol = this->curSol();
411 const auto origPriVars = curSol[globalI];
412 const auto origVolVars = curVolVars;
413
414 // element solution container to be deflected
415 auto elemSol = elementSolution(element, curSol, gridGeometry);
416
417 NumEqVector partialDeriv;
418
419 // derivatives in the neighbors with respect to the current elements
420 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
421 {
422 // reset derivatives of element dof with respect to itself
423 partialDeriv = 0.0;
424
425 auto evalStorage = [&](Scalar priVar)
426 {
427 // update the volume variables and calculate
428 // the residual with the deflected primary variables
429 elemSol[0][pvIdx] = priVar;
430 curVolVars.update(elemSol, this->problem(), element, scv);
431 return this->evalLocalStorageResidual()[0];
432 };
433
434 // for non-ghosts compute the derivative numerically
435 if (!this->elementIsGhost())
436 {
437 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
438 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
439 NumericDifferentiation::partialDerivative(evalStorage, elemSol[0][pvIdx], partialDeriv, storageResidual,
440 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
441 }
442
443 // for ghost elements we assemble a 1.0 where the primary variable and zero everywhere else
444 // as we always solve for a delta of the solution with respect to the initial solution this
445 // results in a delta of zero for ghosts
446 else partialDeriv[pvIdx] = 1.0;
447
448 // restore the original state of the scv's volume variables
449 curVolVars = origVolVars;
450
451 // restore the current element solution
452 elemSol[0][pvIdx] = origPriVars[pvIdx];
453
454 // add the current partial derivatives to the global jacobian matrix
455 // only diagonal entries for explicit jacobians
456 if constexpr (Problem::enableInternalDirichletConstraints())
457 {
458 // check if own element has internal Dirichlet constraint
459 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
460 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
461
462 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
463 {
464 if (internalDirichletConstraints[eqIdx])
465 {
466 residual[eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
467 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
468 }
469 else
470 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
471 }
472 }
473 else
474 {
475 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
476 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
477 }
478 }
479
480 // return the original residual
481 return residual;
482 }
483 };
484
485 /*!
486 * \ingroup Assembly
487 * \ingroup CCDiscretization
488 * \brief Cell-centered scheme local assembler using analytic (hand-coded) differentiation and implicit time discretization
489 */
490 template<class TypeTag, class Assembler>
491 9328248 class CCLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true>
492 : public CCLocalAssemblerBase<TypeTag, Assembler,
493 CCLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true>, true>
494 {
495 using ThisType = CCLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true>;
496 using ParentType = CCLocalAssemblerBase<TypeTag, Assembler, ThisType, true>;
497 using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>;
498 using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>;
499 using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
500 using Problem = typename GridVariables::GridVolumeVariables::Problem;
501
502 enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() };
503
504 public:
505
2/8
✓ Branch 2 taken 8793950 times.
✗ Branch 3 not taken.
✓ Branch 6 taken 23328 times.
✗ Branch 7 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
8817278 using ParentType::ParentType;
506
507 /*!
508 * \brief Computes the derivatives with respect to the given element and adds them
509 * to the global matrix.
510 *
511 * \return The element residual at the current solution.
512 */
513 8315950 NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix& A, const GridVariables& gridVariables)
514 {
515 // treat ghost separately, we always want zero update for ghosts
516
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 8315950 times.
8315950 if (this->elementIsGhost())
517 {
518 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
519 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
520 A[globalI][globalI][pvIdx][pvIdx] = 1.0;
521
522 // return zero residual
523 return NumEqVector(0.0);
524 }
525
526 // assemble the undeflected residual
527 8315950 auto residual = this->evalLocalResidual()[0];
528
529 // get some aliases for convenience
530 8315950 const auto& problem = this->problem();
531 8315950 const auto& element = this->element();
532 8315950 const auto& fvGeometry = this->fvGeometry();
533 8315950 const auto& curElemVolVars = this->curElemVolVars();
534 8315950 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
535
536 // get reference to the element's current vol vars
537 24947850 const auto globalI = this->assembler().gridGeometry().elementMapper().index(element);
538
1/2
✓ Branch 0 taken 8304980 times.
✗ Branch 1 not taken.
8315950 const auto& scv = fvGeometry.scv(globalI);
539 8315950 const auto& volVars = curElemVolVars[scv];
540
541 // if the problem is instationary, add derivative of storage term
542
2/2
✓ Branch 0 taken 8261916 times.
✓ Branch 1 taken 54034 times.
8315950 if (!this->assembler().isStationaryProblem())
543 24785748 this->localResidual().addStorageDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv);
544
545 // add source term derivatives
546 24947850 this->localResidual().addSourceDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv);
547
548 // add flux derivatives for each scvf
549
6/6
✓ Branch 0 taken 16934916 times.
✓ Branch 1 taken 8315950 times.
✓ Branch 2 taken 16934916 times.
✓ Branch 3 taken 8315950 times.
✓ Branch 4 taken 36088 times.
✓ Branch 5 taken 852 times.
33566816 for (const auto& scvf : scvfs(fvGeometry))
550 {
551 // inner faces
552
2/2
✓ Branch 0 taken 16882338 times.
✓ Branch 1 taken 52578 times.
16934916 if (!scvf.boundary())
553 33764676 this->localResidual().addFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
554
555 // boundary faces
556 else
557 {
558
2/2
✓ Branch 0 taken 4806 times.
✓ Branch 1 taken 758 times.
52578 const auto& bcTypes = problem.boundaryTypes(element, scvf);
559
560 // add Dirichlet boundary flux derivatives
561
6/8
✓ Branch 0 taken 6908 times.
✓ Branch 1 taken 45670 times.
✓ Branch 2 taken 6908 times.
✓ Branch 3 taken 45670 times.
✓ Branch 4 taken 6908 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 6908 times.
✗ Branch 7 not taken.
105156 if (bcTypes.hasDirichlet() && !bcTypes.hasNeumann())
562 13816 this->localResidual().addCCDirichletFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
563
564 // add Robin ("solution dependent Neumann") boundary flux derivatives
565
4/8
✓ Branch 0 taken 45670 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 45670 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 45670 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 45670 times.
91340 else if (bcTypes.hasNeumann() && !bcTypes.hasDirichlet())
566 76188 this->localResidual().addRobinFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
567
568 else
569 DUNE_THROW(Dune::NotImplemented, "Mixed boundary conditions. Use pure boundary conditions by converting Dirichlet BCs to Robin BCs");
570 }
571 }
572
573 if constexpr (Problem::enableInternalDirichletConstraints())
574 {
575 // check if own element has internal Dirichlet constraint
576 200 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
577 200 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
578
579
2/2
✓ Branch 0 taken 100 times.
✓ Branch 1 taken 100 times.
200 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
580 {
581
2/2
✓ Branch 0 taken 100 times.
✓ Branch 1 taken 100 times.
200 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
582 {
583
4/4
✓ Branch 0 taken 99 times.
✓ Branch 1 taken 1 times.
✓ Branch 2 taken 99 times.
✓ Branch 3 taken 1 times.
200 if (internalDirichletConstraints[eqIdx])
584 {
585 1 residual[eqIdx] = volVars.priVars()[eqIdx] - dirichletValues[eqIdx];
586 2 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
587
588 // inner faces
589
4/4
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
6 for (const auto& scvf : scvfs(fvGeometry))
590
2/2
✓ Branch 0 taken 3 times.
✓ Branch 1 taken 1 times.
4 if (!scvf.boundary())
591
3/6
✗ Branch 0 not taken.
✓ Branch 1 taken 3 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 3 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 3 times.
15 A[globalI][fvGeometry.scv(scvf.outsideScvIdx()).dofIndex()][eqIdx][pvIdx] = 0.0;
592 }
593 }
594 }
595 }
596
597 // return element residual
598 8315950 return residual;
599 }
600 };
601
602 /*!
603 * \ingroup Assembly
604 * \ingroup CCDiscretization
605 * \brief Cell-centered scheme local assembler using analytic (hand-coded) differentiation and explicit time discretization
606 */
607 template<class TypeTag, class Assembler>
608 3586800 class CCLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/false>
609 : public CCLocalAssemblerBase<TypeTag, Assembler,
610 CCLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false>, false>
611 {
612 using ThisType = CCLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false>;
613 using ParentType = CCLocalAssemblerBase<TypeTag, Assembler, ThisType, false>;
614 using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>;
615 using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>;
616 using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
617 using Problem = typename GridVariables::GridVolumeVariables::Problem;
618
619 enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() };
620
621 public:
622
2/4
✓ Branch 2 taken 1826800 times.
✗ Branch 3 not taken.
✓ Branch 6 taken 5000 times.
✗ Branch 7 not taken.
1831800 using ParentType::ParentType;
623
624 /*!
625 * \brief Computes the derivatives with respect to the given element and adds them
626 * to the global matrix.
627 *
628 * \return The element residual at the current solution.
629 */
630 1331800 NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix& A, const GridVariables& gridVariables)
631 {
632 // treat ghost separately, we always want zero update for ghosts
633
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 1331800 times.
1331800 if (this->elementIsGhost())
634 {
635 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
636 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
637 A[globalI][globalI][pvIdx][pvIdx] = 1.0;
638
639 // return zero residual
640 return NumEqVector(0.0);
641 }
642
643 // assemble the undeflected residual
644 1331800 const auto residual = this->evalLocalResidual()[0];
645
646 // get reference to the element's current vol vars
647 3995400 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
648
2/4
✓ Branch 0 taken 76800 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 76800 times.
✗ Branch 3 not taken.
2663600 const auto& scv = this->fvGeometry().scv(globalI);
649 2663600 const auto& volVars = this->curElemVolVars()[scv];
650
651 // add hand-code derivative of storage term
652 6659000 this->localResidual().addStorageDerivatives(A[globalI][globalI], this->problem(), this->element(), this->fvGeometry(), volVars, scv);
653
654 if constexpr (Problem::enableInternalDirichletConstraints())
655 {
656 // check if own element has internal Dirichlet constraint
657 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
658 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
659
660 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
661 {
662 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
663 {
664 if (internalDirichletConstraints[eqIdx])
665 {
666 residual[eqIdx] = volVars.priVars()[eqIdx] - dirichletValues[eqIdx];
667 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
668 }
669 }
670 }
671 }
672
673 // return the original residual
674 1331800 return residual;
675 }
676 };
677
678 } // end namespace Dumux
679
680 #endif
681