version 3.11-dev
fem/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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_FE_INTEGRATION_POINT_DATA_HH
13#define DUMUX_DISCRETIZATION_FE_INTEGRATION_POINT_DATA_HH
14
15#include <vector>
16
17namespace Dumux {
18
19template<class GlobalPosition, class LocalBasis>
21{
22 using LocalPosition = typename LocalBasis::Traits::DomainType;
23 using RangeType = typename LocalBasis::Traits::RangeType;
24 using JacobianType = typename LocalBasis::Traits::JacobianType;
25
26 using ShapeValues = std::vector<RangeType>;
27 using ShapeJacobians = std::vector<JacobianType>;
28 using Gradients = std::vector<GlobalPosition>;
29
30public:
31 // The default constructor
33
34 // The constructor
35 template<class Geometry>
36 FEInterpolationPointData(const Geometry& geometry,
37 const LocalPosition& local,
38 const GlobalPosition& global,
39 const LocalBasis& localBasis)
40 : local_(local),
41 global_(global)
42 {
43 auto numLocalDofs = localBasis.size();
44
45 // set the shape values
46 shapeValues_.resize(numLocalDofs);
47 localBasis.evaluateFunction(local, shapeValues_);
48
49 // the local shape function gradients
50 shapeJacobians_.resize(numLocalDofs);
51 localBasis.evaluateJacobian(local, shapeJacobians_);
52
53 // the global shape function gradients
54 const auto jacInvT = geometry.jacobianInverseTransposed(local);
55 gradients_.resize(numLocalDofs, GlobalPosition(0.0));
56 for (unsigned int i = 0; i < numLocalDofs; ++i)
57 jacInvT.umv(shapeJacobians_[i][0], gradients_[i]);
58 }
59
60 // The constructor
61 template<class Geometry>
62 FEInterpolationPointData(const Geometry& geometry,
63 const LocalPosition& local,
64 const LocalBasis& localBasis)
65 : FEInterpolationPointData(geometry, local, geometry.global(local), localBasis)
66 { }
67
69 const ShapeValues& shapeValues() const
70 { return shapeValues_; }
71
73 const RangeType& shapeValue(int i) const
74 { return shapeValues_[i]; }
75
77 const ShapeJacobians& shapeJacobians() const
78 { return shapeJacobians_; }
79
81 const GlobalPosition& gradN(int i) const
82 { return gradients_[i]; }
83
85 const LocalPosition& local() const
86 { return local_; }
87
89 const GlobalPosition& global() const
90 { return global_; }
91
92private:
93 LocalPosition local_;
94 GlobalPosition global_;
95
96 ShapeValues shapeValues_;
97 ShapeJacobians shapeJacobians_;
98 Gradients gradients_;
99};
100
105template<class BaseClass, class BoundaryFlag, class LocalIndex>
106class FEFaceInterpolationPointData : public BaseClass
107{
108public:
109 using GlobalPosition = std::remove_cvref_t<decltype(std::declval<BaseClass>().global())>;
110 using LocalPosition = std::remove_cvref_t<decltype(std::declval<BaseClass>().local())>;
111
112 template<class... Args>
113 FEFaceInterpolationPointData(GlobalPosition&& n, BoundaryFlag&& bFlag, LocalIndex index, Args&&... args)
114 : BaseClass(std::forward<Args>(args)...), normal_(std::move(n)), boundaryFlag_(std::move(bFlag)), index_(index) {}
115
116 template<class... Args>
117 FEFaceInterpolationPointData(const GlobalPosition& n, const BoundaryFlag& bFlag, LocalIndex index, Args&&... args)
118 : BaseClass(std::forward<Args>(args)...), normal_(n), boundaryFlag_(bFlag), index_(index) {}
119
122 { return normal_; }
123
126 { return boundaryFlag_.get(); }
127
129 LocalIndex intersectionIndex() const
130 { return index_; }
131
132private:
133 const GlobalPosition& normal_;
134 BoundaryFlag boundaryFlag_;
135 LocalIndex index_;
136};
137
138} // end namespace Dumux
139
140#endif
Boundary flag to store e.g. in sub control volume faces.
Definition: boundaryflag.hh:58
std::size_t value_type
Definition: boundaryflag.hh:28
value_type get() const
Definition: boundaryflag.hh:41
An interpolation point related to an intersection.
Definition: fem/interpolationpointdata.hh:107
FEFaceInterpolationPointData(GlobalPosition &&n, BoundaryFlag &&bFlag, LocalIndex index, Args &&... args)
Definition: fem/interpolationpointdata.hh:113
std::remove_cvref_t< decltype(std::declval< BaseClass >().local())> LocalPosition
Definition: fem/interpolationpointdata.hh:110
FEFaceInterpolationPointData(const GlobalPosition &n, const BoundaryFlag &bFlag, LocalIndex index, Args &&... args)
Definition: fem/interpolationpointdata.hh:117
std::remove_cvref_t< decltype(std::declval< BaseClass >().global())> GlobalPosition
Definition: fem/interpolationpointdata.hh:109
const GlobalPosition & unitOuterNormal() const
The unit outer normal vector at the quadrature point.
Definition: fem/interpolationpointdata.hh:121
BoundaryFlag::value_type boundaryFlag() const
Return the boundary flag.
Definition: fem/interpolationpointdata.hh:125
LocalIndex intersectionIndex() const
The local index of an intersection (default is intersection.indexInInside())
Definition: fem/interpolationpointdata.hh:129
Definition: fem/interpolationpointdata.hh:21
const ShapeValues & shapeValues() const
The shape values at the quadrature point.
Definition: fem/interpolationpointdata.hh:69
const LocalPosition & local() const
The local position of the quadrature point.
Definition: fem/interpolationpointdata.hh:85
const ShapeJacobians & shapeJacobians() const
The shape value gradients at the quadrature point.
Definition: fem/interpolationpointdata.hh:77
const GlobalPosition & global() const
The global position of the quadrature point.
Definition: fem/interpolationpointdata.hh:89
const RangeType & shapeValue(int i) const
The shape value of a local dof at the quadrature point.
Definition: fem/interpolationpointdata.hh:73
FEInterpolationPointData(const Geometry &geometry, const LocalPosition &local, const GlobalPosition &global, const LocalBasis &localBasis)
Definition: fem/interpolationpointdata.hh:36
FEInterpolationPointData(const Geometry &geometry, const LocalPosition &local, const LocalBasis &localBasis)
Definition: fem/interpolationpointdata.hh:62
const GlobalPosition & gradN(int i) const
The shape value gradient of a local dof at the quadrature point.
Definition: fem/interpolationpointdata.hh:81
Definition: adapt.hh:17