24#ifndef DUMUX_FV_LOCAL_ASSEMBLER_BASE_HH
25#define DUMUX_FV_LOCAL_ASSEMBLER_BASE_HH
27#include <dune/common/reservedvector.hh>
28#include <dune/grid/common/gridenums.hh>
29#include <dune/istl/matrixindexset.hh>
46template<
class TypeTag,
class Assembler,
class Implementation,
bool useImplicitAssembly>
54 using SolutionVector =
typename Assembler::ResidualType;
57 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
58 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
63 using Element =
typename GridView::template Codim<0>::Entity;
75 const SolutionVector&
curSol)
92 const SolutionVector&
curSol,
114 {
return useImplicitAssembly; }
123 if (this->
assembler().isStationaryProblem())
124 DUNE_THROW(Dune::InvalidStateException,
"Using explicit jacobian assembler with stationary local residual");
168 return localResidual_.evalFluxAndSource(element_, fvGeometry_, elemVolVars, elemFluxVarsCache_, elemBcTypes_);
178 return localResidual_.evalStorage(element_, fvGeometry_, prevElemVolVars_, curElemVolVars_);
190 const auto& prevSol = this->
assembler().prevSol();
203 if (!this->
assembler().isStationaryProblem())
217 template<
typename ApplyFunction,
class P = Problem,
typename std::enable_if_t<P::enableInternalDirichletConstra
ints(),
int> = 0>
222 for (
const auto& scvI : scvs(this->
fvGeometry()))
224 const auto internalDirichletConstraints = this->
problem().hasInternalDirichletConstraint(this->
element(), scvI);
225 if (internalDirichletConstraints.any())
227 const auto dirichletValues = this->
problem().internalDirichlet(this->
element(), scvI);
229 for (
int eqIdx = 0; eqIdx < internalDirichletConstraints.size(); ++eqIdx)
230 if (internalDirichletConstraints[eqIdx])
231 applyDirichlet(scvI, dirichletValues, eqIdx, eqIdx);
236 template<
typename ApplyFunction,
class P = Problem,
typename std::enable_if_t<!P::enableInternalDirichletConstra
ints(),
int> = 0>
242 {
return assembler_.problem(); }
246 {
return assembler_; }
254 {
return elementIsGhost_; }
262 {
return fvGeometry_; }
266 {
return curElemVolVars_; }
270 {
return prevElemVolVars_; }
274 {
return elemFluxVarsCache_; }
278 {
return localResidual_; }
282 {
return elemBcTypes_; }
286 {
return fvGeometry_; }
290 {
return curElemVolVars_; }
294 {
return prevElemVolVars_; }
298 {
return elemFluxVarsCache_; }
302 {
return elemBcTypes_; }
306 {
return localResidual_; }
310 {
return *
static_cast<Implementation*
>(
this); }
313 {
return *
static_cast<const Implementation*
>(
this); }
315 template<class T = TypeTag, typename std::enable_if_t<!getPropValue<T, Properties::EnableGridVolumeVariablesCache>(),
int> = 0>
316 VolumeVariables&
getVolVarAccess(GridVolumeVariables& gridVolVars, ElementVolumeVariables& elemVolVars,
const SubControlVolume& scv)
317 {
return elemVolVars[scv]; }
319 template<class T = TypeTag, typename std::enable_if_t<getPropValue<T, Properties::EnableGridVolumeVariablesCache>(),
int> = 0>
320 VolumeVariables&
getVolVarAccess(GridVolumeVariables& gridVolVars, ElementVolumeVariables& elemVolVars,
const SubControlVolume& scv)
321 {
return gridVolVars.volVars(scv); }
325 const Assembler& assembler_;
326 const Element& element_;
327 const SolutionVector& curSol_;
329 FVElementGeometry fvGeometry_;
330 ElementVolumeVariables curElemVolVars_;
331 ElementVolumeVariables prevElemVolVars_;
332 ElementFluxVariablesCache elemFluxVarsCache_;
333 ElementBoundaryTypes elemBcTypes_;
336 bool elementIsGhost_;
An enum class to define various differentiation methods available in order to compute the derivatives...
A arithmetic block vector type based on DUNE's reserved vector.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
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
const ElementVolumeVariables & prevElemVolVars() const
The element volume variables of the provious time step.
Definition: fvlocalassemblerbase.hh:293
void bindLocalViews()
Convenience function bind and prepare all relevant variables required for the evaluation of the local...
Definition: fvlocalassemblerbase.hh:185
ElementResidualVector evalLocalResidual(const ElementVolumeVariables &elemVolVars) const
Evaluates the complete local residual for the current element.
Definition: fvlocalassemblerbase.hh:137
ElementVolumeVariables & curElemVolVars()
The current element volume variables.
Definition: fvlocalassemblerbase.hh:265
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: fvlocalassemblerbase.hh:68
ElementBoundaryTypes & elemBcTypes()
The element's boundary types.
Definition: fvlocalassemblerbase.hh:281
Implementation & asImp_()
Definition: fvlocalassemblerbase.hh:309
ElementVolumeVariables & prevElemVolVars()
The element volume variables of the provious time step.
Definition: fvlocalassemblerbase.hh:269
ElementResidualVector evalLocalResidual() const
Convenience function to evaluate the complete local residual for the current element....
Definition: fvlocalassemblerbase.hh:120
const Problem & problem() const
The problem.
Definition: fvlocalassemblerbase.hh:241
const ElementBoundaryTypes & elemBcTypes() const
The element's boundary types.
Definition: fvlocalassemblerbase.hh:301
const ElementFluxVariablesCache & elemFluxVarsCache() const
The element flux variables cache.
Definition: fvlocalassemblerbase.hh:297
const Implementation & asImp_() const
Definition: fvlocalassemblerbase.hh:312
FVLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementVolumeVariables &prevElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const LocalResidual &localResidual, const bool elementIsGhost)
The constructor. General version explicitly expecting each argument.
Definition: fvlocalassemblerbase.hh:90
const ElementVolumeVariables & curElemVolVars() const
The current element volume variables.
Definition: fvlocalassemblerbase.hh:289
FVLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol)
The constructor. Delegates to the general constructor.
Definition: fvlocalassemblerbase.hh:73
const LocalResidual & localResidual() const
The local residual for the current element.
Definition: fvlocalassemblerbase.hh:305
FVElementGeometry & fvGeometry()
The global finite volume geometry.
Definition: fvlocalassemblerbase.hh:261
ElementResidualVector evalLocalFluxAndSourceResidual(const ElementVolumeVariables &elemVolVars) const
Evaluates the flux and source terms (i.e, the terms without a time derivative) of the local residual ...
Definition: fvlocalassemblerbase.hh:166
const Assembler & assembler() const
The assembler.
Definition: fvlocalassemblerbase.hh:245
ElementFluxVariablesCache & elemFluxVarsCache()
The element flux variables cache.
Definition: fvlocalassemblerbase.hh:273
bool elementIsGhost() const
Returns if element is a ghost entity.
Definition: fvlocalassemblerbase.hh:253
LocalResidual & localResidual()
The local residual for the current element.
Definition: fvlocalassemblerbase.hh:277
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: fvlocalassemblerbase.hh:113
ElementResidualVector evalLocalFluxAndSourceResidual() const
Convenience function to evaluate the flux and source terms (i.e, the terms without a time derivative)...
Definition: fvlocalassemblerbase.hh:154
void enforceInternalDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforces Dirichlet constraints if enabled in the problem.
Definition: fvlocalassemblerbase.hh:218
const Element & element() const
The current element.
Definition: fvlocalassemblerbase.hh:249
VolumeVariables & getVolVarAccess(GridVolumeVariables &gridVolVars, ElementVolumeVariables &elemVolVars, const SubControlVolume &scv)
Definition: fvlocalassemblerbase.hh:316
const FVElementGeometry & fvGeometry() const
The finite volume geometry.
Definition: fvlocalassemblerbase.hh:285
GetPropType< TypeTag, Properties::LocalResidual > LocalResidual
Definition: fvlocalassemblerbase.hh:67
ElementResidualVector evalLocalStorageResidual() const
Convenience function to evaluate storage term (i.e, the term with a time derivative) of the local res...
Definition: fvlocalassemblerbase.hh:176
const SolutionVector & curSol() const
The current solution.
Definition: fvlocalassemblerbase.hh:257
Declares all properties used in Dumux.