43 get_QGaussLobatto_points(
const unsigned int degree)
52 return std::vector<Point<1>>();
60template <
int dim,
int spacedim>
77template <
int dim,
int spacedim>
81 Polynomials::generate_complete_Lagrange_basis(points.get_points())),
93template <
int dim,
int spacedim>
101 std::ostringstream namebuf;
102 bool equidistant =
true;
103 std::vector<double> points(this->degree + 1);
108 std::vector<unsigned int> lexicographic =
110 for (
unsigned int j = 0; j <= this->degree; ++j)
111 points[j] = this->unit_support_points[lexicographic[j]][0];
114 for (
unsigned int j = 0; j <= this->degree; ++j)
115 if (std::fabs(points[j] -
static_cast<double>(j) / this->degree) > 1e-15)
121 if (equidistant ==
true)
123 if (this->degree > 2)
125 <<
">(QIterated(QTrapezoid()," << this->degree <<
"))";
128 << this->degree <<
")";
134 bool gauss_lobatto =
true;
135 for (
unsigned int j = 0; j <= this->degree; ++j)
136 if (points[j] != points_gl.point(j)[0])
138 gauss_lobatto =
false;
141 if (gauss_lobatto ==
true)
143 << this->degree <<
")";
146 <<
">(QUnknownNodes(" << this->degree <<
"))";
148 return namebuf.str();
153template <
int dim,
int spacedim>
157 std::vector<double> &nodal_values)
const
160 this->get_unit_support_points().size());
164 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
168 nodal_values[i] = support_point_values[i](0);
174template <
int dim,
int spacedim>
175std::unique_ptr<FiniteElement<dim, spacedim>>
178 return std::make_unique<FE_Q<dim, spacedim>>(*this);
183template <
int dim,
int spacedim>
187 const unsigned int codim)
const
206 if (this->degree < fe_q_other->degree)
208 else if (this->degree == fe_q_other->degree)
216 if (this->degree < fe_p_other->degree)
218 else if (this->degree == fe_p_other->degree)
226 if (this->degree < fe_wp_other->degree)
228 else if (this->degree == fe_wp_other->degree)
236 if (this->degree < fe_pp_other->degree)
238 else if (this->degree == fe_pp_other->degree)
246 if (fe_nothing->is_dominating())
257 if (this->degree == 1)
259 if (fe_hermite_other->degree > 1)
264 else if (this->degree >= fe_hermite_other->degree)
void initialize(const std::vector< Point< 1 > > &support_points_1d)
FE_Q(const unsigned int p)
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 std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
virtual std::string get_name() const override
const unsigned int degree
const std::vector< Point< dim > > & get_points() const
const std::vector< unsigned int > & get_numbering_inverse() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertDimension(dim1, dim2)
#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)