16 #include <deal.II/base/array_view.h>
17 #include <deal.II/base/memory_consumption.h>
18 #include <deal.II/base/qprojector.h>
19 #include <deal.II/base/quadrature.h>
20 #include <deal.II/base/signaling_nan.h>
21 #include <deal.II/base/std_cxx14/memory.h>
22 #include <deal.II/base/tensor.h>
24 #include <deal.II/dofs/dof_accessor.h>
26 #include <deal.II/fe/fe_values.h>
27 #include <deal.II/fe/mapping_cartesian.h>
29 #include <deal.II/grid/tria.h>
30 #include <deal.II/grid/tria_iterator.h>
32 #include <deal.II/lac/full_matrix.h>
38 DEAL_II_NAMESPACE_OPEN
41 template <
int dim,
int spacedim>
46 template <
int dim,
int spacedim>
50 , volume_element(
numbers::signaling_nan<double>())
51 , quadrature_points(q.get_points())
56 template <
int dim,
int spacedim>
68 template <
int dim,
int spacedim>
77 template <
int dim,
int spacedim>
80 const UpdateFlags in)
const
88 if (out & update_boundary_forms)
89 out |= update_normal_vectors;
96 template <
int dim,
int spacedim>
97 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
101 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
102 std_cxx14::make_unique<InternalData>(q);
103 auto &data = dynamic_cast<InternalData &>(*data_ptr);
115 template <
int dim,
int spacedim>
116 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
118 const UpdateFlags update_flags,
121 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
122 std_cxx14::make_unique<InternalData>(
124 auto &data = dynamic_cast<InternalData &>(*data_ptr);
133 data.update_each = update_flags;
140 template <
int dim,
int spacedim>
141 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
143 const UpdateFlags update_flags,
146 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
147 std_cxx14::make_unique<InternalData>(
149 auto &data = dynamic_cast<InternalData &>(*data_ptr);
158 data.update_each = update_flags;
165 template <
int dim,
int spacedim>
169 const unsigned int face_no,
170 const unsigned int sub_no,
197 cell->face(face_no)->has_children() ||
203 !cell->face(face_no)->has_children() ||
204 (sub_no < cell->face(face_no)->n_children()),
205 ExcIndexRange(sub_no, 0, cell->face(face_no)->n_children()));
249 if (update_flags & update_quadrature_points)
258 cell->face_orientation(face_no),
259 cell->face_flip(face_no),
260 cell->face_rotation(face_no),
261 quadrature_points.size()) :
266 cell->face_orientation(face_no),
267 cell->face_flip(face_no),
268 cell->face_rotation(face_no),
269 quadrature_points.size(),
270 cell->subface_case(face_no))));
272 for (
unsigned int i = 0; i < quadrature_points.size(); ++i)
274 quadrature_points[i] = start;
275 for (
unsigned int d = 0; d < dim; ++d)
276 quadrature_points[i](d) +=
289 if (update_flags & update_normal_vectors)
299 std::fill(normal_vectors.begin(),
300 normal_vectors.end(),
312 std::fill(normal_vectors.begin(),
313 normal_vectors.end(),
327 std::fill(normal_vectors.begin(),
328 normal_vectors.end(),
341 template <
int dim,
int spacedim>
353 Assert(dynamic_cast<const InternalData *>(&internal_data) !=
nullptr,
355 const InternalData &data = static_cast<const InternalData &>(internal_data);
357 std::vector<Tensor<1, dim>> dummy;
369 if (data.update_each & (update_JxW_values | update_volume_elements))
372 double J = data.cell_extents[0];
373 for (
unsigned int d = 1; d < dim; ++d)
374 J *= data.cell_extents[d];
375 data.volume_element = J;
376 if (data.update_each & update_JxW_values)
377 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
382 if (data.update_each & update_jacobians)
384 for (
unsigned int i = 0; i < output_data.
jacobians.size(); ++i)
387 for (
unsigned int j = 0; j < dim; ++j)
388 output_data.
jacobians[i][j][j] = data.cell_extents[j];
392 if (data.update_each & update_jacobian_grads)
394 for (
unsigned int i = 0; i < output_data.
jacobian_grads.size(); ++i)
397 if (data.update_each & update_jacobian_pushed_forward_grads)
399 for (
unsigned int i = 0;
406 if (data.update_each & update_jacobian_2nd_derivatives)
413 if (data.update_each & update_jacobian_pushed_forward_2nd_derivatives)
415 for (
unsigned int i = 0;
421 if (data.update_each & update_jacobian_3rd_derivatives)
428 if (data.update_each & update_jacobian_pushed_forward_3rd_derivatives)
430 for (
unsigned int i = 0;
438 if (data.update_each & update_inverse_jacobians)
443 for (
unsigned int j = 0; j < dim; ++j)
447 return cell_similarity;
452 template <
int dim,
int spacedim>
456 const unsigned int face_no,
466 Assert(dynamic_cast<const InternalData *>(&internal_data) !=
nullptr,
468 const InternalData &data = static_cast<const InternalData &>(internal_data);
481 for (
unsigned int d = 0;
d < dim; ++
d)
483 J *= data.cell_extents[
d];
485 if (data.update_each & update_JxW_values)
486 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
489 if (data.update_each & update_boundary_forms)
490 for (
unsigned int i = 0; i < output_data.
boundary_forms.size(); ++i)
493 if (data.update_each & update_volume_elements)
495 J = data.cell_extents[0];
496 for (
unsigned int d = 1;
d < dim; ++
d)
497 J *= data.cell_extents[d];
498 data.volume_element = J;
501 if (data.update_each & update_jacobians)
502 for (
unsigned int i = 0; i < output_data.
jacobians.size(); ++i)
505 for (
unsigned int d = 0;
d < dim; ++
d)
506 output_data.
jacobians[i][d][d] = data.cell_extents[d];
509 if (data.update_each & update_jacobian_grads)
510 for (
unsigned int i = 0; i < output_data.
jacobian_grads.size(); ++i)
513 if (data.update_each & update_jacobian_pushed_forward_grads)
514 for (
unsigned int i = 0;
519 if (data.update_each & update_jacobian_2nd_derivatives)
525 if (data.update_each & update_jacobian_pushed_forward_2nd_derivatives)
526 for (
unsigned int i = 0;
532 if (data.update_each & update_jacobian_3rd_derivatives)
538 if (data.update_each & update_jacobian_pushed_forward_3rd_derivatives)
539 for (
unsigned int i = 0;
545 if (data.update_each & update_inverse_jacobians)
549 for (
unsigned int d = 0;
d < dim; ++
d)
556 template <
int dim,
int spacedim>
560 const unsigned int face_no,
561 const unsigned int subface_no,
569 Assert(dynamic_cast<const InternalData *>(&internal_data) !=
nullptr,
571 const InternalData &data = static_cast<const InternalData &>(internal_data);
584 for (
unsigned int d = 0;
d < dim; ++
d)
586 J *= data.cell_extents[
d];
588 if (data.update_each & update_JxW_values)
594 const unsigned int n_subfaces =
595 cell->face(face_no)->has_children() ?
596 cell->face(face_no)->n_children() :
598 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
602 if (data.update_each & update_boundary_forms)
603 for (
unsigned int i = 0; i < output_data.
boundary_forms.size(); ++i)
606 if (data.update_each & update_volume_elements)
608 J = data.cell_extents[0];
609 for (
unsigned int d = 1;
d < dim; ++
d)
610 J *= data.cell_extents[d];
611 data.volume_element = J;
614 if (data.update_each & update_jacobians)
615 for (
unsigned int i = 0; i < output_data.
jacobians.size(); ++i)
618 for (
unsigned int d = 0;
d < dim; ++
d)
619 output_data.
jacobians[i][d][d] = data.cell_extents[d];
622 if (data.update_each & update_jacobian_grads)
623 for (
unsigned int i = 0; i < output_data.
jacobian_grads.size(); ++i)
626 if (data.update_each & update_jacobian_pushed_forward_grads)
627 for (
unsigned int i = 0;
632 if (data.update_each & update_jacobian_2nd_derivatives)
638 if (data.update_each & update_jacobian_pushed_forward_2nd_derivatives)
639 for (
unsigned int i = 0;
645 if (data.update_each & update_jacobian_3rd_derivatives)
651 if (data.update_each & update_jacobian_pushed_forward_3rd_derivatives)
652 for (
unsigned int i = 0;
658 if (data.update_each & update_inverse_jacobians)
662 for (
unsigned int d = 0;
d < dim; ++
d)
669 template <
int dim,
int spacedim>
678 Assert(dynamic_cast<const InternalData *>(&mapping_data) !=
nullptr,
680 const InternalData &data = static_cast<const InternalData &>(mapping_data);
682 switch (mapping_type)
688 "update_covariant_transformation"));
690 for (
unsigned int i = 0; i < output.size(); ++i)
691 for (
unsigned int d = 0; d < dim; ++d)
700 "update_contravariant_transformation"));
702 for (
unsigned int i = 0; i < output.size(); ++i)
703 for (
unsigned int d = 0; d < dim; ++d)
711 "update_contravariant_transformation"));
714 "update_volume_elements"));
716 for (
unsigned int i = 0; i < output.size(); ++i)
717 for (
unsigned int d = 0; d < dim; ++d)
729 template <
int dim,
int spacedim>
738 Assert(dynamic_cast<const InternalData *>(&mapping_data) !=
nullptr,
740 const InternalData &data = static_cast<const InternalData &>(mapping_data);
742 switch (mapping_type)
748 "update_covariant_transformation"));
750 for (
unsigned int i = 0; i < output.size(); ++i)
751 for (
unsigned int d1 = 0; d1 < dim; ++d1)
752 for (
unsigned int d2 = 0; d2 < dim; ++d2)
753 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2];
761 "update_contravariant_transformation"));
763 for (
unsigned int i = 0; i < output.size(); ++i)
764 for (
unsigned int d1 = 0; d1 < dim; ++d1)
765 for (
unsigned int d2 = 0; d2 < dim; ++d2)
766 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2];
774 "update_covariant_transformation"));
776 for (
unsigned int i = 0; i < output.size(); ++i)
777 for (
unsigned int d1 = 0; d1 < dim; ++d1)
778 for (
unsigned int d2 = 0; d2 < dim; ++d2)
779 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2] /
788 "update_contravariant_transformation"));
790 for (
unsigned int i = 0; i < output.size(); ++i)
791 for (
unsigned int d1 = 0; d1 < dim; ++d1)
792 for (
unsigned int d2 = 0; d2 < dim; ++d2)
793 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
802 "update_contravariant_transformation"));
805 "update_volume_elements"));
807 for (
unsigned int i = 0; i < output.size(); ++i)
808 for (
unsigned int d1 = 0; d1 < dim; ++d1)
809 for (
unsigned int d2 = 0; d2 < dim; ++d2)
810 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
819 "update_contravariant_transformation"));
822 "update_volume_elements"));
824 for (
unsigned int i = 0; i < output.size(); ++i)
825 for (
unsigned int d1 = 0; d1 < dim; ++d1)
826 for (
unsigned int d2 = 0; d2 < dim; ++d2)
827 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
839 template <
int dim,
int spacedim>
848 Assert(dynamic_cast<const InternalData *>(&mapping_data) !=
nullptr,
850 const InternalData &data = static_cast<const InternalData &>(mapping_data);
852 switch (mapping_type)
858 "update_covariant_transformation"));
860 for (
unsigned int i = 0; i < output.size(); ++i)
861 for (
unsigned int d1 = 0; d1 < dim; ++d1)
862 for (
unsigned int d2 = 0; d2 < dim; ++d2)
863 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2];
871 "update_contravariant_transformation"));
873 for (
unsigned int i = 0; i < output.size(); ++i)
874 for (
unsigned int d1 = 0; d1 < dim; ++d1)
875 for (
unsigned int d2 = 0; d2 < dim; ++d2)
876 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2];
884 "update_covariant_transformation"));
886 for (
unsigned int i = 0; i < output.size(); ++i)
887 for (
unsigned int d1 = 0; d1 < dim; ++d1)
888 for (
unsigned int d2 = 0; d2 < dim; ++d2)
889 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2] /
898 "update_contravariant_transformation"));
900 for (
unsigned int i = 0; i < output.size(); ++i)
901 for (
unsigned int d1 = 0; d1 < dim; ++d1)
902 for (
unsigned int d2 = 0; d2 < dim; ++d2)
903 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
912 "update_contravariant_transformation"));
915 "update_volume_elements"));
917 for (
unsigned int i = 0; i < output.size(); ++i)
918 for (
unsigned int d1 = 0; d1 < dim; ++d1)
919 for (
unsigned int d2 = 0; d2 < dim; ++d2)
920 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
929 "update_contravariant_transformation"));
932 "update_volume_elements"));
934 for (
unsigned int i = 0; i < output.size(); ++i)
935 for (
unsigned int d1 = 0; d1 < dim; ++d1)
936 for (
unsigned int d2 = 0; d2 < dim; ++d2)
937 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
948 template <
int dim,
int spacedim>
957 Assert(dynamic_cast<const InternalData *>(&mapping_data) !=
nullptr,
959 const InternalData &data = static_cast<const InternalData &>(mapping_data);
961 switch (mapping_type)
967 "update_covariant_transformation"));
969 for (
unsigned int q = 0; q < output.size(); ++q)
970 for (
unsigned int i = 0; i < spacedim; ++i)
971 for (
unsigned int j = 0; j < spacedim; ++j)
972 for (
unsigned int k = 0; k < spacedim; ++k)
974 output[q][i][j][k] = input[q][i][j][k] /
985 template <
int dim,
int spacedim>
994 Assert(dynamic_cast<const InternalData *>(&mapping_data) !=
nullptr,
996 const InternalData &data = static_cast<const InternalData &>(mapping_data);
998 switch (mapping_type)
1004 "update_covariant_transformation"));
1007 "update_contravariant_transformation"));
1009 for (
unsigned int q = 0; q < output.size(); ++q)
1010 for (
unsigned int i = 0; i < spacedim; ++i)
1011 for (
unsigned int j = 0; j < spacedim; ++j)
1012 for (
unsigned int k = 0; k < spacedim; ++k)
1014 output[q][i][j][k] =
1025 "update_covariant_transformation"));
1027 for (
unsigned int q = 0; q < output.size(); ++q)
1028 for (
unsigned int i = 0; i < spacedim; ++i)
1029 for (
unsigned int j = 0; j < spacedim; ++j)
1030 for (
unsigned int k = 0; k < spacedim; ++k)
1032 output[q][i][j][k] =
1044 "update_covariant_transformation"));
1047 "update_contravariant_transformation"));
1050 "update_volume_elements"));
1052 for (
unsigned int q = 0; q < output.size(); ++q)
1053 for (
unsigned int i = 0; i < spacedim; ++i)
1054 for (
unsigned int j = 0; j < spacedim; ++j)
1055 for (
unsigned int k = 0; k < spacedim; ++k)
1057 output[q][i][j][k] =
1072 template <
int dim,
int spacedim>
1083 length[0] = cell->vertex(1)(0) - start(0);
1086 length[0] = cell->vertex(1)(0) - start(0);
1087 length[1] = cell->vertex(2)(1) - start(1);
1090 length[0] = cell->vertex(1)(0) - start(0);
1091 length[1] = cell->vertex(2)(1) - start(1);
1092 length[2] = cell->vertex(4)(2) - start(2);
1099 for (
unsigned int d = 0; d < dim; ++d)
1100 p_real(d) += length[d] * p(d);
1107 template <
int dim,
int spacedim>
1113 if (dim != spacedim)
1122 real(0) /= cell->vertex(1)(0) - start(0);
1125 real(0) /= cell->vertex(1)(0) - start(0);
1126 real(1) /= cell->vertex(2)(1) - start(1);
1129 real(0) /= cell->vertex(1)(0) - start(0);
1130 real(1) /= cell->vertex(2)(1) - start(1);
1131 real(2) /= cell->vertex(4)(2) - start(2);
1140 template <
int dim,
int spacedim>
1141 std::unique_ptr<Mapping<dim, spacedim>>
1144 return std_cxx14::make_unique<MappingCartesian<dim, spacedim>>(*this);
1150 #include "mapping_cartesian.inst"
1153 DEAL_II_NAMESPACE_CLOSE