Dune Core Modules (2.9.1)

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 (C) 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/localinterpolation.hh>
20#include <dune/localfunctions/common/localkey.hh>
21
22namespace Dune { namespace Impl
23{
30 template<class D, class R, unsigned int dim>
31 class CrouzeixRaviartLocalBasis
32 {
33 public:
34 using Traits = LocalBasisTraits<D,dim,FieldVector<D,dim>,R,1,FieldVector<R,1>,FieldMatrix<R,1,dim> >;
35
38 static constexpr unsigned int size ()
39 {
40 return dim+1;
41 }
42
44 void evaluateFunction(const typename Traits::DomainType& x,
45 std::vector<typename Traits::RangeType>& out) const
46 {
47 out.resize(size());
48
49 std::fill(out.begin(), out.end()-1, 1.0);
50 out.back() = 1.0-dim;
51
52 for (unsigned int i=0; i<dim; i++)
53 {
54 out[i] -= dim * x[dim-i-1];
55 out.back() += dim*x[i];
56 }
57 }
58
64 void evaluateJacobian(const typename Traits::DomainType& x,
65 std::vector<typename Traits::JacobianType>& out) const
66 {
67 out.resize(size());
68
69 for (unsigned i=0; i<dim; i++)
70 for (unsigned j=0; j<dim; j++)
71 out[i][0][j] = (i==(dim-1-j)) ? -(double)dim : 0;
72
73 std::fill(out.back()[0].begin(), out.back()[0].end(), dim);
74 }
75
82 void partial(const std::array<unsigned int,dim>& order,
83 const typename Traits::DomainType& in,
84 std::vector<typename Traits::RangeType>& out) const
85 {
86 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
87
88 out.resize(size());
89
90 if (totalOrder == 0) {
91 evaluateFunction(in, out);
92 return;
93 }
94
95 if (totalOrder==1)
96 {
97 auto direction = std::find(order.begin(), order.end(), 1)-order.begin();
98
99 for (unsigned int i=0; i<dim; i++)
100 out[i] = (i==(dim-1-direction)) ? -(double)dim : 0.0;
101
102 out.back()[0] = dim;
103 }
104 else // all higher order derivatives are zero
105 std::fill(out.begin(), out.end(), 0);
106 }
107
109 static constexpr unsigned int order ()
110 {
111 return 1;
112 }
113 };
114
119 template<unsigned int dim>
120 class CrouzeixRaviartLocalCoefficients
121 {
122 public:
124 CrouzeixRaviartLocalCoefficients ()
125 : localKeys_(size())
126 {
127 for (std::size_t i=0; i<size(); i++)
128 localKeys_[i] = LocalKey(i,1,0);
129 }
130
132 static constexpr std::size_t size ()
133 {
134 return dim+1;
135 }
136
138 const LocalKey& localKey (std::size_t i) const
139 {
140 return localKeys_[i];
141 }
142
143 private:
144 std::vector<LocalKey> localKeys_;
145 };
146
151 template<class LocalBasis>
152 class CrouzeixRaviartLocalInterpolation
153 {
154 public:
155
163 template<typename F, typename C>
164 void interpolate (const F& ff, std::vector<C>& out) const
165 {
166 constexpr auto dim = LocalBasis::Traits::dimDomain;
167
168 auto&& f = Impl::makeFunctionWithCallOperator<typename LocalBasis::Traits::DomainType>(ff);
169
170 out.resize(LocalBasis::size());
171
172 // Evaluate at the facet midpoints
174 for (int i=0; i<refSimplex.size(1); i++)
175 out[i] = f(refSimplex.template geometry<1>(i).center());
176 }
177 };
178
179} } // namespace Dune::Impl
180
181namespace Dune
182{
189 template<class D, class R, int dim>
191 {
192 public:
196 Impl::CrouzeixRaviartLocalCoefficients<dim>,
197 Impl::CrouzeixRaviartLocalInterpolation<Impl::CrouzeixRaviartLocalBasis<D,R,dim> > >;
198
201 const typename Traits::LocalBasisType& localBasis() const
202 {
203 return basis_;
204 }
205
209 {
210 return coefficients_;
211 }
212
216 {
217 return interpolation_;
218 }
219
221 static constexpr std::size_t size()
222 {
223 return dim+1;
224 }
225
228 static constexpr GeometryType type()
229 {
230 return GeometryTypes::simplex(dim);
231 }
232
233 private:
234 Impl::CrouzeixRaviartLocalBasis<D,R,dim> basis_;
235 Impl::CrouzeixRaviartLocalCoefficients<dim> coefficients_;
236 Impl::CrouzeixRaviartLocalInterpolation<Impl::CrouzeixRaviartLocalBasis<D,R,dim> > interpolation_;
237 };
238
239} // namespace Dune
240
241#endif // DUNE_LOCALFUNCTIONS_CROUZEIXRAVIART_HH
Crouzeix-Raviart finite element.
Definition: crouzeixraviart.hh:191
static constexpr GeometryType type()
The reference element that the local finite element is defined on.
Definition: crouzeixraviart.hh:228
const Traits::LocalCoefficientsType & localCoefficients() const
Returns the assignment of the degrees of freedom to the element subentities.
Definition: crouzeixraviart.hh:208
const Traits::LocalInterpolationType & localInterpolation() const
Returns object that evaluates degrees of freedom.
Definition: crouzeixraviart.hh:215
const Traits::LocalBasisType & localBasis() const
Returns the local basis, i.e., the set of shape functions.
Definition: crouzeixraviart.hh:201
static constexpr std::size_t size()
The number of shape functions.
Definition: crouzeixraviart.hh:221
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:126
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:464
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:291
Dune namespace.
Definition: alignedallocator.hh:13
static const ReferenceElement & simplex()
get simplex reference elements
Definition: referenceelements.hh:204
D DomainType
domain type
Definition: localbasis.hh:42
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 (Nov 21, 23:30, 2024)