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>
57template<std::
size_t id,
class TypeTag,
class Assembler,
class Implementation,
bool implicit>
65 using SolutionVector =
typename Assembler::SolutionVector;
70 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
71 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
72 using ElementFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
73 using Scalar =
typename GridVariables::Scalar;
75 using GridGeometry =
typename GridVariables::GridGeometry;
76 using FVElementGeometry =
typename GridGeometry::LocalView;
77 using SubControlVolume =
typename GridGeometry::SubControlVolume;
78 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
79 using GridView =
typename GridGeometry::GridView;
80 using Element =
typename GridView::template Codim<0>::Entity;
82 using CouplingManager =
typename Assembler::CouplingManager;
85 static_assert(!Problem::enableInternalDirichletConstraints(),
"Internal Dirichlet constraints are currently not implemented for cc-methods!");
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 *= -scv.volume()*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>
288 using FVElementGeometry =
typename GridGeometry::LocalView;
289 using GridView =
typename GridGeometry::GridView;
290 using Element =
typename GridView::template Codim<0>::Entity;
293 enum { dim = GridView::dimension };
297 static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize;
298 static constexpr auto domainI = Dune::index_constant<id>();
301 using ParentType::ParentType;
309 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
321 const auto& gridGeometry =
fvGeometry.gridGeometry();
326 const auto globalI = gridGeometry.elementMapper().index(
element);
327 const auto& connectivityMap = gridGeometry.connectivityMap();
328 const auto numNeighbors = connectivityMap[globalI].size();
331 Dune::ReservedVector<Element, maxElementStencilSize> neighborElements;
332 neighborElements.resize(numNeighbors);
336 Residuals origResiduals(numNeighbors + 1); origResiduals = 0.0;
342 for (
const auto& dataJ : connectivityMap[globalI])
344 neighborElements[j-1] = gridGeometry.element(dataJ.globalJ);
345 for (
const auto scvfIdx : dataJ.scvfsJ)
358 const auto origPriVars =
curSol[globalI];
359 const auto origVolVars = curVolVars;
366 Residuals partialDerivs(numNeighbors + 1);
368 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
377 auto evalResiduals = [&](Scalar priVar)
379 Residuals partialDerivsTmp(numNeighbors + 1);
380 partialDerivsTmp = 0.0;
382 elemSol[0][pvIdx] = priVar;
383 this->
couplingManager().updateCouplingContext(domainI, *
this, domainI, globalI, elemSol[0], pvIdx);
386 if (enableGridFluxVarsCache)
393 for (std::size_t k = 0; k < connectivityMap[globalI].size(); ++k)
394 for (
auto scvfIdx : connectivityMap[globalI][k].scvfsJ)
397 return partialDerivsTmp;
404 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
407 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
410 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
414 for (
const auto& dataJ : connectivityMap[globalI])
415 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j++][eqIdx];
419 curVolVars = origVolVars;
422 elemSol[0][pvIdx] = origPriVars[pvIdx];
425 this->
couplingManager().updateCouplingContext(domainI, *
this, domainI, globalI, elemSol[0], pvIdx);
434 if (enableGridFluxVarsCache)
438 typename ParentType::LocalResidual::ElementResidualVector origElementResidual({origResiduals[0]});
439 this->
couplingManager().evalAdditionalDomainDerivatives(domainI, *
this, origElementResidual, A, gridVariables);
442 return origResiduals[0];
449 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
451 const LocalResidualValues& res, GridVariables& gridVariables)
460 const auto& gridGeometry =
fvGeometry.gridGeometry();
465 const auto globalI = gridGeometry.elementMapper().index(
element);
467 const auto& curSolJ = this->
curSol()[domainJ];
473 auto updateCoupledVariables = [&] ()
477 if (enableGridFluxVarsCache)
479 if (enableGridVolVarsCache)
481 this->
couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
492 if (enableGridVolVarsCache)
499 for (
const auto& globalJ : stencil)
502 const auto origPriVarsJ = curSolJ[globalJ];
503 auto priVarsJ = origPriVarsJ;
505 const auto origResidual = this->
couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ)[0];
507 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
509 auto evalCouplingResidual = [&](Scalar priVar)
512 priVarsJ[pvIdx] = priVar;
513 this->
couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
514 updateCoupledVariables();
515 return this->
couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ)[0];
519 LocalResidualValues partialDeriv(0.0);
520 const auto& paramGroup = this->
assembler().problem(domainJ).paramGroup();
522 static const auto epsCoupl = this->
couplingManager().numericEpsilon(domainJ, paramGroup);
524 epsCoupl(priVarsJ[pvIdx], pvIdx), numDiffMethod);
527 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
528 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
531 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
534 this->
couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
540 updateCoupledVariables();
551template<std::
size_t id,
class TypeTag,
class Assembler>
554 SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>, false>
563 static constexpr auto domainI = Dune::index_constant<id>();
566 using ParentType::ParentType;
578 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
581 if (this->
assembler().isStationaryProblem())
582 DUNE_THROW(Dune::InvalidStateException,
"Using explicit jacobian assembler with stationary local residual");
597 const auto& gridGeometry =
fvGeometry.gridGeometry();
601 const auto globalI = gridGeometry.elementMapper().index(
element);
608 const auto origPriVars =
curSol[globalI];
609 const auto origVolVars = curVolVars;
615 LocalResidualValues partialDeriv;
616 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
621 auto evalStorage = [&](Scalar priVar)
625 elemSol[0][pvIdx] = priVar;
636 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
642 else partialDeriv[pvIdx] = 1.0;
646 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
647 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
650 curVolVars = origVolVars;
653 elemSol[0][pvIdx] = origPriVars[pvIdx];
666 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
668 const LocalResidualValues& res, GridVariables& gridVariables)
677template<std::
size_t id,
class TypeTag,
class Assembler>
680 SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::analytic, true>, true>
687 using Element =
typename GridView::template Codim<0>::Entity;
690 enum { dim = GridView::dimension };
692 static constexpr auto domainI = Dune::index_constant<id>();
695 using ParentType::ParentType;
703 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
714 const auto globalI = this->
assembler().gridGeometry(domainI).elementMapper().index(
element);
719 if (!this->
assembler().isStationaryProblem())
729 if (!scvf.boundary())
738 if (bcTypes.hasDirichlet() && !bcTypes.hasNeumann())
742 else if (bcTypes.hasNeumann() && !bcTypes.hasDirichlet())
746 DUNE_THROW(Dune::NotImplemented,
"Mixed boundary conditions. Use pure boundary conditions by converting Dirichlet BCs to Robin BCs");
760 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
762 const LocalResidualValues& res, GridVariables& gridVariables)
771 const auto& gridGeometry =
fvGeometry.gridGeometry();
776 const auto globalI = gridGeometry.elementMapper().index(
element);
779 for (
const auto globalJ : stencil)
781 const auto& elementJ = this->
assembler().gridGeometry(domainJ).element(globalJ);
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 base class for all local assemblers.
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.
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
@ analytic
Definition diffmethod.hh:38
@ 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:375
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
Definition common/pdesolver.hh:35
ElementVolumeVariables & curElemVolVars()
Definition fvlocalassemblerbase.hh:263
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition fvlocalassemblerbase.hh:68
Implementation & asImp_()
Definition fvlocalassemblerbase.hh:307
ElementVolumeVariables & prevElemVolVars()
Definition fvlocalassemblerbase.hh:267
ElementResidualVector evalLocalResidual() const
Definition fvlocalassemblerbase.hh:120
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
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
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
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()
Definition subdomaincclocalassembler.hh:249
static constexpr auto domainId
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
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:450
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:310
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:579
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:667
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:761
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:704
Declares all properties used in Dumux.