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>
102 const QuantityCallback& quantityCallback)
const
103 -> std::remove_reference_t<
decltype(quantityCallback.operator()(0))>
105 using RawReturnType =
decltype(quantityCallback.operator()(0));
106 using ReturnType = std::remove_const_t<std::remove_reference_t<RawReturnType>>;
108 Scalar interiorDistance;
109 Scalar exteriorDistance;
110 computeDistances_(interiorDistance, exteriorDistance, elemCtx, fapIdx);
112 const auto& face = elemCtx.stencil(0).interiorFace(fapIdx);
113 const auto i = face.interiorIndex();
114 const auto j = face.exteriorIndex();
115 const auto focusDofIdx = elemCtx.focusDofIndex();
119 if (i == focusDofIdx) {
120 value = quantityCallback(i) * interiorDistance;
122 value = getValue(quantityCallback(i)) * interiorDistance;
125 if (j == focusDofIdx) {
126 value += quantityCallback(j) * exteriorDistance;
129 value += getValue(quantityCallback(j)) * exteriorDistance;
132 value /= interiorDistance + exteriorDistance;
148 template <
class QuantityCallback>
151 const QuantityCallback& quantityCallback)
const
152 -> std::remove_reference_t<
decltype(quantityCallback.operator()(0))>
154 using RawReturnType =
decltype(quantityCallback.operator()(0));
155 using ReturnType = std::remove_const_t<std::remove_reference_t<RawReturnType>>;
157 Scalar interiorDistance;
158 Scalar exteriorDistance;
159 computeDistances_(interiorDistance, exteriorDistance, elemCtx, fapIdx);
161 const auto& face = elemCtx.stencil(0).interiorFace(fapIdx);
162 const auto i = face.interiorIndex();
163 const auto j = face.exteriorIndex();
164 const auto focusDofIdx = elemCtx.focusDofIndex();
168 if (i == focusDofIdx) {
169 value = quantityCallback(i);
170 for (
int k = 0; k < value.size(); ++k) {
171 value[k] *= interiorDistance;
175 const auto& dofVal = getValue(quantityCallback(i));
176 for (
int k = 0; k < dofVal.size(); ++k) {
177 value[k] = getValue(dofVal[k]) * interiorDistance;
181 if (j == focusDofIdx) {
182 const auto& dofVal = quantityCallback(j);
183 for (
int k = 0; k < dofVal.size(); ++k) {
184 value[k] += dofVal[k] * exteriorDistance;
188 const auto& dofVal = quantityCallback(j);
189 for (
int k = 0; k < dofVal.size(); ++k) {
190 value[k] += getValue(dofVal[k]) * exteriorDistance;
194 const Scalar totDistance = interiorDistance + exteriorDistance;
195 for (
int k = 0; k < value.size(); ++k) {
196 value[k] /= totDistance;
214 template <
class QuantityCallback>
216 const ElementContext& elemCtx,
218 const QuantityCallback& quantityCallback)
const
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();
232 deltay = getValue(quantityCallback(j)) - quantityCallback(i);
234 else if (j == focusIdx) {
235 deltay = quantityCallback(j) - getValue(quantityCallback(i));
238 deltay = getValue(quantityCallback(j)) - getValue(quantityCallback(i));
241 Scalar distSquared = 0.0;
242 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
243 const Scalar tmp = exteriorPos[dimIdx] - interiorPos[dimIdx];
244 distSquared += tmp * tmp;
251 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
252 const Scalar tmp = exteriorPos[dimIdx] - interiorPos[dimIdx];
253 quantityGrad[dimIdx] = deltay *( tmp / distSquared);
271 template <
class QuantityCallback>
274 const QuantityCallback& quantityCallback)
275 ->
decltype(quantityCallback.boundaryValue())
276 {
return quantityCallback.boundaryValue(); }
293 template <
class QuantityCallback>
295 const ElementContext& elemCtx,
297 const QuantityCallback& quantityCallback)
const
299 const auto& stencil = elemCtx.stencil(0);
300 const auto& face = stencil.boundaryFace(faceIdx);
303 if (face.interiorIndex() == elemCtx.focusDofIndex()) {
304 deltay = quantityCallback.boundaryValue() - quantityCallback(face.interiorIndex());
307 deltay = getValue(quantityCallback.boundaryValue()) -
308 getValue(quantityCallback(face.interiorIndex()));
311 const auto& boundaryFacePos = face.integrationPos();
312 const auto& interiorPos = stencil.subControlVolume(face.interiorIndex()).center();
314 Scalar distSquared = 0;
315 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
316 const Scalar tmp = boundaryFacePos[dimIdx] - interiorPos[dimIdx];
317 distSquared += tmp * tmp;
325 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
326 const Scalar tmp = boundaryFacePos[dimIdx] - interiorPos[dimIdx];
327 quantityGrad[dimIdx] = deltay * (tmp / distSquared);
332 void computeDistances_(Scalar& interiorDistance,
333 Scalar& exteriorDistance,
334 const ElementContext& elemCtx,
335 unsigned fapIdx)
const
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();
349 interiorDistance = 0.0;
350 exteriorDistance = 0.0;
351 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
352 interiorDistance += (interiorPos[dimIdx] - integrationPos[dimIdx]) * normal[dimIdx];
353 exteriorDistance += (exteriorPos[dimIdx] - integrationPos[dimIdx]) * normal[dimIdx];
356 interiorDistance = std::sqrt(std::abs(interiorDistance));
357 exteriorDistance = std::sqrt(std::abs(exteriorDistance));
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.
Definition: blackoilboundaryratevector.hh:39
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