version 3.11-dev
cvfe/interpolationpointdata.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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_CVFE_IP_DATA_HH
13#define DUMUX_CVFE_IP_DATA_HH
14
15#include <type_traits>
16#include <utility>
17#include <vector>
18
19#include <dumux/common/typetraits/localdofs_.hh>
21#include <dumux/common/concepts/ipdata_.hh>
22
23namespace Dumux::CVFE {
24
29template<class LocalPos, class GlobalPos>
31{
32public:
33 using LocalPosition = LocalPos;
34 using GlobalPosition = GlobalPos;
35
37 : local_(std::move(localPos)), global_(std::move(pos)) {}
39 : local_(localPos), global_(pos) {}
40
42 const GlobalPosition& global() const
43 { return global_; }
44
46 const LocalPosition& local() const
47 { return local_; }
48
49private:
50 LocalPosition local_;
51 GlobalPosition global_;
52};
53
58template<class LocalPosition, class GlobalPosition, class LocalIndex>
59class LocalDofInterpolationPointData : public InterpolationPointData<LocalPosition, GlobalPosition>
60{
62public:
63 LocalDofInterpolationPointData(LocalPosition&& localPos, GlobalPosition&& pos, LocalIndex index)
64 : ParentType(localPos, pos), localDofIndex_(index) {}
65 LocalDofInterpolationPointData(const LocalPosition& localPos, const GlobalPosition& pos, LocalIndex index)
66 : ParentType(localPos, pos), localDofIndex_(index) {}
67
69 LocalIndex localDofIndex() const
70 { return localDofIndex_; }
71
72private:
73 LocalIndex localDofIndex_;
74};
75
80template<class LocalMapping, class GlobalPos>
82{
83public:
84 using LocalPosition = std::invoke_result_t<LocalMapping, const GlobalPos&>;
85 using GlobalPosition = GlobalPos;
86
87 InterpolationPointDataLocalMapping(LocalMapping&& mapping, GlobalPosition&& pos) : localMapping_(std::move(mapping)), global_(std::move(pos)) {}
88 InterpolationPointDataLocalMapping(LocalMapping&& mapping, const GlobalPosition& pos) : localMapping_(std::move(mapping)), global_(pos) {}
89
91 const GlobalPosition& global() const
92 { return global_; }
93
95 const LocalPosition local() const
96 { return localMapping_(global_); }
97
98private:
99 LocalMapping localMapping_;
100 GlobalPosition global_;
101};
102
107template<class BaseClass, class LocalIndex>
108class FaceInterpolationPointData : public BaseClass
109{
110public:
111 using GlobalPosition = std::remove_cvref_t<decltype(std::declval<BaseClass>().global())>;
112 using LocalPosition = std::remove_cvref_t<decltype(std::declval<BaseClass>().local())>;
113
114 template<class... Args>
115 FaceInterpolationPointData(GlobalPosition&& n, LocalIndex index, Args&&... args)
116 : BaseClass(std::forward<Args>(args)...), normal_(std::move(n)), scvfIndex_(index) {}
117
118 template<class... Args>
119 FaceInterpolationPointData(const GlobalPosition& n, LocalIndex index, Args&&... args)
120 : BaseClass(std::forward<Args>(args)...), normal_(n), scvfIndex_(index) {}
121
124 { return normal_; }
125
127 LocalIndex scvfIndex() const
128 { return scvfIndex_; }
129
130private:
131 GlobalPosition normal_;
132 LocalIndex scvfIndex_;
133};
134
139template<class IpData>
141{
142public:
143 template<class... Args>
144 IndexedQuadratureInterpolationPointData(std::size_t qpIdx, Args&&... args)
145 : IpData(std::forward<Args>(args)...), qpIndex_(qpIdx) {}
146
148 std::size_t qpIndex() const
149 { return qpIndex_; }
150
151private:
152 std::size_t qpIndex_;
153};
154
159template< class GridGeometry >
161{
162 using GridView = typename GridGeometry::GridView;
163 using Element = typename GridView::template Codim<0>::Entity;
164 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
165
166 using LocalBasis = typename GridGeometry::FeCache::FiniteElementType::Traits::LocalBasisType;
167 using RangeType = typename LocalBasis::Traits::RangeType;
168 using JacobianType = typename LocalBasis::Traits::JacobianType;
169
170 using Gradients = std::vector<GlobalPosition>;
171
172 using FVElementGeometry = typename GridGeometry::LocalView;
173
174public:
176 static constexpr bool isSolDependent = false;
177
179 template< class Problem, class ElementVariables >
180 void update(const Problem& problem,
181 const Element& element,
182 const FVElementGeometry& fvGeometry,
183 const ElementVariables& elemVars,
184 const GlobalPosition& globalPos)
185 {
186 update_(fvGeometry, fvGeometry.elementGeometry().local(globalPos));
187 }
188
190 template< class Problem, class ElementVariables, Concept::IpData IpData >
191 void update(const Problem& problem,
192 const Element& element,
193 const FVElementGeometry& fvGeometry,
194 const ElementVariables& elemVars,
195 const IpData& ipData)
196 {
197 update_(fvGeometry, ipData.local());
198 }
199
201 const std::vector<JacobianType>& shapeJacobian() const { return shapeJacobian_; }
203 const std::vector<RangeType>& shapeValues() const { return shapeValues_; }
205 const GlobalPosition& gradN(unsigned int localIdx) const { return gradN_[localIdx]; }
206
207private:
209 template<class LocalPosition>
210 void update_(const FVElementGeometry& fvGeometry,
211 const LocalPosition& localPos)
212 {
213 const auto& geometry = fvGeometry.elementGeometry();
214 const auto& localBasis = fvGeometry.feLocalBasis();
215
216 // evaluate shape functions and gradients at the interpolation point
217 const auto jacInvT = geometry.jacobianInverseTransposed(localPos);
218 localBasis.evaluateJacobian(localPos, shapeJacobian_);
219 localBasis.evaluateFunction(localPos, shapeValues_); // shape values for rho
220
221 // compute the gradN for every local dof
222 gradN_.resize(Detail::LocalDofs::numLocalDofs(fvGeometry));
223 for (const auto& localDof: localDofs(fvGeometry))
224 jacInvT.mv(shapeJacobian_[localDof.index()][0], gradN_[localDof.index()]);
225 }
226
227 Gradients gradN_;
228 std::vector<JacobianType> shapeJacobian_;
229 std::vector<RangeType> shapeValues_;
230};
231
232} // end namespace Dumux::CVFE
233
234#endif
An interpolation point related to a face of an element.
Definition: cvfe/interpolationpointdata.hh:109
LocalIndex scvfIndex() const
The local index of an scvf.
Definition: cvfe/interpolationpointdata.hh:127
FaceInterpolationPointData(const GlobalPosition &n, LocalIndex index, Args &&... args)
Definition: cvfe/interpolationpointdata.hh:119
FaceInterpolationPointData(GlobalPosition &&n, LocalIndex index, Args &&... args)
Definition: cvfe/interpolationpointdata.hh:115
const GlobalPosition & unitOuterNormal() const
The unit outer normal vector at the quadrature point.
Definition: cvfe/interpolationpointdata.hh:123
std::remove_cvref_t< decltype(std::declval< BaseClass >().local())> LocalPosition
Definition: cvfe/interpolationpointdata.hh:112
std::remove_cvref_t< decltype(std::declval< BaseClass >().global())> GlobalPosition
Definition: cvfe/interpolationpointdata.hh:111
Wraps interpolation point data and adds a quadrature point index for use in quadrature loops.
Definition: cvfe/interpolationpointdata.hh:141
std::size_t qpIndex() const
The quadrature point index.
Definition: cvfe/interpolationpointdata.hh:148
IndexedQuadratureInterpolationPointData(std::size_t qpIdx, Args &&... args)
Definition: cvfe/interpolationpointdata.hh:144
An interpolation point related to an element that includes global and local positions.
Definition: cvfe/interpolationpointdata.hh:31
const LocalPosition & local() const
The local position of the quadrature point.
Definition: cvfe/interpolationpointdata.hh:46
GlobalPos GlobalPosition
Definition: cvfe/interpolationpointdata.hh:34
const GlobalPosition & global() const
The global position of the quadrature point.
Definition: cvfe/interpolationpointdata.hh:42
InterpolationPointData(LocalPosition &&localPos, GlobalPosition &&pos)
Definition: cvfe/interpolationpointdata.hh:36
LocalPos LocalPosition
Definition: cvfe/interpolationpointdata.hh:33
InterpolationPointData(const LocalPosition &localPos, const GlobalPosition &pos)
Definition: cvfe/interpolationpointdata.hh:38
An interpolation point related to a global position of an element, giving its local positions by a ma...
Definition: cvfe/interpolationpointdata.hh:82
InterpolationPointDataLocalMapping(LocalMapping &&mapping, GlobalPosition &&pos)
Definition: cvfe/interpolationpointdata.hh:87
const LocalPosition local() const
The local position of the quadrature point.
Definition: cvfe/interpolationpointdata.hh:95
InterpolationPointDataLocalMapping(LocalMapping &&mapping, const GlobalPosition &pos)
Definition: cvfe/interpolationpointdata.hh:88
std::invoke_result_t< LocalMapping, const GlobalPos & > LocalPosition
Definition: cvfe/interpolationpointdata.hh:84
const GlobalPosition & global() const
The global position of the quadrature point.
Definition: cvfe/interpolationpointdata.hh:91
GlobalPos GlobalPosition
Definition: cvfe/interpolationpointdata.hh:85
Interpolation point data related to a local basis.
Definition: cvfe/interpolationpointdata.hh:161
static constexpr bool isSolDependent
whether the cache needs an update when the solution changes
Definition: cvfe/interpolationpointdata.hh:176
const GlobalPosition & gradN(unsigned int localIdx) const
returns the shape function gradients in global coordinates at the interpolation point
Definition: cvfe/interpolationpointdata.hh:205
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const GlobalPosition &globalPos)
update the cache for a given global position
Definition: cvfe/interpolationpointdata.hh:180
const std::vector< JacobianType > & shapeJacobian() const
returns the shape function gradients in local coordinates at the interpolation point
Definition: cvfe/interpolationpointdata.hh:201
const std::vector< RangeType > & shapeValues() const
returns the shape function values at the interpolation point
Definition: cvfe/interpolationpointdata.hh:203
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const IpData &ipData)
update the cache for interpolation point data
Definition: cvfe/interpolationpointdata.hh:191
An interpolation point related to a localDof of an element, giving its global and local positions.
Definition: cvfe/interpolationpointdata.hh:60
LocalDofInterpolationPointData(const LocalPosition &localPos, const GlobalPosition &pos, LocalIndex index)
Definition: cvfe/interpolationpointdata.hh:65
LocalDofInterpolationPointData(LocalPosition &&localPos, GlobalPosition &&pos, LocalIndex index)
Definition: cvfe/interpolationpointdata.hh:63
LocalIndex localDofIndex() const
The local index of the corresponding dof.
Definition: cvfe/interpolationpointdata.hh:69
Class representing dofs on elements for control-volume finite element schemes.
Definition: cvfe/interpolationpointdata.hh:23
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition: localdof.hh:50