31 std::vector<unsigned int>
34 std::vector<unsigned int> dpo(dim + 1);
61 std::vector<Point<dim>>
70 const Point<dim> centroid = reference_cell.template barycenter<dim>();
80 points.push_back(centroid);
87 for (
const auto &face_no : reference_cell.face_indices())
90 for (
const auto face_vertex_no :
91 reference_cell.face_reference_cell(0).vertex_indices())
93 const auto vertex_no =
94 reference_cell.face_to_cell_vertices(
100 reference_cell.template vertex<dim>(vertex_no);
104 reference_cell.face_reference_cell(0).n_vertices();
105 points.push_back(midpoint);
108 points.push_back(centroid);
121 std::vector<Point<0>>
124 std::vector<Point<0>> points;
125 points.emplace_back();
136 const Point<dim> centroid = reference_cell.template barycenter<dim>();
138 auto M = [](
const unsigned int d) {
158 for (
const auto &vertex : reference_cell.vertex_indices())
159 c_bubble = c_bubble * M(vertex);
160 c_bubble = c_bubble / c_bubble.
value(centroid);
162 std::vector<BarycentricPolynomial<dim>> bubble_functions;
165 bubble_functions.push_back(c_bubble);
172 for (
const auto &face_no : reference_cell.face_indices())
175 for (
const auto face_vertex_no :
176 reference_cell.face_reference_cell(0).vertex_indices())
177 vertices.push_back(reference_cell.face_to_cell_vertices(
185 bubble_functions.push_back(b -
186 b.value(centroid) * c_bubble);
189 bubble_functions.push_back(c_bubble);
194 const std::vector<Point<dim>> support_points =
196 const std::vector<Point<dim>> bubble_support_points(
197 support_points.begin() + fe_p.n(), support_points.end());
198 Assert(bubble_support_points.size() == bubble_functions.size(),
200 const unsigned int n_bubbles = bubble_support_points.size();
203 std::vector<BarycentricPolynomial<dim>> lump_polys;
204 for (
unsigned int i = 0; i < fe_p.n(); ++i)
208 for (
unsigned int j = 0; j < n_bubbles; ++j)
211 p.
value(bubble_support_points[j]) * bubble_functions[j];
214 lump_polys.push_back(p);
217 for (
auto &p : bubble_functions)
218 lump_polys.push_back(std::move(p));
223 for (
const auto &p : lump_polys)
227 for (
unsigned int d = 0; d < dim; ++d)
264template <
int dim,
int spacedim>
266 const unsigned int degree)
275 , approximation_degree(degree)
280template <
int dim,
int spacedim>
285 "(" + std::to_string(approximation_degree) +
")";
290template <
int dim,
int spacedim>
291std::unique_ptr<FiniteElement<dim, spacedim>>
294 return std::make_unique<FE_SimplexP_Bubbles<dim, spacedim>>(*this);
298#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
ReferenceCell reference_cell() const
const std::vector< Point< dim > > & get_unit_support_points() const
static constexpr unsigned char default_combined_face_orientation()
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
#define DEAL_II_NOT_IMPLEMENTED()
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)
constexpr const ReferenceCell & get_simplex()
std::string dim_string(const int dim, const int spacedim)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)