35 template <
int dim,
int spacedim>
45 get_riaf_vector(degree))
48 ExcMessage(
"This element can only be used for polynomial degrees " 49 "greater than zero"));
55 for (
unsigned int d = 0;
d < dim; ++
d)
63 template <
int dim,
int spacedim>
78 ExcMessage(
"This element can only be used for polynomial degrees " 85 for (
unsigned int d = 0;
d < dim; ++
d)
93 template <
int dim,
int spacedim>
101 std::ostringstream namebuf;
103 const unsigned int n_points = this->
degree + 1;
104 std::vector<double> points(n_points);
108 unsigned int index = 0;
113 if ((dim > 1) ? (unit_support_points[j](1) == 0 &&
114 ((dim > 2) ? unit_support_points[j](2) == 0 :
true)) :
118 points[index] = unit_support_points[j](0);
120 points[n_points - 1] = unit_support_points[j](0);
122 points[index - 1] = unit_support_points[j](0);
128 Assert(index == n_points || (dim == 1 && index == n_points + 1),
130 "Could not decode support points in one coordinate direction."));
133 for (
unsigned int j = 0; j < n_points; j++)
144 <<
">(QIterated(QTrapez()," << this->
degree <<
"))";
147 << this->degree <<
")";
154 for (
unsigned int j = 0; j < n_points; j++)
155 if (points[j] != points_gl.
point(j)(0))
165 <<
">(QUnknownNodes(" << this->degree <<
"))";
167 return namebuf.str();
172 template <
int dim,
int spacedim>
173 std::unique_ptr<FiniteElement<dim, spacedim>>
176 return std_cxx14::make_unique<FE_Q_DG0<dim, spacedim>>(*this);
181 template <
int dim,
int spacedim>
185 std::vector<double> & nodal_dofs)
const 198 const std::pair<unsigned int, unsigned int> index =
200 nodal_dofs[i] = support_point_values[i](index.first);
204 nodal_dofs[nodal_dofs.size() - 1] = 0.;
209 template <
int dim,
int spacedim>
219 (x_source_fe.
get_name().find(
"FE_Q_DG0<") == 0) ||
220 (dynamic_cast<const FEQDG0 *>(&x_source_fe) !=
nullptr),
235 template <
int dim,
int spacedim>
239 std::vector<bool> riaf(Utilities::fixed_power<dim>(deg + 1) + 1,
false);
240 riaf[riaf.size() - 1] =
true;
246 template <
int dim,
int spacedim>
247 std::vector<unsigned int>
250 std::vector<unsigned int> dpo(dim + 1, 1
U);
251 for (
unsigned int i = 1; i < dpo.size(); ++i)
252 dpo[i] = dpo[i - 1] * (deg - 1);
260 template <
int dim,
int spacedim>
263 const unsigned int shape_index,
264 const unsigned int face_index)
const 276 template <
int dim,
int spacedim>
277 std::pair<Table<2, bool>, std::vector<unsigned int>>
284 constant_modes(0, i) =
true;
287 constant_modes(1, this->dofs_per_cell - 1) =
true;
289 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
290 constant_modes, std::vector<unsigned int>(2, 0));
295 template <
int dim,
int spacedim>
299 const unsigned int codim)
const 318 if (this->degree < fe_dg0_other->
degree)
320 else if (this->degree == fe_dg0_other->degree)
328 if (fe_nothing->is_dominating())
343 #include "fe_q_dg0.inst"
#define AssertDimension(dim1, dim2)
const std::vector< Point< dim > > & get_points() const
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
const unsigned int degree
const Point< dim > & point(const unsigned int i) const
#define AssertThrow(cond, exc)
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
std::vector< Point< dim > > unit_support_points
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define DEAL_II_NAMESPACE_CLOSE
virtual std::string get_name() const =0
Expression fabs(const Expression &x)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
std::string dim_string(const int dim, const int spacedim)
unsigned int size() const
const unsigned int dofs_per_cell
static std::vector< bool > get_riaf_vector(const unsigned int degree)
void initialize(const std::vector< Point< 1 >> &support_points_1d)
FE_Q_DG0(const unsigned int p)
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
unsigned int n_components() const
#define DEAL_II_NAMESPACE_OPEN
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
static ::ExceptionBase & ExcNotImplemented()
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::string get_name() 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
void L2(Vector< number > &result, const FEValuesBase< dim > &fe, const std::vector< double > &input, const double factor=1.)
std::vector< Polynomial< double > > generate_complete_Lagrange_basis(const std::vector< Point< 1 >> &points)
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)