26#ifndef DUMUX_MULTIDOMAIN_CC_LOCAL_ASSEMBLER_HH
27#define DUMUX_MULTIDOMAIN_CC_LOCAL_ASSEMBLER_HH
29#include <dune/common/indices.hh>
30#include <dune/common/reservedvector.hh>
31#include <dune/common/hybridutilities.hh>
32#include <dune/grid/common/gridenums.hh>
33#include <dune/istl/matrixindexset.hh>
58template<std::
size_t id,
class TypeTag,
class Assembler,
class Implementation,
bool implicit>
66 using SolutionVector =
typename Assembler::SolutionVector;
71 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
72 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
73 using ElementFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
74 using Scalar =
typename GridVariables::Scalar;
76 using GridGeometry =
typename GridVariables::GridGeometry;
77 using FVElementGeometry =
typename GridGeometry::LocalView;
78 using SubControlVolume =
typename GridGeometry::SubControlVolume;
79 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
81 using GridView =
typename GridGeometry::GridView;
82 using Element =
typename GridView::template Codim<0>::Entity;
84 using CouplingManager =
typename Assembler::CouplingManager;
89 static constexpr auto domainId =
typename Dune::index_constant<id>();
93 using ParentType::ParentType;
98 const SolutionVector&
curSol,
116 template<
class JacobianMatrixRow,
class Gr
idVariablesTuple>
119 this->
asImp_().bindLocalViews();
122 const auto globalI = this->
fvGeometry().gridGeometry().elementMapper().index(this->
element());
123 res[globalI] = this->
asImp_().assembleJacobianAndResidualImplInverse(jacRow[
domainId], *std::get<domainId>(gridVariables));
126 using namespace Dune::Hybrid;
127 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](
auto&& i)
138 template<std::size_t otherId,
class JacRow,
class GridVariables,
139 typename std::enable_if_t<(otherId == id),
int> = 0>
141 const LocalResidualValues& res, GridVariables& gridVariables)
147 template<std::size_t otherId,
class JacRow,
class GridVariables,
148 typename std::enable_if_t<(otherId != id),
int> = 0>
150 const LocalResidualValues& res, GridVariables& gridVariables)
152 this->
asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
160 this->
asImp_().bindLocalViews();
161 const auto globalI = this->
fvGeometry().gridGeometry().elementMapper().index(this->
element());
171 ElementResidualVector residual(this->
fvGeometry().numScv());
177 const auto& curVolVars = elemVolVars[scv];
179 source *= -Extrusion::volume(scv)*curVolVars.extrusionFactor();
180 residual[scv.indexInElement()] = std::move(source);
204 const SubControlVolumeFace& scvf)
const
230 if (!this->
assembler().isStationaryProblem())
250 {
return couplingManager_; }
267template<std::
size_t id,
class TypeTag,
class Assembler, DiffMethod DM = DiffMethod::numeric,
bool implicit = true>
276template<std::
size_t id,
class TypeTag,
class Assembler>
289 using FVElementGeometry =
typename GridGeometry::LocalView;
290 using GridView =
typename GridGeometry::GridView;
291 using Element =
typename GridView::template Codim<0>::Entity;
294 enum { dim = GridView::dimension };
296 static constexpr bool enableGridFluxVarsCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
297 static constexpr bool enableGridVolVarsCache = getPropValue<TypeTag, Properties::EnableGridVolumeVariablesCache>();
298 static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize;
299 static constexpr auto domainI = Dune::index_constant<id>();
310 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
320 const auto& element = this->element();
321 const auto& fvGeometry = this->fvGeometry();
322 const auto& gridGeometry = fvGeometry.gridGeometry();
323 auto&& curElemVolVars = this->curElemVolVars();
324 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
327 const auto globalI = gridGeometry.elementMapper().index(element);
328 const auto& connectivityMap = gridGeometry.connectivityMap();
329 const auto numNeighbors = connectivityMap[globalI].size();
332 Dune::ReservedVector<Element, maxElementStencilSize> neighborElements;
333 neighborElements.resize(numNeighbors);
337 Residuals origResiduals(numNeighbors + 1); origResiduals = 0.0;
338 origResiduals[0] = this->evalLocalResidual()[0];
343 const auto evalNeighborFlux = [&] (
const auto& neighbor,
const auto& scvf)
345 if (neighbor.partitionType() == Dune::GhostEntity)
346 return LocalResidualValues(0.0);
348 return this->evalFluxResidual(neighbor, scvf);
354 for (
const auto& dataJ : connectivityMap[globalI])
356 neighborElements[j-1] = gridGeometry.element(dataJ.globalJ);
357 for (
const auto scvfIdx : dataJ.scvfsJ)
358 origResiduals[j] += evalNeighborFlux(neighborElements[j-1], fvGeometry.scvf(scvfIdx));
364 const auto& scv = fvGeometry.scv(globalI);
365 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
369 const auto& curSol = this->curSol()[domainI];
370 const auto origPriVars = curSol[globalI];
371 const auto origVolVars = curVolVars;
374 auto elemSol = elementSolution<FVElementGeometry>(origPriVars);
378 Residuals partialDerivs(numNeighbors + 1);
380 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
384 auto evalResiduals = [&](Scalar priVar)
386 Residuals partialDerivsTmp(numNeighbors + 1);
387 partialDerivsTmp = 0.0;
389 elemSol[0][pvIdx] = priVar;
390 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, globalI, elemSol[0], pvIdx);
391 curVolVars.update(elemSol, this->problem(), element, scv);
392 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
393 if (enableGridFluxVarsCache)
394 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
397 partialDerivsTmp[0] = this->evalLocalResidual()[0];
400 for (std::size_t k = 0; k < connectivityMap[globalI].size(); ++k)
401 for (
auto scvfIdx : connectivityMap[globalI][k].scvfsJ)
402 partialDerivsTmp[k+1] += evalNeighborFlux(neighborElements[k], fvGeometry.scvf(scvfIdx));
404 return partialDerivsTmp;
408 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
411 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
416 if (this->elementIsGhost())
418 partialDerivs[0] = 0.0;
419 partialDerivs[0][pvIdx] = 1.0;
423 if constexpr (Problem::enableInternalDirichletConstraints())
426 const auto internalDirichletConstraintsOwnElement = this->problem().hasInternalDirichletConstraint(this->element(), scv);
427 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
429 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
431 if (internalDirichletConstraintsOwnElement[eqIdx])
433 origResiduals[0][eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
434 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
437 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
442 for (
const auto& dataJ : connectivityMap[globalI])
444 const auto& neighborElement = neighborElements[j-1];
445 const auto& neighborScv = fvGeometry.scv(dataJ.globalJ);
446 const auto internalDirichletConstraintsNeighbor = this->problem().hasInternalDirichletConstraint(neighborElement, neighborScv);
448 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
450 if (internalDirichletConstraintsNeighbor[eqIdx])
451 A[dataJ.globalJ][globalI][eqIdx][pvIdx] = 0.0;
453 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j][eqIdx];
461 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
464 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
468 for (
const auto& dataJ : connectivityMap[globalI])
469 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j++][eqIdx];
474 curVolVars = origVolVars;
477 elemSol[0][pvIdx] = origPriVars[pvIdx];
480 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, globalI, elemSol[0], pvIdx);
489 if (enableGridFluxVarsCache)
490 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
493 typename ParentType::LocalResidual::ElementResidualVector origElementResidual({origResiduals[0]});
494 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *
this, origElementResidual, A, gridVariables);
497 return origResiduals[0];
504 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
506 const LocalResidualValues& res, GridVariables& gridVariables)
513 const auto& element = this->element();
514 const auto& fvGeometry = this->fvGeometry();
515 const auto& gridGeometry = fvGeometry.gridGeometry();
516 auto&& curElemVolVars = this->curElemVolVars();
517 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
520 const auto globalI = gridGeometry.elementMapper().index(element);
521 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
522 const auto& curSolJ = this->curSol()[domainJ];
525 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
528 auto updateCoupledVariables = [&] ()
532 if (enableGridFluxVarsCache)
534 if (enableGridVolVarsCache)
536 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
537 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), elemFluxVarsCache);
541 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, gridVariables.gridFluxVarsCache());
542 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, elemFluxVarsCache);
547 if (enableGridVolVarsCache)
548 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), elemFluxVarsCache);
550 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, elemFluxVarsCache);
554 for (
const auto& globalJ : stencil)
557 const auto origPriVarsJ = curSolJ[globalJ];
558 auto priVarsJ = origPriVarsJ;
560 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ)[0];
562 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
564 auto evalCouplingResidual = [&](Scalar priVar)
567 priVarsJ[pvIdx] = priVar;
568 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
569 updateCoupledVariables();
570 return this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ)[0];
574 LocalResidualValues partialDeriv(0.0);
575 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
576 static const int numDiffMethod = getParamFromGroup<int>(paramGroup,
"Assembly.NumericDifferenceMethod");
577 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
579 epsCoupl(priVarsJ[pvIdx], pvIdx), numDiffMethod);
582 if constexpr (Problem::enableInternalDirichletConstraints())
584 const auto& scv = fvGeometry.scv(globalI);
585 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
586 if (internalDirichletConstraints.none())
588 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
589 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
593 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
595 if (internalDirichletConstraints[eqIdx])
596 A[globalI][globalJ][eqIdx][pvIdx] = 0.0;
598 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
604 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
605 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
609 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
612 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
618 updateCoupledVariables();
629template<std::
size_t id,
class TypeTag,
class Assembler>
632 SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>, false>
642 static constexpr auto domainI = Dune::index_constant<id>();
657 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
660 if (this->assembler().isStationaryProblem())
661 DUNE_THROW(Dune::InvalidStateException,
"Using explicit jacobian assembler with stationary local residual");
664 const auto residual = this->evalLocalResidual()[0];
665 const auto storageResidual = this->evalLocalStorageResidual();
674 const auto& element = this->element();
675 const auto& fvGeometry = this->fvGeometry();
676 const auto& gridGeometry = fvGeometry.gridGeometry();
677 auto&& curElemVolVars = this->curElemVolVars();
680 const auto globalI = gridGeometry.elementMapper().index(element);
681 const auto& scv = fvGeometry.scv(globalI);
682 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
686 const auto& curSol = this->curSol()[domainI];
687 const auto origPriVars = curSol[globalI];
688 const auto origVolVars = curVolVars;
694 LocalResidualValues partialDeriv;
695 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
700 auto evalStorage = [&](Scalar priVar)
704 elemSol[0][pvIdx] = priVar;
705 curVolVars.update(elemSol, this->problem(), element, scv);
706 return this->evalLocalStorageResidual();
710 if (!this->elementIsGhost())
713 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
715 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
721 else partialDeriv[pvIdx] = 1.0;
725 if constexpr (Problem::enableInternalDirichletConstraints())
728 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
729 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
731 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
733 if (internalDirichletConstraints[eqIdx])
735 residual[eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
736 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
739 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
744 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
745 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
749 curVolVars = origVolVars;
752 elemSol[0][pvIdx] = origPriVars[pvIdx];
765 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
767 const LocalResidualValues& res, GridVariables& gridVariables)
776template<std::
size_t id,
class TypeTag,
class Assembler>
779 SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::analytic, true>, true>
786 using Element =
typename GridView::template Codim<0>::Entity;
790 enum { dim = GridView::dimension };
792 static constexpr auto domainI = Dune::index_constant<id>();
803 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
807 if (this->elementIsGhost())
809 const auto globalI = this->assembler().gridGeometry(domainI).elementMapper().index(this->element());
810 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
811 A[globalI][globalI][pvIdx][pvIdx] = 1.0;
814 return LocalResidualValues(0.0);
818 const auto& problem = this->problem();
819 const auto& element = this->element();
820 const auto& fvGeometry = this->fvGeometry();
821 const auto& curElemVolVars = this->curElemVolVars();
822 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
825 const auto globalI = this->assembler().gridGeometry(domainI).elementMapper().index(element);
826 const auto& scv = fvGeometry.scv(globalI);
827 const auto& volVars = curElemVolVars[scv];
830 if (!this->assembler().isStationaryProblem())
831 this->localResidual().addStorageDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv);
834 this->localResidual().addSourceDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv);
837 for (
const auto& scvf : scvfs(fvGeometry))
840 if (!scvf.boundary())
841 this->localResidual().addFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
846 const auto& bcTypes = problem.boundaryTypes(element, scvf);
849 if (bcTypes.hasDirichlet() && !bcTypes.hasNeumann())
850 this->localResidual().addCCDirichletFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
853 else if (bcTypes.hasNeumann() && !bcTypes.hasDirichlet())
854 this->localResidual().addRobinFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
857 DUNE_THROW(Dune::NotImplemented,
"Mixed boundary conditions. Use pure boundary conditions by converting Dirichlet BCs to Robin BCs");
861 if constexpr (Problem::enableInternalDirichletConstraints())
864 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
865 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
867 auto residual = this->evalLocalResidual()[0];
869 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
871 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
873 if (internalDirichletConstraints[eqIdx])
875 residual[eqIdx] = volVars.priVars()[eqIdx] - dirichletValues[eqIdx];
876 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
879 for (
const auto& scvf : scvfs(fvGeometry))
880 if (!scvf.boundary())
881 A[globalI][fvGeometry.scv(scvf.outsideScvIdx()).dofIndex()][eqIdx][pvIdx] = 0.0;
889 return this->evalLocalResidual()[0];
898 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
900 const LocalResidualValues& res, GridVariables& gridVariables)
907 const auto& element = this->element();
908 const auto& fvGeometry = this->fvGeometry();
909 const auto& gridGeometry = fvGeometry.gridGeometry();
910 auto&& curElemVolVars = this->curElemVolVars();
914 const auto globalI = gridGeometry.elementMapper().index(element);
915 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
917 for (
const auto globalJ : stencil)
919 const auto& elementJ = this->assembler().gridGeometry(domainJ).element(globalJ);
920 this->couplingManager().addCouplingDerivatives(A[globalI][globalJ], domainI, element, fvGeometry, curElemVolVars, domainJ, elementJ);
923 if constexpr (Problem::enableInternalDirichletConstraints())
925 const auto& scv = fvGeometry.scv(globalI);
926 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
927 if (internalDirichletConstraints.any())
929 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
930 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
931 if (internalDirichletConstraints[eqIdx])
932 A[globalI][globalJ][eqIdx][pvIdx] = 0.0;
An enum class to define various differentiation methods available in order to compute the derivatives...
An adapter class for local assemblers using numeric differentiation.
A class for numeric differentiation.
A arithmetic block vector type based on DUNE's reserved vector.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Element solution classes and factory functions.
Helper classes to compute the integration elements.
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
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition: diffmethod.hh:37
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
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
Implementation & asImp_()
Definition: fvlocalassemblerbase.hh:309
ElementVolumeVariables & prevElemVolVars()
The element volume variables of the provious time step.
Definition: fvlocalassemblerbase.hh:269
ElementResidualVector evalLocalResidual() const
Convenience function to evaluate the complete local residual for the current element....
Definition: fvlocalassemblerbase.hh:120
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
LocalResidual & localResidual()
The local residual for the current element.
Definition: fvlocalassemblerbase.hh:277
const Element & element() const
The current element.
Definition: fvlocalassemblerbase.hh:249
const SolutionVector & curSol() const
The current solution.
Definition: fvlocalassemblerbase.hh:257
A helper class for local assemblers using numeric differentiation to determine the epsilon.
Definition: numericepsilon.hh:41
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 arithmetic block vector type based on DUNE's reserved vector.
Definition: reservedblockvector.hh:38
Definition: multidomain/couplingmanager.hh:46
A base class for all multidomain local assemblers.
Definition: subdomaincclocalassembler.hh:60
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: subdomaincclocalassembler.hh:213
SubDomainCCLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
the constructor
Definition: subdomaincclocalassembler.hh:96
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition: subdomaincclocalassembler.hh:189
LocalResidualValues evalFluxResidual(const Element &neighbor, const SubControlVolumeFace &scvf) const
Evaluates the fluxes depending on the chose time discretization scheme.
Definition: subdomaincclocalassembler.hh:203
void assembleJacobianAndResidual(JacobianMatrixRow &jacRow, SubSolutionVector &res, GridVariablesTuple &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: subdomaincclocalassembler.hh:117
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacRow &jacRow, const LocalResidualValues &res, GridVariables &gridVariables)
Assemble the entries in a coupling block of the jacobian. There is no coupling block between a domain...
Definition: subdomaincclocalassembler.hh:140
const Problem & problem() const
return reference to the underlying problem
Definition: subdomaincclocalassembler.hh:245
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: subdomaincclocalassembler.hh:249
static constexpr auto domainId
export the domain id of this sub-domain
Definition: subdomaincclocalassembler.hh:89
GetPropType< TypeTag, Properties::LocalResidual > LocalResidual
the local residual type of this domain
Definition: subdomaincclocalassembler.hh:91
ElementResidualVector evalLocalSourceResidual(const Element &element, const ElementVolumeVariables &elemVolVars) const
Evaluates the local source term for an element and given element volume variables.
Definition: subdomaincclocalassembler.hh:168
LocalResidualValues evalLocalStorageResidual() const
Evaluates the storage terms within the element.
Definition: subdomaincclocalassembler.hh:195
void assembleResidual(SubSolutionVector &res)
Assemble the residual only.
Definition: subdomaincclocalassembler.hh:158
The cell-centered scheme multidomain local assembler.
Definition: subdomaincclocalassembler.hh:268
Cell-centered scheme multidomain local assembler using numeric differentiation and implicit time disc...
Definition: subdomaincclocalassembler.hh:280
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const LocalResidualValues &res, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomaincclocalassembler.hh:505
LocalResidualValues assembleJacobianAndResidualImplInverse(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomaincclocalassembler.hh:311
Cell-centered scheme multidomain local assembler using numeric differentiation and explicit time disc...
Definition: subdomaincclocalassembler.hh:633
LocalResidualValues assembleJacobianAndResidualImplInverse(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomaincclocalassembler.hh:658
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const LocalResidualValues &res, GridVariables &gridVariables)
Computes the coupling derivatives with respect to the given element and adds them to the global matri...
Definition: subdomaincclocalassembler.hh:766
Cell-centered scheme local assembler using analytic differentiation and implicit time discretization.
Definition: subdomaincclocalassembler.hh:780
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const LocalResidualValues &res, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomaincclocalassembler.hh:899
LocalResidualValues assembleJacobianAndResidualImpl(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomaincclocalassembler.hh:804
Declares all properties used in Dumux.