28#ifndef EWOMS_P1FE_GRADIENT_CALCULATOR_HH
29#define EWOMS_P1FE_GRADIENT_CALCULATOR_HH
31#include <dune/common/fvector.hh>
33#include <dune/geometry/type.hh>
35#include <dune/common/version.hh>
39#if HAVE_DUNE_LOCALFUNCTIONS
40#include <dune/localfunctions/lagrange/pqkfactory.hh>
60template<
class TypeTag>
69 enum { dim = GridView::dimension };
74 enum { maxDof = (2 << dim) };
75 enum { maxFap = maxDof };
77 using CoordScalar =
typename GridView::ctype;
78 using DimVector = Dune::FieldVector<Scalar, dim>;
80#if HAVE_DUNE_LOCALFUNCTIONS
83 using LocalBasisTraits =
typename LocalFiniteElement::Traits::LocalBasisType::Traits;
84 using ShapeJacobian =
typename LocalBasisTraits::JacobianType;
95 template <
bool prepareValues = true,
bool prepareGradients = true>
100#if !HAVE_DUNE_LOCALFUNCTIONS
102 throw std::logic_error(
"The dune-localfunctions module is required in order to use"
103 " finite element gradients");
105 const auto& stencil = elemCtx.stencil(
timeIdx);
111 for (
unsigned faceIdx = 0; faceIdx < stencil.numInteriorFaces(); ++faceIdx) {
112 const auto&
localFacePos = stencil.interiorFace(faceIdx).localPos();
128 const auto&
geom = elemCtx.element().geometry();
132 std::size_t numVertices = elemCtx.numDof(
timeIdx);
157 template <
class QuantityCallback>
161 -> std::remove_reference_t<typename QuantityCallback::ResultType>
164#if !HAVE_DUNE_LOCALFUNCTIONS
166 throw std::logic_error(
"The dune-localfunctions module is required in order to use"
167 " finite element gradients");
169 using QuantityConstType = std::remove_reference_t<typename QuantityCallback::ResultType>;
170 using QuantityType = std::remove_const_t<QuantityConstType>;
177 if (std::is_same_v<QuantityType, Scalar> ||
178 elemCtx.focusDofIndex() ==
vertIdx)
206 template <
class QuantityCallback>
210 -> std::remove_reference_t<typename QuantityCallback::ResultType>
213#if !HAVE_DUNE_LOCALFUNCTIONS
215 throw std::logic_error(
"The dune-localfunctions module is required in order to use"
216 " finite element gradients");
218 using QuantityConstType = std::remove_reference_t<typename QuantityCallback::ResultType>;
219 using QuantityType = std::remove_const_t<QuantityConstType>;
221 using RawFieldType =
decltype(std::declval<QuantityType>()[0]);
222 using FieldType = std::remove_const_t<std::remove_reference_t<RawFieldType>>;
230 if (std::is_same_v<QuantityType, Scalar> ||
231 elemCtx.focusDofIndex() ==
vertIdx)
234 for (
unsigned k = 0;
k < tmp.size(); ++
k) {
240 for (
unsigned k = 0;
k < tmp.size(); ++
k) {
266 template <
class QuantityCallback,
class EvalDimVector>
273#if !HAVE_DUNE_LOCALFUNCTIONS
275 throw std::logic_error(
"The dune-localfunctions module is required in order to use"
276 " finite element gradients");
278 using QuantityConstType = std::remove_reference_t<typename QuantityCallback::ResultType>;
279 using QuantityType = std::remove_const_t<QuantityConstType>;
285 if (std::is_same_v<QuantityType, Scalar> ||
286 elemCtx.focusDofIndex() ==
vertIdx)
323 template <
class QuantityCallback>
345 template <
class QuantityCallback,
class EvalDimVector>
347 const ElementContext& elemCtx,
352#if HAVE_DUNE_LOCALFUNCTIONS
358#if HAVE_DUNE_LOCALFUNCTIONS
362 std::array<std::vector<Dune::FieldVector<Scalar, 1>>, maxFap>
p1Value_{};
363 std::array<std::array<DimVector, maxDof>, maxFap>{};
367#if HAVE_DUNE_LOCALFUNCTIONS
368template<
class TypeTag>
369typename P1FeGradientCalculator<TypeTag>::LocalFiniteElementCache
370P1FeGradientCalculator<TypeTag>::feCache_;
This class calculates gradients of arbitrary quantities at flux integration points using the two-poin...
Definition fvbasegradientcalculator.hh:52
auto calculateScalarValue(const ElementContext &elemCtx, unsigned fapIdx, const QuantityCallback &quantityCallback) const -> std::remove_reference_t< decltype(quantityCallback.operator()(0))>
Calculates the value of an arbitrary scalar quantity at any interior flux approximation point.
Definition fvbasegradientcalculator.hh:100
void calculateBoundaryGradient(EvalDimVector &quantityGrad, const ElementContext &elemCtx, unsigned faceIdx, const QuantityCallback &quantityCallback) const
Calculates the gradient of an arbitrary quantity at any flux approximation point on the boundary.
Definition fvbasegradientcalculator.hh:294
auto calculateVectorValue(const ElementContext &elemCtx, unsigned fapIdx, const QuantityCallback &quantityCallback) const -> std::remove_reference_t< decltype(quantityCallback.operator()(0))>
Calculates the value of an arbitrary vectorial quantity at any interior flux approximation point.
Definition fvbasegradientcalculator.hh:149
void calculateGradient(EvalDimVector &quantityGrad, const ElementContext &elemCtx, unsigned fapIdx, const QuantityCallback &quantityCallback) const
Calculates the gradient of an arbitrary quantity at any flux approximation point.
Definition fvbasegradientcalculator.hh:215
auto calculateBoundaryValue(const ElementContext &, unsigned, const QuantityCallback &quantityCallback) -> decltype(quantityCallback.boundaryValue())
Calculates the value of an arbitrary quantity at any flux approximation point on the grid boundary.
Definition fvbasegradientcalculator.hh:272
This class calculates gradients of arbitrary quantities at flux integration points using first order ...
Definition p1fegradientcalculator.hh:62
auto calculateVectorValue(const ElementContext &elemCtx, unsigned fapIdx, const QuantityCallback &quantityCallback) const -> std::remove_reference_t< typename QuantityCallback::ResultType >
Calculates the value of an arbitrary quantity at any interior flux approximation point.
Definition p1fegradientcalculator.hh:207
void calculateBoundaryGradient(EvalDimVector &quantityGrad, const ElementContext &elemCtx, unsigned fapIdx, const QuantityCallback &quantityCallback) const
Calculates the gradient of an arbitrary quantity at any flux approximation point on the boundary.
Definition p1fegradientcalculator.hh:346
auto calculateScalarValue(const ElementContext &elemCtx, unsigned fapIdx, const QuantityCallback &quantityCallback) const -> std::remove_reference_t< typename QuantityCallback::ResultType >
Calculates the value of an arbitrary quantity at any interior flux approximation point.
Definition p1fegradientcalculator.hh:158
auto calculateBoundaryValue(const ElementContext &elemCtx, unsigned fapIdx, const QuantityCallback &quantityCallback) -> decltype(ParentType::calculateBoundaryValue(elemCtx, fapIdx, quantityCallback))
Calculates the value of an arbitrary quantity at any flux approximation point on the grid boundary.
Definition p1fegradientcalculator.hh:324
void prepare(const ElementContext &elemCtx, unsigned timeIdx)
Precomputes the common values to calculate gradients and values of quantities at any flux approximati...
Definition p1fegradientcalculator.hh:96
void calculateGradient(EvalDimVector &quantityGrad, const ElementContext &elemCtx, unsigned fapIdx, const QuantityCallback &quantityCallback) const
Calculates the gradient of an arbitrary quantity at any flux approximation point.
Definition p1fegradientcalculator.hh:267
This class calculates gradients of arbitrary quantities at flux integration points using the two-poin...
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
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:233
Declares the basic properties used by the common infrastructure of the vertex-centered finite volume ...