3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
subdomainstaggeredlocalassembler.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_STAGGERED_LOCAL_ASSEMBLER_HH
27#define DUMUX_MULTIDOMAIN_STAGGERED_LOCAL_ASSEMBLER_HH
28
29#include <dune/common/reservedvector.hh>
30#include <dune/common/indices.hh>
31#include <dune/common/hybridutilities.hh>
32#include <dune/grid/common/gridenums.hh> // for GhostEntity
33
42
43namespace Dumux {
44
56template<std::size_t id, class TypeTag, class Assembler, class Implementation, bool isImplicit = true>
57class SubDomainStaggeredLocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler,Implementation, isImplicit>
58{
60
63 using SolutionVector = typename Assembler::SolutionVector;
64
66 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
67 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
68 using Scalar = typename GridVariables::Scalar;
69
70 using ElementFaceVariables = typename GetPropType<TypeTag, Properties::GridFaceVariables>::LocalView;
71 using CellCenterResidualValue = typename LocalResidual::CellCenterResidualValue;
72 using FaceResidualValue = typename LocalResidual::FaceResidualValue;
73
74 using GridGeometry = typename GridVariables::GridGeometry;
75 using FVElementGeometry = typename GridGeometry::LocalView;
76 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
78 using Element = typename GridView::template Codim<0>::Entity;
79
80 using CouplingManager = typename Assembler::CouplingManager;
81
82 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
83
84public:
85 static constexpr auto domainId = typename Dune::index_constant<id>();
86 static constexpr auto cellCenterId = typename Dune::index_constant<0>();
87 static constexpr auto faceId = typename Dune::index_constant<1>();
88
89 static constexpr auto numEqCellCenter = CellCenterResidualValue::dimension;
90 static constexpr auto faceOffset = numEqCellCenter;
91
92 using ParentType::ParentType;
93
95 const Element& element,
96 const SolutionVector& curSol,
97 CouplingManager& couplingManager)
99 element,
100 curSol,
101 localView(assembler.gridGeometry(domainId)),
102 localView(assembler.gridVariables(domainId).curGridVolVars()),
103 localView(assembler.gridVariables(domainId).prevGridVolVars()),
104 localView(assembler.gridVariables(domainId).gridFluxVarsCache()),
106 (element.partitionType() == Dune::GhostEntity))
107 , curElemFaceVars_(localView(assembler.gridVariables(domainId).curGridFaceVars()))
108 , prevElemFaceVars_(localView(assembler.gridVariables(domainId).prevGridFaceVars()))
109 , couplingManager_(couplingManager)
110 {}
111
116 template<class JacobianMatrixRow, class SubSol, class GridVariablesTuple>
117 void assembleJacobianAndResidual(JacobianMatrixRow& jacRow, SubSol& res, GridVariablesTuple& gridVariables)
118 {
119 this->asImp_().bindLocalViews();
120 this->elemBcTypes().update(problem(), this->element(), this->fvGeometry());
121
122 assembleJacobianAndResidualImpl_(domainId, jacRow, res, gridVariables);
123 }
124
128 template<class SubSol>
129 void assembleResidual(SubSol& res)
130 {
131 this->asImp_().bindLocalViews();
132 this->elemBcTypes().update(problem(), this->element(), this->fvGeometry());
133
134 assembleResidualImpl_(domainId, res);
135 }
136
141 CellCenterResidualValue evalLocalResidualForCellCenter() const
142 {
143 if (!isImplicit)
144 if (this->assembler().isStationaryProblem())
145 DUNE_THROW(Dune::InvalidStateException, "Using explicit jacobian assembler with stationary local residual");
146
147 if (this->elementIsGhost())
148 {
149 return CellCenterResidualValue(0.0);
150 }
151
154 }
155
161 CellCenterResidualValue evalLocalResidualForCellCenter(const ElementVolumeVariables& elemVolVars,
162 const ElementFaceVariables& elemFaceVars) const
163 {
164 auto residual = evalLocalFluxAndSourceResidualForCellCenter(elemVolVars, elemFaceVars);
165
166 if (!this->assembler().isStationaryProblem())
168
169 // handle cells with a fixed Dirichlet value
170 const auto cellCenterGlobalI = problem().gridGeometry().elementMapper().index(this->element());
171 const auto& scvI = this->fvGeometry().scv(cellCenterGlobalI);
172 for (int pvIdx = 0; pvIdx < numEqCellCenter; ++pvIdx)
173 {
174 static constexpr auto offset = numEq - numEqCellCenter;
175 if (this->problem().isDirichletCell(this->element(), this->fvGeometry(), scvI, pvIdx + offset))
176 residual[pvIdx] = this->curSol()[cellCenterId][cellCenterGlobalI][pvIdx] - this->problem().dirichlet(this->element(), scvI)[pvIdx + offset];
177 }
178
179 return residual;
180 }
181
187 CellCenterResidualValue evalLocalFluxAndSourceResidualForCellCenter() const
188 {
191 }
192
200 CellCenterResidualValue evalLocalFluxAndSourceResidualForCellCenter(const ElementVolumeVariables& elemVolVars, const ElementFaceVariables& elemFaceVars) const
201 {
202 return this->localResidual().evalFluxAndSourceForCellCenter(this->element(), this->fvGeometry(), elemVolVars, elemFaceVars, this->elemBcTypes(), this->elemFluxVarsCache());
203 }
204
210 CellCenterResidualValue evalLocalStorageResidualForCellCenter() const
211 {
212 return this->localResidual().evalStorageForCellCenter(this->element(), this->fvGeometry(), this->prevElemVolVars(), this->curElemVolVars());
213 }
214
220 FaceResidualValue evalLocalResidualForFace(const SubControlVolumeFace& scvf) const
221 {
222 if (!isImplicit)
223 if (this->assembler().isStationaryProblem())
224 DUNE_THROW(Dune::InvalidStateException, "Using explicit jacobian assembler with stationary local residual");
225
226 if (this->elementIsGhost())
227 {
228 return FaceResidualValue(0.0);
229 }
230
233 }
234
241 FaceResidualValue evalLocalResidualForFace(const SubControlVolumeFace& scvf,
242 const ElementVolumeVariables& elemVolVars,
243 const ElementFaceVariables& elemFaceVars) const
244 {
245 auto residual = evalLocalFluxAndSourceResidualForFace(scvf, elemVolVars, elemFaceVars);
246
247 if (!this->assembler().isStationaryProblem())
248 residual += evalLocalStorageResidualForFace(scvf);
249
250 this->localResidual().evalDirichletBoundariesForFace(residual, this->problem(), this->element(),
251 this->fvGeometry(), scvf, elemVolVars, elemFaceVars,
252 this->elemBcTypes(), this->elemFluxVarsCache());
253
254 return residual;
255 }
256
263 FaceResidualValue evalLocalFluxAndSourceResidualForFace(const SubControlVolumeFace& scvf) const
264 {
267 }
268
276 FaceResidualValue evalLocalFluxAndSourceResidualForFace(const SubControlVolumeFace& scvf,
277 const ElementVolumeVariables& elemVolVars,
278 const ElementFaceVariables& elemFaceVars) const
279 {
280 return this->localResidual().evalFluxAndSourceForFace(this->element(), this->fvGeometry(), elemVolVars, elemFaceVars, this->elemBcTypes(), this->elemFluxVarsCache(), scvf);
281 }
282
289 FaceResidualValue evalLocalStorageResidualForFace(const SubControlVolumeFace& scvf) const
290 {
291 return this->localResidual().evalStorageForFace(this->element(), this->fvGeometry(), this->prevElemVolVars(), this->curElemVolVars(), this->prevElemFaceVars(), this->curElemFaceVars(), scvf);
292 }
293
294 const Problem& problem() const
295 { return this->assembler().problem(domainId); }
296
298 ElementFaceVariables& curElemFaceVars()
299 { return curElemFaceVars_; }
300
302 ElementFaceVariables& prevElemFaceVars()
303 { return prevElemFaceVars_; }
304
306 const ElementFaceVariables& curElemFaceVars() const
307 { return curElemFaceVars_; }
308
310 const ElementFaceVariables& prevElemFaceVars() const
311 { return prevElemFaceVars_; }
312
313 CouplingManager& couplingManager()
314 { return couplingManager_; }
315
316private:
317
319 template<class SubSol>
320 void assembleResidualImpl_(Dune::index_constant<0>, SubSol& res)
321 {
322 const auto cellCenterGlobalI = problem().gridGeometry().elementMapper().index(this->element());
323 res[cellCenterGlobalI] = this->asImp_().assembleCellCenterResidualImpl();
324 }
325
327 template<class SubSol>
328 void assembleResidualImpl_(Dune::index_constant<1>, SubSol& res)
329 {
330 for (auto&& scvf : scvfs(this->fvGeometry()))
331 res[scvf.dofIndex()] += this->asImp_().assembleFaceResidualImpl(scvf);
332 }
333
335 template<class JacobianMatrixRow, class SubSol, class GridVariablesTuple>
336 auto assembleJacobianAndResidualImpl_(Dune::index_constant<0>, JacobianMatrixRow& jacRow, SubSol& res, GridVariablesTuple& gridVariables)
337 {
338 auto& gridVariablesI = *std::get<domainId>(gridVariables);
339 const auto cellCenterGlobalI = problem().gridGeometry().elementMapper().index(this->element());
340 const auto residual = this->asImp_().assembleCellCenterJacobianAndResidualImpl(jacRow[domainId], gridVariablesI);
341 res[cellCenterGlobalI] = residual;
342
343
344 // for the coupling blocks
345 using namespace Dune::Hybrid;
346 static constexpr auto otherDomainIds = makeIncompleteIntegerSequence<JacobianMatrixRow::size(), domainId>{};
347 forEach(otherDomainIds, [&](auto&& domainJ)
348 {
349 this->asImp_().assembleJacobianCellCenterCoupling(domainJ, jacRow[domainJ], residual, gridVariablesI);
350 });
351
352 // handle cells with a fixed Dirichlet value
353 incorporateDirichletCells_(jacRow);
354 }
355
357 template<class JacobianMatrixRow, class SubSol, class GridVariablesTuple>
358 void assembleJacobianAndResidualImpl_(Dune::index_constant<1>, JacobianMatrixRow& jacRow, SubSol& res, GridVariablesTuple& gridVariables)
359 {
360 auto& gridVariablesI = *std::get<domainId>(gridVariables);
361 const auto residual = this->asImp_().assembleFaceJacobianAndResidualImpl(jacRow[domainId], gridVariablesI);
362
363 for(auto&& scvf : scvfs(this->fvGeometry()))
364 res[scvf.dofIndex()] += residual[scvf.localFaceIdx()];
365
366 // for the coupling blocks
367 using namespace Dune::Hybrid;
368 static constexpr auto otherDomainIds = makeIncompleteIntegerSequence<JacobianMatrixRow::size(), domainId>{};
369 forEach(otherDomainIds, [&](auto&& domainJ)
370 {
371 this->asImp_().assembleJacobianFaceCoupling(domainJ, jacRow[domainJ], residual, gridVariablesI);
372 });
373 }
374
376 template<class JacobianMatrixRow>
377 void incorporateDirichletCells_(JacobianMatrixRow& jacRow)
378 {
379 const auto cellCenterGlobalI = problem().gridGeometry().elementMapper().index(this->element());
380
381 // overwrite the partial derivative with zero in case a fixed Dirichlet BC is used
382 static constexpr auto offset = numEq - numEqCellCenter;
383 for (int eqIdx = 0; eqIdx < numEqCellCenter; ++eqIdx)
384 {
385 if (problem().isDirichletCell(this->element(), this->fvGeometry(), this->fvGeometry().scv(cellCenterGlobalI), eqIdx + offset))
386 {
387 using namespace Dune::Hybrid;
388 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&, domainId = domainId](auto&& i)
389 {
390 auto& ccRowI = jacRow[i][cellCenterGlobalI];
391 for (auto col = ccRowI.begin(); col != ccRowI.end(); ++col)
392 {
393 ccRowI[col.index()][eqIdx] = 0.0;
394 // set the diagonal entry to 1.0
395 if ((i == domainId) && (col.index() == cellCenterGlobalI))
396 ccRowI[col.index()][eqIdx][eqIdx] = 1.0;
397 }
398 });
399 }
400 }
401 }
402
403 ElementFaceVariables curElemFaceVars_;
404 ElementFaceVariables prevElemFaceVars_;
405 CouplingManager& couplingManager_;
406};
407
418template<std::size_t id, class TypeTag, class Assembler, class Implementation>
419class SubDomainStaggeredLocalAssemblerImplicitBase : public SubDomainStaggeredLocalAssemblerBase<id, TypeTag, Assembler, Implementation>
420{
422 static constexpr auto domainId = Dune::index_constant<id>();
423public:
424 using ParentType::ParentType;
425
427 {
428 // get some references for convenience
429 auto& couplingManager = this->couplingManager();
430 const auto& element = this->element();
431 const auto& curSol = this->curSol();
432 auto&& fvGeometry = this->fvGeometry();
433 auto&& curElemVolVars = this->curElemVolVars();
434 auto&& curElemFaceVars = this->curElemFaceVars();
435 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
436
437 // bind the caches
438 couplingManager.bindCouplingContext(domainId, element, this->assembler());
439 fvGeometry.bind(element);
443 if (!this->assembler().isStationaryProblem())
444 {
445 this->prevElemVolVars().bindElement(element, fvGeometry, this->assembler().prevSol());
446 this->prevElemFaceVars().bindElement(element, fvGeometry, this->assembler().prevSol());
447 }
448 }
449};
450
462template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric, bool implicit = true>
464
485template<std::size_t id, class TypeTag, class Assembler>
486class SubDomainStaggeredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>
487: public SubDomainStaggeredLocalAssemblerImplicitBase<id, TypeTag, Assembler,
488 SubDomainStaggeredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, true> >
489{
490 using ThisType = SubDomainStaggeredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>;
494 using CellCenterResidualValue = typename LocalResidual::CellCenterResidualValue;
495 using FaceResidualValue = typename LocalResidual::FaceResidualValue;
496 using Element = typename GetPropType<TypeTag, Properties::GridView>::template Codim<0>::Entity;
498 using ElementFaceVariables = typename GetPropType<TypeTag, Properties::GridFaceVariables>::LocalView;
499 using FaceVariables = typename ElementFaceVariables::FaceVariables;
501 using FVElementGeometry = typename GridGeometry::LocalView;
502 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
504 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
507
508 static constexpr bool enableGridFluxVarsCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
509 static constexpr int maxNeighbors = 4*(2*ModelTraits::dim());
510 static constexpr auto domainI = Dune::index_constant<id>();
511 static constexpr auto cellCenterId = typename Dune::index_constant<0>();
512 static constexpr auto faceId = typename Dune::index_constant<1>();
513
514 static constexpr auto numEq = ModelTraits::numEq();
515 static constexpr auto numEqCellCenter = CellCenterPrimaryVariables::dimension;
516 static constexpr auto numEqFace = FacePrimaryVariables::dimension;
517
518public:
520
521 CellCenterResidualValue assembleCellCenterResidualImpl()
522 {
523 return this->evalLocalResidualForCellCenter();
524 }
525
526 FaceResidualValue assembleFaceResidualImpl(const SubControlVolumeFace& scvf)
527 {
528 return this->evalLocalResidualForFace(scvf);
529 }
530
537 template<class JacobianMatrixDiagBlock, class GridVariables>
538 CellCenterResidualValue assembleCellCenterJacobianAndResidualImpl(JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
539 {
540 assert(domainI == cellCenterId);
541
542 // get some aliases for convenience
543 const auto& element = this->element();
544 const auto& fvGeometry = this->fvGeometry();
545 auto&& curElemVolVars = this->curElemVolVars();
546 const auto& gridGeometry = this->problem().gridGeometry();
547 const auto& curSol = this->curSol()[domainI];
548
549 const auto cellCenterGlobalI = gridGeometry.elementMapper().index(element);
550 const auto origResidual = this->evalLocalResidualForCellCenter();
551
553 // Calculate derivatives of all cell center residuals in the element w.r.t. to other cell center dofs. //
555
556 // lambda to evaluate the derivatives for cell center dofs with respect to neighbor cells
557 auto evaluateCellCenterDerivatives = [&](const std::size_t globalJ)
558 {
559 // get the volVars of the element with respect to which we are going to build the derivative
560 auto&& scvJ = fvGeometry.scv(globalJ);
561 const auto elementJ = fvGeometry.gridGeometry().element(globalJ);
562 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scvJ);
563 const auto origVolVars(curVolVars);
564
565 for (int pvIdx = 0; pvIdx < numEqCellCenter; ++pvIdx)
566 {
567 CellCenterPrimaryVariables cellCenterPriVars = curSol[globalJ];
568 using PrimaryVariables = typename VolumeVariables::PrimaryVariables;
569 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
570
571 constexpr auto offset = numEq - numEqCellCenter;
572
573 auto evalResidual = [&](Scalar priVar)
574 {
575 // update the volume variables
576 priVars[pvIdx + offset] = priVar;
577 auto elemSol = elementSolution<FVElementGeometry>(std::move(priVars));
578 curVolVars.update(elemSol, this->problem(), elementJ, scvJ);
579
580 // update the coupling context
581 cellCenterPriVars[pvIdx] = priVar;
582 this->couplingManager().updateCouplingContext(domainI, *this, domainI, globalJ, cellCenterPriVars, pvIdx);
583
584 // compute element residual
585 return this->evalLocalResidualForCellCenter();
586 };
587
588 // create the vector storing the partial derivatives
589 CellCenterResidualValue partialDeriv(0.0);
590
591 // derive the residuals numerically
592 const auto& paramGroup = this->problem().paramGroup();
593 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
594 static const auto eps = this->couplingManager().numericEpsilon(domainI, paramGroup);
595 NumericDifferentiation::partialDerivative(evalResidual, priVars[pvIdx + offset], partialDeriv, origResidual,
596 eps(priVars[pvIdx + offset], pvIdx), numDiffMethod);
597
598 // update the global jacobian matrix with the current partial derivatives
599 updateGlobalJacobian_(A, cellCenterGlobalI, globalJ, pvIdx, partialDeriv);
600
601 // restore the original volVars
602 curVolVars = origVolVars;
603
604 // restore the undeflected state of the coupling context
605 this->couplingManager().updateCouplingContext(domainI, *this, domainI, globalJ, curSol[globalJ], pvIdx);
606 }
607 };
608
609 // get the list of cell center dofs that have an influence on the cell center resdiual of the current element
610 const auto& connectivityMap = gridGeometry.connectivityMap();
611
612 // evaluate derivatives w.r.t. own dof
613 evaluateCellCenterDerivatives(cellCenterGlobalI);
614
615 // evaluate derivatives w.r.t. all other related cell center dofs
616 for (const auto& globalJ : connectivityMap(cellCenterId, cellCenterId, cellCenterGlobalI))
617 evaluateCellCenterDerivatives(globalJ);
618
619 return origResidual;
620 }
621
628 template<class JacobianMatrixDiagBlock, class GridVariables>
629 auto assembleFaceJacobianAndResidualImpl(JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
630 {
631 assert(domainI == faceId);
632
633 // get some aliases for convenience
634 const auto& problem = this->problem();
635 const auto& element = this->element();
636 const auto& fvGeometry = this->fvGeometry();
637 const auto& gridGeometry = this->problem().gridGeometry();
638 const auto& curSol = this->curSol()[domainI];
639
640 using FaceSolutionVector = GetPropType<TypeTag, Properties::FaceSolutionVector>; // TODO: use reserved vector
641 FaceSolutionVector origResiduals;
642 origResiduals.resize(fvGeometry.numScvf());
643 origResiduals = 0.0;
644
645 // treat the local residua of the face dofs:
646 for (auto&& scvf : scvfs(fvGeometry))
647 origResiduals[scvf.localFaceIdx()] = this->evalLocalResidualForFace(scvf);
648
650 // Calculate derivatives of all face residuals in the element w.r.t. to other face dofs. //
651 // Note that we do an element-wise assembly, therefore this is only the contribution of the //
652 // current element while the contribution of the element on the opposite side of the scvf will //
653 // be added separately. //
655
656 for (auto&& scvf : scvfs(fvGeometry))
657 {
658 // set the actual dof index
659 const auto faceGlobalI = scvf.dofIndex();
660
662 const auto origFaceSolution = FaceSolution(scvf, curSol, gridGeometry);
663
664 // Lambda to evaluate the derivatives for faces
665 auto evaluateFaceDerivatives = [&](const std::size_t globalJ)
666 {
667 // get the faceVars of the face with respect to which we are going to build the derivative
668 auto& faceVars = getFaceVarAccess(gridVariables.curGridFaceVars(), this->curElemFaceVars(), scvf);
669 const auto origFaceVars = faceVars;
670
671 for (int pvIdx = 0; pvIdx < numEqFace; ++pvIdx)
672 {
673 auto faceSolution = origFaceSolution;
674
675 auto evalResidual = [&](Scalar priVar)
676 {
677 // update the volume variables
678 faceSolution[globalJ][pvIdx] = priVar;
679 faceVars.update(faceSolution, problem, element, fvGeometry, scvf);
680
681 // update the coupling context
682 this->couplingManager().updateCouplingContext(domainI, *this, domainI, globalJ, faceSolution[globalJ], pvIdx);
683
684 // compute face residual
685 return this->evalLocalResidualForFace(scvf);
686 };
687
688 // derive the residuals numerically
689 FaceResidualValue partialDeriv(0.0);
690 const auto& paramGroup = problem.paramGroup();
691 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
692 static const auto eps = this->couplingManager().numericEpsilon(domainI, paramGroup);
693 NumericDifferentiation::partialDerivative(evalResidual, faceSolution[globalJ][pvIdx], partialDeriv, origResiduals[scvf.localFaceIdx()],
694 eps(faceSolution[globalJ][pvIdx], pvIdx), numDiffMethod);
695
696 // update the global jacobian matrix with the current partial derivatives
697 updateGlobalJacobian_(A, faceGlobalI, globalJ, pvIdx, partialDeriv);
698
699 // restore the original faceVars
700 faceVars = origFaceVars;
701
702 // restore the undeflected state of the coupling context
703 this->couplingManager().updateCouplingContext(domainI, *this, domainI, globalJ, origFaceSolution[globalJ], pvIdx);
704 }
705 };
706
707 // evaluate derivatives w.r.t. own dof
708 evaluateFaceDerivatives(scvf.dofIndex());
709
710 // get the list of face dofs that have an influence on the resdiual of the current face
711 const auto& connectivityMap = gridGeometry.connectivityMap();
712
713 // evaluate derivatives w.r.t. all other related face dofs
714 for (const auto& globalJ : connectivityMap(faceId, faceId, scvf.index()))
715 evaluateFaceDerivatives(globalJ);
716 }
717
718 return origResiduals;
719 }
720
727 template<class JacobianBlock, class GridVariables>
728 void assembleJacobianCellCenterCoupling(Dune::index_constant<1> domainJ, JacobianBlock& A,
729 const CellCenterResidualValue& origResidual, GridVariables& gridVariables)
730 {
732 // Calculate derivatives of all cell center residuals in the element w.r.t. to all coupled faces dofs //
734
735 // get some aliases for convenience
736 const auto& element = this->element();
737 const auto& fvGeometry = this->fvGeometry();
738 const auto& gridGeometry = this->problem().gridGeometry();
739 const auto& curSol = this->curSol()[domainJ];
740 // build derivatives with for cell center dofs w.r.t. cell center dofs
741 const auto cellCenterGlobalI = gridGeometry.elementMapper().index(element);
742
743 for (const auto& scvfJ : scvfs(fvGeometry))
744 {
745 const auto globalJ = scvfJ.dofIndex();
746
747 // get the faceVars of the face with respect to which we are going to build the derivative
748 auto& faceVars = getFaceVarAccess(gridVariables.curGridFaceVars(), this->curElemFaceVars(), scvfJ);
749 const auto origFaceVars(faceVars);
750
751 for (int pvIdx = 0; pvIdx < numEqFace; ++pvIdx)
752 {
753 auto facePriVars = curSol[globalJ];
754
755 auto evalResidual = [&](Scalar priVar)
756 {
757 // update the face variables
758 facePriVars[pvIdx] = priVar;
759 faceVars.updateOwnFaceOnly(facePriVars);
760
761 // update the coupling context
762 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, facePriVars, pvIdx);
763
764 // compute element residual
765 return this->evalLocalResidualForCellCenter();
766 };
767
768 // create the vector storing the partial derivatives
769 CellCenterResidualValue partialDeriv(0.0);
770
771 // derive the residuals numerically
772 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
773 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
774 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
775 NumericDifferentiation::partialDerivative(evalResidual, facePriVars[pvIdx], partialDeriv, origResidual,
776 epsCoupl(facePriVars[pvIdx], pvIdx), numDiffMethod);
777
778 // update the global jacobian matrix with the current partial derivatives
779 updateGlobalJacobian_(A, cellCenterGlobalI, globalJ, pvIdx, partialDeriv);
780
781 // restore the original faceVars
782 faceVars = origFaceVars;
783
784 // restore the undeflected state of the coupling context
785 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, curSol[globalJ], pvIdx);
786 }
787 }
788 }
789
790 template<std::size_t otherId, class JacobianBlock, class GridVariables>
791 void assembleJacobianCellCenterCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
792 const CellCenterResidualValue& res, GridVariables& gridVariables)
793 {
795 // Calculate derivatives of all cell center residuals in the element w.r.t. all dofs in the coupling stencil. //
797
798 // get some aliases for convenience
799 const auto& element = this->element();
800
801 // get stencil informations
802 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
803
804 if (stencil.empty())
805 return;
806
807 for (const auto globalJ : stencil)
808 {
809 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ);
810 const auto& curSol = this->curSol()[domainJ];
811 const auto origPriVarsJ = curSol[globalJ];
812
813 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
814 {
815 auto evalCouplingResidual = [&](Scalar priVar)
816 {
817 auto deflectedPriVarsJ = origPriVarsJ;
818 deflectedPriVarsJ[pvIdx] = priVar;
819 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, deflectedPriVarsJ, pvIdx);
820 return this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ);
821 };
822
823 // create the vector storing the partial derivatives
824 CellCenterResidualValue partialDeriv(0.0);
825
826 // derive the residuals numerically
827 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
828 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
829 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
830 NumericDifferentiation::partialDerivative(evalCouplingResidual, origPriVarsJ[pvIdx], partialDeriv, origResidual,
831 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
832
833 // update the global stiffness matrix with the current partial derivatives
834 const auto cellCenterGlobalI = this->problem().gridGeometry().elementMapper().index(element);
835 updateGlobalJacobian_(A, cellCenterGlobalI, globalJ, pvIdx, partialDeriv);
836
837 // restore the undeflected state of the coupling context
838 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, origPriVarsJ, pvIdx);
839 }
840 }
841 }
842
849 template<class JacobianBlock, class ElementResidualVector, class GridVariables>
850 void assembleJacobianFaceCoupling(Dune::index_constant<0> domainJ, JacobianBlock& A,
851 const ElementResidualVector& origResiduals, GridVariables& gridVariables)
852 {
854 // Calculate derivatives of all face residuals in the element w.r.t. all coupled cell center dofs //
856
857 // get some aliases for convenience
858 const auto& problem = this->problem();
859 const auto& fvGeometry = this->fvGeometry();
860 const auto& gridGeometry = this->problem().gridGeometry();
861 const auto& connectivityMap = gridGeometry.connectivityMap();
862 const auto& curSol = this->curSol()[domainJ];
863
864 // build derivatives with for cell center dofs w.r.t. cell center dofs
865 for (auto&& scvf : scvfs(fvGeometry))
866 {
867 // set the actual dof index
868 const auto faceGlobalI = scvf.dofIndex();
869
870 // build derivatives with for face dofs w.r.t. cell center dofs
871 for (const auto& globalJ : connectivityMap(faceId, cellCenterId, scvf.index()))
872 {
873 // get the volVars of the element with respect to which we are going to build the derivative
874 auto&& scvJ = fvGeometry.scv(globalJ);
875 const auto elementJ = fvGeometry.gridGeometry().element(globalJ);
876 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), this->curElemVolVars(), scvJ);
877 const auto origVolVars(curVolVars);
878 const auto origCellCenterPriVars = curSol[globalJ];
879
880 for (int pvIdx = 0; pvIdx < numEqCellCenter; ++pvIdx)
881 {
882 using PrimaryVariables = typename VolumeVariables::PrimaryVariables;
883 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(origCellCenterPriVars);
884
885 constexpr auto offset = PrimaryVariables::dimension - CellCenterPrimaryVariables::dimension;
886
887 auto evalResidual = [&](Scalar priVar)
888 {
889 // update the volume variables
890 priVars[pvIdx + offset] = priVar;
891 auto elemSol = elementSolution<FVElementGeometry>(std::move(priVars));
892 curVolVars.update(elemSol, problem, elementJ, scvJ);
893
894 // update the coupling context
895 auto deflectedCellCenterPriVars = origCellCenterPriVars;
896 deflectedCellCenterPriVars[pvIdx] = priVar;
897 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, deflectedCellCenterPriVars, pvIdx);
898
899 // compute face residual
900 return this->evalLocalResidualForFace(scvf);
901 };
902
903 // derive the residuals numerically
904 FaceResidualValue partialDeriv(0.0);
905 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
906 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
907 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
908 NumericDifferentiation::partialDerivative(evalResidual, priVars[pvIdx + offset], partialDeriv, origResiduals[scvf.localFaceIdx()],
909 epsCoupl(priVars[pvIdx + offset], pvIdx), numDiffMethod);
910
911 // update the global jacobian matrix with the current partial derivatives
912 updateGlobalJacobian_(A, faceGlobalI, globalJ, pvIdx, partialDeriv);
913
914 // restore the original volVars
915 curVolVars = origVolVars;
916
917 // restore the undeflected state of the coupling context
918 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, origCellCenterPriVars, pvIdx);
919 }
920 }
921 }
922 }
923
924 template<std::size_t otherId, class JacobianBlock, class ElementResidualVector, class GridVariables>
925 void assembleJacobianFaceCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
926 const ElementResidualVector& res, GridVariables& gridVariables)
927 {
929 // Calculate derivatives of all face residuals in the element w.r.t. all dofs in the coupling stencil. //
931
932 // get some aliases for convenience
933 const auto& fvGeometry = this->fvGeometry();
934 const auto& curSol = this->curSol()[domainJ];
935
936 // build derivatives with for cell center dofs w.r.t. cell center dofs
937 for (auto&& scvf : scvfs(fvGeometry))
938 {
939 // set the actual dof index
940 const auto faceGlobalI = scvf.dofIndex();
941
942 // get stencil informations
943 const auto& stencil = this->couplingManager().couplingStencil(domainI, scvf, domainJ);
944
945 if (stencil.empty())
946 continue;
947
948 // build derivatives with for face dofs w.r.t. cell center dofs
949 for (const auto& globalJ : stencil)
950 {
951 const auto origPriVarsJ = curSol[globalJ];
952 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, scvf, *this, domainJ, globalJ);
953
954 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
955 {
956 auto evalCouplingResidual = [&](Scalar priVar)
957 {
958 auto deflectedPriVars = origPriVarsJ;
959 deflectedPriVars[pvIdx] = priVar;
960 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, deflectedPriVars, pvIdx);
961 return this->couplingManager().evalCouplingResidual(domainI, scvf, *this, domainJ, globalJ);
962 };
963
964 // derive the residuals numerically
965 FaceResidualValue partialDeriv(0.0);
966 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
967 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
968 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
969 NumericDifferentiation::partialDerivative(evalCouplingResidual, origPriVarsJ[pvIdx], partialDeriv, origResidual,
970 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
971
972 // update the global stiffness matrix with the current partial derivatives
973 updateGlobalJacobian_(A, faceGlobalI, globalJ, pvIdx, partialDeriv);
974
975 // restore the undeflected state of the coupling context
976 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, origPriVarsJ, pvIdx);
977 }
978 }
979 }
980 }
981
982 template<class JacobianMatrixDiagBlock, class GridVariables>
983 void evalAdditionalDerivatives(const std::vector<std::size_t>& additionalDofDependencies,
984 JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
985 { }
986
992 template<class SubMatrix, class CCOrFacePrimaryVariables>
993 static void updateGlobalJacobian_(SubMatrix& matrix,
994 const int globalI,
995 const int globalJ,
996 const int pvIdx,
997 const CCOrFacePrimaryVariables& partialDeriv)
998 {
999 for (int eqIdx = 0; eqIdx < partialDeriv.size(); eqIdx++)
1000 {
1001 // A[i][col][eqIdx][pvIdx] is the rate of change of
1002 // the residual of equation 'eqIdx' at dof 'i'
1003 // depending on the primary variable 'pvIdx' at dof
1004 // 'col'.
1005
1006 assert(pvIdx >= 0);
1007 assert(eqIdx < matrix[globalI][globalJ].size());
1008 assert(pvIdx < matrix[globalI][globalJ][eqIdx].size());
1009 matrix[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
1010 }
1011 }
1012
1013private:
1014
1015 template<class T = TypeTag>
1016 static typename std::enable_if<!getPropValue<T, Properties::EnableGridFaceVariablesCache>(), FaceVariables&>::type
1017 getFaceVarAccess(GridFaceVariables& gridFaceVariables, ElementFaceVariables& elemFaceVars, const SubControlVolumeFace& scvf)
1018 { return elemFaceVars[scvf]; }
1019
1020 template<class T = TypeTag>
1021 static typename std::enable_if<getPropValue<T, Properties::EnableGridFaceVariablesCache>(), FaceVariables&>::type
1022 getFaceVarAccess(GridFaceVariables& gridFaceVariables, ElementFaceVariables& elemFaceVars, const SubControlVolumeFace& scvf)
1023 { return gridFaceVariables.faceVars(scvf.index()); }
1024};
1025
1026} // end namespace Dumux
1027
1028#endif
An enum class to define various differentiation methods available in order to compute the derivatives...
Utilities for template meta programming.
A class for numeric differentiation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A arithmetic block vector type based on DUNE's reserved vector.
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
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
typename Detail::ConcatSeq< decltype(std::make_index_sequence< e >{}), e+1, decltype(std::make_index_sequence<(n > e) ?(n - e - 1) :0 >{})>::type makeIncompleteIntegerSequence
Definition: utility.hh:71
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Definition: common/properties/model.hh:34
constexpr double eps
epsilon for checking direction of scvf normals
Definition: test_tpfafvgeometry_nonconforming.cc:44
A base class for all local assemblers.
Definition: fvlocalassemblerbase.hh:48
ElementVolumeVariables & curElemVolVars()
The current element volume variables.
Definition: fvlocalassemblerbase.hh:263
ElementBoundaryTypes & elemBcTypes()
The element's boundary types.
Definition: fvlocalassemblerbase.hh:279
Implementation & asImp_()
Definition: fvlocalassemblerbase.hh:307
ElementVolumeVariables & prevElemVolVars()
The element volume variables of the provious time step.
Definition: fvlocalassemblerbase.hh:267
FVElementGeometry & fvGeometry()
The global finite volume geometry.
Definition: fvlocalassemblerbase.hh:259
const Assembler & assembler() const
The assembler.
Definition: fvlocalassemblerbase.hh:243
ElementFluxVariablesCache & elemFluxVarsCache()
The element flux variables cache.
Definition: fvlocalassemblerbase.hh:271
bool elementIsGhost() const
Returns if element is a ghost entity.
Definition: fvlocalassemblerbase.hh:251
LocalResidual & localResidual()
The local residual for the current element.
Definition: fvlocalassemblerbase.hh:275
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: fvlocalassemblerbase.hh:113
const Element & element() const
The current element.
Definition: fvlocalassemblerbase.hh:247
const SolutionVector & curSol() const
The current solution.
Definition: fvlocalassemblerbase.hh:255
static void partialDerivative(const Function &function, Scalar x0, FunctionEvalType &derivative, const FunctionEvalType &fx0, const int numericDifferenceMethod=1)
Computes the derivative of a function with repect to a function parameter.
Definition: numericdifferentiation.hh:61
A base class for all multidomain local assemblers (staggered)
Definition: subdomainstaggeredlocalassembler.hh:58
FaceResidualValue evalLocalResidualForFace(const SubControlVolumeFace &scvf) const
Convenience function to evaluate the local residual for the current face. Automatically chooses the t...
Definition: subdomainstaggeredlocalassembler.hh:220
const ElementFaceVariables & curElemFaceVars() const
The current element volume variables.
Definition: subdomainstaggeredlocalassembler.hh:306
FaceResidualValue evalLocalFluxAndSourceResidualForFace(const SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars) const
Evaluates the flux and source terms (i.e, the terms without a time derivative) of the local residual ...
Definition: subdomainstaggeredlocalassembler.hh:276
void assembleJacobianAndResidual(JacobianMatrixRow &jacRow, SubSol &res, GridVariablesTuple &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: subdomainstaggeredlocalassembler.hh:117
ElementFaceVariables & prevElemFaceVars()
The element volume variables of the provious time step.
Definition: subdomainstaggeredlocalassembler.hh:302
static constexpr auto domainId
Definition: subdomainstaggeredlocalassembler.hh:85
void assembleResidual(SubSol &res)
Assemble the residual only.
Definition: subdomainstaggeredlocalassembler.hh:129
const Problem & problem() const
Definition: subdomainstaggeredlocalassembler.hh:294
SubDomainStaggeredLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: subdomainstaggeredlocalassembler.hh:94
FaceResidualValue evalLocalFluxAndSourceResidualForFace(const SubControlVolumeFace &scvf) const
Convenience function to evaluate the flux and source terms (i.e, the terms without a time derivative)...
Definition: subdomainstaggeredlocalassembler.hh:263
CellCenterResidualValue evalLocalFluxAndSourceResidualForCellCenter() const
Convenience function to evaluate the flux and source terms (i.e, the terms without a time derivative)...
Definition: subdomainstaggeredlocalassembler.hh:187
static constexpr auto faceId
Definition: subdomainstaggeredlocalassembler.hh:87
static constexpr auto faceOffset
Definition: subdomainstaggeredlocalassembler.hh:90
FaceResidualValue evalLocalResidualForFace(const SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars) const
Evaluates the complete local residual for the current face.
Definition: subdomainstaggeredlocalassembler.hh:241
static constexpr auto numEqCellCenter
Definition: subdomainstaggeredlocalassembler.hh:89
const ElementFaceVariables & prevElemFaceVars() const
The element volume variables of the provious time step.
Definition: subdomainstaggeredlocalassembler.hh:310
CellCenterResidualValue evalLocalResidualForCellCenter() const
Convenience function to evaluate the complete local residual for the current element....
Definition: subdomainstaggeredlocalassembler.hh:141
FaceResidualValue evalLocalStorageResidualForFace(const SubControlVolumeFace &scvf) const
Convenience function to evaluate storage term (i.e, the term with a time derivative) of the local res...
Definition: subdomainstaggeredlocalassembler.hh:289
CouplingManager & couplingManager()
Definition: subdomainstaggeredlocalassembler.hh:313
static constexpr auto cellCenterId
Definition: subdomainstaggeredlocalassembler.hh:86
CellCenterResidualValue evalLocalResidualForCellCenter(const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars) const
Evaluates the complete local residual for the current cell center.
Definition: subdomainstaggeredlocalassembler.hh:161
CellCenterResidualValue evalLocalFluxAndSourceResidualForCellCenter(const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars) const
Evaluates the flux and source terms (i.e, the terms without a time derivative) of the local residual ...
Definition: subdomainstaggeredlocalassembler.hh:200
CellCenterResidualValue evalLocalStorageResidualForCellCenter() const
Convenience function to evaluate storage term (i.e, the term with a time derivative) of the local res...
Definition: subdomainstaggeredlocalassembler.hh:210
ElementFaceVariables & curElemFaceVars()
The current element volume variables.
Definition: subdomainstaggeredlocalassembler.hh:298
A base class for all implicit multidomain local assemblers (staggered)
Definition: subdomainstaggeredlocalassembler.hh:420
void bindLocalViews()
Definition: subdomainstaggeredlocalassembler.hh:426
The staggered multidomain local assembler.
Definition: subdomainstaggeredlocalassembler.hh:463
Staggered scheme local assembler using numeric differentiation and implicit time discretization.
Definition: subdomainstaggeredlocalassembler.hh:489
void assembleJacobianFaceCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const ElementResidualVector &res, GridVariables &gridVariables)
Definition: subdomainstaggeredlocalassembler.hh:925
FaceResidualValue assembleFaceResidualImpl(const SubControlVolumeFace &scvf)
Definition: subdomainstaggeredlocalassembler.hh:526
CellCenterResidualValue assembleCellCenterJacobianAndResidualImpl(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainstaggeredlocalassembler.hh:538
CellCenterResidualValue assembleCellCenterResidualImpl()
Definition: subdomainstaggeredlocalassembler.hh:521
static void updateGlobalJacobian_(SubMatrix &matrix, const int globalI, const int globalJ, const int pvIdx, const CCOrFacePrimaryVariables &partialDeriv)
Updates the current global Jacobian matrix with the partial derivatives of all equations in regard to...
Definition: subdomainstaggeredlocalassembler.hh:993
void assembleJacobianCellCenterCoupling(Dune::index_constant< 1 > domainJ, JacobianBlock &A, const CellCenterResidualValue &origResidual, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainstaggeredlocalassembler.hh:728
void evalAdditionalDerivatives(const std::vector< std::size_t > &additionalDofDependencies, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Definition: subdomainstaggeredlocalassembler.hh:983
auto assembleFaceJacobianAndResidualImpl(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainstaggeredlocalassembler.hh:629
void assembleJacobianFaceCoupling(Dune::index_constant< 0 > domainJ, JacobianBlock &A, const ElementResidualVector &origResiduals, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainstaggeredlocalassembler.hh:850
void assembleJacobianCellCenterCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const CellCenterResidualValue &res, GridVariables &gridVariables)
Definition: subdomainstaggeredlocalassembler.hh:791
Declares all properties used in Dumux.
The local element solution class for staggered methods.