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