fvbasegradientcalculator.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_FV_BASE_GRADIENT_CALCULATOR_HH
29#define EWOMS_FV_BASE_GRADIENT_CALCULATOR_HH
30
31#include "fvbaseproperties.hh"
32
33#include <dune/common/fvector.hh>
34
35namespace Opm {
36template<class TypeTag>
37class EcfvDiscretization;
38
45template<class TypeTag>
47{
53
54 enum { dim = GridView::dimension };
55 enum { dimWorld = GridView::dimensionworld };
56
57 // maximum number of flux approximation points. to calculate this,
58 // we assume that the geometry with the most pointsq is a cube.
59 enum { maxFap = 2 << dim };
60
61 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
62 using EvalDimVector = Dune::FieldVector<Evaluation, dimWorld>;
63
64public:
69 static void registerParameters()
70 { }
71
79 template <bool prepareValues = true, bool prepareGradients = true>
80 void prepare(const ElementContext&, unsigned)
81 { /* noting to do */ }
82
94 template <class QuantityCallback>
95 auto calculateScalarValue(const ElementContext& elemCtx,
96 unsigned fapIdx,
97 const QuantityCallback& quantityCallback) const
98 -> typename std::remove_reference<decltype(quantityCallback.operator()(0))>::type
99 {
100 using RawReturnType = decltype(quantityCallback.operator()(0));
101 using ReturnType = typename std::remove_const<typename std::remove_reference<RawReturnType>::type>::type;
102
103 Scalar interiorDistance;
104 Scalar exteriorDistance;
105 computeDistances_(interiorDistance, exteriorDistance, elemCtx, fapIdx);
106
107 const auto& face = elemCtx.stencil(/*timeIdx=*/0).interiorFace(fapIdx);
108 auto i = face.interiorIndex();
109 auto j = face.exteriorIndex();
110 auto focusDofIdx = elemCtx.focusDofIndex();
111
112 // use the average weighted by distance...
113 ReturnType value;
114 if (i == focusDofIdx)
115 value = quantityCallback(i)*interiorDistance;
116 else
117 value = getValue(quantityCallback(i))*interiorDistance;
118
119 if (j == focusDofIdx)
120 value += quantityCallback(j)*exteriorDistance;
121 else
122 value += getValue(quantityCallback(j))*exteriorDistance;
123
124 value /= interiorDistance + exteriorDistance;
125
126 return value;
127 }
128
140 template <class QuantityCallback>
141 auto calculateVectorValue(const ElementContext& elemCtx,
142 unsigned fapIdx,
143 const QuantityCallback& quantityCallback) const
144 -> typename std::remove_reference<decltype(quantityCallback.operator()(0))>::type
145 {
146 using RawReturnType = decltype(quantityCallback.operator()(0));
147 using ReturnType = typename std::remove_const<typename std::remove_reference<RawReturnType>::type>::type;
148
149 Scalar interiorDistance;
150 Scalar exteriorDistance;
151 computeDistances_(interiorDistance, exteriorDistance, elemCtx, fapIdx);
152
153 const auto& face = elemCtx.stencil(/*timeIdx=*/0).interiorFace(fapIdx);
154 auto i = face.interiorIndex();
155 auto j = face.exteriorIndex();
156 auto focusDofIdx = elemCtx.focusDofIndex();
157
158 // use the average weighted by distance...
159 ReturnType value;
160 if (i == focusDofIdx) {
161 value = quantityCallback(i);
162 for (int k = 0; k < value.size(); ++k)
163 value[k] *= interiorDistance;
164 }
165 else {
166 const auto& dofVal = getValue(quantityCallback(i));
167 for (int k = 0; k < dofVal.size(); ++k)
168 value[k] = getValue(dofVal[k])*interiorDistance;
169 }
170
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;
175 }
176 else {
177 const auto& dofVal = quantityCallback(j);
178 for (int k = 0; k < dofVal.size(); ++k)
179 value[k] += getValue(dofVal[k])*exteriorDistance;
180 }
181
182 Scalar totDistance = interiorDistance + exteriorDistance;
183 for (int k = 0; k < value.size(); ++k)
184 value[k] /= totDistance;
185
186 return value;
187 }
188
200 template <class QuantityCallback>
201 void calculateGradient(EvalDimVector& quantityGrad,
202 const ElementContext& elemCtx,
203 unsigned fapIdx,
204 const QuantityCallback& quantityCallback) const
205 {
206 const auto& stencil = elemCtx.stencil(/*timeIdx=*/0);
207 const auto& face = stencil.interiorFace(fapIdx);
208
209 auto i = face.interiorIndex();
210 auto j = face.exteriorIndex();
211 auto focusIdx = elemCtx.focusDofIndex();
212
213 const auto& interiorPos = stencil.subControlVolume(i).globalPos();
214 const auto& exteriorPos = stencil.subControlVolume(j).globalPos();
215
216 Evaluation deltay;
217 if (i == focusIdx) {
218 deltay =
219 getValue(quantityCallback(j))
220 - quantityCallback(i);
221 }
222 else if (j == focusIdx) {
223 deltay =
224 quantityCallback(j)
225 - getValue(quantityCallback(i));
226 }
227 else
228 deltay =
229 getValue(quantityCallback(j))
230 - getValue(quantityCallback(i));
231
232 Scalar distSquared = 0.0;
233 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
234 Scalar tmp = exteriorPos[dimIdx] - interiorPos[dimIdx];
235 distSquared += tmp*tmp;
236 }
237
238 // divide the gradient by the squared distance between the centers of the
239 // sub-control volumes: the gradient is the normalized directional vector between
240 // the two centers times the ratio of the difference of the values and their
241 // distance, i.e., d/abs(d) * delta y / abs(d) = d*delta y / abs(d)^2.
242 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
243 Scalar tmp = exteriorPos[dimIdx] - interiorPos[dimIdx];
244 quantityGrad[dimIdx] = deltay*(tmp/distSquared);
245 }
246 }
247
262 template <class QuantityCallback>
263 auto calculateBoundaryValue(const ElementContext&,
264 unsigned,
265 const QuantityCallback& quantityCallback)
266 -> decltype(quantityCallback.boundaryValue())
267 { return quantityCallback.boundaryValue(); }
268
283 template <class QuantityCallback>
284 void calculateBoundaryGradient(EvalDimVector& quantityGrad,
285 const ElementContext& elemCtx,
286 unsigned faceIdx,
287 const QuantityCallback& quantityCallback) const
288 {
289 const auto& stencil = elemCtx.stencil(/*timeIdx=*/0);
290 const auto& face = stencil.boundaryFace(faceIdx);
291
292 Evaluation deltay;
293 if (face.interiorIndex() == elemCtx.focusDofIndex())
294 deltay = quantityCallback.boundaryValue() - quantityCallback(face.interiorIndex());
295 else
296 deltay =
297 getValue(quantityCallback.boundaryValue())
298 - getValue(quantityCallback(face.interiorIndex()));
299
300 const auto& boundaryFacePos = face.integrationPos();
301 const auto& interiorPos = stencil.subControlVolume(face.interiorIndex()).center();
302
303 Scalar distSquared = 0;
304 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
305 Scalar tmp = boundaryFacePos[dimIdx] - interiorPos[dimIdx];
306 distSquared += tmp*tmp;
307 }
308
309 // divide the gradient by the squared distance between the center of the
310 // sub-control and the center of the boundary face: the gradient is the
311 // normalized directional vector between the two centers times the ratio of the
312 // difference of the values and their distance, i.e., d/abs(d) * deltay / abs(d)
313 // = d*deltay / abs(d)^2.
314 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
315 Scalar tmp = boundaryFacePos[dimIdx] - interiorPos[dimIdx];
316 quantityGrad[dimIdx] = deltay*(tmp/distSquared);
317 }
318 }
319
320private:
321 void computeDistances_(Scalar& interiorDistance,
322 Scalar& exteriorDistance,
323 const ElementContext& elemCtx,
324 unsigned fapIdx) const
325 {
326 const auto& stencil = elemCtx.stencil(/*timeIdx=*/0);
327 const auto& face = stencil.interiorFace(fapIdx);
328
329 // calculate the distances of the position of the interior and of the exterior
330 // finite volume to the position of the integration point.
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();
337
338 interiorDistance = 0.0;
339 exteriorDistance = 0.0;
340 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
341 interiorDistance +=
342 (interiorPos[dimIdx] - integrationPos[dimIdx])
343 * normal[dimIdx];
344
345 exteriorDistance +=
346 (exteriorPos[dimIdx] - integrationPos[dimIdx])
347 * normal[dimIdx];
348 }
349
350 interiorDistance = std::sqrt(std::abs(interiorDistance));
351 exteriorDistance = std::sqrt(std::abs(exteriorDistance));
352 }
353};
354} // namespace Opm
355
356#endif
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