3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
subdomainpq1bubblelocalassembler.hh
Go to the documentation of this file.
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
26#ifndef DUMUX_MULTIDOMAIN_PQ1BUBBLE_SUBDOMAIN_LOCAL_ASSEMBLER_HH
27#define DUMUX_MULTIDOMAIN_PQ1BUBBLE_SUBDOMAIN_LOCAL_ASSEMBLER_HH
28
29#include <dune/common/indices.hh>
30#include <dune/common/hybridutilities.hh>
31#include <dune/grid/common/gridenums.hh> // for GhostEntity
32
39
40namespace Dumux {
41
53template<std::size_t id, class TypeTag, class Assembler, class Implementation, DiffMethod dm, bool implicit>
54class SubDomainPQ1BubbleLocalAssemblerBase : public PQ1BubbleLocalAssembler<TypeTag, Assembler, dm, implicit, Implementation>
55{
57
59 using SolutionVector = typename Assembler::SolutionVector;
61
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 GridView = typename GridGeometry::GridView;
72 using Element = typename GridView::template Codim<0>::Entity;
73
74 using CouplingManager = typename Assembler::CouplingManager;
75
76 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
77
78public:
80 static constexpr auto domainId = typename Dune::index_constant<id>();
82 using ParentType::ParentType;
84 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
85
86 // the constructor
88 const Assembler& assembler,
89 const Element& element,
90 const SolutionVector& curSol,
91 CouplingManager& couplingManager
92 )
93 : ParentType(assembler,
94 element,
95 curSol,
96 localView(assembler.gridGeometry(domainId)),
97 localView(assembler.gridVariables(domainId).curGridVolVars()),
98 localView(assembler.gridVariables(domainId).prevGridVolVars()),
99 localView(assembler.gridVariables(domainId).gridFluxVarsCache()),
100 assembler.localResidual(domainId),
101 (element.partitionType() == Dune::GhostEntity))
102 , couplingManager_(couplingManager)
103 {}
104
109 template<class JacobianMatrixRow, class GridVariablesTuple>
110 void assembleJacobianAndResidual(JacobianMatrixRow& jacRow, SubSolutionVector& res, GridVariablesTuple& gridVariables)
111 {
112 auto assembleCouplingBlocks = [&](const auto& residual)
113 {
114 // assemble the coupling blocks
115 using namespace Dune::Hybrid;
116 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](auto&& i)
117 {
118 if (i != id)
119 this->assembleJacobianCoupling(i, jacRow, residual, gridVariables);
120 });
121 };
122
123 // the coupled model does not support partial reassembly yet
124 const DefaultPartialReassembler* noReassembler = nullptr;
125 ParentType::assembleJacobianAndResidual(jacRow[domainId], res, *std::get<domainId>(gridVariables), noReassembler, assembleCouplingBlocks);
126 }
127
132 template<std::size_t otherId, class JacRow, class GridVariables,
133 typename std::enable_if_t<(otherId == id), int> = 0>
134 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
135 const ElementResidualVector& res, GridVariables& gridVariables)
136 {}
137
141 template<std::size_t otherId, class JacRow, class GridVariables,
142 typename std::enable_if_t<(otherId != id), int> = 0>
143 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
144 const ElementResidualVector& res, GridVariables& gridVariables)
145 {
146 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
147 }
148
152 ElementResidualVector evalLocalSourceResidual(const Element& element, const ElementVolumeVariables& elemVolVars) const
153 {
154 // initialize the residual vector for all scvs in this element
155 ElementResidualVector residual(this->fvGeometry().numScv());
156
157 // evaluate the volume terms (storage + source terms)
158 // forward to the local residual specialized for the discretization methods
159 for (auto&& scv : scvs(this->fvGeometry()))
160 {
161 const auto& curVolVars = elemVolVars[scv];
162 auto source = this->localResidual().computeSource(problem(), element, this->fvGeometry(), elemVolVars, scv);
163 source *= -scv.volume()*curVolVars.extrusionFactor();
164 residual[scv.localDofIndex()] = std::move(source);
165 }
166
167 return residual;
168 }
169
173 ElementResidualVector evalLocalSourceResidual(const Element& neighbor) const
174 { return this->evalLocalSourceResidual(neighbor, implicit ? this->curElemVolVars() : this->prevElemVolVars()); }
175
180 {
181 // get some references for convenience
182 const auto& element = this->element();
183 const auto& curSol = this->curSol(domainId);
184 auto&& fvGeometry = this->fvGeometry();
185 auto&& curElemVolVars = this->curElemVolVars();
186 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
187
188 // bind the caches
189 couplingManager_.bindCouplingContext(domainId, element, this->assembler());
190 fvGeometry.bind(element);
191
192 if (implicit)
193 {
194 curElemVolVars.bind(element, fvGeometry, curSol);
195 elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars);
196 if (!this->assembler().isStationaryProblem())
197 this->prevElemVolVars().bindElement(element, fvGeometry, this->assembler().prevSol()[domainId]);
198 }
199 else
200 {
201 auto& prevElemVolVars = this->prevElemVolVars();
202 const auto& prevSol = this->assembler().prevSol()[domainId];
203
204 curElemVolVars.bindElement(element, fvGeometry, curSol);
205 prevElemVolVars.bind(element, fvGeometry, prevSol);
206 elemFluxVarsCache.bind(element, fvGeometry, prevElemVolVars);
207 }
208
209 this->elemBcTypes().update(problem(), this->element(), this->fvGeometry());
210 }
211
213 template<std::size_t i = domainId>
214 const Problem& problem(Dune::index_constant<i> dId = domainId) const
215 { return this->assembler().problem(domainId); }
216
218 template<std::size_t i = domainId>
219 const auto& curSol(Dune::index_constant<i> dId = domainId) const
220 { return ParentType::curSol()[dId]; }
221
223 CouplingManager& couplingManager()
224 { return couplingManager_; }
225
226private:
227 CouplingManager& couplingManager_;
228};
229
241template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric, bool implicit = true>
243
250template<std::size_t id, class TypeTag, class Assembler>
251class SubDomainPQ1BubbleLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>
252: public SubDomainPQ1BubbleLocalAssemblerBase<id, TypeTag, Assembler,
254{
255 using ThisType = SubDomainPQ1BubbleLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>;
256 using ParentType = SubDomainPQ1BubbleLocalAssemblerBase<id, TypeTag, Assembler, ThisType, DiffMethod::numeric, /*implicit=*/true>;
259
262 using FVElementGeometry = typename GridGeometry::LocalView;
263 using GridView = typename GridGeometry::GridView;
264 using Element = typename GridView::template Codim<0>::Entity;
265 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
266
268 enum { dim = GridView::dimension };
269
272 static constexpr auto domainI = Dune::index_constant<id>();
273
274public:
278
282 template<class ElemSol>
283 void maybeUpdateCouplingContext(const SubControlVolume& scv, ElemSol& elemSol, const int pvIdx)
284 {
285 this->couplingManager().updateCouplingContext(domainI, *this, domainI, scv.dofIndex(), elemSol[scv.localDofIndex()], pvIdx);
286 }
287
291 template<class JacobianMatrixDiagBlock, class GridVariables>
292 void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector& origResiduals, const JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
293 {
294 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *this, origResiduals, A, gridVariables);
295 }
296
303 template<std::size_t otherId, class JacobianBlock, class GridVariables>
304 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
305 const ElementResidualVector& res, GridVariables& gridVariables)
306 {
308 // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. //
310
311 // get some aliases for convenience
312 const auto& element = this->element();
313 const auto& fvGeometry = this->fvGeometry();
314 auto&& curElemVolVars = this->curElemVolVars();
315 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
316
317 // convenience lambda for call to update self
318 auto updateCoupledVariables = [&] ()
319 {
320 // Update ourself after the context has been modified. Depending on the
321 // type of caching, other objects might have to be updated. All ifs can be optimized away.
322 if constexpr (enableGridFluxVarsCache)
323 {
324 if constexpr (enableGridVolVarsCache)
325 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
326 else
327 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, gridVariables.gridFluxVarsCache());
328 }
329 else
330 {
331 if constexpr (enableGridVolVarsCache)
332 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache);
333 else
334 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache);
335 }
336 };
337
338 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
339 const auto& curSolJ = this->curSol(domainJ);
340 for (const auto globalJ : stencil)
341 {
342 // undeflected privars and privars to be deflected
343 const auto origPriVarsJ = curSolJ[globalJ];
344 auto priVarsJ = origPriVarsJ;
345
346 // the undeflected coupling residual
347 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ);
348
349 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
350 {
351 auto evalCouplingResidual = [&](Scalar priVar)
352 {
353 priVarsJ[pvIdx] = priVar;
354 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
355 updateCoupledVariables();
356 return this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ);
357 };
358
359 // derive the residuals numerically
360 ElementResidualVector partialDerivs(fvGeometry.numScv());
361
362 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
363 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
364 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
365
366 NumericDifferentiation::partialDerivative(evalCouplingResidual, origPriVarsJ[pvIdx], partialDerivs, origResidual,
367 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
368
369 // update the global stiffness matrix with the current partial derivatives
370 for (const auto& scv : scvs(fvGeometry))
371 {
372 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
373 {
374 // A[i][col][eqIdx][pvIdx] is the rate of change of
375 // the residual of equation 'eqIdx' at dof 'i'
376 // depending on the primary variable 'pvIdx' at dof
377 // 'col'.
378 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
379
380 // If the dof is coupled by a Dirichlet condition,
381 // set the derived value only once (i.e. overwrite existing values).
382 if (this->elemBcTypes().hasDirichlet())
383 {
384 const auto bcTypes = this->elemBcTypes().get(fvGeometry, scv);
385 if (bcTypes.isCouplingDirichlet(eqIdx))
386 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][eqIdx];
387 else if (bcTypes.isDirichlet(eqIdx))
388 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
389 }
390
391 // enforce internal Dirichlet constraints
392 if constexpr (Problem::enableInternalDirichletConstraints())
393 {
394 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
395 if (internalDirichletConstraints[eqIdx])
396 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
397 }
398 }
399 }
400
401 // restore the current element solution
402 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
403
404 // restore the undeflected state of the coupling context
405 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
406 }
407
408 // Restore original state of the flux vars cache and/or vol vars.
409 // This has to be done in case they depend on variables of domainJ before
410 // we continue with the numeric derivative w.r.t the next globalJ. Otherwise,
411 // the next "origResidual" will be incorrect.
412 updateCoupledVariables();
413 }
414 }
415};
416
417} // end namespace Dumux
418
419#endif
An enum class to define various differentiation methods available in order to compute the derivatives...
An adapter class for local assemblers using numeric differentiation.
An assembler for Jacobian and residual contribution per element (pq1bubble method)
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A class for numeric differentiation.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
Definition: deprecated.hh:149
Definition: partialreassembler.hh:44
An assembler for Jacobian and residual contribution per element (PQ1Bubble methods)
Definition: pq1bubblelocalassembler.hh:304
static void partialDerivative(const Function &function, Scalar x0, FunctionEvalType &derivative, const FunctionEvalType &fx0, const int numericDifferenceMethod=1)
Computes the derivative of a function with respect to a function parameter.
Definition: numericdifferentiation.hh:61
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:60
A base class for PQ1Bubble local assemblers.
Definition: subdomainpq1bubblelocalassembler.hh:55
const Problem & problem(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: subdomainpq1bubblelocalassembler.hh:214
static constexpr auto domainId
export the domain id of this sub-domain
Definition: subdomainpq1bubblelocalassembler.hh:80
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: subdomainpq1bubblelocalassembler.hh:84
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: subdomainpq1bubblelocalassembler.hh:223
ElementResidualVector evalLocalSourceResidual(const Element &element, const ElementVolumeVariables &elemVolVars) const
Evaluates the local source term for an element and given element volume variables.
Definition: subdomainpq1bubblelocalassembler.hh:152
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacRow &jacRow, const ElementResidualVector &res, GridVariables &gridVariables)
Assemble the entries in a coupling block of the jacobian. There is no coupling block between a domain...
Definition: subdomainpq1bubblelocalassembler.hh:134
void assembleJacobianAndResidual(JacobianMatrixRow &jacRow, SubSolutionVector &res, GridVariablesTuple &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: subdomainpq1bubblelocalassembler.hh:110
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition: subdomainpq1bubblelocalassembler.hh:173
SubDomainPQ1BubbleLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: subdomainpq1bubblelocalassembler.hh:87
const auto & curSol(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: subdomainpq1bubblelocalassembler.hh:219
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: subdomainpq1bubblelocalassembler.hh:179
The PQ1Bubble scheme multidomain local assembler.
Definition: subdomainpq1bubblelocalassembler.hh:242
Control-volume fe staggered scheme multi domain local assembler using numeric differentiation and imp...
Definition: subdomainpq1bubblelocalassembler.hh:254
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, const JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: subdomainpq1bubblelocalassembler.hh:292
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const ElementResidualVector &res, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainpq1bubblelocalassembler.hh:304
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: subdomainpq1bubblelocalassembler.hh:277
void maybeUpdateCouplingContext(const SubControlVolume &scv, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: subdomainpq1bubblelocalassembler.hh:283
Declares all properties used in Dumux.