34template <
int dim,
int spacedim>
43 get_riaf_vector(degree))
46 ExcMessage(
"This element can only be used for polynomial degrees "
47 "greater than zero"));
53 for (
unsigned int d = 0; d < dim; ++d)
61template <
int dim,
int spacedim>
65 Polynomials::generate_complete_Lagrange_basis(points.get_points())),
70 get_riaf_vector(points.size() - 1))
76 ExcMessage(
"This element can only be used for polynomial degrees "
83 for (
unsigned int d = 0; d < dim; ++d)
91template <
int dim,
int spacedim>
99 std::ostringstream namebuf;
101 const unsigned int n_points = this->degree + 1;
102 std::vector<double> points(n_points);
103 const unsigned int dofs_per_cell = this->n_dofs_per_cell();
104 const std::vector<Point<dim>> &unit_support_points =
105 this->unit_support_points;
106 unsigned int index = 0;
109 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
111 if ((dim > 1) ? (unit_support_points[j](1) == 0 &&
112 ((dim > 2) ? unit_support_points[j](2) == 0 :
true)) :
116 points[index] = unit_support_points[j](0);
118 points[n_points - 1] = unit_support_points[j](0);
120 points[index - 1] = unit_support_points[j](0);
126 Assert(index == n_points || (dim == 1 && index == n_points + 1),
128 "Could not decode support points in one coordinate direction."));
131 for (
unsigned int j = 0; j < n_points; ++j)
132 if (std::fabs(points[j] -
static_cast<double>(j) / this->degree) > 1e-15)
140 if (this->degree > 2)
142 <<
">(QIterated(QTrapezoid()," << this->degree <<
"))";
145 << this->degree <<
")";
152 for (
unsigned int j = 0; j < n_points; ++j)
153 if (points[j] != points_gl.
point(j)(0))
160 << this->degree <<
")";
163 <<
">(QUnknownNodes(" << this->degree <<
"))";
165 return namebuf.str();
170template <
int dim,
int spacedim>
171std::unique_ptr<FiniteElement<dim, spacedim>>
174 return std::make_unique<FE_Q_DG0<dim, spacedim>>(*this);
179template <
int dim,
int spacedim>
183 std::vector<double> & nodal_dofs)
const
185 Assert(support_point_values.size() == this->unit_support_points.size(),
187 this->unit_support_points.size()));
188 Assert(nodal_dofs.size() == this->n_dofs_per_cell(),
190 Assert(support_point_values[0].size() == this->n_components(),
192 this->n_components()));
194 for (
unsigned int i = 0; i < this->n_dofs_per_cell() - 1; ++i)
196 const std::pair<unsigned int, unsigned int> index =
197 this->system_to_component_index(i);
198 nodal_dofs[i] = support_point_values[i](index.first);
202 nodal_dofs.back() = 0.;
207template <
int dim,
int spacedim>
217 (x_source_fe.
get_name().find(
"FE_Q_DG0<") == 0) ||
218 (
dynamic_cast<const FEQDG0 *
>(&x_source_fe) !=
nullptr),
221 Assert(interpolation_matrix.
m() == this->n_dofs_per_cell(),
223 this->n_dofs_per_cell()));
229 x_source_fe, interpolation_matrix);
234template <
int dim,
int spacedim>
238 std::vector<bool> riaf(Utilities::fixed_power<dim>(deg + 1) + 1,
false);
245template <
int dim,
int spacedim>
246std::vector<unsigned int>
249 std::vector<unsigned int> dpo(dim + 1, 1U);
250 for (
unsigned int i = 1; i < dpo.size(); ++i)
251 dpo[i] = dpo[i - 1] * (deg - 1);
259template <
int dim,
int spacedim>
262 const unsigned int shape_index,
263 const unsigned int face_index)
const
266 if (shape_index == this->n_dofs_per_cell() - 1)
275template <
int dim,
int spacedim>
276std::pair<Table<2, bool>, std::vector<unsigned int>>
282 for (
unsigned int i = 0; i < this->n_dofs_per_cell() - 1; ++i)
283 constant_modes(0, i) =
true;
286 constant_modes(1, this->n_dofs_per_cell() - 1) =
true;
288 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
289 constant_modes, std::vector<unsigned int>(2, 0));
294template <
int dim,
int spacedim>
298 const unsigned int codim)
const
317 if (this->degree < fe_dg0_other->degree)
319 else if (this->degree == fe_dg0_other->degree)
327 if (fe_nothing->is_dominating())
342#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 Point< dim > & point(const unsigned int i) const
const std::vector< Point< dim > > & get_points() const
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#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)
@ either_element_can_dominate
@ other_element_dominates
@ neither_element_dominates
std::string dim_string(const int dim, const int spacedim)