Line | Branch | Exec | Source |
---|---|---|---|
1 | // -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*- | ||
2 | // vi: set et ts=4 sw=4 sts=4: | ||
3 | // | ||
4 | // SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder | ||
5 | // SPDX-License-Identifier: GPL-3.0-or-later | ||
6 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup Assembly | ||
10 | * \ingroup FaceCenteredStaggeredDiscretization | ||
11 | * \brief An assembler for Jacobian and residual contribution per element (face-centered staggered methods) | ||
12 | */ | ||
13 | #ifndef DUMUX_FC_LOCAL_ASSEMBLER_HH | ||
14 | #define DUMUX_FC_LOCAL_ASSEMBLER_HH | ||
15 | |||
16 | #include <dune/grid/common/gridenums.hh> | ||
17 | |||
18 | #include <dumux/common/properties.hh> | ||
19 | #include <dumux/common/parameters.hh> | ||
20 | #include <dumux/common/numericdifferentiation.hh> | ||
21 | #include <dumux/assembly/numericepsilon.hh> | ||
22 | #include <dumux/assembly/diffmethod.hh> | ||
23 | #include <dumux/assembly/fvlocalassemblerbase.hh> | ||
24 | #include <dumux/assembly/entitycolor.hh> | ||
25 | #include <dumux/assembly/partialreassembler.hh> | ||
26 | #include <dumux/discretization/facecentered/staggered/elementsolution.hh> | ||
27 | |||
28 | namespace Dumux { | ||
29 | |||
30 | namespace Detail { | ||
31 | |||
32 | struct NoOpFunctor | ||
33 | { | ||
34 | template<class... Args> | ||
35 | constexpr void operator()(Args&&...) const {} | ||
36 | }; | ||
37 | |||
38 | template<class T, class Default> | ||
39 | using NonVoidOrDefault_t = std::conditional_t<!std::is_same_v<T, void>, T, Default>; | ||
40 | |||
41 | } // end namespace Detail | ||
42 | |||
43 | /*! | ||
44 | * \ingroup Assembly | ||
45 | * \ingroup FaceCenteredStaggeredDiscretization | ||
46 | * \brief A base class for all local cell-centered assemblers | ||
47 | * \tparam TypeTag The TypeTag | ||
48 | * \tparam Assembler The assembler type | ||
49 | * \tparam Implementation The actual implementation | ||
50 | * \tparam implicit Specifies whether the time discretization is implicit or not (i.e. explicit) | ||
51 | */ | ||
52 | template<class TypeTag, class Assembler, class Implementation, bool implicit> | ||
53 | 2605838 | class FaceCenteredLocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit> | |
54 | { | ||
55 | using ParentType = FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit>; | ||
56 | using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView; | ||
57 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
58 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
59 | using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>; | ||
60 | |||
61 | static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq(); | ||
62 | |||
63 | public: | ||
64 | |||
65 | 2605838 | using ParentType::ParentType; | |
66 | |||
67 | 6728 | void bindLocalViews() | |
68 | { | ||
69 | 6728 | ParentType::bindLocalViews(); | |
70 | 26912 | this->elemBcTypes().update(this->asImp_().problem(), this->element(), this->fvGeometry()); | |
71 | 6728 | } | |
72 | |||
73 | /*! | ||
74 | * \brief Computes the derivatives with respect to the given element and adds them | ||
75 | * to the global matrix. The element residual is written into the right hand side. | ||
76 | */ | ||
77 | template <class ResidualVector, class PartialReassembler = DefaultPartialReassembler, class CouplingFunction = Detail::NoOpFunctor> | ||
78 | 1900838 | void assembleJacobianAndResidual(JacobianMatrix& jac, ResidualVector& res, GridVariables& gridVariables, | |
79 | const PartialReassembler* partialReassembler, | ||
80 | const CouplingFunction& maybeAssembleCouplingBlocks = CouplingFunction{}) | ||
81 | { | ||
82 | static_assert(!std::decay_t<decltype(this->asImp_().problem())>::enableInternalDirichletConstraints(), | ||
83 | "Internal Dirichlet constraints are currently not implemented for face-centered staggered models!"); | ||
84 | |||
85 | 1900838 | this->asImp_().bindLocalViews(); | |
86 | 3801676 | const auto& gridGeometry = this->asImp_().problem().gridGeometry(); | |
87 | 3801676 | const auto eIdxGlobal = gridGeometry.elementMapper().index(this->element()); | |
88 | if (partialReassembler | ||
89 | && partialReassembler->elementColor(eIdxGlobal) == EntityColor::green) | ||
90 | { | ||
91 | const auto residual = this->asImp_().evalLocalResidual(); // forward to the internal implementation | ||
92 | for (const auto& scv : scvs(this->fvGeometry())) | ||
93 | res[scv.dofIndex()] += residual[scv.localDofIndex()]; | ||
94 | |||
95 | // assemble the coupling blocks for coupled models (does nothing if not coupled) | ||
96 | maybeAssembleCouplingBlocks(residual); | ||
97 | } | ||
98 |
1/2✓ Branch 0 taken 1900838 times.
✗ Branch 1 not taken.
|
1900838 | else if (!this->elementIsGhost()) |
99 | { | ||
100 | 1900838 | const auto residual = this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables, partialReassembler); // forward to the internal implementation | |
101 | |||
102 |
2/2✓ Branch 1 taken 1900510 times.
✓ Branch 2 taken 328 times.
|
1900838 | if (this->element().partitionType() == Dune::InteriorEntity) |
103 | { | ||
104 |
6/6✓ Branch 0 taken 124144 times.
✓ Branch 1 taken 31036 times.
✓ Branch 2 taken 7486232 times.
✓ Branch 3 taken 1869474 times.
✓ Branch 4 taken 7486232 times.
✓ Branch 5 taken 1869474 times.
|
18928664 | for (const auto& scv : scvs(this->fvGeometry())) |
105 | 30441504 | res[scv.dofIndex()] += residual[scv.localDofIndex()]; | |
106 | } | ||
107 | else | ||
108 | { | ||
109 | // handle residual and matrix entries for parallel runs | ||
110 |
4/6✗ Branch 0 not taken.
✗ Branch 1 not taken.
✓ Branch 2 taken 1312 times.
✓ Branch 3 taken 328 times.
✓ Branch 4 taken 1312 times.
✓ Branch 5 taken 328 times.
|
1968 | for (const auto& scv : scvs(this->fvGeometry())) |
111 | { | ||
112 | 1312 | const auto& facet = this->element().template subEntity <1> (scv.indexInElement()); | |
113 | // make sure that the residual at border entities is consistent by adding the | ||
114 | // the contribution from the neighboring overlap element's scv | ||
115 |
4/4✓ Branch 0 taken 320 times.
✓ Branch 1 taken 992 times.
✓ Branch 2 taken 320 times.
✓ Branch 3 taken 992 times.
|
2624 | if (facet.partitionType() == Dune::BorderEntity) |
116 | 960 | res[scv.dofIndex()] += residual[scv.localDofIndex()]; | |
117 | |||
118 | // set the matrix entries of all DOFs within the overlap region (except the border DOF) | ||
119 | // to 1.0 and the residual entries to 0.0 | ||
120 | else | ||
121 | { | ||
122 | 992 | const auto idx = scv.dofIndex(); | |
123 |
0/4✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
1984 | jac[idx][idx] = 0.0; |
124 |
4/9✗ Branch 1 not taken.
✓ Branch 2 taken 992 times.
✓ Branch 3 taken 992 times.
✓ Branch 4 taken 992 times.
✓ Branch 5 taken 992 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
2976 | for (int i = 0; i < jac[idx][idx].size(); ++i) |
125 |
0/4✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
1984 | jac[idx][idx][i][i] = 1.0; |
126 | 2976 | res[idx] = 0; | |
127 | } | ||
128 | } | ||
129 | } | ||
130 | |||
131 | // assemble the coupling blocks for coupled models (does nothing if not coupled) | ||
132 | 1894110 | maybeAssembleCouplingBlocks(residual); | |
133 | } | ||
134 | else | ||
135 | ✗ | DUNE_THROW(Dune::NotImplemented, "Ghost elements not supported"); | |
136 | |||
137 | |||
138 | 2113474 | auto applyDirichlet = [&] (const auto& scvI, | |
139 | const auto& dirichletValues, | ||
140 | const auto eqIdx, | ||
141 | 106318 | const auto pvIdx) | |
142 | { | ||
143 | 212636 | res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx]; | |
144 | |||
145 | 106318 | auto& row = jac[scvI.dofIndex()]; | |
146 |
4/6✗ Branch 0 not taken.
✓ Branch 1 taken 745904 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 745904 times.
✓ Branch 4 taken 639586 times.
✓ Branch 5 taken 106318 times.
|
1385490 | for (auto col = row.begin(); col != row.end(); ++col) |
147 | 1279172 | row[col.index()][eqIdx] = 0.0; | |
148 | |||
149 | 212636 | jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0; | |
150 | |||
151 | // if a periodic dof has Dirichlet values also apply the same Dirichlet values to the other dof TODO periodic | ||
152 | // if (this->assembler().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex())) | ||
153 | // { | ||
154 | // const auto periodicDof = this->assembler().gridGeometry().periodicallyMappedDof(scvI.dofIndex()); | ||
155 | // res[periodicDof][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx]; | ||
156 | // const auto end = jac[periodicDof].end(); | ||
157 | // for (auto it = jac[periodicDof].begin(); it != end; ++it) | ||
158 | // (*it) = periodicDof != it.index() ? 0.0 : 1.0; | ||
159 | // } | ||
160 | }; | ||
161 | |||
162 | 3801676 | this->asImp_().enforceDirichletConstraints(applyDirichlet); | |
163 | 1900838 | } | |
164 | |||
165 | /*! | ||
166 | * \brief Computes the derivatives with respect to the given element and adds them | ||
167 | * to the global matrix. | ||
168 | */ | ||
169 | void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables) | ||
170 | { | ||
171 | this->asImp_().bindLocalViews(); | ||
172 | this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation | ||
173 | |||
174 | auto applyDirichlet = [&] (const auto& scvI, | ||
175 | const auto& dirichletValues, | ||
176 | const auto eqIdx, | ||
177 | const auto pvIdx) | ||
178 | { | ||
179 | auto& row = jac[scvI.dofIndex()]; | ||
180 | for (auto col = row.begin(); col != row.end(); ++col) | ||
181 | row[col.index()][eqIdx] = 0.0; | ||
182 | |||
183 | jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0; | ||
184 | }; | ||
185 | |||
186 | this->asImp_().enforceDirichletConstraints(applyDirichlet); | ||
187 | } | ||
188 | |||
189 | /*! | ||
190 | * \brief Assemble the residual only | ||
191 | */ | ||
192 | template<class ResidualVector> | ||
193 | 705000 | void assembleResidual(ResidualVector& res) | |
194 | { | ||
195 | 705000 | this->asImp_().bindLocalViews(); | |
196 |
0/2✗ Branch 0 not taken.
✗ Branch 1 not taken.
|
705000 | const auto residual = this->evalLocalResidual(); |
197 | |||
198 |
4/6✗ Branch 0 not taken.
✗ Branch 1 not taken.
✓ Branch 2 taken 2820000 times.
✓ Branch 3 taken 705000 times.
✓ Branch 4 taken 2820000 times.
✓ Branch 5 taken 705000 times.
|
7050000 | for (const auto& scv : scvs(this->fvGeometry())) |
199 | 11280000 | res[scv.dofIndex()] += residual[scv.localDofIndex()]; | |
200 | |||
201 | 705000 | auto applyDirichlet = [&] (const auto& scvI, | |
202 | const auto& dirichletValues, | ||
203 | const auto eqIdx, | ||
204 | ✗ | const auto pvIdx) | |
205 | { | ||
206 | 27600 | res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[eqIdx] - dirichletValues[pvIdx]; | |
207 | }; | ||
208 | |||
209 | 1410000 | this->asImp_().enforceDirichletConstraints(applyDirichlet); | |
210 | 705000 | } | |
211 | |||
212 | /*! | ||
213 | * \brief Enforce Dirichlet constraints | ||
214 | */ | ||
215 | template<typename ApplyFunction> | ||
216 | ✗ | void enforceDirichletConstraints(const ApplyFunction& applyDirichlet) | |
217 | { | ||
218 | // enforce Dirichlet boundary conditions | ||
219 | 2605838 | this->asImp_().evalDirichletBoundaries(applyDirichlet); | |
220 | // take care of internal Dirichlet constraints (if enabled) | ||
221 | 1900838 | this->asImp_().enforceInternalDirichletConstraints(applyDirichlet); | |
222 | ✗ | } | |
223 | |||
224 | /*! | ||
225 | * \brief Evaluates Dirichlet boundaries | ||
226 | */ | ||
227 | template< typename ApplyDirichletFunctionType > | ||
228 | 5211676 | void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet) | |
229 | { | ||
230 | // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0 | ||
231 | // and set the residual to (privar - dirichletvalue) | ||
232 |
4/4✓ Branch 0 taken 131379 times.
✓ Branch 1 taken 2474459 times.
✓ Branch 2 taken 131379 times.
✓ Branch 3 taken 2474459 times.
|
10423352 | if (this->elemBcTypes().hasDirichlet()) |
233 | { | ||
234 |
6/6✓ Branch 0 taken 41250 times.
✓ Branch 1 taken 3162 times.
✓ Branch 2 taken 1699648 times.
✓ Branch 3 taken 128217 times.
✓ Branch 4 taken 1699648 times.
✓ Branch 5 taken 128217 times.
|
4013636 | for (const auto& scvf : scvfs(this->fvGeometry())) |
235 | { | ||
236 |
4/4✓ Branch 0 taken 664474 times.
✓ Branch 1 taken 1076424 times.
✓ Branch 2 taken 135794 times.
✓ Branch 3 taken 528680 times.
|
3481796 | if (scvf.isFrontal() && scvf.boundary()) |
237 | { | ||
238 | 543176 | const auto bcTypes = this->elemBcTypes()[scvf.localIndex()]; | |
239 |
4/4✓ Branch 0 taken 134118 times.
✓ Branch 1 taken 1676 times.
✓ Branch 2 taken 134118 times.
✓ Branch 3 taken 1676 times.
|
543176 | if (bcTypes.hasDirichlet()) |
240 | { | ||
241 |
3/6✗ Branch 0 not taken.
✓ Branch 1 taken 130856 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 130856 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 130856 times.
|
804708 | const auto& scv = this->fvGeometry().scv(scvf.insideScvIdx()); |
242 | 516072 | const auto dirichletValues = this->asImp_().problem().dirichlet(this->element(), scvf); | |
243 | |||
244 | // set the Dirichlet conditions in residual and jacobian | ||
245 |
3/3✓ Branch 0 taken 123918 times.
✓ Branch 1 taken 134118 times.
✓ Branch 2 taken 10200 times.
|
536472 | for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) |
246 | { | ||
247 |
2/2✓ Branch 0 taken 134118 times.
✓ Branch 1 taken 270366 times.
|
808968 | for (int pvIdx = 0; pvIdx < GridView::dimension; ++pvIdx) |
248 | { | ||
249 |
4/4✓ Branch 0 taken 2300 times.
✓ Branch 1 taken 268066 times.
✓ Branch 2 taken 133908 times.
✓ Branch 3 taken 134158 times.
|
540732 | if (bcTypes.isDirichlet(pvIdx) && pvIdx == scv.dofAxis()) // TODO? |
250 | 267836 | applyDirichlet(scv, dirichletValues, eqIdx, pvIdx); | |
251 | } | ||
252 | } | ||
253 | } | ||
254 | } | ||
255 | } | ||
256 | } | ||
257 | 5211676 | } | |
258 | |||
259 | /*! | ||
260 | * \brief Update the coupling context for coupled models. | ||
261 | * \note This does nothing per default (not a coupled model). | ||
262 | */ | ||
263 | template<class... Args> | ||
264 | ✗ | void maybeUpdateCouplingContext(Args&&...) {} | |
265 | |||
266 | /*! | ||
267 | * \brief Update the additional domain derivatives for coupled models. | ||
268 | * \note This does nothing per default (not a coupled model). | ||
269 | */ | ||
270 | template<class... Args> | ||
271 | ✗ | void maybeEvalAdditionalDomainDerivatives(Args&&...) {} | |
272 | }; | ||
273 | |||
274 | /*! | ||
275 | * \ingroup Assembly | ||
276 | * \ingroup FaceCenteredStaggeredDiscretization | ||
277 | * \brief An assembler for Jacobian and residual contribution per element (Face-centered methods) | ||
278 | * \tparam TypeTag The TypeTag | ||
279 | * \tparam diffMethod The differentiation method to residual compute derivatives | ||
280 | * \tparam implicit Specifies whether the time discretization is implicit or not (i.e. explicit) | ||
281 | * \tparam Implementation The actual implementation, if void this class is the actual implementation | ||
282 | */ | ||
283 | template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, bool implicit = true, class Implementation = void> | ||
284 | class FaceCenteredLocalAssembler; | ||
285 | |||
286 | /*! | ||
287 | * \ingroup Assembly | ||
288 | * \ingroup FaceCenteredStaggeredDiscretization | ||
289 | * \brief Face-centered scheme local assembler using numeric differentiation and implicit time discretization | ||
290 | */ | ||
291 | template<class TypeTag, class Assembler, class Implementation> | ||
292 | 2605838 | class FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true, Implementation> | |
293 | : public FaceCenteredLocalAssemblerBase< | ||
294 | TypeTag, Assembler, | ||
295 | Detail::NonVoidOrDefault_t<Implementation, FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true, Implementation>>, | ||
296 | /*implicit=*/true | ||
297 | > | ||
298 | { | ||
299 | using ThisType = FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true, Implementation>; | ||
300 | using ParentType = FaceCenteredLocalAssemblerBase<TypeTag, Assembler, Detail::NonVoidOrDefault_t<Implementation, ThisType>, true>; | ||
301 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
302 | using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity; | ||
303 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
304 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
305 | using SubControlVolume = typename FVElementGeometry::SubControlVolume; | ||
306 | using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; | ||
307 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
308 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
309 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
310 | |||
311 | static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq(); | ||
312 | static constexpr bool enableGridFluxVarsCache = GetPropType<TypeTag, Properties::GridVariables>::GridFluxVariablesCache::cachingEnabled; | ||
313 | |||
314 | public: | ||
315 | |||
316 | using LocalResidual = typename ParentType::LocalResidual; | ||
317 | using ElementResidualVector = typename LocalResidual::ElementResidualVector; | ||
318 |
2/5✓ Branch 1 taken 2599110 times.
✓ Branch 2 taken 6728 times.
✗ Branch 3 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
2605838 | using ParentType::ParentType; |
319 | |||
320 | /*! | ||
321 | * \brief Computes the derivatives with respect to the given element and adds them | ||
322 | * to the global matrix. | ||
323 | * | ||
324 | * \return The element residual at the current solution. | ||
325 | */ | ||
326 | template <class PartialReassembler = DefaultPartialReassembler> | ||
327 | 1900838 | ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables, | |
328 | const PartialReassembler* partialReassembler = nullptr) | ||
329 | { | ||
330 | // get some aliases for convenience | ||
331 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 160 times.
|
1900838 | const auto& problem = this->asImp_().problem(); |
332 | 1900838 | const auto& element = this->element(); | |
333 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 160 times.
|
1900838 | const auto& fvGeometry = this->fvGeometry(); |
334 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 160 times.
|
1900838 | const auto& curSol = this->asImp_().curSol(); |
335 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 160 times.
|
1900838 | auto&& curElemVolVars = this->curElemVolVars(); |
336 | |||
337 | // get the vector of the actual element residuals | ||
338 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 160 times.
|
1900838 | const auto origResiduals = this->evalLocalResidual(); |
339 | |||
340 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
341 | // // | ||
342 | // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the // | ||
343 | // neighboring elements we do so by computing the derivatives of the fluxes which depend on the // | ||
344 | // actual element. In the actual element we evaluate the derivative of the entire residual. // | ||
345 | // // | ||
346 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
347 | |||
348 | // one residual per element facet | ||
349 | 1900838 | const auto numElementResiduals = fvGeometry.numScv(); | |
350 | |||
351 | // create the vector storing the partial derivatives | ||
352 | 1900838 | ElementResidualVector partialDerivs(numElementResiduals); | |
353 | |||
354 | 62430054 | const auto evalSource = [&](ElementResidualVector& residual, const SubControlVolume& scv) | |
355 | { | ||
356 | 60529216 | this->localResidual().evalSource(residual, problem, element, fvGeometry, curElemVolVars, scv); | |
357 | }; | ||
358 | |||
359 | 64369862 | const auto evalStorage = [&](ElementResidualVector& residual, const SubControlVolume& scv) | |
360 | { | ||
361 | 62469024 | this->localResidual().evalStorage(residual, problem, element, fvGeometry, this->prevElemVolVars(), curElemVolVars, scv); | |
362 | }; | ||
363 | |||
364 | 550834802 | const auto evalFlux = [&](ElementResidualVector& residual, const SubControlVolumeFace& scvf) | |
365 | { | ||
366 |
2/2✓ Branch 0 taken 182977540 times.
✓ Branch 1 taken 672 times.
|
182978212 | if (!scvf.processorBoundary()) |
367 | 548932620 | this->localResidual().evalFlux(residual, problem, element, fvGeometry, curElemVolVars, this->elemBcTypes(), this->elemFluxVarsCache(), scvf); | |
368 | }; | ||
369 | |||
370 | 60529216 | const auto evalDerivative = [&] (const auto& scvI, const auto& scvJ) | |
371 | { | ||
372 | // derivative w.r.t. own DOF | ||
373 |
2/2✓ Branch 0 taken 60529216 times.
✓ Branch 1 taken 60529216 times.
|
121058432 | for (int pvIdx = 0; pvIdx < numEq; pvIdx++) |
374 | { | ||
375 | 60529216 | partialDerivs = 0.0; | |
376 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 60529216 times.
✓ Branch 3 taken 60522976 times.
✗ Branch 4 not taken.
|
60751424 | const auto& otherElement = fvGeometry.gridGeometry().element(scvJ.elementIndex()); |
377 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 60529216 times.
✓ Branch 3 taken 1417544 times.
✗ Branch 4 not taken.
|
60529216 | auto otherElemSol = elementSolution(otherElement, curSol, fvGeometry.gridGeometry()); // TODO allow selective creation of elemsol (for one scv) |
378 |
4/4✓ Branch 0 taken 68 times.
✓ Branch 1 taken 59551244 times.
✓ Branch 2 taken 68 times.
✓ Branch 3 taken 59551244 times.
|
121058432 | auto& curOtherVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scvJ); |
379 | 60529216 | const VolumeVariables origOtherVolVars(curOtherVolVars); | |
380 | |||
381 | 1123128868 | auto evalResiduals = [&](Scalar priVar) | |
382 | { | ||
383 | // update the volume variables and compute element residual | ||
384 | 60529216 | otherElemSol[scvJ.localDofIndex()][pvIdx] = priVar; | |
385 | 60529216 | curOtherVolVars.update(otherElemSol, problem, otherElement, scvJ); | |
386 | 60529216 | this->asImp_().maybeUpdateCouplingContext(scvJ, otherElemSol, pvIdx); | |
387 | |||
388 | 60529216 | ElementResidualVector residual(numElementResiduals); | |
389 | 60529216 | residual = 0; | |
390 | |||
391 | 60529216 | evalSource(residual, scvI); | |
392 | |||
393 |
4/4✓ Branch 0 taken 31234512 times.
✓ Branch 1 taken 29294704 times.
✓ Branch 2 taken 31234512 times.
✓ Branch 3 taken 29090832 times.
|
120854560 | if (!this->assembler().isStationaryProblem()) |
394 | 31234512 | evalStorage(residual, scvI); | |
395 | |||
396 |
4/4✓ Branch 1 taken 182978212 times.
✓ Branch 2 taken 60529216 times.
✓ Branch 3 taken 182978212 times.
✓ Branch 4 taken 60529216 times.
|
243507428 | for (const auto& scvf : scvfs(fvGeometry, scvI)) |
397 | 182978212 | evalFlux(residual, scvf); | |
398 | |||
399 | 60529216 | return residual; | |
400 | }; | ||
401 | |||
402 | // derive the residuals numerically | ||
403 |
7/10✓ Branch 0 taken 72 times.
✓ Branch 1 taken 60529144 times.
✓ Branch 3 taken 60 times.
✓ Branch 4 taken 12 times.
✓ Branch 6 taken 60 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 60 times.
✗ Branch 10 not taken.
✓ Branch 12 taken 60 times.
✗ Branch 13 not taken.
|
60529216 | static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()}; |
404 |
7/10✓ Branch 0 taken 81 times.
✓ Branch 1 taken 60529135 times.
✓ Branch 3 taken 60 times.
✓ Branch 4 taken 21 times.
✓ Branch 6 taken 60 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 60 times.
✗ Branch 10 not taken.
✓ Branch 12 taken 60 times.
✗ Branch 13 not taken.
|
60529216 | static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod"); |
405 |
2/4✓ Branch 1 taken 60529216 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 60529216 times.
✗ Branch 5 not taken.
|
121058432 | NumericDifferentiation::partialDerivative(evalResiduals, otherElemSol[scvJ.localDofIndex()][pvIdx], partialDerivs, origResiduals, |
406 | 121058432 | eps_(otherElemSol[scvJ.localDofIndex()][pvIdx], pvIdx), numDiffMethod); | |
407 | |||
408 | 484238128 | const auto updateJacobian = [&]() | |
409 | { | ||
410 |
2/2✓ Branch 0 taken 60530096 times.
✓ Branch 1 taken 60530096 times.
|
121060192 | for (int eqIdx = 0; eqIdx < numEq; eqIdx++) |
411 | { | ||
412 | // A[i][col][eqIdx][pvIdx] is the rate of change of | ||
413 | // the residual of equation 'eqIdx' at dof 'i' | ||
414 | // depending on the primary variable 'pvIdx' at dof | ||
415 | // 'col'. | ||
416 | 181590288 | A[scvI.dofIndex()][scvJ.dofIndex()][eqIdx][pvIdx] += partialDerivs[scvI.localDofIndex()][eqIdx]; | |
417 | } | ||
418 | }; | ||
419 | |||
420 | using GeometryHelper = typename std::decay_t<decltype(fvGeometry.gridGeometry())>::GeometryHelper; | ||
421 | using LocalIntersectionMapper = typename std::decay_t<decltype(fvGeometry.gridGeometry())>::LocalIntersectionMapper; | ||
422 | 228448 | LocalIntersectionMapper localIsMapper; | |
423 | |||
424 |
4/9✗ Branch 0 not taken.
✓ Branch 1 taken 60529216 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 60307008 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 60529216 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 59978472 times.
✗ Branch 8 not taken.
|
60529216 | const bool isParallel = fvGeometry.gridGeometry().gridView().comm().size() > 1; |
425 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 228448 times.
|
228448 | if (isParallel) |
426 |
0/4✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
60506272 | localIsMapper.update(fvGeometry.gridGeometry().gridView(), element); |
427 | |||
428 |
3/4✓ Branch 1 taken 60529216 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 60299136 times.
✓ Branch 4 taken 1632 times.
|
60529216 | if (element.partitionType() == Dune::InteriorEntity) |
429 |
1/2✓ Branch 1 taken 60527584 times.
✗ Branch 2 not taken.
|
60527584 | updateJacobian(); |
430 | else | ||
431 | { | ||
432 | 1632 | const auto localIdxI = scvI.indexInElement(); | |
433 | 1632 | const auto localIdxJ = scvJ.indexInElement(); | |
434 | |||
435 |
0/4✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
1632 | const auto& facetI = GeometryHelper::facet(localIsMapper.refToRealIdx(localIdxI), element); |
436 | // add contribution of opposite scv lying within the overlap/ghost zone | ||
437 |
4/4✓ Branch 0 taken 640 times.
✓ Branch 1 taken 992 times.
✓ Branch 2 taken 640 times.
✓ Branch 3 taken 992 times.
|
3264 | if (facetI.partitionType() == Dune::BorderEntity && |
438 |
5/6✓ Branch 0 taken 320 times.
✓ Branch 1 taken 320 times.
✓ Branch 2 taken 320 times.
✓ Branch 3 taken 320 times.
✓ Branch 4 taken 320 times.
✗ Branch 5 not taken.
|
1280 | (localIdxJ == GeometryHelper::localOppositeIdx(localIdxI) || scvJ.dofIndex() == scvI.dofIndex())) |
439 |
1/2✓ Branch 1 taken 640 times.
✗ Branch 2 not taken.
|
640 | updateJacobian(); |
440 | } | ||
441 | |||
442 |
5/6✓ Branch 0 taken 102752 times.
✓ Branch 1 taken 60097928 times.
✓ Branch 3 taken 102752 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 101120 times.
✓ Branch 6 taken 1632 times.
|
60200680 | if (isParallel && element.partitionType() == Dune::InteriorEntity) |
443 | { | ||
444 | 101120 | const auto localIdxI = scvI.indexInElement(); | |
445 |
0/4✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
101120 | const auto& facetI = GeometryHelper::facet(localIsMapper.refToRealIdx(localIdxI), element); |
446 |
4/4✓ Branch 0 taken 2544 times.
✓ Branch 1 taken 98576 times.
✓ Branch 2 taken 2544 times.
✓ Branch 3 taken 98576 times.
|
202240 | if (facetI.partitionType() == Dune::BorderEntity) |
447 | { | ||
448 |
4/4✓ Branch 1 taken 7632 times.
✓ Branch 2 taken 2544 times.
✓ Branch 3 taken 7632 times.
✓ Branch 4 taken 2544 times.
|
10176 | for (const auto& scvf : scvfs(fvGeometry, scvI)) |
449 | { | ||
450 |
4/4✓ Branch 0 taken 5088 times.
✓ Branch 1 taken 2544 times.
✓ Branch 2 taken 4976 times.
✓ Branch 3 taken 112 times.
|
7632 | if (scvf.isFrontal() || scvf.boundary()) |
451 | continue; | ||
452 | |||
453 | // parallel scvs TODO drawing | ||
454 |
4/4✓ Branch 0 taken 624 times.
✓ Branch 1 taken 4352 times.
✓ Branch 2 taken 624 times.
✓ Branch 3 taken 4352 times.
|
9952 | if (scvf.outsideScvIdx() == scvJ.index()) |
455 |
1/2✓ Branch 1 taken 624 times.
✗ Branch 2 not taken.
|
624 | updateJacobian(); |
456 | else | ||
457 | { | ||
458 | // normal scvs | ||
459 | 4352 | const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf); | |
460 |
1/2✓ Branch 0 taken 4352 times.
✗ Branch 1 not taken.
|
4352 | if (orthogonalScvf.boundary()) |
461 | continue; | ||
462 | |||
463 |
8/8✓ Branch 0 taken 3728 times.
✓ Branch 1 taken 624 times.
✓ Branch 2 taken 3728 times.
✓ Branch 3 taken 624 times.
✓ Branch 4 taken 624 times.
✓ Branch 5 taken 3104 times.
✓ Branch 6 taken 624 times.
✓ Branch 7 taken 3104 times.
|
8704 | if (orthogonalScvf.insideScvIdx() == scvJ.index() || orthogonalScvf.outsideScvIdx() == scvJ.index()) |
464 |
1/2✓ Branch 1 taken 1248 times.
✗ Branch 2 not taken.
|
1248 | updateJacobian(); |
465 | } | ||
466 | } | ||
467 | } | ||
468 | } | ||
469 | |||
470 | // restore the original state of the scv's volume variables | ||
471 | 60529216 | curOtherVolVars = origOtherVolVars; | |
472 | |||
473 | // restore the original element solution | ||
474 | 121058432 | otherElemSol[scvJ.localDofIndex()][pvIdx] = curSol[scvJ.dofIndex()][pvIdx]; | |
475 | 120825032 | this->asImp_().maybeUpdateCouplingContext(scvJ, otherElemSol, pvIdx); | |
476 | // TODO additional dof dependencies | ||
477 | } | ||
478 | }; | ||
479 | |||
480 | // calculation of the derivatives | ||
481 |
5/5✓ Branch 0 taken 124144 times.
✓ Branch 1 taken 7518580 times.
✓ Branch 2 taken 1869802 times.
✓ Branch 3 taken 7487544 times.
✓ Branch 4 taken 1869802 times.
|
9543562 | for (const auto& scvI : scvs(fvGeometry)) |
482 | { | ||
483 | // derivative w.r.t. own DOFs | ||
484 | 7611688 | evalDerivative(scvI, scvI); | |
485 | |||
486 | // derivative w.r.t. other DOFs | ||
487 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 7611688 times.
|
7611688 | const auto& otherScvIndices = fvGeometry.gridGeometry().connectivityMap()[scvI.index()]; |
488 |
4/4✓ Branch 0 taken 52917528 times.
✓ Branch 1 taken 7611688 times.
✓ Branch 2 taken 52917528 times.
✓ Branch 3 taken 7611688 times.
|
75752592 | for (const auto globalJ : otherScvIndices) |
489 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 52063768 times.
|
52917528 | evalDerivative(scvI, fvGeometry.scv(globalJ)); |
490 | } | ||
491 | |||
492 | // evaluate additional derivatives that might arise from the coupling | ||
493 | 1900838 | this->asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables); | |
494 | |||
495 | 1900838 | return origResiduals; | |
496 | } | ||
497 | }; | ||
498 | |||
499 | |||
500 | /*! | ||
501 | * \ingroup Assembly | ||
502 | * \ingroup FaceCenteredStaggeredDiscretization | ||
503 | * \brief TODO docme | ||
504 | */ | ||
505 | template<class TypeTag, class Assembler, class Implementation> | ||
506 | class FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false, Implementation> | ||
507 | : public FaceCenteredLocalAssemblerBase< | ||
508 | TypeTag, Assembler, | ||
509 | Detail::NonVoidOrDefault_t<Implementation, FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false, Implementation>>, | ||
510 | /*implicit=*/false | ||
511 | > | ||
512 | { | ||
513 | using ThisType = FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false, Implementation>; | ||
514 | using ParentType = FaceCenteredLocalAssemblerBase<TypeTag, Assembler, Detail::NonVoidOrDefault_t<Implementation, ThisType>, false>; | ||
515 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
516 | using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity; | ||
517 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
518 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
519 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
520 | |||
521 | enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() }; | ||
522 | |||
523 | public: | ||
524 | using LocalResidual = typename ParentType::LocalResidual; | ||
525 | using ElementResidualVector = typename LocalResidual::ElementResidualVector; | ||
526 | using ParentType::ParentType; | ||
527 | |||
528 | /*! | ||
529 | * \brief Computes the derivatives with respect to the given element and adds them | ||
530 | * to the global matrix. | ||
531 | * | ||
532 | * \return The element residual at the current solution. | ||
533 | */ | ||
534 | template <class PartialReassembler = DefaultPartialReassembler> | ||
535 | ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables, | ||
536 | const PartialReassembler* partialReassembler = nullptr) | ||
537 | { | ||
538 | if (partialReassembler) | ||
539 | DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization"); | ||
540 | |||
541 | // get some aliases for convenience | ||
542 | const auto& problem = this->asImp_().problem(); | ||
543 | const auto& element = this->element(); | ||
544 | const auto& fvGeometry = this->fvGeometry(); | ||
545 | const auto& curSol = this->asImp_().curSol(); | ||
546 | auto&& curElemVolVars = this->curElemVolVars(); | ||
547 | |||
548 | // get the vector of the actual element residuals | ||
549 | const auto origResiduals = this->evalLocalResidual(); | ||
550 | const auto origStorageResiduals = this->evalLocalStorageResidual(); | ||
551 | |||
552 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
553 | // // | ||
554 | // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the // | ||
555 | // neighboring elements we do so by computing the derivatives of the fluxes which depend on the // | ||
556 | // actual element. In the actual element we evaluate the derivative of the entire residual. // | ||
557 | // // | ||
558 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
559 | |||
560 | // create the element solution | ||
561 | auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry()); | ||
562 | |||
563 | // create the vector storing the partial derivatives | ||
564 | ElementResidualVector partialDerivs(fvGeometry.numScv()); | ||
565 | |||
566 | // calculation of the derivatives | ||
567 | for (const auto& scv : scvs(fvGeometry)) | ||
568 | { | ||
569 | // dof index and corresponding actual pri vars | ||
570 | const auto dofIdx = scv.dofIndex(); | ||
571 | auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv); | ||
572 | const VolumeVariables origVolVars(curVolVars); | ||
573 | |||
574 | // calculate derivatives w.r.t to the privars at the dof at hand | ||
575 | for (int pvIdx = 0; pvIdx < numEq; pvIdx++) | ||
576 | { | ||
577 | partialDerivs = 0.0; | ||
578 | |||
579 | auto evalStorage = [&](Scalar priVar) | ||
580 | { | ||
581 | // auto partialDerivsTmp = partialDerivs; | ||
582 | elemSol[scv.localDofIndex()][pvIdx] = priVar; | ||
583 | curVolVars.update(elemSol, problem, element, scv); | ||
584 | return this->evalLocalStorageResidual(); | ||
585 | }; | ||
586 | |||
587 | // derive the residuals numerically | ||
588 | static const NumericEpsilon<Scalar, numEq> eps_{problem.paramGroup()}; | ||
589 | static const int numDiffMethod = getParamFromGroup<int>(problem.paramGroup(), "Assembly.NumericDifferenceMethod"); | ||
590 | NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origStorageResiduals, | ||
591 | eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod); | ||
592 | |||
593 | // update the global stiffness matrix with the current partial derivatives | ||
594 | for (int eqIdx = 0; eqIdx < numEq; eqIdx++) | ||
595 | { | ||
596 | // A[i][col][eqIdx][pvIdx] is the rate of change of | ||
597 | // the residual of equation 'eqIdx' at dof 'i' | ||
598 | // depending on the primary variable 'pvIdx' at dof | ||
599 | // 'col'. | ||
600 | A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx]; | ||
601 | } | ||
602 | |||
603 | // restore the original state of the scv's volume variables | ||
604 | curVolVars = origVolVars; | ||
605 | |||
606 | // restore the original element solution | ||
607 | elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx]; | ||
608 | } | ||
609 | } | ||
610 | return origResiduals; | ||
611 | } | ||
612 | }; | ||
613 | |||
614 | /*! | ||
615 | * \ingroup Assembly | ||
616 | * \ingroup FaceCenteredStaggeredDiscretization | ||
617 | * \brief TODO docme | ||
618 | */ | ||
619 | template<class TypeTag, class Assembler, class Implementation> | ||
620 | class FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true, Implementation> | ||
621 | : public FaceCenteredLocalAssemblerBase< | ||
622 | TypeTag, Assembler, | ||
623 | FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true, Implementation>, | ||
624 | /*implicit=*/true | ||
625 | > | ||
626 | { | ||
627 | using ThisType = FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true, Implementation>; | ||
628 | using ParentType = FaceCenteredLocalAssemblerBase<TypeTag, Assembler, Detail::NonVoidOrDefault_t<Implementation, ThisType>, true>; | ||
629 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
630 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
631 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
632 | |||
633 | enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() }; | ||
634 | |||
635 | public: | ||
636 | using LocalResidual = typename ParentType::LocalResidual; | ||
637 | using ElementResidualVector = typename LocalResidual::ElementResidualVector; | ||
638 | using ParentType::ParentType; | ||
639 | |||
640 | /*! | ||
641 | * \brief Computes the derivatives with respect to the given element and adds them | ||
642 | * to the global matrix. | ||
643 | * | ||
644 | * \return The element residual at the current solution. | ||
645 | */ | ||
646 | template <class PartialReassembler = DefaultPartialReassembler> | ||
647 | ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables, | ||
648 | const PartialReassembler* partialReassembler = nullptr) | ||
649 | { | ||
650 | if (partialReassembler) | ||
651 | DUNE_THROW(Dune::NotImplemented, "partial reassembly for analytic differentiation"); | ||
652 | |||
653 | // get some aliases for convenience | ||
654 | const auto& element = this->element(); | ||
655 | const auto& fvGeometry = this->fvGeometry(); | ||
656 | const auto& problem = this->asImp_().problem(); | ||
657 | const auto& curElemVolVars = this->curElemVolVars(); | ||
658 | const auto& elemFluxVarsCache = this->elemFluxVarsCache(); | ||
659 | |||
660 | // get the vecor of the actual element residuals | ||
661 | const auto origResiduals = this->evalLocalResidual(); | ||
662 | |||
663 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
664 | // // | ||
665 | // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the // | ||
666 | // neighboring elements we do so by computing the derivatives of the fluxes which depend on the // | ||
667 | // actual element. In the actual element we evaluate the derivative of the entire residual. // | ||
668 | // // | ||
669 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
670 | |||
671 | // calculation of the source and storage derivatives | ||
672 | for (const auto& scv : scvs(fvGeometry)) | ||
673 | { | ||
674 | // dof index and corresponding actual pri vars | ||
675 | const auto dofIdx = scv.dofIndex(); | ||
676 | const auto& volVars = curElemVolVars[scv]; | ||
677 | |||
678 | // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping) | ||
679 | // only if the problem is instationary we add derivative of storage term | ||
680 | // TODO if e.g. porosity depends on all dofs in the element, we would have off-diagonal matrix entries!? | ||
681 | if (!this->assembler().isStationaryProblem()) | ||
682 | this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx], | ||
683 | problem, | ||
684 | element, | ||
685 | fvGeometry, | ||
686 | volVars, | ||
687 | scv); | ||
688 | |||
689 | // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping) | ||
690 | // add source term derivatives | ||
691 | this->localResidual().addSourceDerivatives(A[dofIdx][dofIdx], | ||
692 | problem, | ||
693 | element, | ||
694 | fvGeometry, | ||
695 | volVars, | ||
696 | scv); | ||
697 | } | ||
698 | |||
699 | // localJacobian[scvIdx][otherScvIdx][eqIdx][priVarIdx] of the fluxes | ||
700 | for (const auto& scvf : scvfs(fvGeometry)) | ||
701 | { | ||
702 | if (!scvf.boundary()) | ||
703 | { | ||
704 | // add flux term derivatives | ||
705 | this->localResidual().addFluxDerivatives(A, | ||
706 | problem, | ||
707 | element, | ||
708 | fvGeometry, | ||
709 | curElemVolVars, | ||
710 | elemFluxVarsCache, | ||
711 | scvf); | ||
712 | } | ||
713 | |||
714 | // the boundary gets special treatment to simplify | ||
715 | // for the user | ||
716 | else | ||
717 | { | ||
718 | const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx()); | ||
719 | if (this->elemBcTypes()[insideScv.localDofIndex()].hasNeumann()) | ||
720 | { | ||
721 | // add flux term derivatives | ||
722 | this->localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()], | ||
723 | problem, | ||
724 | element, | ||
725 | fvGeometry, | ||
726 | curElemVolVars, | ||
727 | elemFluxVarsCache, | ||
728 | scvf); | ||
729 | } | ||
730 | } | ||
731 | } | ||
732 | |||
733 | return origResiduals; | ||
734 | } | ||
735 | }; | ||
736 | |||
737 | /*! | ||
738 | * \ingroup Assembly | ||
739 | * \ingroup FaceCenteredStaggeredDiscretization | ||
740 | * \brief TODO docme | ||
741 | */ | ||
742 | template<class TypeTag, class Assembler, class Implementation> | ||
743 | class FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/false, Implementation> | ||
744 | : public FaceCenteredLocalAssemblerBase< | ||
745 | TypeTag, Assembler, | ||
746 | FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false, Implementation>, | ||
747 | /*implicit=*/false | ||
748 | > | ||
749 | { | ||
750 | using ThisType = FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false, Implementation>; | ||
751 | using ParentType = FaceCenteredLocalAssemblerBase<TypeTag, Assembler, Detail::NonVoidOrDefault_t<Implementation, ThisType>, false>; | ||
752 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
753 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
754 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
755 | |||
756 | enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() }; | ||
757 | |||
758 | public: | ||
759 | using LocalResidual = typename ParentType::LocalResidual; | ||
760 | using ElementResidualVector = typename LocalResidual::ElementResidualVector; | ||
761 | using ParentType::ParentType; | ||
762 | |||
763 | /*! | ||
764 | * \brief Computes the derivatives with respect to the given element and adds them | ||
765 | * to the global matrix. | ||
766 | * | ||
767 | * \return The element residual at the current solution. | ||
768 | */ | ||
769 | template <class PartialReassembler = DefaultPartialReassembler> | ||
770 | ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables, | ||
771 | const PartialReassembler* partialReassembler = nullptr) | ||
772 | { | ||
773 | if (partialReassembler) | ||
774 | DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization"); | ||
775 | |||
776 | // get some aliases for convenience | ||
777 | const auto& element = this->element(); | ||
778 | const auto& fvGeometry = this->fvGeometry(); | ||
779 | const auto& problem = this->asImp_().problem(); | ||
780 | const auto& curElemVolVars = this->curElemVolVars(); | ||
781 | |||
782 | // get the vector of the actual element residuals | ||
783 | const auto origResiduals = this->evalLocalResidual(); | ||
784 | |||
785 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
786 | // // | ||
787 | // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the // | ||
788 | // neighboring elements we do so by computing the derivatives of the fluxes which depend on the // | ||
789 | // actual element. In the actual element we evaluate the derivative of the entire residual. // | ||
790 | // // | ||
791 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
792 | |||
793 | // calculation of the source and storage derivatives | ||
794 | for (const auto& scv : scvs(fvGeometry)) | ||
795 | { | ||
796 | // dof index and corresponding actual pri vars | ||
797 | const auto dofIdx = scv.dofIndex(); | ||
798 | const auto& volVars = curElemVolVars[scv]; | ||
799 | |||
800 | // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping) | ||
801 | // only if the problem is instationary we add derivative of storage term | ||
802 | this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx], | ||
803 | problem, | ||
804 | element, | ||
805 | fvGeometry, | ||
806 | volVars, | ||
807 | scv); | ||
808 | } | ||
809 | |||
810 | return origResiduals; | ||
811 | } | ||
812 | }; | ||
813 | |||
814 | } // end namespace Dumux | ||
815 | |||
816 | #endif | ||
817 |