32 std::vector<unsigned int>
35 std::vector<unsigned int> dpo(dim + 1);
62 std::vector<Point<dim>>
81 points.push_back(centroid);
88 const double q13 = 1.0 / 3.0;
89 points.emplace_back(q13, q13, 0.0);
90 points.emplace_back(q13, 0.0, q13);
91 points.emplace_back(0.0, q13, q13);
92 points.emplace_back(q13, q13, q13);
93 points.push_back(centroid);
106 std::vector<Point<0>>
109 std::vector<Point<0>> points;
110 points.emplace_back();
123 auto M = [](
const unsigned int d) {
143 for (
unsigned int d = 0; d < dim + 1; ++d)
144 c_bubble = c_bubble * M(d);
145 c_bubble = c_bubble / c_bubble.
value(centroid);
147 std::vector<BarycentricPolynomial<dim>> bubble_functions;
150 bubble_functions.push_back(c_bubble);
157 auto b0 = 27 * M(0) * M(1) * M(2);
158 bubble_functions.push_back(b0 - b0.value(centroid) * c_bubble);
159 auto b1 = 27 * M(0) * M(1) * M(3);
160 bubble_functions.push_back(b1 - b1.value(centroid) * c_bubble);
161 auto b2 = 27 * M(0) * M(2) * M(3);
162 bubble_functions.push_back(b2 - b2.value(centroid) * c_bubble);
163 auto b3 = 27 * M(1) * M(2) * M(3);
164 bubble_functions.push_back(b3 - b3.value(centroid) * c_bubble);
166 bubble_functions.push_back(c_bubble);
171 const std::vector<Point<dim>> support_points =
172 unit_support_points<dim>(degree);
173 const std::vector<Point<dim>> bubble_support_points(
174 support_points.begin() + fe_p.n(), support_points.end());
175 Assert(bubble_support_points.size() == bubble_functions.size(),
177 const unsigned int n_bubbles = bubble_support_points.size();
180 std::vector<BarycentricPolynomial<dim>> lump_polys;
181 for (
unsigned int i = 0; i < fe_p.n(); ++i)
185 for (
unsigned int j = 0; j < n_bubbles; ++j)
188 p.
value(bubble_support_points[j]) * bubble_functions[j];
191 lump_polys.push_back(p);
194 for (
auto &p : bubble_functions)
195 lump_polys.push_back(std::move(p));
200 for (
const auto &p : lump_polys)
204 for (
unsigned int d = 0; d < dim; ++d)
230 const auto polys = get_basis<dim>(degree);
232 ReferenceCells::get_simplex<dim>(),
241template <
int dim,
int spacedim>
243 const unsigned int degree)
251 , approximation_degree(degree)
256template <
int dim,
int spacedim>
261 "(" + std::to_string(approximation_degree) +
")";
266template <
int dim,
int spacedim>
267std::unique_ptr<FiniteElement<dim, spacedim>>
270 return std::make_unique<FE_SimplexP_Bubbles<dim, spacedim>>(*this);
274#include "fe_simplex_p_bubbles.inst"
Number value(const Point< dim > &point) const
static BarycentricPolynomial< dim, Number > monomial(const unsigned int d)
static BarycentricPolynomials< dim > get_fe_p_basis(const unsigned int degree)
FE_SimplexP_Bubbles(const unsigned int degree)
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
virtual std::string get_name() const override
const unsigned int degree
const std::vector< Point< dim > > & get_unit_support_points() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
BarycentricPolynomials< dim > get_basis(const unsigned int degree)
std::vector< Point< dim > > unit_support_points(const unsigned int degree)
std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
std::vector< Point< 0 > > unit_support_points< 0 >(const unsigned int)
FiniteElementData< dim > get_fe_data(const unsigned int degree)
std::string dim_string(const int dim, const int spacedim)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)