49 const unsigned int start_index2d[6] = {0, 1, 4, 10, 20, 35};
50 const double points2d[35][2] = {
51 {0, 0}, {0, 0}, {1, 0}, {0, 1}, {0, 0},
52 {1, 0}, {0, 1}, {1, 1}, {0.5, 0}, {0, 0.5},
53 {0, 0}, {1, 0}, {0, 1}, {1, 1}, {1. / 3., 0},
54 {2. / 3., 0}, {0, 1. / 3.}, {0, 2. / 3.}, {0.5, 1}, {1, 0.5},
55 {0, 0}, {1, 0}, {0, 1}, {1, 1}, {0.25, 0},
56 {0.5, 0}, {0.75, 0}, {0, 0.25}, {0, 0.5}, {0, 0.75},
57 {1. / 3., 1}, {2. / 3., 1}, {1, 1. / 3.}, {1, 2. / 3.}, {0.5, 0.5}};
65 const unsigned int start_index3d[6] = {0, 1, 5, 15 };
66 const double points3d[35][3] = {{0, 0, 0},
85 generate_unit_points(
const unsigned int, std::vector<
Point<dim>> &);
89 generate_unit_points(
const unsigned int k, std::vector<
Point<1>> &p)
92 const double h = 1. / k;
93 for (
unsigned int i = 0; i < p.size(); ++i)
99 generate_unit_points(
const unsigned int k, std::vector<
Point<2>> &p)
102 Assert(p.size() == start_index2d[k + 1] - start_index2d[k],
104 for (
unsigned int i = 0; i < p.size(); ++i)
106 p[i](0) = points2d[start_index2d[k] + i][0];
107 p[i](1) = points2d[start_index2d[k] + i][1];
113 generate_unit_points(
const unsigned int k, std::vector<
Point<3>> &p)
116 Assert(p.size() == start_index3d[k + 1] - start_index3d[k],
118 for (
unsigned int i = 0; i < p.size(); ++i)
120 p[i](0) = points3d[start_index3d[k] + i][0];
121 p[i](1) = points3d[start_index3d[k] + i][1];
122 p[i](2) = points3d[start_index3d[k] + i][2];
176 std::ostringstream namebuf;
177 namebuf <<
"FE_DGPMonomial<" << dim <<
">(" << this->degree <<
")";
179 return namebuf.str();
185std::unique_ptr<FiniteElement<dim, dim>>
188 return std::make_unique<FE_DGPMonomial<dim>>(*this);
203 if (source_dgp_monomial)
208 const unsigned int m = interpolation_matrix.
m();
209 const unsigned int n = interpolation_matrix.
n();
212 Assert(m == this->n_dofs_per_cell(),
217 const unsigned int min_mn =
218 interpolation_matrix.
m() < interpolation_matrix.
n() ?
219 interpolation_matrix.
m() :
220 interpolation_matrix.
n();
222 for (
unsigned int i = 0; i < min_mn; ++i)
223 interpolation_matrix(i, i) = 1.;
227 std::vector<Point<dim>> unit_points(this->n_dofs_per_cell());
228 internal::FE_DGPMonomial::generate_unit_points(this->degree, unit_points);
233 for (
unsigned int k = 0; k < unit_points.size(); ++k)
234 source_fe_matrix(k, j) = source_fe.
shape_value(j, unit_points[k]);
237 this->n_dofs_per_cell());
238 for (
unsigned int j = 0; j < this->n_dofs_per_cell(); ++j)
239 for (
unsigned int k = 0; k < unit_points.size(); ++k)
241 this->poly_space->compute_value(j, unit_points[k]);
245 this_matrix.
mmult(interpolation_matrix, source_fe_matrix);
265std::vector<unsigned int>
268 std::vector<unsigned int> dpo(dim + 1, 0U);
270 for (
unsigned int i = 1; i < dim; ++i)
272 dpo[dim] *= deg + 1 + i;
284 const unsigned int)
const
292 (void)interpolation_matrix;
298 Assert(interpolation_matrix.
m() == 0,
300 Assert(interpolation_matrix.
n() == 0,
312 const unsigned int)
const
320 (void)interpolation_matrix;
326 Assert(interpolation_matrix.
m() == 0,
328 Assert(interpolation_matrix.
n() == 0,
344std::vector<std::pair<unsigned int, unsigned int>>
351 return std::vector<std::pair<unsigned int, unsigned int>>();
355 return std::vector<std::pair<unsigned int, unsigned int>>();
362std::vector<std::pair<unsigned int, unsigned int>>
369 return std::vector<std::pair<unsigned int, unsigned int>>();
373 return std::vector<std::pair<unsigned int, unsigned int>>();
380std::vector<std::pair<unsigned int, unsigned int>>
382 const unsigned int)
const
387 return std::vector<std::pair<unsigned int, unsigned int>>();
391 return std::vector<std::pair<unsigned int, unsigned int>>();
400 const unsigned int codim)
const
417 if (this->degree < fe_monomial_other->degree)
419 else if (this->degree == fe_monomial_other->degree)
427 if (fe_nothing->is_dominating())
445 const unsigned int face_index)
const
447 return face_index == 1 || (face_index == 0 && this->degree == 0);
455 const unsigned int face_index)
const
457 bool support_on_face =
false;
458 if (face_index == 1 || face_index == 2)
459 support_on_face =
true;
462 auto *
const polynomial_space_p =
465 const std::array<unsigned int, 2> degrees =
466 polynomial_space_p->directional_degrees(shape_index);
468 if ((face_index == 0 && degrees[1] == 0) ||
469 (face_index == 3 && degrees[0] == 0))
470 support_on_face =
true;
472 return support_on_face;
480 const unsigned int face_index)
const
482 bool support_on_face =
false;
483 if (face_index == 1 || face_index == 3 || face_index == 4)
484 support_on_face =
true;
487 auto *
const polynomial_space_p =
490 const std::array<unsigned int, 3> degrees =
491 polynomial_space_p->directional_degrees(shape_index);
493 if ((face_index == 0 && degrees[1] == 0) ||
494 (face_index == 2 && degrees[2] == 0) ||
495 (face_index == 5 && degrees[0] == 0))
496 support_on_face =
true;
498 return support_on_face;
514#include "fe_dgp_monomial.inst"
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim > &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_quad_dof_identities(const FiniteElement< dim > &fe_other, const unsigned int face_no=0) const override
FE_DGPMonomial(const unsigned int p)
virtual void get_subface_interpolation_matrix(const FiniteElement< dim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual void get_interpolation_matrix(const FiniteElement< dim > &source, FullMatrix< double > &matrix) const override
virtual void get_face_interpolation_matrix(const FiniteElement< dim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual std::size_t memory_consumption() const override
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim > &fe_other, const unsigned int codim=0) const override final
virtual bool hp_constraints_are_implemented() const override
virtual std::string get_name() const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim > &fe_other) const override
void initialize_restriction()
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
const std::unique_ptr< ScalarPolynomialsBase< dim > > poly_space
unsigned int n_dofs_per_cell() const
virtual std::string get_name() const =0
std::vector< std::vector< FullMatrix< double > > > restriction
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
std::vector< std::vector< FullMatrix< double > > > prolongation
void mmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define AssertThrow(cond, exc)
@ either_element_can_dominate
@ other_element_dominates
@ neither_element_dominates