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> 29 DEAL_II_NAMESPACE_OPEN
32 template <
int dim,
int spacedim>
40 get_riaf_vector(degree))
43 ExcMessage (
"This element can only be used for polynomial degrees " 44 "greater than zero"));
46 std::vector<Point<1> > support_points_1d(degree+1);
47 for (
unsigned int i=0; i<=
degree; ++i)
48 support_points_1d[i][0] = static_cast<double>(i)/
degree;
54 for (
unsigned int d=0; d<dim; ++d)
62 template <
int dim,
int spacedim>
76 ExcMessage (
"This element can only be used for polynomial degrees " 83 for (
unsigned int d=0; d<dim; ++d)
91 template <
int dim,
int spacedim>
99 std::ostringstream namebuf;
101 const unsigned int n_points = this->
degree +1;
102 std::vector<double> points(n_points);
105 unsigned int index = 0;
110 if ((dim>1) ? (unit_support_points[j](1)==0 &&
111 ((dim>2) ? unit_support_points[j](2)==0:
true)) :
true)
114 points[index] = unit_support_points[j](0);
116 points[n_points-1] = unit_support_points[j](0);
118 points[index-1] = unit_support_points[j](0);
124 Assert (index == n_points || (dim==1 && index == n_points+1),
125 ExcMessage (
"Could not decode support points in one coordinate direction."));
128 for (
unsigned int j=0; j<n_points; j++)
129 if (std::fabs(points[j] - (
double)j/this->
degree) > 1e-15)
136 namebuf <<
"FE_Q_DG0<" 138 <<
">(" << this->
degree <<
")";
145 for (
unsigned int j=0; j<n_points; j++)
146 if (points[j] != points_gl.
point(j)(0))
152 namebuf <<
"FE_Q_DG0<" 154 <<
">(QGaussLobatto(" << this->degree+1 <<
"))";
156 namebuf <<
"FE_Q_DG0<" 158 <<
">(QUnknownNodes(" << this->degree <<
"))";
160 return namebuf.str();
165 template <
int dim,
int spacedim>
174 template <
int dim,
int spacedim>
177 const std::vector<double> &values)
const 180 ExcDimensionMismatch(values.size(),
183 ExcDimensionMismatch(local_dofs.size(),this->
dofs_per_cell));
187 std::copy(values.begin(), values.end(), local_dofs.begin());
190 local_dofs[local_dofs.size()-1] = 0.;
195 template <
int dim,
int spacedim>
199 unsigned int offset)
const 202 ExcDimensionMismatch(values.size(),
205 ExcDimensionMismatch(local_dofs.size(),this->
dofs_per_cell));
207 ExcDimensionMismatch(values[0].size(),offset+this->
n_components()));
211 const std::pair<unsigned int, unsigned int> index
213 local_dofs[i] = values[i](offset+index.first);
217 local_dofs[local_dofs.size()-1] = 0.;
222 template <
int dim,
int spacedim>
225 std::vector<double> &local_dofs,
226 const VectorSlice<
const std::vector<std::vector<double> > > &values)
const 229 ExcDimensionMismatch(values.size(),
232 ExcDimensionMismatch(local_dofs.size(),this->
dofs_per_cell));
234 ExcDimensionMismatch(values.size(), this->
n_components()));
238 const std::pair<unsigned int, unsigned int> index
240 local_dofs[i] = values[index.first][i];
244 local_dofs[local_dofs.size()-1] = 0.;
249 template <
int dim,
int spacedim>
261 (dynamic_cast<const FEQDG0 *>(&x_source_fe) != 0),
263 ExcInterpolationNotImplemented());
266 ExcDimensionMismatch (interpolation_matrix.
m(),
269 ExcDimensionMismatch (interpolation_matrix.
m(),
278 template <
int dim,
int spacedim>
282 std::vector<bool> riaf(Utilities::fixed_power<dim>(deg+1)+1,
false);
283 riaf[riaf.size()-1]=
true;
289 template <
int dim,
int spacedim>
290 std::vector<unsigned int>
293 std::vector<unsigned int> dpo(dim+1, 1U);
294 for (
unsigned int i=1; i<dpo.size(); ++i)
295 dpo[i]=dpo[i-1]*(deg-1);
303 template <
int dim,
int spacedim>
306 const unsigned int face_index)
const 317 template <
int dim,
int spacedim>
318 std::pair<Table<2,bool>, std::vector<unsigned int> >
325 constant_modes(0, i) =
true;
328 constant_modes(1, this->dofs_per_cell-1) =
true;
330 return std::pair<Table<2,bool>, std::vector<unsigned int> >
331 (constant_modes, std::vector<unsigned int> (2, 0));
337 #include "fe_q_dg0.inst" 339 DEAL_II_NAMESPACE_CLOSE
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const
#define AssertDimension(dim1, dim2)
const std::vector< Point< dim > > & get_points() const
::ExceptionBase & ExcMessage(std::string arg1)
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
#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
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)
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
unsigned int n_components() const
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual FiniteElement< dim, spacedim > * clone() const
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)