3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
subdomainfclocalassembler.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_FACECENTERED_LOCAL_ASSEMBLER_HH
27#define DUMUX_MULTIDOMAIN_FACECENTERED_LOCAL_ASSEMBLER_HH
28
29#include <dune/common/indices.hh>
30#include <dune/common/hybridutilities.hh>
31#include <dune/grid/common/gridenums.hh>
32
39
40namespace Dumux {
41
53template<std::size_t id, class TypeTag, class Assembler, class Implementation, DiffMethod dm, bool implicit>
54class SubDomainFaceCenteredLocalAssemblerBase : public FaceCenteredLocalAssembler<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
87 explicit SubDomainFaceCenteredLocalAssemblerBase(const Assembler& assembler,
88 const Element& element,
89 const SolutionVector& curSol,
90 CouplingManager& couplingManager)
91 : ParentType(assembler,
92 element,
93 curSol,
94 localView(assembler.gridGeometry(domainId)),
95 localView(assembler.gridVariables(domainId).curGridVolVars()),
96 localView(assembler.gridVariables(domainId).prevGridVolVars()),
97 localView(assembler.gridVariables(domainId).gridFluxVarsCache()),
98 assembler.localResidual(domainId),
99 (element.partitionType() == Dune::GhostEntity))
100 , couplingManager_(couplingManager)
101 {}
102
107 template<class JacobianMatrixRow, class GridVariablesTuple>
108 void assembleJacobianAndResidual(JacobianMatrixRow& jacRow, SubSolutionVector& res, GridVariablesTuple& gridVariables)
109 {
110 auto assembleCouplingBlocks = [&](const auto& residual)
111 {
112 // assemble the coupling blocks
113 using namespace Dune::Hybrid;
114 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](auto&& i)
115 {
116 if (i != id)
117 this->assembleJacobianCoupling(i, jacRow, residual, gridVariables);
118 });
119 };
120
121 // the coupled model does not support partial reassembly yet
122 const DefaultPartialReassembler* noReassembler = nullptr;
123 ParentType::assembleJacobianAndResidual(jacRow[domainId], res, *std::get<domainId>(gridVariables), noReassembler, assembleCouplingBlocks);
124 }
125
130 template<std::size_t otherId, class JacRow, class GridVariables,
131 typename std::enable_if_t<(otherId == id), int> = 0>
132 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
133 const ElementResidualVector& res, GridVariables& gridVariables)
134 {}
135
139 template<std::size_t otherId, class JacRow, class GridVariables,
140 typename std::enable_if_t<(otherId != id), int> = 0>
141 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
142 const ElementResidualVector& res, GridVariables& gridVariables)
143 {
144 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
145 }
146
150 ElementResidualVector evalLocalSourceResidual(const Element& element, const ElementVolumeVariables& elemVolVars) const
151 {
152 // initialize the residual vector for all scvs in this element
153 ElementResidualVector residual(this->fvGeometry().numScv());
154
155 // evaluate the volume terms (storage + source terms)
156 // forward to the local residual specialized for the discretization methods
157 for (auto&& scv : scvs(this->fvGeometry()))
158 {
159 const auto& curVolVars = elemVolVars[scv];
160 auto source = this->localResidual().computeSource(problem(), element, this->fvGeometry(), elemVolVars, scv);
161 source *= -scv.volume()*curVolVars.extrusionFactor();
162 residual[scv.localDofIndex()] = std::move(source);
163 }
164
165 return residual;
166 }
167
171 ElementResidualVector evalLocalSourceResidual(const Element& neighbor) const
172 { return this->evalLocalSourceResidual(neighbor, implicit ? this->curElemVolVars() : this->prevElemVolVars()); }
173
178 {
179 // get some references for convenience
180 const auto& element = this->element();
181 const auto& curSol = this->curSol(domainId);
182 auto&& fvGeometry = this->fvGeometry();
183 auto&& curElemVolVars = this->curElemVolVars();
184 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
185
186 // bind the caches
187 couplingManager_.bindCouplingContext(domainId, element, this->assembler());
188 fvGeometry.bind(element);
189
190 if (implicit)
191 {
192 curElemVolVars.bind(element, fvGeometry, curSol);
193 elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars);
194 if (!this->assembler().isStationaryProblem())
195 this->prevElemVolVars().bindElement(element, fvGeometry, this->assembler().prevSol()[domainId]);
196 }
197 else
198 {
199 auto& prevElemVolVars = this->prevElemVolVars();
200 const auto& prevSol = this->assembler().prevSol()[domainId];
201
202 curElemVolVars.bindElement(element, fvGeometry, curSol);
203 prevElemVolVars.bind(element, fvGeometry, prevSol);
204 elemFluxVarsCache.bind(element, fvGeometry, prevElemVolVars);
205 }
206
207 this->elemBcTypes().update(problem(), this->element(), this->fvGeometry());
208 }
209
211 template<std::size_t i = domainId>
212 const Problem& problem(Dune::index_constant<i> dId = domainId) const
213 { return this->assembler().problem(domainId); }
214
216 template<std::size_t i = domainId>
217 const auto& curSol(Dune::index_constant<i> dId = domainId) const
218 { return ParentType::curSol()[dId]; }
219
221 CouplingManager& couplingManager()
222 { return couplingManager_; }
223
224private:
225 CouplingManager& couplingManager_;
226};
227
239template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric, bool implicit = true>
241
248template<std::size_t id, class TypeTag, class Assembler>
249class SubDomainFaceCenteredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>
251 id, TypeTag, Assembler,
253 DiffMethod::numeric, /*implicit=*/true
254>
255{
256 using ThisType = SubDomainFaceCenteredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>;
257 using ParentType = SubDomainFaceCenteredLocalAssemblerBase<id, TypeTag, Assembler, ThisType, DiffMethod::numeric, /*implicit=*/true>;
260
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
301 template<std::size_t otherId, class JacobianBlock, class GridVariables>
302 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
303 const ElementResidualVector& res, GridVariables& gridVariables)
304 {
306 // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. //
308
309 // get some aliases for convenience
310 const auto& element = this->element();
311 const auto& fvGeometry = this->fvGeometry();
312 auto&& curElemVolVars = this->curElemVolVars();
313 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
314
315 // the solution vector of the other domain
316 const auto& curSolJ = this->curSol(domainJ);
317
318 // convenience lambda for call to update self
319 auto updateCoupledVariables = [&] ()
320 {
321 // Update ourself after the context has been modified. Depending on the
322 // type of caching, other objects might have to be updated. All ifs can be optimized away.
323 if constexpr (enableGridFluxVarsCache)
324 {
325 if constexpr (enableGridVolVarsCache)
326 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
327 else
328 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, gridVariables.gridFluxVarsCache());
329 }
330 else
331 {
332 if constexpr (enableGridVolVarsCache)
333 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache);
334 else
335 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache);
336 }
337 };
338
339 for (const auto& scv : scvs(fvGeometry))
340 {
341 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, scv, domainJ);
342
343 for (const auto globalJ : stencil)
344 {
345 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, scv, domainJ, globalJ); // TODO is this necessary?
346 // undeflected privars and privars to be deflected
347 const auto origPriVarsJ = curSolJ[globalJ];
348 auto priVarsJ = origPriVarsJ;
349
350 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
351 {
352 auto evalCouplingResidual = [&](Scalar priVar)
353 {
354 priVarsJ[pvIdx] = priVar;
355 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
356 updateCoupledVariables();
357 return this->couplingManager().evalCouplingResidual(domainI, *this, scv, domainJ, globalJ);
358 };
359
360 // derive the residuals numerically
361 ElementResidualVector partialDerivs(element.subEntities(1));
362
363 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
364 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
365 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
366
367 NumericDifferentiation::partialDerivative(evalCouplingResidual, origPriVarsJ[pvIdx], partialDerivs, origResidual,
368 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
369
370 // update the global stiffness matrix with the current partial derivatives
371 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
372 {
373 // A[i][col][eqIdx][pvIdx] is the rate of change of
374 // the residual of equation 'eqIdx' at dof 'i'
375 // depending on the primary variable 'pvIdx' at dof
376 // 'col'.
377 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
378 }
379
380 // handle Dirichlet boundary conditions
381 // TODO internal constraints
382 if (scv.boundary() && this->elemBcTypes().hasDirichlet())
383 {
384 const auto bcTypes = this->elemBcTypes()[fvGeometry.frontalScvfOnBoundary(scv).localIndex()];
385 if (bcTypes.hasDirichlet())
386 {
387 // If the dof is coupled by a Dirichlet condition,
388 // set the derived value only once (i.e. overwrite existing values).
389 // For other dofs, add the contribution of the partial derivative.
390 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
391 {
392 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
393 {
394 if (bcTypes.isCouplingDirichlet(eqIdx))
395 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][pvIdx];
396 else if (bcTypes.isDirichlet(eqIdx))
397 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
398 }
399 }
400 }
401 }
402
403 // restore the current element solution
404 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
405
406 // restore the undeflected state of the coupling context
407 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
408 }
409
410 // Restore original state of the flux vars cache and/or vol vars.
411 // This has to be done in case they depend on variables of domainJ before
412 // we continue with the numeric derivative w.r.t the next globalJ. Otherwise,
413 // the next "origResidual" will be incorrect.
414 updateCoupledVariables();
415 }
416 }
417 }
418};
419
426template<std::size_t id, class TypeTag, class Assembler>
427class SubDomainFaceCenteredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>
429 id, TypeTag, Assembler,
430 SubDomainFaceCenteredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>,
431 DiffMethod::numeric, /*implicit=*/false
432>
433{
434 using ThisType = SubDomainFaceCenteredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>;
435 using ParentType = SubDomainFaceCenteredLocalAssemblerBase<id, TypeTag, Assembler, ThisType, DiffMethod::numeric, /*implicit=*/false>;
437 using FVElementGeometry = typename GridGeometry::LocalView;
438 using GridView = typename GridGeometry::GridView;
439 using Element = typename GridView::template Codim<0>::Entity;
440 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
441 static constexpr auto domainI = Dune::index_constant<id>();
442public:
446
450 template<class ElemSol>
451 void maybeUpdateCouplingContext(const SubControlVolume& scv, ElemSol& elemSol, const int pvIdx)
452 {
453 this->couplingManager().updateCouplingContext(domainI, *this, domainI, scv.dofIndex(), elemSol[scv.localDofIndex()], pvIdx);
454 }
455
459 template<class JacobianMatrixDiagBlock, class GridVariables>
460 void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector& origResiduals, const JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
461 {}
462
467 template<std::size_t otherId, class JacobianBlock, class GridVariables>
468 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
469 const ElementResidualVector& res, GridVariables& gridVariables)
470 {}
471};
472
473} // end namespace Dumux
474
475#endif
An assembler for Jacobian and residual contribution per element (face-centered staggered methods)
An enum class to define various differentiation methods available in order to compute the derivatives...
An adapter class for local assemblers using numeric differentiation.
A class for numeric differentiation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition: diffmethod.hh:37
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
An assembler for Jacobian and residual contribution per element (Face-centered methods)
Definition: fclocalassembler.hh:296
Definition: partialreassembler.hh:44
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 all face-centered staggered local assemblers.
Definition: subdomainfclocalassembler.hh:55
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition: subdomainfclocalassembler.hh:171
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: subdomainfclocalassembler.hh:221
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: subdomainfclocalassembler.hh:108
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: subdomainfclocalassembler.hh:177
const auto & curSol(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: subdomainfclocalassembler.hh:217
const Problem & problem(Dune::index_constant< i > dId=domainId) const
return reference to the underlying problem
Definition: subdomainfclocalassembler.hh:212
ElementResidualVector evalLocalSourceResidual(const Element &element, const ElementVolumeVariables &elemVolVars) const
Evaluates the local source term for an element and given element volume variables.
Definition: subdomainfclocalassembler.hh:150
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: subdomainfclocalassembler.hh:132
SubDomainFaceCenteredLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: subdomainfclocalassembler.hh:87
static constexpr auto domainId
export the domain id of this sub-domain
Definition: subdomainfclocalassembler.hh:80
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: subdomainfclocalassembler.hh:84
The face-centered staggered scheme multidomain local assembler.
Definition: subdomainfclocalassembler.hh:240
Face-centered staggered scheme multi domain local assembler using numeric differentiation and implici...
Definition: subdomainfclocalassembler.hh:255
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: subdomainfclocalassembler.hh:302
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: subdomainfclocalassembler.hh:277
void maybeUpdateCouplingContext(const SubControlVolume &scv, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: subdomainfclocalassembler.hh:283
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, const JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: subdomainfclocalassembler.hh:292
Face-centered staggered scheme multi domain local assembler using numeric differentiation and explici...
Definition: subdomainfclocalassembler.hh:433
typename ParentType::LocalResidual::ElementResidualVector ElementResidualVector
export element residual vector type
Definition: subdomainfclocalassembler.hh:445
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: subdomainfclocalassembler.hh:468
void maybeUpdateCouplingContext(const SubControlVolume &scv, ElemSol &elemSol, const int pvIdx)
Update the coupling context for coupled models.
Definition: subdomainfclocalassembler.hh:451
void maybeEvalAdditionalDomainDerivatives(const ElementResidualVector &origResiduals, const JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Update the additional domain derivatives for coupled models.
Definition: subdomainfclocalassembler.hh:460
Declares all properties used in Dumux.