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