Reference documentation for deal.II version 8.5.1
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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #include <deal.II/base/tensor_product_polynomials_bubbles.h>
18 #include <deal.II/base/exceptions.h>
19 #include <deal.II/base/table.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 Point<dim> &p) const
118 {
119  const unsigned int q_degree = this->polynomials.size()-1;
120  const unsigned int max_q_indices = this->n_tensor_pols;
121  const unsigned int n_bubbles = ((q_degree<=1)?1:dim);
122  (void)n_bubbles;
123  Assert (i<max_q_indices+n_bubbles, ExcInternalError());
124 
125  // treat the regular basis functions
126  if (i<max_q_indices)
128 
129  const unsigned int comp = i - this->n_tensor_pols;
130 
131  double v [dim+1][3];
132  {
133  for (unsigned int c=0; c<dim; ++c)
134  {
135  v[c][0] = 4*p(c)*(1-p(c));
136  v[c][1] = 4*(1-2*p(c));
137  v[c][2] = -8;
138  }
139 
140  double tmp=1.;
141  for (unsigned int i=0; i<q_degree-1; ++i)
142  tmp *= 2*p(comp)-1;
143  v[dim][0] = tmp;
144 
145  if (q_degree>=2)
146  {
147  double tmp = 2*(q_degree-1);
148  for (unsigned int i=0; i<q_degree-2; ++i)
149  tmp *= 2*p(comp)-1;
150  v[dim][1] = tmp;
151  }
152  else
153  v[dim][1] = 0.;
154 
155  if (q_degree>=3)
156  {
157  double tmp=4*(q_degree-2)*(q_degree-1);
158  for (unsigned int i=0; i<q_degree-3; ++i)
159  tmp *= 2*p(comp)-1;
160  v[dim][2] = tmp;
161  }
162  else
163  v[dim][2] = 0.;
164  }
165 
166  //calculate (\partial_j \partial_k \psi) * monomial
167  Tensor<2,dim> grad_grad_1;
168  for (unsigned int d1=0; d1<dim; ++d1)
169  for (unsigned int d2=0; d2<dim; ++d2)
170  {
171  grad_grad_1[d1][d2] = v[dim][0];
172  for (unsigned int x=0; x<dim; ++x)
173  {
174  unsigned int derivative=0;
175  if (d1==x || d2==x)
176  {
177  if (d1==d2)
178  derivative=2;
179  else
180  derivative=1;
181  }
182  grad_grad_1[d1][d2] *= v[x][derivative];
183  }
184  }
185 
186  //calculate (\partial_j \psi) *(\partial_k monomial)
187  // and (\partial_k \psi) *(\partial_j monomial)
188  Tensor<2,dim> grad_grad_2;
189  Tensor<2,dim> grad_grad_3;
190  for (unsigned int d=0; d<dim; ++d)
191  {
192  grad_grad_2[d][comp] = v[dim][1];
193  grad_grad_3[comp][d] = v[dim][1];
194  for (unsigned int x=0; x<dim; ++x)
195  {
196  grad_grad_2[d][comp] *= v[x][d==x];
197  grad_grad_3[comp][d] *= v[x][d==x];
198  }
199  }
200 
201  //calculate \psi *(\partial j \partial_k monomial) and sum
202  Tensor<2,dim> grad_grad;
203  double psi_value = 1.;
204  for (unsigned int x=0; x<dim; ++x)
205  psi_value *= v[x][0];
206 
207  for (unsigned int d1=0; d1<dim; ++d1)
208  for (unsigned int d2=0; d2<dim; ++d2)
209  grad_grad[d1][d2] = grad_grad_1[d1][d2]
210  +grad_grad_2[d1][d2]
211  +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 compute (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 (grad_grads.size()==max_q_indices+n_bubbles || grad_grads.size()==0,
236  ExcDimensionMismatch2(grad_grads.size(), max_q_indices+n_bubbles, 0));
237  Assert (third_derivatives.size()==max_q_indices+n_bubbles || third_derivatives.size()==0,
238  ExcDimensionMismatch2(third_derivatives.size(), max_q_indices+n_bubbles, 0));
239  Assert (fourth_derivatives.size()==max_q_indices+n_bubbles || fourth_derivatives.size()==0,
240  ExcDimensionMismatch2(fourth_derivatives.size(), max_q_indices+n_bubbles, 0));
241 
242  bool do_values = false, do_grads = false, do_grad_grads = false;
243  bool do_3rd_derivatives = false, do_4th_derivatives = false;
244  if (values.empty() == false)
245  {
246  values.resize(this->n_tensor_pols);
247  do_values = true;
248  }
249  if (grads.empty() == false)
250  {
251  grads.resize(this->n_tensor_pols);
252  do_grads = true;
253  }
254  if (grad_grads.empty() == false)
255  {
256  grad_grads.resize(this->n_tensor_pols);
257  do_grad_grads = true;
258  }
259  if (third_derivatives.empty() == false)
260  {
261  third_derivatives.resize(this->n_tensor_pols);
262  do_3rd_derivatives = true;
263  }
264  if (fourth_derivatives.empty() == false)
265  {
266  fourth_derivatives.resize(this->n_tensor_pols);
267  do_4th_derivatives = true;
268  }
269 
270  this->TensorProductPolynomials<dim>::compute(p, values, grads, grad_grads, third_derivatives, fourth_derivatives);
271 
272  for (unsigned int i=this->n_tensor_pols; i<this->n_tensor_pols+n_bubbles; ++i)
273  {
274  if (do_values)
275  values.push_back(compute_value(i,p));
276  if (do_grads)
277  grads.push_back(compute_grad(i,p));
278  if (do_grad_grads)
279  grad_grads.push_back(compute_grad_grad(i,p));
280  if (do_3rd_derivatives)
281  third_derivatives.push_back(compute_derivative<3>(i,p));
282  if (do_4th_derivatives)
283  fourth_derivatives.push_back(compute_derivative<4>(i,p));
284  }
285 }
286 
287 
288 /* ------------------- explicit instantiations -------------- */
292 
293 DEAL_II_NAMESPACE_CLOSE
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
static ::ExceptionBase & ExcDimensionMismatch2(int arg1, int arg2, int arg3)
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:313
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
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const
double compute_value(const unsigned int i, const Point< dim > &p) const
static ::ExceptionBase & ExcNotImplemented()
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const
double compute_value(const unsigned int i, const Point< dim > &p) const
static ::ExceptionBase & ExcInternalError()