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:
95 template <bool prepareValues = true, bool prepareGradients = true>
96 void prepare([[maybe_unused]] const ElementContext& elemCtx,
97 [[maybe_unused]] unsigned timeIdx)
98 {
99 if constexpr (getPropValue<TypeTag, Properties::UseP1FiniteElementGradients>()) {
100#if !HAVE_DUNE_LOCALFUNCTIONS
101 // The dune-localfunctions module is required for P1 finite element gradients
102 throw std::logic_error("The dune-localfunctions module is required in order to use"
103 " finite element gradients");
104#else
105 const auto& stencil = elemCtx.stencil(timeIdx);
106
107 const LocalFiniteElement& localFE = feCache_.get(elemCtx.element().type());
108 localFiniteElement_ = &localFE;
109
110 // loop over all face centeres
111 for (unsigned faceIdx = 0; faceIdx < stencil.numInteriorFaces(); ++faceIdx) {
112 const auto& localFacePos = stencil.interiorFace(faceIdx).localPos();
113
114 // Evaluate the P1 shape functions and their gradients at all
115 // flux approximation points.
116 if (prepareValues) {
117 localFE.localBasis().evaluateFunction(localFacePos, p1Value_[faceIdx]);
118 }
119
120 if (prepareGradients) {
121 // first, get the shape function's gradient in local coordinates
122 std::vector<ShapeJacobian> localGradient;
123 localFE.localBasis().evaluateJacobian(localFacePos, localGradient);
124
125 // convert to a gradient in global space by
126 // multiplying with the inverse transposed jacobian of
127 // the position
128 const auto& geom = elemCtx.element().geometry();
129 const auto& jacInvT =
130 geom.jacobianInverseTransposed(localFacePos);
131
132 std::size_t numVertices = elemCtx.numDof(timeIdx);
133 for (unsigned vertIdx = 0; vertIdx < numVertices; vertIdx++) {
134 jacInvT.mv(/*xVector=*/localGradient[vertIdx][0],
135 /*destVector=*/p1Gradient_[faceIdx][vertIdx]);
136 }
137 }
138 }
139#endif
140 }
141 else {
142 ParentType::template prepare<prepareValues, prepareGradients>(elemCtx, timeIdx);
143 }
144 }
145
157 template <class QuantityCallback>
158 auto calculateScalarValue([[maybe_unused]] const ElementContext& elemCtx,
159 [[maybe_unused]] unsigned fapIdx,
160 [[maybe_unused]] const QuantityCallback& quantityCallback) const
161 -> std::remove_reference_t<typename QuantityCallback::ResultType>
162 {
163 if constexpr (getPropValue<TypeTag, Properties::UseP1FiniteElementGradients>()) {
164#if !HAVE_DUNE_LOCALFUNCTIONS
165 // The dune-localfunctions module is required for P1 finite element gradients
166 throw std::logic_error("The dune-localfunctions module is required in order to use"
167 " finite element gradients");
168#else
169 using QuantityConstType = std::remove_reference_t<typename QuantityCallback::ResultType>;
170 using QuantityType = std::remove_const_t<QuantityConstType>;
171 using Toolbox = MathToolbox<QuantityType>;
172
173 // If the user does not want to use two-point gradients, we
174 // use P1 finite element gradients..
175 QuantityType value(0.0);
176 for (unsigned vertIdx = 0; vertIdx < elemCtx.numDof(/*timeIdx=*/0); ++vertIdx) {
177 if (std::is_same_v<QuantityType, Scalar> ||
178 elemCtx.focusDofIndex() == vertIdx)
179 {
180 value += quantityCallback(vertIdx)*p1Value_[fapIdx][vertIdx];
181 }
182 else {
183 value += Toolbox::value(quantityCallback(vertIdx))*p1Value_[fapIdx][vertIdx];
184 }
185 }
186
187 return value;
188#endif
189 }
190 else {
191 return ParentType::calculateScalarValue(elemCtx, fapIdx, quantityCallback);
192 }
193 }
194
206 template <class QuantityCallback>
207 auto calculateVectorValue([[maybe_unused]] const ElementContext& elemCtx,
208 [[maybe_unused]] unsigned fapIdx,
209 [[maybe_unused]] const QuantityCallback& quantityCallback) const
210 -> std::remove_reference_t<typename QuantityCallback::ResultType>
211 {
212 if constexpr (getPropValue<TypeTag, Properties::UseP1FiniteElementGradients>()) {
213#if !HAVE_DUNE_LOCALFUNCTIONS
214 // The dune-localfunctions module is required for P1 finite element gradients
215 throw std::logic_error("The dune-localfunctions module is required in order to use"
216 " finite element gradients");
217#else
218 using QuantityConstType = std::remove_reference_t<typename QuantityCallback::ResultType>;
219 using QuantityType = std::remove_const_t<QuantityConstType>;
220
221 using RawFieldType = decltype(std::declval<QuantityType>()[0]);
222 using FieldType = std::remove_const_t<std::remove_reference_t<RawFieldType>>;
223
224 using Toolbox = MathToolbox<FieldType>;
225
226 // If the user does not want to use two-point gradients, we
227 // use P1 finite element gradients..
228 QuantityType value(0.0);
229 for (unsigned vertIdx = 0; vertIdx < elemCtx.numDof(/*timeIdx=*/0); ++vertIdx) {
230 if (std::is_same_v<QuantityType, Scalar> ||
231 elemCtx.focusDofIndex() == vertIdx)
232 {
233 const auto& tmp = quantityCallback(vertIdx);
234 for (unsigned k = 0; k < tmp.size(); ++k) {
235 value[k] += tmp[k]*p1Value_[fapIdx][vertIdx];
236 }
237 }
238 else {
239 const auto& tmp = quantityCallback(vertIdx);
240 for (unsigned k = 0; k < tmp.size(); ++k) {
241 value[k] += Toolbox::value(tmp[k])*p1Value_[fapIdx][vertIdx];
242 }
243 }
244 }
245
246 return value;
247#endif
248 }
249 else {
250 return ParentType::calculateVectorValue(elemCtx, fapIdx, quantityCallback);
251 }
252 }
253
266 template <class QuantityCallback, class EvalDimVector>
267 void calculateGradient([[maybe_unused]] EvalDimVector& quantityGrad,
268 [[maybe_unused]] const ElementContext& elemCtx,
269 [[maybe_unused]] unsigned fapIdx,
270 [[maybe_unused]] const QuantityCallback& quantityCallback) const
271 {
272 if constexpr (getPropValue<TypeTag, Properties::UseP1FiniteElementGradients>()) {
273#if !HAVE_DUNE_LOCALFUNCTIONS
274 // The dune-localfunctions module is required for P1 finite element gradients
275 throw std::logic_error("The dune-localfunctions module is required in order to use"
276 " finite element gradients");
277#else
278 using QuantityConstType = std::remove_reference_t<typename QuantityCallback::ResultType>;
279 using QuantityType = std::remove_const_t<QuantityConstType>;
280
281 // If the user does not want two-point gradients, we use P1 finite element
282 // gradients...
283 quantityGrad = 0.0;
284 for (unsigned vertIdx = 0; vertIdx < elemCtx.numDof(/*timeIdx=*/0); ++vertIdx) {
285 if (std::is_same_v<QuantityType, Scalar> ||
286 elemCtx.focusDofIndex() == vertIdx)
287 {
288 const auto& dofVal = quantityCallback(vertIdx);
289 const auto& tmp = p1Gradient_[fapIdx][vertIdx];
290 for (int dimIdx = 0; dimIdx < dim; ++dimIdx) {
291 quantityGrad[dimIdx] += dofVal * tmp[dimIdx];
292 }
293 }
294 else {
295 const auto& dofVal = quantityCallback(vertIdx);
296 const auto& tmp = p1Gradient_[fapIdx][vertIdx];
297 for (int dimIdx = 0; dimIdx < dim; ++dimIdx) {
298 quantityGrad[dimIdx] += scalarValue(dofVal) * tmp[dimIdx];
299 }
300 }
301 }
302#endif
303 }
304 else {
305 ParentType::calculateGradient(quantityGrad, elemCtx, fapIdx, quantityCallback);
306 }
307 }
308
323 template <class QuantityCallback>
324 auto calculateBoundaryValue(const ElementContext& elemCtx,
325 unsigned fapIdx,
326 const QuantityCallback& quantityCallback)
327 -> decltype(ParentType::calculateBoundaryValue(elemCtx, fapIdx, quantityCallback))
328 { return ParentType::calculateBoundaryValue(elemCtx, fapIdx, quantityCallback); }
329
345 template <class QuantityCallback, class EvalDimVector>
346 void calculateBoundaryGradient(EvalDimVector& quantityGrad,
347 const ElementContext& elemCtx,
348 unsigned fapIdx,
349 const QuantityCallback& quantityCallback) const
350 { ParentType::calculateBoundaryGradient(quantityGrad, elemCtx, fapIdx, quantityCallback); }
351
352#if HAVE_DUNE_LOCALFUNCTIONS
353 static LocalFiniteElementCache& localFiniteElementCache()
354 { return feCache_; }
355#endif
356
357private:
358#if HAVE_DUNE_LOCALFUNCTIONS
359 static LocalFiniteElementCache feCache_;
360
361 const LocalFiniteElement* localFiniteElement_{nullptr};
362 std::array<std::vector<Dune::FieldVector<Scalar, 1>>, maxFap> p1Value_{};
363 std::array<std::array<DimVector, maxDof>, maxFap>{};
364#endif // HAVE_DUNE_LOCALFUNCTIONS
365};
366
367#if HAVE_DUNE_LOCALFUNCTIONS
368template<class TypeTag>
369typename P1FeGradientCalculator<TypeTag>::LocalFiniteElementCache
370P1FeGradientCalculator<TypeTag>::feCache_;
371#endif
372
373} // namespace Opm
374
375#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: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
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
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:207
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:346
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:158
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:324
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:96
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:267
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 ...