17 #include <deal.II/base/tensor_product_polynomials_bubbles.h> 18 #include <deal.II/base/exceptions.h> 19 #include <deal.II/base/table.h> 21 DEAL_II_NAMESPACE_OPEN
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);
43 const unsigned int comp = i - this->n_tensor_pols;
47 for (
unsigned int j=0; j<dim; ++j)
48 value*=4*p(j)*(1-p(j));
50 for (
unsigned int i=0; i<q_degree-1; ++i)
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);
82 const unsigned int comp = i - this->n_tensor_pols;
85 for (
unsigned int d=0; d<dim ; ++d)
89 for (
unsigned j=0; j<dim; ++j)
90 grad[d] *= (d==j ? 4*(1-2*p(j)) : 4*p(j)*(1-p(j)));
92 for (
unsigned int i=0; i<q_degree-1; ++i)
100 for (
unsigned int j=0; j < dim; ++j)
101 value*=4*p(j)*(1-p(j));
103 double tmp=value*2*(q_degree-1);
104 for (
unsigned int i=0; i<q_degree-2; ++i)
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);
129 const unsigned int comp = i - this->n_tensor_pols;
133 for (
unsigned int c=0; c<dim; ++c)
135 v[c][0] = 4*p(c)*(1-p(c));
136 v[c][1] = 4*(1-2*p(c));
141 for (
unsigned int i=0; i<q_degree-1; ++i)
147 double tmp = 2*(q_degree-1);
148 for (
unsigned int i=0; i<q_degree-2; ++i)
157 double tmp=4*(q_degree-2)*(q_degree-1);
158 for (
unsigned int i=0; i<q_degree-3; ++i)
168 for (
unsigned int d1=0; d1<dim; ++d1)
169 for (
unsigned int d2=0; d2<dim; ++d2)
171 grad_grad_1[d1][d2] = v[dim][0];
172 for (
unsigned int x=0; x<dim; ++x)
174 unsigned int derivative=0;
182 grad_grad_1[d1][d2] *= v[x][derivative];
190 for (
unsigned int d=0; d<dim; ++d)
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)
196 grad_grad_2[d][comp] *= v[x][d==x];
197 grad_grad_3[comp][d] *= v[x][d==x];
203 double psi_value = 1.;
204 for (
unsigned int x=0; x<dim; ++x)
205 psi_value *= v[x][0];
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]
211 +grad_grad_3[d1][d2];
212 grad_grad[comp][comp]+=psi_value*v[dim][2];
221 std::vector<double> &values,
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,
233 Assert (grads.size()==max_q_indices+n_bubbles || grads.size()==0,
235 Assert (grad_grads.size()==max_q_indices+n_bubbles || grad_grads.size()==0,
237 Assert (third_derivatives.size()==max_q_indices+n_bubbles || third_derivatives.size()==0,
239 Assert (fourth_derivatives.size()==max_q_indices+n_bubbles || fourth_derivatives.size()==0,
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)
246 values.resize(this->n_tensor_pols);
249 if (grads.empty() ==
false)
251 grads.resize(this->n_tensor_pols);
254 if (grad_grads.empty() ==
false)
256 grad_grads.resize(this->n_tensor_pols);
257 do_grad_grads =
true;
259 if (third_derivatives.empty() ==
false)
261 third_derivatives.resize(this->n_tensor_pols);
262 do_3rd_derivatives =
true;
264 if (fourth_derivatives.empty() ==
false)
266 fourth_derivatives.resize(this->n_tensor_pols);
267 do_4th_derivatives =
true;
272 for (
unsigned int i=this->n_tensor_pols; i<this->n_tensor_pols+n_bubbles; ++i)
275 values.push_back(compute_value(i,p));
277 grads.push_back(compute_grad(i,p));
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));
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)
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()