56 get_embedding_computation_tolerance(
const unsigned int p)
63 return 1.e-15 * std::exp(std::pow(p, 1.075));
80 std::vector<
bool>(dim, true)))
107 deallog <<
"Face Embedding" << std::endl;
111 for (
unsigned int i = 0; i < GeometryInfo<dim>::max_children_per_face; ++i)
114 FETools::compute_face_embedding_matrices<dim, double>(
119 internal::FE_Nedelec::get_embedding_computation_tolerance(order));
134 for (
unsigned int i = 0; i < GeometryInfo<2>::max_children_per_face;
139 face_embeddings[i](j, k);
150 unsigned int target_row = 0;
152 for (
unsigned int i = 0; i < 2; ++i)
153 for (
unsigned int j = this->
degree; j < 2 * this->
degree;
157 face_embeddings[2 * i](j, k);
159 for (
unsigned int i = 0; i < 2; ++i)
160 for (
unsigned int j = 3 * this->degree;
161 j < GeometryInfo<3>::lines_per_face * this->
degree;
165 face_embeddings[i](j, k);
167 for (
unsigned int i = 0; i < 2; ++i)
168 for (
unsigned int j = 0; j < 2; ++j)
169 for (
unsigned int k = i * this->degree;
170 k < (i + 1) * this->degree;
174 face_embeddings[i + 2 * j](k,
l);
176 for (
unsigned int i = 0; i < 2; ++i)
177 for (
unsigned int j = 0; j < 2; ++j)
178 for (
unsigned int k = (i + 2) * this->
degree;
179 k < (i + 3) * this->degree;
183 face_embeddings[2 * i + j](k,
l);
185 for (
unsigned int i = 0; i < GeometryInfo<3>::max_children_per_face;
187 for (
unsigned int j =
193 face_embeddings[i](j, k);
216 std::ostringstream namebuf;
217 namebuf <<
"FE_Nedelec<" << dim <<
">(" << this->degree - 1 <<
")";
219 return namebuf.str();
224 std::unique_ptr<FiniteElement<dim, dim>>
227 return std_cxx14::make_unique<FE_Nedelec<dim>>(*this);
258 const std::vector<Polynomials::Polynomial<double>> &lobatto_polynomials =
260 std::vector<Polynomials::Polynomial<double>> lobatto_polynomials_grad(order +
263 for (
unsigned int i = 0; i < lobatto_polynomials_grad.size(); ++i)
264 lobatto_polynomials_grad[i] = lobatto_polynomials[i + 1].derivative();
268 const QGauss<dim - 1> reference_edge_quadrature(order + 1);
269 const unsigned int n_edge_points = reference_edge_quadrature.size();
270 const unsigned int n_boundary_points =
275 this->generalized_face_support_points.resize(n_edge_points);
278 for (
unsigned int q_point = 0; q_point < n_edge_points; ++q_point)
279 this->generalized_face_support_points[q_point] =
280 reference_edge_quadrature.point(q_point);
288 const unsigned int n_interior_points = quadrature.
size();
290 this->generalized_support_points.resize(n_boundary_points +
292 boundary_weights.reinit(n_edge_points, order);
294 for (
unsigned int q_point = 0; q_point < n_edge_points; ++q_point)
296 for (
unsigned int line = 0; line < GeometryInfo<dim>::lines_per_cell;
298 this->generalized_support_points[line * n_edge_points + q_point] =
300 line,
true,
false,
false, n_edge_points) +
303 for (
unsigned int i = 0; i < order; ++i)
304 boundary_weights(q_point, i) =
305 reference_edge_quadrature.weight(q_point) *
306 lobatto_polynomials_grad[i + 1].value(
307 this->generalized_face_support_points[q_point](0));
310 for (
unsigned int q_point = 0; q_point < n_interior_points; ++q_point)
311 this->generalized_support_points[q_point + n_boundary_points] =
312 quadrature.
point(q_point);
319 this->generalized_support_points.resize(n_boundary_points);
321 for (
unsigned int line = 0; line < GeometryInfo<dim>::lines_per_cell;
323 for (
unsigned int q_point = 0; q_point < n_edge_points; ++q_point)
324 this->generalized_support_points[line * n_edge_points + q_point] =
326 line,
true,
false,
false, n_edge_points) +
340 const std::vector<Polynomials::Polynomial<double>> &lobatto_polynomials =
342 std::vector<Polynomials::Polynomial<double>> lobatto_polynomials_grad(order +
345 for (
unsigned int i = 0; i < lobatto_polynomials_grad.size(); ++i)
346 lobatto_polynomials_grad[i] = lobatto_polynomials[i + 1].derivative();
350 const QGauss<1> reference_edge_quadrature(order + 1);
351 const unsigned int n_edge_points = reference_edge_quadrature.
size();
360 const QGauss<dim - 1> reference_face_quadrature(order + 1);
361 const unsigned int n_face_points = reference_face_quadrature.size();
362 const unsigned int n_boundary_points =
366 const unsigned int n_interior_points = quadrature.
size();
368 boundary_weights.reinit(n_edge_points + n_face_points,
369 2 * (order + 1) * order);
370 this->generalized_face_support_points.resize(4 * n_edge_points +
372 this->generalized_support_points.resize(n_boundary_points +
376 for (
unsigned int q_point = 0; q_point < n_edge_points; ++q_point)
378 for (
unsigned int line = 0;
381 this->generalized_face_support_points[line * n_edge_points +
383 edge_quadrature.point(
385 line,
true,
false,
false, n_edge_points) +
388 for (
unsigned int i = 0; i < 2; ++i)
389 for (
unsigned int j = 0; j < 2; ++j)
391 this->generalized_support_points[q_point +
392 (i + 4 * j) * n_edge_points] =
394 this->generalized_support_points[q_point + (i + 4 * j + 2) *
397 this->generalized_support_points[q_point + (i + 2 * (j + 4)) *
402 for (
unsigned int i = 0; i < order; ++i)
403 boundary_weights(q_point, i) =
404 reference_edge_quadrature.
weight(q_point) *
405 lobatto_polynomials_grad[i + 1].value(
406 this->generalized_face_support_points[q_point](1));
410 for (
unsigned int q_point = 0; q_point < n_face_points; ++q_point)
412 this->generalized_face_support_points[q_point + 4 * n_edge_points] =
413 reference_face_quadrature.point(q_point);
415 for (
unsigned int i = 0; i <= order; ++i)
416 for (
unsigned int j = 0; j < order; ++j)
418 boundary_weights(q_point + n_edge_points, 2 * (i * order + j)) =
419 reference_face_quadrature.weight(q_point) *
420 lobatto_polynomials_grad[i].value(
421 this->generalized_face_support_points[q_point +
424 lobatto_polynomials[j + 2].value(
425 this->generalized_face_support_points[q_point +
428 boundary_weights(q_point + n_edge_points,
429 2 * (i * order + j) + 1) =
430 reference_face_quadrature.weight(q_point) *
431 lobatto_polynomials_grad[i].value(
432 this->generalized_face_support_points[q_point +
435 lobatto_polynomials[j + 2].value(
436 this->generalized_face_support_points[q_point +
446 for (
unsigned int q_point = 0; q_point < n_face_points; ++q_point)
448 this->generalized_support_points[face * n_face_points + q_point +
452 face,
true,
false,
false, n_face_points) +
457 for (
unsigned int q_point = 0; q_point < n_interior_points; ++q_point)
458 this->generalized_support_points[q_point + n_boundary_points] =
459 quadrature.
point(q_point);
464 this->generalized_face_support_points.resize(4 * n_edge_points);
465 this->generalized_support_points.resize(
468 for (
unsigned int q_point = 0; q_point < n_edge_points; ++q_point)
470 for (
unsigned int line = 0;
473 this->generalized_face_support_points[line * n_edge_points +
475 edge_quadrature.point(
477 line,
true,
false,
false, n_edge_points) +
480 for (
unsigned int i = 0; i < 2; ++i)
481 for (
unsigned int j = 0; j < 2; ++j)
483 this->generalized_support_points[q_point +
484 (i + 4 * j) * n_edge_points] =
486 this->generalized_support_points[q_point + (i + 4 * j + 2) *
489 this->generalized_support_points[q_point + (i + 2 * (j + 4)) *
506 for (
unsigned int i = 0; i < GeometryInfo<1>::max_children_per_cell; ++i)
507 this->restriction[0][i].
reinit(0, 0);
523 const QGauss<1> edge_quadrature(2 * this->degree);
524 const std::vector<Point<1>> &edge_quadrature_points =
526 const unsigned int n_edge_quadrature_points = edge_quadrature.
size();
537 for (
unsigned int dof = 0; dof < this->dofs_per_cell; ++dof)
538 for (
unsigned int q_point = 0; q_point < n_edge_quadrature_points;
541 const double weight = 2.0 * edge_quadrature.
weight(q_point);
543 if (edge_quadrature_points[q_point](0) < 0.5)
546 0.0, 2.0 * edge_quadrature_points[q_point](0));
548 this->restriction[index][0](0, dof) +=
550 this->shape_value_component(dof, quadrature_point, 1);
551 quadrature_point(0) = 1.0;
552 this->restriction[index][1](this->degree, dof) +=
554 this->shape_value_component(dof, quadrature_point, 1);
555 quadrature_point(0) = quadrature_point(1);
556 quadrature_point(1) = 0.0;
557 this->restriction[index][0](2 * this->degree, dof) +=
559 this->shape_value_component(dof, quadrature_point, 0);
560 quadrature_point(1) = 1.0;
561 this->restriction[index][2](3 * this->degree, dof) +=
563 this->shape_value_component(dof, quadrature_point, 0);
569 0.0, 2.0 * edge_quadrature_points[q_point](0) - 1.0);
571 this->restriction[index][2](0, dof) +=
573 this->shape_value_component(dof, quadrature_point, 1);
574 quadrature_point(0) = 1.0;
575 this->restriction[index][3](this->degree, dof) +=
577 this->shape_value_component(dof, quadrature_point, 1);
578 quadrature_point(0) = quadrature_point(1);
579 quadrature_point(1) = 0.0;
580 this->restriction[index][1](2 * this->degree, dof) +=
582 this->shape_value_component(dof, quadrature_point, 0);
583 quadrature_point(1) = 1.0;
584 this->restriction[index][3](3 * this->degree, dof) +=
586 this->shape_value_component(dof, quadrature_point, 0);
594 if (this->degree > 1)
596 const unsigned int deg = this->degree - 1;
597 const std::vector<Polynomials::Polynomial<double>>
598 &legendre_polynomials =
604 n_edge_quadrature_points);
606 for (
unsigned int q_point = 0;
607 q_point < n_edge_quadrature_points;
610 const double weight =
611 std::sqrt(edge_quadrature.
weight(q_point));
613 for (
unsigned int i = 0; i < deg; ++i)
614 assembling_matrix(i, q_point) =
615 weight * legendre_polynomials[i + 1].value(
616 edge_quadrature_points[q_point](0));
621 assembling_matrix.
mTmult(system_matrix, assembling_matrix);
622 system_matrix_inv.
invert(system_matrix);
629 for (
unsigned int dof = 0; dof < this->dofs_per_cell; ++dof)
630 for (
unsigned int i = 0; i < 2; ++i)
634 for (
unsigned int q_point = 0;
635 q_point < n_edge_quadrature_points;
638 const double weight = edge_quadrature.
weight(q_point);
640 i, edge_quadrature_points[q_point](0));
642 edge_quadrature_points[q_point](0), i);
644 if (edge_quadrature_points[q_point](0) < 0.5)
647 i, 2.0 * edge_quadrature_points[q_point](0));
651 (2.0 * this->shape_value_component(
652 dof, quadrature_point_2, 1) -
653 this->restriction[index][i](i * this->degree,
655 this->shape_value_component(i * this->degree,
660 this->restriction[index][i + 2](i * this->degree,
662 this->shape_value_component(i * this->degree,
666 2.0 * edge_quadrature_points[q_point](0), i);
669 (2.0 * this->shape_value_component(
670 dof, quadrature_point_2, 0) -
671 this->restriction[index][2 * i]((i + 2) *
674 this->shape_value_component((i + 2) *
680 this->restriction[index][2 * i + 1](
681 (i + 2) * this->degree, dof) *
682 this->shape_value_component(
683 (i + 2) * this->degree, quadrature_point_1, 0);
690 this->restriction[index][i](i * this->degree,
692 this->shape_value_component(i * this->degree,
698 2.0 * edge_quadrature_points[q_point](0) - 1.0);
702 (2.0 * this->shape_value_component(
703 dof, quadrature_point_2, 1) -
704 this->restriction[index][i + 2](i * this->degree,
706 this->shape_value_component(i * this->degree,
711 this->restriction[index][2 * i]((i + 2) *
714 this->shape_value_component(
715 (i + 2) * this->degree, quadrature_point_1, 0);
717 2.0 * edge_quadrature_points[q_point](0) - 1.0,
721 (2.0 * this->shape_value_component(
722 dof, quadrature_point_2, 0) -
723 this->restriction[index][2 * i + 1](
724 (i + 2) * this->degree, dof) *
725 this->shape_value_component((i + 2) *
731 for (
unsigned int j = 0; j < this->degree - 1; ++j)
734 legendre_polynomials[j + 1].value(
735 edge_quadrature_points[q_point](0));
737 for (
unsigned int k = 0; k < tmp.
size(); ++k)
738 system_rhs(j, k) += tmp(k) * L_j;
742 system_matrix_inv.
mmult(solution, system_rhs);
744 for (
unsigned int j = 0; j < this->degree - 1; ++j)
745 for (
unsigned int k = 0; k < 2; ++k)
747 if (std::abs(solution(j, k)) > 1
e-14)
748 this->restriction[index][i + 2 * k](
749 i * this->degree + j + 1, dof) = solution(j, k);
751 if (std::abs(solution(j, k + 2)) > 1
e-14)
752 this->restriction[index][2 * i + k](
753 (i + 2) * this->degree + j + 1, dof) =
761 const std::vector<Polynomials::Polynomial<double>>
762 &lobatto_polynomials =
764 const unsigned int n_boundary_dofs =
766 const unsigned int n_quadrature_points = quadrature.
size();
771 n_quadrature_points);
773 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
776 const double weight = std::sqrt(quadrature.
weight(q_point));
778 for (
unsigned int i = 0; i < this->degree; ++i)
781 weight * legendre_polynomials[i].value(
784 for (
unsigned int j = 0; j < this->degree - 1; ++j)
785 assembling_matrix(i * (this->degree - 1) + j,
787 L_i * lobatto_polynomials[j + 2].value(
793 assembling_matrix.
m());
795 assembling_matrix.
mTmult(system_matrix, assembling_matrix);
796 system_matrix_inv.
reinit(system_matrix.
m(), system_matrix.
m());
797 system_matrix_inv.
invert(system_matrix);
800 solution.
reinit(system_matrix_inv.
m(), 8);
801 system_rhs.
reinit(system_matrix_inv.
m(), 8);
804 for (
unsigned int dof = 0; dof < this->dofs_per_cell; ++dof)
808 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
821 tmp(0) += 2.0 * this->shape_value_component(
822 dof, quadrature_point, 0);
823 tmp(1) += 2.0 * this->shape_value_component(
824 dof, quadrature_point, 1);
833 tmp(4) += 2.0 * this->shape_value_component(
834 dof, quadrature_point, 0);
835 tmp(5) += 2.0 * this->shape_value_component(
836 dof, quadrature_point, 1);
847 2.0 * this->shape_value_component(dof,
851 2.0 * this->shape_value_component(dof,
863 2.0 * this->shape_value_component(dof,
867 2.0 * this->shape_value_component(dof,
872 for (
unsigned int i = 0; i < 2; ++i)
873 for (
unsigned int j = 0; j < this->degree; ++j)
876 this->restriction[index][i](j + 2 * this->degree,
878 this->shape_value_component(
879 j + 2 * this->degree,
883 this->restriction[index][i](i * this->degree + j,
885 this->shape_value_component(
886 i * this->degree + j,
889 tmp(2 * (i + 2)) -= this->restriction[index][i + 2](
890 j + 3 * this->degree, dof) *
891 this->shape_value_component(
892 j + 3 * this->degree,
895 tmp(2 * i + 5) -= this->restriction[index][i + 2](
896 i * this->degree + j, dof) *
897 this->shape_value_component(
898 i * this->degree + j,
903 tmp *= quadrature.
weight(q_point);
905 for (
unsigned int i = 0; i < this->degree; ++i)
907 const double L_i_0 = legendre_polynomials[i].value(
909 const double L_i_1 = legendre_polynomials[i].value(
912 for (
unsigned int j = 0; j < this->degree - 1; ++j)
915 L_i_0 * lobatto_polynomials[j + 2].value(
918 L_i_1 * lobatto_polynomials[j + 2].value(
921 for (
unsigned int k = 0; k < 4; ++k)
923 system_rhs(i * (this->degree - 1) + j,
924 2 * k) += tmp(2 * k) * l_j_0;
925 system_rhs(i * (this->degree - 1) + j,
927 tmp(2 * k + 1) * l_j_1;
933 system_matrix_inv.
mmult(solution, system_rhs);
935 for (
unsigned int i = 0; i < this->degree; ++i)
936 for (
unsigned int j = 0; j < this->degree - 1; ++j)
937 for (
unsigned int k = 0; k < 4; ++k)
939 if (std::abs(solution(i * (this->degree - 1) + j,
941 this->restriction[index][k](i * (this->degree - 1) +
944 solution(i * (this->degree - 1) + j, 2 * k);
946 if (std::abs(solution(i * (this->degree - 1) + j,
948 this->restriction[index][k](
949 i + (this->degree - 1 + j) * this->degree +
952 solution(i * (this->degree - 1) + j, 2 * k + 1);
966 for (
unsigned int dof = 0; dof < this->dofs_per_cell; ++dof)
967 for (
unsigned int q_point = 0; q_point < n_edge_quadrature_points;
970 const double weight = 2.0 * edge_quadrature.
weight(q_point);
972 if (edge_quadrature_points[q_point](0) < 0.5)
973 for (
unsigned int i = 0; i < 2; ++i)
974 for (
unsigned int j = 0; j < 2; ++j)
977 i, 2.0 * edge_quadrature_points[q_point](0), j);
979 this->restriction[index][i + 4 * j]((i + 4 * j) *
983 this->shape_value_component(dof, quadrature_point, 1);
985 Point<dim>(2.0 * edge_quadrature_points[q_point](0),
988 this->restriction[index][2 * (i + 2 * j)](
989 (i + 4 * j + 2) * this->degree, dof) +=
991 this->shape_value_component(dof, quadrature_point, 0);
995 2.0 * edge_quadrature_points[q_point](0));
996 this->restriction[index][i + 2 * j]((i + 2 * (j + 4)) *
1000 this->shape_value_component(dof, quadrature_point, 2);
1004 for (
unsigned int i = 0; i < 2; ++i)
1005 for (
unsigned int j = 0; j < 2; ++j)
1008 i, 2.0 * edge_quadrature_points[q_point](0) - 1.0, j);
1010 this->restriction[index][i + 4 * j + 2]((i + 4 * j) *
1014 this->shape_value_component(dof, quadrature_point, 1);
1016 2.0 * edge_quadrature_points[q_point](0) - 1.0, i, j);
1017 this->restriction[index][2 * (i + 2 * j) + 1](
1018 (i + 4 * j + 2) * this->degree, dof) +=
1020 this->shape_value_component(dof, quadrature_point, 0);
1022 i, j, 2.0 * edge_quadrature_points[q_point](0) - 1.0);
1023 this->restriction[index][i + 2 * (j + 2)](
1024 (i + 2 * (j + 4)) * this->degree, dof) +=
1026 this->shape_value_component(dof, quadrature_point, 2);
1034 if (this->degree > 1)
1036 const unsigned int deg = this->degree - 1;
1037 const std::vector<Polynomials::Polynomial<double>>
1038 &legendre_polynomials =
1044 n_edge_quadrature_points);
1046 for (
unsigned int q_point = 0;
1047 q_point < n_edge_quadrature_points;
1050 const double weight =
1051 std::sqrt(edge_quadrature.
weight(q_point));
1053 for (
unsigned int i = 0; i < deg; ++i)
1054 assembling_matrix(i, q_point) =
1055 weight * legendre_polynomials[i + 1].value(
1056 edge_quadrature_points[q_point](0));
1061 assembling_matrix.
mTmult(system_matrix, assembling_matrix);
1062 system_matrix_inv.
invert(system_matrix);
1069 for (
unsigned int i = 0; i < 2; ++i)
1070 for (
unsigned int j = 0; j < 2; ++j)
1071 for (
unsigned int dof = 0; dof < this->dofs_per_cell; ++dof)
1075 for (
unsigned int q_point = 0;
1076 q_point < n_edge_quadrature_points;
1079 const double weight = edge_quadrature.
weight(q_point);
1081 i, edge_quadrature_points[q_point](0), j);
1083 edge_quadrature_points[q_point](0), i, j);
1085 i, j, edge_quadrature_points[q_point](0));
1087 if (edge_quadrature_points[q_point](0) < 0.5)
1090 i, 2.0 * edge_quadrature_points[q_point](0), j);
1093 weight * (2.0 * this->shape_value_component(
1094 dof, quadrature_point_3, 1) -
1095 this->restriction[index][i + 4 * j](
1096 (i + 4 * j) * this->degree, dof) *
1097 this->shape_value_component(
1098 (i + 4 * j) * this->degree,
1103 this->restriction[index][i + 4 * j + 2](
1104 (i + 4 * j) * this->degree, dof) *
1105 this->shape_value_component((i + 4 * j) *
1110 2.0 * edge_quadrature_points[q_point](0), i, j);
1113 (2.0 * this->shape_value_component(
1114 dof, quadrature_point_3, 0) -
1115 this->restriction[index][2 * (i + 2 * j)](
1116 (i + 4 * j + 2) * this->degree, dof) *
1117 this->shape_value_component(
1118 (i + 4 * j + 2) * this->degree,
1123 this->restriction[index][2 * (i + 2 * j) + 1](
1124 (i + 4 * j + 2) * this->degree, dof) *
1125 this->shape_value_component((i + 4 * j + 2) *
1130 i, j, 2.0 * edge_quadrature_points[q_point](0));
1133 (2.0 * this->shape_value_component(
1134 dof, quadrature_point_3, 2) -
1135 this->restriction[index][i + 2 * j](
1136 (i + 2 * (j + 4)) * this->degree, dof) *
1137 this->shape_value_component(
1138 (i + 2 * (j + 4)) * this->degree,
1143 this->restriction[index][i + 2 * (j + 2)](
1144 (i + 2 * (j + 4)) * this->degree, dof) *
1145 this->shape_value_component((i + 2 * (j + 4)) *
1155 this->restriction[index][i + 4 * j](
1156 (i + 4 * j) * this->degree, dof) *
1157 this->shape_value_component((i + 4 * j) *
1164 2.0 * edge_quadrature_points[q_point](0) - 1.0,
1168 (2.0 * this->shape_value_component(
1169 dof, quadrature_point_3, 1) -
1170 this->restriction[index][i + 4 * j + 2](
1171 (i + 4 * j) * this->degree, dof) *
1172 this->shape_value_component(
1173 (i + 4 * j) * this->degree,
1178 this->restriction[index][2 * (i + 2 * j)](
1179 (i + 4 * j + 2) * this->degree, dof) *
1180 this->shape_value_component((i + 4 * j + 2) *
1185 2.0 * edge_quadrature_points[q_point](0) - 1.0,
1190 (2.0 * this->shape_value_component(
1191 dof, quadrature_point_3, 0) -
1192 this->restriction[index][2 * (i + 2 * j) + 1](
1193 (i + 4 * j + 2) * this->degree, dof) *
1194 this->shape_value_component(
1195 (i + 4 * j + 2) * this->degree,
1200 this->restriction[index][i + 2 * j](
1201 (i + 2 * (j + 4)) * this->degree, dof) *
1202 this->shape_value_component((i + 2 * (j + 4)) *
1209 2.0 * edge_quadrature_points[q_point](0) - 1.0);
1212 (2.0 * this->shape_value_component(
1213 dof, quadrature_point_3, 2) -
1214 this->restriction[index][i + 2 * (j + 2)](
1215 (i + 2 * (j + 4)) * this->degree, dof) *
1216 this->shape_value_component(
1217 (i + 2 * (j + 4)) * this->degree,
1222 for (
unsigned int k = 0; k < deg; ++k)
1225 legendre_polynomials[k + 1].value(
1226 edge_quadrature_points[q_point](0));
1228 for (
unsigned int l = 0;
l < tmp.
size(); ++
l)
1229 system_rhs(k,
l) += tmp(
l) * L_k;
1233 system_matrix_inv.
mmult(solution, system_rhs);
1235 for (
unsigned int k = 0; k < 2; ++k)
1236 for (
unsigned int l = 0;
l < deg; ++
l)
1238 if (std::abs(solution(
l, k)) > 1
e-14)
1239 this->restriction[index][i + 2 * (2 * j + k)](
1240 (i + 4 * j) * this->degree +
l + 1, dof) =
1243 if (std::abs(solution(
l, k + 2)) > 1
e-14)
1244 this->restriction[index][2 * (i + 2 * j) + k](
1245 (i + 4 * j + 2) * this->degree +
l + 1, dof) =
1248 if (std::abs(solution(
l, k + 4)) > 1
e-14)
1249 this->restriction[index][i + 2 * (j + 2 * k)](
1250 (i + 2 * (j + 4)) * this->degree +
l + 1, dof) =
1255 const QGauss<2> face_quadrature(2 * this->degree);
1256 const std::vector<Point<2>> &face_quadrature_points =
1258 const std::vector<Polynomials::Polynomial<double>>
1259 &lobatto_polynomials =
1261 const unsigned int n_edge_dofs =
1263 const unsigned int n_face_quadrature_points =
1264 face_quadrature.
size();
1268 n_face_quadrature_points);
1270 for (
unsigned int q_point = 0;
1271 q_point < n_face_quadrature_points;
1274 const double weight =
1275 std::sqrt(face_quadrature.
weight(q_point));
1277 for (
unsigned int i = 0; i <= deg; ++i)
1280 weight * legendre_polynomials[i].value(
1281 face_quadrature_points[q_point](0));
1283 for (
unsigned int j = 0; j < deg; ++j)
1284 assembling_matrix(i * deg + j, q_point) =
1285 L_i * lobatto_polynomials[j + 2].value(
1286 face_quadrature_points[q_point](1));
1291 assembling_matrix.
m());
1293 assembling_matrix.
mTmult(system_matrix, assembling_matrix);
1294 system_matrix_inv.
reinit(system_matrix.
m(), system_matrix.
m());
1295 system_matrix_inv.
invert(system_matrix);
1298 solution.
reinit(system_matrix_inv.
m(), 24);
1299 system_rhs.
reinit(system_matrix_inv.
m(), 24);
1302 for (
unsigned int i = 0; i < 2; ++i)
1303 for (
unsigned int dof = 0; dof < this->dofs_per_cell; ++dof)
1307 for (
unsigned int q_point = 0;
1308 q_point < n_face_quadrature_points;
1313 if (face_quadrature_points[q_point](0) < 0.5)
1315 if (face_quadrature_points[q_point](1) < 0.5)
1319 2.0 * face_quadrature_points[q_point](0),
1320 2.0 * face_quadrature_points[q_point](1));
1322 tmp(0) += 2.0 * this->shape_value_component(
1323 dof, quadrature_point_0, 1);
1324 tmp(1) += 2.0 * this->shape_value_component(
1325 dof, quadrature_point_0, 2);
1327 2.0 * face_quadrature_points[q_point](0),
1329 2.0 * face_quadrature_points[q_point](1));
1330 tmp(8) += 2.0 * this->shape_value_component(
1331 dof, quadrature_point_0, 2);
1332 tmp(9) += 2.0 * this->shape_value_component(
1333 dof, quadrature_point_0, 0);
1335 2.0 * face_quadrature_points[q_point](0),
1336 2.0 * face_quadrature_points[q_point](1),
1338 tmp(16) += 2.0 * this->shape_value_component(
1339 dof, quadrature_point_0, 0);
1340 tmp(17) += 2.0 * this->shape_value_component(
1341 dof, quadrature_point_0, 1);
1348 2.0 * face_quadrature_points[q_point](0),
1349 2.0 * face_quadrature_points[q_point](1) -
1352 tmp(2) += 2.0 * this->shape_value_component(
1353 dof, quadrature_point_0, 1);
1354 tmp(3) += 2.0 * this->shape_value_component(
1355 dof, quadrature_point_0, 2);
1357 2.0 * face_quadrature_points[q_point](0),
1359 2.0 * face_quadrature_points[q_point](1) -
1361 tmp(10) += 2.0 * this->shape_value_component(
1362 dof, quadrature_point_0, 2);
1363 tmp(11) += 2.0 * this->shape_value_component(
1364 dof, quadrature_point_0, 0);
1366 2.0 * face_quadrature_points[q_point](0),
1367 2.0 * face_quadrature_points[q_point](1) -
1370 tmp(18) += 2.0 * this->shape_value_component(
1371 dof, quadrature_point_0, 0);
1372 tmp(19) += 2.0 * this->shape_value_component(
1373 dof, quadrature_point_0, 1);
1377 else if (face_quadrature_points[q_point](1) < 0.5)
1381 2.0 * face_quadrature_points[q_point](0) - 1.0,
1382 2.0 * face_quadrature_points[q_point](1));
1384 tmp(4) += 2.0 * this->shape_value_component(
1385 dof, quadrature_point_0, 1);
1386 tmp(5) += 2.0 * this->shape_value_component(
1387 dof, quadrature_point_0, 2);
1389 2.0 * face_quadrature_points[q_point](0) - 1.0,
1391 2.0 * face_quadrature_points[q_point](1));
1392 tmp(12) += 2.0 * this->shape_value_component(
1393 dof, quadrature_point_0, 2);
1394 tmp(13) += 2.0 * this->shape_value_component(
1395 dof, quadrature_point_0, 0);
1397 2.0 * face_quadrature_points[q_point](0) - 1.0,
1398 2.0 * face_quadrature_points[q_point](1),
1400 tmp(20) += 2.0 * this->shape_value_component(
1401 dof, quadrature_point_0, 0);
1402 tmp(21) += 2.0 * this->shape_value_component(
1403 dof, quadrature_point_0, 1);
1410 2.0 * face_quadrature_points[q_point](0) - 1.0,
1411 2.0 * face_quadrature_points[q_point](1) - 1.0);
1413 tmp(6) += 2.0 * this->shape_value_component(
1414 dof, quadrature_point_0, 1);
1415 tmp(7) += 2.0 * this->shape_value_component(
1416 dof, quadrature_point_0, 2);
1418 2.0 * face_quadrature_points[q_point](0) - 1.0,
1420 2.0 * face_quadrature_points[q_point](1) - 1.0);
1421 tmp(14) += 2.0 * this->shape_value_component(
1422 dof, quadrature_point_0, 2);
1423 tmp(15) += 2.0 * this->shape_value_component(
1424 dof, quadrature_point_0, 0);
1426 2.0 * face_quadrature_points[q_point](0) - 1.0,
1427 2.0 * face_quadrature_points[q_point](1) - 1.0,
1429 tmp(22) += 2.0 * this->shape_value_component(
1430 dof, quadrature_point_0, 0);
1431 tmp(23) += 2.0 * this->shape_value_component(
1432 dof, quadrature_point_0, 1);
1437 face_quadrature_points[q_point](0),
1438 face_quadrature_points[q_point](1));
1440 face_quadrature_points[q_point](0),
1442 face_quadrature_points[q_point](1));
1444 face_quadrature_points[q_point](0),
1445 face_quadrature_points[q_point](1),
1448 for (
unsigned int j = 0; j < 2; ++j)
1449 for (
unsigned int k = 0; k < 2; ++k)
1450 for (
unsigned int l = 0;
l <= deg; ++
l)
1452 tmp(2 * (j + 2 * k)) -=
1453 this->restriction[index][i + 2 * (2 * j + k)](
1454 (i + 4 * j) * this->degree +
l, dof) *
1455 this->shape_value_component(
1456 (i + 4 * j) * this->degree +
l,
1459 tmp(2 * (j + 2 * k) + 1) -=
1460 this->restriction[index][i + 2 * (2 * j + k)](
1461 (i + 2 * (k + 4)) * this->degree +
l, dof) *
1462 this->shape_value_component(
1463 (i + 2 * (k + 4)) * this->degree +
l,
1466 tmp(2 * (j + 2 * (k + 2))) -=
1467 this->restriction[index][2 * (i + 2 * j) + k](
1468 (2 * (i + 4) + k) * this->degree +
l, dof) *
1469 this->shape_value_component(
1470 (2 * (i + 4) + k) * this->degree +
l,
1473 tmp(2 * (j + 2 * k) + 9) -=
1474 this->restriction[index][2 * (i + 2 * j) + k](
1475 (i + 4 * j + 2) * this->degree +
l, dof) *
1476 this->shape_value_component(
1477 (i + 4 * j + 2) * this->degree +
l,
1480 tmp(2 * (j + 2 * (k + 4))) -=
1481 this->restriction[index][2 * (2 * i + j) + k](
1482 (4 * i + j + 2) * this->degree +
l, dof) *
1483 this->shape_value_component(
1484 (4 * i + j + 2) * this->degree +
l,
1487 tmp(2 * (j + 2 * k) + 17) -=
1488 this->restriction[index][2 * (2 * i + j) + k](
1489 (4 * i + k) * this->degree +
l, dof) *
1490 this->shape_value_component(
1491 (4 * i + k) * this->degree +
l,
1496 tmp *= face_quadrature.
weight(q_point);
1498 for (
unsigned int j = 0; j <= deg; ++j)
1500 const double L_j_0 = legendre_polynomials[j].value(
1501 face_quadrature_points[q_point](0));
1502 const double L_j_1 = legendre_polynomials[j].value(
1503 face_quadrature_points[q_point](1));
1505 for (
unsigned int k = 0; k < deg; ++k)
1507 const double l_k_0 =
1508 L_j_0 * lobatto_polynomials[k + 2].value(
1509 face_quadrature_points[q_point](1));
1510 const double l_k_1 =
1511 L_j_1 * lobatto_polynomials[k + 2].value(
1512 face_quadrature_points[q_point](0));
1514 for (
unsigned int l = 0;
l < 4; ++
l)
1516 system_rhs(j * deg + k, 2 *
l) +=
1518 system_rhs(j * deg + k, 2 *
l + 1) +=
1519 tmp(2 *
l + 1) * l_k_1;
1520 system_rhs(j * deg + k, 2 * (
l + 4)) +=
1521 tmp(2 * (
l + 4)) * l_k_1;
1522 system_rhs(j * deg + k, 2 *
l + 9) +=
1523 tmp(2 *
l + 9) * l_k_0;
1524 system_rhs(j * deg + k, 2 * (
l + 8)) +=
1525 tmp(2 * (
l + 8)) * l_k_0;
1526 system_rhs(j * deg + k, 2 *
l + 17) +=
1527 tmp(2 *
l + 17) * l_k_1;
1533 system_matrix_inv.
mmult(solution, system_rhs);
1535 for (
unsigned int j = 0; j < 2; ++j)
1536 for (
unsigned int k = 0; k < 2; ++k)
1537 for (
unsigned int l = 0;
l <= deg; ++
l)
1538 for (
unsigned int m = 0; m < deg; ++m)
1540 if (std::abs(solution(
l * deg + m,
1541 2 * (j + 2 * k))) > 1
e-14)
1542 this->restriction[index][i + 2 * (2 * j + k)](
1543 (2 * i * this->degree +
l) * deg + m +
1545 dof) = solution(
l * deg + m, 2 * (j + 2 * k));
1547 if (std::abs(solution(
l * deg + m,
1548 2 * (j + 2 * k) + 1)) >
1550 this->restriction[index][i + 2 * (2 * j + k)](
1551 ((2 * i + 1) * deg + m) * this->degree +
l +
1554 solution(
l * deg + m, 2 * (j + 2 * k) + 1);
1556 if (std::abs(solution(
l * deg + m,
1557 2 * (j + 2 * (k + 2)))) >
1559 this->restriction[index][2 * (i + 2 * j) + k](
1560 (2 * (i + 2) * this->degree +
l) * deg + m +
1563 solution(
l * deg + m, 2 * (j + 2 * (k + 2)));
1565 if (std::abs(solution(
l * deg + m,
1566 2 * (j + 2 * k) + 9)) >
1568 this->restriction[index][2 * (i + 2 * j) + k](
1569 ((2 * i + 5) * deg + m) * this->degree +
l +
1572 solution(
l * deg + m, 2 * (j + 2 * k) + 9);
1574 if (std::abs(solution(
l * deg + m,
1575 2 * (j + 2 * (k + 4)))) >
1577 this->restriction[index][2 * (2 * i + j) + k](
1578 (2 * (i + 4) * this->degree +
l) * deg + m +
1581 solution(
l * deg + m, 2 * (j + 2 * (k + 4)));
1583 if (std::abs(solution(
l * deg + m,
1584 2 * (j + 2 * k) + 17)) >
1586 this->restriction[index][2 * (2 * i + j) + k](
1587 ((2 * i + 9) * deg + m) * this->degree +
l +
1590 solution(
l * deg + m, 2 * (j + 2 * k) + 17);
1597 const unsigned int n_boundary_dofs =
1600 const unsigned int n_quadrature_points = quadrature.
size();
1604 n_quadrature_points);
1606 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
1609 const double weight = std::sqrt(quadrature.
weight(q_point));
1611 for (
unsigned int i = 0; i <= deg; ++i)
1614 weight * legendre_polynomials[i].value(
1617 for (
unsigned int j = 0; j < deg; ++j)
1620 L_i * lobatto_polynomials[j + 2].value(
1623 for (
unsigned int k = 0; k < deg; ++k)
1624 assembling_matrix((i * deg + j) * deg + k,
1626 l_j * lobatto_polynomials[k + 2].value(
1633 assembling_matrix.
m());
1635 assembling_matrix.
mTmult(system_matrix, assembling_matrix);
1636 system_matrix_inv.
reinit(system_matrix.
m(), system_matrix.
m());
1637 system_matrix_inv.
invert(system_matrix);
1640 solution.
reinit(system_matrix_inv.
m(), 24);
1641 system_rhs.
reinit(system_matrix_inv.
m(), 24);
1644 for (
unsigned int dof = 0; dof < this->dofs_per_cell; ++dof)
1648 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
1664 tmp(0) += 2.0 * this->shape_value_component(
1665 dof, quadrature_point, 0);
1666 tmp(1) += 2.0 * this->shape_value_component(
1667 dof, quadrature_point, 1);
1668 tmp(2) += 2.0 * this->shape_value_component(
1669 dof, quadrature_point, 2);
1679 tmp(3) += 2.0 * this->shape_value_component(
1680 dof, quadrature_point, 0);
1681 tmp(4) += 2.0 * this->shape_value_component(
1682 dof, quadrature_point, 1);
1683 tmp(5) += 2.0 * this->shape_value_component(
1684 dof, quadrature_point, 2);
1695 tmp(6) += 2.0 * this->shape_value_component(
1696 dof, quadrature_point, 0);
1697 tmp(7) += 2.0 * this->shape_value_component(
1698 dof, quadrature_point, 1);
1699 tmp(8) += 2.0 * this->shape_value_component(
1700 dof, quadrature_point, 2);
1710 tmp(9) += 2.0 * this->shape_value_component(
1711 dof, quadrature_point, 0);
1712 tmp(10) += 2.0 * this->shape_value_component(
1713 dof, quadrature_point, 1);
1714 tmp(11) += 2.0 * this->shape_value_component(
1715 dof, quadrature_point, 2);
1728 tmp(12) += 2.0 * this->shape_value_component(
1729 dof, quadrature_point, 0);
1730 tmp(13) += 2.0 * this->shape_value_component(
1731 dof, quadrature_point, 1);
1732 tmp(14) += 2.0 * this->shape_value_component(
1733 dof, quadrature_point, 2);
1743 tmp(15) += 2.0 * this->shape_value_component(
1744 dof, quadrature_point, 0);
1745 tmp(16) += 2.0 * this->shape_value_component(
1746 dof, quadrature_point, 1);
1747 tmp(17) += 2.0 * this->shape_value_component(
1748 dof, quadrature_point, 2);
1760 2.0 * this->shape_value_component(dof,
1764 2.0 * this->shape_value_component(dof,
1768 2.0 * this->shape_value_component(dof,
1781 2.0 * this->shape_value_component(dof,
1785 2.0 * this->shape_value_component(dof,
1789 2.0 * this->shape_value_component(dof,
1794 for (
unsigned int i = 0; i < 2; ++i)
1795 for (
unsigned int j = 0; j < 2; ++j)
1796 for (
unsigned int k = 0; k < 2; ++k)
1797 for (
unsigned int l = 0;
l <= deg; ++
l)
1799 tmp(3 * (i + 2 * (j + 2 * k))) -=
1800 this->restriction[index][2 * (2 * i + j) + k](
1801 (4 * i + j + 2) * this->degree +
l, dof) *
1802 this->shape_value_component(
1803 (4 * i + j + 2) * this->degree +
l,
1806 tmp(3 * (i + 2 * (j + 2 * k)) + 1) -=
1807 this->restriction[index][2 * (2 * i + j) + k](
1808 (4 * i + k) * this->degree +
l, dof) *
1809 this->shape_value_component(
1810 (4 * i + k) * this->degree +
l,
1813 tmp(3 * (i + 2 * (j + 2 * k)) + 2) -=
1814 this->restriction[index][2 * (2 * i + j) + k](
1815 (2 * (j + 4) + k) * this->degree +
l, dof) *
1816 this->shape_value_component(
1817 (2 * (j + 4) + k) * this->degree +
l,
1821 for (
unsigned int m = 0; m < deg; ++m)
1823 tmp(3 * (i + 2 * (j + 2 * k))) -=
1824 this->restriction[index][2 * (2 * i + j) +
1826 ((2 * j + 5) * deg + m) * this->degree +
1829 this->shape_value_component(
1830 ((2 * j + 5) * deg + m) * this->degree +
1834 tmp(3 * (i + 2 * (j + 2 * k))) -=
1835 this->restriction[index][2 * (2 * i + j) +
1837 (2 * (i + 4) * this->degree +
l) * deg +
1840 this->shape_value_component(
1841 (2 * (i + 4) * this->degree +
l) * deg +
1845 tmp(3 * (i + 2 * (j + 2 * k)) + 1) -=
1846 this->restriction[index][2 * (2 * i + j) +
1848 (2 * k * this->degree +
l) * deg + m +
1851 this->shape_value_component(
1852 (2 * k * this->degree +
l) * deg + m +
1856 tmp(3 * (i + 2 * (j + 2 * k)) + 1) -=
1857 this->restriction[index][2 * (2 * i + j) +
1859 ((2 * i + 9) * deg + m) * this->degree +
1862 this->shape_value_component(
1863 ((2 * i + 9) * deg + m) * this->degree +
1867 tmp(3 * (i + 2 * (j + 2 * k)) + 2) -=
1868 this->restriction[index][2 * (2 * i + j) +
1870 ((2 * k + 1) * deg + m) * this->degree +
1873 this->shape_value_component(
1874 ((2 * k + 1) * deg + m) * this->degree +
1878 tmp(3 * (i + 2 * (j + 2 * k)) + 2) -=
1879 this->restriction[index][2 * (2 * i + j) +
1881 (2 * (j + 2) * this->degree +
l) * deg +
1884 this->shape_value_component(
1885 (2 * (j + 2) * this->degree +
l) * deg +
1892 tmp *= quadrature.
weight(q_point);
1894 for (
unsigned int i = 0; i <= deg; ++i)
1896 const double L_i_0 = legendre_polynomials[i].value(
1898 const double L_i_1 = legendre_polynomials[i].value(
1900 const double L_i_2 = legendre_polynomials[i].value(
1903 for (
unsigned int j = 0; j < deg; ++j)
1905 const double l_j_0 =
1906 L_i_0 * lobatto_polynomials[j + 2].value(
1908 const double l_j_1 =
1909 L_i_1 * lobatto_polynomials[j + 2].value(
1911 const double l_j_2 =
1912 L_i_2 * lobatto_polynomials[j + 2].value(
1915 for (
unsigned int k = 0; k < deg; ++k)
1917 const double l_k_0 =
1918 l_j_0 * lobatto_polynomials[k + 2].value(
1920 const double l_k_1 =
1921 l_j_1 * lobatto_polynomials[k + 2].value(
1923 const double l_k_2 =
1924 l_j_2 * lobatto_polynomials[k + 2].value(
1927 for (
unsigned int l = 0;
l < 8; ++
l)
1929 system_rhs((i * deg + j) * deg + k,
1930 3 *
l) += tmp(3 *
l) * l_k_0;
1931 system_rhs((i * deg + j) * deg + k,
1933 tmp(3 *
l + 1) * l_k_1;
1934 system_rhs((i * deg + j) * deg + k,
1936 tmp(3 *
l + 2) * l_k_2;
1943 system_matrix_inv.
mmult(solution, system_rhs);
1945 for (
unsigned int i = 0; i < 2; ++i)
1946 for (
unsigned int j = 0; j < 2; ++j)
1947 for (
unsigned int k = 0; k < 2; ++k)
1948 for (
unsigned int l = 0;
l <= deg; ++
l)
1949 for (
unsigned int m = 0; m < deg; ++m)
1950 for (
unsigned int n = 0; n < deg; ++n)
1953 solution((
l * deg + m) * deg + n,
1954 3 * (i + 2 * (j + 2 * k)))) >
1956 this->restriction[index][2 * (2 * i + j) + k](
1957 (
l * deg + m) * deg + n + n_boundary_dofs,
1958 dof) = solution((
l * deg + m) * deg + n,
1959 3 * (i + 2 * (j + 2 * k)));
1962 solution((
l * deg + m) * deg + n,
1963 3 * (i + 2 * (j + 2 * k)) + 1)) >
1965 this->restriction[index][2 * (2 * i + j) + k](
1966 (
l + (m + deg) * this->degree) * deg + n +
1969 solution((
l * deg + m) * deg + n,
1970 3 * (i + 2 * (j + 2 * k)) + 1);
1973 solution((
l * deg + m) * deg + n,
1974 3 * (i + 2 * (j + 2 * k)) + 2)) >
1976 this->restriction[index][2 * (2 * i + j) + k](
1978 ((m + 2 * deg) * deg + n) * this->degree +
1981 solution((
l * deg + m) * deg + n,
1982 3 * (i + 2 * (j + 2 * k)) + 2);
1998 std::vector<unsigned int>
2001 std::vector<unsigned int> dpo(dim + 1);
2010 dpo[1] = degree + 1;
2011 dpo[2] = 2 * degree * (degree + 1);
2014 dpo[3] = 3 * degree * degree * (degree + 1);
2036 const unsigned int face_index)
const
2041 const unsigned int deg = this->degree - 1;
2048 if (!((shape_index > deg) && (shape_index < 2 * this->degree)))
2055 if ((shape_index > deg) &&
2064 if (shape_index < 3 * this->degree)
2071 if (!((shape_index >= 2 * this->degree) &&
2072 (shape_index < 3 * this->degree)))
2089 if (((shape_index > deg) && (shape_index < 2 * this->degree)) ||
2090 ((shape_index >= 5 * this->degree) &&
2091 (shape_index < 6 * this->degree)) ||
2092 ((shape_index >= 9 * this->degree) &&
2093 (shape_index < 10 * this->degree)) ||
2094 ((shape_index >= 11 * this->degree) &&
2120 if (((shape_index > deg) && (shape_index < 4 * this->degree)) ||
2121 ((shape_index >= 5 * this->degree) &&
2122 (shape_index < 8 * this->degree)) ||
2123 ((shape_index >= 9 * this->degree) &&
2124 (shape_index < 10 * this->degree)) ||
2125 ((shape_index >= 11 * this->degree) &&
2151 if ((shape_index < 3 * this->degree) ||
2152 ((shape_index >= 4 * this->degree) &&
2153 (shape_index < 7 * this->degree)) ||
2154 ((shape_index >= 8 * this->degree) &&
2155 (shape_index < 10 * this->degree)) ||
2179 if ((shape_index < 2 * this->degree) ||
2180 ((shape_index >= 3 * this->degree) &&
2181 (shape_index < 6 * this->degree)) ||
2182 ((shape_index >= 7 * this->degree) &&
2183 (shape_index < 8 * this->degree)) ||
2184 ((shape_index >= 10 * this->degree) &&
2210 if ((shape_index < 4 * this->degree) ||
2211 ((shape_index >= 8 * this->degree) &&
2232 if (((shape_index >= 4 * this->degree) &&
2275 const unsigned int codim)
const
2285 if (this->degree < fe_nedelec_other->degree)
2287 else if (this->degree == fe_nedelec_other->degree)
2295 if (fe_nothing->is_dominating())
2316 std::vector<std::pair<unsigned int, unsigned int>>
2321 return std::vector<std::pair<unsigned int, unsigned int>>();
2325 std::vector<std::pair<unsigned int, unsigned int>>
2340 std::vector<std::pair<unsigned int, unsigned int>> identities;
2342 for (
unsigned int i = 0;
2343 i <
std::min(fe_nedelec_other->degree, this->degree);
2345 identities.emplace_back(i, i);
2356 return std::vector<std::pair<unsigned int, unsigned int>>();
2362 return std::vector<std::pair<unsigned int, unsigned int>>();
2367 std::vector<std::pair<unsigned int, unsigned int>>
2382 const unsigned int p = fe_nedelec_other->degree;
2383 const unsigned int q = this->degree;
2384 const unsigned int p_min =
std::min(p, q);
2385 std::vector<std::pair<unsigned int, unsigned int>> identities;
2387 for (
unsigned int i = 0; i < p_min; ++i)
2388 for (
unsigned int j = 0; j < p_min - 1; ++j)
2390 identities.emplace_back(i * (q - 1) + j, i * (p - 1) + j);
2391 identities.emplace_back(i + (j + q - 1) * q, i + (j + p - 1) * p);
2403 return std::vector<std::pair<unsigned int, unsigned int>>();
2409 return std::vector<std::pair<unsigned int, unsigned int>>();
2434 Assert(interpolation_matrix.
n() == this->dofs_per_face,
2455 interpolation_matrix = 0;
2459 for (
unsigned int i = 0; i < this->degree; ++i)
2460 interpolation_matrix(i, i) = 1.0;
2470 const unsigned int p = source_fe.
degree;
2471 const unsigned int q = this->degree;
2473 for (
unsigned int i = 0; i < q; ++i)
2475 for (
unsigned int j = 1; j < GeometryInfo<dim>::lines_per_face; ++j)
2476 interpolation_matrix(j * p + i, j * q + i) = 1.0;
2478 for (
unsigned int j = 0; j < q - 1; ++j)
2483 i * (q - 1) + j) = 1.0;
2487 (j + q - 1) * q) = 1.0;
2526 const unsigned int subface,
2537 Assert(interpolation_matrix.
n() == this->dofs_per_face,
2558 interpolation_matrix = 0.0;
2562 const std::vector<Point<1>> &edge_quadrature_points =
2564 const unsigned int n_edge_quadrature_points = edge_quadrature.
size();
2570 for (
unsigned int dof = 0; dof < this->dofs_per_face; ++dof)
2571 for (
unsigned int q_point = 0; q_point < n_edge_quadrature_points;
2575 0.0, 0.5 * (edge_quadrature_points[q_point](0) + subface));
2577 interpolation_matrix(0, dof) +=
2578 0.5 * edge_quadrature.
weight(q_point) *
2579 this->shape_value_component(dof, quadrature_point, 1);
2582 if (source_fe.
degree > 1)
2584 const std::vector<Polynomials::Polynomial<double>>
2585 &legendre_polynomials =
2593 n_edge_quadrature_points);
2595 for (
unsigned int q_point = 0;
2596 q_point < n_edge_quadrature_points;
2599 const double weight =
2600 std::sqrt(edge_quadrature.
weight(q_point));
2602 for (
unsigned int i = 0; i < source_fe.
degree - 1; ++i)
2603 assembling_matrix(i, q_point) =
2604 weight * legendre_polynomials[i + 1].value(
2605 edge_quadrature_points[q_point](0));
2611 assembling_matrix.
mTmult(system_matrix, assembling_matrix);
2612 system_matrix_inv.
invert(system_matrix);
2618 for (
unsigned int dof = 0; dof < this->dofs_per_face; ++dof)
2622 for (
unsigned int q_point = 0;
2623 q_point < n_edge_quadrature_points;
2628 0.5 * (edge_quadrature_points[q_point](0) + subface));
2630 0.0, edge_quadrature_points[q_point](0));
2632 edge_quadrature.
weight(q_point) *
2633 (0.5 * this->shape_value_component(dof,
2636 interpolation_matrix(0, dof) *
2641 for (
unsigned int i = 0; i < source_fe.
degree - 1; ++i)
2643 tmp * legendre_polynomials[i + 1].value(
2644 edge_quadrature_points[q_point](0));
2647 system_matrix_inv.
vmult(solution, system_rhs);
2649 for (
unsigned int i = 0; i < source_fe.
degree - 1; ++i)
2650 if (std::abs(solution(i)) > 1
e-14)
2651 interpolation_matrix(i + 1, dof) = solution(i);
2660 const double shifts[4][2] = {{0.0, 0.0},
2665 for (
unsigned int dof = 0; dof < this->dofs_per_face; ++dof)
2666 for (
unsigned int q_point = 0; q_point < n_edge_quadrature_points;
2669 const double weight = 0.5 * edge_quadrature.
weight(q_point);
2671 for (
unsigned int i = 0; i < 2; ++i)
2674 0.5 * (i + shifts[subface][0]),
2675 0.5 * (edge_quadrature_points[q_point](0) +
2676 shifts[subface][1]),
2679 interpolation_matrix(i * source_fe.
degree, dof) +=
2681 this->shape_value_component(
2682 this->face_to_cell_index(dof, 4), quadrature_point, 1);
2684 Point<dim>(0.5 * (edge_quadrature_points[q_point](0) +
2685 shifts[subface][0]),
2686 0.5 * (i + shifts[subface][1]),
2688 interpolation_matrix((i + 2) * source_fe.
degree, dof) +=
2690 this->shape_value_component(
2691 this->face_to_cell_index(dof, 4), quadrature_point, 0);
2695 if (source_fe.
degree > 1)
2697 const std::vector<Polynomials::Polynomial<double>>
2698 &legendre_polynomials =
2706 n_edge_quadrature_points);
2708 for (
unsigned int q_point = 0;
2709 q_point < n_edge_quadrature_points;
2712 const double weight =
2713 std::sqrt(edge_quadrature.
weight(q_point));
2715 for (
unsigned int i = 0; i < source_fe.
degree - 1; ++i)
2716 assembling_matrix(i, q_point) =
2717 weight * legendre_polynomials[i + 1].value(
2718 edge_quadrature_points[q_point](0));
2724 assembling_matrix.
mTmult(system_matrix, assembling_matrix);
2725 system_matrix_inv.
invert(system_matrix);
2734 for (
unsigned int dof = 0; dof < this->dofs_per_face; ++dof)
2738 for (
unsigned int q_point = 0;
2739 q_point < n_edge_quadrature_points;
2742 const double weight = edge_quadrature.
weight(q_point);
2744 for (
unsigned int i = 0; i < 2; ++i)
2747 0.5 * (i + shifts[subface][0]),
2748 0.5 * (edge_quadrature_points[q_point](0) +
2749 shifts[subface][1]),
2752 i, edge_quadrature_points[q_point](0), 0.0);
2756 (0.5 * this->shape_value_component(
2757 this->face_to_cell_index(dof, 4),
2760 interpolation_matrix(i * source_fe.
degree, dof) *
2762 i * source_fe.
degree, quadrature_point_1, 1));
2763 quadrature_point_0 =
2765 (edge_quadrature_points[q_point](0) +
2766 shifts[subface][0]),
2767 0.5 * (i + shifts[subface][1]),
2769 quadrature_point_1 =
2770 Point<dim>(edge_quadrature_points[q_point](0),
2775 (0.5 * this->shape_value_component(
2776 this->face_to_cell_index(dof, 4),
2779 interpolation_matrix((i + 2) * source_fe.
degree,
2782 (i + 2) * source_fe.
degree,
2787 for (
unsigned int i = 0; i < source_fe.
degree - 1; ++i)
2789 const double L_i = legendre_polynomials[i + 1].value(
2790 edge_quadrature_points[q_point](0));
2792 for (
unsigned int j = 0;
2793 j < GeometryInfo<dim>::lines_per_face;
2795 system_rhs(i, j) += tmp(j) * L_i;
2799 system_matrix_inv.
mmult(solution, system_rhs);
2801 for (
unsigned int i = 0;
2802 i < GeometryInfo<dim>::lines_per_face;
2804 for (
unsigned int j = 0; j < source_fe.
degree - 1; ++j)
2805 if (std::abs(solution(j, i)) > 1
e-14)
2806 interpolation_matrix(i * source_fe.
degree + j + 1,
2807 dof) = solution(j, i);
2813 const std::vector<Polynomials::Polynomial<double>>
2814 &lobatto_polynomials =
2817 const unsigned int n_boundary_dofs =
2819 const unsigned int n_quadrature_points = quadrature.
size();
2824 n_quadrature_points);
2826 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
2829 const double weight = std::sqrt(quadrature.
weight(q_point));
2831 for (
unsigned int i = 0; i < source_fe.
degree; ++i)
2834 weight * legendre_polynomials[i].value(
2837 for (
unsigned int j = 0; j < source_fe.
degree - 1; ++j)
2838 assembling_matrix(i * (source_fe.
degree - 1) + j,
2840 L_i * lobatto_polynomials[j + 2].value(
2846 assembling_matrix.
m());
2848 assembling_matrix.
mTmult(system_matrix, assembling_matrix);
2849 system_matrix_inv.
reinit(system_matrix.
m(), system_matrix.
m());
2850 system_matrix_inv.
invert(system_matrix);
2853 solution.
reinit(system_matrix_inv.
m(), 2);
2854 system_rhs.
reinit(system_matrix_inv.
m(), 2);
2857 for (
unsigned int dof = 0; dof < this->dofs_per_face; ++dof)
2861 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
2870 tmp(0) = 0.5 * this->shape_value_component(
2871 this->face_to_cell_index(dof, 4),
2874 tmp(1) = 0.5 * this->shape_value_component(
2875 this->face_to_cell_index(dof, 4),
2883 for (
unsigned int i = 0; i < 2; ++i)
2884 for (
unsigned int j = 0; j < source_fe.
degree; ++j)
2886 tmp(0) -= interpolation_matrix(
2887 (i + 2) * source_fe.
degree + j, dof) *
2889 (i + 2) * source_fe.
degree + j,
2893 interpolation_matrix(i * source_fe.
degree + j,
2896 i * source_fe.
degree + j, quadrature_point, 1);
2899 tmp *= quadrature.
weight(q_point);
2901 for (
unsigned int i = 0; i < source_fe.
degree; ++i)
2903 const double L_i_0 = legendre_polynomials[i].value(
2905 const double L_i_1 = legendre_polynomials[i].value(
2908 for (
unsigned int j = 0; j < source_fe.
degree - 1;
2911 system_rhs(i * (source_fe.
degree - 1) + j, 0) +=
2913 lobatto_polynomials[j + 2].value(
2915 system_rhs(i * (source_fe.
degree - 1) + j, 1) +=
2917 lobatto_polynomials[j + 2].value(
2923 system_matrix_inv.
mmult(solution, system_rhs);
2925 for (
unsigned int i = 0; i < source_fe.
degree; ++i)
2926 for (
unsigned int j = 0; j < source_fe.
degree - 1; ++j)
2928 if (std::abs(solution(i * (source_fe.
degree - 1) + j,
2930 interpolation_matrix(i * (source_fe.
degree - 1) + j +
2933 solution(i * (source_fe.
degree - 1) + j, 0);
2935 if (std::abs(solution(i * (source_fe.
degree - 1) + j,
2937 interpolation_matrix(
2940 dof) = solution(i * (source_fe.
degree - 1) + j, 1);
2956 const unsigned int child,
2963 "Prolongation matrices are only available for refined cells!"));
2967 if (this->prolongation[refinement_case - 1][child].n() == 0)
2969 std::lock_guard<std::mutex> lock(this->mutex);
2972 if (this->prolongation[refinement_case - 1][child].n() ==
2973 this->dofs_per_cell)
2974 return this->prolongation[refinement_case - 1][child];
2985 #ifdef DEBUG_NEDELEC
2986 deallog <<
"Embedding" << std::endl;
2994 internal::FE_Nedelec::get_embedding_computation_tolerance(
2996 #ifdef DEBUG_NEDELEC
2997 deallog <<
"Restriction" << std::endl;
3005 return this->prolongation[refinement_case - 1][child];
3011 const unsigned int child,
3018 "Restriction matrices are only available for refined cells!"));
3023 if (this->restriction[refinement_case - 1][child].n() == 0)
3025 std::lock_guard<std::mutex> lock(this->mutex);
3028 if (this->restriction[refinement_case - 1][child].n() ==
3029 this->dofs_per_cell)
3030 return this->restriction[refinement_case - 1][child];
3041 #ifdef DEBUG_NEDELEC
3042 deallog <<
"Embedding" << std::endl;
3050 internal::FE_Nedelec::get_embedding_computation_tolerance(
3052 #ifdef DEBUG_NEDELEC
3053 deallog <<
"Restriction" << std::endl;
3061 return this->restriction[refinement_case - 1][child];
3075 std::vector<double> & nodal_values)
const
3077 const unsigned int deg = this->degree - 1;
3078 Assert(support_point_values.size() == this->generalized_support_points.size(),
3080 this->generalized_support_points.size()));
3084 Assert(nodal_values.size() == this->dofs_per_cell,
3086 std::fill(nodal_values.begin(), nodal_values.end(), 0.0);
3094 const QGauss<dim - 1> reference_edge_quadrature(this->degree);
3095 const unsigned int n_edge_points = reference_edge_quadrature.size();
3097 for (
unsigned int i = 0; i < 2; ++i)
3098 for (
unsigned int j = 0; j < 2; ++j)
3100 for (
unsigned int q_point = 0; q_point < n_edge_points;
3102 nodal_values[(i + 2 * j) * this->degree] +=
3103 reference_edge_quadrature.weight(q_point) *
3104 support_point_values[q_point + (i + 2 * j) * n_edge_points]
3110 if (std::abs(nodal_values[(i + 2 * j) * this->degree]) < 1
e-14)
3111 nodal_values[(i + 2 * j) * this->degree] = 0.0;
3124 if (this->degree - 1 > 1)
3129 const std::vector<Polynomials::Polynomial<double>>
3130 &lobatto_polynomials =
3134 std::vector<Polynomials::Polynomial<double>>
3135 lobatto_polynomials_grad(this->degree);
3137 for (
unsigned int i = 0; i < lobatto_polynomials_grad.size(); ++i)
3138 lobatto_polynomials_grad[i] =
3139 lobatto_polynomials[i + 1].derivative();
3144 for (
unsigned int i = 0; i < system_matrix.
m(); ++i)
3145 for (
unsigned int j = 0; j < system_matrix.
n(); ++j)
3146 for (
unsigned int q_point = 0; q_point < n_edge_points;
3148 system_matrix(i, j) +=
3149 boundary_weights(q_point, j) *
3150 lobatto_polynomials_grad[i + 1].value(
3151 this->generalized_face_support_points[q_point](0));
3156 system_matrix_inv.
invert(system_matrix);
3163 for (
unsigned int line = 0;
3164 line < GeometryInfo<dim>::lines_per_cell;
3170 for (
unsigned int q_point = 0; q_point < n_edge_points;
3174 support_point_values[line * n_edge_points + q_point]
3175 [line_coordinate[line]] -
3176 nodal_values[line * this->degree] *
3177 this->shape_value_component(
3178 line * this->degree,
3179 this->generalized_support_points[line *
3182 line_coordinate[line]);
3184 for (
unsigned int i = 0; i < system_rhs.
size(); ++i)
3185 system_rhs(i) += boundary_weights(q_point, i) * tmp;
3188 system_matrix_inv.
vmult(solution, system_rhs);
3194 for (
unsigned int i = 0; i < solution.
size(); ++i)
3195 if (std::abs(solution(i)) > 1
e-14)
3196 nodal_values[line * this->degree + i + 1] = solution(i);
3208 const QGauss<dim> reference_quadrature(this->degree);
3209 const unsigned int n_interior_points =
3210 reference_quadrature.
size();
3211 const std::vector<Polynomials::Polynomial<double>>
3212 &legendre_polynomials =
3216 system_matrix.
reinit((this->degree - 1) * this->degree,
3217 (this->degree - 1) * this->degree);
3220 for (
unsigned int i = 0; i < this->degree; ++i)
3221 for (
unsigned int j = 0; j < this->degree - 1; ++j)
3222 for (
unsigned int k = 0; k < this->degree; ++k)
3223 for (
unsigned int l = 0;
l < this->degree - 1; ++
l)
3224 for (
unsigned int q_point = 0;
3225 q_point < n_interior_points;
3227 system_matrix(i * (this->degree - 1) + j,
3228 k * (this->degree - 1) +
l) +=
3229 reference_quadrature.
weight(q_point) *
3230 legendre_polynomials[i].value(
3231 this->generalized_support_points
3233 n_edge_points](0)) *
3234 lobatto_polynomials[j + 2].value(
3235 this->generalized_support_points
3237 n_edge_points](1)) *
3238 lobatto_polynomials_grad[k].value(
3239 this->generalized_support_points
3241 n_edge_points](0)) *
3242 lobatto_polynomials[
l + 2].value(
3243 this->generalized_support_points
3247 system_matrix_inv.
reinit(system_matrix.
m(), system_matrix.
m());
3248 system_matrix_inv.
invert(system_matrix);
3252 system_rhs.
reinit(system_matrix_inv.
m());
3255 for (
unsigned int q_point = 0; q_point < n_interior_points;
3259 support_point_values[q_point +
3263 for (
unsigned int i = 0; i < 2; ++i)
3264 for (
unsigned int j = 0; j <= deg; ++j)
3265 tmp -= nodal_values[(i + 2) * this->degree + j] *
3266 this->shape_value_component(
3267 (i + 2) * this->degree + j,
3268 this->generalized_support_points
3273 for (
unsigned int i = 0; i <= deg; ++i)
3274 for (
unsigned int j = 0; j < deg; ++j)
3275 system_rhs(i * deg + j) +=
3276 reference_quadrature.
weight(q_point) * tmp *
3277 lobatto_polynomials_grad[i].value(
3278 this->generalized_support_points
3280 n_edge_points](0)) *
3281 lobatto_polynomials[j + 2].value(
3282 this->generalized_support_points
3287 solution.
reinit(system_matrix.
m());
3288 system_matrix_inv.
vmult(solution, system_rhs);
3294 for (
unsigned int i = 0; i <= deg; ++i)
3295 for (
unsigned int j = 0; j < deg; ++j)
3296 if (std::abs(solution(i * deg + j)) > 1
e-14)
3299 solution(i * deg + j);
3306 for (
unsigned int q_point = 0; q_point < n_interior_points;
3310 support_point_values[q_point +
3314 for (
unsigned int i = 0; i < 2; ++i)
3315 for (
unsigned int j = 0; j <= deg; ++j)
3316 tmp -= nodal_values[i * this->degree + j] *
3317 this->shape_value_component(
3318 i * this->degree + j,
3319 this->generalized_support_points
3324 for (
unsigned int i = 0; i <= deg; ++i)
3325 for (
unsigned int j = 0; j < deg; ++j)
3326 system_rhs(i * deg + j) +=
3327 reference_quadrature.
weight(q_point) * tmp *
3328 lobatto_polynomials_grad[i].value(
3329 this->generalized_support_points
3331 n_edge_points](1)) *
3332 lobatto_polynomials[j + 2].value(
3333 this->generalized_support_points
3338 system_matrix_inv.
vmult(solution, system_rhs);
3344 for (
unsigned int i = 0; i <= deg; ++i)
3345 for (
unsigned int j = 0; j < deg; ++j)
3346 if (std::abs(solution(i * deg + j)) > 1
e-14)
3349 this->degree] = solution(i * deg + j);
3359 const QGauss<1> reference_edge_quadrature(this->degree);
3360 const unsigned int n_edge_points = reference_edge_quadrature.
size();
3362 for (
unsigned int q_point = 0; q_point < n_edge_points; ++q_point)
3364 for (
unsigned int i = 0; i < 4; ++i)
3365 nodal_values[(i + 8) * this->degree] +=
3366 reference_edge_quadrature.
weight(q_point) *
3367 support_point_values[q_point + (i + 8) * n_edge_points][2];
3369 for (
unsigned int i = 0; i < 2; ++i)
3370 for (
unsigned int j = 0; j < 2; ++j)
3371 for (
unsigned int k = 0; k < 2; ++k)
3372 nodal_values[(i + 2 * (2 * j + k)) * this->degree] +=
3373 reference_edge_quadrature.
weight(q_point) *
3374 support_point_values[q_point + (i + 2 * (2 * j + k)) *
3375 n_edge_points][1 - k];
3382 for (
unsigned int i = 0; i < 4; ++i)
3383 if (std::abs(nodal_values[(i + 8) * this->degree]) < 1
e-14)
3384 nodal_values[(i + 8) * this->degree] = 0.0;
3386 for (
unsigned int i = 0; i < 2; ++i)
3387 for (
unsigned int j = 0; j < 2; ++j)
3388 for (
unsigned int k = 0; k < 2; ++k)
3390 nodal_values[(i + 2 * (2 * j + k)) * this->degree]) <
3392 nodal_values[(i + 2 * (2 * j + k)) * this->degree] = 0.0;
3403 if (this->degree > 1)
3408 const std::vector<Polynomials::Polynomial<double>>
3409 &lobatto_polynomials =
3413 std::vector<Polynomials::Polynomial<double>>
3414 lobatto_polynomials_grad(this->degree);
3416 for (
unsigned int i = 0; i < lobatto_polynomials_grad.size(); ++i)
3417 lobatto_polynomials_grad[i] =
3418 lobatto_polynomials[i + 1].derivative();
3423 for (
unsigned int i = 0; i < system_matrix.
m(); ++i)
3424 for (
unsigned int j = 0; j < system_matrix.
n(); ++j)
3425 for (
unsigned int q_point = 0; q_point < n_edge_points;
3427 system_matrix(i, j) +=
3428 boundary_weights(q_point, j) *
3429 lobatto_polynomials_grad[i + 1].value(
3430 this->generalized_face_support_points[q_point](1));
3435 system_matrix_inv.
invert(system_matrix);
3439 1, 1, 0, 0, 1, 1, 0, 0, 2, 2, 2, 2};
3443 for (
unsigned int line = 0;
3444 line < GeometryInfo<dim>::lines_per_cell;
3450 for (
unsigned int q_point = 0; q_point < this->degree;
3454 support_point_values[line * this->degree + q_point]
3455 [line_coordinate[line]] -
3456 nodal_values[line * this->degree] *
3457 this->shape_value_component(
3458 line * this->degree,
3460 ->generalized_support_points[line * this->degree +
3462 line_coordinate[line]);
3464 for (
unsigned int i = 0; i < system_rhs.
size(); ++i)
3465 system_rhs(i) += boundary_weights(q_point, i) * tmp;
3468 system_matrix_inv.
vmult(solution, system_rhs);
3474 for (
unsigned int i = 0; i < solution.
size(); ++i)
3475 if (std::abs(solution(i)) > 1
e-14)
3476 nodal_values[line * this->degree + i + 1] = solution(i);
3487 const std::vector<Polynomials::Polynomial<double>>
3488 &legendre_polynomials =
3491 const unsigned int n_face_points = n_edge_points * n_edge_points;
3493 system_matrix.
reinit((this->degree - 1) * this->degree,
3494 (this->degree - 1) * this->degree);
3497 for (
unsigned int i = 0; i < this->degree; ++i)
3498 for (
unsigned int j = 0; j < this->degree - 1; ++j)
3499 for (
unsigned int k = 0; k < this->degree; ++k)
3500 for (
unsigned int l = 0;
l < this->degree - 1; ++
l)
3501 for (
unsigned int q_point = 0; q_point < n_face_points;
3503 system_matrix(i * (this->degree - 1) + j,
3504 k * (this->degree - 1) +
l) +=
3505 boundary_weights(q_point + n_edge_points,
3506 2 * (k * (this->degree - 1) +
l)) *
3507 legendre_polynomials[i].value(
3508 this->generalized_face_support_points
3509 [q_point + 4 * n_edge_points](0)) *
3510 lobatto_polynomials[j + 2].value(
3511 this->generalized_face_support_points
3512 [q_point + 4 * n_edge_points](1));
3514 system_matrix_inv.
reinit(system_matrix.
m(), system_matrix.
m());
3515 system_matrix_inv.
invert(system_matrix);
3516 solution.
reinit(system_matrix.
m());
3517 system_rhs.
reinit(system_matrix.
m());
3521 {1, 2}, {1, 2}, {2, 0}, {2, 0}, {0, 1}, {0, 1}};
3538 for (
unsigned int q_point = 0; q_point < n_face_points;
3542 support_point_values[q_point +
3545 [face_coordinates[face][0]];
3547 for (
unsigned int i = 0; i < 2; ++i)
3548 for (
unsigned int j = 0; j <= deg; ++j)
3550 nodal_values[edge_indices[face][i] * this->degree +
3552 this->shape_value_component(
3553 edge_indices[face][i] * this->degree + j,
3554 this->generalized_support_points
3557 face_coordinates[face][0]);
3559 for (
unsigned int i = 0; i <= deg; ++i)
3560 for (
unsigned int j = 0; j < deg; ++j)
3561 system_rhs(i * deg + j) +=
3562 boundary_weights(q_point + n_edge_points,
3563 2 * (i * deg + j)) *
3567 system_matrix_inv.
vmult(solution, system_rhs);
3573 for (
unsigned int i = 0; i <= deg; ++i)
3574 for (
unsigned int j = 0; j < deg; ++j)
3575 if (std::abs(solution(i * deg + j)) > 1
e-14)
3576 nodal_values[(2 * face * this->degree + i +
3580 solution(i * deg + j);
3587 for (
unsigned int q_point = 0; q_point < n_face_points;
3591 support_point_values[q_point +
3594 [face_coordinates[face][1]];
3596 for (
unsigned int i = 2;
3597 i < GeometryInfo<dim>::lines_per_face;
3599 for (
unsigned int j = 0; j <= deg; ++j)
3601 nodal_values[edge_indices[face][i] * this->degree +
3603 this->shape_value_component(
3604 edge_indices[face][i] * this->degree + j,
3605 this->generalized_support_points
3608 face_coordinates[face][1]);
3610 for (
unsigned int i = 0; i <= deg; ++i)
3611 for (
unsigned int j = 0; j < deg; ++j)
3612 system_rhs(i * deg + j) +=
3613 boundary_weights(q_point + n_edge_points,
3614 2 * (i * deg + j) + 1) *
3618 system_matrix_inv.
vmult(solution, system_rhs);
3624 for (
unsigned int i = 0; i <= deg; ++i)
3625 for (
unsigned int j = 0; j < deg; ++j)
3626 if (std::abs(solution(i * deg + j)) > 1
e-14)
3627 nodal_values[((2 * face + 1) * deg + j +
3630 i] = solution(i * deg + j);
3638 const QGauss<dim> reference_quadrature(this->degree);
3639 const unsigned int n_interior_points =
3640 reference_quadrature.
size();
3644 system_matrix.
reinit(this->degree * deg * deg,
3645 this->degree * deg * deg);
3648 for (
unsigned int i = 0; i <= deg; ++i)
3649 for (
unsigned int j = 0; j < deg; ++j)
3650 for (
unsigned int k = 0; k < deg; ++k)
3651 for (
unsigned int l = 0;
l <= deg; ++
l)
3652 for (
unsigned int m = 0; m < deg; ++m)
3653 for (
unsigned int n = 0; n < deg; ++n)
3654 for (
unsigned int q_point = 0;
3655 q_point < n_interior_points;
3657 system_matrix((i * deg + j) * deg + k,
3658 (
l * deg + m) * deg + n) +=
3659 reference_quadrature.
weight(q_point) *
3660 legendre_polynomials[i].value(
3661 this->generalized_support_points
3666 n_face_points](0)) *
3667 lobatto_polynomials[j + 2].value(
3668 this->generalized_support_points
3673 n_face_points](1)) *
3674 lobatto_polynomials[k + 2].value(
3675 this->generalized_support_points
3680 n_face_points](2)) *
3681 lobatto_polynomials_grad[
l].value(
3682 this->generalized_support_points
3687 n_face_points](0)) *
3688 lobatto_polynomials[m + 2].value(
3689 this->generalized_support_points
3694 n_face_points](1)) *
3695 lobatto_polynomials[n + 2].value(
3696 this->generalized_support_points
3703 system_matrix_inv.
reinit(system_matrix.
m(), system_matrix.
m());
3704 system_matrix_inv.
invert(system_matrix);
3706 system_rhs.
reinit(system_matrix.
m());
3709 for (
unsigned int q_point = 0; q_point < n_interior_points;
3713 support_point_values[q_point +
3719 for (
unsigned int i = 0; i <= deg; ++i)
3721 for (
unsigned int j = 0; j < 2; ++j)
3722 for (
unsigned int k = 0; k < 2; ++k)
3724 nodal_values[i + (j + 4 * k + 2) * this->degree] *
3725 this->shape_value_component(
3726 i + (j + 4 * k + 2) * this->degree,
3727 this->generalized_support_points
3735 for (
unsigned int j = 0; j < deg; ++j)
3736 for (
unsigned int k = 0; k < 4; ++k)
3738 nodal_values[(i + 2 * (k + 2) * this->degree +
3743 this->shape_value_component(
3744 (i + 2 * (k + 2) * this->degree +
3748 this->generalized_support_points
3757 for (
unsigned int i = 0; i <= deg; ++i)
3758 for (
unsigned int j = 0; j < deg; ++j)
3759 for (
unsigned int k = 0; k < deg; ++k)
3760 system_rhs((i * deg + j) * deg + k) +=
3761 reference_quadrature.
weight(q_point) * tmp *
3762 lobatto_polynomials_grad[i].value(
3763 this->generalized_support_points
3768 n_face_points](0)) *
3769 lobatto_polynomials[j + 2].value(
3770 this->generalized_support_points
3775 n_face_points](1)) *
3776 lobatto_polynomials[k + 2].value(
3777 this->generalized_support_points
3786 system_matrix_inv.
vmult(solution, system_rhs);
3792 for (
unsigned int i = 0; i <= deg; ++i)
3793 for (
unsigned int j = 0; j < deg; ++j)
3794 for (
unsigned int k = 0; k < deg; ++k)
3795 if (std::abs(solution((i * deg + j) * deg + k)) > 1
e-14)
3802 solution((i * deg + j) * deg + k);
3807 for (
unsigned int q_point = 0; q_point < n_interior_points;
3811 support_point_values[q_point +
3817 for (
unsigned int i = 0; i <= deg; ++i)
3818 for (
unsigned int j = 0; j < 2; ++j)
3820 for (
unsigned int k = 0; k < 2; ++k)
3821 tmp -= nodal_values[i + (4 * j + k) * this->degree] *
3822 this->shape_value_component(
3823 i + (4 * j + k) * this->degree,
3824 this->generalized_support_points
3832 for (
unsigned int k = 0; k < deg; ++k)
3834 nodal_values[(i + 2 * j * this->degree +
3839 this->shape_value_component(
3840 (i + 2 * j * this->degree +
3844 this->generalized_support_points
3852 ((2 * j + 9) * deg + k +
3855 this->shape_value_component(
3856 i + ((2 * j + 9) * deg + k +
3859 this->generalized_support_points
3868 for (
unsigned int i = 0; i <= deg; ++i)
3869 for (
unsigned int j = 0; j < deg; ++j)
3870 for (
unsigned int k = 0; k < deg; ++k)
3871 system_rhs((i * deg + j) * deg + k) +=
3872 reference_quadrature.
weight(q_point) * tmp *
3873 lobatto_polynomials_grad[i].value(
3874 this->generalized_support_points
3879 n_face_points](1)) *
3880 lobatto_polynomials[j + 2].value(
3881 this->generalized_support_points
3886 n_face_points](0)) *
3887 lobatto_polynomials[k + 2].value(
3888 this->generalized_support_points
3896 system_matrix_inv.
vmult(solution, system_rhs);
3902 for (
unsigned int i = 0; i <= deg; ++i)
3903 for (
unsigned int j = 0; j < deg; ++j)
3904 for (
unsigned int k = 0; k < deg; ++k)
3905 if (std::abs(solution((i * deg + j) * deg + k)) > 1
e-14)
3906 nodal_values[((i + this->degree +
3913 solution((i * deg + j) * deg + k);
3918 for (
unsigned int q_point = 0; q_point < n_interior_points;
3922 support_point_values[q_point +
3928 for (
unsigned int i = 0; i <= deg; ++i)
3929 for (
unsigned int j = 0; j < 4; ++j)
3931 tmp -= nodal_values[i + (j + 8) * this->degree] *
3932 this->shape_value_component(
3933 i + (j + 8) * this->degree,
3934 this->generalized_support_points
3942 for (
unsigned int k = 0; k < deg; ++k)
3945 ((2 * j + 1) * deg + k +
3948 this->shape_value_component(
3949 i + ((2 * j + 1) * deg + k +
3952 this->generalized_support_points
3961 for (
unsigned int i = 0; i <= deg; ++i)
3962 for (
unsigned int j = 0; j < deg; ++j)
3963 for (
unsigned int k = 0; k < deg; ++k)
3964 system_rhs((i * deg + j) * deg + k) +=
3965 reference_quadrature.
weight(q_point) * tmp *
3966 lobatto_polynomials_grad[i].value(
3967 this->generalized_support_points
3972 n_face_points](2)) *
3973 lobatto_polynomials[j + 2].
value(
3974 this->generalized_support_points
3979 n_face_points](0)) *
3980 lobatto_polynomials[k + 2].
value(
3981 this->generalized_support_points
3989 system_matrix_inv.
vmult(solution, system_rhs);
3995 for (
unsigned int i = 0; i <= deg; ++i)
3996 for (
unsigned int j = 0; j < deg; ++j)
3997 for (
unsigned int k = 0; k < deg; ++k)
3998 if (std::abs(solution((i * deg + j) * deg + k)) > 1
e-14)
4004 this->degree] = solution((i * deg + j) * deg + k);
4018 std::pair<Table<2, bool>, std::vector<unsigned int>>
4022 for (
unsigned int d = 0;
d < dim; ++
d)
4023 for (
unsigned int i = 0; i < this->dofs_per_cell; ++i)
4024 constant_modes(
d, i) =
true;
4025 std::vector<unsigned int> components;
4026 for (
unsigned int d = 0;
d < dim; ++
d)
4027 components.push_back(
d);
4028 return std::pair<Table<2, bool>, std::vector<unsigned int>>(constant_modes,
4046 #include "fe_nedelec.inst"