dune-localfunctions  2.5.0
raviartthomas02dlocalbasis.hh
Go to the documentation of this file.
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 
11 
12 namespace Dune
13 {
22  template<class D, class R>
24  {
25  public:
26  typedef LocalBasisTraits<D,2,Dune::FieldVector<D,2>,R,2,Dune::FieldVector<R,2>,
27  Dune::FieldMatrix<R,2,2> > Traits;
28 
31  {
32  sign0 = sign1 = sign2 = 1.0;
33  }
34 
36  RT02DLocalBasis (unsigned int s)
37  {
38  sign0 = sign1 = sign2 = 1.0;
39  if (s&1) sign0 = -1.0;
40  if (s&2) sign1 = -1.0;
41  if (s&4) sign2 = -1.0;
42  }
43 
45  unsigned int size () const
46  {
47  return 3;
48  }
49 
51  inline void evaluateFunction (const typename Traits::DomainType& in,
52  std::vector<typename Traits::RangeType>& out) const
53  {
54  out.resize(3);
55  out[0][0] = sign0*in[0]; out[0][1]=sign0*(in[1]-D(1));
56  out[1][0] = sign1*(in[0]-D(1)); out[1][1]=sign1*in[1];
57  out[2][0] = sign2*in[0]; out[2][1]=sign2*in[1];
58  }
59 
61  inline void
62  evaluateJacobian (const typename Traits::DomainType& in, // position
63  std::vector<typename Traits::JacobianType>& out) const // return value
64  {
65  out.resize(3);
66  out[0][0][0] = sign0; out[0][0][1] = 0;
67  out[0][1][0] = 0; out[0][1][1] = sign0;
68  out[1][0][0] = sign1; out[1][0][1] = 0;
69  out[1][1][0] = 0; out[1][1][1] = sign1;
70  out[2][0][0] = sign2; out[2][0][1] = 0;
71  out[2][1][0] = 0; out[2][1][1] = sign2;
72  }
73 
75  void partial (const std::array<unsigned int, 2>& order,
76  const typename Traits::DomainType& in, // position
77  std::vector<typename Traits::RangeType>& out) const // return value
78  {
79  auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
80  if (totalOrder == 0) {
81  evaluateFunction(in, out);
82  } else if (totalOrder == 1) {
83  auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
84  out.resize(size());
85 
86  out[0][direction] = sign0;
87  out[0][1-direction] = 0;
88  out[1][direction] = sign1;
89  out[1][1-direction] = 0;
90  out[2][direction] = sign2;
91  out[2][1-direction] = 0;
92  } else {
93  out.resize(size());
94  for (std::size_t i = 0; i < size(); ++i)
95  for (std::size_t j = 0; j < 2; ++j)
96  out[i][j] = 0;
97  }
98 
99  }
100 
102  unsigned int order () const
103  {
104  return 1;
105  }
106 
107  private:
108  R sign0, sign1, sign2;
109  };
110 }
111 #endif
unsigned int order() const
Polynomial order of the shape functions.
Definition: raviartthomas02dlocalbasis.hh:102
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:15
D DomainType
domain type
Definition: localbasis.hh:49
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: raviartthomas02dlocalbasis.hh:62
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas02dlocalbasis.hh:51
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition: raviartthomas02dlocalbasis.hh:27
unsigned int size() const
number of shape functions
Definition: raviartthomas02dlocalbasis.hh:45
Lowest order Raviart-Thomas shape functions on the reference triangle.
Definition: raviartthomas02dlocalbasis.hh:23
RT02DLocalBasis()
Standard constructor.
Definition: raviartthomas02dlocalbasis.hh:30
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
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:75
RT02DLocalBasis(unsigned int s)
Make set numer s, where 0<=s<8.
Definition: raviartthomas02dlocalbasis.hh:36