17 #include <deal.II/base/geometry_info.h> 18 #include <deal.II/base/polynomials_rannacher_turek.h> 20 DEAL_II_NAMESPACE_OPEN
39 return (0.75 - 2.5 * p(0) + 1.5 * p(1) +
40 1.5 * (p(0) * p(0) - p(1) * p(1)));
44 return (-0.25 - 0.5 * p(0) + 1.5 * p(1) +
45 1.5 * (p(0) * p(0) - p(1) * p(1)));
49 return (0.75 + 1.5 * p(0) - 2.5 * p(1) -
50 1.5 * (p(0) * p(0) - p(1) * p(1)));
54 return (-0.25 + 1.5 * p(0) - 0.5 * p(1) -
55 1.5 * (p(0) * p(0) - p(1) * p(1)));
73 grad[0] = -2.5 + 3 * p(0);
74 grad[1] = 1.5 - 3 * p(1);
78 grad[0] = -0.5 + 3.0 * p(0);
79 grad[1] = 1.5 - 3.0 * p(1);
83 grad[0] = 1.5 - 3.0 * p(0);
84 grad[1] = -2.5 + 3.0 * p(1);
88 grad[0] = 1.5 - 3.0 * p(0);
89 grad[1] = -0.5 + 3.0 * p(1);
104 const unsigned int i,
114 grad_grad[1][1] = -3;
121 grad_grad[1][1] = -3;
125 grad_grad[0][0] = -3;
132 grad_grad[0][0] = -3;
146 std::vector<double> & values,
153 Assert(values.size() == n_pols || values.size() == 0,
155 Assert(grads.size() == n_pols || grads.size() == 0,
157 Assert(grad_grads.size() == n_pols || grad_grads.size() == 0,
159 Assert(third_derivatives.size() == n_pols || third_derivatives.size() == 0,
161 Assert(fourth_derivatives.size() == n_pols || fourth_derivatives.size() == 0,
164 for (
unsigned int i = 0; i < n_pols; ++i)
166 if (values.size() != 0)
168 values[i] = compute_value(i, unit_point);
170 if (grads.size() != 0)
172 grads[i] = compute_grad(i, unit_point);
174 if (grad_grads.size() != 0)
176 grad_grads[i] = compute_grad_grad(i, unit_point);
178 if (third_derivatives.size() != 0)
180 third_derivatives[i] = compute_derivative<3>(i, unit_point);
182 if (fourth_derivatives.size() != 0)
184 fourth_derivatives[i] = compute_derivative<4>(i, unit_point);
191 #include "polynomials_rannacher_turek.inst" 193 DEAL_II_NAMESPACE_CLOSE
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const
PolynomialsRannacherTurek()
double compute_value(const unsigned int i, const Point< dim > &p) const
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const
static ::ExceptionBase & ExcNotImplemented()
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