Reference documentation for deal.II version 9.1.0-pre
tensor_product_polynomials_bubbles.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2012 - 2017 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 
17 #include <deal.II/base/exceptions.h>
18 #include <deal.II/base/table.h>
19 #include <deal.II/base/tensor_product_polynomials_bubbles.h>
20 
21 DEAL_II_NAMESPACE_OPEN
22 
23 
24 
25 /* ------------------- TensorProductPolynomialsBubbles -------------- */
26 
27 
28 template <int dim>
29 double
31  const Point<dim> & p) const
32 {
33  const unsigned int q_degree = this->polynomials.size() - 1;
34  const unsigned int max_q_indices = this->n_tensor_pols;
35  const unsigned int n_bubbles = ((q_degree <= 1) ? 1 : dim);
36  (void)n_bubbles;
37  Assert(i < max_q_indices + n_bubbles, ExcInternalError());
38 
39  // treat the regular basis functions
40  if (i < max_q_indices)
42 
43  const unsigned int comp = i - this->n_tensor_pols;
44 
45  // compute \prod_{i=1}^d 4*(1-x_i^2)(p)
46  double value = 1.;
47  for (unsigned int j = 0; j < dim; ++j)
48  value *= 4 * p(j) * (1 - p(j));
49  // and multiply with (2x_i-1)^{r-1}
50  for (unsigned int i = 0; i < q_degree - 1; ++i)
51  value *= 2 * p(comp) - 1;
52  return value;
53 }
54 
55 
56 
57 template <>
58 double
60  const Point<0> &) const
61 {
62  Assert(false, ExcNotImplemented());
63  return 0.;
64 }
65 
66 
67 template <int dim>
70  const Point<dim> & p) const
71 {
72  const unsigned int q_degree = this->polynomials.size() - 1;
73  const unsigned int max_q_indices = this->n_tensor_pols;
74  const unsigned int n_bubbles = ((q_degree <= 1) ? 1 : dim);
75  (void)n_bubbles;
76  Assert(i < max_q_indices + n_bubbles, ExcInternalError());
77 
78  // treat the regular basis functions
79  if (i < max_q_indices)
81 
82  const unsigned int comp = i - this->n_tensor_pols;
83  Tensor<1, dim> grad;
84 
85  for (unsigned int d = 0; d < dim; ++d)
86  {
87  grad[d] = 1.;
88  // compute grad(4*\prod_{i=1}^d (x_i(1-x_i)))(p)
89  for (unsigned j = 0; j < dim; ++j)
90  grad[d] *= (d == j ? 4 * (1 - 2 * p(j)) : 4 * p(j) * (1 - p(j)));
91  // and multiply with (2*x_i-1)^{r-1}
92  for (unsigned int i = 0; i < q_degree - 1; ++i)
93  grad[d] *= 2 * p(comp) - 1;
94  }
95 
96  if (q_degree >= 2)
97  {
98  // add \prod_{i=1}^d 4*(x_i(1-x_i))(p)
99  double value = 1.;
100  for (unsigned int j = 0; j < dim; ++j)
101  value *= 4 * p(j) * (1 - p(j));
102  // and multiply with grad(2*x_i-1)^{r-1}
103  double tmp = value * 2 * (q_degree - 1);
104  for (unsigned int i = 0; i < q_degree - 2; ++i)
105  tmp *= 2 * p(comp) - 1;
106  grad[comp] += tmp;
107  }
108 
109  return grad;
110 }
111 
112 
113 
114 template <int dim>
117  const unsigned int i,
118  const Point<dim> & p) const
119 {
120  const unsigned int q_degree = this->polynomials.size() - 1;
121  const unsigned int max_q_indices = this->n_tensor_pols;
122  const unsigned int n_bubbles = ((q_degree <= 1) ? 1 : dim);
123  (void)n_bubbles;
124  Assert(i < max_q_indices + n_bubbles, ExcInternalError());
125 
126  // treat the regular basis functions
127  if (i < max_q_indices)
129 
130  const unsigned int comp = i - this->n_tensor_pols;
131 
132  double v[dim + 1][3];
133  {
134  for (unsigned int c = 0; c < dim; ++c)
135  {
136  v[c][0] = 4 * p(c) * (1 - p(c));
137  v[c][1] = 4 * (1 - 2 * p(c));
138  v[c][2] = -8;
139  }
140 
141  double tmp = 1.;
142  for (unsigned int i = 0; i < q_degree - 1; ++i)
143  tmp *= 2 * p(comp) - 1;
144  v[dim][0] = tmp;
145 
146  if (q_degree >= 2)
147  {
148  double tmp = 2 * (q_degree - 1);
149  for (unsigned int i = 0; i < q_degree - 2; ++i)
150  tmp *= 2 * p(comp) - 1;
151  v[dim][1] = tmp;
152  }
153  else
154  v[dim][1] = 0.;
155 
156  if (q_degree >= 3)
157  {
158  double tmp = 4 * (q_degree - 2) * (q_degree - 1);
159  for (unsigned int i = 0; i < q_degree - 3; ++i)
160  tmp *= 2 * p(comp) - 1;
161  v[dim][2] = tmp;
162  }
163  else
164  v[dim][2] = 0.;
165  }
166 
167  // calculate (\partial_j \partial_k \psi) * monomial
168  Tensor<2, dim> grad_grad_1;
169  for (unsigned int d1 = 0; d1 < dim; ++d1)
170  for (unsigned int d2 = 0; d2 < dim; ++d2)
171  {
172  grad_grad_1[d1][d2] = v[dim][0];
173  for (unsigned int x = 0; x < dim; ++x)
174  {
175  unsigned int derivative = 0;
176  if (d1 == x || d2 == x)
177  {
178  if (d1 == d2)
179  derivative = 2;
180  else
181  derivative = 1;
182  }
183  grad_grad_1[d1][d2] *= v[x][derivative];
184  }
185  }
186 
187  // calculate (\partial_j \psi) *(\partial_k monomial)
188  // and (\partial_k \psi) *(\partial_j monomial)
189  Tensor<2, dim> grad_grad_2;
190  Tensor<2, dim> grad_grad_3;
191  for (unsigned int d = 0; d < dim; ++d)
192  {
193  grad_grad_2[d][comp] = v[dim][1];
194  grad_grad_3[comp][d] = v[dim][1];
195  for (unsigned int x = 0; x < dim; ++x)
196  {
197  grad_grad_2[d][comp] *= v[x][d == x];
198  grad_grad_3[comp][d] *= v[x][d == x];
199  }
200  }
201 
202  // calculate \psi *(\partial j \partial_k monomial) and sum
203  Tensor<2, dim> grad_grad;
204  double psi_value = 1.;
205  for (unsigned int x = 0; x < dim; ++x)
206  psi_value *= v[x][0];
207 
208  for (unsigned int d1 = 0; d1 < dim; ++d1)
209  for (unsigned int d2 = 0; d2 < dim; ++d2)
210  grad_grad[d1][d2] =
211  grad_grad_1[d1][d2] + grad_grad_2[d1][d2] + grad_grad_3[d1][d2];
212  grad_grad[comp][comp] += psi_value * v[dim][2];
213 
214  return grad_grad;
215 }
216 
217 template <int dim>
218 void
220  const Point<dim> & p,
221  std::vector<double> & values,
222  std::vector<Tensor<1, dim>> &grads,
223  std::vector<Tensor<2, dim>> &grad_grads,
224  std::vector<Tensor<3, dim>> &third_derivatives,
225  std::vector<Tensor<4, dim>> &fourth_derivatives) const
226 {
227  const unsigned int q_degree = this->polynomials.size() - 1;
228  const unsigned int max_q_indices = this->n_tensor_pols;
229  (void)max_q_indices;
230  const unsigned int n_bubbles = ((q_degree <= 1) ? 1 : dim);
231  Assert(values.size() == max_q_indices + n_bubbles || values.size() == 0,
232  ExcDimensionMismatch2(values.size(), max_q_indices + n_bubbles, 0));
233  Assert(grads.size() == max_q_indices + n_bubbles || grads.size() == 0,
234  ExcDimensionMismatch2(grads.size(), max_q_indices + n_bubbles, 0));
235  Assert(
236  grad_grads.size() == max_q_indices + n_bubbles || grad_grads.size() == 0,
237  ExcDimensionMismatch2(grad_grads.size(), max_q_indices + n_bubbles, 0));
238  Assert(third_derivatives.size() == max_q_indices + n_bubbles ||
239  third_derivatives.size() == 0,
240  ExcDimensionMismatch2(third_derivatives.size(),
241  max_q_indices + n_bubbles,
242  0));
243  Assert(fourth_derivatives.size() == max_q_indices + n_bubbles ||
244  fourth_derivatives.size() == 0,
245  ExcDimensionMismatch2(fourth_derivatives.size(),
246  max_q_indices + n_bubbles,
247  0));
248 
249  bool do_values = false, do_grads = false, do_grad_grads = false;
250  bool do_3rd_derivatives = false, do_4th_derivatives = false;
251  if (values.empty() == false)
252  {
253  values.resize(this->n_tensor_pols);
254  do_values = true;
255  }
256  if (grads.empty() == false)
257  {
258  grads.resize(this->n_tensor_pols);
259  do_grads = true;
260  }
261  if (grad_grads.empty() == false)
262  {
263  grad_grads.resize(this->n_tensor_pols);
264  do_grad_grads = true;
265  }
266  if (third_derivatives.empty() == false)
267  {
268  third_derivatives.resize(this->n_tensor_pols);
269  do_3rd_derivatives = true;
270  }
271  if (fourth_derivatives.empty() == false)
272  {
273  fourth_derivatives.resize(this->n_tensor_pols);
274  do_4th_derivatives = true;
275  }
276 
278  p, values, grads, grad_grads, third_derivatives, fourth_derivatives);
279 
280  for (unsigned int i = this->n_tensor_pols;
281  i < this->n_tensor_pols + n_bubbles;
282  ++i)
283  {
284  if (do_values)
285  values.push_back(compute_value(i, p));
286  if (do_grads)
287  grads.push_back(compute_grad(i, p));
288  if (do_grad_grads)
289  grad_grads.push_back(compute_grad_grad(i, p));
290  if (do_3rd_derivatives)
291  third_derivatives.push_back(compute_derivative<3>(i, p));
292  if (do_4th_derivatives)
293  fourth_derivatives.push_back(compute_derivative<4>(i, p));
294  }
295 }
296 
297 
298 /* ------------------- explicit instantiations -------------- */
302 
303 DEAL_II_NAMESPACE_CLOSE
double compute_value(const unsigned int i, const Point< dim > &p) const
double compute_value(const unsigned int i, const Point< dim > &p) const
static::ExceptionBase & ExcDimensionMismatch2(int arg1, int arg2, int arg3)
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const
void compute(const Point< dim > &unit_point, std::vector< double > &values, std::vector< Tensor< 1, dim >> &grads, std::vector< Tensor< 2, dim >> &grad_grads, std::vector< Tensor< 3, dim >> &third_derivatives, std::vector< Tensor< 4, dim >> &fourth_derivatives) const
#define Assert(cond, exc)
Definition: exceptions.h:1227
void compute(const Point< dim > &unit_point, std::vector< double > &values, std::vector< Tensor< 1, dim >> &grads, std::vector< Tensor< 2, dim >> &grad_grads, std::vector< Tensor< 3, dim >> &third_derivatives, std::vector< Tensor< 4, dim >> &fourth_derivatives) const
static::ExceptionBase & ExcNotImplemented()
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const
static::ExceptionBase & ExcInternalError()