17 #include <deal.II/base/quadrature.h> 18 #include <deal.II/base/qprojector.h> 19 #include <deal.II/base/template_constraints.h> 20 #include <deal.II/fe/fe_q_dg0.h> 21 #include <deal.II/fe/fe_nothing.h> 22 #include <deal.II/base/quadrature_lib.h> 23 #include <deal.II/dofs/dof_accessor.h> 28 #include <deal.II/base/std_cxx14/memory.h> 30 DEAL_II_NAMESPACE_OPEN
33 template <
int dim,
int spacedim>
41 get_riaf_vector(degree))
44 ExcMessage (
"This element can only be used for polynomial degrees " 45 "greater than zero"));
51 for (
unsigned int d=0; d<dim; ++d)
59 template <
int dim,
int spacedim>
67 get_riaf_vector(points.size()-1))
73 ExcMessage (
"This element can only be used for polynomial degrees " 80 for (
unsigned int d=0; d<dim; ++d)
88 template <
int dim,
int spacedim>
96 std::ostringstream namebuf;
98 const unsigned int n_points = this->degree +1;
99 std::vector<double> points(n_points);
100 const unsigned int dofs_per_cell = this->dofs_per_cell;
101 const std::vector<Point<dim> > &unit_support_points = this->unit_support_points;
102 unsigned int index = 0;
105 for (
unsigned int j=0; j<dofs_per_cell; j++)
107 if ((dim>1) ? (unit_support_points[j](1)==0 &&
108 ((dim>2) ? unit_support_points[j](2)==0:
true)) :
true)
111 points[index] = unit_support_points[j](0);
113 points[n_points-1] = unit_support_points[j](0);
115 points[index-1] = unit_support_points[j](0);
121 Assert (index == n_points || (dim==1 && index == n_points+1),
122 ExcMessage (
"Could not decode support points in one coordinate direction."));
125 for (
unsigned int j=0; j<n_points; j++)
126 if (std::fabs(points[j] - (
double)j/this->degree) > 1e-15)
134 if (this->degree > 2)
135 namebuf <<
"FE_Q_DG0<" 137 <<
">(QIterated(QTrapez()," << this->degree <<
"))";
139 namebuf <<
"FE_Q_DG0<" 141 <<
">(" << this->degree <<
")";
149 for (
unsigned int j=0; j<n_points; j++)
150 if (points[j] != points_gl.
point(j)(0))
156 namebuf <<
"FE_Q_DG0<" 158 <<
">(" << this->degree <<
")";
160 namebuf <<
"FE_Q_DG0<" 162 <<
">(QUnknownNodes(" << this->degree <<
"))";
164 return namebuf.str();
169 template <
int dim,
int spacedim>
170 std::unique_ptr<FiniteElement<dim,spacedim> >
173 return std_cxx14::make_unique<FE_Q_DG0<dim,spacedim>>(*this);
178 template <
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->dofs_per_cell,
192 for (
unsigned int i=0; i<this->dofs_per_cell-1; ++i)
194 const std::pair<unsigned int, unsigned int> index
195 = this->system_to_component_index(i);
196 nodal_dofs[i] = support_point_values[i](index.first);
200 nodal_dofs[nodal_dofs.size()-1] = 0.;
205 template <
int dim,
int spacedim>
216 (dynamic_cast<const FEQDG0 *>(&x_source_fe) !=
nullptr),
218 ExcInterpolationNotImplemented()));
220 Assert (interpolation_matrix.
m() == this->dofs_per_cell,
222 this->dofs_per_cell));
228 get_interpolation_matrix(x_source_fe, interpolation_matrix);
233 template <
int dim,
int spacedim>
237 std::vector<bool> riaf(Utilities::fixed_power<dim>(deg+1)+1,
false);
238 riaf[riaf.size()-1]=
true;
244 template <
int dim,
int spacedim>
245 std::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);
258 template <
int dim,
int spacedim>
261 const unsigned int face_index)
const 264 if (shape_index == this->dofs_per_cell-1)
272 template <
int dim,
int spacedim>
273 std::pair<Table<2,bool>, std::vector<unsigned int> >
279 for (
unsigned int i=0; i<this->dofs_per_cell-1; ++i)
280 constant_modes(0, i) =
true;
283 constant_modes(1, this->dofs_per_cell-1) =
true;
285 return std::pair<Table<2,bool>, std::vector<unsigned int> >
286 (constant_modes, std::vector<unsigned int> (2, 0));
292 #include "fe_q_dg0.inst" 294 DEAL_II_NAMESPACE_CLOSE
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
#define AssertDimension(dim1, dim2)
const std::vector< Point< dim > > & get_points() const
const unsigned int degree
const Point< dim > & point(const unsigned int i) const
#define AssertThrow(cond, exc)
std::vector< Point< dim > > unit_support_points
virtual std::string get_name() const
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
#define Assert(cond, exc)
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
virtual std::string get_name() const =0
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)
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)