28#ifndef EWOMS_FV_BASE_GRADIENT_CALCULATOR_HH
29#define EWOMS_FV_BASE_GRADIENT_CALCULATOR_HH
31#include <dune/common/fvector.hh>
41template<
class TypeTag>
42class EcfvDiscretization;
50template<
class TypeTag>
59 enum { dim = GridView::dimension };
60 enum { dimWorld = GridView::dimensionworld };
64 enum { maxFap = 2 << dim };
66 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
67 using EvalDimVector = Dune::FieldVector<Evaluation, dimWorld>;
84 template <
bool prepareValues = true,
bool prepareGradients = true>
85 void prepare(
const ElementContext&,
unsigned)
99 template <
class QuantityCallback>
106 using ReturnType = std::remove_const_t<std::remove_reference_t<RawReturnType>>;
112 const auto&
face = elemCtx.stencil(0).interiorFace(
fapIdx);
113 const auto i =
face.interiorIndex();
114 const auto j =
face.exteriorIndex();
148 template <
class QuantityCallback>
155 using ReturnType = std::remove_const_t<std::remove_reference_t<RawReturnType>>;
161 const auto&
face = elemCtx.stencil(0).interiorFace(
fapIdx);
162 const auto i =
face.interiorIndex();
163 const auto j =
face.exteriorIndex();
170 for (
int k = 0;
k < value.size(); ++
k) {
176 for (
int k = 0;
k <
dofVal.size(); ++
k) {
183 for (
int k = 0;
k <
dofVal.size(); ++
k) {
189 for (
int k = 0;
k <
dofVal.size(); ++
k) {
195 for (
int k = 0;
k < value.size(); ++
k) {
214 template <
class QuantityCallback>
216 const ElementContext& elemCtx,
220 const auto& stencil = elemCtx.stencil(0);
221 const auto&
face = stencil.interiorFace(
fapIdx);
223 const auto i =
face.interiorIndex();
224 const auto j =
face.exteriorIndex();
225 const auto focusIdx = elemCtx.focusDofIndex();
227 const auto&
interiorPos = stencil.subControlVolume(i).globalPos();
228 const auto&
exteriorPos = stencil.subControlVolume(j).globalPos();
271 template <
class QuantityCallback>
293 template <
class QuantityCallback>
295 const ElementContext& elemCtx,
299 const auto& stencil = elemCtx.stencil(0);
300 const auto&
face = stencil.boundaryFace(faceIdx);
303 if (
face.interiorIndex() == elemCtx.focusDofIndex()) {
312 const auto&
interiorPos = stencil.subControlVolume(
face.interiorIndex()).center();
334 const ElementContext& elemCtx,
337 const auto& stencil = elemCtx.stencil(0);
338 const auto&
face = stencil.interiorFace(
fapIdx);
342 const auto& normal =
face.normal();
343 const auto i =
face.interiorIndex();
344 const auto j =
face.exteriorIndex();
345 const auto&
interiorPos = stencil.subControlVolume(i).globalPos();
346 const auto&
exteriorPos = stencil.subControlVolume(j).globalPos();
347 const auto& integrationPos =
face.integrationPos();
Defines a type tags and some fundamental properties all models.
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
void prepare(const ElementContext &, unsigned)
Precomputes the common values to calculate gradients and values of quantities at every interior flux ...
Definition fvbasegradientcalculator.hh:85
static void registerParameters()
Register all run-time parameters for the gradient calculator of the base class of the discretization.
Definition fvbasegradientcalculator.hh:74
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
Declare the properties used by the infrastructure code of the finite volume discretizations.
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