Reference documentation for deal.II version 8.4.2
tensor_product_polynomials_bubbles.cc
1 // ---------------------------------------------------------------------
2 // @f$Id@f$
3 //
4 // Copyright (C) 2012 - 2015 by the deal.II authors
5 //
6 // This file is part of the deal.II library.
7 //
8 // The deal.II library is free software; you can use it, redistribute
9 // it, and/or modify it under the terms of the GNU Lesser General
10 // Public License as published by the Free Software Foundation; either
11 // version 2.1 of the License, or (at your option) any later version.
12 // The full text of the license can be found in the file LICENSE at
13 // the top level of the deal.II distribution.
14 //
15 // ---------------------------------------------------------------------
16 
17 
18 #include <deal.II/base/tensor_product_polynomials_bubbles.h>
19 #include <deal.II/base/exceptions.h>
20 #include <deal.II/base/table.h>
21 
22 DEAL_II_NAMESPACE_OPEN
23 
24 
25 
26 /* ------------------- TensorProductPolynomialsBubbles -------------- */
27 
28 
29 template <int dim>
30 double
32  const Point<dim> &p) const
33 {
34  const unsigned int q_degree = this->polynomials.size()-1;
35  const unsigned int max_q_indices = this->n_tensor_pols;
36  const unsigned int n_bubbles = ((q_degree<=1)?1:dim);
37  (void)n_bubbles;
38  Assert (i<max_q_indices+n_bubbles, ExcInternalError());
39 
40  // treat the regular basis functions
41  if (i<max_q_indices)
43 
44  const unsigned int comp = i - this->n_tensor_pols;
45 
46  //compute \prod_{i=1}^d 4*(1-x_i^2)(p)
47  double value=1.;
48  for (unsigned int j=0; j<dim; ++j)
49  value*=4*p(j)*(1-p(j));
50  // and multiply with (2x_i-1)^{r-1}
51  for (unsigned int i=0; i<q_degree-1; ++i)
52  value*=2*p(comp)-1;
53  return value;
54 }
55 
56 
57 
58 template <>
59 double
61  const Point<0> &) const
62 {
63  Assert (false, ExcNotImplemented());
64  return 0.;
65 }
66 
67 
68 template <int dim>
71  const Point<dim> &p) const
72 {
73  const unsigned int q_degree = this->polynomials.size()-1;
74  const unsigned int max_q_indices = this->n_tensor_pols;
75  const unsigned int n_bubbles = ((q_degree<=1)?1:dim);
76  (void)n_bubbles;
77  Assert (i<max_q_indices+n_bubbles, ExcInternalError());
78 
79  // treat the regular basis functions
80  if (i<max_q_indices)
82 
83  const unsigned int comp = i - this->n_tensor_pols;
84  Tensor<1,dim> grad;
85 
86  for (unsigned int d=0; d<dim ; ++d)
87  {
88  grad[d] = 1.;
89  //compute grad(4*\prod_{i=1}^d (x_i(1-x_i)))(p)
90  for (unsigned j=0; j<dim; ++j)
91  grad[d] *= (d==j ? 4*(1-2*p(j)) : 4*p(j)*(1-p(j)));
92  // and multiply with (2*x_i-1)^{r-1}
93  for (unsigned int i=0; i<q_degree-1; ++i)
94  grad[d]*=2*p(comp)-1;
95  }
96 
97  if (q_degree>=2)
98  {
99  //add \prod_{i=1}^d 4*(x_i(1-x_i))(p)
100  double value=1.;
101  for (unsigned int j=0; j < dim; ++j)
102  value*=4*p(j)*(1-p(j));
103  //and multiply with grad(2*x_i-1)^{r-1}
104  double tmp=value*2*(q_degree-1);
105  for (unsigned int i=0; i<q_degree-2; ++i)
106  tmp*=2*p(comp)-1;
107  grad[comp]+=tmp;
108  }
109 
110  return grad;
111 }
112 
113 
114 
115 template <int dim>
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] = grad_grad_1[d1][d2]
211  +grad_grad_2[d1][d2]
212  +grad_grad_3[d1][d2];
213  grad_grad[comp][comp]+=psi_value*v[dim][2];
214 
215  return grad_grad;
216 }
217 
218 template <int dim>
219 void
221 compute (const Point<dim> &p,
222  std::vector<double> &values,
223  std::vector<Tensor<1,dim> > &grads,
224  std::vector<Tensor<2,dim> > &grad_grads,
225  std::vector<Tensor<3,dim> > &third_derivatives,
226  std::vector<Tensor<4,dim> > &fourth_derivatives) const
227 {
228  const unsigned int q_degree = this->polynomials.size()-1;
229  const unsigned int max_q_indices = this->n_tensor_pols;
230  (void) max_q_indices;
231  const unsigned int n_bubbles = ((q_degree<=1)?1:dim);
232  Assert (values.size()==max_q_indices+n_bubbles || values.size()==0,
233  ExcDimensionMismatch2(values.size(), max_q_indices+n_bubbles, 0));
234  Assert (grads.size()==max_q_indices+n_bubbles || grads.size()==0,
235  ExcDimensionMismatch2(grads.size(), max_q_indices+n_bubbles, 0));
236  Assert (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 || third_derivatives.size()==0,
239  ExcDimensionMismatch2(third_derivatives.size(), max_q_indices+n_bubbles, 0));
240  Assert (fourth_derivatives.size()==max_q_indices+n_bubbles || fourth_derivatives.size()==0,
241  ExcDimensionMismatch2(fourth_derivatives.size(), max_q_indices+n_bubbles, 0));
242 
243  bool do_values = false, do_grads = false, do_grad_grads = false;
244  bool do_3rd_derivatives = false, do_4th_derivatives = false;
245  if (values.empty() == false)
246  {
247  values.resize(this->n_tensor_pols);
248  do_values = true;
249  }
250  if (grads.empty() == false)
251  {
252  grads.resize(this->n_tensor_pols);
253  do_grads = true;
254  }
255  if (grad_grads.empty() == false)
256  {
257  grad_grads.resize(this->n_tensor_pols);
258  do_grad_grads = true;
259  }
260  if (third_derivatives.empty() == false)
261  {
262  third_derivatives.resize(this->n_tensor_pols);
263  do_3rd_derivatives = true;
264  }
265  if (fourth_derivatives.empty() == false)
266  {
267  fourth_derivatives.resize(this->n_tensor_pols);
268  do_4th_derivatives = true;
269  }
270 
271  this->TensorProductPolynomials<dim>::compute(p, values, grads, grad_grads, third_derivatives, fourth_derivatives);
272 
273  for (unsigned int i=this->n_tensor_pols; i<this->n_tensor_pols+n_bubbles; ++i)
274  {
275  if (do_values)
276  values.push_back(compute_value(i,p));
277  if (do_grads)
278  grads.push_back(compute_grad(i,p));
279  if (do_grad_grads)
280  grad_grads.push_back(compute_grad_grad(i,p));
281  if (do_3rd_derivatives)
282  third_derivatives.push_back(compute_derivative<3>(i,p));
283  if (do_4th_derivatives)
284  fourth_derivatives.push_back(compute_derivative<4>(i,p));
285  }
286 }
287 
288 
289 /* ------------------- explicit instantiations -------------- */
293 
294 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
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:294
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
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