28#ifndef EWOMS_FV_BASE_GRADIENT_CALCULATOR_HH
29#define EWOMS_FV_BASE_GRADIENT_CALCULATOR_HH
33#include <dune/common/fvector.hh>
36template<
class TypeTag>
37class EcfvDiscretization;
45template<
class TypeTag>
54 enum { dim = GridView::dimension };
55 enum { dimWorld = GridView::dimensionworld };
59 enum { maxFap = 2 << dim };
61 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
62 using EvalDimVector = Dune::FieldVector<Evaluation, dimWorld>;
79 template <
bool prepareValues = true,
bool prepareGradients = true>
80 void prepare(
const ElementContext&,
unsigned)
94 template <
class QuantityCallback>
97 const QuantityCallback& quantityCallback)
const
98 ->
typename std::remove_reference<
decltype(quantityCallback.operator()(0))>::type
100 using RawReturnType =
decltype(quantityCallback.operator()(0));
101 using ReturnType =
typename std::remove_const<typename std::remove_reference<RawReturnType>::type>::type;
103 Scalar interiorDistance;
104 Scalar exteriorDistance;
105 computeDistances_(interiorDistance, exteriorDistance, elemCtx, fapIdx);
107 const auto& face = elemCtx.stencil(0).interiorFace(fapIdx);
108 auto i = face.interiorIndex();
109 auto j = face.exteriorIndex();
110 auto focusDofIdx = elemCtx.focusDofIndex();
114 if (i == focusDofIdx)
115 value = quantityCallback(i)*interiorDistance;
117 value = getValue(quantityCallback(i))*interiorDistance;
119 if (j == focusDofIdx)
120 value += quantityCallback(j)*exteriorDistance;
122 value += getValue(quantityCallback(j))*exteriorDistance;
124 value /= interiorDistance + exteriorDistance;
140 template <
class QuantityCallback>
143 const QuantityCallback& quantityCallback)
const
144 ->
typename std::remove_reference<
decltype(quantityCallback.operator()(0))>::type
146 using RawReturnType =
decltype(quantityCallback.operator()(0));
147 using ReturnType =
typename std::remove_const<typename std::remove_reference<RawReturnType>::type>::type;
149 Scalar interiorDistance;
150 Scalar exteriorDistance;
151 computeDistances_(interiorDistance, exteriorDistance, elemCtx, fapIdx);
153 const auto& face = elemCtx.stencil(0).interiorFace(fapIdx);
154 auto i = face.interiorIndex();
155 auto j = face.exteriorIndex();
156 auto focusDofIdx = elemCtx.focusDofIndex();
160 if (i == focusDofIdx) {
161 value = quantityCallback(i);
162 for (
int k = 0; k < value.size(); ++k)
163 value[k] *= interiorDistance;
166 const auto& dofVal = getValue(quantityCallback(i));
167 for (
int k = 0; k < dofVal.size(); ++k)
168 value[k] = getValue(dofVal[k])*interiorDistance;
171 if (j == focusDofIdx) {
172 const auto& dofVal = quantityCallback(j);
173 for (
int k = 0; k < dofVal.size(); ++k)
174 value[k] += dofVal[k]*exteriorDistance;
177 const auto& dofVal = quantityCallback(j);
178 for (
int k = 0; k < dofVal.size(); ++k)
179 value[k] += getValue(dofVal[k])*exteriorDistance;
182 Scalar totDistance = interiorDistance + exteriorDistance;
183 for (
int k = 0; k < value.size(); ++k)
184 value[k] /= totDistance;
200 template <
class QuantityCallback>
202 const ElementContext& elemCtx,
204 const QuantityCallback& quantityCallback)
const
206 const auto& stencil = elemCtx.stencil(0);
207 const auto& face = stencil.interiorFace(fapIdx);
209 auto i = face.interiorIndex();
210 auto j = face.exteriorIndex();
211 auto focusIdx = elemCtx.focusDofIndex();
213 const auto& interiorPos = stencil.subControlVolume(i).globalPos();
214 const auto& exteriorPos = stencil.subControlVolume(j).globalPos();
219 getValue(quantityCallback(j))
220 - quantityCallback(i);
222 else if (j == focusIdx) {
225 - getValue(quantityCallback(i));
229 getValue(quantityCallback(j))
230 - getValue(quantityCallback(i));
232 Scalar distSquared = 0.0;
233 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
234 Scalar tmp = exteriorPos[dimIdx] - interiorPos[dimIdx];
235 distSquared += tmp*tmp;
242 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
243 Scalar tmp = exteriorPos[dimIdx] - interiorPos[dimIdx];
244 quantityGrad[dimIdx] = deltay*(tmp/distSquared);
262 template <
class QuantityCallback>
265 const QuantityCallback& quantityCallback)
266 ->
decltype(quantityCallback.boundaryValue())
267 {
return quantityCallback.boundaryValue(); }
283 template <
class QuantityCallback>
285 const ElementContext& elemCtx,
287 const QuantityCallback& quantityCallback)
const
289 const auto& stencil = elemCtx.stencil(0);
290 const auto& face = stencil.boundaryFace(faceIdx);
293 if (face.interiorIndex() == elemCtx.focusDofIndex())
294 deltay = quantityCallback.boundaryValue() - quantityCallback(face.interiorIndex());
297 getValue(quantityCallback.boundaryValue())
298 - getValue(quantityCallback(face.interiorIndex()));
300 const auto& boundaryFacePos = face.integrationPos();
301 const auto& interiorPos = stencil.subControlVolume(face.interiorIndex()).center();
303 Scalar distSquared = 0;
304 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
305 Scalar tmp = boundaryFacePos[dimIdx] - interiorPos[dimIdx];
306 distSquared += tmp*tmp;
314 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
315 Scalar tmp = boundaryFacePos[dimIdx] - interiorPos[dimIdx];
316 quantityGrad[dimIdx] = deltay*(tmp/distSquared);
321 void computeDistances_(Scalar& interiorDistance,
322 Scalar& exteriorDistance,
323 const ElementContext& elemCtx,
324 unsigned fapIdx)
const
326 const auto& stencil = elemCtx.stencil(0);
327 const auto& face = stencil.interiorFace(fapIdx);
331 const auto& normal = face.normal();
332 auto i = face.interiorIndex();
333 auto j = face.exteriorIndex();
334 const auto& interiorPos = stencil.subControlVolume(i).globalPos();
335 const auto& exteriorPos = stencil.subControlVolume(j).globalPos();
336 const auto& integrationPos = face.integrationPos();
338 interiorDistance = 0.0;
339 exteriorDistance = 0.0;
340 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
342 (interiorPos[dimIdx] - integrationPos[dimIdx])
346 (exteriorPos[dimIdx] - integrationPos[dimIdx])
350 interiorDistance = std::sqrt(std::abs(interiorDistance));
351 exteriorDistance = std::sqrt(std::abs(exteriorDistance));
This class calculates gradients of arbitrary quantities at flux integration points using the two-poin...
Definition: fvbasegradientcalculator.hh:47
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:284
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:201
void prepare(const ElementContext &, unsigned)
Precomputes the common values to calculate gradients and values of quantities at every interior flux ...
Definition: fvbasegradientcalculator.hh:80
auto calculateScalarValue(const ElementContext &elemCtx, unsigned fapIdx, const QuantityCallback &quantityCallback) const -> typename std::remove_reference< decltype(quantityCallback.operator()(0))>::type
Calculates the value of an arbitrary scalar quantity at any interior flux approximation point.
Definition: fvbasegradientcalculator.hh:95
static void registerParameters()
Register all run-time parameters for the gradient calculator of the base class of the discretization.
Definition: fvbasegradientcalculator.hh:69
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:263
auto calculateVectorValue(const ElementContext &elemCtx, unsigned fapIdx, const QuantityCallback &quantityCallback) const -> typename std::remove_reference< decltype(quantityCallback.operator()(0))>::type
Calculates the value of an arbitrary vectorial quantity at any interior flux approximation point.
Definition: fvbasegradientcalculator.hh:141
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilboundaryratevector.hh:37
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:235