33template <
int dim,
int spacedim>
42 get_riaf_vector(degree))
45 ExcMessage(
"This element can only be used for polynomial degrees "
46 "greater than zero"));
52 for (
unsigned int d = 0; d < dim; ++d)
60template <
int dim,
int spacedim>
64 Polynomials::generate_complete_Lagrange_basis(points.get_points())),
69 get_riaf_vector(points.size() - 1))
75 ExcMessage(
"This element can only be used for polynomial degrees "
82 for (
unsigned int d = 0; d < dim; ++d)
90template <
int dim,
int spacedim>
98 std::ostringstream namebuf;
100 const unsigned int n_points = this->degree + 1;
101 std::vector<double> points(n_points);
102 const unsigned int dofs_per_cell = this->n_dofs_per_cell();
103 const std::vector<Point<dim>> &unit_support_points =
104 this->unit_support_points;
105 unsigned int index = 0;
108 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
110 if ((dim > 1) ? (unit_support_points[j][1] == 0 &&
111 ((dim > 2) ? unit_support_points[j][2] == 0 :
true)) :
115 points[index] = unit_support_points[j][0];
117 points[n_points - 1] = unit_support_points[j][0];
119 points[index - 1] = unit_support_points[j][0];
125 Assert(index == n_points || (dim == 1 && index == n_points + 1),
127 "Could not decode support points in one coordinate direction."));
130 for (
unsigned int j = 0; j < n_points; ++j)
131 if (std::fabs(points[j] -
static_cast<double>(j) / this->degree) > 1e-15)
139 if (this->degree > 2)
141 <<
">(QIterated(QTrapezoid()," << this->degree <<
"))";
144 << this->degree <<
")";
151 for (
unsigned int j = 0; j < n_points; ++j)
152 if (points[j] != points_gl.point(j)[0])
159 << this->degree <<
")";
162 <<
">(QUnknownNodes(" << this->degree <<
"))";
164 return namebuf.str();
169template <
int dim,
int spacedim>
170std::unique_ptr<FiniteElement<dim, spacedim>>
173 return std::make_unique<FE_Q_DG0<dim, spacedim>>(*this);
178template <
int dim,
int spacedim>
182 std::vector<double> &nodal_dofs)
const
184 Assert(support_point_values.size() == this->unit_support_points.size(),
186 this->unit_support_points.size()));
187 Assert(nodal_dofs.size() == this->n_dofs_per_cell(),
189 Assert(support_point_values[0].size() == this->n_components(),
191 this->n_components()));
193 for (
unsigned int i = 0; i < this->n_dofs_per_cell() - 1; ++i)
195 const std::pair<unsigned int, unsigned int> index =
196 this->system_to_component_index(i);
197 nodal_dofs[i] = support_point_values[i](index.first);
201 nodal_dofs.back() = 0.;
206template <
int dim,
int spacedim>
216 (x_source_fe.
get_name().find(
"FE_Q_DG0<") == 0) ||
217 (
dynamic_cast<const FEQDG0 *
>(&x_source_fe) !=
nullptr),
220 Assert(interpolation_matrix.
m() == this->n_dofs_per_cell(),
222 this->n_dofs_per_cell()));
228 x_source_fe, interpolation_matrix);
233template <
int dim,
int spacedim>
244template <
int dim,
int spacedim>
245std::vector<unsigned int>
248 std::vector<unsigned int> dpo(dim + 1, 1U);
249 for (
unsigned int i = 1; i < dpo.size(); ++i)
250 dpo[i] = dpo[i - 1] * (deg - 1);
258template <
int dim,
int spacedim>
261 const unsigned int shape_index,
262 const unsigned int face_index)
const
265 if (shape_index == this->n_dofs_per_cell() - 1)
274template <
int dim,
int spacedim>
275std::pair<Table<2, bool>, std::vector<unsigned int>>
281 for (
unsigned int i = 0; i < this->n_dofs_per_cell() - 1; ++i)
282 constant_modes(0, i) =
true;
285 constant_modes(1, this->n_dofs_per_cell() - 1) =
true;
287 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
288 constant_modes, std::vector<unsigned int>(2, 0));
293template <
int dim,
int spacedim>
297 const unsigned int codim)
const
316 if (this->degree < fe_dg0_other->degree)
318 else if (this->degree == fe_dg0_other->degree)
326 if (fe_nothing->is_dominating())
341#include "fe_q_dg0.inst"
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const override
void initialize(const std::vector< Point< 1 > > &support_points_1d)
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const override
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
FE_Q_DG0(const unsigned int p)
virtual std::string get_name() const override
static std::vector< bool > get_riaf_vector(const unsigned int degree)
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const override
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
const unsigned int degree
unsigned int n_dofs_per_cell() const
virtual std::string get_name() const =0
std::vector< Point< dim > > unit_support_points
const std::vector< Point< dim > > & get_points() const
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
#define DEAL_II_NOT_IMPLEMENTED()
@ either_element_can_dominate
@ other_element_dominates
@ neither_element_dominates
std::string dim_string(const int dim, const int spacedim)
constexpr T fixed_power(const T t)