Dune Core Modules (2.7.1)

raviartthomas02dlocalbasis.hh
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #ifndef DUNE_RT0TRIANGLELOCALBASIS_HH
4 #define DUNE_RT0TRIANGLELOCALBASIS_HH
5 
6 #include <numeric>
7 
8 #include <dune/common/fmatrix.hh>
9 
10 #include <dune/localfunctions/common/localbasis.hh>
11 
12 namespace Dune
13 {
22  template<class D, class R>
24  {
25  public:
28 
30  RT02DLocalBasis (std::bitset<3> s = 0)
31  {
32  for (int i=0; i<3; i++)
33  sign_[i] = s[i] ? -1.0 : 1.0;
34  }
35 
37  unsigned int size () const
38  {
39  return 3;
40  }
41 
43  inline void evaluateFunction (const typename Traits::DomainType& in,
44  std::vector<typename Traits::RangeType>& out) const
45  {
46  out.resize(3);
47  out[0] = {sign_[0]*in[0], sign_[0]*(in[1]-D(1))};
48  out[1] = {sign_[1]*(in[0]-D(1)), sign_[1]*in[1]};
49  out[2] = {sign_[2]*in[0], sign_[2]*in[1]};
50  }
51 
53  inline void
54  evaluateJacobian (const typename Traits::DomainType& in, // position
55  std::vector<typename Traits::JacobianType>& out) const // return value
56  {
57  out.resize(3);
58  for (int i=0; i<3; i++)
59  {
60  out[i][0] = {sign_[i], 0};
61  out[i][1] = { 0, sign_[i]};
62  }
63  }
64 
66  void partial (const std::array<unsigned int, 2>& order,
67  const typename Traits::DomainType& in, // position
68  std::vector<typename Traits::RangeType>& out) const // return value
69  {
70  auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
71  if (totalOrder == 0) {
72  evaluateFunction(in, out);
73  } else if (totalOrder == 1) {
74  auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
75  out.resize(size());
76 
77  for (int i=0; i<3; i++)
78  {
79  out[i][direction] = sign_[i];
80  out[i][1-direction] = 0;
81  }
82  } else {
83  out.resize(size());
84  for (std::size_t i = 0; i < size(); ++i)
85  for (std::size_t j = 0; j < 2; ++j)
86  out[i][j] = 0;
87  }
88 
89  }
90 
92  unsigned int order () const
93  {
94  return 1;
95  }
96 
97  private:
98 
99  // Signs of the edge normals
100  std::array<R,3> sign_;
101  };
102 }
103 #endif
A dense n x m matrix.
Definition: fmatrix.hh:69
vector space out of a tensor product of fields.
Definition: fvector.hh:96
Lowest order Raviart-Thomas shape functions on the reference triangle.
Definition: raviartthomas02dlocalbasis.hh:24
unsigned int order() const
Polynomial order of the shape functions.
Definition: raviartthomas02dlocalbasis.hh:92
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: raviartthomas02dlocalbasis.hh:54
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas02dlocalbasis.hh:43
unsigned int size() const
number of shape functions
Definition: raviartthomas02dlocalbasis.hh:37
RT02DLocalBasis(std::bitset< 3 > s=0)
Make set number s, where 0 <= s < 8.
Definition: raviartthomas02dlocalbasis.hh:30
void partial(const std::array< unsigned int, 2 > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition: raviartthomas02dlocalbasis.hh:66
Implements a matrix constructed from a given type representing a field and compile-time given number ...
T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:290
Dune namespace.
Definition: alignedallocator.hh:14
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:32
D DomainType
domain type
Definition: localbasis.hh:43
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 16, 22:29, 2024)