17 #include <deal.II/base/quadrature.h> 18 #include <deal.II/grid/tria.h> 19 #include <deal.II/grid/tria_iterator.h> 20 #include <deal.II/dofs/dof_accessor.h> 21 #include <deal.II/fe/fe.h> 22 #include <deal.II/fe/mapping.h> 23 #include <deal.II/fe/fe_dgp_nonparametric.h> 24 #include <deal.II/fe/fe_values.h> 27 #include <deal.II/base/std_cxx14/memory.h> 30 DEAL_II_NAMESPACE_OPEN
32 template <
int dim,
int spacedim>
42 std::vector<bool>(1,true))),
43 polynomial_space (
Polynomials::Legendre::generate_complete_basis(degree))
47 ref_case<RefinementCase<dim>::isotropic_refinement+1; ++ref_case)
56 for (
unsigned int i=0; i<nc; ++i)
58 this->
prolongation[ref_case-1][i].reinit (n_dofs, n_dofs);
61 for (
unsigned int j=0; j<n_dofs; ++j)
97 template <
int dim,
int spacedim>
108 std::ostringstream namebuf;
109 namebuf <<
"FE_DGPNonparametric<" 111 <<
">(" << this->degree <<
")";
113 return namebuf.str();
118 template <
int dim,
int spacedim>
119 std::unique_ptr<FiniteElement<dim,spacedim> >
122 return std_cxx14::make_unique<FE_DGPNonparametric<dim,spacedim>>(*this);
127 template <
int dim,
int spacedim>
141 template <
int dim,
int spacedim>
145 const unsigned int component)
const 158 template <
int dim,
int spacedim>
171 template <
int dim,
int spacedim>
175 const unsigned int component)
const 188 template <
int dim,
int spacedim>
202 template <
int dim,
int spacedim>
206 const unsigned int component)
const 223 template <
int dim,
int spacedim>
224 std::vector<unsigned int>
227 std::vector<unsigned int> dpo(dim+1, static_cast<unsigned int>(0));
229 for (
unsigned int i=1; i<dim; ++i)
239 template <
int dim,
int spacedim>
257 template <
int dim,
int spacedim>
258 std::unique_ptr<typename FiniteElement<dim,spacedim>::InternalDataBase>
266 auto data = std_cxx14::make_unique<typename FiniteElement<dim,spacedim>::InternalDataBase>();
267 data->update_each = requires_update_flags(update_flags);
272 return std::move(data);
281 template <
int dim,
int spacedim>
289 const ::internal::FEValuesImplementation::MappingRelatedData<dim, spacedim> &mapping_data,
295 const unsigned int n_q_points = mapping_data.quadrature_points.size();
300 std::vector<Tensor<3,dim> > empty_vector_of_3rd_order_tensors;
301 std::vector<Tensor<4,dim> > empty_vector_of_4th_order_tensors;
304 for (
unsigned int i=0; i<n_q_points; ++i)
306 polynomial_space.compute(mapping_data.quadrature_points[i],
307 values, grads, grad_grads,
308 empty_vector_of_3rd_order_tensors,
309 empty_vector_of_4th_order_tensors);
312 for (
unsigned int k=0; k<this->dofs_per_cell; ++k)
316 for (
unsigned int k=0; k<this->dofs_per_cell; ++k)
320 for (
unsigned int k=0; k<this->dofs_per_cell; ++k)
327 template <
int dim,
int spacedim>
335 const ::internal::FEValuesImplementation::MappingRelatedData<dim, spacedim> &mapping_data,
341 const unsigned int n_q_points = mapping_data.quadrature_points.size();
346 std::vector<Tensor<3,dim> > empty_vector_of_3rd_order_tensors;
347 std::vector<Tensor<4,dim> > empty_vector_of_4th_order_tensors;
350 for (
unsigned int i=0; i<n_q_points; ++i)
352 polynomial_space.compute(mapping_data.quadrature_points[i],
353 values, grads, grad_grads,
354 empty_vector_of_3rd_order_tensors,
355 empty_vector_of_4th_order_tensors);
358 for (
unsigned int k=0; k<this->dofs_per_cell; ++k)
362 for (
unsigned int k=0; k<this->dofs_per_cell; ++k)
366 for (
unsigned int k=0; k<this->dofs_per_cell; ++k)
373 template <
int dim,
int spacedim>
382 const ::internal::FEValuesImplementation::MappingRelatedData<dim, spacedim> &mapping_data,
388 const unsigned int n_q_points = mapping_data.quadrature_points.size();
393 std::vector<Tensor<3,dim> > empty_vector_of_3rd_order_tensors;
394 std::vector<Tensor<4,dim> > empty_vector_of_4th_order_tensors;
397 for (
unsigned int i=0; i<n_q_points; ++i)
399 polynomial_space.compute(mapping_data.quadrature_points[i],
400 values, grads, grad_grads,
401 empty_vector_of_3rd_order_tensors,
402 empty_vector_of_4th_order_tensors);
405 for (
unsigned int k=0; k<this->dofs_per_cell; ++k)
409 for (
unsigned int k=0; k<this->dofs_per_cell; ++k)
413 for (
unsigned int k=0; k<this->dofs_per_cell; ++k)
420 template <
int dim,
int spacedim>
432 (void)interpolation_matrix;
437 ExcInterpolationNotImplemented()));
439 Assert (interpolation_matrix.
m() == 0,
442 Assert (interpolation_matrix.
n() == 0,
449 template <
int dim,
int spacedim>
462 (void)interpolation_matrix;
468 Assert (interpolation_matrix.
m() == 0,
471 Assert (interpolation_matrix.
n() == 0,
478 template <
int dim,
int spacedim>
487 template <
int dim,
int spacedim>
488 std::vector<std::pair<unsigned int, unsigned int> >
496 std::vector<std::pair<unsigned int, unsigned int> > ();
500 return std::vector<std::pair<unsigned int, unsigned int> > ();
506 template <
int dim,
int spacedim>
507 std::vector<std::pair<unsigned int, unsigned int> >
515 std::vector<std::pair<unsigned int, unsigned int> > ();
519 return std::vector<std::pair<unsigned int, unsigned int> > ();
525 template <
int dim,
int spacedim>
526 std::vector<std::pair<unsigned int, unsigned int> >
534 std::vector<std::pair<unsigned int, unsigned int> > ();
538 return std::vector<std::pair<unsigned int, unsigned int> > ();
544 template <
int dim,
int spacedim>
561 template <
int dim,
int spacedim>
564 const unsigned int)
const 571 template <
int dim,
int spacedim>
581 template <
int dim,
int spacedim>
591 #include "fe_dgp_nonparametric.inst" 594 DEAL_II_NAMESPACE_CLOSE
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
friend class FE_DGPNonparametric
virtual std::size_t memory_consumption() const
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
Transformed quadrature points.
virtual bool hp_constraints_are_implemented() const
#define AssertThrow(cond, exc)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
unsigned int get_degree() const
std::vector< std::vector< FullMatrix< double > > > prolongation
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
virtual std::string get_name() const =0
Second derivatives of shape functions.
std::string dim_string(const int dim, const int spacedim)
const unsigned int dofs_per_cell
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
Shape function gradients.
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const
static ::ExceptionBase & ExcNotImplemented()
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
static ::ExceptionBase & ExcInternalError()
virtual std::string get_name() const