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 CCDiscretization | ||
11 | * \ingroup MultiDomain | ||
12 | * \brief A multidomain local assembler for Jacobian and residual contribution per element (cell-centered methods) | ||
13 | */ | ||
14 | #ifndef DUMUX_MULTIDOMAIN_CC_LOCAL_ASSEMBLER_HH | ||
15 | #define DUMUX_MULTIDOMAIN_CC_LOCAL_ASSEMBLER_HH | ||
16 | |||
17 | #include <dune/common/indices.hh> | ||
18 | #include <dune/common/reservedvector.hh> | ||
19 | #include <dune/common/hybridutilities.hh> | ||
20 | #include <dune/grid/common/gridenums.hh> // for GhostEntity | ||
21 | #include <dune/istl/matrixindexset.hh> | ||
22 | |||
23 | #include <dumux/common/reservedblockvector.hh> | ||
24 | #include <dumux/common/properties.hh> | ||
25 | #include <dumux/common/parameters.hh> | ||
26 | #include <dumux/common/numericdifferentiation.hh> | ||
27 | #include <dumux/common/numeqvector.hh> | ||
28 | #include <dumux/discretization/elementsolution.hh> | ||
29 | #include <dumux/discretization/extrusion.hh> | ||
30 | #include <dumux/assembly/numericepsilon.hh> | ||
31 | #include <dumux/assembly/diffmethod.hh> | ||
32 | #include <dumux/assembly/fvlocalassemblerbase.hh> | ||
33 | |||
34 | namespace Dumux { | ||
35 | |||
36 | /*! | ||
37 | * \ingroup Assembly | ||
38 | * \ingroup CCDiscretization | ||
39 | * \ingroup MultiDomain | ||
40 | * \brief A base class for all multidomain local assemblers | ||
41 | * \tparam id the id of the sub domain | ||
42 | * \tparam TypeTag the TypeTag | ||
43 | * \tparam Assembler the assembler type | ||
44 | * \tparam Implementation the actual assembler implementation | ||
45 | * \tparam implicit Specifies whether the time discretization is implicit or not not (i.e. explicit) | ||
46 | */ | ||
47 | template<std::size_t id, class TypeTag, class Assembler, class Implementation, bool implicit> | ||
48 | 35520443 | class SubDomainCCLocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler,Implementation, implicit> | |
49 | { | ||
50 | using ParentType = FVLocalAssemblerBase<TypeTag, Assembler,Implementation, implicit>; | ||
51 | |||
52 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
53 | using LocalResidualValues = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; | ||
54 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
55 | using SolutionVector = typename Assembler::SolutionVector; | ||
56 | using ElementBoundaryTypes = GetPropType<TypeTag, Properties::ElementBoundaryTypes>; | ||
57 | |||
58 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
59 | using GridVolumeVariables = typename GridVariables::GridVolumeVariables; | ||
60 | using ElementVolumeVariables = typename GridVolumeVariables::LocalView; | ||
61 | using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView; | ||
62 | using Scalar = typename GridVariables::Scalar; | ||
63 | |||
64 | using GridGeometry = typename GridVariables::GridGeometry; | ||
65 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
66 | using SubControlVolume = typename GridGeometry::SubControlVolume; | ||
67 | using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace; | ||
68 | using Extrusion = Extrusion_t<GridGeometry>; | ||
69 | using GridView = typename GridGeometry::GridView; | ||
70 | using Element = typename GridView::template Codim<0>::Entity; | ||
71 | |||
72 | using CouplingManager = typename Assembler::CouplingManager; | ||
73 | using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector; | ||
74 | |||
75 | public: | ||
76 | //! export the domain id of this sub-domain | ||
77 | static constexpr auto domainId = typename Dune::index_constant<id>(); | ||
78 | //! the local residual type of this domain | ||
79 | using LocalResidual = typename ParentType::LocalResidual; | ||
80 | //! pull up constructor of parent class | ||
81 | using ParentType::ParentType; | ||
82 | |||
83 | //! the constructor | ||
84 | 36009438 | explicit SubDomainCCLocalAssemblerBase(const Assembler& assembler, | |
85 | const Element& element, | ||
86 | const SolutionVector& curSol, | ||
87 | CouplingManager& couplingManager) | ||
88 | : ParentType(assembler, | ||
89 | element, | ||
90 | curSol, | ||
91 | localView(assembler.gridGeometry(domainId)), | ||
92 |
1/2✓ Branch 1 taken 19500881 times.
✗ Branch 2 not taken.
|
36009438 | localView(assembler.gridVariables(domainId).curGridVolVars()), |
93 |
1/2✓ Branch 1 taken 19500881 times.
✗ Branch 2 not taken.
|
36009438 | localView(assembler.gridVariables(domainId).prevGridVolVars()), |
94 |
3/6✓ Branch 1 taken 19500881 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 16023098 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 16023098 times.
✗ Branch 8 not taken.
|
94966766 | localView(assembler.gridVariables(domainId).gridFluxVarsCache()), |
95 | assembler.localResidual(domainId), | ||
96 | 68297550 | (element.partitionType() == Dune::GhostEntity)) | |
97 |
30/33✓ Branch 0 taken 419186 times.
✓ Branch 1 taken 1886765 times.
✓ Branch 2 taken 17194930 times.
✓ Branch 3 taken 419186 times.
✓ Branch 4 taken 1886765 times.
✓ Branch 5 taken 17194930 times.
✓ Branch 6 taken 419186 times.
✓ Branch 7 taken 1886765 times.
✓ Branch 8 taken 17194930 times.
✓ Branch 9 taken 419186 times.
✓ Branch 10 taken 1886765 times.
✓ Branch 11 taken 17194930 times.
✓ Branch 12 taken 419186 times.
✓ Branch 13 taken 1886765 times.
✓ Branch 14 taken 17194930 times.
✓ Branch 15 taken 419186 times.
✓ Branch 16 taken 1886765 times.
✓ Branch 17 taken 17194930 times.
✓ Branch 18 taken 419186 times.
✓ Branch 19 taken 1886765 times.
✓ Branch 20 taken 17194930 times.
✓ Branch 21 taken 419186 times.
✓ Branch 22 taken 1886765 times.
✓ Branch 23 taken 17194930 times.
✓ Branch 24 taken 419186 times.
✓ Branch 25 taken 1886765 times.
✓ Branch 26 taken 3372993 times.
✗ Branch 27 not taken.
✓ Branch 28 taken 108326 times.
✓ Branch 29 taken 3369457 times.
✗ Branch 30 not taken.
✓ Branch 31 taken 108326 times.
✗ Branch 32 not taken.
|
70436694 | , couplingManager_(couplingManager) |
98 | 36009438 | {} | |
99 | |||
100 | /*! | ||
101 | * \brief Computes the derivatives with respect to the given element and adds them | ||
102 | * to the global matrix. The element residual is written into the right hand side. | ||
103 | */ | ||
104 | template<class JacobianMatrixRow, class SubResidualVector, class GridVariablesTuple> | ||
105 | 21940231 | void assembleJacobianAndResidual(JacobianMatrixRow& jacRow, SubResidualVector& res, GridVariablesTuple& gridVariables) | |
106 | { | ||
107 | 21940231 | this->asImp_().bindLocalViews(); | |
108 | |||
109 | // for the diagonal jacobian block | ||
110 | 65820693 | const auto globalI = this->fvGeometry().gridGeometry().elementMapper().index(this->element()); | |
111 |
1/2✓ Branch 5 taken 1142879 times.
✗ Branch 6 not taken.
|
87760924 | res[globalI] = this->asImp_().assembleJacobianAndResidualImplInverse(jacRow[domainId], *std::get<domainId>(gridVariables)); |
112 | |||
113 | // for the coupling blocks | ||
114 | using namespace Dune::Hybrid; | ||
115 |
1/2✓ Branch 1 taken 1142879 times.
✗ Branch 2 not taken.
|
39671153 | forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](auto&& i) |
116 | { | ||
117 | if (i != id) | ||
118 | 19262163 | this->assembleJacobianCoupling(i, jacRow, res[globalI], gridVariables); | |
119 | }); | ||
120 | 21940231 | } | |
121 | |||
122 | /*! | ||
123 | * \brief Assemble the entries in a coupling block of the jacobian. | ||
124 | * There is no coupling block between a domain and itself. | ||
125 | */ | ||
126 | template<std::size_t otherId, class JacRow, class GridVariables, | ||
127 | typename std::enable_if_t<(otherId == id), int> = 0> | ||
128 | void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow, | ||
129 | const LocalResidualValues& res, GridVariables& gridVariables) | ||
130 | {} | ||
131 | |||
132 | /*! | ||
133 | * \brief Assemble the entries in a coupling block of the jacobian. | ||
134 | */ | ||
135 | template<std::size_t otherId, class JacRow, class GridVariables, | ||
136 | typename std::enable_if_t<(otherId != id), int> = 0> | ||
137 | ✗ | void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow, | |
138 | const LocalResidualValues& res, GridVariables& gridVariables) | ||
139 | { | ||
140 | 44600652 | this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables)); | |
141 | ✗ | } | |
142 | |||
143 | /*! | ||
144 | * \brief Assemble the residual only | ||
145 | */ | ||
146 | template<class SubResidualVector> | ||
147 | 14069207 | void assembleResidual(SubResidualVector& res) | |
148 | { | ||
149 | 14069207 | this->asImp_().bindLocalViews(); | |
150 | 42207621 | const auto globalI = this->fvGeometry().gridGeometry().elementMapper().index(this->element()); | |
151 |
0/2✗ Branch 0 not taken.
✗ Branch 1 not taken.
|
14069207 | res[globalI] = this->evalLocalResidual()[0]; // forward to the internal implementation |
152 | |||
153 | if constexpr (Problem::enableInternalDirichletConstraints()) | ||
154 | { | ||
155 | 705000 | const auto applyDirichlet = [&] (const auto& scvI, | |
156 | const auto& dirichletValues, | ||
157 | const auto eqIdx, | ||
158 | ✗ | const auto pvIdx) | |
159 | { | ||
160 | 220 | res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx]; | |
161 | }; | ||
162 | |||
163 | 705000 | this->asImp_().enforceInternalDirichletConstraints(applyDirichlet); | |
164 | } | ||
165 | 14069205 | } | |
166 | |||
167 | /*! | ||
168 | * \brief Evaluates the local source term for an element and given element volume variables | ||
169 | */ | ||
170 | 1640806 | ElementResidualVector evalLocalSourceResidual(const Element& element, const ElementVolumeVariables& elemVolVars) const | |
171 | { | ||
172 | // initialize the residual vector for all scvs in this element | ||
173 | 4922418 | ElementResidualVector residual(this->fvGeometry().numScv()); | |
174 | |||
175 | // evaluate the volume terms (storage + source terms) | ||
176 | // forward to the local residual specialized for the discretization methods | ||
177 |
4/4✓ Branch 0 taken 821243 times.
✓ Branch 1 taken 821243 times.
✓ Branch 2 taken 821243 times.
✓ Branch 3 taken 821243 times.
|
6563224 | for (auto&& scv : scvs(this->fvGeometry())) |
178 | { | ||
179 | 1640806 | const auto& curVolVars = elemVolVars[scv]; | |
180 |
1/2✓ Branch 0 taken 76979 times.
✗ Branch 1 not taken.
|
1640806 | auto source = this->localResidual().computeSource(problem(), element, this->fvGeometry(), elemVolVars, scv); |
181 | 3281612 | source *= -Extrusion::volume(this->fvGeometry(), scv)*curVolVars.extrusionFactor(); | |
182 | 3281612 | residual[scv.indexInElement()] = std::move(source); | |
183 | } | ||
184 | |||
185 | 1640806 | return residual; | |
186 | } | ||
187 | |||
188 | /*! | ||
189 | * \brief Evaluates the local source term depending on time discretization scheme | ||
190 | */ | ||
191 | ElementResidualVector evalLocalSourceResidual(const Element& neighbor) const | ||
192 | 1642486 | { return this->evalLocalSourceResidual(neighbor, implicit ? this->curElemVolVars() : this->prevElemVolVars()); } | |
193 | |||
194 | /*! | ||
195 | * \brief Evaluates the storage terms within the element | ||
196 | */ | ||
197 | LocalResidualValues evalLocalStorageResidual() const | ||
198 | { | ||
199 | 12168000 | return this->localResidual().evalStorage(this->element(), this->fvGeometry(), this->prevElemVolVars(), this->curElemVolVars())[0]; | |
200 | } | ||
201 | |||
202 | /*! | ||
203 | * \brief Evaluates the fluxes depending on the chose time discretization scheme | ||
204 | */ | ||
205 | LocalResidualValues evalFluxResidual(const Element& neighbor, | ||
206 | const SubControlVolumeFace& scvf) const | ||
207 | { | ||
208 | 225045848 | const auto& elemVolVars = implicit ? this->curElemVolVars() : this->prevElemVolVars(); | |
209 | 337568772 | return this->localResidual().evalFlux(problem(), neighbor, this->fvGeometry(), elemVolVars, this->elemFluxVarsCache(), scvf); | |
210 | } | ||
211 | |||
212 | /*! | ||
213 | * \brief Prepares all local views necessary for local assembly. | ||
214 | */ | ||
215 | 36009438 | void bindLocalViews() | |
216 | { | ||
217 | // get some references for convenience | ||
218 | 36009438 | const auto& element = this->element(); | |
219 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 173500 times.
|
36009438 | const auto& curSol = this->curSol()[domainId]; |
220 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 173500 times.
|
36009438 | auto&& fvGeometry = this->fvGeometry(); |
221 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 173500 times.
|
36009438 | auto&& curElemVolVars = this->curElemVolVars(); |
222 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 173500 times.
|
36009438 | auto&& elemFluxVarsCache = this->elemFluxVarsCache(); |
223 | |||
224 | // bind the caches | ||
225 | 36009438 | couplingManager_.bindCouplingContext(domainId, element, this->assembler()); | |
226 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 173500 times.
|
36009438 | fvGeometry.bind(element); |
227 | |||
228 | if (implicit) | ||
229 | { | ||
230 | 31953438 | curElemVolVars.bind(element, fvGeometry, curSol); | |
231 | 31953438 | elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars); | |
232 |
4/4✓ Branch 0 taken 621328 times.
✓ Branch 1 taken 828455 times.
✓ Branch 2 taken 621328 times.
✓ Branch 3 taken 828455 times.
|
34428212 | if (!this->assembler().isStationaryProblem()) |
233 | 4234000 | this->prevElemVolVars().bindElement(element, fvGeometry, this->assembler().prevSol()[domainId]); | |
234 | } | ||
235 | else | ||
236 | { | ||
237 | 4056000 | auto& prevElemVolVars = this->prevElemVolVars(); | |
238 | 8112000 | const auto& prevSol = this->assembler().prevSol()[domainId]; | |
239 | |||
240 | 4056000 | curElemVolVars.bindElement(element, fvGeometry, curSol); | |
241 | 4056000 | prevElemVolVars.bind(element, fvGeometry, prevSol); | |
242 | 4056000 | elemFluxVarsCache.bind(element, fvGeometry, prevElemVolVars); | |
243 | } | ||
244 | 36009438 | } | |
245 | |||
246 | //! return reference to the underlying problem | ||
247 | const Problem& problem() const | ||
248 |
9/17✓ Branch 1 taken 1680 times.
✓ Branch 2 taken 705000 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 75300 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 6 times.
✓ Branch 13 taken 13 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 6 times.
✓ Branch 16 taken 12 times.
✗ Branch 17 not taken.
|
114049317 | { return this->assembler().problem(domainId); } |
249 | |||
250 | //! return reference to the coupling manager | ||
251 | ✗ | CouplingManager& couplingManager() | |
252 | ✗ | { return couplingManager_; } | |
253 | |||
254 | private: | ||
255 | CouplingManager& couplingManager_; //!< the coupling manager | ||
256 | }; | ||
257 | |||
258 | /*! | ||
259 | * \ingroup Assembly | ||
260 | * \ingroup CCDiscretization | ||
261 | * \ingroup MultiDomain | ||
262 | * \brief The cell-centered scheme multidomain local assembler | ||
263 | * \tparam id the id of the sub domain | ||
264 | * \tparam TypeTag the TypeTag | ||
265 | * \tparam Assembler the assembler type | ||
266 | * \tparam DM the numeric differentiation method | ||
267 | * \tparam implicit whether the assembler is explicit or implicit in time | ||
268 | */ | ||
269 | template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric, bool implicit = true> | ||
270 | class SubDomainCCLocalAssembler; | ||
271 | |||
272 | /*! | ||
273 | * \ingroup Assembly | ||
274 | * \ingroup CCDiscretization | ||
275 | * \ingroup MultiDomain | ||
276 | * \brief Cell-centered scheme multidomain local assembler using numeric differentiation and implicit time discretization | ||
277 | */ | ||
278 | template<std::size_t id, class TypeTag, class Assembler> | ||
279 | 33492442 | class SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true> | |
280 | : public SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler, | ||
281 | SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, true>, true> | ||
282 | { | ||
283 | using ThisType = SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>; | ||
284 | using ParentType = SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler, ThisType, /*implicit=*/true>; | ||
285 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
286 | |||
287 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
288 | using LocalResidualValues = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; | ||
289 | |||
290 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
291 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
292 | using GridView = typename GridGeometry::GridView; | ||
293 | using Element = typename GridView::template Codim<0>::Entity; | ||
294 | |||
295 | enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() }; | ||
296 | enum { dim = GridView::dimension }; | ||
297 | |||
298 | static constexpr bool enableGridFluxVarsCache = GetPropType<TypeTag, Properties::GridVariables>::GridFluxVariablesCache::cachingEnabled; | ||
299 | static constexpr bool enableGridVolVarsCache = GetPropType<TypeTag, Properties::GridVariables>::GridVolumeVariables::cachingEnabled; | ||
300 | static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize; | ||
301 | static constexpr auto domainI = Dune::index_constant<id>(); | ||
302 | |||
303 | public: | ||
304 | 17472881 | using ParentType::ParentType; | |
305 | |||
306 | /*! | ||
307 | * \brief Computes the derivatives with respect to the given element and adds them | ||
308 | * to the global matrix. | ||
309 | * | ||
310 | * \return The element residual at the current solution. | ||
311 | */ | ||
312 | template<class JacobianMatrixDiagBlock, class GridVariables> | ||
313 | 17884231 | LocalResidualValues assembleJacobianAndResidualImplInverse(JacobianMatrixDiagBlock& A, GridVariables& gridVariables) | |
314 | { | ||
315 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
316 | // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the // | ||
317 | // neighboring elements we do so by computing the derivatives of the fluxes which depend on the // | ||
318 | // actual element. In the actual element we evaluate the derivative of the entire residual. // | ||
319 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
320 | |||
321 | // get some aliases for convenience | ||
322 | 17884231 | const auto& element = this->element(); | |
323 | 17884231 | const auto& fvGeometry = this->fvGeometry(); | |
324 | 17884231 | const auto& gridGeometry = fvGeometry.gridGeometry(); | |
325 | 17884231 | auto&& curElemVolVars = this->curElemVolVars(); | |
326 | 17884231 | auto&& elemFluxVarsCache = this->elemFluxVarsCache(); | |
327 | |||
328 | // get stencil information | ||
329 | 35768462 | const auto globalI = gridGeometry.elementMapper().index(element); | |
330 | 17884231 | const auto& connectivityMap = gridGeometry.connectivityMap(); | |
331 | 35768462 | const auto numNeighbors = connectivityMap[globalI].size(); | |
332 | |||
333 | // container to store the neighboring elements | ||
334 | 17884231 | Dune::ReservedVector<Element, maxElementStencilSize> neighborElements; | |
335 | 17884231 | neighborElements.resize(numNeighbors); | |
336 | |||
337 | // assemble the undeflected residual | ||
338 | using Residuals = ReservedBlockVector<LocalResidualValues, maxElementStencilSize>; | ||
339 | 35768462 | Residuals origResiduals(numNeighbors + 1); origResiduals = 0.0; | |
340 |
1/3✗ Branch 0 not taken.
✓ Branch 1 taken 1125317 times.
✗ Branch 2 not taken.
|
20011079 | origResiduals[0] = this->evalLocalResidual()[0]; |
341 | |||
342 | // lambda for convenient evaluation of the fluxes across scvfs in the neighbors | ||
343 | // if the neighbor is a ghost we don't want to add anything to their residual | ||
344 | // so we return 0 and omit computing the flux | ||
345 | 168614351 | const auto evalNeighborFlux = [&] (const auto& neighbor, const auto& scvf) | |
346 | { | ||
347 |
4/12✓ Branch 0 taken 959064 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 50616388 times.
✓ Branch 3 taken 53295168 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 3876268 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
112511484 | if (neighbor.partitionType() == Dune::GhostEntity) |
348 | ✗ | return LocalResidualValues(0.0); | |
349 | else | ||
350 | 112528644 | return this->evalFluxResidual(neighbor, scvf); | |
351 | }; | ||
352 | |||
353 | // get the elements in which we need to evaluate the fluxes | ||
354 | // and calculate these in the undeflected state | ||
355 | 17884231 | unsigned int j = 1; | |
356 |
4/4✓ Branch 0 taken 48534620 times.
✓ Branch 1 taken 10007284 times.
✓ Branch 2 taken 48534620 times.
✓ Branch 3 taken 10007284 times.
|
249373268 | for (const auto& dataJ : connectivityMap[globalI]) |
357 | { | ||
358 |
1/2✓ Branch 1 taken 2684894 times.
✗ Branch 2 not taken.
|
95981760 | neighborElements[j-1] = gridGeometry.element(dataJ.globalJ); |
359 |
4/4✓ Branch 0 taken 51651416 times.
✓ Branch 1 taken 48534620 times.
✓ Branch 2 taken 4688708 times.
✓ Branch 3 taken 1571912 times.
|
371283696 | for (const auto scvfIdx : dataJ.scvfsJ) |
360 |
5/8✓ Branch 1 taken 13632 times.
✓ Branch 2 taken 5788058 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 13632 times.
✓ Branch 5 taken 5788058 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 13632 times.
✗ Branch 8 not taken.
|
261837784 | origResiduals[j] += evalNeighborFlux(neighborElements[j-1], fvGeometry.scvf(scvfIdx)); |
361 | |||
362 | 88918172 | ++j; | |
363 | } | ||
364 | |||
365 | // reference to the element's scv (needed later) and corresponding vol vars | ||
366 |
1/2✓ Branch 0 taken 1449783 times.
✗ Branch 1 not taken.
|
17884231 | const auto& scv = fvGeometry.scv(globalI); |
367 | 35768462 | auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv); | |
368 | |||
369 | // save a copy of the original privars and vol vars in order | ||
370 | // to restore the original solution after deflection | ||
371 | 17884231 | const auto& curSol = this->curSol()[domainI]; | |
372 | 17884231 | const auto origPriVars = curSol[globalI]; | |
373 | 17884231 | const auto origVolVars = curVolVars; | |
374 | |||
375 | // element solution container to be deflected | ||
376 | 17884231 | auto elemSol = elementSolution<FVElementGeometry>(origPriVars); | |
377 | |||
378 | // derivatives in the neighbors with respect to the current elements | ||
379 | // in index 0 we save the derivative of the element residual with respect to it's own dofs | ||
380 | 17884231 | Residuals partialDerivs(numNeighbors + 1); | |
381 | |||
382 |
2/2✓ Branch 0 taken 11298288 times.
✓ Branch 1 taken 10007284 times.
|
37704758 | for (int pvIdx = 0; pvIdx < numEq; ++pvIdx) |
383 | { | ||
384 | 19820527 | partialDerivs = 0.0; | |
385 | |||
386 | 65773631 | auto evalResiduals = [&](Scalar priVar) | |
387 | { | ||
388 | 11298288 | Residuals partialDerivsTmp(numNeighbors + 1); | |
389 | 11298288 | partialDerivsTmp = 0.0; | |
390 | // update the volume variables and the flux var cache | ||
391 | 13558332 | elemSol[0][pvIdx] = priVar; | |
392 | 22596576 | this->couplingManager().updateCouplingContext(domainI, *this, domainI, globalI, elemSol[0], pvIdx); | |
393 | 22596576 | curVolVars.update(elemSol, this->problem(), element, scv); | |
394 | 62837981 | elemFluxVarsCache.update(element, fvGeometry, curElemVolVars); | |
395 | if (enableGridFluxVarsCache) | ||
396 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 118184 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 217464 times.
|
9938725 | gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars); |
397 | |||
398 | // calculate the residual with the deflected primary variables | ||
399 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 1053360 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 435864 times.
|
14276736 | partialDerivsTmp[0] = this->evalLocalResidual()[0]; |
400 | |||
401 | // calculate the fluxes in the neighbors with the deflected primary variables | ||
402 |
18/18✓ Branch 0 taken 14255232 times.
✓ Branch 1 taken 4111561 times.
✓ Branch 2 taken 14255232 times.
✓ Branch 3 taken 4111561 times.
✓ Branch 4 taken 14255232 times.
✓ Branch 5 taken 4111561 times.
✓ Branch 6 taken 40156100 times.
✓ Branch 7 taken 7180946 times.
✓ Branch 8 taken 40156100 times.
✓ Branch 9 taken 7180946 times.
✓ Branch 10 taken 40156100 times.
✓ Branch 11 taken 7180946 times.
✓ Branch 12 taken 20892 times.
✓ Branch 13 taken 5781 times.
✓ Branch 14 taken 20892 times.
✓ Branch 15 taken 5781 times.
✓ Branch 16 taken 20892 times.
✓ Branch 17 taken 5781 times.
|
65730512 | for (std::size_t k = 0; k < connectivityMap[globalI].size(); ++k) |
403 |
6/6✓ Branch 0 taken 14255232 times.
✓ Branch 1 taken 14255232 times.
✓ Branch 2 taken 46595384 times.
✓ Branch 3 taken 40156100 times.
✓ Branch 4 taken 9677784 times.
✓ Branch 5 taken 3238500 times.
|
342689520 | for (auto scvfIdx : connectivityMap[globalI][k].scvfsJ) |
404 | 165885732 | partialDerivsTmp[k+1] += evalNeighborFlux(neighborElements[k], fvGeometry.scvf(scvfIdx)); | |
405 | |||
406 | 11298288 | return partialDerivsTmp; | |
407 | }; | ||
408 | |||
409 | // derive the residuals numerically | ||
410 |
7/10✓ Branch 0 taken 345 times.
✓ Branch 1 taken 11297943 times.
✓ Branch 3 taken 215 times.
✓ Branch 4 taken 130 times.
✓ Branch 6 taken 215 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 215 times.
✗ Branch 10 not taken.
✓ Branch 12 taken 215 times.
✗ Branch 13 not taken.
|
19820527 | static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod"); |
411 |
7/10✓ Branch 0 taken 474 times.
✓ Branch 1 taken 11297814 times.
✓ Branch 3 taken 215 times.
✓ Branch 4 taken 259 times.
✓ Branch 6 taken 215 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 215 times.
✗ Branch 10 not taken.
✓ Branch 12 taken 215 times.
✗ Branch 13 not taken.
|
19820527 | static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()}; |
412 |
2/4✓ Branch 1 taken 11298288 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 2308812 times.
✗ Branch 5 not taken.
|
23256683 | NumericDifferentiation::partialDerivative(evalResiduals, elemSol[0][pvIdx], partialDerivs, origResiduals, |
413 | 23256683 | eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod); | |
414 | |||
415 | // Correct derivative for ghost elements, i.e. set a 1 for the derivative w.r.t. the | ||
416 | // current primary variable and a 0 elsewhere. As we always solve for a delta of the | ||
417 | // solution with respect to the initial one, this results in a delta of 0 for ghosts. | ||
418 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 11298288 times.
|
19820527 | if (this->elementIsGhost()) |
419 | { | ||
420 | ✗ | partialDerivs[0] = 0.0; | |
421 | ✗ | partialDerivs[0][pvIdx] = 1.0; | |
422 | } | ||
423 | |||
424 | // restore the original state of the scv's volume variables | ||
425 | 19820527 | curVolVars = origVolVars; | |
426 | |||
427 | // restore the current element solution | ||
428 |
5/6✓ Branch 0 taken 167 times.
✓ Branch 1 taken 244773 times.
✓ Branch 2 taken 167 times.
✓ Branch 3 taken 231333 times.
✓ Branch 4 taken 13440 times.
✗ Branch 5 not taken.
|
23256683 | elemSol[0][pvIdx] = origPriVars[pvIdx]; |
429 | |||
430 | // restore the undeflected state of the coupling context | ||
431 |
2/3✓ Branch 0 taken 705216 times.
✓ Branch 1 taken 1199012 times.
✗ Branch 2 not taken.
|
19820527 | this->couplingManager().updateCouplingContext(domainI, *this, domainI, globalI, elemSol[0], pvIdx); |
432 | |||
433 | // add the current partial derivatives to the global jacobian matrix | ||
434 | if constexpr (Problem::enableInternalDirichletConstraints()) | ||
435 | { | ||
436 | // check if own element has internal Dirichlet constraint | ||
437 |
5/6✓ Branch 0 taken 705216 times.
✓ Branch 1 taken 439060 times.
✓ Branch 2 taken 705216 times.
✓ Branch 3 taken 431988 times.
✓ Branch 4 taken 7072 times.
✗ Branch 5 not taken.
|
3216072 | const auto internalDirichletConstraintsOwnElement = this->problem().hasInternalDirichletConstraint(this->element(), scv); |
438 | 3216272 | const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv); | |
439 | |||
440 |
3/3✓ Branch 0 taken 2473336 times.
✓ Branch 1 taken 1608136 times.
✓ Branch 2 taken 200 times.
|
4081672 | for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) |
441 | { | ||
442 |
4/4✓ Branch 0 taken 6560 times.
✓ Branch 1 taken 2466976 times.
✓ Branch 2 taken 6560 times.
✓ Branch 3 taken 1482976 times.
|
3963072 | if (internalDirichletConstraintsOwnElement[eqIdx]) |
443 | { | ||
444 |
8/8✓ Branch 0 taken 92 times.
✓ Branch 1 taken 6468 times.
✓ Branch 2 taken 92 times.
✓ Branch 3 taken 75 times.
✓ Branch 4 taken 92 times.
✓ Branch 5 taken 75 times.
✓ Branch 6 taken 92 times.
✓ Branch 7 taken 75 times.
|
7061 | origResiduals[0][eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx]; |
445 |
4/6✓ Branch 0 taken 92 times.
✓ Branch 1 taken 6468 times.
✓ Branch 3 taken 6560 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 6560 times.
✗ Branch 7 not taken.
|
6744 | A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0; |
446 | } | ||
447 | else | ||
448 |
4/8✓ Branch 1 taken 2466976 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 2466976 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 2466976 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1492733 times.
✗ Branch 11 not taken.
|
8893661 | A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx]; |
449 | } | ||
450 | |||
451 | // off-diagonal entries | ||
452 | 1608136 | j = 1; | |
453 |
4/4✓ Branch 0 taken 6293486 times.
✓ Branch 1 taken 1608136 times.
✓ Branch 2 taken 6293486 times.
✓ Branch 3 taken 1608136 times.
|
19019516 | for (const auto& dataJ : connectivityMap[globalI]) |
454 | { | ||
455 |
2/3✓ Branch 0 taken 3006232 times.
✓ Branch 1 taken 1718424 times.
✗ Branch 2 not taken.
|
4990326 | const auto& neighborElement = neighborElements[j-1]; |
456 |
2/3✓ Branch 0 taken 3006232 times.
✓ Branch 1 taken 1718424 times.
✗ Branch 2 not taken.
|
4990326 | const auto& neighborScv = fvGeometry.scv(dataJ.globalJ); |
457 |
5/6✓ Branch 0 taken 2792832 times.
✓ Branch 1 taken 1718424 times.
✓ Branch 2 taken 2792832 times.
✓ Branch 3 taken 1691160 times.
✓ Branch 4 taken 27264 times.
✗ Branch 5 not taken.
|
9714982 | const auto internalDirichletConstraintsNeighbor = this->problem().hasInternalDirichletConstraint(neighborElement, neighborScv); |
458 | |||
459 |
3/3✓ Branch 0 taken 9359496 times.
✓ Branch 1 taken 6293486 times.
✓ Branch 2 taken 265670 times.
|
15918652 | for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) |
460 | { | ||
461 |
4/4✓ Branch 0 taken 19130 times.
✓ Branch 1 taken 9606036 times.
✓ Branch 2 taken 19130 times.
✓ Branch 3 taken 5836596 times.
|
15480892 | if (internalDirichletConstraintsNeighbor[eqIdx]) |
462 |
2/4✓ Branch 1 taken 19130 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 19130 times.
✗ Branch 5 not taken.
|
38260 | A[dataJ.globalJ][globalI][eqIdx][pvIdx] = 0.0; |
463 | else | ||
464 |
4/8✓ Branch 1 taken 9606036 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 9606036 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 9606036 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 5749226 times.
✗ Branch 11 not taken.
|
34567334 | A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j][eqIdx]; |
465 | } | ||
466 | |||
467 | 6293486 | ++j; | |
468 | } | ||
469 | } | ||
470 | else // no internal Dirichlet constraints specified | ||
471 | { | ||
472 |
2/2✓ Branch 0 taken 11953152 times.
✓ Branch 1 taken 9690152 times.
|
40304846 | for (int eqIdx = 0; eqIdx < numEq; eqIdx++) |
473 | { | ||
474 | // the diagonal entries | ||
475 |
4/8✓ Branch 1 taken 11953152 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 11953152 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 11953152 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 3944312 times.
✗ Branch 11 not taken.
|
72966085 | A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx]; |
476 | |||
477 | // off-diagonal entries | ||
478 | 22092455 | j = 1; | |
479 |
4/4✓ Branch 0 taken 27192832 times.
✓ Branch 1 taken 44048858 times.
✓ Branch 2 taken 27192832 times.
✓ Branch 3 taken 44048858 times.
|
311932032 | for (const auto& dataJ : connectivityMap[globalI]) |
480 |
4/8✓ Branch 1 taken 59288538 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 59288538 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 59288538 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 19183992 times.
✗ Branch 11 not taken.
|
369130998 | A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j++][eqIdx]; |
481 | } | ||
482 | } | ||
483 | } | ||
484 | |||
485 | // restore original state of the flux vars cache in case of global caching. | ||
486 | // This has to be done in order to guarantee that everything is in an undeflected | ||
487 | // state before the assembly of another element is called. In the case of local caching | ||
488 | // this is obsolete because the elemFluxVarsCache used here goes out of scope after this. | ||
489 | // We only have to do this for the last primary variable, for all others the flux var cache | ||
490 | // is updated with the correct element volume variables before residual evaluations | ||
491 | if (enableGridFluxVarsCache) | ||
492 | 15761057 | gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars); | |
493 | |||
494 | // compute potential additional derivatives causes by the coupling | ||
495 | 35768462 | typename ParentType::LocalResidual::ElementResidualVector origElementResidual({origResiduals[0]}); | |
496 | 17884231 | this->couplingManager().evalAdditionalDomainDerivatives(domainI, *this, origElementResidual, A, gridVariables); | |
497 | |||
498 | // return the original residual | ||
499 | 36741200 | return origResiduals[0]; | |
500 | } | ||
501 | |||
502 | /*! | ||
503 | * \brief Computes the derivatives with respect to the given element and adds them | ||
504 | * to the global matrix. | ||
505 | */ | ||
506 | template<std::size_t otherId, class JacobianBlock, class GridVariables> | ||
507 | 20564411 | void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A, | |
508 | const LocalResidualValues& res, GridVariables& gridVariables) | ||
509 | { | ||
510 | //////////////////////////////////////////////////////////////////////////////////////////////////////// | ||
511 | // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. // | ||
512 | //////////////////////////////////////////////////////////////////////////////////////////////////////// | ||
513 | |||
514 | // get some aliases for convenience | ||
515 | 20564411 | const auto& element = this->element(); | |
516 | 20564411 | const auto& fvGeometry = this->fvGeometry(); | |
517 | 20564411 | const auto& gridGeometry = fvGeometry.gridGeometry(); | |
518 | 20564411 | auto&& curElemVolVars = this->curElemVolVars(); | |
519 | 20564411 | auto&& elemFluxVarsCache = this->elemFluxVarsCache(); | |
520 | |||
521 | // get stencil information | ||
522 | 41128822 | const auto globalI = gridGeometry.elementMapper().index(element); | |
523 | 20564411 | const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ); | |
524 | 20564411 | const auto& curSolJ = this->curSol()[domainJ]; | |
525 | |||
526 | // make sure the flux variables cache does not contain any artifacts from prior differentiation | ||
527 | 20564411 | elemFluxVarsCache.update(element, fvGeometry, curElemVolVars); | |
528 | |||
529 | // convenience lambda for call to update self | ||
530 | 51210678 | auto updateCoupledVariables = [&] () | |
531 | { | ||
532 | // Update ourself after the context has been modified. Depending on the | ||
533 | // type of caching, other objects might have to be updated. All ifs can be optimized away. | ||
534 | if (enableGridFluxVarsCache) | ||
535 | { | ||
536 | if (enableGridVolVarsCache) | ||
537 | { | ||
538 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1089200 times.
|
34583726 | this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache()); |
539 | 56272531 | this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache); | |
540 | } | ||
541 | else | ||
542 | { | ||
543 | 668372 | this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, gridVariables.gridFluxVarsCache()); | |
544 | this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache); | ||
545 | } | ||
546 | } | ||
547 | else | ||
548 | { | ||
549 | if (enableGridVolVarsCache) | ||
550 | this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache); | ||
551 | else | ||
552 | 1141136 | this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache); | |
553 | } | ||
554 | }; | ||
555 | |||
556 |
4/4✓ Branch 0 taken 12382593 times.
✓ Branch 1 taken 11150163 times.
✓ Branch 2 taken 12382593 times.
✓ Branch 3 taken 11150163 times.
|
79546161 | for (const auto& globalJ : stencil) |
557 | { | ||
558 | // undeflected privars and privars to be deflected | ||
559 |
1/3✗ Branch 0 not taken.
✓ Branch 1 taken 12382593 times.
✗ Branch 2 not taken.
|
17585648 | const auto origPriVarsJ = curSolJ[globalJ]; |
560 | 17585648 | auto priVarsJ = origPriVarsJ; | |
561 | |||
562 |
1/3✗ Branch 0 not taken.
✓ Branch 1 taken 12382593 times.
✗ Branch 2 not taken.
|
18828672 | const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ)[0]; |
563 | |||
564 |
2/2✓ Branch 0 taken 13190879 times.
✓ Branch 1 taken 12382593 times.
|
36242092 | for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx) |
565 | { | ||
566 | 58408098 | auto evalCouplingResidual = [&](Scalar priVar) | |
567 | { | ||
568 | // update the volume variables and the flux var cache | ||
569 | 24815252 | priVarsJ[pvIdx] = priVar; | |
570 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1207072 times.
|
13190879 | this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx); |
571 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1089200 times.
|
13126339 | updateCoupledVariables(); |
572 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 1219532 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 7020 times.
|
13190879 | return this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ)[0]; |
573 | }; | ||
574 | |||
575 | // derive the residuals numerically | ||
576 |
2/2✓ Branch 0 taken 233 times.
✓ Branch 1 taken 11173362 times.
|
18656444 | LocalResidualValues partialDeriv(0.0); |
577 |
4/4✓ Branch 0 taken 280 times.
✓ Branch 1 taken 13190599 times.
✓ Branch 2 taken 280 times.
✓ Branch 3 taken 13190599 times.
|
37312888 | const auto& paramGroup = this->assembler().problem(domainJ).paramGroup(); |
578 |
5/6✓ Branch 0 taken 280 times.
✓ Branch 1 taken 13190599 times.
✓ Branch 3 taken 263 times.
✓ Branch 4 taken 17 times.
✓ Branch 6 taken 263 times.
✗ Branch 7 not taken.
|
18656444 | static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod"); |
579 |
5/6✓ Branch 0 taken 290 times.
✓ Branch 1 taken 13190589 times.
✓ Branch 3 taken 263 times.
✓ Branch 4 taken 27 times.
✓ Branch 6 taken 263 times.
✗ Branch 7 not taken.
|
18656444 | static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup); |
580 |
2/4✓ Branch 1 taken 13190879 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1535338 times.
✗ Branch 5 not taken.
|
20689666 | NumericDifferentiation::partialDerivative(evalCouplingResidual, priVarsJ[pvIdx], partialDeriv, origResidual, |
581 | 20689666 | epsCoupl(priVarsJ[pvIdx], pvIdx), numDiffMethod); | |
582 | |||
583 | // add the current partial derivatives to the global jacobian matrix | ||
584 | if constexpr (Problem::enableInternalDirichletConstraints()) | ||
585 | { | ||
586 |
2/3✓ Branch 0 taken 2834640 times.
✓ Branch 1 taken 1358864 times.
✗ Branch 2 not taken.
|
5018624 | const auto& scv = fvGeometry.scv(globalI); |
587 |
5/6✓ Branch 0 taken 2820496 times.
✓ Branch 1 taken 1263008 times.
✓ Branch 2 taken 2820496 times.
✓ Branch 3 taken 1234720 times.
✓ Branch 4 taken 28288 times.
✗ Branch 5 not taken.
|
10037248 | const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv); |
588 |
4/4✓ Branch 0 taken 4328772 times.
✓ Branch 1 taken 25852 times.
✓ Branch 2 taken 4328772 times.
✓ Branch 3 taken 25852 times.
|
10037248 | if (internalDirichletConstraints.none()) |
589 | { | ||
590 |
2/2✓ Branch 0 taken 6406004 times.
✓ Branch 1 taken 4992772 times.
|
11398776 | for (int eqIdx = 0; eqIdx < numEq; eqIdx++) |
591 |
3/6✓ Branch 1 taken 6406004 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 6406004 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 2509332 times.
✗ Branch 8 not taken.
|
15321340 | A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx]; |
592 | } | ||
593 | else | ||
594 | { | ||
595 |
2/2✓ Branch 0 taken 26220 times.
✓ Branch 1 taken 25852 times.
|
52072 | for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) |
596 | { | ||
597 |
4/4✓ Branch 0 taken 25216 times.
✓ Branch 1 taken 1004 times.
✓ Branch 2 taken 25216 times.
✓ Branch 3 taken 368 times.
|
51804 | if (internalDirichletConstraints[eqIdx]) |
598 |
2/4✓ Branch 1 taken 25852 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 25852 times.
✗ Branch 5 not taken.
|
51704 | A[globalI][globalJ][eqIdx][pvIdx] = 0.0; |
599 | else | ||
600 |
3/6✓ Branch 1 taken 368 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 368 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 368 times.
✗ Branch 8 not taken.
|
1104 | A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx]; |
601 | } | ||
602 | } | ||
603 | } | ||
604 | else | ||
605 | { | ||
606 |
2/3✓ Branch 0 taken 9176819 times.
✓ Branch 1 taken 8172255 times.
✗ Branch 2 not taken.
|
28864768 | for (int eqIdx = 0; eqIdx < numEq; eqIdx++) |
607 |
3/6✓ Branch 1 taken 9176819 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 9176819 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1925448 times.
✗ Branch 8 not taken.
|
33464792 | A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx]; |
608 | } | ||
609 | |||
610 | // restore the current element solution | ||
611 |
2/4✓ Branch 1 taken 125888 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 125888 times.
✗ Branch 5 not taken.
|
20689666 | priVarsJ[pvIdx] = origPriVarsJ[pvIdx]; |
612 | |||
613 | // restore the undeflected state of the coupling context | ||
614 |
1/2✓ Branch 1 taken 160965 times.
✗ Branch 2 not taken.
|
36968294 | this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx); |
615 | } | ||
616 | |||
617 | // restore original state of the flux vars cache and/or vol vars in case of global caching. | ||
618 | // This has to be done in order to guarantee that the undeflected residual computation done | ||
619 | // above is correct when jumping to the next coupled dof of domainJ. | ||
620 |
1/2✓ Branch 1 taken 35692 times.
✗ Branch 2 not taken.
|
17585648 | updateCoupledVariables(); |
621 | } | ||
622 | 20564411 | } | |
623 | }; | ||
624 | |||
625 | /*! | ||
626 | * \ingroup Assembly | ||
627 | * \ingroup CCDiscretization | ||
628 | * \ingroup MultiDomain | ||
629 | * \brief Cell-centered scheme multidomain local assembler using numeric differentiation and explicit time discretization | ||
630 | */ | ||
631 | template<std::size_t id, class TypeTag, class Assembler> | ||
632 | 2028000 | class SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false> | |
633 | : public SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler, | ||
634 | SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>, false> | ||
635 | { | ||
636 | using ThisType = SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>; | ||
637 | using ParentType = SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler, ThisType, /*implicit=*/false>; | ||
638 | |||
639 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
640 | using LocalResidualValues = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; | ||
641 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
642 | |||
643 | static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq(); | ||
644 | static constexpr auto domainI = Dune::index_constant<id>(); | ||
645 | |||
646 | public: | ||
647 | 2028000 | using ParentType::ParentType; | |
648 | |||
649 | /*! | ||
650 | * \brief Computes the derivatives with respect to the given element and adds them | ||
651 | * to the global matrix. | ||
652 | * | ||
653 | * \note In an explicit scheme, only the storage terms need to be differentiated. | ||
654 | * Thus, this can be done as in the uncoupled case as the coupling can only | ||
655 | * enter sources or fluxes. | ||
656 | * | ||
657 | * \return The element residual at the current solution. | ||
658 | */ | ||
659 | template<class JacobianMatrixDiagBlock, class GridVariables> | ||
660 | 4056000 | LocalResidualValues assembleJacobianAndResidualImplInverse(JacobianMatrixDiagBlock& A, GridVariables& gridVariables) | |
661 | { | ||
662 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 2028000 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2028000 times.
|
8112000 | if (this->assembler().isStationaryProblem()) |
663 | ✗ | DUNE_THROW(Dune::InvalidStateException, "Using explicit jacobian assembler with stationary local residual"); | |
664 | |||
665 | // assemble the undeflected residual | ||
666 | 4056000 | const auto residual = this->evalLocalResidual()[0]; | |
667 | 4056000 | const auto storageResidual = this->evalLocalStorageResidual(); | |
668 | |||
669 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
670 | // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the // | ||
671 | // neighboring elements all derivatives are zero. For the assembled element only the storage // | ||
672 | // derivatives are non-zero. // | ||
673 | ////////////////////////////////////////////////////////////////////////////////////////////////// | ||
674 | |||
675 | // get some aliases for convenience | ||
676 | 4056000 | const auto& element = this->element(); | |
677 | 4056000 | const auto& fvGeometry = this->fvGeometry(); | |
678 | 4056000 | const auto& gridGeometry = fvGeometry.gridGeometry(); | |
679 | 4056000 | auto&& curElemVolVars = this->curElemVolVars(); | |
680 | |||
681 | // reference to the element's scv (needed later) and corresponding vol vars | ||
682 | 8112000 | const auto globalI = gridGeometry.elementMapper().index(element); | |
683 |
1/2✓ Branch 0 taken 2028000 times.
✗ Branch 1 not taken.
|
4056000 | const auto& scv = fvGeometry.scv(globalI); |
684 | 8112000 | auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv); | |
685 | |||
686 | // save a copy of the original privars and vol vars in order | ||
687 | // to restore the original solution after deflection | ||
688 | 4056000 | const auto& curSol = this->curSol()[domainI]; | |
689 | 4056000 | const auto origPriVars = curSol[globalI]; | |
690 | 4056000 | const auto origVolVars = curVolVars; | |
691 | |||
692 | // element solution container to be deflected | ||
693 | 4056000 | auto elemSol = elementSolution(element, curSol, gridGeometry); | |
694 | |||
695 | // derivatives in the neighbors with respect to the current elements | ||
696 | 4056000 | LocalResidualValues partialDeriv; | |
697 |
2/2✓ Branch 0 taken 2028000 times.
✓ Branch 1 taken 2028000 times.
|
8112000 | for (int pvIdx = 0; pvIdx < numEq; ++pvIdx) |
698 | { | ||
699 | // reset derivatives of element dof with respect to itself | ||
700 |
1/2✓ Branch 0 taken 2028000 times.
✗ Branch 1 not taken.
|
4056000 | partialDeriv = 0.0; |
701 | |||
702 | 12168000 | auto evalStorage = [&](Scalar priVar) | |
703 | { | ||
704 | // update the volume variables and calculate | ||
705 | // the residual with the deflected primary variables | ||
706 | 2028000 | elemSol[0][pvIdx] = priVar; | |
707 | 4056000 | curVolVars.update(elemSol, this->problem(), element, scv); | |
708 | 2028000 | return this->evalLocalStorageResidual(); | |
709 | }; | ||
710 | |||
711 | // for non-ghosts compute the derivative numerically | ||
712 |
1/2✓ Branch 0 taken 2028000 times.
✗ Branch 1 not taken.
|
4056000 | if (!this->elementIsGhost()) |
713 | { | ||
714 |
6/10✓ Branch 0 taken 4 times.
✓ Branch 1 taken 2027996 times.
✓ Branch 3 taken 4 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 4 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 4 times.
✗ Branch 10 not taken.
✓ Branch 12 taken 4 times.
✗ Branch 13 not taken.
|
4056000 | static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()}; |
715 |
6/10✓ Branch 0 taken 4 times.
✓ Branch 1 taken 2027996 times.
✓ Branch 3 taken 4 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 4 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 4 times.
✗ Branch 10 not taken.
✓ Branch 12 taken 4 times.
✗ Branch 13 not taken.
|
4056000 | static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod"); |
716 |
1/2✓ Branch 1 taken 2028000 times.
✗ Branch 2 not taken.
|
4056000 | NumericDifferentiation::partialDerivative(evalStorage, elemSol[0][pvIdx], partialDeriv, storageResidual, |
717 | 4056000 | eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod); | |
718 | } | ||
719 | |||
720 | // for ghost elements we assemble a 1.0 where the primary variable and zero everywhere else | ||
721 | // as we always solve for a delta of the solution with respect to the initial solution this | ||
722 | // results in a delta of zero for ghosts | ||
723 | ✗ | else partialDeriv[pvIdx] = 1.0; | |
724 | |||
725 | // restore the original state of the scv's volume variables | ||
726 | 4056000 | curVolVars = origVolVars; | |
727 | |||
728 | // restore the current element solution | ||
729 | 4056000 | elemSol[0][pvIdx] = origPriVars[pvIdx]; | |
730 | |||
731 | // add the current partial derivatives to the global jacobian matrix | ||
732 | // only diagonal entries for explicit jacobians | ||
733 | if constexpr (Problem::enableInternalDirichletConstraints()) | ||
734 | { | ||
735 | // check if own element has internal Dirichlet constraint | ||
736 | const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv); | ||
737 | const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv); | ||
738 | |||
739 | for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) | ||
740 | { | ||
741 | if (internalDirichletConstraints[eqIdx]) | ||
742 | { | ||
743 | residual[eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx]; | ||
744 | A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0; | ||
745 | } | ||
746 | else | ||
747 | A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx]; | ||
748 | } | ||
749 | } | ||
750 | else | ||
751 | { | ||
752 |
2/2✓ Branch 0 taken 2028000 times.
✓ Branch 1 taken 2028000 times.
|
8112000 | for (int eqIdx = 0; eqIdx < numEq; eqIdx++) |
753 |
2/4✓ Branch 1 taken 2028000 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 2028000 times.
✗ Branch 5 not taken.
|
8112000 | A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx]; |
754 | } | ||
755 | } | ||
756 | |||
757 | // return the original residual | ||
758 | 4056000 | return residual; | |
759 | } | ||
760 | |||
761 | /*! | ||
762 | * \brief Computes the coupling derivatives with respect to the given element and adds them | ||
763 | * to the global matrix. | ||
764 | * \note Since the coupling can only enter sources or fluxes and these are evaluated on | ||
765 | * the old time level (explicit scheme), the coupling blocks are empty. | ||
766 | */ | ||
767 | template<std::size_t otherId, class JacobianBlock, class GridVariables> | ||
768 | ✗ | void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A, | |
769 | const LocalResidualValues& res, GridVariables& gridVariables) | ||
770 | ✗ | {} | |
771 | }; | ||
772 | |||
773 | /*! | ||
774 | * \ingroup Assembly | ||
775 | * \ingroup CCDiscretization | ||
776 | * \ingroup MultiDomain | ||
777 | * \brief Cell-centered scheme local assembler using analytic differentiation and implicit time discretization | ||
778 | */ | ||
779 | template<std::size_t id, class TypeTag, class Assembler> | ||
780 | class SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true> | ||
781 | : public SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler, | ||
782 | SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::analytic, true>, true> | ||
783 | { | ||
784 | using ThisType = SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true>; | ||
785 | using ParentType = SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler, ThisType, /*implicit=*/true>; | ||
786 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
787 | using LocalResidualValues = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; | ||
788 | using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView; | ||
789 | using Element = typename GridView::template Codim<0>::Entity; | ||
790 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
791 | |||
792 | enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() }; | ||
793 | enum { dim = GridView::dimension }; | ||
794 | |||
795 | static constexpr auto domainI = Dune::index_constant<id>(); | ||
796 | |||
797 | public: | ||
798 | using ParentType::ParentType; | ||
799 | |||
800 | /*! | ||
801 | * \brief Computes the derivatives with respect to the given element and adds them | ||
802 | * to the global matrix. | ||
803 | * | ||
804 | * \return The element residual at the current solution. | ||
805 | */ | ||
806 | template<class JacobianMatrixDiagBlock, class GridVariables> | ||
807 | LocalResidualValues assembleJacobianAndResidualImpl(JacobianMatrixDiagBlock& A, GridVariables& gridVariables) | ||
808 | { | ||
809 | // treat ghost separately, we always want zero update for ghosts | ||
810 | if (this->elementIsGhost()) | ||
811 | { | ||
812 | const auto globalI = this->assembler().gridGeometry(domainI).elementMapper().index(this->element()); | ||
813 | for (int pvIdx = 0; pvIdx < numEq; ++pvIdx) | ||
814 | A[globalI][globalI][pvIdx][pvIdx] = 1.0; | ||
815 | |||
816 | // return zero residual | ||
817 | return LocalResidualValues(0.0); | ||
818 | } | ||
819 | |||
820 | // get some aliases for convenience | ||
821 | const auto& problem = this->problem(); | ||
822 | const auto& element = this->element(); | ||
823 | const auto& fvGeometry = this->fvGeometry(); | ||
824 | const auto& curElemVolVars = this->curElemVolVars(); | ||
825 | const auto& elemFluxVarsCache = this->elemFluxVarsCache(); | ||
826 | |||
827 | // get reference to the element's current vol vars | ||
828 | const auto globalI = this->assembler().gridGeometry(domainI).elementMapper().index(element); | ||
829 | const auto& scv = fvGeometry.scv(globalI); | ||
830 | const auto& volVars = curElemVolVars[scv]; | ||
831 | |||
832 | // if the problem is instationary, add derivative of storage term | ||
833 | if (!this->assembler().isStationaryProblem()) | ||
834 | this->localResidual().addStorageDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv); | ||
835 | |||
836 | // add source term derivatives | ||
837 | this->localResidual().addSourceDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv); | ||
838 | |||
839 | // add flux derivatives for each scvf | ||
840 | for (const auto& scvf : scvfs(fvGeometry)) | ||
841 | { | ||
842 | // inner faces | ||
843 | if (!scvf.boundary()) | ||
844 | this->localResidual().addFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf); | ||
845 | |||
846 | // boundary faces | ||
847 | else | ||
848 | { | ||
849 | const auto& bcTypes = problem.boundaryTypes(element, scvf); | ||
850 | |||
851 | // add Dirichlet boundary flux derivatives | ||
852 | if (bcTypes.hasDirichlet() && !bcTypes.hasNeumann()) | ||
853 | this->localResidual().addCCDirichletFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf); | ||
854 | |||
855 | // add Robin ("solution dependent Neumann") boundary flux derivatives | ||
856 | else if (bcTypes.hasNeumann() && !bcTypes.hasDirichlet()) | ||
857 | this->localResidual().addRobinFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf); | ||
858 | |||
859 | else | ||
860 | DUNE_THROW(Dune::NotImplemented, "Mixed boundary conditions. Use pure boundary conditions by converting Dirichlet BCs to Robin BCs"); | ||
861 | } | ||
862 | } | ||
863 | |||
864 | if constexpr (Problem::enableInternalDirichletConstraints()) | ||
865 | { | ||
866 | // check if own element has internal Dirichlet constraint | ||
867 | const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv); | ||
868 | const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv); | ||
869 | |||
870 | auto residual = this->evalLocalResidual()[0]; | ||
871 | |||
872 | for (int pvIdx = 0; pvIdx < numEq; ++pvIdx) | ||
873 | { | ||
874 | for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) | ||
875 | { | ||
876 | if (internalDirichletConstraints[eqIdx]) | ||
877 | { | ||
878 | residual[eqIdx] = volVars.priVars()[eqIdx] - dirichletValues[eqIdx]; | ||
879 | A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0; | ||
880 | |||
881 | // inner faces | ||
882 | for (const auto& scvf : scvfs(fvGeometry)) | ||
883 | if (!scvf.boundary()) | ||
884 | A[globalI][fvGeometry.scv(scvf.outsideScvIdx()).dofIndex()][eqIdx][pvIdx] = 0.0; | ||
885 | } | ||
886 | } | ||
887 | } | ||
888 | return residual; | ||
889 | } | ||
890 | else | ||
891 | // return element residual | ||
892 | return this->evalLocalResidual()[0]; | ||
893 | } | ||
894 | |||
895 | /*! | ||
896 | * \brief Computes the derivatives with respect to the given element and adds them | ||
897 | * to the global matrix. | ||
898 | */ | ||
899 | template<std::size_t otherId, class JacobianBlock, class GridVariables> | ||
900 | void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A, | ||
901 | const LocalResidualValues& res, GridVariables& gridVariables) | ||
902 | { | ||
903 | //////////////////////////////////////////////////////////////////////////////////////////////////////// | ||
904 | // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. // | ||
905 | //////////////////////////////////////////////////////////////////////////////////////////////////////// | ||
906 | |||
907 | // get some aliases for convenience | ||
908 | const auto& element = this->element(); | ||
909 | const auto& fvGeometry = this->fvGeometry(); | ||
910 | const auto& gridGeometry = fvGeometry.gridGeometry(); | ||
911 | auto&& curElemVolVars = this->curElemVolVars(); | ||
912 | // auto&& elemFluxVarsCache = this->elemFluxVarsCache(); | ||
913 | |||
914 | // get stencil information | ||
915 | const auto globalI = gridGeometry.elementMapper().index(element); | ||
916 | const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ); | ||
917 | |||
918 | for (const auto globalJ : stencil) | ||
919 | { | ||
920 | const auto& elementJ = this->assembler().gridGeometry(domainJ).element(globalJ); | ||
921 | this->couplingManager().addCouplingDerivatives(A[globalI][globalJ], domainI, element, fvGeometry, curElemVolVars, domainJ, elementJ); | ||
922 | |||
923 | // add the current partial derivatives to the global jacobian matrix | ||
924 | if constexpr (Problem::enableInternalDirichletConstraints()) | ||
925 | { | ||
926 | const auto& scv = fvGeometry.scv(globalI); | ||
927 | const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv); | ||
928 | if (internalDirichletConstraints.any()) | ||
929 | { | ||
930 | for (int pvIdx = 0; pvIdx < numEq; ++pvIdx) | ||
931 | for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) | ||
932 | if (internalDirichletConstraints[eqIdx]) | ||
933 | A[globalI][globalJ][eqIdx][pvIdx] = 0.0; | ||
934 | } | ||
935 | } | ||
936 | } | ||
937 | } | ||
938 | }; | ||
939 | |||
940 | } // end namespace Dumux | ||
941 | |||
942 | #endif | ||
943 |