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 Assembly | ||
10 | * \ingroup CVFEDiscretization | ||
11 | * \brief An assembler for Jacobian and residual contribution per element (CVFE methods) | ||
12 | */ | ||
13 | #ifndef DUMUX_CVFE_LOCAL_ASSEMBLER_HH | ||
14 | #define DUMUX_CVFE_LOCAL_ASSEMBLER_HH | ||
15 | |||
16 | #include <dune/common/exceptions.hh> | ||
17 | #include <dune/common/hybridutilities.hh> | ||
18 | #include <dune/common/reservedvector.hh> | ||
19 | #include <dune/grid/common/gridenums.hh> | ||
20 | #include <dune/istl/matrixindexset.hh> | ||
21 | #include <dune/istl/bvector.hh> | ||
22 | |||
23 | #include <dumux/common/properties.hh> | ||
24 | #include <dumux/common/parameters.hh> | ||
25 | #include <dumux/common/numericdifferentiation.hh> | ||
26 | |||
27 | #include <dumux/assembly/numericepsilon.hh> | ||
28 | #include <dumux/assembly/diffmethod.hh> | ||
29 | #include <dumux/assembly/fvlocalassemblerbase.hh> | ||
30 | #include <dumux/assembly/partialreassembler.hh> | ||
31 | #include <dumux/assembly/entitycolor.hh> | ||
32 | |||
33 | #include <dumux/discretization/cvfe/elementsolution.hh> | ||
34 | |||
35 | #include "volvardeflectionhelper_.hh" | ||
36 | |||
37 | namespace Dumux { | ||
38 | |||
39 | #ifndef DOXYGEN | ||
40 | namespace Detail::CVFE { | ||
41 | |||
42 | struct NoOperator | ||
43 | { | ||
44 | template<class... Args> | ||
45 | constexpr void operator()(Args&&...) const {} | ||
46 | }; | ||
47 | |||
48 | template<class X, class Y> | ||
49 | using Impl = std::conditional_t<!std::is_same_v<X, void>, X, Y>; | ||
50 | |||
51 | } // end namespace Detail | ||
52 | #endif // DOXYGEN | ||
53 | |||
54 | /*! | ||
55 | * \ingroup Assembly | ||
56 | * \ingroup CVFEDiscretization | ||
57 | * \brief A base class for all local CVFE assemblers | ||
58 | * \tparam TypeTag The TypeTag | ||
59 | * \tparam Assembler The assembler type | ||
60 | * \tparam Implementation The actual implementation | ||
61 | * \tparam implicit Specifies whether the time discretization is implicit or not not (i.e. explicit) | ||
62 | */ | ||
63 | template<class TypeTag, class Assembler, class Implementation, bool implicit> | ||
64 | 20200727 | class CVFELocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit> | |
65 | { | ||
66 | using ParentType = FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit>; | ||
67 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
68 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
69 | using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView; | ||
70 | using SolutionVector = typename Assembler::SolutionVector; | ||
71 | |||
72 | static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq(); | ||
73 | static constexpr int dim = GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension; | ||
74 | |||
75 | public: | ||
76 | |||
77 | 16662088 | using ParentType::ParentType; | |
78 | |||
79 | 13569764 | void bindLocalViews() | |
80 | { | ||
81 | 13569764 | ParentType::bindLocalViews(); | |
82 | 13569764 | this->elemBcTypes().update(this->asImp_().problem(), this->element(), this->fvGeometry()); | |
83 | 13569764 | } | |
84 | |||
85 | |||
86 | /*! | ||
87 | * \brief Computes the derivatives with respect to the given element and adds them | ||
88 | * to the global matrix. The element residual is written into the right hand side. | ||
89 | */ | ||
90 | template <class ResidualVector, class PartialReassembler = DefaultPartialReassembler, class CouplingFunction = Detail::CVFE::NoOperator> | ||
91 | 16632986 | void assembleJacobianAndResidual(JacobianMatrix& jac, ResidualVector& res, GridVariables& gridVariables, | |
92 | const PartialReassembler* partialReassembler = nullptr, | ||
93 | const CouplingFunction& maybeAssembleCouplingBlocks = {}) | ||
94 | { | ||
95 | 16632986 | this->asImp_().bindLocalViews(); | |
96 | 16632986 | const auto eIdxGlobal = this->asImp_().problem().gridGeometry().elementMapper().index(this->element()); | |
97 | if (partialReassembler | ||
98 |
4/4✓ Branch 0 taken 452312 times.
✓ Branch 1 taken 10605623 times.
✓ Branch 2 taken 101622 times.
✓ Branch 3 taken 350690 times.
|
12150199 | && partialReassembler->elementColor(eIdxGlobal) == EntityColor::green) |
99 | { | ||
100 | 101622 | const auto residual = this->asImp_().evalLocalResidual(); // forward to the internal implementation | |
101 |
2/2✓ Branch 0 taken 406488 times.
✓ Branch 1 taken 101622 times.
|
508110 | for (const auto& scv : scvs(this->fvGeometry())) |
102 | 812976 | res[scv.dofIndex()] += residual[scv.localDofIndex()]; | |
103 | |||
104 | // assemble the coupling blocks for coupled models (does nothing if not coupled) | ||
105 | maybeAssembleCouplingBlocks(residual); | ||
106 | } | ||
107 |
2/2✓ Branch 0 taken 13491491 times.
✓ Branch 1 taken 69504 times.
|
16531364 | else if (!this->elementIsGhost()) |
108 | { | ||
109 | 16461860 | const auto residual = this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables, partialReassembler); // forward to the internal implementation | |
110 |
2/2✓ Branch 0 taken 45463904 times.
✓ Branch 1 taken 13491491 times.
|
70409936 | for (const auto& scv : scvs(this->fvGeometry())) |
111 | 91603184 | res[scv.dofIndex()] += residual[scv.localDofIndex()]; | |
112 | |||
113 | // assemble the coupling blocks for coupled models (does nothing if not coupled) | ||
114 | 4327926 | maybeAssembleCouplingBlocks(residual); | |
115 | } | ||
116 | else | ||
117 | { | ||
118 | // Treatment of ghost elements | ||
119 | 69504 | assert(this->elementIsGhost()); | |
120 | |||
121 | // handle dofs per codimension | ||
122 | 69504 | const auto& gridGeometry = this->asImp_().problem().gridGeometry(); | |
123 | 419904 | Dune::Hybrid::forEach(std::make_integer_sequence<int, dim>{}, [&](auto d) | |
124 | { | ||
125 | 139008 | constexpr int codim = dim - d; | |
126 | 139008 | const auto& localCoeffs = gridGeometry.feCache().get(this->element().type()).localCoefficients(); | |
127 |
6/14✓ Branch 0 taken 272592 times.
✓ Branch 1 taken 72696 times.
✓ Branch 2 taken 273948 times.
✓ Branch 3 taken 68148 times.
✓ Branch 4 taken 4548 times.
✓ Branch 5 taken 1356 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 6 not taken.
✗ Branch 9 not taken.
|
828060 | for (int idx = 0; idx < localCoeffs.size(); ++idx) |
128 | { | ||
129 |
4/14✓ Branch 0 taken 272592 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 4548 times.
✓ Branch 3 taken 272592 times.
✓ Branch 4 taken 4548 times.
✗ Branch 5 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 6 not taken.
✗ Branch 9 not taken.
|
554280 | const auto& localKey = localCoeffs.localKey(idx); |
130 | |||
131 | // skip if we are not handling this codim right now | ||
132 |
4/12✓ Branch 0 taken 272592 times.
✓ Branch 1 taken 4548 times.
✓ Branch 2 taken 4548 times.
✓ Branch 3 taken 272592 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
554280 | if (localKey.codim() != codim) |
133 | 414792 | continue; | |
134 | |||
135 | // do not change the non-ghost entities | ||
136 |
4/18✓ Branch 1 taken 1356 times.
✓ Branch 2 taken 3192 times.
✓ Branch 4 taken 133416 times.
✓ Branch 5 taken 133416 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 0 not taken.
✗ Branch 3 not taken.
✗ Branch 6 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✗ Branch 9 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
277140 | auto entity = this->element().template subEntity<codim>(localKey.subEntity()); |
137 |
7/32✓ Branch 0 taken 3192 times.
✓ Branch 1 taken 1356 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 133416 times.
✓ Branch 4 taken 133416 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 5760 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2880 times.
✓ Branch 9 taken 2880 times.
✗ 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 26 not taken.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
|
282900 | if (entity.partitionType() == Dune::InteriorEntity || entity.partitionType() == Dune::BorderEntity) |
138 | 2880 | continue; | |
139 | |||
140 | // WARNING: this only works if the mapping from codim+subEntity to | ||
141 | // global dofIndex is unique (on dof per entity of this codim). | ||
142 | // For more general mappings, we should use a proper local-global mapping here. | ||
143 | // For example through dune-functions. | ||
144 |
0/11✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 6 not taken.
✗ Branch 3 not taken.
✗ Branch 7 not taken.
|
139488 | const auto dofIndex = gridGeometry.dofMapper().index(entity); |
145 | |||
146 | // this might be a vector-valued dof | ||
147 | using BlockType = typename JacobianMatrix::block_type; | ||
148 |
2/11✓ Branch 1 taken 3192 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 133416 times.
✗ Branch 5 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 9 not taken.
✗ Branch 12 not taken.
✗ Branch 6 not taken.
|
139488 | BlockType &J = jac[dofIndex][dofIndex]; |
149 |
4/12✓ Branch 0 taken 6384 times.
✓ Branch 1 taken 3192 times.
✓ Branch 2 taken 332712 times.
✓ Branch 3 taken 136296 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
478584 | for (int j = 0; j < BlockType::rows; ++j) |
150 | 339096 | J[j][j] = 1.0; | |
151 | |||
152 | // set residual for the ghost dof | ||
153 |
0/8✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
271776 | res[dofIndex] = 0; |
154 | } | ||
155 | }); | ||
156 | } | ||
157 | |||
158 | 19988708 | auto applyDirichlet = [&] (const auto& scvI, | |
159 | const auto& dirichletValues, | ||
160 | const auto eqIdx, | ||
161 | const auto pvIdx) | ||
162 | { | ||
163 | 1677861 | res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx]; | |
164 | |||
165 | 1677861 | auto& row = jac[scvI.dofIndex()]; | |
166 |
9/12✗ Branch 0 not taken.
✓ Branch 1 taken 12473180 times.
✓ Branch 2 taken 10958589 times.
✓ Branch 3 taken 1514591 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 1497657 times.
✓ Branch 6 taken 1335331 times.
✓ Branch 7 taken 162326 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 5656 times.
✓ Branch 10 taken 4712 times.
✓ Branch 11 taken 944 times.
|
13976493 | for (auto col = row.begin(); col != row.end(); ++col) |
167 | 23189509 | row[col.index()][eqIdx] = 0.0; | |
168 | |||
169 | 1677861 | jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0; | |
170 | |||
171 | // if a periodic dof has Dirichlet values also apply the same Dirichlet values to the other dof | ||
172 | 2529417 | if (this->asImp_().problem().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex())) | |
173 | { | ||
174 | 4 | const auto periodicDof = this->asImp_().problem().gridGeometry().periodicallyMappedDof(scvI.dofIndex()); | |
175 | 4 | res[periodicDof][eqIdx] = this->asImp_().curSol()[periodicDof][pvIdx] - dirichletValues[pvIdx]; | |
176 | |||
177 | 4 | auto& rowP = jac[periodicDof]; | |
178 |
3/12✗ Branch 0 not taken.
✓ Branch 1 taken 30 times.
✓ Branch 2 taken 26 times.
✓ Branch 3 taken 4 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
30 | for (auto col = rowP.begin(); col != rowP.end(); ++col) |
179 | 26 | rowP[col.index()][eqIdx] = 0.0; | |
180 | |||
181 | 4 | rowP[periodicDof][eqIdx][pvIdx] = 1.0; | |
182 | } | ||
183 | }; | ||
184 | |||
185 | 33265972 | this->asImp_().enforceDirichletConstraints(applyDirichlet); | |
186 | 16632986 | } | |
187 | |||
188 | /*! | ||
189 | * \brief Computes the derivatives with respect to the given element and adds them | ||
190 | * to the global matrix. | ||
191 | */ | ||
192 | 1136 | void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables) | |
193 | { | ||
194 | 1136 | this->asImp_().bindLocalViews(); | |
195 | 1136 | this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation | |
196 | |||
197 | 1816 | auto applyDirichlet = [&] (const auto& scvI, | |
198 | const auto& dirichletValues, | ||
199 | const auto eqIdx, | ||
200 | const auto pvIdx) | ||
201 | { | ||
202 | 336 | auto& row = jac[scvI.dofIndex()]; | |
203 |
3/8✗ Branch 0 not taken.
✓ Branch 1 taken 2416 times.
✓ Branch 2 taken 2080 times.
✓ Branch 3 taken 336 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
2416 | for (auto col = row.begin(); col != row.end(); ++col) |
204 | 2080 | row[col.index()][eqIdx] = 0.0; | |
205 | |||
206 | 336 | jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0; | |
207 | }; | ||
208 | |||
209 | 2272 | this->asImp_().enforceDirichletConstraints(applyDirichlet); | |
210 | 1136 | } | |
211 | |||
212 | /*! | ||
213 | * \brief Assemble the residual only | ||
214 | */ | ||
215 | template <class ResidualVector> | ||
216 | 2820927 | void assembleResidual(ResidualVector& res) | |
217 | { | ||
218 | 2820927 | this->asImp_().bindLocalViews(); | |
219 | 2820927 | const auto residual = this->evalLocalResidual(); | |
220 | |||
221 |
2/2✓ Branch 0 taken 8924562 times.
✓ Branch 1 taken 2298815 times.
|
13747617 | for (const auto& scv : scvs(this->fvGeometry())) |
222 | 18773356 | res[scv.dofIndex()] += residual[scv.localDofIndex()]; | |
223 | |||
224 | 2962729 | auto applyDirichlet = [&] (const auto& scvI, | |
225 | const auto& dirichletValues, | ||
226 | const auto eqIdx, | ||
227 | const auto pvIdx) | ||
228 | { | ||
229 | 141802 | res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx]; | |
230 | }; | ||
231 | |||
232 | 5641854 | this->asImp_().enforceDirichletConstraints(applyDirichlet); | |
233 | 2820927 | } | |
234 | |||
235 | //! Enforce Dirichlet constraints | ||
236 | template<typename ApplyFunction> | ||
237 | 15962568 | void enforceDirichletConstraints(const ApplyFunction& applyDirichlet) | |
238 | { | ||
239 | // enforce Dirichlet boundary conditions | ||
240 | 15962568 | this->asImp_().evalDirichletBoundaries(applyDirichlet); | |
241 | // take care of internal Dirichlet constraints (if enabled) | ||
242 | 15956168 | this->asImp_().enforceInternalDirichletConstraints(applyDirichlet); | |
243 | } | ||
244 | |||
245 | /*! | ||
246 | * \brief Evaluates Dirichlet boundaries | ||
247 | */ | ||
248 | template< typename ApplyDirichletFunctionType > | ||
249 | 31907652 | void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet) | |
250 | { | ||
251 | // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0 | ||
252 | // and set the residual to (privar - dirichletvalue) | ||
253 |
2/2✓ Branch 0 taken 542862 times.
✓ Branch 1 taken 15419706 times.
|
31907652 | if (this->elemBcTypes().hasDirichlet()) |
254 | { | ||
255 |
2/2✓ Branch 0 taken 1874911 times.
✓ Branch 1 taken 542862 times.
|
4821760 | for (const auto& scvI : scvs(this->fvGeometry())) |
256 | { | ||
257 |
2/2✓ Branch 1 taken 896588 times.
✓ Branch 2 taken 978323 times.
|
3737762 | const auto bcTypes = this->elemBcTypes().get(this->fvGeometry(), scvI); |
258 |
2/2✓ Branch 0 taken 896588 times.
✓ Branch 1 taken 978323 times.
|
3737762 | if (bcTypes.hasDirichlet()) |
259 | { | ||
260 |
2/3✓ Branch 0 taken 180773 times.
✓ Branch 1 taken 285285 times.
✗ Branch 2 not taken.
|
1787248 | const auto dirichletValues = this->asImp_().problem().dirichlet(this->element(), scvI); |
261 | |||
262 | // set the Dirichlet conditions in residual and jacobian | ||
263 |
3/3✓ Branch 0 taken 1792169 times.
✓ Branch 1 taken 931436 times.
✓ Branch 2 taken 17424 times.
|
5470202 | for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) |
264 | { | ||
265 |
3/4✓ Branch 0 taken 24454 times.
✓ Branch 1 taken 1819987 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 24454 times.
|
3682954 | if (bcTypes.isDirichlet(eqIdx)) |
266 | { | ||
267 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1819987 times.
|
3634046 | const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx); |
268 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1819987 times.
|
3634046 | assert(0 <= pvIdx && pvIdx < numEq); |
269 |
1/2✓ Branch 1 taken 311512 times.
✗ Branch 2 not taken.
|
3634046 | applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx); |
270 | } | ||
271 | } | ||
272 | } | ||
273 | } | ||
274 | } | ||
275 | 31907652 | } | |
276 | |||
277 | /*! | ||
278 | * \brief Update the coupling context for coupled models. | ||
279 | * \note This does nothing per default (not a coupled model). | ||
280 | */ | ||
281 | template<class... Args> | ||
282 | void maybeUpdateCouplingContext(Args&&...) {} | ||
283 | |||
284 | /*! | ||
285 | * \brief Update the additional domain derivatives for coupled models. | ||
286 | * \note This does nothing per default (not a coupled model). | ||
287 | */ | ||
288 | template<class... Args> | ||
289 | void maybeEvalAdditionalDomainDerivatives(Args&&...) {} | ||
290 | }; | ||
291 | |||
292 | /*! | ||
293 | * \ingroup Assembly | ||
294 | * \ingroup CVFEDiscretization | ||
295 | * \brief An assembler for Jacobian and residual contribution per element (CVFE methods) | ||
296 | * \tparam TypeTag The TypeTag | ||
297 | * \tparam diffMethod The differentiation method to residual compute derivatives | ||
298 | * \tparam implicit Specifies whether the time discretization is implicit or not not (i.e. explicit) | ||
299 | * \tparam Implementation via CRTP | ||
300 | */ | ||
301 | template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, bool implicit = true, class Implementation = void> | ||
302 | class CVFELocalAssembler; | ||
303 | |||
304 | /*! | ||
305 | * \ingroup Assembly | ||
306 | * \ingroup CVFEDiscretization | ||
307 | * \brief Control volume finite element local assembler using numeric differentiation and implicit time discretization | ||
308 | */ | ||
309 | template<class TypeTag, class Assembler, class Implementation> | ||
310 |
4/16✓ Branch 0 taken 2202672 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 503303 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✓ Branch 8 taken 261128 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✓ Branch 12 taken 261056 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
16917130 | class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true, Implementation> |
311 | : public CVFELocalAssemblerBase<TypeTag, Assembler, | ||
312 | Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true, Implementation>>, | ||
313 | true> | ||
314 | { | ||
315 | using ThisType = CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true, Implementation>; | ||
316 | using ParentType = CVFELocalAssemblerBase<TypeTag, Assembler, Detail::CVFE::Impl<Implementation, ThisType>, true>; | ||
317 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
318 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
319 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
320 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
321 | |||
322 | static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq(); | ||
323 | static constexpr int dim = GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension; | ||
324 | |||
325 | static constexpr bool enableGridFluxVarsCache | ||
326 | = GridVariables::GridFluxVariablesCache::cachingEnabled; | ||
327 | static constexpr bool solutionDependentFluxVarsCache | ||
328 | = GridVariables::GridFluxVariablesCache::FluxVariablesCache::isSolDependent; | ||
329 | |||
330 | public: | ||
331 | |||
332 | using LocalResidual = typename ParentType::LocalResidual; | ||
333 | using ElementResidualVector = typename LocalResidual::ElementResidualVector; | ||
334 |
6/12✓ Branch 2 taken 9725537 times.
✗ Branch 3 not taken.
✓ Branch 6 taken 566220 times.
✗ Branch 7 not taken.
✓ Branch 10 taken 72 times.
✗ Branch 11 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
✓ Branch 1 taken 931296 times.
✓ Branch 5 taken 268160 times.
✓ Branch 9 taken 7708 times.
✗ Branch 4 not taken.
|
13682112 | using ParentType::ParentType; |
335 | |||
336 | /*! | ||
337 | * \brief Computes the derivatives with respect to the given element and adds them | ||
338 | * to the global matrix. | ||
339 | * | ||
340 | * \return The element residual at the current solution. | ||
341 | */ | ||
342 | template <class PartialReassembler = DefaultPartialReassembler> | ||
343 | 12227870 | ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables, | |
344 | const PartialReassembler* partialReassembler = nullptr) | ||
345 | { | ||
346 | // get some aliases for convenience | ||
347 | 12227870 | const auto& element = this->element(); | |
348 | 12227870 | const auto& fvGeometry = this->fvGeometry(); | |
349 | 12227870 | const auto& curSol = this->asImp_().curSol(); | |
350 | |||
351 | 12227870 | auto&& curElemVolVars = this->curElemVolVars(); | |
352 | 12227870 | auto&& elemFluxVarsCache = this->elemFluxVarsCache(); | |
353 | |||
354 | // get the vector of the actual element residuals | ||
355 | 12227870 | const auto origResiduals = this->evalLocalResidual(); | |
356 | |||
357 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
358 | // // | ||
359 | // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the // | ||
360 | // neighboring elements we do so by computing the derivatives of the fluxes which depend on the // | ||
361 | // actual element. In the actual element we evaluate the derivative of the entire residual. // | ||
362 | // // | ||
363 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
364 | |||
365 | // if all volvars in the stencil have to be updated or if it's enough to only update the | ||
366 | // volVars for the scv whose associated dof has been deflected | ||
367 |
5/6✓ Branch 0 taken 172 times.
✓ Branch 1 taken 10271329 times.
✓ Branch 3 taken 163 times.
✓ Branch 4 taken 9 times.
✓ Branch 6 taken 163 times.
✗ Branch 7 not taken.
|
12227870 | static const bool updateAllVolVars = getParamFromGroup<bool>( |
368 |
1/2✓ Branch 1 taken 163 times.
✗ Branch 2 not taken.
|
452 | this->asImp_().problem().paramGroup(), "Assembly.BoxVolVarsDependOnAllElementDofs", false |
369 | ); | ||
370 | |||
371 | // create the element solution | ||
372 | 13942736 | auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry()); | |
373 | |||
374 | // create the vector storing the partial derivatives | ||
375 | 12227870 | ElementResidualVector partialDerivs(fvGeometry.numScv()); | |
376 | |||
377 | 12227870 | Detail::VolVarsDeflectionHelper deflectionHelper( | |
378 | 28658638 | [&] (const auto& scv) -> VolumeVariables& { | |
379 | 234852218 | return this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv); | |
380 | }, | ||
381 | fvGeometry, | ||
382 | updateAllVolVars | ||
383 | ); | ||
384 | |||
385 | // calculation of the derivatives | ||
386 |
4/4✓ Branch 0 taken 22054684 times.
✓ Branch 1 taken 14188362 times.
✓ Branch 2 taken 26177315 times.
✓ Branch 3 taken 6148870 times.
|
53963088 | for (const auto& scv : scvs(fvGeometry)) |
387 | { | ||
388 | // dof index and corresponding actual pri vars | ||
389 |
1/2✓ Branch 0 taken 22054684 times.
✗ Branch 1 not taken.
|
41735218 | const auto dofIdx = scv.dofIndex(); |
390 | 41735218 | deflectionHelper.setCurrent(scv); | |
391 | |||
392 | // calculate derivatives w.r.t to the privars at the dof at hand | ||
393 |
2/2✓ Branch 0 taken 84562816 times.
✓ Branch 1 taken 36243046 times.
|
134019494 | for (int pvIdx = 0; pvIdx < numEq; pvIdx++) |
394 | { | ||
395 | 92284276 | partialDerivs = 0.0; | |
396 | |||
397 | 185890536 | auto evalResiduals = [&](Scalar priVar) | |
398 | { | ||
399 | // update the volume variables and compute element residual | ||
400 | 93606260 | elemSol[scv.localDofIndex()][pvIdx] = priVar; | |
401 | 93606260 | deflectionHelper.deflect(elemSol, scv, this->asImp_().problem()); | |
402 | if constexpr (solutionDependentFluxVarsCache) | ||
403 | { | ||
404 | 3127718 | elemFluxVarsCache.update(element, fvGeometry, curElemVolVars); | |
405 | if constexpr (enableGridFluxVarsCache) | ||
406 | gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars); | ||
407 | } | ||
408 | 93606260 | this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx); | |
409 | 93606260 | return this->evalLocalResidual(); | |
410 | }; | ||
411 | |||
412 | // derive the residuals numerically | ||
413 |
5/6✓ Branch 0 taken 209 times.
✓ Branch 1 taken 84562607 times.
✓ Branch 3 taken 163 times.
✓ Branch 4 taken 46 times.
✓ Branch 6 taken 163 times.
✗ Branch 7 not taken.
|
92284276 | static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()}; |
414 |
5/6✓ Branch 0 taken 235 times.
✓ Branch 1 taken 84562581 times.
✓ Branch 3 taken 163 times.
✓ Branch 4 taken 72 times.
✓ Branch 6 taken 163 times.
✗ Branch 7 not taken.
|
92284276 | static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod"); |
415 |
1/2✓ Branch 2 taken 84562816 times.
✗ Branch 3 not taken.
|
92284276 | NumericDifferentiation::partialDerivative(evalResiduals, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals, |
416 | 92284276 | eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod); | |
417 | |||
418 | // update the global stiffness matrix with the current partial derivatives | ||
419 |
2/2✓ Branch 0 taken 312887950 times.
✓ Branch 1 taken 84562816 times.
|
431994646 | for (const auto& scvJ : scvs(fvGeometry)) |
420 | { | ||
421 | // don't add derivatives for green dofs | ||
422 | if (!partialReassembler | ||
423 |
4/4✓ Branch 0 taken 11222080 times.
✓ Branch 1 taken 284227662 times.
✓ Branch 2 taken 10849608 times.
✓ Branch 3 taken 372472 times.
|
299818798 | || partialReassembler->dofColor(scvJ.dofIndex()) != EntityColor::green) |
424 | { | ||
425 |
2/2✓ Branch 0 taken 821068658 times.
✓ Branch 1 taken 312515478 times.
|
1208715008 | for (int eqIdx = 0; eqIdx < numEq; eqIdx++) |
426 | { | ||
427 | // A[i][col][eqIdx][pvIdx] is the rate of change of | ||
428 | // the residual of equation 'eqIdx' at dof 'i' | ||
429 | // depending on the primary variable 'pvIdx' at dof | ||
430 | // 'col'. | ||
431 |
1/2✓ Branch 1 taken 821068658 times.
✗ Branch 2 not taken.
|
869377110 | A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx]; |
432 | } | ||
433 | } | ||
434 | } | ||
435 | |||
436 | // restore the original state of the scv's volume variables | ||
437 | 92284276 | deflectionHelper.restore(scv); | |
438 | |||
439 | // restore the original element solution | ||
440 |
1/2✓ Branch 1 taken 3982 times.
✗ Branch 2 not taken.
|
92284276 | elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx]; |
441 |
1/2✓ Branch 1 taken 3982 times.
✗ Branch 2 not taken.
|
92341316 | this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx); |
442 | } | ||
443 | } | ||
444 | |||
445 | // restore original state of the flux vars cache in case of global caching. | ||
446 | // In the case of local caching this is obsolete because the elemFluxVarsCache used here goes out of scope after this. | ||
447 | if constexpr (enableGridFluxVarsCache) | ||
448 | 3850564 | gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars); | |
449 | |||
450 | // evaluate additional derivatives that might arise from the coupling (no-op if not coupled) | ||
451 | 10858843 | this->asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables); | |
452 | |||
453 | 12227870 | return origResiduals; | |
454 | } | ||
455 | |||
456 | }; // implicit CVFEAssembler with numeric Jacobian | ||
457 | |||
458 | /*! | ||
459 | * \ingroup Assembly | ||
460 | * \ingroup CVFEDiscretization | ||
461 | * \brief Control volume finite element local assembler using numeric differentiation and explicit time discretization | ||
462 | */ | ||
463 | template<class TypeTag, class Assembler, class Implementation> | ||
464 | 1014000 | class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false, Implementation> | |
465 | : public CVFELocalAssemblerBase<TypeTag, Assembler, | ||
466 | Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false, Implementation>>, | ||
467 | false> | ||
468 | { | ||
469 | using ThisType = CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false, Implementation>; | ||
470 | using ParentType = CVFELocalAssemblerBase<TypeTag, Assembler, Detail::CVFE::Impl<Implementation, ThisType>, false>; | ||
471 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
472 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
473 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
474 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
475 | |||
476 | static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq(); | ||
477 | static constexpr int dim = GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension; | ||
478 | |||
479 | public: | ||
480 | |||
481 | using LocalResidual = typename ParentType::LocalResidual; | ||
482 | using ElementResidualVector = typename LocalResidual::ElementResidualVector; | ||
483 |
2/4✓ Branch 1 taken 50000 times.
✗ Branch 2 not taken.
✓ Branch 5 taken 964000 times.
✗ Branch 6 not taken.
|
1014000 | using ParentType::ParentType; |
484 | |||
485 | /*! | ||
486 | * \brief Computes the derivatives with respect to the given element and adds them | ||
487 | * to the global matrix. | ||
488 | * | ||
489 | * \return The element residual at the current solution. | ||
490 | */ | ||
491 | template <class PartialReassembler = DefaultPartialReassembler> | ||
492 | 1014000 | ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables, | |
493 | const PartialReassembler* partialReassembler = nullptr) | ||
494 | { | ||
495 | 1014000 | if (partialReassembler) | |
496 | ✗ | DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization"); | |
497 | |||
498 | // get some aliases for convenience | ||
499 | 1014000 | const auto& element = this->element(); | |
500 | 1014000 | const auto& fvGeometry = this->fvGeometry(); | |
501 | 1014000 | const auto& curSol = this->asImp_().curSol(); | |
502 | 1014000 | auto&& curElemVolVars = this->curElemVolVars(); | |
503 | |||
504 | // get the vecor of the actual element residuals | ||
505 | 1014000 | const auto origResiduals = this->evalLocalResidual(); | |
506 | 1014000 | const auto origStorageResiduals = this->evalLocalStorageResidual(); | |
507 | |||
508 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
509 | // // | ||
510 | // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the // | ||
511 | // neighboring elements we do so by computing the derivatives of the fluxes which depend on the // | ||
512 | // actual element. In the actual element we evaluate the derivative of the entire residual. // | ||
513 | // // | ||
514 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
515 | |||
516 | // create the element solution | ||
517 | 2028000 | auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry()); | |
518 | |||
519 | // create the vector storing the partial derivatives | ||
520 | 2028000 | ElementResidualVector partialDerivs(fvGeometry.numScv()); | |
521 | |||
522 | // calculation of the derivatives | ||
523 | 4006000 | for (const auto& scv : scvs(fvGeometry)) | |
524 | { | ||
525 | // dof index and corresponding actual pri vars | ||
526 | 2992000 | const auto dofIdx = scv.dofIndex(); | |
527 | 2992000 | auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv); | |
528 | 2992000 | const VolumeVariables origVolVars(curVolVars); | |
529 | |||
530 | // calculate derivatives w.r.t to the privars at the dof at hand | ||
531 | 5984000 | for (int pvIdx = 0; pvIdx < numEq; pvIdx++) | |
532 | { | ||
533 | 2992000 | partialDerivs = 0.0; | |
534 | |||
535 | 5984000 | auto evalStorage = [&](Scalar priVar) | |
536 | { | ||
537 | // auto partialDerivsTmp = partialDerivs; | ||
538 | 2992000 | elemSol[scv.localDofIndex()][pvIdx] = priVar; | |
539 | 2992000 | curVolVars.update(elemSol, this->asImp_().problem(), element, scv); | |
540 | 2992000 | return this->evalLocalStorageResidual(); | |
541 | }; | ||
542 | |||
543 | // derive the residuals numerically | ||
544 | 2992000 | static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()}; | |
545 | 2992000 | static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod"); | |
546 | 2992000 | NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origStorageResiduals, | |
547 | 2992000 | eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod); | |
548 | |||
549 | // update the global stiffness matrix with the current partial derivatives | ||
550 | 5984000 | for (int eqIdx = 0; eqIdx < numEq; eqIdx++) | |
551 | { | ||
552 | // A[i][col][eqIdx][pvIdx] is the rate of change of | ||
553 | // the residual of equation 'eqIdx' at dof 'i' | ||
554 | // depending on the primary variable 'pvIdx' at dof | ||
555 | // 'col'. | ||
556 | 2992000 | A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx]; | |
557 | } | ||
558 | |||
559 | // restore the original state of the scv's volume variables | ||
560 | 2992000 | curVolVars = origVolVars; | |
561 | |||
562 | // restore the original element solution | ||
563 | 2992000 | elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx]; | |
564 | // TODO additional dof dependencies | ||
565 | } | ||
566 | } | ||
567 | |||
568 | 1014000 | return origResiduals; | |
569 | } | ||
570 | }; // explicit CVFEAssembler with numeric Jacobian | ||
571 | |||
572 | /*! | ||
573 | * \ingroup Assembly | ||
574 | * \ingroup CVFEDiscretization | ||
575 | * \brief Control volume finite element local assembler using analytic differentiation and implicit time discretization | ||
576 | */ | ||
577 | template<class TypeTag, class Assembler, class Implementation> | ||
578 |
2/8✓ Branch 0 taken 500000 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 5000 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1259597 | class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true, Implementation> |
579 | : public CVFELocalAssemblerBase<TypeTag, Assembler, | ||
580 | Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true, Implementation>>, | ||
581 | true> | ||
582 | { | ||
583 | using ThisType = CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true, Implementation>; | ||
584 | using ParentType = CVFELocalAssemblerBase<TypeTag, Assembler, Detail::CVFE::Impl<Implementation, ThisType>, true>; | ||
585 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
586 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
587 | |||
588 | public: | ||
589 | |||
590 | using LocalResidual = typename ParentType::LocalResidual; | ||
591 | using ElementResidualVector = typename LocalResidual::ElementResidualVector; | ||
592 |
2/4✓ Branch 2 taken 729133 times.
✗ Branch 3 not taken.
✓ Branch 6 taken 25464 times.
✗ Branch 7 not taken.
|
754597 | using ParentType::ParentType; |
593 | |||
594 | /*! | ||
595 | * \brief Computes the derivatives with respect to the given element and adds them | ||
596 | * to the global matrix. | ||
597 | * | ||
598 | * \return The element residual at the current solution. | ||
599 | */ | ||
600 | template <class PartialReassembler = DefaultPartialReassembler> | ||
601 | 253869 | ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables, | |
602 | const PartialReassembler* partialReassembler = nullptr) | ||
603 | { | ||
604 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 253869 times.
|
253869 | if (partialReassembler) |
605 | ✗ | DUNE_THROW(Dune::NotImplemented, "partial reassembly for analytic differentiation"); | |
606 | |||
607 | // get some aliases for convenience | ||
608 | 253869 | const auto& element = this->element(); | |
609 | 253869 | const auto& fvGeometry = this->fvGeometry(); | |
610 | 253869 | const auto& problem = this->asImp_().problem(); | |
611 | 253869 | const auto& curElemVolVars = this->curElemVolVars(); | |
612 | 253869 | const auto& elemFluxVarsCache = this->elemFluxVarsCache(); | |
613 | |||
614 | // get the vecor of the actual element residuals | ||
615 | 253869 | const auto origResiduals = this->evalLocalResidual(); | |
616 | |||
617 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
618 | // // | ||
619 | // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the // | ||
620 | // neighboring elements we do so by computing the derivatives of the fluxes which depend on the // | ||
621 | // actual element. In the actual element we evaluate the derivative of the entire residual. // | ||
622 | // // | ||
623 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
624 | |||
625 | // calculation of the source and storage derivatives | ||
626 |
4/4✓ Branch 0 taken 611000 times.
✓ Branch 1 taken 134266 times.
✓ Branch 2 taken 745266 times.
✓ Branch 3 taken 253869 times.
|
999135 | for (const auto& scv : scvs(fvGeometry)) |
627 | { | ||
628 | // dof index and corresponding actual pri vars | ||
629 |
2/2✓ Branch 0 taken 611000 times.
✓ Branch 1 taken 134266 times.
|
745266 | const auto dofIdx = scv.dofIndex(); |
630 | 745266 | const auto& volVars = curElemVolVars[scv]; | |
631 | |||
632 | // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping) | ||
633 | // only if the problem is instationary we add derivative of storage term | ||
634 | // TODO if e.g. porosity depends on all dofs in the element, we would have off-diagonal matrix entries!? | ||
635 |
2/2✓ Branch 0 taken 611000 times.
✓ Branch 1 taken 134266 times.
|
745266 | if (!this->assembler().isStationaryProblem()) |
636 | 611000 | this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx], | |
637 | problem, | ||
638 | element, | ||
639 | fvGeometry, | ||
640 | volVars, | ||
641 | scv); | ||
642 | |||
643 | // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping) | ||
644 | // add source term derivatives | ||
645 | 745266 | this->localResidual().addSourceDerivatives(A[dofIdx][dofIdx], | |
646 | problem, | ||
647 | element, | ||
648 | fvGeometry, | ||
649 | volVars, | ||
650 | scv); | ||
651 | } | ||
652 | |||
653 | // localJacobian[scvIdx][otherScvIdx][eqIdx][priVarIdx] of the fluxes | ||
654 |
4/4✓ Branch 0 taken 610161 times.
✓ Branch 1 taken 26963 times.
✓ Branch 2 taken 635914 times.
✓ Branch 3 taken 253264 times.
|
890388 | for (const auto& scvf : scvfs(fvGeometry)) |
655 | { | ||
656 |
2/2✓ Branch 0 taken 609556 times.
✓ Branch 1 taken 26358 times.
|
635914 | if (!scvf.boundary()) |
657 | { | ||
658 | // add flux term derivatives | ||
659 | 610161 | this->localResidual().addFluxDerivatives(A, | |
660 | problem, | ||
661 | element, | ||
662 | fvGeometry, | ||
663 | curElemVolVars, | ||
664 | elemFluxVarsCache, | ||
665 | scvf); | ||
666 | } | ||
667 | |||
668 | // the boundary gets special treatment to simplify | ||
669 | // for the user | ||
670 | else | ||
671 | { | ||
672 | 26358 | const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx()); | |
673 | 26358 | if (this->elemBcTypes().get(fvGeometry, insideScv).hasNeumann()) | |
674 | { | ||
675 | // add flux term derivatives | ||
676 | 26358 | this->localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()], | |
677 | problem, | ||
678 | element, | ||
679 | fvGeometry, | ||
680 | curElemVolVars, | ||
681 | elemFluxVarsCache, | ||
682 | scvf); | ||
683 | } | ||
684 | } | ||
685 | } | ||
686 | |||
687 | 253869 | return origResiduals; | |
688 | } | ||
689 | |||
690 | }; // implicit CVFEAssembler with analytic Jacobian | ||
691 | |||
692 | /*! | ||
693 | * \ingroup Assembly | ||
694 | * \ingroup CVFEDiscretization | ||
695 | * \brief Control volume finite element local assembler using analytic differentiation and explicit time discretization | ||
696 | */ | ||
697 | template<class TypeTag, class Assembler, class Implementation> | ||
698 |
2/8✓ Branch 0 taken 500000 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 5000 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1010000 | class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/false, Implementation> |
699 | : public CVFELocalAssemblerBase<TypeTag, Assembler, | ||
700 | Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false, Implementation>>, | ||
701 | false> | ||
702 | { | ||
703 | using ThisType = CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false, Implementation>; | ||
704 | using ParentType = CVFELocalAssemblerBase<TypeTag, Assembler, Detail::CVFE::Impl<Implementation, ThisType>, false>; | ||
705 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
706 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
707 | |||
708 | public: | ||
709 | |||
710 | using LocalResidual = typename ParentType::LocalResidual; | ||
711 | using ElementResidualVector = typename LocalResidual::ElementResidualVector; | ||
712 |
2/4✓ Branch 2 taken 500000 times.
✗ Branch 3 not taken.
✓ Branch 6 taken 5000 times.
✗ Branch 7 not taken.
|
505000 | using ParentType::ParentType; |
713 | |||
714 | /*! | ||
715 | * \brief Computes the derivatives with respect to the given element and adds them | ||
716 | * to the global matrix. | ||
717 | * | ||
718 | * \return The element residual at the current solution. | ||
719 | */ | ||
720 | template <class PartialReassembler = DefaultPartialReassembler> | ||
721 | 5000 | ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables, | |
722 | const PartialReassembler* partialReassembler = nullptr) | ||
723 | { | ||
724 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5000 times.
|
5000 | if (partialReassembler) |
725 | ✗ | DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization"); | |
726 | |||
727 | // get some aliases for convenience | ||
728 | 5000 | const auto& element = this->element(); | |
729 | 5000 | const auto& fvGeometry = this->fvGeometry(); | |
730 | 5000 | const auto& problem = this->asImp_().problem(); | |
731 | 5000 | const auto& curElemVolVars = this->curElemVolVars(); | |
732 | |||
733 | // get the vecor of the actual element residuals | ||
734 | 5000 | const auto origResiduals = this->evalLocalResidual(); | |
735 | |||
736 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
737 | // // | ||
738 | // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the // | ||
739 | // neighboring elements we do so by computing the derivatives of the fluxes which depend on the // | ||
740 | // actual element. In the actual element we evaluate the derivative of the entire residual. // | ||
741 | // // | ||
742 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
743 | |||
744 | // calculation of the source and storage derivatives | ||
745 |
2/2✓ Branch 1 taken 20000 times.
✓ Branch 2 taken 5000 times.
|
25000 | for (const auto& scv : scvs(fvGeometry)) |
746 | { | ||
747 | // dof index and corresponding actual pri vars | ||
748 | 20000 | const auto dofIdx = scv.dofIndex(); | |
749 | 20000 | const auto& volVars = curElemVolVars[scv]; | |
750 | |||
751 | // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping) | ||
752 | // only if the problem is instationary we add derivative of storage term | ||
753 | 20000 | this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx], | |
754 | problem, | ||
755 | element, | ||
756 | fvGeometry, | ||
757 | volVars, | ||
758 | scv); | ||
759 | } | ||
760 | |||
761 | 5000 | return origResiduals; | |
762 | } | ||
763 | |||
764 | }; // explicit CVFEAssembler with analytic Jacobian | ||
765 | |||
766 | } // end namespace Dumux | ||
767 | |||
768 | #endif | ||
769 |