DUNE PDELab (git)

crouzeixraviart.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3// SPDX-FileCopyrightInfo: Copyright © DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5#ifndef DUNE_LOCALFUNCTIONS_CROUZEIXRAVIART_HH
6#define DUNE_LOCALFUNCTIONS_CROUZEIXRAVIART_HH
7
8#include <array>
9#include <numeric>
10
13
14#include <dune/geometry/type.hh>
15#include <dune/geometry/referenceelements.hh>
16
17#include <dune/localfunctions/common/localbasis.hh>
18#include <dune/localfunctions/common/localfiniteelementtraits.hh>
19#include <dune/localfunctions/common/localkey.hh>
20
21namespace Dune { namespace Impl
22{
29 template<class D, class R, unsigned int dim>
30 class CrouzeixRaviartLocalBasis
31 {
32 public:
33 using Traits = LocalBasisTraits<D,dim,FieldVector<D,dim>,R,1,FieldVector<R,1>,FieldMatrix<R,1,dim> >;
34
37 static constexpr unsigned int size ()
38 {
39 return dim+1;
40 }
41
43 void evaluateFunction(const typename Traits::DomainType& x,
44 std::vector<typename Traits::RangeType>& out) const
45 {
46 out.resize(size());
47
48 std::fill(out.begin(), out.end()-1, 1.0);
49 out.back() = 1.0-dim;
50
51 for (unsigned int i=0; i<dim; i++)
52 {
53 out[i] -= dim * x[dim-i-1];
54 out.back() += dim*x[i];
55 }
56 }
57
63 void evaluateJacobian(const typename Traits::DomainType& x,
64 std::vector<typename Traits::JacobianType>& out) const
65 {
66 out.resize(size());
67
68 for (unsigned i=0; i<dim; i++)
69 for (unsigned j=0; j<dim; j++)
70 out[i][0][j] = (i==(dim-1-j)) ? -(double)dim : 0;
71
72 std::fill(out.back()[0].begin(), out.back()[0].end(), dim);
73 }
74
81 void partial(const std::array<unsigned int,dim>& order,
82 const typename Traits::DomainType& in,
83 std::vector<typename Traits::RangeType>& out) const
84 {
85 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
86
87 out.resize(size());
88
89 if (totalOrder == 0) {
90 evaluateFunction(in, out);
91 return;
92 }
93
94 if (totalOrder==1)
95 {
96 auto direction = std::find(order.begin(), order.end(), 1)-order.begin();
97
98 for (unsigned int i=0; i<dim; i++)
99 out[i] = (i==(dim-1-direction)) ? -(double)dim : 0.0;
100
101 out.back()[0] = dim;
102 }
103 else // all higher order derivatives are zero
104 std::fill(out.begin(), out.end(), 0);
105 }
106
108 static constexpr unsigned int order ()
109 {
110 return 1;
111 }
112 };
113
118 template<unsigned int dim>
119 class CrouzeixRaviartLocalCoefficients
120 {
121 public:
123 CrouzeixRaviartLocalCoefficients ()
124 : localKeys_(size())
125 {
126 for (std::size_t i=0; i<size(); i++)
127 localKeys_[i] = LocalKey(i,1,0);
128 }
129
131 static constexpr std::size_t size ()
132 {
133 return dim+1;
134 }
135
137 const LocalKey& localKey (std::size_t i) const
138 {
139 return localKeys_[i];
140 }
141
142 private:
143 std::vector<LocalKey> localKeys_;
144 };
145
150 template<class LocalBasis>
151 class CrouzeixRaviartLocalInterpolation
152 {
153 public:
154
162 template<typename F, typename C>
163 void interpolate (const F& f, std::vector<C>& out) const
164 {
165 constexpr auto dim = LocalBasis::Traits::dimDomain;
166
167 out.resize(LocalBasis::size());
168
169 // Evaluate at the facet midpoints
171 for (int i=0; i<refSimplex.size(1); i++)
172 out[i] = f(refSimplex.template geometry<1>(i).center());
173 }
174 };
175
176} } // namespace Dune::Impl
177
178namespace Dune
179{
186 template<class D, class R, int dim>
188 {
189 public:
193 Impl::CrouzeixRaviartLocalCoefficients<dim>,
194 Impl::CrouzeixRaviartLocalInterpolation<Impl::CrouzeixRaviartLocalBasis<D,R,dim> > >;
195
198 const typename Traits::LocalBasisType& localBasis() const
199 {
200 return basis_;
201 }
202
206 {
207 return coefficients_;
208 }
209
213 {
214 return interpolation_;
215 }
216
218 static constexpr std::size_t size()
219 {
220 return dim+1;
221 }
222
225 static constexpr GeometryType type()
226 {
227 return GeometryTypes::simplex(dim);
228 }
229
230 private:
231 Impl::CrouzeixRaviartLocalBasis<D,R,dim> basis_;
232 Impl::CrouzeixRaviartLocalCoefficients<dim> coefficients_;
233 Impl::CrouzeixRaviartLocalInterpolation<Impl::CrouzeixRaviartLocalBasis<D,R,dim> > interpolation_;
234 };
235
236} // namespace Dune
237
238#endif // DUNE_LOCALFUNCTIONS_CROUZEIXRAVIART_HH
Crouzeix-Raviart finite element.
Definition: crouzeixraviart.hh:188
static constexpr GeometryType type()
The reference element that the local finite element is defined on.
Definition: crouzeixraviart.hh:225
const Traits::LocalCoefficientsType & localCoefficients() const
Returns the assignment of the degrees of freedom to the element subentities.
Definition: crouzeixraviart.hh:205
const Traits::LocalInterpolationType & localInterpolation() const
Returns object that evaluates degrees of freedom.
Definition: crouzeixraviart.hh:212
const Traits::LocalBasisType & localBasis() const
Returns the local basis, i.e., the set of shape functions.
Definition: crouzeixraviart.hh:198
static constexpr std::size_t size()
The number of shape functions.
Definition: crouzeixraviart.hh:218
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:114
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
constexpr GeometryType simplex(unsigned int dim)
Returns a GeometryType representing a simplex of dimension dim.
Definition: type.hh:453
void interpolate(const F &f, const GFS &gfs, XG &xg)
interpolation from a given grid function
Definition: interpolate.hh:177
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:279
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
static const ReferenceElement & simplex()
get simplex reference elements
Definition: referenceelements.hh:162
D DomainType
domain type
Definition: localbasis.hh:43
traits helper struct
Definition: localfiniteelementtraits.hh:13
LB LocalBasisType
Definition: localfiniteelementtraits.hh:16
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:20
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:24
A unique label for each type of element that can occur in a grid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)