36template <
int dim,
int spacedim>
38 const unsigned int degree)
52 , polynomial_space(
Polynomials::Legendre::generate_complete_basis(degree))
56 ref_case < RefinementCase<dim>::isotropic_refinement + 1;
65 const unsigned int nc =
67 for (
unsigned int i = 0; i < nc; ++i)
69 this->
prolongation[ref_case - 1][i].reinit(n_dofs, n_dofs);
72 for (
unsigned int j = 0; j < n_dofs; ++j)
108template <
int dim,
int spacedim>
119 std::ostringstream namebuf;
121 <<
">(" << this->degree <<
")";
123 return namebuf.str();
128template <
int dim,
int spacedim>
129std::unique_ptr<FiniteElement<dim, spacedim>>
132 return std::make_unique<FE_DGPNonparametric<dim, spacedim>>(*this);
137template <
int dim,
int spacedim>
152template <
int dim,
int spacedim>
155 const unsigned int i,
157 const unsigned int component)
const
171template <
int dim,
int spacedim>
185template <
int dim,
int spacedim>
188 const unsigned int i,
190 const unsigned int component)
const
204template <
int dim,
int spacedim>
219template <
int dim,
int spacedim>
222 const unsigned int i,
224 const unsigned int component)
const
242template <
int dim,
int spacedim>
243std::vector<unsigned int>
246 std::vector<unsigned int> dpo(dim + 1,
static_cast<unsigned int>(0));
248 for (
unsigned int i = 1; i < dim; ++i)
250 dpo[dim] *= deg + 1 + i;
258template <
int dim,
int spacedim>
277template <
int dim,
int spacedim>
278std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
289 std::make_unique<typename FiniteElement<dim, spacedim>::InternalDataBase>();
290 data_ptr->update_each = requires_update_flags(update_flags);
304template <
int dim,
int spacedim>
324 std::vector<double> values(
326 std::vector<Tensor<1, dim>> grads(
328 std::vector<Tensor<2, dim>> grad_grads(
330 std::vector<Tensor<3, dim>> empty_vector_of_3rd_order_tensors;
331 std::vector<Tensor<4, dim>> empty_vector_of_4th_order_tensors;
334 for (
unsigned int i = 0; i < n_q_points; ++i)
340 empty_vector_of_3rd_order_tensors,
341 empty_vector_of_4th_order_tensors);
344 for (
unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
345 output_data.shape_values[k][i] = values[k];
348 for (
unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
349 output_data.shape_gradients[k][i] = grads[k];
352 for (
unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
353 output_data.shape_hessians[k][i] = grad_grads[k];
359template <
int dim,
int spacedim>
379 std::vector<double> values(
381 std::vector<Tensor<1, dim>> grads(
383 std::vector<Tensor<2, dim>> grad_grads(
385 std::vector<Tensor<3, dim>> empty_vector_of_3rd_order_tensors;
386 std::vector<Tensor<4, dim>> empty_vector_of_4th_order_tensors;
389 for (
unsigned int i = 0; i < n_q_points; ++i)
395 empty_vector_of_3rd_order_tensors,
396 empty_vector_of_4th_order_tensors);
399 for (
unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
400 output_data.shape_values[k][i] = values[k];
403 for (
unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
404 output_data.shape_gradients[k][i] = grads[k];
407 for (
unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
408 output_data.shape_hessians[k][i] = grad_grads[k];
414template <
int dim,
int spacedim>
435 std::vector<double> values(
437 std::vector<Tensor<1, dim>> grads(
439 std::vector<Tensor<2, dim>> grad_grads(
441 std::vector<Tensor<3, dim>> empty_vector_of_3rd_order_tensors;
442 std::vector<Tensor<4, dim>> empty_vector_of_4th_order_tensors;
445 for (
unsigned int i = 0; i < n_q_points; ++i)
451 empty_vector_of_3rd_order_tensors,
452 empty_vector_of_4th_order_tensors);
455 for (
unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
456 output_data.shape_values[k][i] = values[k];
459 for (
unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
460 output_data.shape_gradients[k][i] = grads[k];
463 for (
unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
464 output_data.shape_hessians[k][i] = grad_grads[k];
470template <
int dim,
int spacedim>
475 const unsigned int)
const
483 (void)interpolation_matrix;
485 (x_source_fe.
get_name().find(
"FE_DGPNonparametric<") == 0) ||
490 Assert(interpolation_matrix.
m() == 0,
492 Assert(interpolation_matrix.
n() == 0,
498template <
int dim,
int spacedim>
504 const unsigned int)
const
512 (void)interpolation_matrix;
514 (x_source_fe.
get_name().find(
"FE_DGPNonparametric<") == 0) ||
519 Assert(interpolation_matrix.
m() == 0,
521 Assert(interpolation_matrix.
n() == 0,
527template <
int dim,
int spacedim>
536template <
int dim,
int spacedim>
537std::vector<std::pair<unsigned int, unsigned int>>
545 return std::vector<std::pair<unsigned int, unsigned int>>();
549 return std::vector<std::pair<unsigned int, unsigned int>>();
555template <
int dim,
int spacedim>
556std::vector<std::pair<unsigned int, unsigned int>>
564 return std::vector<std::pair<unsigned int, unsigned int>>();
568 return std::vector<std::pair<unsigned int, unsigned int>>();
574template <
int dim,
int spacedim>
575std::vector<std::pair<unsigned int, unsigned int>>
578 const unsigned int)
const
584 return std::vector<std::pair<unsigned int, unsigned int>>();
588 return std::vector<std::pair<unsigned int, unsigned int>>();
594template <
int dim,
int spacedim>
598 const unsigned int codim)
const
615 if (this->degree < fe_nonparametric_other->degree)
617 else if (this->degree == fe_nonparametric_other->degree)
625 if (fe_nothing->is_dominating())
640template <
int dim,
int spacedim>
644 const unsigned int)
const
651template <
int dim,
int spacedim>
661template <
int dim,
int spacedim>
671#include "fe_dgp_nonparametric.inst"
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 override
unsigned int get_degree() const
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const override
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
virtual bool hp_constraints_are_implemented() const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other, const unsigned int face_no=0) const override
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual std::string get_name() const override
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual void fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
friend class FE_DGPNonparametric
virtual std::size_t memory_consumption() const override
unsigned int n_dofs_per_cell() const
virtual std::string get_name() const =0
std::vector< std::vector< FullMatrix< double > > > prolongation
Abstract base class for mapping classes.
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define AssertThrow(cond, exc)
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
@ either_element_can_dominate
@ other_element_dominates
@ neither_element_dominates
std::string dim_string(const int dim, const int spacedim)
static unsigned int n_children(const RefinementCase< dim > &refinement_case)