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 ShapeGradients = std::vector<GlobalPosition>;
28
29public:
30 // The default constructor
32
33 // The constructor
34 template<class Geometry>
35 FEInterpolationPointData(const Geometry& geometry,
36 const LocalPosition& local,
37 const LocalBasis& localBasis)
38 : local_(local),
39 global_(geometry.global(local))
40 {
41 auto numLocalDofs = localBasis.size();
42
43 // set the shape values
44 shapeValues_.resize(numLocalDofs);
45 localBasis.evaluateFunction(local, shapeValues_);
46
47 // the local shape function gradients
48 std::vector<JacobianType> shapeGrads(numLocalDofs);
49 localBasis.evaluateJacobian(local, shapeGrads);
50
51 // the global shape function gradients
52 const auto jacInvT = geometry.jacobianInverseTransposed(local);
53 shapeGradients_.resize(numLocalDofs, GlobalPosition(0.0));
54 for (unsigned int i = 0; i < numLocalDofs; ++i)
55 jacInvT.umv(shapeGrads[i][0], shapeGradients_[i]);
56 }
57
59 const ShapeValues& shapeValues() const
60 { return shapeValues_; }
61
63 const RangeType& shapeValue(int i) const
64 { return shapeValues_[i]; }
65
67 const ShapeGradients& shapeGradients() const
68 { return shapeGradients_; }
69
71 const GlobalPosition& gradN(int i) const
72 { return shapeGradients_[i]; }
73
75 const LocalPosition& local() const
76 { return local_; }
77
79 const GlobalPosition& global() const
80 { return global_; }
81
82private:
83 LocalPosition local_;
84 GlobalPosition global_;
85
86 ShapeValues shapeValues_;
87 ShapeGradients shapeGradients_;
88};
89
94template<class GlobalPosition, class LocalBasis, class BoundaryFlag>
95class FEFaceInterpolationPointData : public FEInterpolationPointData<GlobalPosition, LocalBasis>
96{
98 using LocalPosition = typename LocalBasis::Traits::DomainType;
99public:
100 // The default constructor
102
103 // The constructor
104 template<class Geometry>
105 FEFaceInterpolationPointData(const Geometry& geometry,
106 const LocalPosition& local,
107 const LocalBasis& localBasis,
108 const GlobalPosition& n,
109 const BoundaryFlag& bFlag)
110 : ParentType(geometry, local, localBasis), normal_(n), boundaryFlag_(bFlag)
111 {}
112
114 const GlobalPosition& unitOuterNormal() const
115 { return normal_; }
116
119 { return boundaryFlag_.get(); }
120
121private:
122 const GlobalPosition& normal_;
123 BoundaryFlag boundaryFlag_;
124};
125
126
127} // end namespace Dumux
128
129#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
Interpolation point data related to a face of an element.
Definition: fem/interpolationpointdata.hh:96
BoundaryFlag::value_type boundaryFlag() const
Return the boundary flag.
Definition: fem/interpolationpointdata.hh:118
FEFaceInterpolationPointData(const Geometry &geometry, const LocalPosition &local, const LocalBasis &localBasis, const GlobalPosition &n, const BoundaryFlag &bFlag)
Definition: fem/interpolationpointdata.hh:105
const GlobalPosition & unitOuterNormal() const
The unit outer normal vector at the quadrature point.
Definition: fem/interpolationpointdata.hh:114
Definition: fem/interpolationpointdata.hh:21
const ShapeValues & shapeValues() const
The shape values at the quadrature point.
Definition: fem/interpolationpointdata.hh:59
const LocalPosition & local() const
The local position of the quadrature point.
Definition: fem/interpolationpointdata.hh:75
const GlobalPosition & global() const
The global position of the quadrature point.
Definition: fem/interpolationpointdata.hh:79
const RangeType & shapeValue(int i) const
The shape value of a local dof at the quadrature point.
Definition: fem/interpolationpointdata.hh:63
FEInterpolationPointData(const Geometry &geometry, const LocalPosition &local, const LocalBasis &localBasis)
Definition: fem/interpolationpointdata.hh:35
const ShapeGradients & shapeGradients() const
The shape value gradients at the quadrature point.
Definition: fem/interpolationpointdata.hh:67
const GlobalPosition & gradN(int i) const
The shape value gradient of a local dof at the quadrature point.
Definition: fem/interpolationpointdata.hh:71
Definition: adapt.hh:17