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 <dune/common/fvector.hh>
32
35
36#include <cmath>
37#include <type_traits>
38
39namespace Opm {
40
41template<class TypeTag>
42class EcfvDiscretization;
43
50template<class TypeTag>
52{
58
59 enum { dim = GridView::dimension };
60 enum { dimWorld = GridView::dimensionworld };
61
62 // maximum number of flux approximation points. to calculate this,
63 // we assume that the geometry with the most points is a cube.
64 enum { maxFap = 2 << dim };
65
66 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
67 using EvalDimVector = Dune::FieldVector<Evaluation, dimWorld>;
68
69public:
74 static void registerParameters()
75 {}
76
84 template <bool prepareValues = true, bool prepareGradients = true>
85 void prepare(const ElementContext&, unsigned)
86 { /* nothing to do */ }
87
99 template <class QuantityCallback>
100 auto calculateScalarValue(const ElementContext& elemCtx,
101 unsigned fapIdx,
102 const QuantityCallback& quantityCallback) const
103 -> std::remove_reference_t<decltype(quantityCallback.operator()(0))>
104 {
105 using RawReturnType = decltype(quantityCallback.operator()(0));
106 using ReturnType = std::remove_const_t<std::remove_reference_t<RawReturnType>>;
107
108 Scalar interiorDistance;
109 Scalar exteriorDistance;
110 computeDistances_(interiorDistance, exteriorDistance, elemCtx, fapIdx);
111
112 const auto& face = elemCtx.stencil(/*timeIdx=*/0).interiorFace(fapIdx);
113 const auto i = face.interiorIndex();
114 const auto j = face.exteriorIndex();
115 const auto focusDofIdx = elemCtx.focusDofIndex();
116
117 // use the average weighted by distance...
118 ReturnType value;
119 if (i == focusDofIdx) {
120 value = quantityCallback(i) * interiorDistance;
121 } else {
122 value = getValue(quantityCallback(i)) * interiorDistance;
123 }
124
125 if (j == focusDofIdx) {
126 value += quantityCallback(j) * exteriorDistance;
127 }
128 else {
129 value += getValue(quantityCallback(j)) * exteriorDistance;
130 }
131
132 value /= interiorDistance + exteriorDistance;
133
134 return value;
135 }
136
148 template <class QuantityCallback>
149 auto calculateVectorValue(const ElementContext& elemCtx,
150 unsigned fapIdx,
151 const QuantityCallback& quantityCallback) const
152 -> std::remove_reference_t<decltype(quantityCallback.operator()(0))>
153 {
154 using RawReturnType = decltype(quantityCallback.operator()(0));
155 using ReturnType = std::remove_const_t<std::remove_reference_t<RawReturnType>>;
156
157 Scalar interiorDistance;
158 Scalar exteriorDistance;
159 computeDistances_(interiorDistance, exteriorDistance, elemCtx, fapIdx);
160
161 const auto& face = elemCtx.stencil(/*timeIdx=*/0).interiorFace(fapIdx);
162 const auto i = face.interiorIndex();
163 const auto j = face.exteriorIndex();
164 const auto focusDofIdx = elemCtx.focusDofIndex();
165
166 // use the average weighted by distance...
167 ReturnType value;
168 if (i == focusDofIdx) {
169 value = quantityCallback(i);
170 for (int k = 0; k < value.size(); ++k) {
171 value[k] *= interiorDistance;
172 }
173 }
174 else {
175 const auto& dofVal = getValue(quantityCallback(i));
176 for (int k = 0; k < dofVal.size(); ++k) {
177 value[k] = getValue(dofVal[k]) * interiorDistance;
178 }
179 }
180
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;
185 }
186 }
187 else {
188 const auto& dofVal = quantityCallback(j);
189 for (int k = 0; k < dofVal.size(); ++k) {
190 value[k] += getValue(dofVal[k]) * exteriorDistance;
191 }
192 }
193
194 const Scalar totDistance = interiorDistance + exteriorDistance;
195 for (int k = 0; k < value.size(); ++k) {
196 value[k] /= totDistance;
197 }
198
199 return value;
200 }
201
213 template <class QuantityCallback>
214 void calculateGradient(EvalDimVector& quantityGrad,
215 const ElementContext& elemCtx,
216 unsigned fapIdx,
217 const QuantityCallback& quantityCallback) const
218 {
219 const auto& stencil = elemCtx.stencil(/*timeIdx=*/0);
220 const auto& face = stencil.interiorFace(fapIdx);
221
222 const auto i = face.interiorIndex();
223 const auto j = face.exteriorIndex();
224 const auto focusIdx = elemCtx.focusDofIndex();
225
226 const auto& interiorPos = stencil.subControlVolume(i).globalPos();
227 const auto& exteriorPos = stencil.subControlVolume(j).globalPos();
228
229 Evaluation deltay;
230 if (i == focusIdx) {
231 deltay = getValue(quantityCallback(j)) - quantityCallback(i);
232 }
233 else if (j == focusIdx) {
234 deltay = quantityCallback(j) - getValue(quantityCallback(i));
235 }
236 else {
237 deltay = getValue(quantityCallback(j)) - getValue(quantityCallback(i));
238 }
239
240 Scalar distSquared = 0.0;
241 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
242 const Scalar tmp = exteriorPos[dimIdx] - interiorPos[dimIdx];
243 distSquared += tmp * tmp;
244 }
245
246 // divide the gradient by the squared distance between the centers of the
247 // sub-control volumes: the gradient is the normalized directional vector between
248 // the two centers times the ratio of the difference of the values and their
249 // distance, i.e., d/abs(d) * delta y / abs(d) = d*delta y / abs(d)^2.
250 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
251 const Scalar tmp = exteriorPos[dimIdx] - interiorPos[dimIdx];
252 quantityGrad[dimIdx] = deltay *( tmp / distSquared);
253 }
254 }
255
270 template <class QuantityCallback>
271 auto calculateBoundaryValue(const ElementContext&,
272 unsigned,
273 const QuantityCallback& quantityCallback)
274 -> decltype(quantityCallback.boundaryValue())
275 { return quantityCallback.boundaryValue(); }
276
291 template <class QuantityCallback>
292 void calculateBoundaryGradient(EvalDimVector& quantityGrad,
293 const ElementContext& elemCtx,
294 unsigned faceIdx,
295 const QuantityCallback& quantityCallback) const
296 {
297 const auto& stencil = elemCtx.stencil(/*timeIdx=*/0);
298 const auto& face = stencil.boundaryFace(faceIdx);
299
300 Evaluation deltay;
301 if (face.interiorIndex() == elemCtx.focusDofIndex()) {
302 deltay = quantityCallback.boundaryValue() - quantityCallback(face.interiorIndex());
303 }
304 else {
305 deltay = getValue(quantityCallback.boundaryValue()) -
306 getValue(quantityCallback(face.interiorIndex()));
307 }
308
309 const auto& boundaryFacePos = face.integrationPos();
310 const auto& interiorPos = stencil.subControlVolume(face.interiorIndex()).center();
311
312 Scalar distSquared = 0;
313 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
314 const Scalar tmp = boundaryFacePos[dimIdx] - interiorPos[dimIdx];
315 distSquared += tmp * tmp;
316 }
317
318 // divide the gradient by the squared distance between the center of the
319 // sub-control and the center of the boundary face: the gradient is the
320 // normalized directional vector between the two centers times the ratio of the
321 // difference of the values and their distance, i.e., d/abs(d) * deltay / abs(d)
322 // = d*deltay / abs(d)^2.
323 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
324 const Scalar tmp = boundaryFacePos[dimIdx] - interiorPos[dimIdx];
325 quantityGrad[dimIdx] = deltay * (tmp / distSquared);
326 }
327 }
328
329private:
330 void computeDistances_(Scalar& interiorDistance,
331 Scalar& exteriorDistance,
332 const ElementContext& elemCtx,
333 unsigned fapIdx) const
334 {
335 const auto& stencil = elemCtx.stencil(/*timeIdx=*/0);
336 const auto& face = stencil.interiorFace(fapIdx);
337
338 // calculate the distances of the position of the interior and of the exterior
339 // finite volume to the position of the integration point.
340 const auto& normal = face.normal();
341 const auto i = face.interiorIndex();
342 const auto j = face.exteriorIndex();
343 const auto& interiorPos = stencil.subControlVolume(i).globalPos();
344 const auto& exteriorPos = stencil.subControlVolume(j).globalPos();
345 const auto& integrationPos = face.integrationPos();
346
347 interiorDistance = 0.0;
348 exteriorDistance = 0.0;
349 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
350 interiorDistance += (interiorPos[dimIdx] - integrationPos[dimIdx]) * normal[dimIdx];
351 exteriorDistance += (exteriorPos[dimIdx] - integrationPos[dimIdx]) * normal[dimIdx];
352 }
353
354 interiorDistance = std::sqrt(std::abs(interiorDistance));
355 exteriorDistance = std::sqrt(std::abs(exteriorDistance));
356 }
357};
358
359} // namespace Opm
360
361#endif
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:292
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:214
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:271
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