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>
59template<std::
size_t id,
class TypeTag,
class Assembler,
class Implementation,
bool implicit>
67 using SolutionVector =
typename Assembler::SolutionVector;
72 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
73 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
74 using ElementFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
75 using Scalar =
typename GridVariables::Scalar;
77 using GridGeometry =
typename GridVariables::GridGeometry;
78 using FVElementGeometry =
typename GridGeometry::LocalView;
79 using SubControlVolume =
typename GridGeometry::SubControlVolume;
80 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
82 using GridView =
typename GridGeometry::GridView;
83 using Element =
typename GridView::template Codim<0>::Entity;
85 using CouplingManager =
typename Assembler::CouplingManager;
90 static constexpr auto domainId =
typename Dune::index_constant<id>();
94 using ParentType::ParentType;
99 const SolutionVector&
curSol,
117 template<
class JacobianMatrixRow,
class Gr
idVariablesTuple>
120 this->
asImp_().bindLocalViews();
123 const auto globalI = this->
fvGeometry().gridGeometry().elementMapper().index(this->
element());
124 res[globalI] = this->
asImp_().assembleJacobianAndResidualImplInverse(jacRow[
domainId], *std::get<domainId>(gridVariables));
127 using namespace Dune::Hybrid;
128 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](
auto&& i)
139 template<std::size_t otherId,
class JacRow,
class GridVariables,
140 typename std::enable_if_t<(otherId == id),
int> = 0>
142 const LocalResidualValues& res, GridVariables& gridVariables)
148 template<std::size_t otherId,
class JacRow,
class GridVariables,
149 typename std::enable_if_t<(otherId != id),
int> = 0>
151 const LocalResidualValues& res, GridVariables& gridVariables)
153 this->
asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
161 this->
asImp_().bindLocalViews();
162 const auto globalI = this->
fvGeometry().gridGeometry().elementMapper().index(this->
element());
165 if constexpr (Problem::enableInternalDirichletConstraints())
167 const auto applyDirichlet = [&] (
const auto& scvI,
168 const auto& dirichletValues,
172 res[scvI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
175 this->
asImp_().enforceInternalDirichletConstraints(applyDirichlet);
185 ElementResidualVector residual(this->
fvGeometry().numScv());
191 const auto& curVolVars = elemVolVars[scv];
194 residual[scv.indexInElement()] = std::move(source);
218 const SubControlVolumeFace& scvf)
const
244 if (!this->
assembler().isStationaryProblem())
264 {
return couplingManager_; }
281template<std::
size_t id,
class TypeTag,
class Assembler, DiffMethod DM = DiffMethod::numeric,
bool implicit = true>
290template<std::
size_t id,
class TypeTag,
class Assembler>
303 using FVElementGeometry =
typename GridGeometry::LocalView;
304 using GridView =
typename GridGeometry::GridView;
305 using Element =
typename GridView::template Codim<0>::Entity;
308 enum { dim = GridView::dimension };
312 static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize;
313 static constexpr auto domainI = Dune::index_constant<id>();
324 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
334 const auto& element = this->element();
335 const auto& fvGeometry = this->fvGeometry();
336 const auto& gridGeometry = fvGeometry.gridGeometry();
337 auto&& curElemVolVars = this->curElemVolVars();
338 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
341 const auto globalI = gridGeometry.elementMapper().index(element);
342 const auto& connectivityMap = gridGeometry.connectivityMap();
343 const auto numNeighbors = connectivityMap[globalI].size();
346 Dune::ReservedVector<Element, maxElementStencilSize> neighborElements;
347 neighborElements.resize(numNeighbors);
351 Residuals origResiduals(numNeighbors + 1); origResiduals = 0.0;
352 origResiduals[0] = this->evalLocalResidual()[0];
357 const auto evalNeighborFlux = [&] (
const auto& neighbor,
const auto& scvf)
359 if (neighbor.partitionType() == Dune::GhostEntity)
360 return LocalResidualValues(0.0);
362 return this->evalFluxResidual(neighbor, scvf);
368 for (
const auto& dataJ : connectivityMap[globalI])
370 neighborElements[j-1] = gridGeometry.element(dataJ.globalJ);
371 for (
const auto scvfIdx : dataJ.scvfsJ)
372 origResiduals[j] += evalNeighborFlux(neighborElements[j-1], fvGeometry.scvf(scvfIdx));
378 const auto& scv = fvGeometry.scv(globalI);
379 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
383 const auto& curSol = this->curSol()[domainI];
384 const auto origPriVars = curSol[globalI];
385 const auto origVolVars = curVolVars;
388 auto elemSol = elementSolution<FVElementGeometry>(origPriVars);
392 Residuals partialDerivs(numNeighbors + 1);
394 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
398 auto evalResiduals = [&](Scalar priVar)
400 Residuals partialDerivsTmp(numNeighbors + 1);
401 partialDerivsTmp = 0.0;
403 elemSol[0][pvIdx] = priVar;
404 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, globalI, elemSol[0], pvIdx);
405 curVolVars.update(elemSol, this->problem(), element, scv);
406 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
407 if (enableGridFluxVarsCache)
408 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
411 partialDerivsTmp[0] = this->evalLocalResidual()[0];
414 for (std::size_t k = 0; k < connectivityMap[globalI].size(); ++k)
415 for (
auto scvfIdx : connectivityMap[globalI][k].scvfsJ)
416 partialDerivsTmp[k+1] += evalNeighborFlux(neighborElements[k], fvGeometry.scvf(scvfIdx));
418 return partialDerivsTmp;
422 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
425 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
430 if (this->elementIsGhost())
432 partialDerivs[0] = 0.0;
433 partialDerivs[0][pvIdx] = 1.0;
437 curVolVars = origVolVars;
440 elemSol[0][pvIdx] = origPriVars[pvIdx];
443 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, globalI, elemSol[0], pvIdx);
446 if constexpr (Problem::enableInternalDirichletConstraints())
449 const auto internalDirichletConstraintsOwnElement = this->problem().hasInternalDirichletConstraint(this->element(), scv);
450 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
452 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
454 if (internalDirichletConstraintsOwnElement[eqIdx])
456 origResiduals[0][eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
457 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
460 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
465 for (
const auto& dataJ : connectivityMap[globalI])
467 const auto& neighborElement = neighborElements[j-1];
468 const auto& neighborScv = fvGeometry.scv(dataJ.globalJ);
469 const auto internalDirichletConstraintsNeighbor = this->problem().hasInternalDirichletConstraint(neighborElement, neighborScv);
471 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
473 if (internalDirichletConstraintsNeighbor[eqIdx])
474 A[dataJ.globalJ][globalI][eqIdx][pvIdx] = 0.0;
476 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j][eqIdx];
484 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
487 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
491 for (
const auto& dataJ : connectivityMap[globalI])
492 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j++][eqIdx];
503 if (enableGridFluxVarsCache)
504 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
507 typename ParentType::LocalResidual::ElementResidualVector origElementResidual({origResiduals[0]});
508 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *
this, origElementResidual, A, gridVariables);
511 return origResiduals[0];
518 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
520 const LocalResidualValues& res, GridVariables& gridVariables)
527 const auto& element = this->element();
528 const auto& fvGeometry = this->fvGeometry();
529 const auto& gridGeometry = fvGeometry.gridGeometry();
530 auto&& curElemVolVars = this->curElemVolVars();
531 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
534 const auto globalI = gridGeometry.elementMapper().index(element);
535 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
536 const auto& curSolJ = this->curSol()[domainJ];
539 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
542 auto updateCoupledVariables = [&] ()
546 if (enableGridFluxVarsCache)
548 if (enableGridVolVarsCache)
550 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
551 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), elemFluxVarsCache);
555 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, gridVariables.gridFluxVarsCache());
556 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, elemFluxVarsCache);
561 if (enableGridVolVarsCache)
562 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), elemFluxVarsCache);
564 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, elemFluxVarsCache);
568 for (
const auto& globalJ : stencil)
571 const auto origPriVarsJ = curSolJ[globalJ];
572 auto priVarsJ = origPriVarsJ;
574 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ)[0];
576 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
578 auto evalCouplingResidual = [&](Scalar priVar)
581 priVarsJ[pvIdx] = priVar;
582 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
583 updateCoupledVariables();
584 return this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ)[0];
588 LocalResidualValues partialDeriv(0.0);
589 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
590 static const int numDiffMethod = getParamFromGroup<int>(paramGroup,
"Assembly.NumericDifferenceMethod");
591 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
593 epsCoupl(priVarsJ[pvIdx], pvIdx), numDiffMethod);
596 if constexpr (Problem::enableInternalDirichletConstraints())
598 const auto& scv = fvGeometry.scv(globalI);
599 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
600 if (internalDirichletConstraints.none())
602 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
603 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
607 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
609 if (internalDirichletConstraints[eqIdx])
610 A[globalI][globalJ][eqIdx][pvIdx] = 0.0;
612 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
618 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
619 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
623 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
626 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
632 updateCoupledVariables();
643template<std::
size_t id,
class TypeTag,
class Assembler>
646 SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>, false>
656 static constexpr auto domainI = Dune::index_constant<id>();
671 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
674 if (this->assembler().isStationaryProblem())
675 DUNE_THROW(Dune::InvalidStateException,
"Using explicit jacobian assembler with stationary local residual");
678 const auto residual = this->evalLocalResidual()[0];
679 const auto storageResidual = this->evalLocalStorageResidual();
688 const auto& element = this->element();
689 const auto& fvGeometry = this->fvGeometry();
690 const auto& gridGeometry = fvGeometry.gridGeometry();
691 auto&& curElemVolVars = this->curElemVolVars();
694 const auto globalI = gridGeometry.elementMapper().index(element);
695 const auto& scv = fvGeometry.scv(globalI);
696 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
700 const auto& curSol = this->curSol()[domainI];
701 const auto origPriVars = curSol[globalI];
702 const auto origVolVars = curVolVars;
708 LocalResidualValues partialDeriv;
709 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
714 auto evalStorage = [&](Scalar priVar)
718 elemSol[0][pvIdx] = priVar;
719 curVolVars.update(elemSol, this->problem(), element, scv);
720 return this->evalLocalStorageResidual();
724 if (!this->elementIsGhost())
727 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
729 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
735 else partialDeriv[pvIdx] = 1.0;
738 curVolVars = origVolVars;
741 elemSol[0][pvIdx] = origPriVars[pvIdx];
745 if constexpr (Problem::enableInternalDirichletConstraints())
748 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
749 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
751 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
753 if (internalDirichletConstraints[eqIdx])
755 residual[eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
756 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
759 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
764 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
765 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
779 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
781 const LocalResidualValues& res, GridVariables& gridVariables)
791template<std::
size_t id,
class TypeTag,
class Assembler>
794 SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::analytic, true>, true>
801 using Element =
typename GridView::template Codim<0>::Entity;
805 enum { dim = GridView::dimension };
807 static constexpr auto domainI = Dune::index_constant<id>();
818 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
822 if (this->elementIsGhost())
824 const auto globalI = this->assembler().gridGeometry(domainI).elementMapper().index(this->element());
825 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
826 A[globalI][globalI][pvIdx][pvIdx] = 1.0;
829 return LocalResidualValues(0.0);
833 const auto& problem = this->problem();
834 const auto& element = this->element();
835 const auto& fvGeometry = this->fvGeometry();
836 const auto& curElemVolVars = this->curElemVolVars();
837 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
840 const auto globalI = this->assembler().gridGeometry(domainI).elementMapper().index(element);
841 const auto& scv = fvGeometry.scv(globalI);
842 const auto& volVars = curElemVolVars[scv];
845 if (!this->assembler().isStationaryProblem())
846 this->localResidual().addStorageDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv);
849 this->localResidual().addSourceDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv);
852 for (
const auto& scvf : scvfs(fvGeometry))
855 if (!scvf.boundary())
856 this->localResidual().addFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
861 const auto& bcTypes = problem.boundaryTypes(element, scvf);
864 if (bcTypes.hasDirichlet() && !bcTypes.hasNeumann())
865 this->localResidual().addCCDirichletFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
868 else if (bcTypes.hasNeumann() && !bcTypes.hasDirichlet())
869 this->localResidual().addRobinFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
872 DUNE_THROW(Dune::NotImplemented,
"Mixed boundary conditions. Use pure boundary conditions by converting Dirichlet BCs to Robin BCs");
876 if constexpr (Problem::enableInternalDirichletConstraints())
879 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
880 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
882 auto residual = this->evalLocalResidual()[0];
884 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
886 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
888 if (internalDirichletConstraints[eqIdx])
890 residual[eqIdx] = volVars.priVars()[eqIdx] - dirichletValues[eqIdx];
891 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
894 for (
const auto& scvf : scvfs(fvGeometry))
895 if (!scvf.boundary())
896 A[globalI][fvGeometry.scv(scvf.outsideScvIdx()).dofIndex()][eqIdx][pvIdx] = 0.0;
904 return this->evalLocalResidual()[0];
911 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
913 const LocalResidualValues& res, GridVariables& gridVariables)
920 const auto& element = this->element();
921 const auto& fvGeometry = this->fvGeometry();
922 const auto& gridGeometry = fvGeometry.gridGeometry();
923 auto&& curElemVolVars = this->curElemVolVars();
927 const auto globalI = gridGeometry.elementMapper().index(element);
928 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
930 for (
const auto globalJ : stencil)
932 const auto& elementJ = this->assembler().gridGeometry(domainJ).element(globalJ);
933 this->couplingManager().addCouplingDerivatives(A[globalI][globalJ], domainI, element, fvGeometry, curElemVolVars, domainJ, elementJ);
936 if constexpr (Problem::enableInternalDirichletConstraints())
938 const auto& scv = fvGeometry.scv(globalI);
939 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
940 if (internalDirichletConstraints.any())
942 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
943 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
944 if (internalDirichletConstraints[eqIdx])
945 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 helper to deduce a vector with the same size as numbers of equations.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A class for numeric differentiation.
A arithmetic block vector type based on DUNE's reserved vector.
Element solution classes and factory functions.
Helper classes to compute the integration elements.
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
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==DiscretizationMethods::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition: diffmethod.hh:37
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:46
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
Definition: deprecated.hh:149
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 respect to a function parameter.
Definition: numericdifferentiation.hh:61
A arithmetic block vector type based on DUNE's reserved vector.
Definition: reservedblockvector.hh:38
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:60
A base class for all multidomain local assemblers.
Definition: subdomaincclocalassembler.hh:61
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: subdomaincclocalassembler.hh:227
SubDomainCCLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
the constructor
Definition: subdomaincclocalassembler.hh:97
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition: subdomaincclocalassembler.hh:203
LocalResidualValues evalFluxResidual(const Element &neighbor, const SubControlVolumeFace &scvf) const
Evaluates the fluxes depending on the chose time discretization scheme.
Definition: subdomaincclocalassembler.hh:217
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:118
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:141
const Problem & problem() const
return reference to the underlying problem
Definition: subdomaincclocalassembler.hh:259
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: subdomaincclocalassembler.hh:263
static constexpr auto domainId
export the domain id of this sub-domain
Definition: subdomaincclocalassembler.hh:90
GetPropType< TypeTag, Properties::LocalResidual > LocalResidual
the local residual type of this domain
Definition: subdomaincclocalassembler.hh:92
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:182
LocalResidualValues evalLocalStorageResidual() const
Evaluates the storage terms within the element.
Definition: subdomaincclocalassembler.hh:209
void assembleResidual(SubSolutionVector &res)
Assemble the residual only.
Definition: subdomaincclocalassembler.hh:159
The cell-centered scheme multidomain local assembler.
Definition: subdomaincclocalassembler.hh:282
Cell-centered scheme multidomain local assembler using numeric differentiation and implicit time disc...
Definition: subdomaincclocalassembler.hh:294
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:519
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:325
Cell-centered scheme multidomain local assembler using numeric differentiation and explicit time disc...
Definition: subdomaincclocalassembler.hh:647
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:672
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:780
Cell-centered scheme local assembler using analytic differentiation and implicit time discretization.
Definition: subdomaincclocalassembler.hh:795
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:912
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:819
Declares all properties used in Dumux.