Dune Core Modules (2.6.0)

raviartthomas1cube2dlocalbasis.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_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE2D_LOCALBASIS_HH
4 #define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE2D_LOCALBASIS_HH
5 
6 #include <numeric>
7 #include <vector>
8 
9 #include <dune/common/fmatrix.hh>
10 
11 #include "../../common/localbasis.hh"
12 
13 namespace Dune
14 {
24  template<class D, class R>
26  {
27 
28  public:
31 
34  {
35  sign0 = sign1 = sign2 = sign3 = 1.0;
36  }
37 
43  RT1Cube2DLocalBasis (unsigned int s)
44  {
45  sign0 = sign1 = sign2 = sign3 = 1.0;
46  if (s & 1)
47  {
48  sign0 = -1.0;
49  }
50  if (s & 2)
51  {
52  sign1 = -1.0;
53  }
54  if (s & 4)
55  {
56  sign2 = -1.0;
57  }
58  if (s & 8)
59  {
60  sign3 = -1.0;
61  }
62  }
63 
65  unsigned int size () const
66  {
67  return 12;
68  }
69 
76  inline void evaluateFunction (const typename Traits::DomainType& in,
77  std::vector<typename Traits::RangeType>& out) const
78  {
79  out.resize(12);
80 
81  out[0][0] = sign0*(-1.0 + 4.0*in[0]-3*in[0]*in[0]);
82  out[0][1] = 0.0;
83  out[1][0] = 3.0 - 12.0*in[0] - 6.0*in[1] + 24.0*in[0]*in[1]+9*in[0]*in[0] - 18.0*in[0]*in[0]*in[1];
84  out[1][1] = 0.0;
85  out[2][0] = sign1*(-2.0*in[0] + 3.0*in[0]*in[0]);
86  out[2][1] = 0.0;
87  out[3][0] = -6.0*in[0] + 12.0*in[0]*in[1] + 9.0*in[0]*in[0] - 18.0*in[0]*in[0]*in[1];
88  out[3][1] = 0.0;
89  out[4][0] = 0.0;
90  out[4][1] = sign2*(-1.0 + 4.0*in[1] - 3.0*in[1]*in[1]);
91  out[5][0] = 0.0;
92  out[5][1] = -3.0 + 6.0*in[0] + 12.0*in[1] - 24.0*in[0]*in[1] - 9.0*in[1]*in[1] + 18.0*in[0]*in[1]*in[1];
93  out[6][0] = 0.0;
94  out[6][1] = sign3*(-2.0*in[1] + 3.0*in[1]*in[1]);
95  out[7][0] = 0.0;
96  out[7][1] = 6.0*in[1] - 12.0*in[0]*in[1] - 9.0*in[1]*in[1] + 18.0*in[0]*in[1]*in[1];
97  out[8][0] = 24.0*in[0] - 36.0*in[0]*in[1] - 24.0*in[0]*in[0] + 36.0*in[0]*in[0]*in[1];
98  out[8][1] = 0.0;
99  out[9][0] = 0.0;
100  out[9][1] = 24.0*in[1] - 36.0*in[0]*in[1] - 24.0*in[1]*in[1] + 36.0*in[0]*in[1]*in[1];
101  out[10][0] = -36.0*in[0] + 72.0*in[0]*in[1] + 36.0*in[0]*in[0] - 72.0*in[0]*in[0]*in[1];
102  out[10][1] = 0.0;
103  out[11][0] = 0.0;
104  out[11][1] = -36.0*in[1] + 72.0*in[0]*in[1] + 36*in[1]*in[1] - 72.0*in[0]*in[1]*in[1];
105  }
106 
113  inline void evaluateJacobian (const typename Traits::DomainType& in,
114  std::vector<typename Traits::JacobianType>& out) const
115  {
116  out.resize(12);
117 
118  out[0][0][0] = sign0*(4.0 - 6.0*in[0]);
119  out[0][0][1] = 0.0;
120  out[0][1][0] = 0.0;
121  out[0][1][1] = 0.0;
122 
123  out[1][0][0] = -12.0 + 24.0*in[1] + 18.0*in[0] - 36.0*in[0]*in[1];
124  out[1][0][1] = -6 + 24.0*in[0] - 18.0*in[0]*in[0];
125  out[1][1][0] = 0.0;
126  out[1][1][1] = 0.0;
127 
128  out[2][0][0] = sign1*(-2.0 + 6.0*in[0]);
129  out[2][0][1] = 0.0;
130  out[2][1][0] = 0.0;
131  out[2][1][1] = 0.0;
132 
133  out[3][0][0] = -6.0 + 12.0*in[1] + 18.0*in[0] - 36.0*in[0]*in[1];
134  out[3][0][1] = 12.0*in[0] - 18.0*in[0]*in[0];
135  out[3][1][0] = 0.0;
136  out[3][1][1] = 0.0;
137 
138  out[4][0][0] = 0.0;
139  out[4][0][1] = 0.0;
140  out[4][1][0] = 0.0;
141  out[4][1][1] = sign2*(4.0 - 6.0*in[1]);
142 
143  out[5][0][0] = 0.0;
144  out[5][0][1] = 0.0;
145  out[5][1][0] = 6.0 - 24.0*in[1] + 18.0*in[1]*in[1];
146  out[5][1][1] = 12.0 - 24.0*in[0] - 18.0*in[1] + 36.0*in[0]*in[1];
147 
148  out[6][0][0] = 0.0;
149  out[6][0][1] = 0.0;
150  out[6][1][0] = 0.0;
151  out[6][1][1] = sign3*(-2.0 + 6.0*in[1]);
152 
153  out[7][0][0] = 0.0;
154  out[7][0][1] = 0.0;
155  out[7][1][0] = -12.0*in[1] + 18.0*in[1]*in[1];
156  out[7][1][1] = 6.0 - 12.0*in[0] - 18.0*in[1] + 36.0*in[1]*in[0];
157 
158  out[8][0][0] = 24.0 - 36.0*in[1] - 48.0*in[0] + 72.0*in[0]*in[1];
159  out[8][0][1] = -36.0*in[0] + 36.0*in[0]*in[0];
160  out[8][1][0] = 0.0;
161  out[8][1][1] = 0.0;
162 
163  out[9][0][0] = 0.0;
164  out[9][0][1] = 0.0;
165  out[9][1][0] = -36.0*in[1] + 36.0*in[1]*in[1];
166  out[9][1][1] = 24.0 - 36.0*in[0] - 48.0*in[1] + 72.0*in[0]*in[1];
167 
168  out[10][0][0] = -36.0 + 72.0*in[1] + 72.0*in[0] - 144.0*in[0]*in[1];
169  out[10][0][1] = 72.0*in[0] - 72.0*in[0]*in[0];
170  out[10][1][0] = 0.0;
171  out[10][1][1] = 0.0;
172 
173  out[11][0][0] = 0.0;
174  out[11][0][1] = 0.0;
175  out[11][1][0] = 72.0*in[1] - 72.0*in[1]*in[1];
176  out[11][1][1] = -36.0 + 72.0*in[0] + 72.0*in[1] - 144.0*in[0]*in[1];
177  }
178 
180  void partial (const std::array<unsigned int, 2>& order,
181  const typename Traits::DomainType& in, // position
182  std::vector<typename Traits::RangeType>& out) const // return value
183  {
184  auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
185  if (totalOrder == 0) {
186  evaluateFunction(in, out);
187  } else {
188  DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
189  }
190  }
191 
193  unsigned int order () const
194  {
195  return 3;
196  }
197 
198  private:
199  R sign0, sign1, sign2, sign3;
200  };
201 }
202 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE2D_LOCALBASIS_HH
A dense n x m matrix.
Definition: fmatrix.hh:68
vector space out of a tensor product of fields.
Definition: fvector.hh:93
Default exception for dummy implementations.
Definition: exceptions.hh:261
First order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas1cube2dlocalbasis.hh:26
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: raviartthomas1cube2dlocalbasis.hh:113
unsigned int order() const
Polynomial order of the shape functions.
Definition: raviartthomas1cube2dlocalbasis.hh:193
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas1cube2dlocalbasis.hh:76
unsigned int size() const
number of shape functions
Definition: raviartthomas1cube2dlocalbasis.hh:65
RT1Cube2DLocalBasis()
Standard constructor.
Definition: raviartthomas1cube2dlocalbasis.hh:33
RT1Cube2DLocalBasis(unsigned int s)
Make set number s, where 0 <= s < 16.
Definition: raviartthomas1cube2dlocalbasis.hh:43
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: raviartthomas1cube2dlocalbasis.hh:180
Implements a matrix constructed from a given type representing a field and compile-time given number ...
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:331
Dune namespace.
Definition: alignedallocator.hh:10
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 7, 22:32, 2024)