48 get_QGaussLobatto_points(
const unsigned int degree)
53 return std::vector<Point<1>>(1,
Point<1>(0.5));
61 template <
int dim,
int spacedim>
76 std::vector<
bool>(1, true)))
93 template <
int dim,
int spacedim>
100 polynomials.size() - 1,
103 polynomials.size() - 1),
105 polynomials.size() - 1)
111 polynomials.size() - 1)
113 std::vector<
bool>(1, true)))
125 template <
int dim,
int spacedim>
133 std::ostringstream namebuf;
135 << this->degree <<
")";
136 return namebuf.str();
141 template <
int dim,
int spacedim>
145 std::vector<double> & nodal_values)
const
148 this->get_unit_support_points().size());
152 for (
unsigned int i = 0; i < this->dofs_per_cell; ++i)
156 nodal_values[i] = support_point_values[i](0);
162 template <
int dim,
int spacedim>
163 std::unique_ptr<FiniteElement<dim, spacedim>>
166 return std_cxx14::make_unique<FE_DGQ<dim, spacedim>>(*this);
175 template <
int dim,
int spacedim>
176 std::vector<unsigned int>
179 std::vector<unsigned int> dpo(dim + 1, 0
U);
181 for (
unsigned int i = 1; i < dim; ++i)
188 template <
int dim,
int spacedim>
191 const char direction)
const
193 const unsigned int n = this->degree + 1;
195 for (
unsigned int i = 1; i < dim; ++i)
204 for (
unsigned int i = n; i > 0;)
214 for (
unsigned int iz = 0; iz < ((dim > 2) ? n : 1); ++iz)
215 for (
unsigned int j = 0; j < n; ++j)
216 for (
unsigned int i = 0; i < n; ++i)
218 unsigned int k = n * i - j + n - 1 + n * n * iz;
225 for (
unsigned int iz = 0; iz < ((dim > 2) ? n : 1); ++iz)
226 for (
unsigned int iy = 0; iy < n; ++iy)
227 for (
unsigned int ix = 0; ix < n; ++ix)
229 unsigned int k = n * ix - iy + n - 1 + n * n * iz;
237 for (
unsigned int iz = 0; iz < n; ++iz)
238 for (
unsigned int iy = 0; iy < n; ++iy)
239 for (
unsigned int ix = 0; ix < n; ++ix)
241 unsigned int k = n * (n * iy - iz + n - 1) + ix;
249 for (
unsigned int iz = 0; iz < n; ++iz)
250 for (
unsigned int iy = 0; iy < n; ++iy)
251 for (
unsigned int ix = 0; ix < n; ++ix)
253 unsigned int k = n * (n * iy - iz + n - 1) + ix;
265 template <
int dim,
int spacedim>
277 typename FE::ExcInterpolationNotImplemented());
284 Assert(interpolation_matrix.
m() == this->dofs_per_cell,
296 this->dofs_per_cell);
300 for (
unsigned int j = 0; j < this->dofs_per_cell; ++j)
305 const Point<dim> p = this->unit_support_points[j];
306 for (
unsigned int i = 0; i < this->dofs_per_cell; ++i)
307 cell_interpolation(j, i) = this->poly_space.compute_value(i, p);
318 cell_interpolation.
mmult(interpolation_matrix, source_interpolation);
321 for (
unsigned int i = 0; i < this->dofs_per_cell; ++i)
323 if (
std::fabs(interpolation_matrix(i, j)) < 1
e-15)
324 interpolation_matrix(i, j) = 0.;
331 for (
unsigned int i = 0; i < this->dofs_per_cell; ++i)
335 sum += interpolation_matrix(i, j);
344 template <
int dim,
int spacedim>
356 (void)interpolation_matrix;
360 typename FE::ExcInterpolationNotImplemented());
362 Assert(interpolation_matrix.
m() == 0,
364 Assert(interpolation_matrix.
n() == 0,
370 template <
int dim,
int spacedim>
383 (void)interpolation_matrix;
387 typename FE::ExcInterpolationNotImplemented());
389 Assert(interpolation_matrix.
m() == 0,
391 Assert(interpolation_matrix.
n() == 0,
397 template <
int dim,
int spacedim>
400 const unsigned int child,
407 "Prolongation matrices are only available for refined cells!"));
411 if (this->prolongation[refinement_case - 1][child].n() == 0)
413 std::lock_guard<std::mutex> lock(this->mutex);
416 if (this->prolongation[refinement_case - 1][child].n() ==
418 return this->prolongation[refinement_case - 1][child];
426 std::vector<std::vector<FullMatrix<double>>> isotropic_matrices(
428 isotropic_matrices.back().resize(
440 isotropic_matrices.back());
467 return this->prolongation[refinement_case - 1][child];
472 template <
int dim,
int spacedim>
475 const unsigned int child,
482 "Restriction matrices are only available for refined cells!"));
486 if (this->restriction[refinement_case - 1][child].n() == 0)
488 std::lock_guard<std::mutex> lock(this->mutex);
491 if (this->restriction[refinement_case - 1][child].n() ==
493 return this->restriction[refinement_case - 1][child];
501 std::vector<std::vector<FullMatrix<double>>> isotropic_matrices(
503 isotropic_matrices.back().resize(
514 this_nonconst.
restriction[refinement_case - 1].swap(
515 isotropic_matrices.back());
542 return this->restriction[refinement_case - 1][child];
547 template <
int dim,
int spacedim>
556 template <
int dim,
int spacedim>
557 std::vector<std::pair<unsigned int, unsigned int>>
565 return std::vector<std::pair<unsigned int, unsigned int>>();
570 template <
int dim,
int spacedim>
571 std::vector<std::pair<unsigned int, unsigned int>>
579 return std::vector<std::pair<unsigned int, unsigned int>>();
584 template <
int dim,
int spacedim>
585 std::vector<std::pair<unsigned int, unsigned int>>
593 return std::vector<std::pair<unsigned int, unsigned int>>();
598 template <
int dim,
int spacedim>
602 const unsigned int codim)
const
622 if (this->degree < fe_dgq_other->degree)
624 else if (this->degree == fe_dgq_other->degree)
632 if (this->degree < fe_q_other->degree)
634 else if (this->degree == fe_q_other->degree)
642 if (this->degree < fe_bernstein_other->degree)
644 else if (this->degree == fe_bernstein_other->degree)
652 if (this->degree < fe_bubbles_other->degree)
654 else if (this->degree == fe_bubbles_other->degree)
662 if (this->degree < fe_dg0_other->degree)
664 else if (this->degree == fe_dg0_other->degree)
672 if (this->degree < fe_q_iso_q1_other->degree)
674 else if (this->degree == fe_q_iso_q1_other->degree)
682 if (this->degree < fe_hierarchical_other->degree)
684 else if (this->degree == fe_hierarchical_other->degree)
692 if (fe_nothing->is_dominating())
707 template <
int dim,
int spacedim>
710 const unsigned int face_index)
const
715 unsigned int n = this->degree + 1;
720 if (this->unit_support_points.empty())
726 bool support_points_on_boundary =
true;
727 for (
unsigned int d = 0;
d < dim; ++
d)
728 if (std::abs(this->unit_support_points[0][
d]) > 1
e-13)
729 support_points_on_boundary =
false;
730 for (
unsigned int d = 0;
d < dim; ++
d)
731 if (std::abs(this->unit_support_points.back()[
d] - 1.) > 1
e-13)
732 support_points_on_boundary =
false;
733 if (support_points_on_boundary ==
false)
736 unsigned int n2 = n * n;
748 return (((shape_index == 0) && (face_index == 0)) ||
749 ((shape_index == this->degree) && (face_index == 1)));
754 if (face_index == 0 && (shape_index % n) == 0)
756 if (face_index == 1 && (shape_index % n) == this->degree)
758 if (face_index == 2 && shape_index < n)
760 if (face_index == 3 && shape_index >= this->dofs_per_cell - n)
767 const unsigned int in2 = shape_index % n2;
770 if (face_index == 0 && (shape_index % n) == 0)
773 if (face_index == 1 && (shape_index % n) == n - 1)
776 if (face_index == 2 && in2 < n)
779 if (face_index == 3 && in2 >= n2 - n)
782 if (face_index == 4 && shape_index < n2)
785 if (face_index == 5 && shape_index >= this->dofs_per_cell - n2)
798 template <
int dim,
int spacedim>
799 std::pair<Table<2, bool>, std::vector<unsigned int>>
803 constant_modes.
fill(
true);
804 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
805 constant_modes, std::vector<unsigned int>(1, 0));
812 template <
int dim,
int spacedim>
825 template <
int dim,
int spacedim>
831 std::ostringstream namebuf;
832 bool equidistant =
true;
833 std::vector<double> points(this->degree + 1);
835 std::vector<unsigned int> lexicographic =
836 this->poly_space.get_numbering_inverse();
837 for (
unsigned int j = 0; j <= this->degree; j++)
838 points[j] = this->unit_support_points[lexicographic[j]][0];
841 for (
unsigned int j = 0; j <= this->degree; j++)
842 if (std::abs(points[j] -
static_cast<double>(j) / this->degree) > 1
e-15)
847 if (this->degree == 0 && std::abs(points[0] - 0.5) < 1
e-15)
850 if (equidistant ==
true)
852 if (this->degree > 2)
853 namebuf <<
"FE_DGQArbitraryNodes<"
855 <<
">(QIterated(QTrapez()," << this->degree <<
"))";
858 << this->degree <<
")";
859 return namebuf.str();
864 bool gauss_lobatto =
true;
865 for (
unsigned int j = 0; j <= this->degree; j++)
866 if (points[j] != points_gl.
point(j)(0))
868 gauss_lobatto =
false;
872 if (gauss_lobatto ==
true)
875 << this->degree <<
")";
876 return namebuf.str();
880 const QGauss<1> points_g(this->degree + 1);
882 for (
unsigned int j = 0; j <= this->degree; j++)
883 if (points[j] != points_g.
point(j)(0))
892 <<
">(QGauss(" << this->degree + 1 <<
"))";
893 return namebuf.str();
898 bool gauss_log =
true;
899 for (
unsigned int j = 0; j <= this->degree; j++)
900 if (points[j] != points_glog.
point(j)(0))
906 if (gauss_log ==
true)
909 <<
">(QGaussLog(" << this->degree + 1 <<
"))";
910 return namebuf.str();
915 <<
">(QUnknownNodes(" << this->degree + 1 <<
"))";
916 return namebuf.str();
921 template <
int dim,
int spacedim>
926 std::vector<double> & nodal_values)
const
929 this->get_unit_support_points().size());
933 for (
unsigned int i = 0; i < this->dofs_per_cell; ++i)
937 nodal_values[i] = support_point_values[i](0);
943 template <
int dim,
int spacedim>
944 std::unique_ptr<FiniteElement<dim, spacedim>>
948 std::vector<Point<1>> qpoints(this->degree + 1);
949 std::vector<unsigned int> lexicographic =
950 this->poly_space.get_numbering_inverse();
951 for (
unsigned int i = 0; i <= this->degree; ++i)
952 qpoints[i] =
Point<1>(this->unit_support_points[lexicographic[i]][0]);
955 return std_cxx14::make_unique<FE_DGQArbitraryNodes<dim, spacedim>>(
963 template <
int dim,
int spacedim>
966 Polynomials::Legendre::generate_complete_basis(degree))
971 template <
int dim,
int spacedim>
972 std::pair<Table<2, bool>, std::vector<unsigned int>>
978 constant_modes(0, 0) =
true;
979 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
980 constant_modes, std::vector<unsigned int>(1, 0));
985 template <
int dim,
int spacedim>
995 template <
int dim,
int spacedim>
996 std::unique_ptr<FiniteElement<dim, spacedim>>
999 return std_cxx14::make_unique<FE_DGQLegendre<dim, spacedim>>(this->degree);
1006 template <
int dim,
int spacedim>
1009 Polynomials::HermiteLikeInterpolation::generate_complete_basis(degree))
1014 template <
int dim,
int spacedim>
1024 template <
int dim,
int spacedim>
1025 std::unique_ptr<FiniteElement<dim, spacedim>>
1028 return std_cxx14::make_unique<FE_DGQHermite<dim, spacedim>>(this->degree);
1034 #include "fe_dgq.inst"