Line | Branch | Exec | Source |
---|---|---|---|
1 | // -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*- | ||
2 | // vi: set et ts=4 sw=4 sts=4: | ||
3 | // | ||
4 | // SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder | ||
5 | // SPDX-License-Identifier: GPL-3.0-or-later | ||
6 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup Experimental | ||
10 | * \ingroup Assembly | ||
11 | * \ingroup CVFEDiscretization | ||
12 | * \ingroup MultiDomain | ||
13 | * \brief An assembler for Jacobian and residual contribution per element (CVFE methods) for multidomain problems | ||
14 | */ | ||
15 | #ifndef DUMUX_EXPERIMENTAL_SUBDOMAIN_CVFE_LOCAL_ASSEMBLER_HH | ||
16 | #define DUMUX_EXPERIMENTAL_SUBDOMAIN_CVFE_LOCAL_ASSEMBLER_HH | ||
17 | |||
18 | #include <type_traits> | ||
19 | |||
20 | #include <dune/common/reservedvector.hh> | ||
21 | #include <dune/common/indices.hh> | ||
22 | #include <dune/common/hybridutilities.hh> | ||
23 | #include <dune/grid/common/gridenums.hh> // for GhostEntity | ||
24 | #include <dune/istl/matrixindexset.hh> | ||
25 | |||
26 | #include <dumux/common/typetraits/localdofs_.hh> | ||
27 | #include <dumux/common/reservedblockvector.hh> | ||
28 | #include <dumux/common/properties.hh> | ||
29 | #include <dumux/common/parameters.hh> | ||
30 | #include <dumux/common/numericdifferentiation.hh> | ||
31 | #include <dumux/common/numeqvector.hh> | ||
32 | #include <dumux/assembly/numericepsilon.hh> | ||
33 | #include <dumux/assembly/diffmethod.hh> | ||
34 | #include <dumux/discretization/extrusion.hh> | ||
35 | |||
36 | #include <dumux/experimental/assembly/cvfelocalassembler.hh> | ||
37 | |||
38 | namespace Dumux::Experimental { | ||
39 | |||
40 | /*! | ||
41 | * \ingroup Experimental | ||
42 | * \ingroup Assembly | ||
43 | * \ingroup CVFEDiscretization | ||
44 | * \ingroup MultiDomain | ||
45 | * \brief A base class for all CVFE subdomain local assemblers | ||
46 | * \tparam id the id of the sub domain | ||
47 | * \tparam TypeTag the TypeTag | ||
48 | * \tparam Assembler the assembler type | ||
49 | * \tparam Implementation the actual implementation type | ||
50 | */ | ||
51 | template<std::size_t id, class TypeTag, class Assembler, class Implementation, DiffMethod dm> | ||
52 | 9472 | class SubDomainCVFELocalAssemblerBase : public Experimental::CVFELocalAssembler<TypeTag, Assembler, dm, Implementation> | |
53 | { | ||
54 | using ParentType = Experimental::CVFELocalAssembler<TypeTag, Assembler, dm, Implementation>; | ||
55 | |||
56 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
57 | using LocalResidualValues = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; | ||
58 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
59 | using SolutionVector = typename Assembler::SolutionVector; | ||
60 | using ElementBoundaryTypes = GetPropType<TypeTag, Properties::ElementBoundaryTypes>; | ||
61 | |||
62 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
63 | using GridVolumeVariables = typename GridVariables::GridVolumeVariables; | ||
64 | using ElementVolumeVariables = typename GridVolumeVariables::LocalView; | ||
65 | using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView; | ||
66 | using Scalar = typename GridVariables::Scalar; | ||
67 | |||
68 | using GridGeometry = typename GridVariables::GridGeometry; | ||
69 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
70 | using SubControlVolume = typename GridGeometry::SubControlVolume; | ||
71 | using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace; | ||
72 | using Extrusion = Extrusion_t<GridGeometry>; | ||
73 | using GridView = typename GridGeometry::GridView; | ||
74 | using Element = typename GridView::template Codim<0>::Entity; | ||
75 | |||
76 | using CouplingManager = typename Assembler::CouplingManager; | ||
77 | |||
78 | static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq(); | ||
79 | |||
80 | public: | ||
81 | //! export the domain id of this sub-domain | ||
82 | static constexpr auto domainId = typename Dune::index_constant<id>(); | ||
83 | //! pull up constructor of parent class | ||
84 | using ParentType::ParentType; | ||
85 | //! export element residual vector type | ||
86 | using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector; | ||
87 | |||
88 | // the constructor | ||
89 | 4736 | explicit SubDomainCVFELocalAssemblerBase( | |
90 | const Assembler& assembler, | ||
91 | const Element& element, | ||
92 | const SolutionVector& curSol, | ||
93 | CouplingManager& couplingManager | ||
94 | ) | ||
95 | : ParentType(assembler, | ||
96 | element, | ||
97 | curSol, | ||
98 | localView(assembler.gridGeometry(domainId)), | ||
99 |
1/2✓ Branch 1 taken 4736 times.
✗ Branch 2 not taken.
|
4736 | localView(assembler.gridVariables(domainId).curGridVolVars()), |
100 |
1/2✓ Branch 1 taken 4736 times.
✗ Branch 2 not taken.
|
4736 | localView(assembler.gridVariables(domainId).gridFluxVarsCache()), |
101 |
1/2✓ Branch 1 taken 4736 times.
✗ Branch 2 not taken.
|
4736 | assembler.localResidual(domainId), |
102 | 4736 | (element.partitionType() == Dune::GhostEntity), | |
103 | 4736 | assembler.isImplicit()) | |
104 |
1/2✓ Branch 3 taken 4736 times.
✗ Branch 4 not taken.
|
9472 | , couplingManager_(couplingManager) |
105 | 4736 | {} | |
106 | |||
107 | /*! | ||
108 | * \brief Computes the derivatives with respect to the given element and adds them | ||
109 | * to the global matrix. The element residual is written into the right hand side. | ||
110 | */ | ||
111 | template<class JacobianMatrixRow, class SubResidualVector, class GridVariablesTuple, class StageParams> | ||
112 | 3456 | void assembleJacobianAndResidual(JacobianMatrixRow& jacRow, SubResidualVector& res, GridVariablesTuple& gridVariables, | |
113 | const StageParams& stageParams, SubResidualVector& temporal, SubResidualVector& spatial, | ||
114 | SubResidualVector& constrainedDofs) | ||
115 | { | ||
116 | 3456 | auto assembleCouplingBlocks = [&](const auto& residual) | |
117 | { | ||
118 | // assemble the coupling blocks | ||
119 | using namespace Dune::Hybrid; | ||
120 | 3456 | forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](auto&& i) | |
121 | { | ||
122 | if constexpr (std::decay_t<decltype(i)>{} != id) | ||
123 | 3456 | this->assembleJacobianCoupling(i, jacRow, residual, gridVariables); | |
124 | }); | ||
125 | }; | ||
126 | |||
127 | // the coupled model does not support partial reassembly yet | ||
128 | 3456 | const DefaultPartialReassembler* noReassembler = nullptr; | |
129 | 3456 | ParentType::assembleJacobianAndResidual( | |
130 | 3456 | jacRow[domainId], res, | |
131 |
1/2✓ Branch 1 taken 3456 times.
✗ Branch 2 not taken.
|
3456 | *std::get<domainId>(gridVariables), |
132 | stageParams, temporal, spatial, | ||
133 | constrainedDofs, | ||
134 | noReassembler, | ||
135 | assembleCouplingBlocks | ||
136 | ); | ||
137 | } | ||
138 | |||
139 | /*! | ||
140 | * \brief Assemble the entries in a coupling block of the jacobian. | ||
141 | * There is no coupling block between a domain and itself. | ||
142 | */ | ||
143 | template<std::size_t otherId, class JacRow, class GridVariables> | ||
144 | 3456 | void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow, | |
145 | const ElementResidualVector& res, GridVariables& gridVariables) | ||
146 | { | ||
147 | if constexpr (id != otherId) | ||
148 | 3456 | this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables)); | |
149 | } | ||
150 | |||
151 | /*! | ||
152 | * \brief Prepares all local views necessary for local assembly. | ||
153 | */ | ||
154 | 4736 | void bindLocalViews() | |
155 | { | ||
156 | // get some references for convenience | ||
157 | 4736 | const auto& element = this->element(); | |
158 | 4736 | const auto& curSol = this->curSol(domainId); | |
159 | 4736 | auto&& fvGeometry = this->fvGeometry(); | |
160 | 4736 | auto&& curElemVolVars = this->curElemVolVars(); | |
161 | 4736 | auto&& elemFluxVarsCache = this->elemFluxVarsCache(); | |
162 | |||
163 | // bind the caches | ||
164 | 4736 | couplingManager_.bindCouplingContext(domainId, element, this->assembler()); | |
165 | 4736 | fvGeometry.bind(element); | |
166 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4736 times.
|
4736 | if (std::abs(this->localResidual().spatialWeight()) < 1e-6) |
167 | ✗ | curElemVolVars.bindElement(element, fvGeometry, curSol); | |
168 | else | ||
169 | { | ||
170 | 4736 | curElemVolVars.bind(element, fvGeometry, curSol); | |
171 | 4736 | elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars); | |
172 | } | ||
173 | |||
174 | 4736 | this->elemBcTypes().update(problem(), this->element(), this->fvGeometry()); | |
175 | 4736 | } | |
176 | |||
177 | //! return reference to the underlying problem | ||
178 | template<std::size_t i = domainId> | ||
179 | 139094 | const Problem& problem(Dune::index_constant<i> dId = domainId) const | |
180 |
3/10✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
|
155078 | { return this->assembler().problem(domainId); } |
181 | |||
182 | //! return reference to the underlying problem | ||
183 | template<std::size_t i = domainId> | ||
184 | 11648 | const auto& curSol(Dune::index_constant<i> dId = domainId) const | |
185 |
2/2✓ Branch 1 taken 2 times.
✓ Branch 2 taken 3454 times.
|
8192 | { return ParentType::curSol()[dId]; } |
186 | |||
187 | //! return reference to the coupling manager | ||
188 | 193536 | CouplingManager& couplingManager() | |
189 | 203906 | { return couplingManager_; } | |
190 | |||
191 | private: | ||
192 | CouplingManager& couplingManager_; //!< the coupling manager | ||
193 | }; | ||
194 | |||
195 | /*! | ||
196 | * \ingroup Experimental | ||
197 | * \ingroup Assembly | ||
198 | * \ingroup CVFEDiscretization | ||
199 | * \ingroup MultiDomain | ||
200 | * \brief The CVFE scheme multidomain local assembler | ||
201 | * \tparam id the id of the sub domain | ||
202 | * \tparam TypeTag the TypeTag | ||
203 | * \tparam Assembler the assembler type | ||
204 | * \tparam DM the numeric differentiation method | ||
205 | */ | ||
206 | template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric> | ||
207 | class SubDomainCVFELocalAssembler; | ||
208 | |||
209 | /*! | ||
210 | * \ingroup Experimental | ||
211 | * \ingroup Assembly | ||
212 | * \ingroup CVFEDiscretization | ||
213 | * \ingroup MultiDomain | ||
214 | * \brief CVFE scheme multi domain local assembler using numeric differentiation | ||
215 | */ | ||
216 | template<std::size_t id, class TypeTag, class Assembler> | ||
217 | 4736 | class SubDomainCVFELocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric> | |
218 | : public SubDomainCVFELocalAssemblerBase<id, TypeTag, Assembler, | ||
219 | SubDomainCVFELocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric>, DiffMethod::numeric> | ||
220 | { | ||
221 | using ThisType = SubDomainCVFELocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric>; | ||
222 | using ParentType = SubDomainCVFELocalAssemblerBase<id, TypeTag, Assembler, ThisType, DiffMethod::numeric>; | ||
223 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
224 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
225 | |||
226 | using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView; | ||
227 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
228 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
229 | using SubControlVolume = typename GridGeometry::SubControlVolume; | ||
230 | using GridView = typename GridGeometry::GridView; | ||
231 | using Element = typename GridView::template Codim<0>::Entity; | ||
232 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
233 | |||
234 | static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq(); | ||
235 | static constexpr int dim = GridView::dimension; | ||
236 | |||
237 | static constexpr bool enableGridFluxVarsCache = GetPropType<TypeTag, Properties::GridVariables>::GridFluxVariablesCache::cachingEnabled; | ||
238 | static constexpr bool enableGridVolVarsCache = GetPropType<TypeTag, Properties::GridVariables>::GridVolumeVariables::cachingEnabled; | ||
239 | static constexpr auto domainI = Dune::index_constant<id>(); | ||
240 | |||
241 | public: | ||
242 | 4736 | using ParentType::ParentType; | |
243 | //! export element residual vector type | ||
244 | using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector; | ||
245 | |||
246 | /*! | ||
247 | * \brief Update the coupling context for coupled models. | ||
248 | */ | ||
249 | template<class ScvOrLocalDof, class ElemSol> | ||
250 | 165888 | void maybeUpdateCouplingContext(const ScvOrLocalDof& scvOrLocalDof, ElemSol& elemSol, const int pvIdx) | |
251 | { | ||
252 |
1/2✓ Branch 1 taken 165888 times.
✗ Branch 2 not taken.
|
165888 | if (this->assembler().isImplicit()) |
253 | 165888 | this->couplingManager().updateCouplingContext(domainI, *this, domainI, scvOrLocalDof.dofIndex(), elemSol[Dumux::Detail::LocalDofs::index(scvOrLocalDof)], pvIdx); | |
254 | 165888 | } | |
255 | |||
256 | /*! | ||
257 | * \brief Update the additional domain derivatives for coupled models. | ||
258 | */ | ||
259 | template<class JacobianMatrixDiagBlock, class GridVariables> | ||
260 | 3456 | void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector& origResiduals, JacobianMatrixDiagBlock& A, GridVariables& gridVariables) | |
261 | { | ||
262 | 3456 | if (this->assembler().isImplicit()) | |
263 | this->couplingManager().evalAdditionalDomainDerivatives(domainI, *this, origResiduals, A, gridVariables); | ||
264 | } | ||
265 | |||
266 | /*! | ||
267 | * \brief Computes the derivatives with respect to the given element and adds them | ||
268 | * to the global matrix. | ||
269 | */ | ||
270 | template<std::size_t otherId, class JacobianBlock, class GridVariables> | ||
271 | 3456 | void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A, | |
272 | const ElementResidualVector& res, GridVariables& gridVariables) | ||
273 | { | ||
274 |
1/2✗ Branch 1 not taken.
✓ Branch 2 taken 3456 times.
|
3456 | if (!this->assembler().isImplicit()) |
275 | ✗ | return; | |
276 | |||
277 | //////////////////////////////////////////////////////////////////////////////////////////////////////// | ||
278 | // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. // | ||
279 | //////////////////////////////////////////////////////////////////////////////////////////////////////// | ||
280 | |||
281 | // get some aliases for convenience | ||
282 | 3456 | const auto& element = this->element(); | |
283 | 3456 | const auto& fvGeometry = this->fvGeometry(); | |
284 | 3456 | auto&& curElemVolVars = this->curElemVolVars(); | |
285 | 3456 | auto&& elemFluxVarsCache = this->elemFluxVarsCache(); | |
286 | |||
287 | // convenience lambda for call to update self | ||
288 | 24192 | auto updateCoupledVariables = [&] () | |
289 | { | ||
290 | // Update ourself after the context has been modified. Depending on the | ||
291 | // type of caching, other objects might have to be updated. All ifs can be optimized away. | ||
292 | if constexpr (enableGridFluxVarsCache) | ||
293 | { | ||
294 | if constexpr (enableGridVolVarsCache) | ||
295 | this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache()); | ||
296 | else | ||
297 | this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, gridVariables.gridFluxVarsCache()); | ||
298 | } | ||
299 | else | ||
300 | { | ||
301 | if constexpr (enableGridVolVarsCache) | ||
302 | this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache); | ||
303 | else | ||
304 | 10368 | this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache); | |
305 | } | ||
306 | }; | ||
307 | |||
308 | // get element stencil information | ||
309 | 3456 | const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ); | |
310 | 3456 | const auto& curSolJ = this->curSol(domainJ); | |
311 |
2/2✓ Branch 0 taken 3456 times.
✓ Branch 1 taken 3456 times.
|
6912 | for (const auto globalJ : stencil) |
312 | { | ||
313 | // undeflected privars and privars to be deflected | ||
314 | 3456 | const auto origPriVarsJ = curSolJ[globalJ]; | |
315 | 3456 | auto priVarsJ = origPriVarsJ; | |
316 | |||
317 | // the undeflected coupling residual | ||
318 | 3456 | const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ); | |
319 | |||
320 |
2/2✓ Branch 0 taken 6912 times.
✓ Branch 1 taken 3456 times.
|
10368 | for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx) |
321 | { | ||
322 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6910 times.
|
6912 | auto evalCouplingResidual = [&](Scalar priVar) |
323 | { | ||
324 | 6912 | priVarsJ[pvIdx] = priVar; | |
325 | 6912 | this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx); | |
326 | 6912 | updateCoupledVariables(); | |
327 | 6912 | return this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ); | |
328 | }; | ||
329 | |||
330 | // derive the residuals numerically | ||
331 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6910 times.
|
6912 | ElementResidualVector partialDerivs(fvGeometry.numScv()); |
332 | |||
333 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6910 times.
|
6912 | const auto& paramGroup = this->assembler().problem(domainJ).paramGroup(); |
334 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6910 times.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
|
6912 | static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod"); |
335 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6910 times.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
|
6912 | static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup); |
336 | |||
337 | 6912 | NumericDifferentiation::partialDerivative(evalCouplingResidual, origPriVarsJ[pvIdx], partialDerivs, origResidual, | |
338 | 6912 | epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod); | |
339 | |||
340 | // update the global stiffness matrix with the current partial derivatives | ||
341 |
2/2✓ Branch 0 taken 55296 times.
✓ Branch 1 taken 6912 times.
|
62208 | for (const auto& scv : scvs(fvGeometry)) |
342 | { | ||
343 |
2/2✓ Branch 0 taken 165888 times.
✓ Branch 1 taken 55296 times.
|
221184 | for (int eqIdx = 0; eqIdx < numEq; eqIdx++) |
344 | { | ||
345 | // A[i][col][eqIdx][pvIdx] is the rate of change of | ||
346 | // the residual of equation 'eqIdx' at dof 'i' | ||
347 | // depending on the primary variable 'pvIdx' at dof | ||
348 | // 'col'. | ||
349 |
2/2✓ Branch 1 taken 145152 times.
✓ Branch 2 taken 20736 times.
|
165888 | A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx]; |
350 | |||
351 | // If the dof is coupled by a Dirichlet condition, | ||
352 | // set the derived value only once (i.e. overwrite existing values). | ||
353 |
2/2✓ Branch 0 taken 145152 times.
✓ Branch 1 taken 20736 times.
|
165888 | if (this->elemBcTypes().hasDirichlet()) |
354 | { | ||
355 | 145152 | const auto bcTypes = this->elemBcTypes().get(fvGeometry, scv); | |
356 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 145152 times.
|
145152 | if (bcTypes.isCouplingDirichlet(eqIdx)) |
357 | ✗ | A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][eqIdx]; | |
358 |
2/2✓ Branch 0 taken 95904 times.
✓ Branch 1 taken 49248 times.
|
241056 | else if (bcTypes.isDirichlet(eqIdx)) |
359 | 95904 | A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0; | |
360 | } | ||
361 | |||
362 | // enforce internal Dirichlet constraints | ||
363 | if constexpr (Problem::enableInternalDirichletConstraints()) | ||
364 | { | ||
365 | const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv); | ||
366 | if (internalDirichletConstraints[eqIdx]) | ||
367 | A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0; | ||
368 | } | ||
369 | |||
370 | } | ||
371 | } | ||
372 | |||
373 | // restore the current element solution | ||
374 | 6912 | priVarsJ[pvIdx] = origPriVarsJ[pvIdx]; | |
375 | |||
376 | // restore the undeflected state of the coupling context | ||
377 | 6912 | this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx); | |
378 | } | ||
379 | |||
380 | // Restore original state of the flux vars cache and/or vol vars. | ||
381 | // This has to be done in case they depend on variables of domainJ before | ||
382 | // we continue with the numeric derivative w.r.t the next globalJ. Otherwise, | ||
383 | // the next "origResidual" will be incorrect. | ||
384 | 3456 | updateCoupledVariables(); | |
385 | } | ||
386 | } | ||
387 | }; | ||
388 | |||
389 | } // end namespace Dumux | ||
390 | |||
391 | #endif | ||
392 |