Dune Core Modules (2.9.0)

nedelec1stkindsimplex.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_NEDELEC_NEDELEC1STKINDSIMPLEX_HH
6 #define DUNE_LOCALFUNCTIONS_NEDELEC_NEDELEC1STKINDSIMPLEX_HH
7 
8 #include <numeric>
9 
10 #include <dune/common/fmatrix.hh>
11 #include <dune/common/fvector.hh>
12 
13 #include <dune/geometry/referenceelements.hh>
14 #include <dune/geometry/type.hh>
15 
16 #include <dune/localfunctions/common/localbasis.hh>
17 #include <dune/localfunctions/common/localfiniteelementtraits.hh>
18 #include <dune/localfunctions/common/localinterpolation.hh> // For deprecated makeFunctionWithCallOperator
19 #include <dune/localfunctions/common/localkey.hh>
20 
21 namespace Dune
22 {
23 namespace Impl
24 {
35  template<class D, class R, int dim, int k>
36  class Nedelec1stKindSimplexLocalBasis
37  {
38  // Number of edges of the reference simplex
39  constexpr static std::size_t numberOfEdges = dim*(dim+1)/2;
40 
41  public:
42  using Traits = LocalBasisTraits<D,dim,FieldVector<D,dim>,
43  R,dim,FieldVector<R,dim>,
44  FieldMatrix<R,dim,dim> >;
45 
52  Nedelec1stKindSimplexLocalBasis()
53  {
54  std::fill(edgeOrientation_.begin(), edgeOrientation_.end(), 1.0);
55  }
56 
59  Nedelec1stKindSimplexLocalBasis(std::bitset<numberOfEdges> edgeOrientation)
60  : Nedelec1stKindSimplexLocalBasis()
61  {
62  for (std::size_t i=0; i<edgeOrientation_.size(); i++)
63  edgeOrientation_[i] *= edgeOrientation[i] ? -1.0 : 1.0;
64  }
65 
67  static constexpr unsigned int size()
68  {
69  static_assert(dim==2 || dim==3, "Nedelec shape functions are implemented only for 2d and 3d simplices.");
70  if (dim==2)
71  return k * (k+2);
72  if (dim==3)
73  return k * (k+2) * (k+3) / 2;
74  }
75 
81  void evaluateFunction(const typename Traits::DomainType& in,
82  std::vector<typename Traits::RangeType>& out) const
83  {
84  static_assert(k==1, "Evaluating Nédélec shape functions is implemented only for first order.");
85  out.resize(size());
86 
87  if (dim==2)
88  {
89  // First-order Nédélec shape functions on a triangle are of the form
90  //
91  // (a1, a2) + b(-x2, x1)^T, a_1, a_2, b \in R
92  out[0] = {D(1) - in[1], in[0]};
93  out[1] = {in[1], -in[0]+D(1)};
94  out[2] = {-in[1], in[0]};
95  }
96 
97  if constexpr (dim==3)
98  {
99  // First-order Nédélec shape functions on a tetrahedron are of the form
100  //
101  // a + b \times x, a, b \in R^3
102  //
103  // The following coefficients create the six basis vectors
104  // that are dual to the edge degrees of freedom:
105  //
106  // a[0] = { 1, 0, 0} b[0] = { 0, -1, 1}
107  // a[1] = { 0, 1, 0} b[1] = { 1, 0, -1}
108  // a[2] = { 0, 0, 0} b[2] = { 0, 0, 1}
109  // a[3] = { 0, 0, 1} b[3] = {-1, 1, 0}
110  // a[4] = { 0, 0, 0} b[4] = { 0, -1, 0}
111  // a[5] = { 0, 0, 0} b[5] = { 1, 0, 0}
112  //
113  // The following implementation uses these values, and simply
114  // skips all the zeros.
115 
116  out[0] = { 1 - in[1] - in[2], in[0] , in[0] };
117  out[1] = { in[1] , 1 - in[0] - in[2], in[1]};
118  out[2] = { - in[1] , in[0] , 0 };
119  out[3] = { in[2], in[2], 1 - in[0] - in[1]};
120  out[4] = { -in[2], 0 , in[0] };
121  out[5] = { 0 , -in[2], in[1]};
122  }
123 
124  for (std::size_t i=0; i<out.size(); i++)
125  out[i] *= edgeOrientation_[i];
126  }
127 
133  void evaluateJacobian(const typename Traits::DomainType& in,
134  std::vector<typename Traits::JacobianType>& out) const
135  {
136  out.resize(size());
137  if (dim==2)
138  {
139  out[0][0] = { 0, -1};
140  out[0][1] = { 1, 0};
141 
142  out[1][0] = { 0, 1};
143  out[1][1] = {-1, 0};
144 
145  out[2][0] = { 0, -1};
146  out[2][1] = { 1, 0};
147  }
148  if (dim==3)
149  {
150  out[0][0] = { 0,-1,-1};
151  out[0][1] = { 1, 0, 0};
152  out[0][2] = { 1, 0, 0};
153 
154  out[1][0] = { 0, 1, 0};
155  out[1][1] = {-1, 0, -1};
156  out[1][2] = { 0, 1, 0};
157 
158  out[2][0] = { 0, -1, 0};
159  out[2][1] = { 1, 0, 0};
160  out[2][2] = { 0, 0, 0};
161 
162  out[3][0] = { 0, 0, 1};
163  out[3][1] = { 0, 0, 1};
164  out[3][2] = {-1, -1, 0};
165 
166  out[4][0] = { 0, 0, -1};
167  out[4][1] = { 0, 0, 0};
168  out[4][2] = { 1, 0, 0};
169 
170  out[5][0] = { 0, 0, 0};
171  out[5][1] = { 0, 0, -1};
172  out[5][2] = { 0, 1, 0};
173  }
174 
175  for (std::size_t i=0; i<out.size(); i++)
176  out[i] *= edgeOrientation_[i];
177 
178  }
179 
186  void partial(const std::array<unsigned int, dim>& order,
187  const typename Traits::DomainType& in,
188  std::vector<typename Traits::RangeType>& out) const
189  {
190  auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
191  if (totalOrder == 0) {
192  evaluateFunction(in, out);
193  } else if (totalOrder == 1) {
194  auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
195  out.resize(size());
196 
197  if (dim==2)
198  {
199  if (direction==0)
200  {
201  out[0] = {0, 1};
202  out[1] = {0, -1};
203  out[2] = {0, 1};
204  }
205  else
206  {
207  out[0] = {-1, 0};
208  out[1] = { 1, 0};
209  out[2] = {-1, 0};
210  }
211  }
212 
213  if (dim==3)
214  {
215  switch (direction)
216  {
217  case 0:
218  out[0] = { 0, 1, 1};
219  out[1] = { 0,-1, 0};
220  out[2] = { 0, 1, 0};
221  out[3] = { 0, 0,-1};
222  out[4] = { 0, 0, 1};
223  out[5] = { 0, 0, 0};
224  break;
225 
226  case 1:
227  out[0] = {-1, 0, 0};
228  out[1] = { 1, 0, 1};
229  out[2] = {-1, 0, 0};
230  out[3] = { 0, 0,-1};
231  out[4] = { 0, 0, 0};
232  out[5] = { 0, 0, 1};
233  break;
234 
235  case 2:
236  out[0] = {-1, 0, 0};
237  out[1] = { 0,-1, 0};
238  out[2] = { 0, 0, 0};
239  out[3] = { 1, 1, 0};
240  out[4] = {-1, 0, 0};
241  out[5] = { 0,-1, 0};
242  break;
243  }
244  }
245 
246  for (std::size_t i=0; i<out.size(); i++)
247  out[i] *= edgeOrientation_[i];
248 
249  } else {
250  out.resize(size());
251  for (std::size_t i = 0; i < size(); ++i)
252  for (std::size_t j = 0; j < dim; ++j)
253  out[i][j] = 0;
254  }
255 
256  }
257 
259  unsigned int order() const
260  {
261  return k;
262  }
263 
264  private:
265 
266  // Orientations of the simplex edges
267  std::array<R,numberOfEdges> edgeOrientation_;
268  };
269 
270 
275  template <int dim, int k>
276  class Nedelec1stKindSimplexLocalCoefficients
277  {
278  public:
280  Nedelec1stKindSimplexLocalCoefficients ()
281  : localKey_(size())
282  {
283  static_assert(k==1, "Only first-order Nédélec local coefficients are implemented.");
284  // Assign all degrees of freedom to edges
285  // TODO: This is correct only for first-order Nédélec elements
286  for (std::size_t i=0; i<size(); i++)
287  localKey_[i] = LocalKey(i,dim-1,0);
288  }
289 
291  std::size_t size() const
292  {
293  static_assert(dim==2 || dim==3, "Nédélec shape functions are implemented only for 2d and 3d simplices.");
294  return (dim==2) ? k * (k+2)
295  : k * (k+2) * (k+3) / 2;
296  }
297 
300  const LocalKey& localKey (std::size_t i) const
301  {
302  return localKey_[i];
303  }
304 
305  private:
306  std::vector<LocalKey> localKey_;
307  };
308 
313  template<class LB>
314  class Nedelec1stKindSimplexLocalInterpolation
315  {
316  static constexpr auto dim = LB::Traits::dimDomain;
317  static constexpr auto size = LB::size();
318 
319  // Number of edges of the reference simplex
320  constexpr static std::size_t numberOfEdges = dim*(dim+1)/2;
321 
322  public:
323 
325  Nedelec1stKindSimplexLocalInterpolation (std::bitset<numberOfEdges> s = 0)
326  {
327  auto refElement = Dune::referenceElement<double,dim>(GeometryTypes::simplex(dim));
328 
329  for (std::size_t i=0; i<numberOfEdges; i++)
330  m_[i] = refElement.position(i,dim-1);
331 
332  for (std::size_t i=0; i<numberOfEdges; i++)
333  {
334  auto vertexIterator = refElement.subEntities(i,dim-1,dim).begin();
335  auto v0 = *vertexIterator;
336  auto v1 = *(++vertexIterator);
337  // By default, edges point from the vertex with the smaller index
338  // to the vertex with the larger index.
339  if (v0>v1)
340  std::swap(v0,v1);
341  edge_[i] = refElement.position(v1,dim) - refElement.position(v0,dim);
342  edge_[i] *= (s[i]) ? -1.0 : 1.0;
343  }
344  }
345 
351  template<typename F, typename C>
352  void interpolate (const F& ff, std::vector<C>& out) const
353  {
354  out.resize(size);
355  auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
356 
357  for (std::size_t i=0; i<size; i++)
358  {
359  auto y = f(m_[i]);
360  out[i] = 0.0;
361  for (int j=0; j<dim; j++)
362  out[i] += y[j]*edge_[i][j];
363  }
364  }
365 
366  private:
367  // Edge midpoints of the reference simplex
368  std::array<typename LB::Traits::DomainType, numberOfEdges> m_;
369  // Edges of the reference simplex
370  std::array<typename LB::Traits::DomainType, numberOfEdges> edge_;
371  };
372 
373 }
374 
375 
401  template<class D, class R, int dim, int k>
403  {
404  public:
406  Impl::Nedelec1stKindSimplexLocalCoefficients<dim,k>,
407  Impl::Nedelec1stKindSimplexLocalInterpolation<Impl::Nedelec1stKindSimplexLocalBasis<D,R,dim,k> > >;
408 
409  static_assert(dim==2 || dim==3, "Nedelec elements are only implemented for 2d and 3d elements.");
410  static_assert(k==1, "Nedelec elements of the first kind are currently only implemented for order k==1.");
411 
415 
421  Nedelec1stKindSimplexLocalFiniteElement (std::bitset<dim*(dim+1)/2> s) :
422  basis_(s),
423  interpolation_(s)
424  {}
425 
426  const typename Traits::LocalBasisType& localBasis () const
427  {
428  return basis_;
429  }
430 
431  const typename Traits::LocalCoefficientsType& localCoefficients () const
432  {
433  return coefficients_;
434  }
435 
436  const typename Traits::LocalInterpolationType& localInterpolation () const
437  {
438  return interpolation_;
439  }
440 
441  static constexpr unsigned int size ()
442  {
443  return Traits::LocalBasisType::size();
444  }
445 
446  static constexpr GeometryType type ()
447  {
448  return GeometryTypes::simplex(dim);
449  }
450 
451  private:
452  typename Traits::LocalBasisType basis_;
453  typename Traits::LocalCoefficientsType coefficients_;
454  typename Traits::LocalInterpolationType interpolation_;
455  };
456 
457 }
458 
459 #endif
Nédélec elements of the first kind for simplex elements.
Definition: nedelec1stkindsimplex.hh:403
Nedelec1stKindSimplexLocalFiniteElement()=default
Default constructor.
Nedelec1stKindSimplexLocalFiniteElement(std::bitset< dim *(dim+1)/2 > s)
Constructor with explicitly given edge orientations.
Definition: nedelec1stkindsimplex.hh:421
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:132
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:463
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:291
Dune namespace.
Definition: alignedallocator.hh:13
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.80.0 (Apr 26, 22:29, 2024)