Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
polynomials_rt_bubbles.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2018 - 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 
18 
22 
23 #include <iomanip>
24 #include <iostream>
25 
27 
28 
29 template <int dim>
31  : TensorPolynomialsBase<dim>(k, n_polynomials(k))
32  , raviart_thomas_space(k - 1)
33  , monomials(k + 2)
34 {
35  Assert(dim >= 2, ExcImpossibleInDim(dim));
36 
37  for (unsigned int i = 0; i < monomials.size(); ++i)
39 }
40 
41 
42 
43 template <int dim>
44 void
46  const Point<dim> & unit_point,
47  std::vector<Tensor<1, dim>> &values,
48  std::vector<Tensor<2, dim>> &grads,
49  std::vector<Tensor<3, dim>> &grad_grads,
50  std::vector<Tensor<4, dim>> &third_derivatives,
51  std::vector<Tensor<5, dim>> &fourth_derivatives) const
52 {
53  Assert(values.size() == this->n() || values.size() == 0,
54  ExcDimensionMismatch(values.size(), this->n()));
55  Assert(grads.size() == this->n() || grads.size() == 0,
56  ExcDimensionMismatch(grads.size(), this->n()));
57  Assert(grad_grads.size() == this->n() || grad_grads.size() == 0,
58  ExcDimensionMismatch(grad_grads.size(), this->n()));
59  Assert(third_derivatives.size() == this->n() || third_derivatives.size() == 0,
60  ExcDimensionMismatch(third_derivatives.size(), this->n()));
61  Assert(fourth_derivatives.size() == this->n() ||
62  fourth_derivatives.size() == 0,
63  ExcDimensionMismatch(fourth_derivatives.size(), this->n()));
64 
65  // Third and fourth derivatives are not implemented
66  (void)third_derivatives;
67  Assert(third_derivatives.size() == 0, ExcNotImplemented());
68  (void)fourth_derivatives;
69  Assert(fourth_derivatives.size() == 0, ExcNotImplemented());
70 
71  const unsigned int n_sub = raviart_thomas_space.n();
72  const unsigned int my_degree = this->degree();
73 
74  // Guard access to the scratch arrays in the following block
75  // using a mutex to make sure they are not used by multiple threads
76  // at once
77  {
78  static std::mutex mutex;
79  std::lock_guard<std::mutex> lock(mutex);
80 
81  static std::vector<Tensor<1, dim>> p_values;
82  static std::vector<Tensor<2, dim>> p_grads;
83  static std::vector<Tensor<3, dim>> p_grad_grads;
84  static std::vector<Tensor<4, dim>> p_third_derivatives;
85  static std::vector<Tensor<5, dim>> p_fourth_derivatives;
86 
87  p_values.resize((values.size() == 0) ? 0 : n_sub);
88  p_grads.resize((grads.size() == 0) ? 0 : n_sub);
89  p_grad_grads.resize((grad_grads.size() == 0) ? 0 : n_sub);
90 
91  // This is the Raviart-Thomas part of the space
92  raviart_thomas_space.evaluate(unit_point,
93  p_values,
94  p_grads,
95  p_grad_grads,
96  p_third_derivatives,
97  p_fourth_derivatives);
98  for (unsigned int i = 0; i < p_values.size(); ++i)
99  values[i] = p_values[i];
100  for (unsigned int i = 0; i < p_grads.size(); ++i)
101  grads[i] = p_grads[i];
102  for (unsigned int i = 0; i < p_grad_grads.size(); ++i)
103  grad_grads[i] = p_grad_grads[i];
104  }
105 
106  // Next we compute the polynomials and derivatives
107  // of the curl part of the space
108  const unsigned int n_derivatives = 3;
109  double monoval_plus[dim][n_derivatives + 1];
110  double monoval[dim][n_derivatives + 1];
111 
112  double monoval_i[dim][n_derivatives + 1];
113  double monoval_j[dim][n_derivatives + 1];
114  double monoval_jplus[dim][n_derivatives + 1];
115 
116  unsigned int start = n_sub;
117 
118  if (dim == 2)
119  {
120  // In 2d the curl part of the space is spanned by the vectors
121  // of two types. The first one is
122  // [ x^i * [y^(k+1)]' ]
123  // [ -[x^i]' * y^(k+1) ]
124  // The second one can be obtained from the first by a cyclic
125  // rotation of the coordinates.
126  // monoval_i = x^i,
127  // monoval_plus = x^(k+1)
128  for (unsigned int d = 0; d < dim; ++d)
129  monomials[my_degree + 1].value(unit_point(d),
130  n_derivatives,
131  monoval_plus[d]);
132 
133  for (unsigned int i = 0; i <= my_degree; ++i, ++start)
134  {
135  for (unsigned int d = 0; d < dim; ++d)
136  monomials[i].value(unit_point(d), n_derivatives, monoval_i[d]);
137 
138  if (values.size() != 0)
139  {
140  values[start][0] = monoval_i[0][0] * monoval_plus[1][1];
141  values[start][1] = -monoval_i[0][1] * monoval_plus[1][0];
142 
143  values[start + my_degree + 1][0] =
144  -monoval_plus[0][0] * monoval_i[1][1];
145  values[start + my_degree + 1][1] =
146  monoval_plus[0][1] * monoval_i[1][0];
147  }
148 
149  if (grads.size() != 0)
150  {
151  grads[start][0][0] = monoval_i[0][1] * monoval_plus[1][1];
152  grads[start][0][1] = monoval_i[0][0] * monoval_plus[1][2];
153  grads[start][1][0] = -monoval_i[0][2] * monoval_plus[1][0];
154  grads[start][1][1] = -monoval_i[0][1] * monoval_plus[1][1];
155 
156  grads[start + my_degree + 1][0][0] =
157  -monoval_plus[0][1] * monoval_i[1][1];
158  grads[start + my_degree + 1][0][1] =
159  -monoval_plus[0][0] * monoval_i[1][2];
160  grads[start + my_degree + 1][1][0] =
161  monoval_plus[0][2] * monoval_i[1][0];
162  grads[start + my_degree + 1][1][1] =
163  monoval_plus[0][1] * monoval_i[1][1];
164  }
165 
166  if (grad_grads.size() != 0)
167  {
168  grad_grads[start][0][0][0] = monoval_i[0][2] * monoval_plus[1][1];
169  grad_grads[start][0][0][1] = monoval_i[0][1] * monoval_plus[1][2];
170  grad_grads[start][0][1][0] = monoval_i[0][1] * monoval_plus[1][2];
171  grad_grads[start][0][1][1] = monoval_i[0][0] * monoval_plus[1][3];
172  grad_grads[start][1][0][0] =
173  -monoval_i[0][3] * monoval_plus[1][0];
174  grad_grads[start][1][0][1] =
175  -monoval_i[0][2] * monoval_plus[1][1];
176  grad_grads[start][1][1][0] =
177  -monoval_i[0][2] * monoval_plus[1][1];
178  grad_grads[start][1][1][1] =
179  -monoval_i[0][1] * monoval_plus[1][2];
180 
181  grad_grads[start + my_degree + 1][0][0][0] =
182  -monoval_plus[0][2] * monoval_i[1][1];
183  grad_grads[start + my_degree + 1][0][0][1] =
184  -monoval_plus[0][1] * monoval_i[1][2];
185  grad_grads[start + my_degree + 1][0][1][0] =
186  -monoval_plus[0][1] * monoval_i[1][2];
187  grad_grads[start + my_degree + 1][0][1][1] =
188  -monoval_plus[0][0] * monoval_i[1][3];
189  grad_grads[start + my_degree + 1][1][0][0] =
190  monoval_plus[0][3] * monoval_i[1][0];
191  grad_grads[start + my_degree + 1][1][0][1] =
192  monoval_plus[0][2] * monoval_i[1][1];
193  grad_grads[start + my_degree + 1][1][1][0] =
194  monoval_plus[0][2] * monoval_i[1][1];
195  grad_grads[start + my_degree + 1][1][1][1] =
196  monoval_plus[0][1] * monoval_i[1][2];
197  }
198  }
199  Assert(start == this->n() - my_degree - 1, ExcInternalError());
200  }
201  else if (dim == 3)
202  {
203  // In 3d the first type of basis vector is
204  // [ x^i * y^j * z^k * (j+k+2) ]
205  // [ -[x^i]' * y^(j+1) * z^k ]
206  // [ -[x^i]' * y^j * z^(k+1) ],
207  // For the second type of basis vector y and z
208  // are swapped. Then for each of these,
209  // two more are obtained by the cyclic rotation
210  // of the coordinates.
211  // monoval = x^k, monoval_plus = x^(k+1)
212  // monoval_* = x^*, monoval_jplus = x^(j+1)
213  for (unsigned int d = 0; d < dim; ++d)
214  {
215  monomials[my_degree + 1].value(unit_point(d),
216  n_derivatives,
217  monoval_plus[d]);
218  monomials[my_degree].value(unit_point(d), n_derivatives, monoval[d]);
219  }
220 
221  const unsigned int n_curls = (my_degree + 1) * (2 * my_degree + 1);
222  // Span of @f$\tilde{B}@f$
223  for (unsigned int i = 0; i <= my_degree; ++i)
224  {
225  for (unsigned int d = 0; d < dim; ++d)
226  monomials[i].value(unit_point(d), n_derivatives, monoval_i[d]);
227 
228  for (unsigned int j = 0; j <= my_degree; ++j)
229  {
230  for (unsigned int d = 0; d < dim; ++d)
231  {
232  monomials[j].value(unit_point(d),
233  n_derivatives,
234  monoval_j[d]);
235  monomials[j + 1].value(unit_point(d),
236  n_derivatives,
237  monoval_jplus[d]);
238  }
239 
240  if (values.size() != 0)
241  {
242  values[start][0] = monoval_i[0][0] * monoval_j[1][0] *
243  monoval[2][0] *
244  static_cast<double>(j + my_degree + 2);
245  values[start][1] =
246  -monoval_i[0][1] * monoval_jplus[1][0] * monoval[2][0];
247  values[start][2] =
248  -monoval_i[0][1] * monoval_j[1][0] * monoval_plus[2][0];
249 
250  values[start + n_curls][0] =
251  -monoval_jplus[0][0] * monoval_i[1][1] * monoval[2][0];
252  values[start + n_curls][1] =
253  monoval_j[0][0] * monoval_i[1][0] * monoval[2][0] *
254  static_cast<double>(j + my_degree + 2);
255  values[start + n_curls][2] =
256  -monoval_j[0][0] * monoval_i[1][1] * monoval_plus[2][0];
257 
258  values[start + 2 * n_curls][0] =
259  -monoval_jplus[0][0] * monoval[1][0] * monoval_i[2][1];
260  values[start + 2 * n_curls][1] =
261  -monoval_j[0][0] * monoval_plus[1][0] * monoval_i[2][1];
262  values[start + 2 * n_curls][2] =
263  monoval_j[0][0] * monoval[1][0] * monoval_i[2][0] *
264  static_cast<double>(j + my_degree + 2);
265 
266  // Only unique triples of powers (i j k)
267  // and (i k j) are allowed, 0 <= i,j <= k
268  if (j != my_degree)
269  {
270  values[start + 1][0] =
271  monoval_i[0][0] * monoval[1][0] * monoval_j[2][0] *
272  static_cast<double>(j + my_degree + 2);
273  values[start + 1][1] =
274  -monoval_i[0][1] * monoval_plus[1][0] * monoval_j[2][0];
275  values[start + 1][2] =
276  -monoval_i[0][1] * monoval[1][0] * monoval_jplus[2][0];
277 
278  values[start + n_curls + 1][0] =
279  -monoval_plus[0][0] * monoval_i[1][1] * monoval_j[2][0];
280  values[start + n_curls + 1][1] =
281  monoval[0][0] * monoval_i[1][0] * monoval_j[2][0] *
282  static_cast<double>(j + my_degree + 2);
283  values[start + n_curls + 1][2] =
284  -monoval[0][0] * monoval_i[1][1] * monoval_jplus[2][0];
285 
286  values[start + 2 * n_curls + 1][0] =
287  -monoval_plus[0][0] * monoval_j[1][0] * monoval_i[2][1];
288  values[start + 2 * n_curls + 1][1] =
289  -monoval[0][0] * monoval_jplus[1][0] * monoval_i[2][1];
290  values[start + 2 * n_curls + 1][2] =
291  monoval[0][0] * monoval_j[1][0] * monoval_i[2][0] *
292  static_cast<double>(j + my_degree + 2);
293  }
294  }
295 
296  if (grads.size() != 0)
297  {
298  grads[start][0][0] = monoval_i[0][1] * monoval_j[1][0] *
299  monoval[2][0] *
300  static_cast<double>(j + my_degree + 2);
301  grads[start][0][1] = monoval_i[0][0] * monoval_j[1][1] *
302  monoval[2][0] *
303  static_cast<double>(j + my_degree + 2);
304  grads[start][0][2] = monoval_i[0][0] * monoval_j[1][0] *
305  monoval[2][1] *
306  static_cast<double>(j + my_degree + 2);
307  grads[start][1][0] =
308  -monoval_i[0][2] * monoval_jplus[1][0] * monoval[2][0];
309  grads[start][1][1] =
310  -monoval_i[0][1] * monoval_jplus[1][1] * monoval[2][0];
311  grads[start][1][2] =
312  -monoval_i[0][1] * monoval_jplus[1][0] * monoval[2][1];
313  grads[start][2][0] =
314  -monoval_i[0][2] * monoval_j[1][0] * monoval_plus[2][0];
315  grads[start][2][1] =
316  -monoval_i[0][1] * monoval_j[1][1] * monoval_plus[2][0];
317  grads[start][2][2] =
318  -monoval_i[0][1] * monoval_j[1][0] * monoval_plus[2][1];
319 
320  grads[start + n_curls][0][0] =
321  -monoval_jplus[0][1] * monoval_i[1][1] * monoval[2][0];
322  grads[start + n_curls][0][1] =
323  -monoval_jplus[0][0] * monoval_i[1][2] * monoval[2][0];
324  grads[start + n_curls][0][2] =
325  -monoval_jplus[0][0] * monoval_i[1][1] * monoval[2][1];
326  grads[start + n_curls][1][0] =
327  monoval_j[0][1] * monoval_i[1][0] * monoval[2][0] *
328  static_cast<double>(j + my_degree + 2);
329  grads[start + n_curls][1][1] =
330  monoval_j[0][0] * monoval_i[1][1] * monoval[2][0] *
331  static_cast<double>(j + my_degree + 2);
332  grads[start + n_curls][1][2] =
333  monoval_j[0][0] * monoval_i[1][0] * monoval[2][1] *
334  static_cast<double>(j + my_degree + 2);
335  grads[start + n_curls][2][0] =
336  -monoval_j[0][1] * monoval_i[1][1] * monoval_plus[2][0];
337  grads[start + n_curls][2][1] =
338  -monoval_j[0][0] * monoval_i[1][2] * monoval_plus[2][0];
339  grads[start + n_curls][2][2] =
340  -monoval_j[0][0] * monoval_i[1][1] * monoval_plus[2][1];
341 
342  grads[start + 2 * n_curls][0][0] =
343  -monoval_jplus[0][1] * monoval[1][0] * monoval_i[2][1];
344  grads[start + 2 * n_curls][0][1] =
345  -monoval_jplus[0][0] * monoval[1][1] * monoval_i[2][1];
346  grads[start + 2 * n_curls][0][2] =
347  -monoval_jplus[0][0] * monoval[1][0] * monoval_i[2][2];
348  grads[start + 2 * n_curls][1][0] =
349  -monoval_j[0][1] * monoval_plus[1][0] * monoval_i[2][1];
350  grads[start + 2 * n_curls][1][1] =
351  -monoval_j[0][0] * monoval_plus[1][1] * monoval_i[2][1];
352  grads[start + 2 * n_curls][1][2] =
353  -monoval_j[0][0] * monoval_plus[1][0] * monoval_i[2][2];
354  grads[start + 2 * n_curls][2][0] =
355  monoval_j[0][1] * monoval[1][0] * monoval_i[2][0] *
356  static_cast<double>(j + my_degree + 2);
357  grads[start + 2 * n_curls][2][1] =
358  monoval_j[0][0] * monoval[1][1] * monoval_i[2][0] *
359  static_cast<double>(j + my_degree + 2);
360  grads[start + 2 * n_curls][2][2] =
361  monoval_j[0][0] * monoval[1][0] * monoval_i[2][1] *
362  static_cast<double>(j + my_degree + 2);
363 
364  if (j != my_degree)
365  {
366  grads[start + 1][0][0] =
367  monoval_i[0][1] * monoval[1][0] * monoval_j[2][0] *
368  static_cast<double>(j + my_degree + 2);
369  grads[start + 1][0][1] =
370  monoval_i[0][0] * monoval[1][1] * monoval_j[2][0] *
371  static_cast<double>(j + my_degree + 2);
372  grads[start + 1][0][2] =
373  monoval_i[0][0] * monoval[1][0] * monoval_j[2][1] *
374  static_cast<double>(j + my_degree + 2);
375  grads[start + 1][1][0] =
376  -monoval_i[0][2] * monoval_plus[1][0] * monoval_j[2][0];
377  grads[start + 1][1][1] =
378  -monoval_i[0][1] * monoval_plus[1][1] * monoval_j[2][0];
379  grads[start + 1][1][2] =
380  -monoval_i[0][1] * monoval_plus[1][0] * monoval_j[2][1];
381  grads[start + 1][2][0] =
382  -monoval_i[0][2] * monoval[1][0] * monoval_jplus[2][0];
383  grads[start + 1][2][1] =
384  -monoval_i[0][1] * monoval[1][1] * monoval_jplus[2][0];
385  grads[start + 1][2][2] =
386  -monoval_i[0][1] * monoval[1][0] * monoval_jplus[2][1];
387 
388  grads[start + n_curls + 1][0][0] =
389  -monoval_plus[0][1] * monoval_i[1][1] * monoval_j[2][0];
390  grads[start + n_curls + 1][0][1] =
391  -monoval_plus[0][0] * monoval_i[1][2] * monoval_j[2][0];
392  grads[start + n_curls + 1][0][2] =
393  -monoval_plus[0][0] * monoval_i[1][1] * monoval_j[2][1];
394  grads[start + n_curls + 1][1][0] =
395  monoval[0][1] * monoval_i[1][0] * monoval_j[2][0] *
396  static_cast<double>(j + my_degree + 2);
397  grads[start + n_curls + 1][1][1] =
398  monoval[0][0] * monoval_i[1][1] * monoval_j[2][0] *
399  static_cast<double>(j + my_degree + 2);
400  grads[start + n_curls + 1][1][2] =
401  monoval[0][0] * monoval_i[1][0] * monoval_j[2][1] *
402  static_cast<double>(j + my_degree + 2);
403  grads[start + n_curls + 1][2][0] =
404  -monoval[0][1] * monoval_i[1][1] * monoval_jplus[2][0];
405  grads[start + n_curls + 1][2][1] =
406  -monoval[0][0] * monoval_i[1][2] * monoval_jplus[2][0];
407  grads[start + n_curls + 1][2][2] =
408  -monoval[0][0] * monoval_i[1][1] * monoval_jplus[2][1];
409 
410  grads[start + 2 * n_curls + 1][0][0] =
411  -monoval_plus[0][1] * monoval_j[1][0] * monoval_i[2][1];
412  grads[start + 2 * n_curls + 1][0][1] =
413  -monoval_plus[0][0] * monoval_j[1][1] * monoval_i[2][1];
414  grads[start + 2 * n_curls + 1][0][2] =
415  -monoval_plus[0][0] * monoval_j[1][0] * monoval_i[2][2];
416  grads[start + 2 * n_curls + 1][1][0] =
417  -monoval[0][1] * monoval_jplus[1][0] * monoval_i[2][1];
418  grads[start + 2 * n_curls + 1][1][1] =
419  -monoval[0][0] * monoval_jplus[1][1] * monoval_i[2][1];
420  grads[start + 2 * n_curls + 1][1][2] =
421  -monoval[0][0] * monoval_jplus[1][0] * monoval_i[2][2];
422  grads[start + 2 * n_curls + 1][2][0] =
423  monoval[0][1] * monoval_j[1][0] * monoval_i[2][0] *
424  static_cast<double>(j + my_degree + 2);
425  grads[start + 2 * n_curls + 1][2][1] =
426  monoval[0][0] * monoval_j[1][1] * monoval_i[2][0] *
427  static_cast<double>(j + my_degree + 2);
428  grads[start + 2 * n_curls + 1][2][2] =
429  monoval[0][0] * monoval_j[1][0] * monoval_i[2][1] *
430  static_cast<double>(j + my_degree + 2);
431  }
432  }
433 
434  if (grad_grads.size() != 0)
435  {
436  grad_grads[start][0][0][0] =
437  monoval_i[0][2] * monoval_j[1][0] * monoval[2][0] *
438  static_cast<double>(j + my_degree + 2);
439  grad_grads[start][0][0][1] =
440  monoval_i[0][1] * monoval_j[1][1] * monoval[2][0] *
441  static_cast<double>(j + my_degree + 2);
442  grad_grads[start][0][0][2] =
443  monoval_i[0][1] * monoval_j[1][0] * monoval[2][1] *
444  static_cast<double>(j + my_degree + 2);
445  grad_grads[start][0][1][0] =
446  monoval_i[0][1] * monoval_j[1][1] * monoval[2][0] *
447  static_cast<double>(j + my_degree + 2);
448  grad_grads[start][0][1][1] =
449  monoval_i[0][0] * monoval_j[1][2] * monoval[2][0] *
450  static_cast<double>(j + my_degree + 2);
451  grad_grads[start][0][1][2] =
452  monoval_i[0][0] * monoval_j[1][1] * monoval[2][1] *
453  static_cast<double>(j + my_degree + 2);
454  grad_grads[start][0][2][0] =
455  monoval_i[0][1] * monoval_j[1][0] * monoval[2][1] *
456  static_cast<double>(j + my_degree + 2);
457  grad_grads[start][0][2][1] =
458  monoval_i[0][0] * monoval_j[1][1] * monoval[2][1] *
459  static_cast<double>(j + my_degree + 2);
460  grad_grads[start][0][2][2] =
461  monoval_i[0][0] * monoval_j[1][0] * monoval[2][2] *
462  static_cast<double>(j + my_degree + 2);
463  grad_grads[start][1][0][0] =
464  -monoval_i[0][3] * monoval_jplus[1][0] * monoval[2][0];
465  grad_grads[start][1][0][1] =
466  -monoval_i[0][2] * monoval_jplus[1][1] * monoval[2][0];
467  grad_grads[start][1][0][2] =
468  -monoval_i[0][2] * monoval_jplus[1][0] * monoval[2][1];
469  grad_grads[start][1][1][0] =
470  -monoval_i[0][2] * monoval_jplus[1][1] * monoval[2][0];
471  grad_grads[start][1][1][1] =
472  -monoval_i[0][1] * monoval_jplus[1][2] * monoval[2][0];
473  grad_grads[start][1][1][2] =
474  -monoval_i[0][1] * monoval_jplus[1][1] * monoval[2][1];
475  grad_grads[start][1][2][0] =
476  -monoval_i[0][2] * monoval_jplus[1][0] * monoval[2][1];
477  grad_grads[start][1][2][1] =
478  -monoval_i[0][1] * monoval_jplus[1][1] * monoval[2][1];
479  grad_grads[start][1][2][2] =
480  -monoval_i[0][1] * monoval_jplus[1][0] * monoval[2][2];
481  grad_grads[start][2][0][0] =
482  -monoval_i[0][3] * monoval_j[1][0] * monoval_plus[2][0];
483  grad_grads[start][2][0][1] =
484  -monoval_i[0][2] * monoval_j[1][1] * monoval_plus[2][0];
485  grad_grads[start][2][0][2] =
486  -monoval_i[0][2] * monoval_j[1][0] * monoval_plus[2][1];
487  grad_grads[start][2][1][0] =
488  -monoval_i[0][2] * monoval_j[1][1] * monoval_plus[2][0];
489  grad_grads[start][2][1][1] =
490  -monoval_i[0][1] * monoval_j[1][2] * monoval_plus[2][0];
491  grad_grads[start][2][1][2] =
492  -monoval_i[0][1] * monoval_j[1][1] * monoval_plus[2][1];
493  grad_grads[start][2][2][0] =
494  -monoval_i[0][2] * monoval_j[1][0] * monoval_plus[2][1];
495  grad_grads[start][2][2][1] =
496  -monoval_i[0][1] * monoval_j[1][1] * monoval_plus[2][1];
497  grad_grads[start][2][2][2] =
498  -monoval_i[0][1] * monoval_j[1][0] * monoval_plus[2][2];
499 
500  grad_grads[start + n_curls][0][0][0] =
501  -monoval_jplus[0][2] * monoval_i[1][1] * monoval[2][0];
502  grad_grads[start + n_curls][0][0][1] =
503  -monoval_jplus[0][1] * monoval_i[1][2] * monoval[2][0];
504  grad_grads[start + n_curls][0][0][2] =
505  -monoval_jplus[0][1] * monoval_i[1][1] * monoval[2][1];
506  grad_grads[start + n_curls][0][1][0] =
507  -monoval_jplus[0][1] * monoval_i[1][2] * monoval[2][0];
508  grad_grads[start + n_curls][0][1][1] =
509  -monoval_jplus[0][0] * monoval_i[1][3] * monoval[2][0];
510  grad_grads[start + n_curls][0][1][2] =
511  -monoval_jplus[0][0] * monoval_i[1][2] * monoval[2][1];
512  grad_grads[start + n_curls][0][2][0] =
513  -monoval_jplus[0][1] * monoval_i[1][1] * monoval[2][1];
514  grad_grads[start + n_curls][0][2][1] =
515  -monoval_jplus[0][0] * monoval_i[1][2] * monoval[2][1];
516  grad_grads[start + n_curls][0][2][2] =
517  -monoval_jplus[0][0] * monoval_i[1][1] * monoval[2][2];
518  grad_grads[start + n_curls][1][0][0] =
519  monoval_j[0][2] * monoval_i[1][0] * monoval[2][0] *
520  static_cast<double>(j + my_degree + 2);
521  grad_grads[start + n_curls][1][0][1] =
522  monoval_j[0][1] * monoval_i[1][1] * monoval[2][0] *
523  static_cast<double>(j + my_degree + 2);
524  grad_grads[start + n_curls][1][0][2] =
525  monoval_j[0][1] * monoval_i[1][0] * monoval[2][1] *
526  static_cast<double>(j + my_degree + 2);
527  grad_grads[start + n_curls][1][1][0] =
528  monoval_j[0][1] * monoval_i[1][1] * monoval[2][0] *
529  static_cast<double>(j + my_degree + 2);
530  grad_grads[start + n_curls][1][1][1] =
531  monoval_j[0][0] * monoval_i[1][2] * monoval[2][0] *
532  static_cast<double>(j + my_degree + 2);
533  grad_grads[start + n_curls][1][1][2] =
534  monoval_j[0][0] * monoval_i[1][1] * monoval[2][1] *
535  static_cast<double>(j + my_degree + 2);
536  grad_grads[start + n_curls][1][2][0] =
537  monoval_j[0][1] * monoval_i[1][0] * monoval[2][1] *
538  static_cast<double>(j + my_degree + 2);
539  grad_grads[start + n_curls][1][2][1] =
540  monoval_j[0][0] * monoval_i[1][1] * monoval[2][1] *
541  static_cast<double>(j + my_degree + 2);
542  grad_grads[start + n_curls][1][2][2] =
543  monoval_j[0][0] * monoval_i[1][0] * monoval[2][2] *
544  static_cast<double>(j + my_degree + 2);
545  grad_grads[start + n_curls][2][0][0] =
546  -monoval_j[0][2] * monoval_i[1][1] * monoval_plus[2][0];
547  grad_grads[start + n_curls][2][0][1] =
548  -monoval_j[0][1] * monoval_i[1][2] * monoval_plus[2][0];
549  grad_grads[start + n_curls][2][0][2] =
550  -monoval_j[0][1] * monoval_i[1][1] * monoval_plus[2][1];
551  grad_grads[start + n_curls][2][1][0] =
552  -monoval_j[0][1] * monoval_i[1][2] * monoval_plus[2][0];
553  grad_grads[start + n_curls][2][1][1] =
554  -monoval_j[0][0] * monoval_i[1][3] * monoval_plus[2][0];
555  grad_grads[start + n_curls][2][1][2] =
556  -monoval_j[0][0] * monoval_i[1][2] * monoval_plus[2][1];
557  grad_grads[start + n_curls][2][2][0] =
558  -monoval_j[0][1] * monoval_i[1][1] * monoval_plus[2][1];
559  grad_grads[start + n_curls][2][2][1] =
560  -monoval_j[0][0] * monoval_i[1][2] * monoval_plus[2][1];
561  grad_grads[start + n_curls][2][2][2] =
562  -monoval_j[0][0] * monoval_i[1][1] * monoval_plus[2][2];
563 
564  grad_grads[start + 2 * n_curls][0][0][0] =
565  -monoval_jplus[0][2] * monoval[1][0] * monoval_i[2][1];
566  grad_grads[start + 2 * n_curls][0][0][1] =
567  -monoval_jplus[0][1] * monoval[1][1] * monoval_i[2][1];
568  grad_grads[start + 2 * n_curls][0][0][2] =
569  -monoval_jplus[0][1] * monoval[1][0] * monoval_i[2][2];
570  grad_grads[start + 2 * n_curls][0][1][0] =
571  -monoval_jplus[0][1] * monoval[1][1] * monoval_i[2][1];
572  grad_grads[start + 2 * n_curls][0][1][1] =
573  -monoval_jplus[0][0] * monoval[1][2] * monoval_i[2][1];
574  grad_grads[start + 2 * n_curls][0][1][2] =
575  -monoval_jplus[0][0] * monoval[1][1] * monoval_i[2][2];
576  grad_grads[start + 2 * n_curls][0][2][0] =
577  -monoval_jplus[0][1] * monoval[1][0] * monoval_i[2][2];
578  grad_grads[start + 2 * n_curls][0][2][1] =
579  -monoval_jplus[0][0] * monoval[1][1] * monoval_i[2][2];
580  grad_grads[start + 2 * n_curls][0][2][2] =
581  -monoval_jplus[0][0] * monoval[1][0] * monoval_i[2][3];
582  grad_grads[start + 2 * n_curls][1][0][0] =
583  -monoval_j[0][2] * monoval_plus[1][0] * monoval_i[2][1];
584  grad_grads[start + 2 * n_curls][1][0][1] =
585  -monoval_j[0][1] * monoval_plus[1][1] * monoval_i[2][1];
586  grad_grads[start + 2 * n_curls][1][0][2] =
587  -monoval_j[0][1] * monoval_plus[1][0] * monoval_i[2][2];
588  grad_grads[start + 2 * n_curls][1][1][0] =
589  -monoval_j[0][1] * monoval_plus[1][1] * monoval_i[2][1];
590  grad_grads[start + 2 * n_curls][1][1][1] =
591  -monoval_j[0][0] * monoval_plus[1][2] * monoval_i[2][1];
592  grad_grads[start + 2 * n_curls][1][1][2] =
593  -monoval_j[0][0] * monoval_plus[1][1] * monoval_i[2][2];
594  grad_grads[start + 2 * n_curls][1][2][0] =
595  -monoval_j[0][1] * monoval_plus[1][0] * monoval_i[2][2];
596  grad_grads[start + 2 * n_curls][1][2][1] =
597  -monoval_j[0][0] * monoval_plus[1][1] * monoval_i[2][2];
598  grad_grads[start + 2 * n_curls][1][2][2] =
599  -monoval_j[0][0] * monoval_plus[1][0] * monoval_i[2][3];
600  grad_grads[start + 2 * n_curls][2][0][0] =
601  monoval_j[0][2] * monoval[1][0] * monoval_i[2][0] *
602  static_cast<double>(j + my_degree + 2);
603  grad_grads[start + 2 * n_curls][2][0][1] =
604  monoval_j[0][1] * monoval[1][1] * monoval_i[2][0] *
605  static_cast<double>(j + my_degree + 2);
606  grad_grads[start + 2 * n_curls][2][0][2] =
607  monoval_j[0][1] * monoval[1][0] * monoval_i[2][1] *
608  static_cast<double>(j + my_degree + 2);
609  grad_grads[start + 2 * n_curls][2][1][0] =
610  monoval_j[0][1] * monoval[1][1] * monoval_i[2][0] *
611  static_cast<double>(j + my_degree + 2);
612  grad_grads[start + 2 * n_curls][2][1][1] =
613  monoval_j[0][0] * monoval[1][2] * monoval_i[2][0] *
614  static_cast<double>(j + my_degree + 2);
615  grad_grads[start + 2 * n_curls][2][1][2] =
616  monoval_j[0][0] * monoval[1][1] * monoval_i[2][1] *
617  static_cast<double>(j + my_degree + 2);
618  grad_grads[start + 2 * n_curls][2][2][0] =
619  monoval_j[0][1] * monoval[1][0] * monoval_i[2][1] *
620  static_cast<double>(j + my_degree + 2);
621  grad_grads[start + 2 * n_curls][2][2][1] =
622  monoval_j[0][0] * monoval[1][1] * monoval_i[2][1] *
623  static_cast<double>(j + my_degree + 2);
624  grad_grads[start + 2 * n_curls][2][2][2] =
625  monoval_j[0][0] * monoval[1][0] * monoval_i[2][2] *
626  static_cast<double>(j + my_degree + 2);
627 
628  if (j != my_degree)
629  {
630  grad_grads[start + 1][0][0][0] =
631  monoval_i[0][2] * monoval[1][0] * monoval_j[2][0] *
632  static_cast<double>(j + my_degree + 2);
633  grad_grads[start + 1][0][0][1] =
634  monoval_i[0][1] * monoval[1][1] * monoval_j[2][0] *
635  static_cast<double>(j + my_degree + 2);
636  grad_grads[start + 1][0][0][2] =
637  monoval_i[0][1] * monoval[1][0] * monoval_j[2][1] *
638  static_cast<double>(j + my_degree + 2);
639  grad_grads[start + 1][0][1][0] =
640  monoval_i[0][1] * monoval[1][1] * monoval_j[2][0] *
641  static_cast<double>(j + my_degree + 2);
642  grad_grads[start + 1][0][1][1] =
643  monoval_i[0][0] * monoval[1][2] * monoval_j[2][0] *
644  static_cast<double>(j + my_degree + 2);
645  grad_grads[start + 1][0][1][2] =
646  monoval_i[0][0] * monoval[1][1] * monoval_j[2][1] *
647  static_cast<double>(j + my_degree + 2);
648  grad_grads[start + 1][0][2][0] =
649  monoval_i[0][1] * monoval[1][0] * monoval_j[2][1] *
650  static_cast<double>(j + my_degree + 2);
651  grad_grads[start + 1][0][2][1] =
652  monoval_i[0][0] * monoval[1][1] * monoval_j[2][1] *
653  static_cast<double>(j + my_degree + 2);
654  grad_grads[start + 1][0][2][2] =
655  monoval_i[0][0] * monoval[1][0] * monoval_j[2][2] *
656  static_cast<double>(j + my_degree + 2);
657  grad_grads[start + 1][1][0][0] =
658  -monoval_i[0][3] * monoval_plus[1][0] * monoval_j[2][0];
659  grad_grads[start + 1][1][0][1] =
660  -monoval_i[0][2] * monoval_plus[1][1] * monoval_j[2][0];
661  grad_grads[start + 1][1][0][2] =
662  -monoval_i[0][2] * monoval_plus[1][0] * monoval_j[2][1];
663  grad_grads[start + 1][1][1][0] =
664  -monoval_i[0][2] * monoval_plus[1][1] * monoval_j[2][0];
665  grad_grads[start + 1][1][1][1] =
666  -monoval_i[0][1] * monoval_plus[1][2] * monoval_j[2][0];
667  grad_grads[start + 1][1][1][2] =
668  -monoval_i[0][1] * monoval_plus[1][1] * monoval_j[2][1];
669  grad_grads[start + 1][1][2][0] =
670  -monoval_i[0][2] * monoval_plus[1][0] * monoval_j[2][1];
671  grad_grads[start + 1][1][2][1] =
672  -monoval_i[0][1] * monoval_plus[1][1] * monoval_j[2][1];
673  grad_grads[start + 1][1][2][2] =
674  -monoval_i[0][1] * monoval_plus[1][0] * monoval_j[2][2];
675  grad_grads[start + 1][2][0][0] =
676  -monoval_i[0][3] * monoval[1][0] * monoval_jplus[2][0];
677  grad_grads[start + 1][2][0][1] =
678  -monoval_i[0][2] * monoval[1][1] * monoval_jplus[2][0];
679  grad_grads[start + 1][2][0][2] =
680  -monoval_i[0][2] * monoval[1][0] * monoval_jplus[2][1];
681  grad_grads[start + 1][2][1][0] =
682  -monoval_i[0][2] * monoval[1][1] * monoval_jplus[2][0];
683  grad_grads[start + 1][2][1][1] =
684  -monoval_i[0][1] * monoval[1][2] * monoval_jplus[2][0];
685  grad_grads[start + 1][2][1][2] =
686  -monoval_i[0][1] * monoval[1][1] * monoval_jplus[2][1];
687  grad_grads[start + 1][2][2][0] =
688  -monoval_i[0][2] * monoval[1][0] * monoval_jplus[2][1];
689  grad_grads[start + 1][2][2][1] =
690  -monoval_i[0][1] * monoval[1][1] * monoval_jplus[2][1];
691  grad_grads[start + 1][2][2][2] =
692  -monoval_i[0][1] * monoval[1][0] * monoval_jplus[2][2];
693 
694  grad_grads[start + n_curls + 1][0][0][0] =
695  -monoval_plus[0][2] * monoval_i[1][1] * monoval_j[2][0];
696  grad_grads[start + n_curls + 1][0][0][1] =
697  -monoval_plus[0][1] * monoval_i[1][2] * monoval_j[2][0];
698  grad_grads[start + n_curls + 1][0][0][2] =
699  -monoval_plus[0][1] * monoval_i[1][1] * monoval_j[2][1];
700  grad_grads[start + n_curls + 1][0][1][0] =
701  -monoval_plus[0][1] * monoval_i[1][2] * monoval_j[2][0];
702  grad_grads[start + n_curls + 1][0][1][1] =
703  -monoval_plus[0][0] * monoval_i[1][3] * monoval_j[2][0];
704  grad_grads[start + n_curls + 1][0][1][2] =
705  -monoval_plus[0][0] * monoval_i[1][2] * monoval_j[2][1];
706  grad_grads[start + n_curls + 1][0][2][0] =
707  -monoval_plus[0][1] * monoval_i[1][1] * monoval_j[2][1];
708  grad_grads[start + n_curls + 1][0][2][1] =
709  -monoval_plus[0][0] * monoval_i[1][2] * monoval_j[2][1];
710  grad_grads[start + n_curls + 1][0][2][2] =
711  -monoval_plus[0][0] * monoval_i[1][1] * monoval_j[2][2];
712  grad_grads[start + n_curls + 1][1][0][0] =
713  monoval[0][2] * monoval_i[1][0] * monoval_j[2][0] *
714  static_cast<double>(j + my_degree + 2);
715  grad_grads[start + n_curls + 1][1][0][1] =
716  monoval[0][1] * monoval_i[1][1] * monoval_j[2][0] *
717  static_cast<double>(j + my_degree + 2);
718  grad_grads[start + n_curls + 1][1][0][2] =
719  monoval[0][1] * monoval_i[1][0] * monoval_j[2][1] *
720  static_cast<double>(j + my_degree + 2);
721  grad_grads[start + n_curls + 1][1][1][0] =
722  monoval[0][1] * monoval_i[1][1] * monoval_j[2][0] *
723  static_cast<double>(j + my_degree + 2);
724  grad_grads[start + n_curls + 1][1][1][1] =
725  monoval[0][0] * monoval_i[1][2] * monoval_j[2][0] *
726  static_cast<double>(j + my_degree + 2);
727  grad_grads[start + n_curls + 1][1][1][2] =
728  monoval[0][0] * monoval_i[1][1] * monoval_j[2][1] *
729  static_cast<double>(j + my_degree + 2);
730  grad_grads[start + n_curls + 1][1][2][0] =
731  monoval[0][1] * monoval_i[1][0] * monoval_j[2][1] *
732  static_cast<double>(j + my_degree + 2);
733  grad_grads[start + n_curls + 1][1][2][1] =
734  monoval[0][0] * monoval_i[1][1] * monoval_j[2][1] *
735  static_cast<double>(j + my_degree + 2);
736  grad_grads[start + n_curls + 1][1][2][2] =
737  monoval[0][0] * monoval_i[1][0] * monoval_j[2][2] *
738  static_cast<double>(j + my_degree + 2);
739  grad_grads[start + n_curls + 1][2][0][0] =
740  -monoval[0][2] * monoval_i[1][1] * monoval_jplus[2][0];
741  grad_grads[start + n_curls + 1][2][0][1] =
742  -monoval[0][1] * monoval_i[1][2] * monoval_jplus[2][0];
743  grad_grads[start + n_curls + 1][2][0][2] =
744  -monoval[0][1] * monoval_i[1][1] * monoval_jplus[2][1];
745  grad_grads[start + n_curls + 1][2][1][0] =
746  -monoval[0][1] * monoval_i[1][2] * monoval_jplus[2][0];
747  grad_grads[start + n_curls + 1][2][1][1] =
748  -monoval[0][0] * monoval_i[1][3] * monoval_jplus[2][0];
749  grad_grads[start + n_curls + 1][2][1][2] =
750  -monoval[0][0] * monoval_i[1][2] * monoval_jplus[2][1];
751  grad_grads[start + n_curls + 1][2][2][0] =
752  -monoval[0][1] * monoval_i[1][1] * monoval_jplus[2][1];
753  grad_grads[start + n_curls + 1][2][2][1] =
754  -monoval[0][0] * monoval_i[1][2] * monoval_jplus[2][1];
755  grad_grads[start + n_curls + 1][2][2][2] =
756  -monoval[0][0] * monoval_i[1][1] * monoval_jplus[2][2];
757 
758  grad_grads[start + 2 * n_curls + 1][0][0][0] =
759  -monoval_plus[0][2] * monoval_j[1][0] * monoval_i[2][1];
760  grad_grads[start + 2 * n_curls + 1][0][0][1] =
761  -monoval_plus[0][1] * monoval_j[1][1] * monoval_i[2][1];
762  grad_grads[start + 2 * n_curls + 1][0][0][2] =
763  -monoval_plus[0][1] * monoval_j[1][0] * monoval_i[2][2];
764  grad_grads[start + 2 * n_curls + 1][0][1][0] =
765  -monoval_plus[0][1] * monoval_j[1][1] * monoval_i[2][1];
766  grad_grads[start + 2 * n_curls + 1][0][1][1] =
767  -monoval_plus[0][0] * monoval_j[1][2] * monoval_i[2][1];
768  grad_grads[start + 2 * n_curls + 1][0][1][2] =
769  -monoval_plus[0][0] * monoval_j[1][1] * monoval_i[2][2];
770  grad_grads[start + 2 * n_curls + 1][0][2][0] =
771  -monoval_plus[0][1] * monoval_j[1][0] * monoval_i[2][2];
772  grad_grads[start + 2 * n_curls + 1][0][2][1] =
773  -monoval_plus[0][0] * monoval_j[1][1] * monoval_i[2][2];
774  grad_grads[start + 2 * n_curls + 1][0][2][2] =
775  -monoval_plus[0][0] * monoval_j[1][0] * monoval_i[2][3];
776  grad_grads[start + 2 * n_curls + 1][1][0][0] =
777  -monoval[0][2] * monoval_jplus[1][0] * monoval_i[2][1];
778  grad_grads[start + 2 * n_curls + 1][1][0][1] =
779  -monoval[0][1] * monoval_jplus[1][1] * monoval_i[2][1];
780  grad_grads[start + 2 * n_curls + 1][1][0][2] =
781  -monoval[0][1] * monoval_jplus[1][0] * monoval_i[2][2];
782  grad_grads[start + 2 * n_curls + 1][1][1][0] =
783  -monoval[0][1] * monoval_jplus[1][1] * monoval_i[2][1];
784  grad_grads[start + 2 * n_curls + 1][1][1][1] =
785  -monoval[0][0] * monoval_jplus[1][2] * monoval_i[2][1];
786  grad_grads[start + 2 * n_curls + 1][1][1][2] =
787  -monoval[0][0] * monoval_jplus[1][1] * monoval_i[2][2];
788  grad_grads[start + 2 * n_curls + 1][1][2][0] =
789  -monoval[0][1] * monoval_jplus[1][0] * monoval_i[2][2];
790  grad_grads[start + 2 * n_curls + 1][1][2][1] =
791  -monoval[0][0] * monoval_jplus[1][1] * monoval_i[2][2];
792  grad_grads[start + 2 * n_curls + 1][1][2][2] =
793  -monoval[0][0] * monoval_jplus[1][0] * monoval_i[2][3];
794  grad_grads[start + 2 * n_curls + 1][2][0][0] =
795  monoval[0][2] * monoval_j[1][0] * monoval_i[2][0] *
796  static_cast<double>(j + my_degree + 2);
797  grad_grads[start + 2 * n_curls + 1][2][0][1] =
798  monoval[0][1] * monoval_j[1][1] * monoval_i[2][0] *
799  static_cast<double>(j + my_degree + 2);
800  grad_grads[start + 2 * n_curls + 1][2][0][2] =
801  monoval[0][1] * monoval_j[1][0] * monoval_i[2][1] *
802  static_cast<double>(j + my_degree + 2);
803  grad_grads[start + 2 * n_curls + 1][2][1][0] =
804  monoval[0][1] * monoval_j[1][1] * monoval_i[2][0] *
805  static_cast<double>(j + my_degree + 2);
806  grad_grads[start + 2 * n_curls + 1][2][1][1] =
807  monoval[0][0] * monoval_j[1][2] * monoval_i[2][0] *
808  static_cast<double>(j + my_degree + 2);
809  grad_grads[start + 2 * n_curls + 1][2][1][2] =
810  monoval[0][0] * monoval_j[1][1] * monoval_i[2][1] *
811  static_cast<double>(j + my_degree + 2);
812  grad_grads[start + 2 * n_curls + 1][2][2][0] =
813  monoval[0][1] * monoval_j[1][0] * monoval_i[2][1] *
814  static_cast<double>(j + my_degree + 2);
815  grad_grads[start + 2 * n_curls + 1][2][2][1] =
816  monoval[0][0] * monoval_j[1][1] * monoval_i[2][1] *
817  static_cast<double>(j + my_degree + 2);
818  grad_grads[start + 2 * n_curls + 1][2][2][2] =
819  monoval[0][0] * monoval_j[1][0] * monoval_i[2][2] *
820  static_cast<double>(j + my_degree + 2);
821  }
822  }
823 
824  if (j == my_degree)
825  start += 1;
826  else
827  start += 2;
828  }
829  }
830  Assert(start == this->n() - 2 * n_curls, ExcInternalError());
831  }
832 }
833 
834 
835 
836 template <int dim>
837 unsigned int
839 {
840  if (dim == 1 || dim == 2 || dim == 3)
841  return dim * Utilities::fixed_power<dim>(k + 1);
842 
843  Assert(false, ExcNotImplemented());
844  return 0;
845 }
846 
847 
848 template <int dim>
849 std::unique_ptr<TensorPolynomialsBase<dim>>
851 {
852  return std_cxx14::make_unique<PolynomialsRT_Bubbles<dim>>(*this);
853 }
854 
855 
856 template class PolynomialsRT_Bubbles<1>;
857 template class PolynomialsRT_Bubbles<2>;
858 template class PolynomialsRT_Bubbles<3>;
859 
860 
PolynomialsRT_Bubbles
Definition: polynomials_rt_bubbles.h:88
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
Physics::Elasticity::Kinematics::d
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
polynomials_rt_bubbles.h
quadrature_lib.h
thread_management.h
PolynomialsRT_Bubbles::monomials
std::vector< Polynomials::Polynomial< double > > monomials
Definition: polynomials_rt_bubbles.h:147
Tensor< 1, dim >
PolynomialsRT_Bubbles::n_polynomials
static unsigned int n_polynomials(const unsigned int degree)
Definition: polynomials_rt_bubbles.cc:838
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
TensorPolynomialsBase
Definition: tensor_polynomials_base.h:62
value
static const bool value
Definition: dof_tools_constraints.cc:433
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
PolynomialsRT_Bubbles::evaluate
void evaluate(const Point< dim > &unit_point, std::vector< Tensor< 1, dim >> &values, std::vector< Tensor< 2, dim >> &grads, std::vector< Tensor< 3, dim >> &grad_grads, std::vector< Tensor< 4, dim >> &third_derivatives, std::vector< Tensor< 5, dim >> &fourth_derivatives) const override
Definition: polynomials_rt_bubbles.cc:45
PolynomialsRT_Bubbles::clone
virtual std::unique_ptr< TensorPolynomialsBase< dim > > clone() const override
Definition: polynomials_rt_bubbles.cc:850
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
StandardExceptions::ExcImpossibleInDim
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
Polynomials::Monomial
Definition: polynomial.h:299
StandardExceptions::ExcDimensionMismatch
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Point< dim >
PolynomialsRT_Bubbles::PolynomialsRT_Bubbles
PolynomialsRT_Bubbles(const unsigned int k)
Definition: polynomials_rt_bubbles.cc:30
memory.h
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359