Dune Core Modules (2.9.0)

raviartthomas03dlocalbasis.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_RAVIARTTHOMAS_RAVIARTTHOMAS03D_RAVIARTTHOMAS03DLOCALBASIS_HH
6 #define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS_RAVIARTTHOMAS03D_RAVIARTTHOMAS03DLOCALBASIS_HH
7 
8 #include <numeric>
9 
10 #include <dune/common/fmatrix.hh>
11 
12 #include <dune/localfunctions/common/localbasis.hh>
13 
14 namespace Dune
15 {
24  template<class D, class R>
26  {
27  public:
30 
32  RT03DLocalBasis (std::bitset<4> s = 0)
33  {
34  for (int i=0; i<4; i++)
35  sign_[i] = s[i] ? -1.0 : 1.0;
36  }
37 
39  unsigned int size () const
40  {
41  return 4;
42  }
43 
45  inline void evaluateFunction (const typename Traits::DomainType& in,
46  std::vector<typename Traits::RangeType>& out) const
47  {
48  out.resize(4);
49  auto c = std::sqrt(2.0);
50  out[0] = {sign_[0]*c* in[0], sign_[0]*c* in[1], sign_[0]*c*(in[2]-D(1))};
51  out[1] = {sign_[1]*c* in[0], sign_[1]*c*(in[1]-D(1)), sign_[1]*c* in[2] };
52  out[2] = {sign_[2]*c*(in[0]-D(1)), sign_[2]*c* in[1], sign_[2]*c* in[2] };
53  out[3] = {sign_[3]*c* in[0], sign_[3]*c* in[1], sign_[3]*c* in[2] };
54  }
55 
57  inline void
58  evaluateJacobian (const typename Traits::DomainType& in, // position
59  std::vector<typename Traits::JacobianType>& out) const // return value
60  {
61  out.resize(4);
62  for (int i=0; i<4; i++)
63  {
64  auto c = std::sqrt(2.0);
65  out[i][0] = {c*sign_[i], 0, 0};
66  out[i][1] = { 0,c*sign_[i], 0};
67  out[i][2] = { 0, 0,c*sign_[i]};
68  }
69  }
70 
72  void partial (const std::array<unsigned int, 3>& order,
73  const typename Traits::DomainType& in, // position
74  std::vector<typename Traits::RangeType>& out) const // return value
75  {
76  auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
77  if (totalOrder == 0) {
78  evaluateFunction(in, out);
79  } else if (totalOrder == 1) {
80  auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
81  out.resize(size());
82 
83  for (int i=0; i<size(); i++)
84  {
85  out[i][direction] = sign_[i]* std::sqrt(2.0) ;
86  out[i][(direction+1)%3] = 0;
87  out[i][(direction+2)%3] = 0;
88  }
89  } else {
90  out.resize(size());
91  for (std::size_t i = 0; i < size(); ++i)
92  for (std::size_t j = 0; j < 3; ++j)
93  out[i][j] = 0;
94  }
95 
96  }
97 
99  unsigned int order () const
100  {
101  return 1;
102  }
103 
104  private:
105 
106  // Signs of the face normals
107  std::array<R,4> sign_;
108  };
109 }
110 #endif
A dense n x m matrix.
Definition: fmatrix.hh:117
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Lowest order Raviart-Thomas shape functions on the reference tetrahedron.
Definition: raviartthomas03dlocalbasis.hh:26
RT03DLocalBasis(std::bitset< 4 > s=0)
Make set number s, where 0 <= s < 16.
Definition: raviartthomas03dlocalbasis.hh:32
void partial(const std::array< unsigned int, 3 > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition: raviartthomas03dlocalbasis.hh:72
unsigned int order() const
Polynomial order of the shape functions.
Definition: raviartthomas03dlocalbasis.hh:99
unsigned int size() const
number of shape functions
Definition: raviartthomas03dlocalbasis.hh:39
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: raviartthomas03dlocalbasis.hh:58
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas03dlocalbasis.hh:45
Implements a matrix constructed from a given type representing a field and compile-time given number ...
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:291
Dune namespace.
Definition: alignedallocator.hh:13
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:34
D DomainType
domain type
Definition: localbasis.hh:42
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 2, 22:35, 2024)