26#ifndef DUMUX_MULTIDOMAIN_STAGGERED_LOCAL_ASSEMBLER_HH
27#define DUMUX_MULTIDOMAIN_STAGGERED_LOCAL_ASSEMBLER_HH
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>
56template<std::
size_t id,
class TypeTag,
class Assembler,
class Implementation,
bool isImplicit = true>
63 using SolutionVector =
typename Assembler::SolutionVector;
66 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
67 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
68 using Scalar =
typename GridVariables::Scalar;
71 using CellCenterResidualValue =
typename LocalResidual::CellCenterResidualValue;
72 using FaceResidualValue =
typename LocalResidual::FaceResidualValue;
74 using GridGeometry =
typename GridVariables::GridGeometry;
75 using FVElementGeometry =
typename GridGeometry::LocalView;
76 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
77 using GridView =
typename GridGeometry::GridView;
78 using Element =
typename GridView::template Codim<0>::Entity;
80 using CouplingManager =
typename Assembler::CouplingManager;
85 static constexpr auto domainId =
typename Dune::index_constant<id>();
86 static constexpr auto cellCenterId = GridGeometry::cellCenterIdx();
87 static constexpr auto faceId = GridGeometry::faceIdx();
92 using ParentType::ParentType;
96 const SolutionVector&
curSol,
116 template<
class JacobianMatrixRow,
class SubSol,
class Gr
idVariablesTuple>
119 this->
asImp_().bindLocalViews();
122 assembleJacobianAndResidualImpl_(
domainId, jacRow, res, gridVariables);
128 template<
class SubSol>
131 this->
asImp_().bindLocalViews();
136 const auto cellCenterGlobalI =
problem().gridGeometry().elementMapper().index(this->
element());
137 res[cellCenterGlobalI] = this->
asImp_().assembleCellCenterResidualImpl();
142 res[scvf.dofIndex()] += this->
asImp_().assembleFaceResidualImpl(scvf);
153 if (this->
assembler().isStationaryProblem())
154 DUNE_THROW(Dune::InvalidStateException,
"Using explicit jacobian assembler with stationary local residual");
158 return CellCenterResidualValue(0.0);
171 const ElementFaceVariables& elemFaceVars)
const
175 if (!this->
assembler().isStationaryProblem())
179 const auto cellCenterGlobalI =
problem().gridGeometry().elementMapper().index(this->
element());
180 const auto& scvI = this->
fvGeometry().scv(cellCenterGlobalI);
232 if (this->
assembler().isStationaryProblem())
233 DUNE_THROW(Dune::InvalidStateException,
"Using explicit jacobian assembler with stationary local residual");
237 return FaceResidualValue(0.0);
251 const ElementVolumeVariables& elemVolVars,
252 const ElementFaceVariables& elemFaceVars)
const
256 if (!this->
assembler().isStationaryProblem())
260 this->
fvGeometry(), scvf, elemVolVars, elemFaceVars,
286 const ElementVolumeVariables& elemVolVars,
287 const ElementFaceVariables& elemFaceVars)
const
308 {
return curElemFaceVars_; }
312 {
return prevElemFaceVars_; }
316 {
return curElemFaceVars_; }
320 {
return prevElemFaceVars_; }
323 {
return couplingManager_; }
328 template<
class JacobianMatrixRow,
class SubSol,
class Gr
idVariablesTuple>
329 auto assembleJacobianAndResidualImpl_(Dune::index_constant<cellCenterId>, JacobianMatrixRow& jacRow, SubSol& res, GridVariablesTuple& gridVariables)
331 auto& gridVariablesI = *std::get<domainId>(gridVariables);
332 const auto cellCenterGlobalI =
problem().gridGeometry().elementMapper().index(this->
element());
333 const auto residual = this->
asImp_().assembleCellCenterJacobianAndResidualImpl(jacRow[
domainId], gridVariablesI);
334 res[cellCenterGlobalI] = residual;
338 using namespace Dune::Hybrid;
340 forEach(otherDomainIds, [&](
auto&& domainJ)
342 this->
asImp_().assembleJacobianCellCenterCoupling(domainJ, jacRow[domainJ], residual, gridVariablesI);
346 incorporateDirichletCells_(jacRow);
350 template<
class JacobianMatrixRow,
class SubSol,
class Gr
idVariablesTuple>
351 void assembleJacobianAndResidualImpl_(Dune::index_constant<faceId>, JacobianMatrixRow& jacRow, SubSol& res, GridVariablesTuple& gridVariables)
353 auto& gridVariablesI = *std::get<domainId>(gridVariables);
354 const auto residual = this->
asImp_().assembleFaceJacobianAndResidualImpl(jacRow[
domainId], gridVariablesI);
357 res[scvf.dofIndex()] += residual[scvf.localFaceIdx()];
360 using namespace Dune::Hybrid;
362 forEach(otherDomainIds, [&](
auto&& domainJ)
364 this->
asImp_().assembleJacobianFaceCoupling(domainJ, jacRow[domainJ], residual, gridVariablesI);
369 template<
class JacobianMatrixRow>
370 void incorporateDirichletCells_(JacobianMatrixRow& jacRow)
372 const auto cellCenterGlobalI =
problem().gridGeometry().elementMapper().index(this->
element());
380 using namespace Dune::Hybrid;
381 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&,
domainId =
domainId](
auto&& i)
383 auto& ccRowI = jacRow[i][cellCenterGlobalI];
384 for (
auto col = ccRowI.begin(); col != ccRowI.end(); ++col)
386 ccRowI[col.index()][eqIdx] = 0.0;
388 if ((i == domainId) && (col.index() == cellCenterGlobalI))
389 ccRowI[col.index()][eqIdx][eqIdx] = 1.0;
396 ElementFaceVariables curElemFaceVars_;
397 ElementFaceVariables prevElemFaceVars_;
398 CouplingManager& couplingManager_;
411template<std::
size_t id,
class TypeTag,
class Assembler,
class Implementation>
415 static constexpr auto domainId = Dune::index_constant<id>();
417 using ParentType::ParentType;
436 if (!this->
assembler().isStationaryProblem())
455template<std::
size_t id,
class TypeTag,
class Assembler, DiffMethod DM = DiffMethod::numeric,
bool implicit = true>
478template<std::
size_t id,
class TypeTag,
class Assembler>
481 SubDomainStaggeredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, true> >
487 using CellCenterResidualValue =
typename LocalResidual::CellCenterResidualValue;
488 using FaceResidualValue =
typename LocalResidual::FaceResidualValue;
492 using FaceVariables =
typename ElementFaceVariables::FaceVariables;
494 using FVElementGeometry =
typename GridGeometry::LocalView;
495 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
501 static constexpr bool enableGridFluxVarsCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
502 static constexpr int maxNeighbors = 4*(2*ModelTraits::dim());
503 static constexpr auto domainI = Dune::index_constant<id>();
504 static constexpr auto cellCenterId = GridGeometry::cellCenterIdx();
505 static constexpr auto faceId = GridGeometry::faceIdx();
507 static constexpr auto numEq = ModelTraits::numEq();
508 static constexpr auto numEqCellCenter = CellCenterPrimaryVariables::dimension;
509 static constexpr auto numEqFace = FacePrimaryVariables::dimension;
516 return this->evalLocalResidualForCellCenter();
521 return this->evalLocalResidualForFace(scvf);
530 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
533 assert(domainI == cellCenterId);
536 const auto& element = this->element();
537 const auto& fvGeometry = this->fvGeometry();
538 auto&& curElemVolVars = this->curElemVolVars();
539 const auto& gridGeometry = this->problem().gridGeometry();
540 const auto& curSol = this->curSol()[domainI];
542 const auto cellCenterGlobalI = gridGeometry.elementMapper().index(element);
543 const auto origResidual = this->evalLocalResidualForCellCenter();
550 auto evaluateCellCenterDerivatives = [&](
const std::size_t globalJ)
553 auto&& scvJ = fvGeometry.scv(globalJ);
554 const auto elementJ = fvGeometry.gridGeometry().element(globalJ);
555 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scvJ);
556 const auto origVolVars(curVolVars);
558 for (
int pvIdx = 0; pvIdx < numEqCellCenter; ++pvIdx)
560 CellCenterPrimaryVariables cellCenterPriVars = curSol[globalJ];
561 using PrimaryVariables =
typename VolumeVariables::PrimaryVariables;
562 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
564 constexpr auto offset = numEq - numEqCellCenter;
566 auto evalResidual = [&](Scalar priVar)
569 priVars[pvIdx + offset] = priVar;
570 auto elemSol = elementSolution<FVElementGeometry>(std::move(priVars));
571 curVolVars.update(elemSol, this->problem(), elementJ, scvJ);
574 cellCenterPriVars[pvIdx] = priVar;
575 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, globalJ, cellCenterPriVars, pvIdx);
578 return this->evalLocalResidualForCellCenter();
582 CellCenterResidualValue partialDeriv(0.0);
585 const auto& paramGroup = this->problem().paramGroup();
586 static const int numDiffMethod = getParamFromGroup<int>(paramGroup,
"Assembly.NumericDifferenceMethod");
587 static const auto eps = this->couplingManager().numericEpsilon(domainI, paramGroup);
589 eps(priVars[pvIdx + offset], pvIdx), numDiffMethod);
592 updateGlobalJacobian_(A, cellCenterGlobalI, globalJ, pvIdx, partialDeriv);
595 curVolVars = origVolVars;
598 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, globalJ, curSol[globalJ], pvIdx);
603 const auto& connectivityMap = gridGeometry.connectivityMap();
606 evaluateCellCenterDerivatives(cellCenterGlobalI);
609 for (
const auto& globalJ : connectivityMap(cellCenterId, cellCenterId, cellCenterGlobalI))
610 evaluateCellCenterDerivatives(globalJ);
621 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
624 assert(domainI == faceId);
627 const auto& problem = this->problem();
628 const auto& element = this->element();
629 const auto& fvGeometry = this->fvGeometry();
630 const auto& gridGeometry = this->problem().gridGeometry();
631 const auto& curSol = this->curSol()[domainI];
634 FaceSolutionVector origResiduals;
635 origResiduals.resize(fvGeometry.numScvf());
639 for (
auto&& scvf : scvfs(fvGeometry))
640 origResiduals[scvf.localFaceIdx()] = this->evalLocalResidualForFace(scvf);
649 for (
auto&& scvf : scvfs(fvGeometry))
652 const auto faceGlobalI = scvf.dofIndex();
655 const auto origFaceSolution = FaceSolution(scvf, curSol, gridGeometry);
658 auto evaluateFaceDerivatives = [&](
const std::size_t globalJ)
661 auto& faceVars = getFaceVarAccess_(gridVariables.curGridFaceVars(), this->curElemFaceVars(), scvf);
662 const auto origFaceVars = faceVars;
664 for (
int pvIdx = 0; pvIdx < numEqFace; ++pvIdx)
666 auto faceSolution = origFaceSolution;
668 auto evalResidual = [&](Scalar priVar)
671 faceSolution[globalJ][pvIdx] = priVar;
672 faceVars.update(faceSolution, problem, element, fvGeometry, scvf);
675 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, globalJ, faceSolution[globalJ], pvIdx);
678 return this->evalLocalResidualForFace(scvf);
682 FaceResidualValue partialDeriv(0.0);
683 const auto& paramGroup = problem.paramGroup();
684 static const int numDiffMethod = getParamFromGroup<int>(paramGroup,
"Assembly.NumericDifferenceMethod");
685 static const auto eps = this->couplingManager().numericEpsilon(domainI, paramGroup);
687 eps(faceSolution[globalJ][pvIdx], pvIdx), numDiffMethod);
690 updateGlobalJacobian_(A, faceGlobalI, globalJ, pvIdx, partialDeriv);
693 faceVars = origFaceVars;
696 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, globalJ, origFaceSolution[globalJ], pvIdx);
701 evaluateFaceDerivatives(scvf.dofIndex());
704 const auto& connectivityMap = gridGeometry.connectivityMap();
707 for (
const auto& globalJ : connectivityMap(faceId, faceId, scvf.index()))
708 evaluateFaceDerivatives(globalJ);
711 return origResiduals;
720 template<
class JacobianBlock,
class Gr
idVariables>
722 const CellCenterResidualValue& origResidual, GridVariables& gridVariables)
729 const auto& element = this->element();
730 const auto& fvGeometry = this->fvGeometry();
731 const auto& gridGeometry = this->problem().gridGeometry();
732 const auto& curSol = this->curSol()[domainJ];
734 const auto cellCenterGlobalI = gridGeometry.elementMapper().index(element);
736 for (
const auto& scvfJ : scvfs(fvGeometry))
738 const auto globalJ = scvfJ.dofIndex();
741 auto& faceVars = getFaceVarAccess_(gridVariables.curGridFaceVars(), this->curElemFaceVars(), scvfJ);
742 const auto origFaceVars(faceVars);
744 for (
int pvIdx = 0; pvIdx < numEqFace; ++pvIdx)
746 auto facePriVars = curSol[globalJ];
748 auto evalResidual = [&](Scalar priVar)
751 facePriVars[pvIdx] = priVar;
752 faceVars.updateOwnFaceOnly(facePriVars);
755 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, facePriVars, pvIdx);
758 return this->evalLocalResidualForCellCenter();
762 CellCenterResidualValue partialDeriv(0.0);
765 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
766 static const int numDiffMethod = getParamFromGroup<int>(paramGroup,
"Assembly.NumericDifferenceMethod");
767 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
769 epsCoupl(facePriVars[pvIdx], pvIdx), numDiffMethod);
772 updateGlobalJacobian_(A, cellCenterGlobalI, globalJ, pvIdx, partialDeriv);
775 faceVars = origFaceVars;
778 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, curSol[globalJ], pvIdx);
783 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
785 const CellCenterResidualValue& res, GridVariables& gridVariables)
792 const auto& element = this->element();
795 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
800 for (
const auto globalJ : stencil)
802 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ);
803 const auto& curSol = this->curSol()[domainJ];
804 const auto origPriVarsJ = curSol[globalJ];
806 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
808 auto evalCouplingResidual = [&](Scalar priVar)
810 auto deflectedPriVarsJ = origPriVarsJ;
811 deflectedPriVarsJ[pvIdx] = priVar;
812 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, deflectedPriVarsJ, pvIdx);
813 return this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ);
817 CellCenterResidualValue partialDeriv(0.0);
820 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
821 static const int numDiffMethod = getParamFromGroup<int>(paramGroup,
"Assembly.NumericDifferenceMethod");
822 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
824 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
827 const auto cellCenterGlobalI = this->problem().gridGeometry().elementMapper().index(element);
828 updateGlobalJacobian_(A, cellCenterGlobalI, globalJ, pvIdx, partialDeriv);
831 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, origPriVarsJ, pvIdx);
842 template<
class JacobianBlock,
class ElementRes
idualVector,
class Gr
idVariables>
851 const auto& problem = this->problem();
852 const auto& fvGeometry = this->fvGeometry();
853 const auto& gridGeometry = this->problem().gridGeometry();
854 const auto& connectivityMap = gridGeometry.connectivityMap();
855 const auto& curSol = this->curSol()[domainJ];
858 for (
auto&& scvf : scvfs(fvGeometry))
861 const auto faceGlobalI = scvf.dofIndex();
864 for (
const auto& globalJ : connectivityMap(faceId, cellCenterId, scvf.index()))
867 auto&& scvJ = fvGeometry.scv(globalJ);
868 const auto elementJ = fvGeometry.gridGeometry().element(globalJ);
869 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), this->curElemVolVars(), scvJ);
870 const auto origVolVars(curVolVars);
871 const auto origCellCenterPriVars = curSol[globalJ];
873 for (
int pvIdx = 0; pvIdx < numEqCellCenter; ++pvIdx)
875 using PrimaryVariables =
typename VolumeVariables::PrimaryVariables;
876 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(origCellCenterPriVars);
878 constexpr auto offset = PrimaryVariables::dimension - CellCenterPrimaryVariables::dimension;
880 auto evalResidual = [&](Scalar priVar)
883 priVars[pvIdx + offset] = priVar;
884 auto elemSol = elementSolution<FVElementGeometry>(std::move(priVars));
885 curVolVars.update(elemSol, problem, elementJ, scvJ);
888 auto deflectedCellCenterPriVars = origCellCenterPriVars;
889 deflectedCellCenterPriVars[pvIdx] = priVar;
890 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, deflectedCellCenterPriVars, pvIdx);
893 return this->evalLocalResidualForFace(scvf);
897 FaceResidualValue partialDeriv(0.0);
898 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
899 static const int numDiffMethod = getParamFromGroup<int>(paramGroup,
"Assembly.NumericDifferenceMethod");
900 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
902 epsCoupl(priVars[pvIdx + offset], pvIdx), numDiffMethod);
905 updateGlobalJacobian_(A, faceGlobalI, globalJ, pvIdx, partialDeriv);
908 curVolVars = origVolVars;
911 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, origCellCenterPriVars, pvIdx);
917 template<std::
size_t otherId,
class JacobianBlock,
class ElementRes
idualVector,
class Gr
idVariables>
926 const auto& fvGeometry = this->fvGeometry();
927 const auto& curSol = this->curSol()[domainJ];
930 for (
auto&& scvf : scvfs(fvGeometry))
933 const auto faceGlobalI = scvf.dofIndex();
936 const auto& stencil = this->couplingManager().couplingStencil(domainI, scvf, domainJ);
942 for (
const auto& globalJ : stencil)
944 const auto origPriVarsJ = curSol[globalJ];
945 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, scvf, *
this, domainJ, globalJ);
947 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
949 auto evalCouplingResidual = [&](Scalar priVar)
951 auto deflectedPriVars = origPriVarsJ;
952 deflectedPriVars[pvIdx] = priVar;
953 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, deflectedPriVars, pvIdx);
954 return this->couplingManager().evalCouplingResidual(domainI, scvf, *
this, domainJ, globalJ);
958 FaceResidualValue partialDeriv(0.0);
959 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
960 static const int numDiffMethod = getParamFromGroup<int>(paramGroup,
"Assembly.NumericDifferenceMethod");
961 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
963 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
966 updateGlobalJacobian_(A, faceGlobalI, globalJ, pvIdx, partialDeriv);
969 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, origPriVarsJ, pvIdx);
975 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
977 JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
985 template<
class SubMatrix,
class CCOrFacePrimaryVariables>
990 const CCOrFacePrimaryVariables& partialDeriv)
992 for (
int eqIdx = 0; eqIdx < partialDeriv.size(); eqIdx++)
1000 assert(eqIdx < matrix[globalI][globalJ].size());
1001 assert(pvIdx < matrix[globalI][globalJ][eqIdx].size());
1002 matrix[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
1008 FaceVariables& getFaceVarAccess_(GridFaceVariables& gridFaceVariables, ElementFaceVariables& elemFaceVars,
const SubControlVolumeFace& scvf)
1010 if constexpr (getPropValue<TypeTag, Properties::EnableGridFaceVariablesCache>())
1011 return gridFaceVariables.faceVars(scvf.index());
1013 return elemFaceVars[scvf];
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.
A arithmetic block vector type based on DUNE's reserved vector.
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
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/pdesolver.hh:35
A base class for all local assemblers.
Definition: fvlocalassemblerbase.hh:48
ElementVolumeVariables & curElemVolVars()
The current element volume variables.
Definition: fvlocalassemblerbase.hh:265
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: fvlocalassemblerbase.hh:68
ElementBoundaryTypes & elemBcTypes()
The element's boundary types.
Definition: fvlocalassemblerbase.hh:281
Implementation & asImp_()
Definition: fvlocalassemblerbase.hh:309
ElementVolumeVariables & prevElemVolVars()
The element volume variables of the provious time step.
Definition: fvlocalassemblerbase.hh:269
FVElementGeometry & fvGeometry()
The global finite volume geometry.
Definition: fvlocalassemblerbase.hh:261
const Assembler & assembler() const
The assembler.
Definition: fvlocalassemblerbase.hh:245
ElementFluxVariablesCache & elemFluxVarsCache()
The element flux variables cache.
Definition: fvlocalassemblerbase.hh:273
bool elementIsGhost() const
Returns if element is a ghost entity.
Definition: fvlocalassemblerbase.hh:253
LocalResidual & localResidual()
The local residual for the current element.
Definition: fvlocalassemblerbase.hh:277
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:249
const SolutionVector & curSol() const
The current solution.
Definition: fvlocalassemblerbase.hh:257
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:229
const ElementFaceVariables & curElemFaceVars() const
The current element volume variables.
Definition: subdomainstaggeredlocalassembler.hh:315
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:285
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:311
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:303
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:272
CellCenterResidualValue evalLocalFluxAndSourceResidualForCellCenter() const
Convenience function to evaluate the flux and source terms (i.e, the terms without a time derivative)...
Definition: subdomainstaggeredlocalassembler.hh:196
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:250
static constexpr auto numEqCellCenter
Definition: subdomainstaggeredlocalassembler.hh:89
const ElementFaceVariables & prevElemFaceVars() const
The element volume variables of the provious time step.
Definition: subdomainstaggeredlocalassembler.hh:319
CellCenterResidualValue evalLocalResidualForCellCenter() const
Convenience function to evaluate the complete local residual for the current element....
Definition: subdomainstaggeredlocalassembler.hh:150
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:298
CouplingManager & couplingManager()
Definition: subdomainstaggeredlocalassembler.hh:322
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:170
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:209
CellCenterResidualValue evalLocalStorageResidualForCellCenter() const
Convenience function to evaluate storage term (i.e, the term with a time derivative) of the local res...
Definition: subdomainstaggeredlocalassembler.hh:219
ElementFaceVariables & curElemFaceVars()
The current element volume variables.
Definition: subdomainstaggeredlocalassembler.hh:307
A base class for all implicit multidomain local assemblers (staggered)
Definition: subdomainstaggeredlocalassembler.hh:413
void bindLocalViews()
Definition: subdomainstaggeredlocalassembler.hh:419
The staggered multidomain local assembler.
Definition: subdomainstaggeredlocalassembler.hh:456
Staggered scheme local assembler using numeric differentiation and implicit time discretization.
Definition: subdomainstaggeredlocalassembler.hh:482
void assembleJacobianFaceCoupling(Dune::index_constant< cellCenterId > 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:843
void assembleJacobianFaceCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const ElementResidualVector &res, GridVariables &gridVariables)
Definition: subdomainstaggeredlocalassembler.hh:918
FaceResidualValue assembleFaceResidualImpl(const SubControlVolumeFace &scvf)
Definition: subdomainstaggeredlocalassembler.hh:519
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:531
void assembleJacobianCellCenterCoupling(Dune::index_constant< faceId > 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:721
CellCenterResidualValue assembleCellCenterResidualImpl()
Definition: subdomainstaggeredlocalassembler.hh:514
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:986
void evalAdditionalDerivatives(const std::vector< std::size_t > &additionalDofDependencies, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Definition: subdomainstaggeredlocalassembler.hh:976
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:622
void assembleJacobianCellCenterCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const CellCenterResidualValue &res, GridVariables &gridVariables)
Definition: subdomainstaggeredlocalassembler.hh:784
Declares all properties used in Dumux.
The local element solution class for staggered methods.