24template <
int dim,
int spacedim>
31 std::vector<
bool>(compute_num_dofs(order), true),
41 for (
unsigned int comp = 0; comp < this->
n_components(); ++comp)
53template <
int dim,
int spacedim>
64template <
int dim,
int spacedim>
69 const unsigned int )
const
78template <
int dim,
int spacedim>
89template <
int dim,
int spacedim>
94 const unsigned int )
const
102template <
int dim,
int spacedim>
113template <
int dim,
int spacedim>
118 const unsigned int )
const
126template <
int dim,
int spacedim>
127std::unique_ptr<typename ::FiniteElement<dim, spacedim>::InternalDataBase>
137 typename ::FiniteElement<dim, spacedim>::InternalDataBase>
138 data_ptr = std::make_unique<InternalData>();
140 data.update_each = requires_update_flags(update_flags);
143 const unsigned int degree(this->degree - 1);
149 const unsigned int n_line_dofs = this->n_dofs_per_line() * lines_per_cell;
152 const unsigned int n_face_dofs = this->n_dofs_per_quad(0) * faces_per_cell;
155 const unsigned int n_q_points = quadrature.
size();
158 data.sigma_imj_values.resize(
160 std::vector<std::vector<double>>(vertices_per_cell,
161 std::vector<double>(vertices_per_cell)));
163 data.sigma_imj_grads.resize(vertices_per_cell,
164 std::vector<std::vector<double>>(
165 vertices_per_cell, std::vector<double>(dim)));
170 data.shape_values.resize(this->n_dofs_per_cell(),
176 data.shape_grads.resize(this->n_dofs_per_cell(),
183 data.shape_hessians.resize(this->n_dofs_per_cell(),
188 std::vector<Point<dim>> p_list(n_q_points);
198 std::vector<std::vector<double>> sigma(
199 n_q_points, std::vector<double>(lines_per_cell));
200 std::vector<std::vector<double>> lambda(
201 n_q_points, std::vector<double>(lines_per_cell));
203 for (
unsigned int q = 0; q < n_q_points; ++q)
205 sigma[q][0] = (1.0 - p_list[q][0]) + (1.0 - p_list[q][1]);
206 sigma[q][1] = p_list[q][0] + (1.0 - p_list[q][1]);
207 sigma[q][2] = (1.0 - p_list[q][0]) + p_list[q][1];
208 sigma[q][3] = p_list[q][0] + p_list[q][1];
210 lambda[q][0] = (1.0 - p_list[q][0]) * (1.0 - p_list[q][1]);
211 lambda[q][1] = p_list[q][0] * (1.0 - p_list[q][1]);
212 lambda[q][2] = (1.0 - p_list[q][0]) * p_list[q][1];
213 lambda[q][3] = p_list[q][0] * p_list[q][1];
214 for (
unsigned int i = 0; i < vertices_per_cell; ++i)
216 for (
unsigned int j = 0; j < vertices_per_cell; ++j)
218 data.sigma_imj_values[q][i][j] =
219 sigma[q][i] - sigma[q][j];
232 {-1, -1}, {1, -1}, {-1, 1}, {1, 1}};
233 int sigma_imj_sign[vertices_per_cell][vertices_per_cell];
234 unsigned int sigma_imj_component[vertices_per_cell]
237 for (
unsigned int i = 0; i < vertices_per_cell; ++i)
239 for (
unsigned int j = 0; j < vertices_per_cell; ++j)
246 sigma_imj_sign[i][j] = (i < j) ? -1 : 1;
247 sigma_imj_sign[i][j] = (i == j) ? 0 : sigma_imj_sign[i][j];
251 sigma_imj_component[i][j] = 0;
252 for (
unsigned int d = 0; d < dim; ++d)
255 sigma_comp_signs[i][d] - sigma_comp_signs[j][d];
260 sigma_imj_component[i][j] = d;
267 data.sigma_imj_grads[i][j][sigma_imj_component[i][j]] =
268 2.0 * sigma_imj_sign[i][j];
276 data.edge_sigma_values.resize(lines_per_cell);
277 data.edge_sigma_grads.resize(lines_per_cell);
278 for (
unsigned int m = 0; m < lines_per_cell; ++m)
280 data.edge_sigma_values[m].resize(n_q_points);
283 data.edge_sigma_grads[m].resize(dim);
293 data.edge_lambda_values.resize(lines_per_cell,
294 std::vector<double>(n_q_points));
295 data.edge_lambda_grads_2d.resize(lines_per_cell,
296 std::vector<double>(dim));
297 for (
unsigned int m = 0; m < lines_per_cell; ++m)
302 const unsigned int e1(
304 const unsigned int e2(
306 for (
unsigned int q = 0; q < n_q_points; ++q)
308 data.edge_sigma_values[m][q] =
309 data.sigma_imj_values[q][e2][e1];
310 data.edge_lambda_values[m][q] = lambda[q][e1] + lambda[q][e2];
313 data.edge_sigma_grads[m][edge_sigma_direction[m]] = -2.0;
315 data.edge_lambda_grads_2d[0] = {-1.0, 0.0};
316 data.edge_lambda_grads_2d[1] = {1.0, 0.0};
317 data.edge_lambda_grads_2d[2] = {0.0, -1.0};
318 data.edge_lambda_grads_2d[3] = {0.0, 1.0};
386 const unsigned int cell_type1_offset = n_line_dofs;
393 const unsigned int cell_type2_offset =
394 cell_type1_offset + degree * degree;
402 const unsigned int cell_type3_offset1 =
403 cell_type2_offset + degree * degree;
404 const unsigned int cell_type3_offset2 = cell_type3_offset1 + degree;
409 std::vector<Point<dim>> cell_points(n_q_points);
410 for (
unsigned int q = 0; q < n_q_points; ++q)
412 for (
unsigned int d = 0; d < dim; ++d)
414 cell_points[q][d] = 2.0 * p_list[q][d] - 1.0;
419 for (
unsigned int q = 0; q < n_q_points; ++q)
431 const unsigned int poly_length =
436 std::vector<std::vector<double>> polyx(
437 degree, std::vector<double>(poly_length));
438 std::vector<std::vector<double>> polyy(
439 degree, std::vector<double>(poly_length));
440 for (
unsigned int i = 0; i < degree; ++i)
445 IntegratedLegendrePolynomials[i + 2].value(
446 cell_points[q][0], polyx[i]);
447 IntegratedLegendrePolynomials[i + 2].value(
448 cell_points[q][1], polyy[i]);
453 for (
unsigned int j = 0; j < degree; ++j)
455 const unsigned int shift_j(j * degree);
456 for (
unsigned int i = 0; i < degree; ++i)
458 const unsigned int shift_ij(i + shift_j);
461 const unsigned int dof_index1(cell_type1_offset +
463 data.shape_values[dof_index1][q][0] =
464 2.0 * polyx[i][1] * polyy[j][0];
465 data.shape_values[dof_index1][q][1] =
466 2.0 * polyx[i][0] * polyy[j][1];
469 const unsigned int dof_index2(cell_type2_offset +
471 data.shape_values[dof_index2][q][0] =
472 data.shape_values[dof_index1][q][0];
473 data.shape_values[dof_index2][q][1] =
474 -1.0 * data.shape_values[dof_index1][q][1];
477 const unsigned int dof_index3_1(cell_type3_offset1 +
479 data.shape_values[dof_index3_1][q][0] = polyy[j][0];
480 data.shape_values[dof_index3_1][q][1] = 0.0;
482 const unsigned int dof_index3_2(cell_type3_offset2 +
484 data.shape_values[dof_index3_2][q][0] = 0.0;
485 data.shape_values[dof_index3_2][q][1] = polyx[j][0];
490 for (
unsigned int j = 0; j < degree; ++j)
492 const unsigned int shift_j(j * degree);
493 for (
unsigned int i = 0; i < degree; ++i)
495 const unsigned int shift_ij(i + shift_j);
498 const unsigned int dof_index1(cell_type1_offset +
500 data.shape_grads[dof_index1][q][0][0] =
501 4.0 * polyx[i][2] * polyy[j][0];
502 data.shape_grads[dof_index1][q][0][1] =
503 4.0 * polyx[i][1] * polyy[j][1];
504 data.shape_grads[dof_index1][q][1][0] =
505 data.shape_grads[dof_index1][q][0][1];
506 data.shape_grads[dof_index1][q][1][1] =
507 4.0 * polyx[i][0] * polyy[j][2];
510 const unsigned int dof_index2(cell_type2_offset +
512 data.shape_grads[dof_index2][q][0][0] =
513 data.shape_grads[dof_index1][q][0][0];
514 data.shape_grads[dof_index2][q][0][1] =
515 data.shape_grads[dof_index1][q][0][1];
516 data.shape_grads[dof_index2][q][1][0] =
517 -1.0 * data.shape_grads[dof_index1][q][1][0];
518 data.shape_grads[dof_index2][q][1][1] =
519 -1.0 * data.shape_grads[dof_index1][q][1][1];
522 const unsigned int dof_index3_1(cell_type3_offset1 +
524 data.shape_grads[dof_index3_1][q][0][0] = 0.0;
525 data.shape_grads[dof_index3_1][q][0][1] =
527 data.shape_grads[dof_index3_1][q][1][0] = 0.0;
528 data.shape_grads[dof_index3_1][q][1][1] = 0.0;
530 const unsigned int dof_index3_2(cell_type3_offset2 +
532 data.shape_grads[dof_index3_2][q][0][0] = 0.0;
533 data.shape_grads[dof_index3_2][q][0][1] = 0.0;
534 data.shape_grads[dof_index3_2][q][1][0] =
536 data.shape_grads[dof_index3_2][q][1][1] = 0.0;
541 for (
unsigned int j = 0; j < degree; ++j)
543 const unsigned int shift_j(j * degree);
544 for (
unsigned int i = 0; i < degree; ++i)
546 const unsigned int shift_ij(i + shift_j);
549 const unsigned int dof_index1(cell_type1_offset +
551 data.shape_hessians[dof_index1][q][0][0][0] =
552 8.0 * polyx[i][3] * polyy[j][0];
553 data.shape_hessians[dof_index1][q][1][0][0] =
554 8.0 * polyx[i][2] * polyy[j][1];
556 data.shape_hessians[dof_index1][q][0][1][0] =
557 data.shape_hessians[dof_index1][q][1][0][0];
558 data.shape_hessians[dof_index1][q][1][1][0] =
559 8.0 * polyx[i][1] * polyy[j][2];
561 data.shape_hessians[dof_index1][q][0][0][1] =
562 data.shape_hessians[dof_index1][q][1][0][0];
563 data.shape_hessians[dof_index1][q][1][0][1] =
564 data.shape_hessians[dof_index1][q][1][1][0];
566 data.shape_hessians[dof_index1][q][0][1][1] =
567 data.shape_hessians[dof_index1][q][1][1][0];
568 data.shape_hessians[dof_index1][q][1][1][1] =
569 8.0 * polyx[i][0] * polyy[j][3];
574 const unsigned int dof_index2(cell_type2_offset +
576 for (
unsigned int d = 0; d < dim; ++d)
578 data.shape_hessians[dof_index2][q][0][0][d] =
579 data.shape_hessians[dof_index1][q][0][0][d];
580 data.shape_hessians[dof_index2][q][0][1][d] =
581 data.shape_hessians[dof_index1][q][0][1][d];
582 data.shape_hessians[dof_index2][q][1][0][d] =
584 data.shape_hessians[dof_index1][q][1][0][d];
585 data.shape_hessians[dof_index2][q][1][1][d] =
587 data.shape_hessians[dof_index1][q][1][1][d];
591 const unsigned int dof_index3_1(cell_type3_offset1 +
593 data.shape_hessians[dof_index3_1][q][0][0][0] = 0.0;
594 data.shape_hessians[dof_index3_1][q][0][0][1] = 0.0;
595 data.shape_hessians[dof_index3_1][q][0][1][0] = 0.0;
596 data.shape_hessians[dof_index3_1][q][0][1][1] =
598 data.shape_hessians[dof_index3_1][q][1][0][0] = 0.0;
599 data.shape_hessians[dof_index3_1][q][1][0][1] = 0.0;
600 data.shape_hessians[dof_index3_1][q][1][1][0] = 0.0;
601 data.shape_hessians[dof_index3_1][q][1][1][1] = 0.0;
603 const unsigned int dof_index3_2(cell_type3_offset2 +
605 data.shape_hessians[dof_index3_2][q][0][0][0] = 0.0;
606 data.shape_hessians[dof_index3_2][q][0][0][1] = 0.0;
607 data.shape_hessians[dof_index3_2][q][0][1][0] = 0.0;
608 data.shape_hessians[dof_index3_2][q][0][1][1] = 0.0;
609 data.shape_hessians[dof_index3_2][q][1][0][0] =
611 data.shape_hessians[dof_index3_2][q][1][0][1] = 0.0;
612 data.shape_hessians[dof_index3_2][q][1][1][0] = 0.0;
613 data.shape_hessians[dof_index3_2][q][1][1][1] = 0.0;
624 std::vector<std::vector<double>> sigma(
625 n_q_points, std::vector<double>(lines_per_cell));
626 std::vector<std::vector<double>> lambda(
627 n_q_points, std::vector<double>(lines_per_cell));
628 for (
unsigned int q = 0; q < n_q_points; ++q)
630 sigma[q][0] = (1.0 - p_list[q][0]) + (1.0 - p_list[q][1]) +
633 p_list[q][0] + (1.0 - p_list[q][1]) + (1 - p_list[q][2]);
635 (1.0 - p_list[q][0]) + p_list[q][1] + (1 - p_list[q][2]);
636 sigma[q][3] = p_list[q][0] + p_list[q][1] + (1 - p_list[q][2]);
638 (1.0 - p_list[q][0]) + (1.0 - p_list[q][1]) + p_list[q][2];
639 sigma[q][5] = p_list[q][0] + (1.0 - p_list[q][1]) + p_list[q][2];
640 sigma[q][6] = (1.0 - p_list[q][0]) + p_list[q][1] + p_list[q][2];
641 sigma[q][7] = p_list[q][0] + p_list[q][1] + p_list[q][2];
643 lambda[q][0] = (1.0 - p_list[q][0]) * (1.0 - p_list[q][1]) *
644 (1.0 - p_list[q][2]);
646 p_list[q][0] * (1.0 - p_list[q][1]) * (1.0 - p_list[q][2]);
648 (1.0 - p_list[q][0]) * p_list[q][1] * (1.0 - p_list[q][2]);
649 lambda[q][3] = p_list[q][0] * p_list[q][1] * (1.0 - p_list[q][2]);
651 (1.0 - p_list[q][0]) * (1.0 - p_list[q][1]) * p_list[q][2];
652 lambda[q][5] = p_list[q][0] * (1.0 - p_list[q][1]) * p_list[q][2];
653 lambda[q][6] = (1.0 - p_list[q][0]) * p_list[q][1] * p_list[q][2];
654 lambda[q][7] = p_list[q][0] * p_list[q][1] * p_list[q][2];
658 for (
unsigned int i = 0; i < vertices_per_cell; ++i)
660 for (
unsigned int j = 0; j < vertices_per_cell; ++j)
662 data.sigma_imj_values[q][i][j] =
663 sigma[q][i] - sigma[q][j];
696 int sigma_imj_sign[vertices_per_cell][vertices_per_cell];
697 unsigned int sigma_imj_component[vertices_per_cell]
700 for (
unsigned int i = 0; i < vertices_per_cell; ++i)
702 for (
unsigned int j = 0; j < vertices_per_cell; ++j)
709 sigma_imj_sign[i][j] = (i < j) ? -1 : 1;
710 sigma_imj_sign[i][j] = (i == j) ? 0 : sigma_imj_sign[i][j];
714 sigma_imj_component[i][j] = 0;
715 for (
unsigned int d = 0; d < dim; ++d)
718 sigma_comp_signs[i][d] - sigma_comp_signs[j][d];
723 sigma_imj_component[i][j] = d;
730 data.sigma_imj_grads[i][j][sigma_imj_component[i][j]] =
731 2.0 * sigma_imj_sign[i][j];
738 data.edge_sigma_values.resize(lines_per_cell);
739 data.edge_lambda_values.resize(lines_per_cell);
740 data.edge_sigma_grads.resize(lines_per_cell);
741 data.edge_lambda_grads_3d.resize(lines_per_cell);
742 data.edge_lambda_gradgrads_3d.resize(lines_per_cell);
743 for (
unsigned int m = 0; m < lines_per_cell; ++m)
745 data.edge_sigma_values[m].resize(n_q_points);
746 data.edge_lambda_values[m].resize(n_q_points);
749 data.edge_sigma_grads[m].resize(dim);
751 data.edge_lambda_grads_3d[m].resize(n_q_points);
752 for (
unsigned int q = 0; q < n_q_points; ++q)
754 data.edge_lambda_grads_3d[m][q].resize(dim);
758 data.edge_lambda_gradgrads_3d[m].resize(dim);
759 for (
unsigned int d = 0; d < dim; ++d)
761 data.edge_lambda_gradgrads_3d[m][d].resize(dim);
768 1, 1, 0, 0, 1, 1, 0, 0, 2, 2, 2, 2};
770 for (
unsigned int m = 0; m < lines_per_cell; ++m)
775 const unsigned int e1(
777 const unsigned int e2(
780 for (
unsigned int q = 0; q < n_q_points; ++q)
782 data.edge_sigma_values[m][q] =
783 data.sigma_imj_values[q][e2][e1];
784 data.edge_lambda_values[m][q] = lambda[q][e1] + lambda[q][e2];
787 data.edge_sigma_grads[m][edge_sigma_direction[m]] = -2.0;
790 for (
unsigned int q = 0; q < n_q_points; ++q)
792 double x(p_list[q][0]);
793 double y(p_list[q][1]);
794 double z(p_list[q][2]);
795 data.edge_lambda_grads_3d[0][q] = {z - 1.0, 0.0, x - 1.0};
796 data.edge_lambda_grads_3d[1][q] = {1.0 - z, 0.0, -x};
797 data.edge_lambda_grads_3d[2][q] = {0.0, z - 1.0, y - 1.0};
798 data.edge_lambda_grads_3d[3][q] = {0.0, 1.0 - z, -y};
799 data.edge_lambda_grads_3d[4][q] = {-z, 0.0, 1.0 - x};
800 data.edge_lambda_grads_3d[5][q] = {z, 0.0, x};
801 data.edge_lambda_grads_3d[6][q] = {0.0, -z, 1.0 - y};
802 data.edge_lambda_grads_3d[7][q] = {0.0, z, y};
803 data.edge_lambda_grads_3d[8][q] = {y - 1.0, x - 1.0, 0.0};
804 data.edge_lambda_grads_3d[9][q] = {1.0 - y, -x, 0.0};
805 data.edge_lambda_grads_3d[10][q] = {-y, 1.0 - x, 0.0};
806 data.edge_lambda_grads_3d[11][q] = {y, x, 0.0};
836 for (
unsigned int m = 0; m < lines_per_cell; ++m)
838 data.edge_lambda_gradgrads_3d[m][edge_lambda_directions[m][0]]
839 [edge_lambda_directions[m][1]] =
841 data.edge_lambda_gradgrads_3d[m][edge_lambda_directions[m][1]]
842 [edge_lambda_directions[m][0]] =
850 data.face_lambda_values.resize(faces_per_cell);
851 data.face_lambda_grads.resize(faces_per_cell);
853 for (
unsigned int m = 0; m < faces_per_cell; ++m)
855 data.face_lambda_values[m].resize(n_q_points);
856 data.face_lambda_grads[m].resize(3);
859 for (
unsigned int q = 0; q < n_q_points; ++q)
861 double x(p_list[q][0]);
862 double y(p_list[q][1]);
863 double z(p_list[q][2]);
864 data.face_lambda_values[0][q] = 1.0 - x;
865 data.face_lambda_values[1][q] = x;
866 data.face_lambda_values[2][q] = 1.0 - y;
867 data.face_lambda_values[3][q] = y;
868 data.face_lambda_values[4][q] = 1.0 - z;
869 data.face_lambda_values[5][q] = z;
872 data.face_lambda_grads[0] = {-1.0, 0.0, 0.0};
873 data.face_lambda_grads[1] = {1.0, 0.0, 0.0};
874 data.face_lambda_grads[2] = {0.0, -1.0, 0.0};
875 data.face_lambda_grads[3] = {0.0, 1.0, 0.0};
876 data.face_lambda_grads[4] = {0.0, 0.0, -1.0};
877 data.face_lambda_grads[5] = {0.0, 0.0, 1.0};
890 const unsigned int cell_type1_offset(n_line_dofs +
903 const unsigned int cell_type2_offset1(
904 cell_type1_offset + degree * degree * degree);
905 const unsigned int cell_type2_offset2(
906 cell_type2_offset1 + degree * degree * degree);
916 const unsigned int cell_type3_offset1(
917 cell_type2_offset2 + degree * degree * degree);
918 const unsigned int cell_type3_offset2(cell_type3_offset1 +
920 const unsigned int cell_type3_offset3(cell_type3_offset2 +
924 std::vector<Point<dim>> cell_points(n_q_points);
925 for (
unsigned int q = 0; q < n_q_points; ++q)
927 for (
unsigned int d = 0; d < dim; ++d)
929 cell_points[q][d] = 2.0 * p_list[q][d] - 1.0;
937 const unsigned int poly_length =
943 for (
unsigned int q = 0; q < n_q_points; ++q)
951 std::vector<std::vector<double>> polyx(
952 degree, std::vector<double>(poly_length));
953 std::vector<std::vector<double>> polyy(
954 degree, std::vector<double>(poly_length));
955 std::vector<std::vector<double>> polyz(
956 degree, std::vector<double>(poly_length));
957 for (
unsigned int i = 0; i < degree; ++i)
960 IntegratedLegendrePolynomials[i + 2].value(
961 cell_points[q][0], polyx[i]);
962 IntegratedLegendrePolynomials[i + 2].value(
963 cell_points[q][1], polyy[i]);
964 IntegratedLegendrePolynomials[i + 2].value(
965 cell_points[q][2], polyz[i]);
970 for (
unsigned int k = 0; k < degree; ++k)
972 const unsigned int shift_k(k * degree * degree);
973 const unsigned int shift_j(
976 for (
unsigned int j = 0; j < degree; ++j)
978 const unsigned int shift_jk(j * degree +
980 for (
unsigned int i = 0; i < degree; ++i)
982 const unsigned int shift_ijk(shift_jk +
986 const unsigned int dof_index1(
987 cell_type1_offset + shift_ijk);
989 data.shape_values[dof_index1][q][0] =
990 2.0 * polyx[i][1] * polyy[j][0] *
992 data.shape_values[dof_index1][q][1] =
993 2.0 * polyx[i][0] * polyy[j][1] *
995 data.shape_values[dof_index1][q][2] =
996 2.0 * polyx[i][0] * polyy[j][0] *
1000 const unsigned int dof_index2_1(
1001 cell_type2_offset1 + shift_ijk);
1002 const unsigned int dof_index2_2(
1003 cell_type2_offset2 + shift_ijk);
1005 data.shape_values[dof_index2_1][q][0] =
1006 data.shape_values[dof_index1][q][0];
1007 data.shape_values[dof_index2_1][q][1] =
1009 data.shape_values[dof_index1][q][1];
1010 data.shape_values[dof_index2_1][q][2] =
1011 data.shape_values[dof_index1][q][2];
1013 data.shape_values[dof_index2_2][q][0] =
1014 data.shape_values[dof_index1][q][0];
1015 data.shape_values[dof_index2_2][q][1] =
1017 data.shape_values[dof_index1][q][1];
1018 data.shape_values[dof_index2_2][q][2] =
1020 data.shape_values[dof_index1][q][2];
1024 const unsigned int shift_ij(
1027 const unsigned int dof_index3_1(
1028 cell_type3_offset1 + shift_ij);
1029 const unsigned int dof_index3_2(
1030 cell_type3_offset2 + shift_ij);
1031 const unsigned int dof_index3_3(
1032 cell_type3_offset3 + shift_ij);
1034 data.shape_values[dof_index3_1][q][0] =
1035 polyy[j][0] * polyz[k][0];
1036 data.shape_values[dof_index3_1][q][1] = 0.0;
1037 data.shape_values[dof_index3_1][q][2] = 0.0;
1039 data.shape_values[dof_index3_2][q][0] = 0.0;
1040 data.shape_values[dof_index3_2][q][1] =
1041 polyx[j][0] * polyz[k][0];
1042 data.shape_values[dof_index3_2][q][2] = 0.0;
1044 data.shape_values[dof_index3_3][q][0] = 0.0;
1045 data.shape_values[dof_index3_3][q][1] = 0.0;
1046 data.shape_values[dof_index3_3][q][2] =
1047 polyx[j][0] * polyy[k][0];
1053 for (
unsigned int k = 0; k < degree; ++k)
1055 const unsigned int shift_k(k * degree * degree);
1056 const unsigned int shift_j(
1059 for (
unsigned int j = 0; j < degree; ++j)
1061 const unsigned int shift_jk(j * degree +
1063 for (
unsigned int i = 0; i < degree; ++i)
1065 const unsigned int shift_ijk(shift_jk +
1069 const unsigned int dof_index1(
1070 cell_type1_offset + shift_ijk);
1072 data.shape_grads[dof_index1][q][0][0] =
1073 4.0 * polyx[i][2] * polyy[j][0] *
1075 data.shape_grads[dof_index1][q][0][1] =
1076 4.0 * polyx[i][1] * polyy[j][1] *
1078 data.shape_grads[dof_index1][q][0][2] =
1079 4.0 * polyx[i][1] * polyy[j][0] *
1082 data.shape_grads[dof_index1][q][1][0] =
1083 data.shape_grads[dof_index1][q][0][1];
1084 data.shape_grads[dof_index1][q][1][1] =
1085 4.0 * polyx[i][0] * polyy[j][2] *
1087 data.shape_grads[dof_index1][q][1][2] =
1088 4.0 * polyx[i][0] * polyy[j][1] *
1091 data.shape_grads[dof_index1][q][2][0] =
1092 data.shape_grads[dof_index1][q][0][2];
1093 data.shape_grads[dof_index1][q][2][1] =
1094 data.shape_grads[dof_index1][q][1][2];
1095 data.shape_grads[dof_index1][q][2][2] =
1096 4.0 * polyx[i][0] * polyy[j][0] *
1100 const unsigned int dof_index2_1(
1101 cell_type2_offset1 + shift_ijk);
1102 const unsigned int dof_index2_2(
1103 cell_type2_offset2 + shift_ijk);
1105 for (
unsigned int d = 0; d < dim; ++d)
1107 data.shape_grads[dof_index2_1][q][0]
1110 .shape_grads[dof_index1][q][0][d];
1111 data.shape_grads[dof_index2_1][q][1]
1115 .shape_grads[dof_index1][q][1][d];
1116 data.shape_grads[dof_index2_1][q][2]
1119 .shape_grads[dof_index1][q][2][d];
1121 data.shape_grads[dof_index2_2][q][0]
1124 .shape_grads[dof_index1][q][0][d];
1125 data.shape_grads[dof_index2_2][q][1]
1129 .shape_grads[dof_index1][q][1][d];
1130 data.shape_grads[dof_index2_2][q][2]
1134 .shape_grads[dof_index1][q][2][d];
1139 const unsigned int shift_ij(
1142 const unsigned int dof_index3_1(
1143 cell_type3_offset1 + shift_ij);
1144 const unsigned int dof_index3_2(
1145 cell_type3_offset2 + shift_ij);
1146 const unsigned int dof_index3_3(
1147 cell_type3_offset3 + shift_ij);
1148 for (
unsigned int d1 = 0; d1 < dim; ++d1)
1150 for (
unsigned int d2 = 0; d2 < dim; ++d2)
1152 data.shape_grads[dof_index3_1][q][d1]
1154 data.shape_grads[dof_index3_2][q][d1]
1156 data.shape_grads[dof_index3_3][q][d1]
1160 data.shape_grads[dof_index3_1][q][0][1] =
1161 2.0 * polyy[j][1] * polyz[k][0];
1162 data.shape_grads[dof_index3_1][q][0][2] =
1163 2.0 * polyy[j][0] * polyz[k][1];
1165 data.shape_grads[dof_index3_2][q][1][0] =
1166 2.0 * polyx[j][1] * polyz[k][0];
1167 data.shape_grads[dof_index3_2][q][1][2] =
1168 2.0 * polyx[j][0] * polyz[k][1];
1170 data.shape_grads[dof_index3_3][q][2][0] =
1171 2.0 * polyx[j][1] * polyy[k][0];
1172 data.shape_grads[dof_index3_3][q][2][1] =
1173 2.0 * polyx[j][0] * polyy[k][1];
1179 for (
unsigned int k = 0; k < degree; ++k)
1181 const unsigned int shift_k(k * degree * degree);
1182 const unsigned int shift_j(
1186 for (
unsigned int j = 0; j < degree; ++j)
1188 const unsigned int shift_jk(j * degree +
1190 for (
unsigned int i = 0; i < degree; ++i)
1192 const unsigned int shift_ijk(shift_jk +
1196 const unsigned int dof_index1(
1197 cell_type1_offset + shift_ijk);
1199 data.shape_hessians[dof_index1][q][0][0]
1201 8.0 * polyx[i][3] * polyy[j][0] *
1203 data.shape_hessians[dof_index1][q][1][0]
1205 8.0 * polyx[i][2] * polyy[j][1] *
1207 data.shape_hessians[dof_index1][q][2][0]
1209 8.0 * polyx[i][2] * polyy[j][0] *
1212 data.shape_hessians[dof_index1][q][0][1]
1214 data.shape_hessians[dof_index1][q][1][0]
1216 data.shape_hessians[dof_index1][q][1][1]
1218 8.0 * polyx[i][1] * polyy[j][2] *
1220 data.shape_hessians[dof_index1][q][2][1]
1222 8.0 * polyx[i][1] * polyy[j][1] *
1225 data.shape_hessians[dof_index1][q][0][2]
1227 data.shape_hessians[dof_index1][q][2][0]
1229 data.shape_hessians[dof_index1][q][1][2]
1231 data.shape_hessians[dof_index1][q][2][1]
1233 data.shape_hessians[dof_index1][q][2][2]
1235 8.0 * polyx[i][1] * polyy[j][0] *
1239 data.shape_hessians[dof_index1][q][0][0]
1241 data.shape_hessians[dof_index1][q][1][0]
1243 data.shape_hessians[dof_index1][q][1][0]
1245 data.shape_hessians[dof_index1][q][1][1]
1247 data.shape_hessians[dof_index1][q][2][0]
1249 data.shape_hessians[dof_index1][q][2][1]
1252 data.shape_hessians[dof_index1][q][0][1]
1254 data.shape_hessians[dof_index1][q][1][1]
1256 data.shape_hessians[dof_index1][q][1][1]
1258 8.0 * polyx[i][0] * polyy[j][3] *
1260 data.shape_hessians[dof_index1][q][2][1]
1262 8.0 * polyx[i][0] * polyy[j][2] *
1265 data.shape_hessians[dof_index1][q][0][2]
1267 data.shape_hessians[dof_index1][q][2][1]
1269 data.shape_hessians[dof_index1][q][1][2]
1271 data.shape_hessians[dof_index1][q][2][1]
1273 data.shape_hessians[dof_index1][q][2][2]
1275 8.0 * polyx[i][0] * polyy[j][1] *
1279 data.shape_hessians[dof_index1][q][0][0]
1281 data.shape_hessians[dof_index1][q][2][0]
1283 data.shape_hessians[dof_index1][q][1][0]
1285 data.shape_hessians[dof_index1][q][2][1]
1287 data.shape_hessians[dof_index1][q][2][0]
1289 data.shape_hessians[dof_index1][q][2][2]
1292 data.shape_hessians[dof_index1][q][0][1]
1294 data.shape_hessians[dof_index1][q][2][1]
1296 data.shape_hessians[dof_index1][q][1][1]
1298 data.shape_hessians[dof_index1][q][2][1]
1300 data.shape_hessians[dof_index1][q][2][1]
1302 data.shape_hessians[dof_index1][q][2][2]
1305 data.shape_hessians[dof_index1][q][0][2]
1307 data.shape_hessians[dof_index1][q][2][2]
1309 data.shape_hessians[dof_index1][q][1][2]
1311 data.shape_hessians[dof_index1][q][2][2]
1313 data.shape_hessians[dof_index1][q][2][2]
1315 8.0 * polyx[i][0] * polyy[j][0] *
1320 const unsigned int dof_index2_1(
1321 cell_type2_offset1 + shift_ijk);
1322 const unsigned int dof_index2_2(
1323 cell_type2_offset2 + shift_ijk);
1325 for (
unsigned int d1 = 0; d1 < dim; ++d1)
1327 for (
unsigned int d2 = 0; d2 < dim;
1331 .shape_hessians[dof_index2_1][q]
1334 .shape_hessians[dof_index1][q]
1337 .shape_hessians[dof_index2_1][q]
1341 .shape_hessians[dof_index1][q]
1344 .shape_hessians[dof_index2_1][q]
1347 .shape_hessians[dof_index1][q]
1351 .shape_hessians[dof_index2_2][q]
1354 .shape_hessians[dof_index1][q]
1357 .shape_hessians[dof_index2_2][q]
1361 .shape_hessians[dof_index1][q]
1364 .shape_hessians[dof_index2_2][q]
1368 .shape_hessians[dof_index1][q]
1375 const unsigned int shift_ij(
1378 const unsigned int dof_index3_1(
1379 cell_type3_offset1 + shift_ij);
1380 const unsigned int dof_index3_2(
1381 cell_type3_offset2 + shift_ij);
1382 const unsigned int dof_index3_3(
1383 cell_type3_offset3 + shift_ij);
1384 for (
unsigned int d1 = 0; d1 < dim; ++d1)
1386 for (
unsigned int d2 = 0; d2 < dim; ++d2)
1388 for (
unsigned int d3 = 0; d3 < dim;
1392 .shape_hessians[dof_index3_1][q]
1396 .shape_hessians[dof_index3_2][q]
1400 .shape_hessians[dof_index3_3][q]
1407 .shape_hessians[dof_index3_1][q][0][1][1] =
1408 4.0 * polyy[j][2] * polyz[k][0];
1410 .shape_hessians[dof_index3_1][q][0][1][2] =
1411 4.0 * polyy[j][1] * polyz[k][1];
1414 .shape_hessians[dof_index3_1][q][0][2][1] =
1416 .shape_hessians[dof_index3_1][q][0][1][2];
1418 .shape_hessians[dof_index3_1][q][0][2][2] =
1419 4.0 * polyy[j][0] * polyz[k][2];
1423 .shape_hessians[dof_index3_2][q][1][0][0] =
1424 4.0 * polyx[j][2] * polyz[k][0];
1426 .shape_hessians[dof_index3_2][q][1][0][2] =
1427 4.0 * polyx[j][1] * polyz[k][1];
1430 .shape_hessians[dof_index3_2][q][1][2][0] =
1432 .shape_hessians[dof_index3_2][q][1][0][2];
1434 .shape_hessians[dof_index3_2][q][1][2][2] =
1435 4.0 * polyx[j][0] * polyz[k][2];
1439 .shape_hessians[dof_index3_3][q][2][0][0] =
1440 4.0 * polyx[j][2] * polyy[k][0];
1442 .shape_hessians[dof_index3_3][q][2][0][1] =
1443 4.0 * polyx[j][1] * polyy[k][1];
1446 .shape_hessians[dof_index3_3][q][2][1][0] =
1448 .shape_hessians[dof_index3_3][q][2][0][1];
1450 .shape_hessians[dof_index3_3][q][2][1][1] =
1451 4.0 * polyx[j][0] * polyy[k][2];
1470template <
int dim,
int spacedim>
1491 const unsigned int n_q_points = quadrature.
size();
1494 fe_data.
shape_values.size() == this->n_dofs_per_cell(),
1496 this->n_dofs_per_cell()));
1502 const unsigned int degree(
1507 const unsigned int vertices_per_line(2);
1511 std::vector<std::vector<unsigned int>> edge_order(
1512 lines_per_cell, std::vector<unsigned int>(vertices_per_line));
1524 std::vector<int> edge_sign(lines_per_cell);
1525 for (
unsigned int m = 0; m < lines_per_cell; ++m)
1527 unsigned int v0_loc =
1529 unsigned int v1_loc =
1531 unsigned int v0_glob = cell->vertex_index(v0_loc);
1532 unsigned int v1_glob = cell->vertex_index(v1_loc);
1534 if (v0_glob > v1_glob)
1537 edge_sign[m] = -1.0;
1581 std::vector<std::vector<double>> edge_sigma_values(
1583 std::vector<std::vector<double>> edge_sigma_grads(
1586 std::vector<std::vector<double>> edge_lambda_values(
1588 std::vector<std::vector<double>> edge_lambda_grads(
1592 for (
unsigned int m = 0; m < lines_per_cell; ++m)
1594 std::transform(edge_sigma_values[m].begin(),
1595 edge_sigma_values[m].end(),
1596 edge_sigma_values[m].begin(),
1597 [&](
const double edge_sigma_value) {
1598 return edge_sign[m] * edge_sigma_value;
1601 std::transform(edge_sigma_grads[m].begin(),
1602 edge_sigma_grads[m].end(),
1603 edge_sigma_grads[m].begin(),
1604 [&](
const double edge_sigma_grad) {
1605 return edge_sign[m] * edge_sigma_grad;
1612 const unsigned int poly_length =
1618 for (
unsigned int m = 0; m < lines_per_cell; ++m)
1620 const unsigned int shift_m(m * this->n_dofs_per_line());
1621 for (
unsigned int q = 0; q < n_q_points; ++q)
1624 std::vector<std::vector<double>> poly(
1625 degree, std::vector<double>(poly_length));
1626 for (
unsigned int i = 1; i < degree + 1; ++i)
1631 IntegratedLegendrePolynomials[i + 1].value(
1632 edge_sigma_values[m][q], poly[i - 1]);
1637 for (
unsigned int d = 0; d < dim; ++d)
1640 0.5 * edge_sigma_grads[m][d] *
1641 edge_lambda_values[m][q];
1644 for (
unsigned int i = 1; i < degree + 1; ++i)
1646 const unsigned int poly_index = i - 1;
1647 const unsigned int dof_index(i + shift_m);
1648 for (
unsigned int d = 0; d < dim; ++d)
1651 edge_sigma_grads[m][d] *
1652 poly[poly_index][1] *
1653 edge_lambda_values[m][q] +
1654 poly[poly_index][0] *
1655 edge_lambda_grads[m][d];
1662 for (
unsigned int d1 = 0; d1 < dim; ++d1)
1664 for (
unsigned int d2 = 0; d2 < dim; ++d2)
1669 0.5 * edge_sigma_grads[m][d1] *
1670 edge_lambda_grads[m][d2];
1674 for (
unsigned int i = 1; i < degree + 1; ++i)
1676 const unsigned int poly_index = i - 1;
1677 const unsigned int dof_index(i + shift_m);
1680 edge_sigma_grads[m][0] *
1681 edge_sigma_grads[m][0] *
1682 edge_lambda_values[m][q] * poly[poly_index][2];
1685 (edge_sigma_grads[m][0] *
1686 edge_lambda_grads[m][1] +
1687 edge_sigma_grads[m][1] *
1688 edge_lambda_grads[m][0]) *
1689 poly[poly_index][1];
1695 edge_sigma_grads[m][1] *
1696 edge_sigma_grads[m][1] *
1697 edge_lambda_values[m][q] * poly[poly_index][2];
1703 for (
unsigned int d1 = 0; d1 < dim; ++d1)
1705 for (
unsigned int d2 = 0; d2 < dim; ++d2)
1707 for (
unsigned int d3 = 0; d3 < dim; ++d3)
1716 for (
unsigned int i = 0; i < degree; ++i)
1718 const unsigned int dof_index(i + 1 + shift_m);
1720 for (
unsigned int d1 = 0; d1 < dim; ++d1)
1722 for (
unsigned int d2 = 0; d2 < dim; ++d2)
1724 for (
unsigned int d3 = 0; d3 < dim; ++d3)
1728 edge_sigma_grads[m][d1] *
1729 edge_sigma_grads[m][d2] *
1730 edge_sigma_grads[m][d3] *
1732 edge_lambda_values[m][q] +
1734 (edge_sigma_grads[m][d1] *
1735 edge_sigma_grads[m][d2] *
1736 edge_lambda_grads[m][d3] +
1737 edge_sigma_grads[m][d3] *
1738 edge_sigma_grads[m][d1] *
1739 edge_lambda_grads[m][d2] +
1740 edge_sigma_grads[m][d3] *
1741 edge_sigma_grads[m][d2] *
1742 edge_lambda_grads[m][d1]);
1760 std::vector<int> edge_sign(lines_per_cell);
1761 for (
unsigned int m = 0; m < lines_per_cell; ++m)
1763 unsigned int v0_loc =
1765 unsigned int v1_loc =
1767 unsigned int v0_glob = cell->vertex_index(v0_loc);
1768 unsigned int v1_glob = cell->vertex_index(v1_loc);
1770 if (v0_glob > v1_glob)
1773 edge_sign[m] = -1.0;
1819 std::vector<std::vector<double>> edge_sigma_values(
1821 std::vector<std::vector<double>> edge_lambda_values(
1823 std::vector<std::vector<double>> edge_sigma_grads(
1825 std::vector<std::vector<std::vector<double>>> edge_lambda_grads(
1827 std::vector<std::vector<std::vector<double>>>
1831 for (
unsigned int m = 0; m < lines_per_cell; ++m)
1833 std::transform(edge_sigma_values[m].begin(),
1834 edge_sigma_values[m].end(),
1835 edge_sigma_values[m].begin(),
1836 [&](
const double edge_sigma_value) {
1837 return edge_sign[m] * edge_sigma_value;
1839 std::transform(edge_sigma_grads[m].begin(),
1840 edge_sigma_grads[m].end(),
1841 edge_sigma_grads[m].begin(),
1842 [&](
const double edge_sigma_grad) {
1843 return edge_sign[m] * edge_sigma_grad;
1851 const unsigned int poly_length =
1856 std::vector<std::vector<double>> poly(
1857 degree, std::vector<double>(poly_length));
1858 for (
unsigned int m = 0; m < lines_per_cell; ++m)
1860 const unsigned int shift_m(m * this->n_dofs_per_line());
1861 for (
unsigned int q = 0; q < n_q_points; ++q)
1867 for (
unsigned int i = 0; i < degree; ++i)
1869 IntegratedLegendrePolynomials[i + 2].value(
1870 edge_sigma_values[m][q], poly[i]);
1876 for (
unsigned int d = 0; d < dim; ++d)
1879 0.5 * edge_sigma_grads[m][d] *
1880 edge_lambda_values[m][q];
1885 for (
unsigned int i = 0; i < degree; ++i)
1887 const unsigned int dof_index(i + 1 + shift_m);
1888 for (
unsigned int d = 0; d < dim; ++d)
1891 edge_sigma_grads[m][d] * poly[i][1] *
1892 edge_lambda_values[m][q] +
1893 poly[i][0] * edge_lambda_grads[m][q][d];
1901 for (
unsigned int d1 = 0; d1 < dim; ++d1)
1903 for (
unsigned int d2 = 0; d2 < dim; ++d2)
1906 0.5 * edge_sigma_grads[m][d1] *
1907 edge_lambda_grads[m][q][d2];
1913 for (
unsigned int i = 0; i < degree; ++i)
1915 const unsigned int dof_index(i + 1 + shift_m);
1917 for (
unsigned int d1 = 0; d1 < dim; ++d1)
1919 for (
unsigned int d2 = 0; d2 < dim; ++d2)
1923 edge_sigma_grads[m][d1] *
1924 edge_sigma_grads[m][d2] *
1926 edge_lambda_values[m][q] +
1927 (edge_sigma_grads[m][d1] *
1928 edge_lambda_grads[m][q][d2] +
1929 edge_sigma_grads[m][d2] *
1930 edge_lambda_grads[m][q][d1]) *
1932 edge_lambda_gradgrads_3d[m][d1]
1943 for (
unsigned int d1 = 0; d1 < dim; ++d1)
1945 for (
unsigned int d2 = 0; d2 < dim; ++d2)
1947 for (
unsigned int d3 = 0; d3 < dim; ++d3)
1951 0.5 * edge_sigma_grads[m][d1] *
1952 edge_lambda_gradgrads_3d[m][d3][d2];
1960 for (
unsigned int i = 0; i < degree; ++i)
1962 const unsigned int dof_index(i + 1 + shift_m);
1964 for (
unsigned int d1 = 0; d1 < dim; ++d1)
1966 for (
unsigned int d2 = 0; d2 < dim; ++d2)
1968 for (
unsigned int d3 = 0; d3 < dim;
1974 edge_sigma_grads[m][d1] *
1975 edge_sigma_grads[m][d2] *
1976 edge_sigma_grads[m][d3] *
1978 edge_lambda_values[m][q] +
1980 (edge_sigma_grads[m][d1] *
1981 edge_sigma_grads[m][d2] *
1982 edge_lambda_grads[m][q]
1984 edge_sigma_grads[m][d3] *
1985 edge_sigma_grads[m][d1] *
1986 edge_lambda_grads[m][q]
1988 edge_sigma_grads[m][d3] *
1989 edge_sigma_grads[m][d2] *
1990 edge_lambda_grads[m][q]
1993 (edge_sigma_grads[m][d1] *
1994 edge_lambda_gradgrads_3d
1996 edge_sigma_grads[m][d2] *
1997 edge_lambda_gradgrads_3d
1999 edge_sigma_grads[m][d3] *
2000 edge_lambda_gradgrads_3d
2022template <
int dim,
int spacedim>
2042 const unsigned int degree(
2053 const unsigned int n_q_points = quadrature.
size();
2056 fe_data.
shape_values.size() == this->n_dofs_per_cell(),
2058 this->n_dofs_per_cell()));
2065 const unsigned int vertices_per_face(
2070 const unsigned int n_line_dofs =
2074 std::vector<std::vector<unsigned int>> face_orientation(
2075 faces_per_cell, std::vector<unsigned int>(vertices_per_face));
2085 {1, 2}, {0, 3}, {0, 3}, {1, 2}};
2087 for (
unsigned int m = 0; m < faces_per_cell; ++m)
2091 unsigned int current_max = 0;
2092 unsigned int current_glob = cell->vertex_index(
2094 for (
unsigned int v = 1; v < vertices_per_face; ++v)
2101 current_glob = cell->vertex_index(
2105 face_orientation[m][0] =
2110 m, vertex_opposite_on_face[current_max]);
2115 m, vertices_adjacent_on_face[current_max][0])) >
2117 m, vertices_adjacent_on_face[current_max][1])))
2119 face_orientation[m][1] =
2121 m, vertices_adjacent_on_face[current_max][0]);
2122 face_orientation[m][3] =
2124 m, vertices_adjacent_on_face[current_max][1]);
2128 face_orientation[m][1] =
2130 m, vertices_adjacent_on_face[current_max][1]);
2131 face_orientation[m][3] =
2133 m, vertices_adjacent_on_face[current_max][0]);
2139 std::vector<std::vector<double>> face_xi_values(
2140 faces_per_cell, std::vector<double>(n_q_points));
2141 std::vector<std::vector<double>> face_xi_grads(
2142 faces_per_cell, std::vector<double>(dim));
2143 std::vector<std::vector<double>> face_eta_values(
2144 faces_per_cell, std::vector<double>(n_q_points));
2145 std::vector<std::vector<double>> face_eta_grads(
2146 faces_per_cell, std::vector<double>(dim));
2148 std::vector<std::vector<double>> face_lambda_values(
2149 faces_per_cell, std::vector<double>(n_q_points));
2150 std::vector<std::vector<double>> face_lambda_grads(
2151 faces_per_cell, std::vector<double>(dim));
2152 for (
unsigned int m = 0; m < faces_per_cell; ++m)
2154 for (
unsigned int q = 0; q < n_q_points; ++q)
2156 face_xi_values[m][q] =
2158 [face_orientation[m][1]];
2159 face_eta_values[m][q] =
2161 [face_orientation[m][3]];
2164 for (
unsigned int d = 0; d < dim; ++d)
2166 face_xi_grads[m][d] =
2168 [face_orientation[m][1]][d];
2169 face_eta_grads[m][d] =
2171 [face_orientation[m][3]][d];
2177 const unsigned int poly_length =
2182 std::vector<std::vector<double>> polyxi(
2183 degree, std::vector<double>(poly_length));
2184 std::vector<std::vector<double>> polyeta(
2185 degree, std::vector<double>(poly_length));
2188 for (
unsigned int m = 0; m < faces_per_cell; ++m)
2191 const unsigned int shift_m(m * this->n_dofs_per_quad(0));
2199 const unsigned int face_type1_offset(n_line_dofs + shift_m);
2209 const unsigned int face_type2_offset(face_type1_offset +
2221 const unsigned int face_type3_offset1(face_type2_offset +
2223 const unsigned int face_type3_offset2(face_type3_offset1 +
2227 for (
unsigned int q = 0; q < n_q_points; ++q)
2236 for (
unsigned int i = 0; i < degree; ++i)
2239 IntegratedLegendrePolynomials[i + 2].value(
2240 face_xi_values[m][q], polyxi[i]);
2241 IntegratedLegendrePolynomials[i + 2].value(
2242 face_eta_values[m][q], polyeta[i]);
2247 for (
unsigned int j = 0; j < degree; ++j)
2249 const unsigned int shift_j(j * degree);
2250 for (
unsigned int i = 0; i < degree; ++i)
2252 const unsigned int shift_ij(shift_j + i);
2254 const unsigned int dof_index1(face_type1_offset +
2256 for (
unsigned int d = 0; d < dim; ++d)
2259 (face_xi_grads[m][d] * polyxi[i][1] *
2261 face_eta_grads[m][d] * polyxi[i][0] *
2263 face_lambda_values[m][q] +
2264 face_lambda_grads[m][d] * polyxi[i][0] *
2268 const unsigned int dof_index2(face_type2_offset +
2270 for (
unsigned int d = 0; d < dim; ++d)
2273 (face_xi_grads[m][d] * polyxi[i][1] *
2275 face_eta_grads[m][d] * polyxi[i][0] *
2277 face_lambda_values[m][q];
2281 const unsigned int dof_index3_1(face_type3_offset1 +
2283 const unsigned int dof_index3_2(face_type3_offset2 +
2285 for (
unsigned int d = 0; d < dim; ++d)
2288 face_xi_grads[m][d] * polyeta[j][0] *
2289 face_lambda_values[m][q];
2292 face_eta_grads[m][d] * polyxi[j][0] *
2293 face_lambda_values[m][q];
2299 for (
unsigned int j = 0; j < degree; ++j)
2301 const unsigned int shift_j(j * degree);
2302 for (
unsigned int i = 0; i < degree; ++i)
2304 const unsigned int shift_ij(shift_j + i);
2306 const unsigned int dof_index1(face_type1_offset +
2308 for (
unsigned int d1 = 0; d1 < dim; ++d1)
2310 for (
unsigned int d2 = 0; d2 < dim; ++d2)
2314 (face_xi_grads[m][d1] *
2315 face_xi_grads[m][d2] * polyxi[i][2] *
2317 (face_xi_grads[m][d1] *
2318 face_eta_grads[m][d2] +
2319 face_xi_grads[m][d2] *
2320 face_eta_grads[m][d1]) *
2321 polyxi[i][1] * polyeta[j][1] +
2322 face_eta_grads[m][d1] *
2323 face_eta_grads[m][d2] *
2324 polyxi[i][0] * polyeta[j][2]) *
2325 face_lambda_values[m][q] +
2326 (face_xi_grads[m][d2] * polyxi[i][1] *
2328 face_eta_grads[m][d2] * polyxi[i][0] *
2330 face_lambda_grads[m][d1] +
2331 (face_xi_grads[m][d1] * polyxi[i][1] *
2333 face_eta_grads[m][d1] * polyxi[i][0] *
2335 face_lambda_grads[m][d2];
2339 const unsigned int dof_index2(face_type2_offset +
2341 for (
unsigned int d1 = 0; d1 < dim; ++d1)
2343 for (
unsigned int d2 = 0; d2 < dim; ++d2)
2347 (face_xi_grads[m][d1] *
2348 face_xi_grads[m][d2] * polyxi[i][2] *
2350 (face_xi_grads[m][d1] *
2351 face_eta_grads[m][d2] -
2352 face_xi_grads[m][d2] *
2353 face_eta_grads[m][d1]) *
2354 polyxi[i][1] * polyeta[j][1] -
2355 face_eta_grads[m][d1] *
2356 face_eta_grads[m][d2] *
2357 polyxi[i][0] * polyeta[j][2]) *
2358 face_lambda_values[m][q] +
2359 (face_xi_grads[m][d1] * polyxi[i][1] *
2361 face_eta_grads[m][d1] * polyxi[i][0] *
2363 face_lambda_grads[m][d2];
2368 const unsigned int dof_index3_1(face_type3_offset1 +
2370 const unsigned int dof_index3_2(face_type3_offset2 +
2372 for (
unsigned int d1 = 0; d1 < dim; ++d1)
2374 for (
unsigned int d2 = 0; d2 < dim; ++d2)
2377 face_xi_grads[m][d1] *
2378 (face_eta_grads[m][d2] * polyeta[j][1] *
2379 face_lambda_values[m][q] +
2380 face_lambda_grads[m][d2] * polyeta[j][0]);
2383 face_eta_grads[m][d1] *
2384 (face_xi_grads[m][d2] * polyxi[j][1] *
2385 face_lambda_values[m][q] +
2386 face_lambda_grads[m][d2] * polyxi[j][0]);
2393 for (
unsigned int j = 0; j < degree; ++j)
2395 const unsigned int shift_j(j * degree);
2396 for (
unsigned int i = 0; i < degree; ++i)
2398 const unsigned int shift_ij(shift_j + i);
2401 const unsigned int dof_index1(face_type1_offset +
2403 for (
unsigned int d1 = 0; d1 < dim; ++d1)
2405 for (
unsigned int d2 = 0; d2 < dim; ++d2)
2407 for (
unsigned int d3 = 0; d3 < dim; ++d3)
2412 face_xi_grads[m][d3] *
2413 (face_eta_grads[m][d1] *
2415 face_eta_grads[m][d2] *
2416 face_lambda_values[m][q] +
2418 face_lambda_grads[m][d2]) +
2420 face_eta_grads[m][d2] *
2421 face_lambda_grads[m][d1]) +
2424 face_eta_grads[m][d1] *
2425 face_eta_grads[m][d2] *
2426 face_eta_grads[m][d3] *
2427 face_lambda_values[m][q] +
2429 (face_eta_grads[m][d1] *
2430 face_eta_grads[m][d2] *
2431 face_lambda_grads[m][d3] +
2432 face_eta_grads[m][d3] *
2433 (face_eta_grads[m][d1] *
2434 face_lambda_grads[m]
2436 face_eta_grads[m][d2] *
2439 (polyxi[i][1] * polyeta[j][1] *
2440 face_eta_grads[m][d3] +
2441 polyxi[i][2] * polyeta[j][0] *
2442 face_xi_grads[m][d3]) *
2443 (face_xi_grads[m][d1] *
2444 face_lambda_grads[m][d2] +
2445 face_xi_grads[m][d2] *
2446 face_lambda_grads[m][d1]) +
2447 face_lambda_grads[m][d3] *
2448 (polyxi[i][2] * polyeta[j][0] *
2449 face_xi_grads[m][d1] *
2450 face_xi_grads[m][d2] +
2451 polyxi[i][1] * polyeta[j][1] *
2452 (face_xi_grads[m][d1] *
2453 face_eta_grads[m][d2] +
2454 face_xi_grads[m][d2] *
2455 face_eta_grads[m][d1])) +
2456 face_lambda_values[m][q] *
2457 (polyxi[i][3] * polyeta[j][0] *
2458 face_xi_grads[m][d1] *
2459 face_xi_grads[m][d2] *
2460 face_xi_grads[m][d3] +
2461 polyxi[i][1] * polyeta[j][2] *
2462 face_eta_grads[m][d3] *
2463 (face_xi_grads[m][d1] *
2464 face_eta_grads[m][d2] +
2465 face_xi_grads[m][d2] *
2466 face_eta_grads[m][d1]) +
2467 polyxi[i][2] * polyeta[j][1] *
2468 (face_xi_grads[m][d3] *
2469 face_xi_grads[m][d2] *
2470 face_eta_grads[m][d1] +
2471 face_xi_grads[m][d1] *
2472 (face_xi_grads[m][d2] *
2473 face_eta_grads[m][d3] +
2474 face_xi_grads[m][d3] *
2475 face_eta_grads[m][d2])));
2481 const unsigned int dof_index2(face_type2_offset +
2483 for (
unsigned int d1 = 0; d1 < dim; ++d1)
2485 for (
unsigned int d2 = 0; d2 < dim; ++d2)
2487 for (
unsigned int d3 = 0; d3 < dim; ++d3)
2491 face_xi_grads[m][d1] *
2492 (polyxi[i][1] * polyeta[j][1] *
2493 (face_eta_grads[m][d2] *
2494 face_lambda_grads[m][d3] +
2495 face_eta_grads[m][d3] *
2496 face_lambda_grads[m][d2]) +
2497 polyxi[i][2] * polyeta[j][0] *
2498 (face_xi_grads[m][d2] *
2499 face_lambda_grads[m][d3] +
2500 face_xi_grads[m][d3] *
2501 face_lambda_grads[m][d2]) +
2502 face_lambda_values[m][q] *
2503 (face_eta_grads[m][d2] *
2506 face_eta_grads[m][d3] +
2509 face_xi_grads[m][d3]) +
2510 face_xi_grads[m][d2] *
2513 face_eta_grads[m][d3] +
2516 face_xi_grads[m][d3]))) -
2518 face_eta_grads[m][d1] *
2519 (face_eta_grads[m][d2] *
2521 face_eta_grads[m][d3] *
2522 face_lambda_values[m][q] +
2524 face_lambda_grads[m][d3]) +
2526 face_eta_grads[m][d3] *
2527 face_lambda_grads[m][d2]) -
2528 face_eta_grads[m][d1] *
2530 face_xi_grads[m][d3] *
2532 face_eta_grads[m][d2] *
2533 face_lambda_values[m][q] +
2535 face_lambda_grads[m][d2]) +
2536 face_xi_grads[m][d2] *
2539 face_eta_grads[m][d3] *
2540 face_lambda_values[m]
2543 face_lambda_grads[m]
2545 polyxi[i][2] * polyeta[j][1] *
2546 face_xi_grads[m][d3] *
2547 face_lambda_values[m][q]));
2553 const unsigned int dof_index3_1(face_type3_offset1 +
2555 const unsigned int dof_index3_2(face_type3_offset2 +
2557 for (
unsigned int d1 = 0; d1 < dim; ++d1)
2559 for (
unsigned int d2 = 0; d2 < dim; ++d2)
2561 for (
unsigned int d3 = 0; d3 < dim; ++d3)
2565 face_xi_grads[m][d1] *
2566 (face_eta_grads[m][d2] *
2568 face_eta_grads[m][d3] *
2569 face_lambda_values[m][q] +
2571 face_lambda_grads[m][d3]) +
2572 face_lambda_grads[m][d2] *
2574 face_eta_grads[m][d3]);
2578 face_eta_grads[m][d1] *
2579 (face_xi_grads[m][d2] *
2581 face_xi_grads[m][d3] *
2582 face_lambda_values[m][q] +
2584 face_lambda_grads[m][d3]) +
2585 face_lambda_grads[m][d2] *
2586 polyxi[j][1] * face_xi_grads[m][d3]);
2600template <
int dim,
int spacedim>
2622 fill_edge_values(cell, quadrature, fe_data);
2623 if (dim == 3 && this->degree > 1)
2625 fill_face_values(cell, quadrature, fe_data);
2629 const unsigned int n_q_points = quadrature.
size();
2632 fe_data.
shape_values.size() == this->n_dofs_per_cell(),
2634 this->n_dofs_per_cell()));
2643 std::vector<Tensor<1, dim>> transformed_shape_values(n_q_points);
2644 for (
unsigned int dof = 0; dof < this->n_dofs_per_cell(); ++dof)
2646 const unsigned int first =
2648 this->get_nonzero_components(dof)
2649 .first_selected_component()];
2655 for (
unsigned int q = 0; q < n_q_points; ++q)
2657 for (
unsigned int d = 0; d < dim; ++d)
2660 transformed_shape_values[q][d];
2670 std::vector<Tensor<2, dim>> input(n_q_points);
2671 std::vector<Tensor<2, dim>> transformed_shape_grads(n_q_points);
2672 for (
unsigned int dof = 0; dof < this->n_dofs_per_cell(); ++dof)
2674 for (
unsigned int q = 0; q < n_q_points; ++q)
2683 const unsigned int first =
2685 this->get_nonzero_components(dof)
2686 .first_selected_component()];
2688 for (
unsigned int q = 0; q < n_q_points; ++q)
2690 for (
unsigned int d1 = 0; d1 < dim; ++d1)
2692 for (
unsigned int d2 = 0; d2 < dim; ++d2)
2694 transformed_shape_grads[q][d1] -=
2701 for (
unsigned int q = 0; q < n_q_points; ++q)
2703 for (
unsigned int d = 0; d < dim; ++d)
2706 transformed_shape_grads[q][d];
2716 std::vector<Tensor<3, dim>> input(n_q_points);
2717 std::vector<Tensor<3, dim>> transformed_shape_hessians(n_q_points);
2718 for (
unsigned int dof = 0; dof < this->n_dofs_per_cell(); ++dof)
2720 for (
unsigned int q = 0; q < n_q_points; ++q)
2729 const unsigned int first =
2731 this->get_nonzero_components(dof)
2732 .first_selected_component()];
2734 for (
unsigned int q = 0; q < n_q_points; ++q)
2736 for (
unsigned int d1 = 0; d1 < dim; ++d1)
2738 for (
unsigned int d2 = 0; d2 < dim; ++d2)
2740 for (
unsigned int d3 = 0; d3 < dim; ++d3)
2742 for (
unsigned int d4 = 0; d4 < dim; ++d4)
2744 transformed_shape_hessians[q][d1][d3][d4] -=
2748 [q][d2][d1][d3][d4]) +
2767 for (
unsigned int q = 0; q < n_q_points; ++q)
2769 for (
unsigned int d = 0; d < dim; ++d)
2772 transformed_shape_hessians[q][d];
2781template <
int dim,
int spacedim>
2785 const unsigned int face_no,
2817 fill_edge_values(cell,
2821 if (dim == 3 && this->degree > 1)
2823 fill_face_values(cell,
2825 this->reference_cell(), quadrature[0]),
2830 const unsigned int n_q_points = quadrature[0].
size();
2834 cell->face_orientation(face_no),
2835 cell->face_flip(face_no),
2836 cell->face_rotation(face_no),
2843 std::vector<Tensor<1, dim>> transformed_shape_values(n_q_points);
2844 for (
unsigned int dof = 0; dof < this->n_dofs_per_cell(); ++dof)
2853 const unsigned int first =
2855 this->get_nonzero_components(dof)
2856 .first_selected_component()];
2858 for (
unsigned int q = 0; q < n_q_points; ++q)
2860 for (
unsigned int d = 0; d < dim; ++d)
2863 transformed_shape_values[q][d];
2872 std::vector<Tensor<2, dim>> input(n_q_points);
2873 std::vector<Tensor<2, dim>> transformed_shape_grads(n_q_points);
2874 for (
unsigned int dof = 0; dof < this->n_dofs_per_cell(); ++dof)
2876 for (
unsigned int q = 0; q < n_q_points; ++q)
2885 const unsigned int first =
2887 this->get_nonzero_components(dof)
2888 .first_selected_component()];
2890 for (
unsigned int q = 0; q < n_q_points; ++q)
2892 for (
unsigned int d1 = 0; d1 < dim; ++d1)
2894 for (
unsigned int d2 = 0; d2 < dim; ++d2)
2896 transformed_shape_grads[q][d1] -=
2903 for (
unsigned int q = 0; q < n_q_points; ++q)
2905 for (
unsigned int d = 0; d < dim; ++d)
2908 transformed_shape_grads[q][d];
2917 std::vector<Tensor<3, dim>> input(n_q_points);
2918 std::vector<Tensor<3, dim>> transformed_shape_hessians(n_q_points);
2919 for (
unsigned int dof = 0; dof < this->n_dofs_per_cell(); ++dof)
2921 for (
unsigned int q = 0; q < n_q_points; ++q)
2929 const unsigned int first =
2931 this->get_nonzero_components(dof)
2932 .first_selected_component()];
2934 for (
unsigned int q = 0; q < n_q_points; ++q)
2936 for (
unsigned int d1 = 0; d1 < dim; ++d1)
2938 for (
unsigned int d2 = 0; d2 < dim; ++d2)
2940 for (
unsigned int d3 = 0; d3 < dim; ++d3)
2942 for (
unsigned int d4 = 0; d4 < dim; ++d4)
2944 transformed_shape_hessians[q][d1][d3][d4] -=
2948 [q][d2][d1][d3][d4]) +
2967 for (
unsigned int q = 0; q < n_q_points; ++q)
2969 for (
unsigned int d = 0; d < dim; ++d)
2972 transformed_shape_hessians[q][d];
2981template <
int dim,
int spacedim>
2985 const unsigned int ,
2986 const unsigned int ,
3001template <
int dim,
int spacedim>
3027template <
int dim,
int spacedim>
3033 std::ostringstream namebuf;
3034 namebuf <<
"FE_NedelecSZ<" << dim <<
">(" << this->degree - 1 <<
")";
3036 return namebuf.str();
3041template <
int dim,
int spacedim>
3042std::unique_ptr<FiniteElement<dim, dim>>
3045 return std::make_unique<FE_NedelecSZ<dim, spacedim>>(*this);
3050template <
int dim,
int spacedim>
3051std::vector<unsigned int>
3060 std::vector<unsigned int> dpo;
3063 dpo.push_back(degree + 1);
3065 dpo.push_back(2 * degree * (degree + 1));
3067 dpo.push_back(3 * degree * degree * (degree + 1));
3074template <
int dim,
int spacedim>
3083 return 2 * (degree + 1) * (degree + 2);
3086 return 3 * (degree + 1) * (degree + 2) * (degree + 2);
3098template <
int dim,
int spacedim>
3103 IntegratedLegendrePolynomials =
3110#include "fe_nedelec_sz.inst"
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
std::vector< std::vector< double > > edge_lambda_grads_2d
std::vector< std::vector< Tensor< 1, dim > > > shape_values
std::vector< std::vector< std::vector< double > > > sigma_imj_values
std::vector< std::vector< double > > edge_sigma_grads
std::vector< std::vector< double > > edge_sigma_values
std::vector< std::vector< DerivativeForm< 2, dim, dim > > > shape_hessians
std::vector< std::vector< std::vector< double > > > sigma_imj_grads
std::vector< std::vector< double > > face_lambda_values
std::vector< std::vector< DerivativeForm< 1, dim, dim > > > shape_grads
std::vector< std::vector< std::vector< double > > > edge_lambda_grads_3d
std::vector< std::vector< double > > edge_lambda_values
std::vector< std::vector< std::vector< double > > > edge_lambda_gradgrads_3d
std::vector< std::vector< double > > face_lambda_grads
unsigned int compute_num_dofs(const unsigned int degree) const
virtual void fill_fe_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, dim > &mapping, const typename Mapping< dim, dim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, dim > &mapping_data, const typename FiniteElement< dim, dim >::InternalDataBase &fedata, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &data) const override
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
virtual void fill_fe_subface_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, dim > &mapping, const typename Mapping< dim, dim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, dim > &mapping_data, const typename FiniteElement< dim, dim >::InternalDataBase &fedata, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &data) const override
FE_NedelecSZ(const unsigned int order)
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
void create_polynomials(const unsigned int degree)
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
virtual void fill_fe_face_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const Mapping< dim, dim > &mapping, const typename Mapping< dim, dim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, dim > &mapping_data, const typename FiniteElement< dim, dim >::InternalDataBase &fedata, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &data) const override
virtual std::string get_name() 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 Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const override
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
void fill_face_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const Quadrature< dim > &quadrature, const InternalData &fedata) const
void fill_edge_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const Quadrature< dim > &quadrature, const InternalData &fedata) const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
unsigned int n_components() const
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > component_to_base_table
static std::vector< Polynomials::Polynomial< double > > generate_complete_basis(const unsigned int degree)
Abstract base class for mapping classes.
virtual void transform(const ArrayView< const Tensor< 1, dim > > &input, const MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim > > &output) const =0
static DataSetDescriptor face(const ReferenceCell &reference_cell, const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)
const std::vector< Point< dim > > & get_points() const
unsigned int size() const
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
@ update_jacobian_pushed_forward_2nd_derivatives
@ update_jacobian_pushed_forward_grads
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_covariant_transformation
Covariant transformation.
@ update_gradients
Shape function gradients.
@ update_default
No update.
@ mapping_covariant_gradient
@ mapping_covariant_hessian
static unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static unsigned int line_to_cell_vertices(const unsigned int line, const unsigned int vertex)