p1fegradientcalculator.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_P1FE_GRADIENT_CALCULATOR_HH
29#define EWOMS_P1FE_GRADIENT_CALCULATOR_HH
30
31#include <dune/common/fvector.hh>
32
33#include <dune/geometry/type.hh>
34
35#include <dune/common/version.hh>
38
39#if HAVE_DUNE_LOCALFUNCTIONS
40#include <dune/localfunctions/lagrange/pqkfactory.hh>
41#include <array>
42#include <cstddef>
43#include <vector>
44#endif // HAVE_DUNE_LOCALFUNCTIONS
45
46#include <stdexcept>
47#include <type_traits>
48
49namespace Opm {
50
60template<class TypeTag>
62{
68
69 enum { dim = GridView::dimension };
70
71 // set the maximum number of degrees of freedom and the maximum
72 // number of flux approximation points per elements. For this, we
73 // assume cubes as the type of element with the most vertices...
74 enum { maxDof = (2 << dim) };
75 enum { maxFap = maxDof };
76
77 using CoordScalar = typename GridView::ctype;
78 using DimVector = Dune::FieldVector<Scalar, dim>;
79
80#if HAVE_DUNE_LOCALFUNCTIONS
81 using LocalFiniteElementCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
82 using LocalFiniteElement = typename LocalFiniteElementCache::FiniteElementType;
83 using LocalBasisTraits = typename LocalFiniteElement::Traits::LocalBasisType::Traits;
84 using ShapeJacobian = typename LocalBasisTraits::JacobianType;
85#endif // HAVE_DUNE_LOCALFUNCTIONS
86
87public:
94 template <bool prepareValues = true, bool prepareGradients = true>
95 void prepare([[maybe_unused]] const ElementContext& elemCtx,
96 [[maybe_unused]] unsigned timeIdx)
97 {
98 if constexpr (getPropValue<TypeTag, Properties::UseP1FiniteElementGradients>()) {
99#if !HAVE_DUNE_LOCALFUNCTIONS
100 // The dune-localfunctions module is required for P1 finite element gradients
101 throw std::logic_error("The dune-localfunctions module is required in order to use"
102 " finite element gradients");
103#else
104 const auto& stencil = elemCtx.stencil(timeIdx);
105
106 const LocalFiniteElement& localFE = feCache_.get(elemCtx.element().type());
107 localFiniteElement_ = &localFE;
108
109 // loop over all face centeres
110 for (unsigned faceIdx = 0; faceIdx < stencil.numInteriorFaces(); ++faceIdx) {
111 const auto& localFacePos = stencil.interiorFace(faceIdx).localPos();
112
113 // Evaluate the P1 shape functions and their gradients at all
114 // flux approximation points.
115 if (prepareValues) {
116 localFE.localBasis().evaluateFunction(localFacePos, p1Value_[faceIdx]);
117 }
118
119 if (prepareGradients) {
120 // first, get the shape function's gradient in local coordinates
121 std::vector<ShapeJacobian> localGradient;
122 localFE.localBasis().evaluateJacobian(localFacePos, localGradient);
123
124 // convert to a gradient in global space by
125 // multiplying with the inverse transposed jacobian of
126 // the position
127 const auto& geom = elemCtx.element().geometry();
128 const auto& jacInvT =
129 geom.jacobianInverseTransposed(localFacePos);
130
131 std::size_t numVertices = elemCtx.numDof(timeIdx);
132 for (unsigned vertIdx = 0; vertIdx < numVertices; vertIdx++) {
133 jacInvT.mv(/*xVector=*/localGradient[vertIdx][0],
134 /*destVector=*/p1Gradient_[faceIdx][vertIdx]);
135 }
136 }
137 }
138#endif
139 }
140 else {
141 ParentType::template prepare<prepareValues, prepareGradients>(elemCtx, timeIdx);
142 }
143 }
144
156 template <class QuantityCallback>
157 auto calculateScalarValue([[maybe_unused]] const ElementContext& elemCtx,
158 [[maybe_unused]] unsigned fapIdx,
159 [[maybe_unused]] const QuantityCallback& quantityCallback) const
160 -> std::remove_reference_t<typename QuantityCallback::ResultType>
161 {
162 if constexpr (getPropValue<TypeTag, Properties::UseP1FiniteElementGradients>()) {
163#if !HAVE_DUNE_LOCALFUNCTIONS
164 // The dune-localfunctions module is required for P1 finite element gradients
165 throw std::logic_error("The dune-localfunctions module is required in order to use"
166 " finite element gradients");
167#else
168 using QuantityConstType = std::remove_reference_t<typename QuantityCallback::ResultType>;
169 using QuantityType = std::remove_const_t<QuantityConstType>;
170 using Toolbox = MathToolbox<QuantityType>;
171
172 // If the user does not want to use two-point gradients, we
173 // use P1 finite element gradients..
174 QuantityType value(0.0);
175 for (unsigned vertIdx = 0; vertIdx < elemCtx.numDof(/*timeIdx=*/0); ++vertIdx) {
176 if (std::is_same_v<QuantityType, Scalar> ||
177 elemCtx.focusDofIndex() == vertIdx)
178 {
179 value += quantityCallback(vertIdx)*p1Value_[fapIdx][vertIdx];
180 }
181 else {
182 value += Toolbox::value(quantityCallback(vertIdx))*p1Value_[fapIdx][vertIdx];
183 }
184 }
185
186 return value;
187#endif
188 }
189 else {
190 return ParentType::calculateScalarValue(elemCtx, fapIdx, quantityCallback);
191 }
192 }
193
205 template <class QuantityCallback>
206 auto calculateVectorValue([[maybe_unused]] const ElementContext& elemCtx,
207 [[maybe_unused]] unsigned fapIdx,
208 [[maybe_unused]] const QuantityCallback& quantityCallback) const
209 -> std::remove_reference_t<typename QuantityCallback::ResultType>
210 {
211 if constexpr (getPropValue<TypeTag, Properties::UseP1FiniteElementGradients>()) {
212#if !HAVE_DUNE_LOCALFUNCTIONS
213 // The dune-localfunctions module is required for P1 finite element gradients
214 throw std::logic_error("The dune-localfunctions module is required in order to use"
215 " finite element gradients");
216#else
217 using QuantityConstType = std::remove_reference_t<typename QuantityCallback::ResultType>;
218 using QuantityType = std::remove_const_t<QuantityConstType>;
219
220 using RawFieldType = decltype(std::declval<QuantityType>()[0]);
221 using FieldType = std::remove_const_t<std::remove_reference_t<RawFieldType>>;
222
223 using Toolbox = MathToolbox<FieldType>;
224
225 // If the user does not want to use two-point gradients, we
226 // use P1 finite element gradients..
227 QuantityType value(0.0);
228 for (unsigned vertIdx = 0; vertIdx < elemCtx.numDof(/*timeIdx=*/0); ++vertIdx) {
229 if (std::is_same_v<QuantityType, Scalar> ||
230 elemCtx.focusDofIndex() == vertIdx)
231 {
232 const auto& tmp = quantityCallback(vertIdx);
233 for (unsigned k = 0; k < tmp.size(); ++k) {
234 value[k] += tmp[k]*p1Value_[fapIdx][vertIdx];
235 }
236 }
237 else {
238 const auto& tmp = quantityCallback(vertIdx);
239 for (unsigned k = 0; k < tmp.size(); ++k) {
240 value[k] += Toolbox::value(tmp[k])*p1Value_[fapIdx][vertIdx];
241 }
242 }
243 }
244
245 return value;
246#endif
247 }
248 else {
249 return ParentType::calculateVectorValue(elemCtx, fapIdx, quantityCallback);
250 }
251 }
252
264 template <class QuantityCallback, class EvalDimVector>
265 void calculateGradient([[maybe_unused]] EvalDimVector& quantityGrad,
266 [[maybe_unused]] const ElementContext& elemCtx,
267 [[maybe_unused]] unsigned fapIdx,
268 [[maybe_unused]] const QuantityCallback& quantityCallback) const
269 {
270 if constexpr (getPropValue<TypeTag, Properties::UseP1FiniteElementGradients>()) {
271#if !HAVE_DUNE_LOCALFUNCTIONS
272 // The dune-localfunctions module is required for P1 finite element gradients
273 throw std::logic_error("The dune-localfunctions module is required in order to use"
274 " finite element gradients");
275#else
276 using QuantityConstType = std::remove_reference_t<typename QuantityCallback::ResultType>;
277 using QuantityType = std::remove_const_t<QuantityConstType>;
278
279 // If the user does not want two-point gradients, we use P1 finite element
280 // gradients...
281 quantityGrad = 0.0;
282 for (unsigned vertIdx = 0; vertIdx < elemCtx.numDof(/*timeIdx=*/0); ++vertIdx) {
283 if (std::is_same_v<QuantityType, Scalar> ||
284 elemCtx.focusDofIndex() == vertIdx)
285 {
286 const auto& dofVal = quantityCallback(vertIdx);
287 const auto& tmp = p1Gradient_[fapIdx][vertIdx];
288 for (int dimIdx = 0; dimIdx < dim; ++dimIdx) {
289 quantityGrad[dimIdx] += dofVal * tmp[dimIdx];
290 }
291 }
292 else {
293 const auto& dofVal = quantityCallback(vertIdx);
294 const auto& tmp = p1Gradient_[fapIdx][vertIdx];
295 for (int dimIdx = 0; dimIdx < dim; ++dimIdx) {
296 quantityGrad[dimIdx] += scalarValue(dofVal) * tmp[dimIdx];
297 }
298 }
299 }
300#endif
301 }
302 else {
303 ParentType::calculateGradient(quantityGrad, elemCtx, fapIdx, quantityCallback);
304 }
305 }
306
321 template <class QuantityCallback>
322 auto calculateBoundaryValue(const ElementContext& elemCtx,
323 unsigned fapIdx,
324 const QuantityCallback& quantityCallback)
325 -> decltype(ParentType::calculateBoundaryValue(elemCtx, fapIdx, quantityCallback))
326 { return ParentType::calculateBoundaryValue(elemCtx, fapIdx, quantityCallback); }
327
342 template <class QuantityCallback, class EvalDimVector>
343 void calculateBoundaryGradient(EvalDimVector& quantityGrad,
344 const ElementContext& elemCtx,
345 unsigned fapIdx,
346 const QuantityCallback& quantityCallback) const
347 { ParentType::calculateBoundaryGradient(quantityGrad, elemCtx, fapIdx, quantityCallback); }
348
349#if HAVE_DUNE_LOCALFUNCTIONS
350 static LocalFiniteElementCache& localFiniteElementCache()
351 { return feCache_; }
352#endif
353
354private:
355#if HAVE_DUNE_LOCALFUNCTIONS
356 static LocalFiniteElementCache feCache_;
357
358 const LocalFiniteElement* localFiniteElement_{nullptr};
359 std::array<std::vector<Dune::FieldVector<Scalar, 1>>, maxFap> p1Value_{};
360 std::array<std::array<DimVector, maxDof>, maxFap>{};
361#endif // HAVE_DUNE_LOCALFUNCTIONS
362};
363
364#if HAVE_DUNE_LOCALFUNCTIONS
365template<class TypeTag>
366typename P1FeGradientCalculator<TypeTag>::LocalFiniteElementCache
367P1FeGradientCalculator<TypeTag>::feCache_;
368#endif
369
370} // namespace Opm
371
372#endif
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
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
This class calculates gradients of arbitrary quantities at flux integration points using first order ...
Definition: p1fegradientcalculator.hh:62
auto calculateVectorValue(const ElementContext &elemCtx, unsigned fapIdx, const QuantityCallback &quantityCallback) const -> std::remove_reference_t< typename QuantityCallback::ResultType >
Calculates the value of an arbitrary quantity at any interior flux approximation point.
Definition: p1fegradientcalculator.hh:206
void calculateBoundaryGradient(EvalDimVector &quantityGrad, const ElementContext &elemCtx, unsigned fapIdx, const QuantityCallback &quantityCallback) const
Calculates the gradient of an arbitrary quantity at any flux approximation point on the boundary.
Definition: p1fegradientcalculator.hh:343
auto calculateScalarValue(const ElementContext &elemCtx, unsigned fapIdx, const QuantityCallback &quantityCallback) const -> std::remove_reference_t< typename QuantityCallback::ResultType >
Calculates the value of an arbitrary quantity at any interior flux approximation point.
Definition: p1fegradientcalculator.hh:157
auto calculateBoundaryValue(const ElementContext &elemCtx, unsigned fapIdx, const QuantityCallback &quantityCallback) -> decltype(ParentType::calculateBoundaryValue(elemCtx, fapIdx, quantityCallback))
Calculates the value of an arbitrary quantity at any flux approximation point on the grid boundary.
Definition: p1fegradientcalculator.hh:322
void prepare(const ElementContext &elemCtx, unsigned timeIdx)
Precomputes the common values to calculate gradients and values of quantities at any flux approximati...
Definition: p1fegradientcalculator.hh:95
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: p1fegradientcalculator.hh:265
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
Declares the basic properties used by the common infrastructure of the vertex-centered finite volume ...