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;
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 =
typename Dune::index_constant<0>();
87 static constexpr auto faceId =
typename Dune::index_constant<1>();
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();
134 assembleResidualImpl_(
domainId, res);
144 if (this->
assembler().isStationaryProblem())
145 DUNE_THROW(Dune::InvalidStateException,
"Using explicit jacobian assembler with stationary local residual");
149 return CellCenterResidualValue(0.0);
162 const ElementFaceVariables& elemFaceVars)
const
166 if (!this->
assembler().isStationaryProblem())
170 const auto cellCenterGlobalI =
problem().gridGeometry().elementMapper().index(this->
element());
171 const auto& scvI = this->
fvGeometry().scv(cellCenterGlobalI);
223 if (this->
assembler().isStationaryProblem())
224 DUNE_THROW(Dune::InvalidStateException,
"Using explicit jacobian assembler with stationary local residual");
228 return FaceResidualValue(0.0);
242 const ElementVolumeVariables& elemVolVars,
243 const ElementFaceVariables& elemFaceVars)
const
247 if (!this->
assembler().isStationaryProblem())
251 this->
fvGeometry(), scvf, elemVolVars, elemFaceVars,
277 const ElementVolumeVariables& elemVolVars,
278 const ElementFaceVariables& elemFaceVars)
const
299 {
return curElemFaceVars_; }
303 {
return prevElemFaceVars_; }
307 {
return curElemFaceVars_; }
311 {
return prevElemFaceVars_; }
314 {
return couplingManager_; }
319 template<
class SubSol>
320 void assembleResidualImpl_(Dune::index_constant<0>, SubSol& res)
322 const auto cellCenterGlobalI =
problem().gridGeometry().elementMapper().index(this->
element());
323 res[cellCenterGlobalI] = this->
asImp_().assembleCellCenterResidualImpl();
327 template<
class SubSol>
328 void assembleResidualImpl_(Dune::index_constant<1>, SubSol& res)
331 res[scvf.dofIndex()] += this->
asImp_().assembleFaceResidualImpl(scvf);
335 template<
class JacobianMatrixRow,
class SubSol,
class Gr
idVariablesTuple>
336 auto assembleJacobianAndResidualImpl_(Dune::index_constant<0>, JacobianMatrixRow& jacRow, SubSol& res, GridVariablesTuple& gridVariables)
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;
345 using namespace Dune::Hybrid;
347 forEach(otherDomainIds, [&](
auto&& domainJ)
349 this->
asImp_().assembleJacobianCellCenterCoupling(domainJ, jacRow[domainJ], residual, gridVariablesI);
353 incorporateDirichletCells_(jacRow);
357 template<
class JacobianMatrixRow,
class SubSol,
class Gr
idVariablesTuple>
358 void assembleJacobianAndResidualImpl_(Dune::index_constant<1>, JacobianMatrixRow& jacRow, SubSol& res, GridVariablesTuple& gridVariables)
360 auto& gridVariablesI = *std::get<domainId>(gridVariables);
361 const auto residual = this->
asImp_().assembleFaceJacobianAndResidualImpl(jacRow[
domainId], gridVariablesI);
364 res[scvf.dofIndex()] += residual[scvf.localFaceIdx()];
367 using namespace Dune::Hybrid;
369 forEach(otherDomainIds, [&](
auto&& domainJ)
371 this->
asImp_().assembleJacobianFaceCoupling(domainJ, jacRow[domainJ], residual, gridVariablesI);
376 template<
class JacobianMatrixRow>
377 void incorporateDirichletCells_(JacobianMatrixRow& jacRow)
379 const auto cellCenterGlobalI =
problem().gridGeometry().elementMapper().index(this->
element());
387 using namespace Dune::Hybrid;
388 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&,
domainId =
domainId](
auto&& i)
390 auto& ccRowI = jacRow[i][cellCenterGlobalI];
391 for (
auto col = ccRowI.begin(); col != ccRowI.end(); ++col)
393 ccRowI[col.index()][eqIdx] = 0.0;
395 if ((i ==
domainId) && (col.index() == cellCenterGlobalI))
396 ccRowI[col.index()][eqIdx][eqIdx] = 1.0;
403 ElementFaceVariables curElemFaceVars_;
404 ElementFaceVariables prevElemFaceVars_;
405 CouplingManager& couplingManager_;
418template<std::
size_t id,
class TypeTag,
class Assembler,
class Implementation>
422 static constexpr auto domainId = Dune::index_constant<id>();
424 using ParentType::ParentType;
443 if (!this->
assembler().isStationaryProblem())
462template<std::
size_t id,
class TypeTag,
class Assembler, DiffMethod DM = DiffMethod::numeric,
bool implicit = true>
485template<std::
size_t id,
class TypeTag,
class Assembler>
488 SubDomainStaggeredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, true> >
494 using CellCenterResidualValue =
typename LocalResidual::CellCenterResidualValue;
495 using FaceResidualValue =
typename LocalResidual::FaceResidualValue;
499 using FaceVariables =
typename ElementFaceVariables::FaceVariables;
501 using FVElementGeometry =
typename GridGeometry::LocalView;
502 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
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>();
514 static constexpr auto numEq = ModelTraits::numEq();
515 static constexpr auto numEqCellCenter = CellCenterPrimaryVariables::dimension;
516 static constexpr auto numEqFace = FacePrimaryVariables::dimension;
519 using ParentType::ParentType;
537 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
540 assert(domainI == cellCenterId);
546 const auto& gridGeometry = this->
problem().gridGeometry();
549 const auto cellCenterGlobalI = gridGeometry.elementMapper().index(
element);
557 auto evaluateCellCenterDerivatives = [&](
const std::size_t globalJ)
561 const auto elementJ =
fvGeometry.gridGeometry().element(globalJ);
563 const auto origVolVars(curVolVars);
565 for (
int pvIdx = 0; pvIdx < numEqCellCenter; ++pvIdx)
567 CellCenterPrimaryVariables cellCenterPriVars =
curSol[globalJ];
568 using PrimaryVariables =
typename VolumeVariables::PrimaryVariables;
571 constexpr auto offset = numEq - numEqCellCenter;
573 auto evalResidual = [&](Scalar priVar)
576 priVars[pvIdx + offset] = priVar;
578 curVolVars.update(elemSol, this->
problem(), elementJ, scvJ);
581 cellCenterPriVars[pvIdx] = priVar;
582 this->
couplingManager().updateCouplingContext(domainI, *
this, domainI, globalJ, cellCenterPriVars, pvIdx);
589 CellCenterResidualValue partialDeriv(0.0);
592 const auto& paramGroup = this->
problem().paramGroup();
594 static const auto eps = this->
couplingManager().numericEpsilon(domainI, paramGroup);
596 eps(priVars[pvIdx + offset], pvIdx), numDiffMethod);
602 curVolVars = origVolVars;
605 this->
couplingManager().updateCouplingContext(domainI, *
this, domainI, globalJ,
curSol[globalJ], pvIdx);
610 const auto& connectivityMap = gridGeometry.connectivityMap();
613 evaluateCellCenterDerivatives(cellCenterGlobalI);
616 for (
const auto& globalJ : connectivityMap(cellCenterId, cellCenterId, cellCenterGlobalI))
617 evaluateCellCenterDerivatives(globalJ);
628 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
631 assert(domainI == faceId);
637 const auto& gridGeometry = this->
problem().gridGeometry();
641 FaceSolutionVector origResiduals;
659 const auto faceGlobalI = scvf.dofIndex();
662 const auto origFaceSolution = FaceSolution(scvf,
curSol, gridGeometry);
665 auto evaluateFaceDerivatives = [&](
const std::size_t globalJ)
668 auto& faceVars = getFaceVarAccess(gridVariables.curGridFaceVars(), this->curElemFaceVars(), scvf);
669 const auto origFaceVars = faceVars;
671 for (
int pvIdx = 0; pvIdx < numEqFace; ++pvIdx)
673 auto faceSolution = origFaceSolution;
675 auto evalResidual = [&](Scalar priVar)
678 faceSolution[globalJ][pvIdx] = priVar;
682 this->
couplingManager().updateCouplingContext(domainI, *
this, domainI, globalJ, faceSolution[globalJ], pvIdx);
689 FaceResidualValue partialDeriv(0.0);
690 const auto& paramGroup =
problem.paramGroup();
692 static const auto eps = this->
couplingManager().numericEpsilon(domainI, paramGroup);
694 eps(faceSolution[globalJ][pvIdx], pvIdx), numDiffMethod);
700 faceVars = origFaceVars;
703 this->
couplingManager().updateCouplingContext(domainI, *
this, domainI, globalJ, origFaceSolution[globalJ], pvIdx);
708 evaluateFaceDerivatives(scvf.dofIndex());
711 const auto& connectivityMap = gridGeometry.connectivityMap();
714 for (
const auto& globalJ : connectivityMap(faceId, faceId, scvf.index()))
715 evaluateFaceDerivatives(globalJ);
718 return origResiduals;
727 template<
class JacobianBlock,
class Gr
idVariables>
729 const CellCenterResidualValue& origResidual, GridVariables& gridVariables)
738 const auto& gridGeometry = this->
problem().gridGeometry();
741 const auto cellCenterGlobalI = gridGeometry.elementMapper().index(
element);
745 const auto globalJ = scvfJ.dofIndex();
748 auto& faceVars = getFaceVarAccess(gridVariables.curGridFaceVars(), this->curElemFaceVars(), scvfJ);
749 const auto origFaceVars(faceVars);
751 for (
int pvIdx = 0; pvIdx < numEqFace; ++pvIdx)
753 auto facePriVars =
curSol[globalJ];
755 auto evalResidual = [&](Scalar priVar)
758 facePriVars[pvIdx] = priVar;
759 faceVars.updateOwnFaceOnly(facePriVars);
762 this->
couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, facePriVars, pvIdx);
769 CellCenterResidualValue partialDeriv(0.0);
772 const auto& paramGroup = this->
assembler().problem(domainJ).paramGroup();
774 static const auto epsCoupl = this->
couplingManager().numericEpsilon(domainJ, paramGroup);
776 epsCoupl(facePriVars[pvIdx], pvIdx), numDiffMethod);
782 faceVars = origFaceVars;
785 this->
couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ,
curSol[globalJ], pvIdx);
790 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
792 const CellCenterResidualValue& res, GridVariables& gridVariables)
807 for (
const auto globalJ : stencil)
809 const auto origResidual = this->
couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ);
811 const auto origPriVarsJ =
curSol[globalJ];
813 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
815 auto evalCouplingResidual = [&](Scalar priVar)
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);
824 CellCenterResidualValue partialDeriv(0.0);
827 const auto& paramGroup = this->
assembler().problem(domainJ).paramGroup();
829 static const auto epsCoupl = this->
couplingManager().numericEpsilon(domainJ, paramGroup);
831 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
834 const auto cellCenterGlobalI = this->
problem().gridGeometry().elementMapper().index(
element);
838 this->
couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, origPriVarsJ, pvIdx);
849 template<
class JacobianBlock,
class ElementRes
idualVector,
class Gr
idVariables>
860 const auto& gridGeometry = this->
problem().gridGeometry();
861 const auto& connectivityMap = gridGeometry.connectivityMap();
868 const auto faceGlobalI = scvf.dofIndex();
871 for (
const auto& globalJ : connectivityMap(faceId, cellCenterId, scvf.index()))
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];
880 for (
int pvIdx = 0; pvIdx < numEqCellCenter; ++pvIdx)
882 using PrimaryVariables =
typename VolumeVariables::PrimaryVariables;
885 constexpr auto offset = PrimaryVariables::dimension - CellCenterPrimaryVariables::dimension;
887 auto evalResidual = [&](Scalar priVar)
890 priVars[pvIdx + offset] = priVar;
892 curVolVars.update(elemSol,
problem, elementJ, scvJ);
895 auto deflectedCellCenterPriVars = origCellCenterPriVars;
896 deflectedCellCenterPriVars[pvIdx] = priVar;
897 this->
couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, deflectedCellCenterPriVars, pvIdx);
904 FaceResidualValue partialDeriv(0.0);
905 const auto& paramGroup = this->
assembler().problem(domainJ).paramGroup();
907 static const auto epsCoupl = this->
couplingManager().numericEpsilon(domainJ, paramGroup);
909 epsCoupl(priVars[pvIdx + offset], pvIdx), numDiffMethod);
915 curVolVars = origVolVars;
918 this->
couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, origCellCenterPriVars, pvIdx);
924 template<std::
size_t otherId,
class JacobianBlock,
class ElementRes
idualVector,
class Gr
idVariables>
940 const auto faceGlobalI = scvf.dofIndex();
943 const auto& stencil = this->
couplingManager().couplingStencil(domainI, scvf, domainJ);
949 for (
const auto& globalJ : stencil)
951 const auto origPriVarsJ =
curSol[globalJ];
952 const auto origResidual = this->
couplingManager().evalCouplingResidual(domainI, scvf, *
this, domainJ, globalJ);
954 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
956 auto evalCouplingResidual = [&](Scalar priVar)
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);
965 FaceResidualValue partialDeriv(0.0);
966 const auto& paramGroup = this->
assembler().problem(domainJ).paramGroup();
968 static const auto epsCoupl = this->
couplingManager().numericEpsilon(domainJ, paramGroup);
970 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
976 this->
couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, origPriVarsJ, pvIdx);
982 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
984 JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
992 template<
class SubMatrix,
class CCOrFacePrimaryVariables>
997 const CCOrFacePrimaryVariables& partialDeriv)
999 for (
int eqIdx = 0; eqIdx < partialDeriv.size(); eqIdx++)
1007 assert(eqIdx < matrix[globalI][globalJ].size());
1008 assert(pvIdx < matrix[globalI][globalJ][eqIdx].size());
1009 matrix[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
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]; }
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()); }
An enum class to define various differentiation methods available in order to compute the derivatives...
A base class for all local assemblers.
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
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethod::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition box/elementsolution.hh:115
PrimaryVariables makePriVarsFromCellCenterPriVars(const CellCenterPrimaryVariables &cellCenterPriVars)
Helper function to create a PrimaryVariables object from CellCenterPrimaryVariables.
Definition staggered/elementsolution.hh:41
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition diffmethod.hh:37
@ numeric
Definition diffmethod.hh:38
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:438
make the local view function available whenever we use the grid geometry
Definition adapt.hh:29
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:153
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
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:70
Definition common/properties/model.hh:34
ElementVolumeVariables & curElemVolVars()
Definition fvlocalassemblerbase.hh:263
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition fvlocalassemblerbase.hh:68
ElementBoundaryTypes & elemBcTypes()
Definition fvlocalassemblerbase.hh:279
Implementation & asImp_()
Definition fvlocalassemblerbase.hh:307
ElementVolumeVariables & prevElemVolVars()
Definition fvlocalassemblerbase.hh:267
FVLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol)
Definition fvlocalassemblerbase.hh:73
FVElementGeometry & fvGeometry()
Definition fvlocalassemblerbase.hh:259
const Assembler & assembler() const
Definition fvlocalassemblerbase.hh:243
ElementFluxVariablesCache & elemFluxVarsCache()
Definition fvlocalassemblerbase.hh:271
bool elementIsGhost() const
Definition fvlocalassemblerbase.hh:251
LocalResidual & localResidual()
Definition fvlocalassemblerbase.hh:275
static constexpr bool isImplicit()
Definition fvlocalassemblerbase.hh:113
const Element & element() const
Definition fvlocalassemblerbase.hh:247
VolumeVariables & getVolVarAccess(GridVolumeVariables &gridVolVars, ElementVolumeVariables &elemVolVars, const SubControlVolume &scv)
Definition fvlocalassemblerbase.hh:314
const SolutionVector & curSol() const
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
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
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.