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 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 | 2761797 | 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 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
61 | |||
62 | static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq(); | ||
63 | |||
64 | public: | ||
65 | |||
66 | 2761797 | using ParentType::ParentType; | |
67 | |||
68 | 6728 | void bindLocalViews() | |
69 | { | ||
70 | 6728 | ParentType::bindLocalViews(); | |
71 | 6728 | this->elemBcTypes().update(this->asImp_().problem(), this->element(), this->fvGeometry()); | |
72 | 6728 | } | |
73 | |||
74 | /*! | ||
75 | * \brief Computes the derivatives with respect to the given element and adds them | ||
76 | * to the global matrix. The element residual is written into the right hand side. | ||
77 | */ | ||
78 | template <class ResidualVector, class PartialReassembler = DefaultPartialReassembler, class CouplingFunction = Detail::NoOpFunctor> | ||
79 | 2056797 | void assembleJacobianAndResidual(JacobianMatrix& jac, ResidualVector& res, GridVariables& gridVariables, | |
80 | const PartialReassembler* partialReassembler, | ||
81 | const CouplingFunction& maybeAssembleCouplingBlocks = CouplingFunction{}) | ||
82 | { | ||
83 | static_assert(!std::decay_t<decltype(this->asImp_().problem())>::enableInternalDirichletConstraints(), | ||
84 | "Internal Dirichlet constraints are currently not implemented for face-centered staggered models!"); | ||
85 | |||
86 | 2056797 | this->asImp_().bindLocalViews(); | |
87 | 2056797 | const auto& gridGeometry = this->asImp_().problem().gridGeometry(); | |
88 | 2056797 | const auto eIdxGlobal = gridGeometry.elementMapper().index(this->element()); | |
89 | if (partialReassembler | ||
90 | && partialReassembler->elementColor(eIdxGlobal) == EntityColor::green) | ||
91 | { | ||
92 | const auto residual = this->asImp_().evalLocalResidual(); // forward to the internal implementation | ||
93 | for (const auto& scv : scvs(this->fvGeometry())) | ||
94 | res[scv.dofIndex()] += residual[scv.localDofIndex()]; | ||
95 | |||
96 | // assemble the coupling blocks for coupled models (does nothing if not coupled) | ||
97 | maybeAssembleCouplingBlocks(residual); | ||
98 | } | ||
99 |
1/2✓ Branch 0 taken 2056797 times.
✗ Branch 1 not taken.
|
2056797 | else if (!this->elementIsGhost()) |
100 | { | ||
101 | 2056797 | const auto residual = this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables, partialReassembler); // forward to the internal implementation | |
102 | |||
103 |
2/2✓ Branch 1 taken 2056469 times.
✓ Branch 2 taken 328 times.
|
2056797 | if (this->element().partitionType() == Dune::InteriorEntity) |
104 | { | ||
105 |
3/3✓ Branch 1 taken 8146762 times.
✓ Branch 2 taken 2025433 times.
✓ Branch 0 taken 124144 times.
|
10296339 | for (const auto& scv : scvs(this->fvGeometry())) |
106 | 8239870 | res[scv.dofIndex()] += residual[scv.localDofIndex()]; | |
107 | } | ||
108 | else | ||
109 | { | ||
110 | // handle residual and matrix entries for parallel runs | ||
111 |
2/3✓ Branch 1 taken 1312 times.
✓ Branch 2 taken 328 times.
✗ Branch 0 not taken.
|
1640 | for (const auto& scv : scvs(this->fvGeometry())) |
112 | { | ||
113 |
2/5✓ Branch 1 taken 992 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 0 taken 320 times.
|
1312 | const auto& facet = this->element().template subEntity <1> (scv.indexInElement()); |
114 | // make sure that the residual at border entities is consistent by adding the | ||
115 | // the contribution from the neighboring overlap element's scv | ||
116 |
2/3✓ Branch 1 taken 992 times.
✗ Branch 2 not taken.
✓ Branch 0 taken 320 times.
|
1312 | if (facet.partitionType() == Dune::BorderEntity) |
117 | 1312 | res[scv.dofIndex()] += residual[scv.localDofIndex()]; | |
118 | |||
119 | // set the matrix entries of all DOFs within the overlap region (except the border DOF) | ||
120 | // to 1.0 and the residual entries to 0.0 | ||
121 | else | ||
122 | { | ||
123 | 992 | const auto idx = scv.dofIndex(); | |
124 |
1/2✓ Branch 1 taken 992 times.
✗ Branch 2 not taken.
|
992 | jac[idx][idx] = 0.0; |
125 |
3/4✓ Branch 1 taken 1984 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 992 times.
✓ Branch 4 taken 992 times.
|
1984 | for (int i = 0; i < jac[idx][idx].size(); ++i) |
126 |
1/2✓ Branch 1 taken 992 times.
✗ Branch 2 not taken.
|
992 | jac[idx][idx][i][i] = 1.0; |
127 | 992 | res[idx] = 0; | |
128 | } | ||
129 | } | ||
130 | } | ||
131 | |||
132 | // assemble the coupling blocks for coupled models (does nothing if not coupled) | ||
133 | 2050069 | maybeAssembleCouplingBlocks(residual); | |
134 | } | ||
135 | else | ||
136 | ✗ | DUNE_THROW(Dune::NotImplemented, "Ghost elements not supported"); | |
137 | |||
138 | |||
139 | 2311185 | auto applyDirichlet = [&] (const auto& scvI, | |
140 | const auto& dirichletValues, | ||
141 | const auto eqIdx, | ||
142 | const auto pvIdx) | ||
143 | { | ||
144 | 127194 | res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx]; | |
145 | |||
146 | 127194 | auto& row = jac[scvI.dofIndex()]; | |
147 |
3/4✗ Branch 0 not taken.
✓ Branch 1 taken 899078 times.
✓ Branch 2 taken 771884 times.
✓ Branch 3 taken 127194 times.
|
899078 | for (auto col = row.begin(); col != row.end(); ++col) |
148 | 771884 | row[col.index()][eqIdx] = 0.0; | |
149 | |||
150 | 127194 | jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0; | |
151 | |||
152 | // if a periodic dof has Dirichlet values also apply the same Dirichlet values to the other dof | ||
153 | 254388 | if (this->asImp_().problem().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex())) | |
154 | { | ||
155 | ✗ | const auto periodicDof = this->asImp_().problem().gridGeometry().periodicallyMappedDof(scvI.dofIndex()); | |
156 | ✗ | res[periodicDof][eqIdx] = this->asImp_().curSol()[periodicDof][pvIdx] - dirichletValues[pvIdx]; | |
157 | |||
158 | ✗ | auto& rowP = jac[periodicDof]; | |
159 | ✗ | for (auto col = rowP.begin(); col != rowP.end(); ++col) | |
160 | ✗ | row[col.index()][eqIdx] = 0.0; | |
161 | |||
162 | ✗ | rowP[periodicDof][eqIdx][pvIdx] = 1.0; | |
163 | } | ||
164 | }; | ||
165 | |||
166 | 4113594 | this->asImp_().enforceDirichletConstraints(applyDirichlet); | |
167 | 2056797 | } | |
168 | |||
169 | /*! | ||
170 | * \brief Computes the derivatives with respect to the given element and adds them | ||
171 | * to the global matrix. | ||
172 | */ | ||
173 | void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables) | ||
174 | { | ||
175 | this->asImp_().bindLocalViews(); | ||
176 | this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation | ||
177 | |||
178 | auto applyDirichlet = [&] (const auto& scvI, | ||
179 | const auto& dirichletValues, | ||
180 | const auto eqIdx, | ||
181 | const auto pvIdx) | ||
182 | { | ||
183 | auto& row = jac[scvI.dofIndex()]; | ||
184 | for (auto col = row.begin(); col != row.end(); ++col) | ||
185 | row[col.index()][eqIdx] = 0.0; | ||
186 | |||
187 | jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0; | ||
188 | }; | ||
189 | |||
190 | this->asImp_().enforceDirichletConstraints(applyDirichlet); | ||
191 | } | ||
192 | |||
193 | /*! | ||
194 | * \brief Assemble the residual only | ||
195 | */ | ||
196 | template<class ResidualVector> | ||
197 | 705000 | void assembleResidual(ResidualVector& res) | |
198 | { | ||
199 | 705000 | this->asImp_().bindLocalViews(); | |
200 | 705000 | const auto residual = this->evalLocalResidual(); | |
201 | |||
202 |
2/3✓ Branch 1 taken 2820000 times.
✓ Branch 2 taken 705000 times.
✗ Branch 0 not taken.
|
3525000 | for (const auto& scv : scvs(this->fvGeometry())) |
203 | 2820000 | res[scv.dofIndex()] += residual[scv.localDofIndex()]; | |
204 | |||
205 | 732600 | auto applyDirichlet = [&] (const auto& scvI, | |
206 | const auto& dirichletValues, | ||
207 | const auto eqIdx, | ||
208 | const auto pvIdx) | ||
209 | { | ||
210 | 27600 | res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx]; | |
211 | }; | ||
212 | |||
213 | 1410000 | this->asImp_().enforceDirichletConstraints(applyDirichlet); | |
214 | 705000 | } | |
215 | |||
216 | /*! | ||
217 | * \brief Enforce Dirichlet constraints | ||
218 | */ | ||
219 | template<typename ApplyFunction> | ||
220 | 2761797 | void enforceDirichletConstraints(const ApplyFunction& applyDirichlet) | |
221 | { | ||
222 | // enforce Dirichlet boundary conditions | ||
223 | 2761797 | this->asImp_().evalDirichletBoundaries(applyDirichlet); | |
224 | // take care of internal Dirichlet constraints (if enabled) | ||
225 | 2056797 | this->asImp_().enforceInternalDirichletConstraints(applyDirichlet); | |
226 | } | ||
227 | |||
228 | /*! | ||
229 | * \brief Evaluates Dirichlet boundaries | ||
230 | */ | ||
231 | template< typename ApplyDirichletFunctionType > | ||
232 | 5523594 | void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet) | |
233 | { | ||
234 | // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0 | ||
235 | // and set the residual to (privar - dirichletvalue) | ||
236 |
2/2✓ Branch 0 taken 148649 times.
✓ Branch 1 taken 2613148 times.
|
5523594 | if (this->elemBcTypes().hasDirichlet()) |
237 | { | ||
238 |
5/5✓ Branch 1 taken 769306 times.
✓ Branch 2 taken 1265962 times.
✓ Branch 3 taken 1971884 times.
✓ Branch 4 taken 145487 times.
✓ Branch 0 taken 15954 times.
|
5805262 | for (const auto& scvf : scvfs(this->fvGeometry())) |
239 | { | ||
240 |
4/4✓ Branch 0 taken 759964 times.
✓ Branch 1 taken 1250008 times.
✓ Branch 2 taken 157776 times.
✓ Branch 3 taken 602188 times.
|
4019944 | if (scvf.isFrontal() && scvf.boundary()) |
241 | { | ||
242 |
2/2✓ Branch 1 taken 154994 times.
✓ Branch 2 taken 2782 times.
|
315552 | const auto bcTypes = this->elemBcTypes()[scvf.localIndex()]; |
243 |
2/2✓ Branch 0 taken 154994 times.
✓ Branch 1 taken 2782 times.
|
315552 | if (bcTypes.hasDirichlet()) |
244 | { | ||
245 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 151732 times.
|
309988 | const auto& scv = this->fvGeometry().scv(scvf.insideScvIdx()); |
246 |
2/3✓ Branch 0 taken 3136 times.
✓ Branch 1 taken 21328 times.
✗ Branch 2 not taken.
|
309988 | const auto dirichletValues = this->asImp_().problem().dirichlet(this->element(), scvf); |
247 | |||
248 | // set the Dirichlet conditions in residual and jacobian | ||
249 |
2/2✓ Branch 0 taken 154994 times.
✓ Branch 1 taken 154994 times.
|
619976 | for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) |
250 | { | ||
251 | static_assert(numEq == 1, "Not yet implemented for more than one vector-valued primary variable"); | ||
252 |
2/2✓ Branch 0 taken 200 times.
✓ Branch 1 taken 154794 times.
|
309988 | const int pvIdx = eqIdx; |
253 | 309988 | const int componentIdx = scv.dofAxis(); | |
254 |
3/4✓ Branch 0 taken 200 times.
✓ Branch 1 taken 154794 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 200 times.
|
309988 | if (bcTypes.isDirichlet(componentIdx)) |
255 |
1/2✓ Branch 1 taken 11920 times.
✗ Branch 2 not taken.
|
309588 | applyDirichlet(scv, std::array<Scalar,1>{{dirichletValues[componentIdx]}}, eqIdx, pvIdx); |
256 | } | ||
257 | } | ||
258 | } | ||
259 | } | ||
260 | } | ||
261 | 5523594 | } | |
262 | |||
263 | /*! | ||
264 | * \brief Update the coupling context for coupled models. | ||
265 | * \note This does nothing per default (not a coupled model). | ||
266 | */ | ||
267 | template<class... Args> | ||
268 | void maybeUpdateCouplingContext(Args&&...) {} | ||
269 | |||
270 | /*! | ||
271 | * \brief Update the additional domain derivatives for coupled models. | ||
272 | * \note This does nothing per default (not a coupled model). | ||
273 | */ | ||
274 | template<class... Args> | ||
275 | void maybeEvalAdditionalDomainDerivatives(Args&&...) {} | ||
276 | }; | ||
277 | |||
278 | /*! | ||
279 | * \ingroup Assembly | ||
280 | * \ingroup FaceCenteredStaggeredDiscretization | ||
281 | * \brief An assembler for Jacobian and residual contribution per element (Face-centered methods) | ||
282 | * \tparam TypeTag The TypeTag | ||
283 | * \tparam diffMethod The differentiation method to residual compute derivatives | ||
284 | * \tparam implicit Specifies whether the time discretization is implicit or not (i.e. explicit) | ||
285 | * \tparam Implementation The actual implementation, if void this class is the actual implementation | ||
286 | */ | ||
287 | template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, bool implicit = true, class Implementation = void> | ||
288 | class FaceCenteredLocalAssembler; | ||
289 | |||
290 | /*! | ||
291 | * \ingroup Assembly | ||
292 | * \ingroup FaceCenteredStaggeredDiscretization | ||
293 | * \brief Face-centered scheme local assembler using numeric differentiation and implicit time discretization | ||
294 | */ | ||
295 | template<class TypeTag, class Assembler, class Implementation> | ||
296 | 2761797 | class FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true, Implementation> | |
297 | : public FaceCenteredLocalAssemblerBase< | ||
298 | TypeTag, Assembler, | ||
299 | Detail::NonVoidOrDefault_t<Implementation, FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true, Implementation>>, | ||
300 | /*implicit=*/true | ||
301 | > | ||
302 | { | ||
303 | using ThisType = FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true, Implementation>; | ||
304 | using ParentType = FaceCenteredLocalAssemblerBase<TypeTag, Assembler, Detail::NonVoidOrDefault_t<Implementation, ThisType>, true>; | ||
305 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
306 | using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity; | ||
307 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
308 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
309 | using SubControlVolume = typename FVElementGeometry::SubControlVolume; | ||
310 | using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; | ||
311 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
312 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
313 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
314 | |||
315 | static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq(); | ||
316 | static constexpr bool enableGridFluxVarsCache = GetPropType<TypeTag, Properties::GridVariables>::GridFluxVariablesCache::cachingEnabled; | ||
317 | |||
318 | public: | ||
319 | |||
320 | using LocalResidual = typename ParentType::LocalResidual; | ||
321 | using ElementResidualVector = typename LocalResidual::ElementResidualVector; | ||
322 |
2/5✓ Branch 2 taken 6728 times.
✗ Branch 3 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✓ Branch 1 taken 2755069 times.
|
2761797 | using ParentType::ParentType; |
323 | |||
324 | /*! | ||
325 | * \brief Computes the derivatives with respect to the given element and adds them | ||
326 | * to the global matrix. | ||
327 | * | ||
328 | * \return The element residual at the current solution. | ||
329 | */ | ||
330 | template <class PartialReassembler = DefaultPartialReassembler> | ||
331 | 2056797 | ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables, | |
332 | const PartialReassembler* partialReassembler = nullptr) | ||
333 | { | ||
334 | // get some aliases for convenience | ||
335 | 2056797 | const auto& problem = this->asImp_().problem(); | |
336 | 2056797 | const auto& element = this->element(); | |
337 | 2056797 | const auto& fvGeometry = this->fvGeometry(); | |
338 | 2056797 | const auto& curSol = this->asImp_().curSol(); | |
339 | 2056797 | auto&& curElemVolVars = this->curElemVolVars(); | |
340 | |||
341 | // get the vector of the actual element residuals | ||
342 | 2056797 | const auto origResiduals = this->evalLocalResidual(); | |
343 | |||
344 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
345 | // // | ||
346 | // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the // | ||
347 | // neighboring elements we do so by computing the derivatives of the fluxes which depend on the // | ||
348 | // actual element. In the actual element we evaluate the derivative of the entire residual. // | ||
349 | // // | ||
350 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
351 | |||
352 | // one residual per element facet | ||
353 | 2056797 | const auto numElementResiduals = fvGeometry.numScv(); | |
354 | |||
355 | // create the vector storing the partial derivatives | ||
356 | 2056797 | ElementResidualVector partialDerivs(numElementResiduals); | |
357 | |||
358 | 68339673 | const auto evalSource = [&](ElementResidualVector& residual, const SubControlVolume& scv) | |
359 | { | ||
360 | 66282876 | this->localResidual().evalSource(residual, problem, element, fvGeometry, curElemVolVars, scv); | |
361 | }; | ||
362 | |||
363 | 72527277 | const auto evalStorage = [&](ElementResidualVector& residual, const SubControlVolume& scv) | |
364 | { | ||
365 | 35235240 | this->localResidual().evalStorage(residual, problem, element, fvGeometry, this->prevElemVolVars(), curElemVolVars, scv); | |
366 | }; | ||
367 | |||
368 |
2/2✓ Branch 0 taken 200866216 times.
✓ Branch 1 taken 3328 times.
|
200869544 | const auto evalFlux = [&](ElementResidualVector& residual, const SubControlVolumeFace& scvf) |
369 | { | ||
370 |
2/2✓ Branch 0 taken 200866216 times.
✓ Branch 1 taken 3328 times.
|
200869544 | if (!scvf.processorBoundary()) |
371 | 200866216 | this->localResidual().evalFlux(residual, problem, element, fvGeometry, curElemVolVars, this->elemBcTypes(), this->elemFluxVarsCache(), scvf); | |
372 | }; | ||
373 | |||
374 | 65562904 | const auto evalDerivative = [&] (const auto& scvI, const auto& scvJ) | |
375 | { | ||
376 | // derivative w.r.t. own DOF | ||
377 |
2/2✓ Branch 0 taken 65531868 times.
✓ Branch 1 taken 65531868 times.
|
131063736 | for (int pvIdx = 0; pvIdx < numEq; pvIdx++) |
378 | { | ||
379 | 65531868 | partialDerivs = 0.0; | |
380 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 65531868 times.
|
65531868 | const auto& otherElement = fvGeometry.gridGeometry().element(scvJ.elementIndex()); |
381 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 64742940 times.
✓ Branch 3 taken 222208 times.
✗ Branch 4 not taken.
|
65531868 | auto otherElemSol = elementSolution(otherElement, curSol, fvGeometry.gridGeometry()); // TODO allow selective creation of elemsol (for one scv) |
382 |
2/2✓ Branch 0 taken 71 times.
✓ Branch 1 taken 64553893 times.
|
65531868 | auto& curOtherVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scvJ); |
383 | 65531868 | const VolumeVariables origOtherVolVars(curOtherVolVars); | |
384 | |||
385 | 131814744 | auto evalResiduals = [&](Scalar priVar) | |
386 | { | ||
387 | // update the volume variables and compute element residual | ||
388 | 66282876 | otherElemSol[scvJ.localDofIndex()][pvIdx] = priVar; | |
389 | 66282876 | curOtherVolVars.update(otherElemSol, problem, otherElement, scvJ); | |
390 | 66282876 | this->asImp_().maybeUpdateCouplingContext(scvJ, otherElemSol, pvIdx); | |
391 | |||
392 | 66282876 | ElementResidualVector residual(numElementResiduals); | |
393 | 66282876 | residual = 0; | |
394 | |||
395 | 66282876 | evalSource(residual, scvI); | |
396 | |||
397 |
2/2✓ Branch 0 taken 35235240 times.
✓ Branch 1 taken 31047636 times.
|
66282876 | if (!this->assembler().isStationaryProblem()) |
398 | 35235240 | evalStorage(residual, scvI); | |
399 | |||
400 |
2/2✓ Branch 2 taken 200869544 times.
✓ Branch 3 taken 66282876 times.
|
468021964 | for (const auto& scvf : scvfs(fvGeometry, scvI)) |
401 | 200869544 | evalFlux(residual, scvf); | |
402 | |||
403 | 66282876 | return residual; | |
404 | }; | ||
405 | |||
406 | // derive the residuals numerically | ||
407 |
5/6✓ Branch 0 taken 75 times.
✓ Branch 1 taken 65531793 times.
✓ Branch 3 taken 71 times.
✓ Branch 4 taken 4 times.
✓ Branch 6 taken 71 times.
✗ Branch 7 not taken.
|
65531868 | static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()}; |
408 |
5/6✓ Branch 0 taken 75 times.
✓ Branch 1 taken 65531793 times.
✓ Branch 3 taken 71 times.
✓ Branch 4 taken 4 times.
✓ Branch 6 taken 71 times.
✗ Branch 7 not taken.
|
65531868 | static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod"); |
409 |
1/2✓ Branch 2 taken 222208 times.
✗ Branch 3 not taken.
|
65531868 | NumericDifferentiation::partialDerivative(evalResiduals, otherElemSol[scvJ.localDofIndex()][pvIdx], partialDerivs, origResiduals, |
410 | 65531868 | eps_(otherElemSol[scvJ.localDofIndex()][pvIdx], pvIdx), numDiffMethod); | |
411 | |||
412 | 196585140 | const auto updateJacobian = [&]() | |
413 | { | ||
414 |
2/2✓ Branch 0 taken 65526636 times.
✓ Branch 1 taken 65526636 times.
|
131053272 | for (int eqIdx = 0; eqIdx < numEq; eqIdx++) |
415 | { | ||
416 | // A[i][col][eqIdx][pvIdx] is the rate of change of | ||
417 | // the residual of equation 'eqIdx' at dof 'i' | ||
418 | // depending on the primary variable 'pvIdx' at dof | ||
419 | // 'col'. | ||
420 | 65526636 | A[scvI.dofIndex()][scvJ.dofIndex()][eqIdx][pvIdx] += partialDerivs[scvI.localDofIndex()][eqIdx]; | |
421 | } | ||
422 | }; | ||
423 | |||
424 |
2/2✓ Branch 1 taken 65524092 times.
✓ Branch 2 taken 7776 times.
|
65531868 | if (element.partitionType() == Dune::InteriorEntity) |
425 |
1/2✓ Branch 1 taken 222208 times.
✗ Branch 2 not taken.
|
65524092 | updateJacobian(); |
426 | else | ||
427 | { | ||
428 | 7776 | const auto localIdxI = scvI.indexInElement(); | |
429 |
2/5✓ Branch 1 taken 5232 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 0 taken 2544 times.
|
7776 | const auto& facetI = element.template subEntity <1> (localIdxI); |
430 | |||
431 |
2/2✓ Branch 0 taken 2544 times.
✓ Branch 1 taken 5232 times.
|
7776 | if (facetI.partitionType() == Dune::BorderEntity) |
432 |
0/2✗ Branch 1 not taken.
✗ Branch 2 not taken.
|
2544 | updateJacobian(); |
433 | ✗ | } | |
434 | |||
435 | // restore the original state of the scv's volume variables | ||
436 | 65531868 | curOtherVolVars = origOtherVolVars; | |
437 | |||
438 | // restore the original element solution | ||
439 | 65531868 | otherElemSol[scvJ.localDofIndex()][pvIdx] = curSol[scvJ.dofIndex()][pvIdx]; | |
440 | 65531868 | this->asImp_().maybeUpdateCouplingContext(scvJ, otherElemSol, pvIdx); | |
441 | // TODO additional dof dependencies | ||
442 | } | ||
443 | }; | ||
444 | |||
445 | // calculation of the derivatives | ||
446 |
3/3✓ Branch 1 taken 8148074 times.
✓ Branch 2 taken 2025761 times.
✓ Branch 0 taken 124144 times.
|
10297979 | for (const auto& scvI : scvs(fvGeometry)) |
447 | { | ||
448 | // derivative w.r.t. own DOFs | ||
449 | 8241182 | evalDerivative(scvI, scvI); | |
450 | |||
451 | // derivative w.r.t. other DOFs | ||
452 | 8241182 | const auto& otherScvIndices = fvGeometry.gridGeometry().connectivityMap()[scvI.index()]; | |
453 |
4/4✓ Branch 0 taken 853760 times.
✓ Branch 1 taken 56561070 times.
✓ Branch 2 taken 56436926 times.
✓ Branch 3 taken 8117038 times.
|
65531868 | for (const auto globalJ : otherScvIndices) |
454 | 57290686 | evalDerivative(scvI, fvGeometry.scv(globalJ)); | |
455 | } | ||
456 | |||
457 | // evaluate additional derivatives that might arise from the coupling | ||
458 | 6728 | this->asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables); | |
459 | |||
460 | 2056797 | return origResiduals; | |
461 | } | ||
462 | }; | ||
463 | |||
464 | |||
465 | /*! | ||
466 | * \ingroup Assembly | ||
467 | * \ingroup FaceCenteredStaggeredDiscretization | ||
468 | * \brief TODO docme | ||
469 | */ | ||
470 | template<class TypeTag, class Assembler, class Implementation> | ||
471 | class FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false, Implementation> | ||
472 | : public FaceCenteredLocalAssemblerBase< | ||
473 | TypeTag, Assembler, | ||
474 | Detail::NonVoidOrDefault_t<Implementation, FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false, Implementation>>, | ||
475 | /*implicit=*/false | ||
476 | > | ||
477 | { | ||
478 | using ThisType = FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false, Implementation>; | ||
479 | using ParentType = FaceCenteredLocalAssemblerBase<TypeTag, Assembler, Detail::NonVoidOrDefault_t<Implementation, ThisType>, false>; | ||
480 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
481 | using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity; | ||
482 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
483 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
484 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
485 | |||
486 | enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() }; | ||
487 | |||
488 | public: | ||
489 | using LocalResidual = typename ParentType::LocalResidual; | ||
490 | using ElementResidualVector = typename LocalResidual::ElementResidualVector; | ||
491 | using ParentType::ParentType; | ||
492 | |||
493 | /*! | ||
494 | * \brief Computes the derivatives with respect to the given element and adds them | ||
495 | * to the global matrix. | ||
496 | * | ||
497 | * \return The element residual at the current solution. | ||
498 | */ | ||
499 | template <class PartialReassembler = DefaultPartialReassembler> | ||
500 | ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables, | ||
501 | const PartialReassembler* partialReassembler = nullptr) | ||
502 | { | ||
503 | if (partialReassembler) | ||
504 | DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization"); | ||
505 | |||
506 | // get some aliases for convenience | ||
507 | const auto& problem = this->asImp_().problem(); | ||
508 | const auto& element = this->element(); | ||
509 | const auto& fvGeometry = this->fvGeometry(); | ||
510 | const auto& curSol = this->asImp_().curSol(); | ||
511 | auto&& curElemVolVars = this->curElemVolVars(); | ||
512 | |||
513 | // get the vector of the actual element residuals | ||
514 | const auto origResiduals = this->evalLocalResidual(); | ||
515 | const auto origStorageResiduals = this->evalLocalStorageResidual(); | ||
516 | |||
517 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
518 | // // | ||
519 | // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the // | ||
520 | // neighboring elements we do so by computing the derivatives of the fluxes which depend on the // | ||
521 | // actual element. In the actual element we evaluate the derivative of the entire residual. // | ||
522 | // // | ||
523 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
524 | |||
525 | // create the element solution | ||
526 | auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry()); | ||
527 | |||
528 | // create the vector storing the partial derivatives | ||
529 | ElementResidualVector partialDerivs(fvGeometry.numScv()); | ||
530 | |||
531 | // calculation of the derivatives | ||
532 | for (const auto& scv : scvs(fvGeometry)) | ||
533 | { | ||
534 | // dof index and corresponding actual pri vars | ||
535 | const auto dofIdx = scv.dofIndex(); | ||
536 | auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv); | ||
537 | const VolumeVariables origVolVars(curVolVars); | ||
538 | |||
539 | // calculate derivatives w.r.t to the privars at the dof at hand | ||
540 | for (int pvIdx = 0; pvIdx < numEq; pvIdx++) | ||
541 | { | ||
542 | partialDerivs = 0.0; | ||
543 | |||
544 | auto evalStorage = [&](Scalar priVar) | ||
545 | { | ||
546 | // auto partialDerivsTmp = partialDerivs; | ||
547 | elemSol[scv.localDofIndex()][pvIdx] = priVar; | ||
548 | curVolVars.update(elemSol, problem, element, scv); | ||
549 | return this->evalLocalStorageResidual(); | ||
550 | }; | ||
551 | |||
552 | // derive the residuals numerically | ||
553 | static const NumericEpsilon<Scalar, numEq> eps_{problem.paramGroup()}; | ||
554 | static const int numDiffMethod = getParamFromGroup<int>(problem.paramGroup(), "Assembly.NumericDifferenceMethod"); | ||
555 | NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origStorageResiduals, | ||
556 | eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod); | ||
557 | |||
558 | // update the global stiffness matrix with the current partial derivatives | ||
559 | for (int eqIdx = 0; eqIdx < numEq; eqIdx++) | ||
560 | { | ||
561 | // A[i][col][eqIdx][pvIdx] is the rate of change of | ||
562 | // the residual of equation 'eqIdx' at dof 'i' | ||
563 | // depending on the primary variable 'pvIdx' at dof | ||
564 | // 'col'. | ||
565 | A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx]; | ||
566 | } | ||
567 | |||
568 | // restore the original state of the scv's volume variables | ||
569 | curVolVars = origVolVars; | ||
570 | |||
571 | // restore the original element solution | ||
572 | elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx]; | ||
573 | } | ||
574 | } | ||
575 | return origResiduals; | ||
576 | } | ||
577 | }; | ||
578 | |||
579 | /*! | ||
580 | * \ingroup Assembly | ||
581 | * \ingroup FaceCenteredStaggeredDiscretization | ||
582 | * \brief TODO docme | ||
583 | */ | ||
584 | template<class TypeTag, class Assembler, class Implementation> | ||
585 | class FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true, Implementation> | ||
586 | : public FaceCenteredLocalAssemblerBase< | ||
587 | TypeTag, Assembler, | ||
588 | FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true, Implementation>, | ||
589 | /*implicit=*/true | ||
590 | > | ||
591 | { | ||
592 | using ThisType = FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true, Implementation>; | ||
593 | using ParentType = FaceCenteredLocalAssemblerBase<TypeTag, Assembler, Detail::NonVoidOrDefault_t<Implementation, ThisType>, true>; | ||
594 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
595 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
596 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
597 | |||
598 | enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() }; | ||
599 | |||
600 | public: | ||
601 | using LocalResidual = typename ParentType::LocalResidual; | ||
602 | using ElementResidualVector = typename LocalResidual::ElementResidualVector; | ||
603 | using ParentType::ParentType; | ||
604 | |||
605 | /*! | ||
606 | * \brief Computes the derivatives with respect to the given element and adds them | ||
607 | * to the global matrix. | ||
608 | * | ||
609 | * \return The element residual at the current solution. | ||
610 | */ | ||
611 | template <class PartialReassembler = DefaultPartialReassembler> | ||
612 | ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables, | ||
613 | const PartialReassembler* partialReassembler = nullptr) | ||
614 | { | ||
615 | if (partialReassembler) | ||
616 | DUNE_THROW(Dune::NotImplemented, "partial reassembly for analytic differentiation"); | ||
617 | |||
618 | // get some aliases for convenience | ||
619 | const auto& element = this->element(); | ||
620 | const auto& fvGeometry = this->fvGeometry(); | ||
621 | const auto& problem = this->asImp_().problem(); | ||
622 | const auto& curElemVolVars = this->curElemVolVars(); | ||
623 | const auto& elemFluxVarsCache = this->elemFluxVarsCache(); | ||
624 | |||
625 | // get the vecor of the actual element residuals | ||
626 | const auto origResiduals = this->evalLocalResidual(); | ||
627 | |||
628 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
629 | // // | ||
630 | // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the // | ||
631 | // neighboring elements we do so by computing the derivatives of the fluxes which depend on the // | ||
632 | // actual element. In the actual element we evaluate the derivative of the entire residual. // | ||
633 | // // | ||
634 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
635 | |||
636 | // calculation of the source and storage derivatives | ||
637 | for (const auto& scv : scvs(fvGeometry)) | ||
638 | { | ||
639 | // dof index and corresponding actual pri vars | ||
640 | const auto dofIdx = scv.dofIndex(); | ||
641 | const auto& volVars = curElemVolVars[scv]; | ||
642 | |||
643 | // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping) | ||
644 | // only if the problem is instationary we add derivative of storage term | ||
645 | // TODO if e.g. porosity depends on all dofs in the element, we would have off-diagonal matrix entries!? | ||
646 | if (!this->assembler().isStationaryProblem()) | ||
647 | this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx], | ||
648 | problem, | ||
649 | element, | ||
650 | fvGeometry, | ||
651 | volVars, | ||
652 | scv); | ||
653 | |||
654 | // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping) | ||
655 | // add source term derivatives | ||
656 | this->localResidual().addSourceDerivatives(A[dofIdx][dofIdx], | ||
657 | problem, | ||
658 | element, | ||
659 | fvGeometry, | ||
660 | volVars, | ||
661 | scv); | ||
662 | } | ||
663 | |||
664 | // localJacobian[scvIdx][otherScvIdx][eqIdx][priVarIdx] of the fluxes | ||
665 | for (const auto& scvf : scvfs(fvGeometry)) | ||
666 | { | ||
667 | if (!scvf.boundary()) | ||
668 | { | ||
669 | // add flux term derivatives | ||
670 | this->localResidual().addFluxDerivatives(A, | ||
671 | problem, | ||
672 | element, | ||
673 | fvGeometry, | ||
674 | curElemVolVars, | ||
675 | elemFluxVarsCache, | ||
676 | scvf); | ||
677 | } | ||
678 | |||
679 | // the boundary gets special treatment to simplify | ||
680 | // for the user | ||
681 | else | ||
682 | { | ||
683 | const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx()); | ||
684 | if (this->elemBcTypes()[insideScv.localDofIndex()].hasNeumann()) | ||
685 | { | ||
686 | // add flux term derivatives | ||
687 | this->localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()], | ||
688 | problem, | ||
689 | element, | ||
690 | fvGeometry, | ||
691 | curElemVolVars, | ||
692 | elemFluxVarsCache, | ||
693 | scvf); | ||
694 | } | ||
695 | } | ||
696 | } | ||
697 | |||
698 | return origResiduals; | ||
699 | } | ||
700 | }; | ||
701 | |||
702 | /*! | ||
703 | * \ingroup Assembly | ||
704 | * \ingroup FaceCenteredStaggeredDiscretization | ||
705 | * \brief TODO docme | ||
706 | */ | ||
707 | template<class TypeTag, class Assembler, class Implementation> | ||
708 | class FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/false, Implementation> | ||
709 | : public FaceCenteredLocalAssemblerBase< | ||
710 | TypeTag, Assembler, | ||
711 | FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false, Implementation>, | ||
712 | /*implicit=*/false | ||
713 | > | ||
714 | { | ||
715 | using ThisType = FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false, Implementation>; | ||
716 | using ParentType = FaceCenteredLocalAssemblerBase<TypeTag, Assembler, Detail::NonVoidOrDefault_t<Implementation, ThisType>, false>; | ||
717 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
718 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
719 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
720 | |||
721 | enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() }; | ||
722 | |||
723 | public: | ||
724 | using LocalResidual = typename ParentType::LocalResidual; | ||
725 | using ElementResidualVector = typename LocalResidual::ElementResidualVector; | ||
726 | using ParentType::ParentType; | ||
727 | |||
728 | /*! | ||
729 | * \brief Computes the derivatives with respect to the given element and adds them | ||
730 | * to the global matrix. | ||
731 | * | ||
732 | * \return The element residual at the current solution. | ||
733 | */ | ||
734 | template <class PartialReassembler = DefaultPartialReassembler> | ||
735 | ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables, | ||
736 | const PartialReassembler* partialReassembler = nullptr) | ||
737 | { | ||
738 | if (partialReassembler) | ||
739 | DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization"); | ||
740 | |||
741 | // get some aliases for convenience | ||
742 | const auto& element = this->element(); | ||
743 | const auto& fvGeometry = this->fvGeometry(); | ||
744 | const auto& problem = this->asImp_().problem(); | ||
745 | const auto& curElemVolVars = this->curElemVolVars(); | ||
746 | |||
747 | // get the vector of the actual element residuals | ||
748 | const auto origResiduals = this->evalLocalResidual(); | ||
749 | |||
750 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
751 | // // | ||
752 | // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the // | ||
753 | // neighboring elements we do so by computing the derivatives of the fluxes which depend on the // | ||
754 | // actual element. In the actual element we evaluate the derivative of the entire residual. // | ||
755 | // // | ||
756 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
757 | |||
758 | // calculation of the source and storage derivatives | ||
759 | for (const auto& scv : scvs(fvGeometry)) | ||
760 | { | ||
761 | // dof index and corresponding actual pri vars | ||
762 | const auto dofIdx = scv.dofIndex(); | ||
763 | const auto& volVars = curElemVolVars[scv]; | ||
764 | |||
765 | // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping) | ||
766 | // only if the problem is instationary we add derivative of storage term | ||
767 | this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx], | ||
768 | problem, | ||
769 | element, | ||
770 | fvGeometry, | ||
771 | volVars, | ||
772 | scv); | ||
773 | } | ||
774 | |||
775 | return origResiduals; | ||
776 | } | ||
777 | }; | ||
778 | |||
779 | } // end namespace Dumux | ||
780 | |||
781 | #endif | ||
782 |