16 #include <deal.II/base/array_view.h> 17 #include <deal.II/base/memory_consumption.h> 18 #include <deal.II/base/polynomial.h> 19 #include <deal.II/base/qprojector.h> 20 #include <deal.II/base/quadrature.h> 21 #include <deal.II/base/quadrature_lib.h> 22 #include <deal.II/base/std_cxx14/memory.h> 23 #include <deal.II/base/tensor_product_polynomials.h> 24 #include <deal.II/base/thread_management.h> 25 #include <deal.II/base/utilities.h> 27 #include <deal.II/dofs/dof_accessor.h> 29 #include <deal.II/fe/fe_q.h> 30 #include <deal.II/fe/fe_system.h> 31 #include <deal.II/fe/fe_tools.h> 32 #include <deal.II/fe/fe_values.h> 33 #include <deal.II/fe/mapping.h> 34 #include <deal.II/fe/mapping_fe_field.h> 35 #include <deal.II/fe/mapping_q1.h> 37 #include <deal.II/grid/tria_iterator.h> 39 #include <deal.II/lac/block_vector.h> 40 #include <deal.II/lac/full_matrix.h> 41 #include <deal.II/lac/la_parallel_block_vector.h> 42 #include <deal.II/lac/la_parallel_vector.h> 43 #include <deal.II/lac/la_vector.h> 44 #include <deal.II/lac/petsc_block_vector.h> 45 #include <deal.II/lac/petsc_vector.h> 46 #include <deal.II/lac/trilinos_parallel_block_vector.h> 47 #include <deal.II/lac/trilinos_tpetra_vector.h> 48 #include <deal.II/lac/trilinos_vector.h> 49 #include <deal.II/lac/vector.h> 51 #include <deal.II/numerics/vector_tools.h> 59 DEAL_II_NAMESPACE_OPEN
62 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
67 , n_shape_functions(fe.dofs_per_cell)
69 , local_dof_indices(fe.dofs_per_cell)
70 , local_dof_values(fe.dofs_per_cell)
75 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
86 template <
int dim,
int spacedim,
typename DoFHandlerType,
typename VectorType>
89 const unsigned int qpoint,
90 const unsigned int shape_nr)
92 Assert(qpoint * n_shape_functions + shape_nr < shape_values.size(),
95 shape_values.size()));
96 return shape_values[qpoint * n_shape_functions + shape_nr];
100 template <
int dim,
int spacedim,
typename DoFHandlerType,
typename VectorType>
103 derivative(
const unsigned int qpoint,
const unsigned int shape_nr)
const 105 Assert(qpoint * n_shape_functions + shape_nr < shape_derivatives.size(),
108 shape_derivatives.size()));
109 return shape_derivatives[qpoint * n_shape_functions + shape_nr];
114 template <
int dim,
int spacedim,
typename DoFHandlerType,
typename VectorType>
117 derivative(
const unsigned int qpoint,
const unsigned int shape_nr)
119 Assert(qpoint * n_shape_functions + shape_nr < shape_derivatives.size(),
122 shape_derivatives.size()));
123 return shape_derivatives[qpoint * n_shape_functions + shape_nr];
127 template <
int dim,
int spacedim,
typename DoFHandlerType,
typename VectorType>
131 const unsigned int shape_nr)
const 133 Assert(qpoint * n_shape_functions + shape_nr <
134 shape_second_derivatives.size(),
137 shape_second_derivatives.size()));
138 return shape_second_derivatives[qpoint * n_shape_functions + shape_nr];
143 template <
int dim,
int spacedim,
typename DoFHandlerType,
typename VectorType>
148 Assert(qpoint * n_shape_functions + shape_nr <
149 shape_second_derivatives.size(),
152 shape_second_derivatives.size()));
153 return shape_second_derivatives[qpoint * n_shape_functions + shape_nr];
157 template <
int dim,
int spacedim,
typename DoFHandlerType,
typename VectorType>
162 Assert(qpoint * n_shape_functions + shape_nr < shape_third_derivatives.size(),
165 shape_third_derivatives.size()));
166 return shape_third_derivatives[qpoint * n_shape_functions + shape_nr];
171 template <
int dim,
int spacedim,
typename DoFHandlerType,
typename VectorType>
176 Assert(qpoint * n_shape_functions + shape_nr < shape_third_derivatives.size(),
179 shape_third_derivatives.size()));
180 return shape_third_derivatives[qpoint * n_shape_functions + shape_nr];
184 template <
int dim,
int spacedim,
typename DoFHandlerType,
typename VectorType>
188 const unsigned int shape_nr)
const 190 Assert(qpoint * n_shape_functions + shape_nr <
191 shape_fourth_derivatives.size(),
194 shape_fourth_derivatives.size()));
195 return shape_fourth_derivatives[qpoint * n_shape_functions + shape_nr];
200 template <
int dim,
int spacedim,
typename DoFHandlerType,
typename VectorType>
205 Assert(qpoint * n_shape_functions + shape_nr <
206 shape_fourth_derivatives.size(),
209 shape_fourth_derivatives.size()));
210 return shape_fourth_derivatives[qpoint * n_shape_functions + shape_nr];
221 get_vertex_quadrature()
224 if (quad.
size() == 0)
227 for (
unsigned int i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
237 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
244 , fe_mask(mask.size() ?
251 get_vertex_quadrature<dim>(),
254 unsigned int size = 0;
255 for (
unsigned int i = 0; i < fe_mask.
size(); ++i)
264 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
267 : euler_vector(mapping.euler_vector)
268 , euler_dof_handler(mapping.euler_dof_handler)
269 , fe_mask(mapping.fe_mask)
270 , fe_to_real(mapping.fe_to_real)
271 , fe_values(mapping.euler_dof_handler->get_fe(),
272 get_vertex_quadrature<dim>(),
278 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
279 inline const double &
281 const unsigned int qpoint,
282 const unsigned int shape_nr)
const 293 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
303 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
319 std::vector<Vector<typename VectorType::value_type>> values(
331 for (
unsigned int i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
332 for (
unsigned int j = 0; j <
fe_to_real.size(); ++j)
341 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
350 const unsigned int n_points = unit_points.size();
352 for (
unsigned int point = 0; point < n_points; ++point)
354 if (data.shape_values.size() != 0)
355 for (
unsigned int i = 0; i < data.n_shape_functions; ++i)
356 data.shape(point, i) = fe->shape_value(i, unit_points[point]);
358 if (data.shape_derivatives.size() != 0)
359 for (
unsigned int i = 0; i < data.n_shape_functions; ++i)
360 data.derivative(point, i) = fe->shape_grad(i, unit_points[point]);
362 if (data.shape_second_derivatives.size() != 0)
363 for (
unsigned int i = 0; i < data.n_shape_functions; ++i)
364 data.second_derivative(point, i) =
365 fe->shape_grad_grad(i, unit_points[point]);
367 if (data.shape_third_derivatives.size() != 0)
368 for (
unsigned int i = 0; i < data.n_shape_functions; ++i)
369 data.third_derivative(point, i) =
370 fe->shape_3rd_derivative(i, unit_points[point]);
372 if (data.shape_fourth_derivatives.size() != 0)
373 for (
unsigned int i = 0; i < data.n_shape_functions; ++i)
374 data.fourth_derivative(point, i) =
375 fe->shape_4th_derivative(i, unit_points[point]);
380 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
392 for (
unsigned int i = 0; i < 5; ++i)
435 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
440 const unsigned int n_original_q_points,
441 InternalData & data)
const 448 const unsigned int n_q_points = q.
size();
453 data.shape_values.resize(data.n_shape_functions * n_q_points);
455 if (data.update_each &
459 data.shape_derivatives.resize(data.n_shape_functions * n_q_points);
462 data.covariant.resize(n_original_q_points);
465 data.contravariant.resize(n_original_q_points);
468 data.volume_elements.resize(n_original_q_points);
470 if (data.update_each &
472 data.shape_second_derivatives.resize(data.n_shape_functions * n_q_points);
476 data.shape_third_derivatives.resize(data.n_shape_functions * n_q_points);
480 data.shape_fourth_derivatives.resize(data.n_shape_functions * n_q_points);
486 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
491 const unsigned int n_original_q_points,
492 InternalData & data)
const 494 compute_data(update_flags, q, n_original_q_points, data);
504 for (
unsigned int i = 0; i < data.unit_tangentials.size(); ++i)
505 data.unit_tangentials[i].resize(n_original_q_points);
511 static const int tangential_orientation[4] = {-1, 1, 1, -1};
512 for (
unsigned int i = 0; i < GeometryInfo<dim>::faces_per_cell;
516 tang[1 - i / 2] = tangential_orientation[i];
517 std::fill(data.unit_tangentials[i].begin(),
518 data.unit_tangentials[i].end(),
524 for (
unsigned int i = 0; i < GeometryInfo<dim>::faces_per_cell;
529 const unsigned int nd =
537 tang1[(nd + 1) % dim] =
542 tang2[(nd + 2) % dim] = 1.;
547 std::fill(data.unit_tangentials[i].begin(),
548 data.unit_tangentials[i].end(),
563 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
564 typename std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
569 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
572 this->compute_data(update_flags, quadrature, quadrature.
size(), data);
579 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
580 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
585 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
589 this->compute_face_data(update_flags, q, quadrature.
size(), data);
595 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
596 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
601 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
605 this->compute_face_data(update_flags, q, quadrature.
size(), data);
614 namespace MappingFEFieldImplementation
627 typename DoFHandlerType>
629 maybe_compute_q_points(
630 const typename ::QProjector<dim>::DataSetDescriptor data_set,
643 for (
unsigned int point = 0; point < quadrature_points.size();
647 const double * shape = &data.shape(point + data_set, 0);
649 for (
unsigned int k = 0; k < data.n_shape_functions; ++k)
654 data.local_dof_values[k] * shape[k];
657 quadrature_points[
point] = result;
672 typename DoFHandlerType>
674 maybe_update_Jacobians(
676 const typename ::QProjector<dim>::DataSetDescriptor data_set,
693 const unsigned int n_q_points = data.contravariant.size();
697 for (
unsigned int point = 0;
point < n_q_points; ++
point)
700 &data.derivative(point + data_set, 0);
704 for (
unsigned int k = 0; k < data.n_shape_functions; ++k)
706 unsigned int comp_k =
710 data.local_dof_values[k] * data_derv[k];
714 for (
unsigned int i = 0; i < spacedim; ++i)
716 data.contravariant[
point][i] = result[i];
726 for (
unsigned int point = 0;
point < data.contravariant.size();
728 data.covariant[point] =
729 (data.contravariant[point]).covariant_form();
736 for (
unsigned int point = 0;
point < data.contravariant.size();
738 data.volume_elements[point] =
739 data.contravariant[point].determinant();
752 typename DoFHandlerType>
754 maybe_update_jacobian_grads(
756 const typename ::QProjector<dim>::DataSetDescriptor data_set,
768 const unsigned int n_q_points = jacobian_grads.size();
772 for (
unsigned int point = 0;
point < n_q_points; ++
point)
775 &data.second_derivative(point + data_set, 0);
779 for (
unsigned int k = 0; k < data.n_shape_functions; ++k)
781 unsigned int comp_k =
784 for (
unsigned int j = 0; j < dim; ++j)
785 for (
unsigned int l = 0;
l < dim; ++
l)
787 (second[k][j][l] * data.local_dof_values[k]);
792 for (
unsigned int i = 0; i < spacedim; ++i)
793 for (
unsigned int j = 0; j < dim; ++j)
794 for (
unsigned int l = 0;
l < dim; ++
l)
795 jacobian_grads[point][i][j][l] = result[i][j][l];
810 typename DoFHandlerType>
812 maybe_update_jacobian_pushed_forward_grads(
814 const typename ::QProjector<dim>::DataSetDescriptor data_set,
826 const unsigned int n_q_points =
827 jacobian_pushed_forward_grads.size();
831 double tmp[spacedim][spacedim][spacedim];
832 for (
unsigned int point = 0;
point < n_q_points; ++
point)
835 &data.second_derivative(point + data_set, 0);
839 for (
unsigned int k = 0; k < data.n_shape_functions; ++k)
841 unsigned int comp_k =
844 for (
unsigned int j = 0; j < dim; ++j)
845 for (
unsigned int l = 0;
l < dim; ++
l)
847 (second[k][j][l] * data.local_dof_values[k]);
851 for (
unsigned int i = 0; i < spacedim; ++i)
852 for (
unsigned int j = 0; j < spacedim; ++j)
853 for (
unsigned int l = 0;
l < dim; ++
l)
856 result[i][0][
l] * data.covariant[
point][j][0];
857 for (
unsigned int jr = 1; jr < dim; ++jr)
859 tmp[i][j][
l] += result[i][jr][
l] *
860 data.covariant[
point][j][jr];
865 for (
unsigned int i = 0; i < spacedim; ++i)
866 for (
unsigned int j = 0; j < spacedim; ++j)
867 for (
unsigned int l = 0;
l < spacedim; ++
l)
869 jacobian_pushed_forward_grads[
point][i][j][
l] =
870 tmp[i][j][0] * data.covariant[
point][
l][0];
871 for (
unsigned int lr = 1; lr < dim; ++lr)
873 jacobian_pushed_forward_grads[
point][i][j][
l] +=
874 tmp[i][j][lr] * data.covariant[
point][
l][lr];
891 typename DoFHandlerType>
893 maybe_update_jacobian_2nd_derivatives(
895 const typename ::QProjector<dim>::DataSetDescriptor data_set,
907 const unsigned int n_q_points = jacobian_2nd_derivatives.size();
911 for (
unsigned int point = 0;
point < n_q_points; ++
point)
914 &data.third_derivative(point + data_set, 0);
918 for (
unsigned int k = 0; k < data.n_shape_functions; ++k)
920 unsigned int comp_k =
923 for (
unsigned int j = 0; j < dim; ++j)
924 for (
unsigned int l = 0;
l < dim; ++
l)
925 for (
unsigned int m = 0; m < dim; ++m)
928 data.local_dof_values[k]);
933 for (
unsigned int i = 0; i < spacedim; ++i)
934 for (
unsigned int j = 0; j < dim; ++j)
935 for (
unsigned int l = 0;
l < dim; ++
l)
936 for (
unsigned int m = 0; m < dim; ++m)
937 jacobian_2nd_derivatives[point][i][j][l][m] =
954 typename DoFHandlerType>
956 maybe_update_jacobian_pushed_forward_2nd_derivatives(
958 const typename ::QProjector<dim>::DataSetDescriptor data_set,
966 &jacobian_pushed_forward_2nd_derivatives)
971 const unsigned int n_q_points =
972 jacobian_pushed_forward_2nd_derivatives.size();
976 double tmp[spacedim][spacedim][spacedim][spacedim];
977 for (
unsigned int point = 0;
point < n_q_points; ++
point)
980 &data.third_derivative(point + data_set, 0);
984 for (
unsigned int k = 0; k < data.n_shape_functions; ++k)
986 unsigned int comp_k =
989 for (
unsigned int j = 0; j < dim; ++j)
990 for (
unsigned int l = 0;
l < dim; ++
l)
991 for (
unsigned int m = 0; m < dim; ++m)
994 data.local_dof_values[k]);
998 for (
unsigned int i = 0; i < spacedim; ++i)
999 for (
unsigned int j = 0; j < spacedim; ++j)
1000 for (
unsigned int l = 0;
l < dim; ++
l)
1001 for (
unsigned int m = 0; m < dim; ++m)
1003 jacobian_pushed_forward_2nd_derivatives
1005 result[i][0][
l][m] *
1006 data.covariant[
point][j][0];
1007 for (
unsigned int jr = 1; jr < dim; ++jr)
1008 jacobian_pushed_forward_2nd_derivatives[point]
1011 result[i][jr][l][m] *
1012 data.covariant[point][j][jr];
1016 for (
unsigned int i = 0; i < spacedim; ++i)
1017 for (
unsigned int j = 0; j < spacedim; ++j)
1018 for (
unsigned int l = 0;
l < spacedim; ++
l)
1019 for (
unsigned int m = 0; m < dim; ++m)
1022 jacobian_pushed_forward_2nd_derivatives[
point]
1025 data.covariant[
point][
l][0];
1026 for (
unsigned int lr = 1; lr < dim; ++lr)
1028 jacobian_pushed_forward_2nd_derivatives
1029 [point][i][j][lr][m] *
1030 data.covariant[point][l][lr];
1034 for (
unsigned int i = 0; i < spacedim; ++i)
1035 for (
unsigned int j = 0; j < spacedim; ++j)
1036 for (
unsigned int l = 0;
l < spacedim; ++
l)
1037 for (
unsigned int m = 0; m < spacedim; ++m)
1039 jacobian_pushed_forward_2nd_derivatives
1041 tmp[i][j][
l][0] * data.covariant[
point][m][0];
1042 for (
unsigned int mr = 1; mr < dim; ++mr)
1043 jacobian_pushed_forward_2nd_derivatives[point]
1047 data.covariant[point][m][mr];
1062 typename VectorType,
1063 typename DoFHandlerType>
1065 maybe_update_jacobian_3rd_derivatives(
1067 const typename ::QProjector<dim>::DataSetDescriptor data_set,
1070 InternalData & data,
1073 const std::vector<unsigned int> &
fe_to_real,
1076 const UpdateFlags update_flags = data.update_each;
1079 const unsigned int n_q_points = jacobian_3rd_derivatives.size();
1083 for (
unsigned int point = 0;
point < n_q_points; ++
point)
1086 &data.fourth_derivative(point + data_set, 0);
1090 for (
unsigned int k = 0; k < data.n_shape_functions; ++k)
1092 unsigned int comp_k =
1094 if (fe_mask[comp_k])
1095 for (
unsigned int j = 0; j < dim; ++j)
1096 for (
unsigned int l = 0;
l < dim; ++
l)
1097 for (
unsigned int m = 0; m < dim; ++m)
1098 for (
unsigned int n = 0; n < dim; ++n)
1100 (fourth[k][j][l][m][n] *
1101 data.local_dof_values[k]);
1107 for (
unsigned int i = 0; i < spacedim; ++i)
1108 for (
unsigned int j = 0; j < dim; ++j)
1109 for (
unsigned int l = 0;
l < dim; ++
l)
1110 for (
unsigned int m = 0; m < dim; ++m)
1111 for (
unsigned int n = 0; n < dim; ++n)
1112 jacobian_3rd_derivatives[point][i][j][l][m][n] =
1113 result[i][j][l][m][n];
1128 typename VectorType,
1129 typename DoFHandlerType>
1131 maybe_update_jacobian_pushed_forward_3rd_derivatives(
1133 const typename ::QProjector<dim>::DataSetDescriptor data_set,
1136 InternalData & data,
1139 const std::vector<unsigned int> &
fe_to_real,
1141 &jacobian_pushed_forward_3rd_derivatives)
1143 const UpdateFlags update_flags = data.update_each;
1146 const unsigned int n_q_points =
1147 jacobian_pushed_forward_3rd_derivatives.size();
1151 double tmp[spacedim][spacedim][spacedim][spacedim][spacedim];
1152 for (
unsigned int point = 0;
point < n_q_points; ++
point)
1155 &data.fourth_derivative(point + data_set, 0);
1159 for (
unsigned int k = 0; k < data.n_shape_functions; ++k)
1161 unsigned int comp_k =
1163 if (fe_mask[comp_k])
1164 for (
unsigned int j = 0; j < dim; ++j)
1165 for (
unsigned int l = 0;
l < dim; ++
l)
1166 for (
unsigned int m = 0; m < dim; ++m)
1167 for (
unsigned int n = 0; n < dim; ++n)
1169 (fourth[k][j][l][m][n] *
1170 data.local_dof_values[k]);
1174 for (
unsigned int i = 0; i < spacedim; ++i)
1175 for (
unsigned int j = 0; j < spacedim; ++j)
1176 for (
unsigned int l = 0;
l < dim; ++
l)
1177 for (
unsigned int m = 0; m < dim; ++m)
1178 for (
unsigned int n = 0; n < dim; ++n)
1180 tmp[i][j][
l][m][n] =
1181 result[i][0][
l][m][n] *
1182 data.covariant[
point][j][0];
1183 for (
unsigned int jr = 1; jr < dim; ++jr)
1184 tmp[i][j][l][m][n] +=
1185 result[i][jr][l][m][n] *
1186 data.covariant[point][j][jr];
1190 for (
unsigned int i = 0; i < spacedim; ++i)
1191 for (
unsigned int j = 0; j < spacedim; ++j)
1192 for (
unsigned int l = 0;
l < spacedim; ++
l)
1193 for (
unsigned int m = 0; m < dim; ++m)
1194 for (
unsigned int n = 0; n < dim; ++n)
1196 jacobian_pushed_forward_3rd_derivatives
1198 tmp[i][j][0][m][n] *
1199 data.covariant[
point][
l][0];
1200 for (
unsigned int lr = 1; lr < dim; ++lr)
1201 jacobian_pushed_forward_3rd_derivatives
1202 [point][i][j][l][m][n] +=
1203 tmp[i][j][lr][m][n] *
1204 data.covariant[point][l][lr];
1208 for (
unsigned int i = 0; i < spacedim; ++i)
1209 for (
unsigned int j = 0; j < spacedim; ++j)
1210 for (
unsigned int l = 0;
l < spacedim; ++
l)
1211 for (
unsigned int m = 0; m < spacedim; ++m)
1212 for (
unsigned int n = 0; n < dim; ++n)
1214 tmp[i][j][
l][m][n] =
1215 jacobian_pushed_forward_3rd_derivatives
1217 data.covariant[
point][m][0];
1218 for (
unsigned int mr = 1; mr < dim; ++mr)
1219 tmp[i][j][l][m][n] +=
1220 jacobian_pushed_forward_3rd_derivatives
1221 [point][i][j][l][mr][n] *
1222 data.covariant[point][m][mr];
1226 for (
unsigned int i = 0; i < spacedim; ++i)
1227 for (
unsigned int j = 0; j < spacedim; ++j)
1228 for (
unsigned int l = 0;
l < spacedim; ++
l)
1229 for (
unsigned int m = 0; m < spacedim; ++m)
1230 for (
unsigned int n = 0; n < spacedim; ++n)
1232 jacobian_pushed_forward_3rd_derivatives
1234 tmp[i][j][
l][m][0] *
1235 data.covariant[
point][n][0];
1236 for (
unsigned int nr = 1; nr < dim; ++nr)
1237 jacobian_pushed_forward_3rd_derivatives
1238 [point][i][j][l][m][n] +=
1239 tmp[i][j][l][m][nr] *
1240 data.covariant[point][n][nr];
1259 typename VectorType,
1260 typename DoFHandlerType>
1262 maybe_compute_face_data(
1263 const ::Mapping<dim, spacedim> &mapping,
1264 const typename ::Triangulation<dim, spacedim>::cell_iterator
1266 const unsigned int face_no,
1267 const unsigned int subface_no,
1268 const std::vector<double> &weights,
1275 const UpdateFlags update_flags = data.update_each;
1279 const unsigned int n_q_points = output_data.
boundary_forms.size();
1288 for (
unsigned int d = 0;
d != dim - 1; ++
d)
1291 data.unit_tangentials.size(),
1294 data.aux[d].size() <=
1296 .unit_tangentials[face_no +
1304 .unit_tangentials[face_no +
1308 make_array_view(data.aux[d]));
1313 if (dim == spacedim)
1315 for (
unsigned int i = 0; i < n_q_points; ++i)
1324 (face_no == 0 ? -1 : +1);
1328 cross_product_2d(data.aux[0][i]);
1332 cross_product_3d(data.aux[0][i], data.aux[1][i]);
1348 for (
unsigned int point = 0;
point < n_q_points; ++
point)
1354 data.contravariant[
point].transpose()[0];
1356 (face_no == 0 ? -1. : +1.) *
1363 data.contravariant[
point].transpose();
1366 cross_product_3d(DX_t[0], DX_t[1]);
1367 cell_normal /= cell_normal.
norm();
1372 cross_product_3d(data.aux[0][point], cell_normal);
1388 const double area_ratio =
1390 cell->subface_case(face_no), subface_no);
1402 for (
unsigned int point = 0;
point < n_q_points; ++
point)
1403 output_data.
jacobians[point] = data.contravariant[point];
1406 for (
unsigned int point = 0;
point < n_q_points; ++
point)
1408 data.covariant[point].transpose();
1420 typename VectorType,
1421 typename DoFHandlerType>
1423 do_fill_fe_face_values(
1424 const ::Mapping<dim, spacedim> &mapping,
1425 const typename ::Triangulation<dim, spacedim>::cell_iterator
1427 const unsigned int face_no,
1428 const unsigned int subface_no,
1429 const typename ::QProjector<dim>::DataSetDescriptor data_set,
1433 InternalData & data,
1436 const std::vector<unsigned int> &
fe_to_real,
1440 maybe_compute_q_points<dim, spacedim, VectorType, DoFHandlerType>(
1448 maybe_update_Jacobians<dim, spacedim, VectorType, DoFHandlerType>(
1451 maybe_update_jacobian_grads<dim, spacedim, VectorType, DoFHandlerType>(
1460 maybe_update_jacobian_pushed_forward_grads<dim,
1472 maybe_update_jacobian_2nd_derivatives<dim,
1484 maybe_update_jacobian_pushed_forward_2nd_derivatives<dim,
1496 maybe_update_jacobian_3rd_derivatives<dim,
1508 maybe_update_jacobian_pushed_forward_3rd_derivatives<dim,
1520 maybe_compute_face_data<dim, spacedim, VectorType, DoFHandlerType>(
1536 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
1548 Assert(dynamic_cast<const InternalData *>(&internal_data) !=
nullptr,
1550 const InternalData &data =
static_cast<const InternalData &
>(internal_data);
1552 const unsigned int n_q_points = quadrature.
size();
1558 internal::MappingFEFieldImplementation::
1559 maybe_compute_q_points<dim, spacedim, VectorType, DoFHandlerType>(
1567 internal::MappingFEFieldImplementation::
1568 maybe_update_Jacobians<dim, spacedim, VectorType, DoFHandlerType>(
1576 const UpdateFlags update_flags = data.update_each;
1577 const std::vector<double> &weights = quadrature.
get_weights();
1593 for (
unsigned int point = 0;
point < n_q_points; ++
point)
1595 if (dim == spacedim)
1597 const double det = data.contravariant[
point].determinant();
1605 1e-12 * Utilities::fixed_power<dim>(
1606 cell->diameter() / std::sqrt(
double(dim))),
1608 cell->center(), det,
point)));
1617 for (
unsigned int i = 0; i < spacedim; ++i)
1618 for (
unsigned int j = 0; j < dim; ++j)
1619 DX_t[j][i] = data.contravariant[point][i][j];
1622 for (
unsigned int i = 0; i < dim; ++i)
1623 for (
unsigned int j = 0; j < dim; ++j)
1624 G[i][j] = DX_t[i] * DX_t[j];
1639 Assert(spacedim - dim == 1,
1641 "There is no cell normal in codim 2."));
1645 cross_product_2d(-DX_t[0]);
1648 cross_product_3d(DX_t[0], DX_t[1]);
1653 if (cell->direction_flag() ==
false)
1666 for (
unsigned int point = 0;
point < n_q_points; ++
point)
1667 output_data.
jacobians[point] = data.contravariant[point];
1675 for (
unsigned int point = 0;
point < n_q_points; ++
point)
1677 data.covariant[point].transpose();
1681 internal::MappingFEFieldImplementation::
1682 maybe_update_jacobian_grads<dim, spacedim, VectorType, DoFHandlerType>(
1693 internal::MappingFEFieldImplementation::
1694 maybe_update_jacobian_pushed_forward_grads<dim,
1707 internal::MappingFEFieldImplementation::maybe_update_jacobian_2nd_derivatives<
1711 DoFHandlerType>(cell_similarity,
1720 internal::MappingFEFieldImplementation::
1721 maybe_update_jacobian_pushed_forward_2nd_derivatives<dim,
1734 internal::MappingFEFieldImplementation::maybe_update_jacobian_3rd_derivatives<
1738 DoFHandlerType>(cell_similarity,
1748 internal::MappingFEFieldImplementation::
1749 maybe_update_jacobian_pushed_forward_3rd_derivatives<dim,
1761 return updated_cell_similarity;
1766 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
1770 const unsigned int face_no,
1778 Assert(dynamic_cast<const InternalData *>(&internal_data) !=
nullptr,
1780 const InternalData &data =
static_cast<const InternalData &
>(internal_data);
1784 internal::MappingFEFieldImplementation::
1785 do_fill_fe_face_values<dim, spacedim, VectorType, DoFHandlerType>(
1791 cell->face_orientation(face_no),
1792 cell->face_flip(face_no),
1793 cell->face_rotation(face_no),
1804 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
1809 const unsigned int face_no,
1810 const unsigned int subface_no,
1818 Assert(dynamic_cast<const InternalData *>(&internal_data) !=
nullptr,
1820 const InternalData &data =
static_cast<const InternalData &
>(internal_data);
1824 internal::MappingFEFieldImplementation::
1825 do_fill_fe_face_values<dim, spacedim, VectorType, DoFHandlerType>(
1832 cell->face_orientation(
1834 cell->face_flip(face_no),
1835 cell->face_rotation(face_no),
1837 cell->subface_case(face_no)),
1849 namespace MappingFEFieldImplementation
1856 typename VectorType,
1857 typename DoFHandlerType>
1868 MappingFEField<dim, spacedim, VectorType, DoFHandlerType>::
1869 InternalData *
>(&mapping_data) !=
nullptr),
1871 const typename ::MappingFEField<dim,
1874 DoFHandlerType>::InternalData
1875 &data =
static_cast< 1877 MappingFEField<dim, spacedim, VectorType, DoFHandlerType>::
1878 InternalData &
>(mapping_data);
1880 switch (mapping_type)
1887 "update_contravariant_transformation"));
1889 for (
unsigned int i = 0; i < output.size(); ++i)
1891 apply_transformation(data.contravariant[i], input[i]);
1901 "update_contravariant_transformation"));
1905 "update_volume_elements"));
1907 for (
unsigned int i = 0; i < output.size(); ++i)
1910 apply_transformation(data.contravariant[i], input[i]);
1911 output[i] /= data.volume_elements[i];
1925 "update_contravariant_transformation"));
1927 for (
unsigned int i = 0; i < output.size(); ++i)
1928 output[i] = apply_transformation(data.covariant[i], input[i]);
1942 typename VectorType,
1943 typename DoFHandlerType>
1945 transform_differential_forms(
1954 MappingFEField<dim, spacedim, VectorType, DoFHandlerType>::
1955 InternalData *
>(&mapping_data) !=
nullptr),
1957 const typename ::MappingFEField<dim,
1960 DoFHandlerType>::InternalData
1961 &data =
static_cast< 1963 MappingFEField<dim, spacedim, VectorType, DoFHandlerType>::
1964 InternalData &
>(mapping_data);
1966 switch (mapping_type)
1973 "update_contravariant_transformation"));
1975 for (
unsigned int i = 0; i < output.size(); ++i)
1976 output[i] = apply_transformation(data.covariant[i], input[i]);
1990 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
2000 internal::MappingFEFieldImplementation::
2001 transform_fields<dim, spacedim, 1, VectorType, DoFHandlerType>(input,
2009 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
2019 internal::MappingFEFieldImplementation::
2020 transform_differential_forms<dim, spacedim, 1, VectorType, DoFHandlerType>(
2021 input, mapping_type, mapping_data, output);
2026 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
2044 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
2053 Assert(dynamic_cast<const InternalData *>(&mapping_data) !=
nullptr,
2057 switch (mapping_type)
2063 "update_covariant_transformation"));
2065 for (
unsigned int q = 0; q < output.
size(); ++q)
2066 for (
unsigned int i = 0; i < spacedim; ++i)
2067 for (
unsigned int j = 0; j < spacedim; ++j)
2068 for (
unsigned int k = 0; k < spacedim; ++k)
2070 output[q][i][j][k] = data.
covariant[q][j][0] *
2073 for (
unsigned int J = 0; J < dim; ++J)
2075 const unsigned int K0 = (0 == J) ? 1 : 0;
2076 for (
unsigned int K = K0; K < dim; ++K)
2077 output[q][i][j][k] += data.
covariant[q][j][J] *
2092 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
2110 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
2121 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> mdata(
2123 Assert(dynamic_cast<InternalData *>(mdata.get()) !=
nullptr,
2132 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
2141 unsigned int comp_i =
2143 if (fe_mask[comp_i])
2153 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
2174 for (
unsigned int d = 0; d < dim; ++d)
2175 initial_p_unit[d] = 0.5;
2188 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> mdata(
2189 get_data(update_flags, point_quadrature));
2190 Assert(dynamic_cast<InternalData *>(mdata.get()) !=
nullptr,
2198 dynamic_cast<InternalData &>(*mdata));
2202 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
2211 const unsigned int n_shapes = mdata.
shape_values.size();
2231 const double eps = 1.e-12 * cell->diameter();
2232 const unsigned int newton_iteration_limit = 20;
2233 unsigned int newton_iteration = 0;
2242 unsigned int comp_k =
2244 if (fe_mask[comp_k])
2245 for (
unsigned int j = 0; j < dim; ++j)
2249 for (
unsigned int j = 0; j < dim; ++j)
2251 f[j] = DF[j] * p_minus_F;
2252 for (
unsigned int l = 0; l < dim; ++l)
2253 df[j][l] = -DF[j] * DF[l];
2259 double step_length = 1;
2272 for (
unsigned int i = 0; i < dim; ++i)
2273 p_unit_trial[i] -= step_length * delta[i];
2283 if (f_trial.
norm() < p_minus_F.
norm())
2285 p_real = p_real_trial;
2286 p_unit = p_unit_trial;
2287 p_minus_F = f_trial;
2290 else if (step_length > 0.05)
2297 if (newton_iteration > newton_iteration_limit)
2314 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
2322 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
2327 return this->fe_mask;
2331 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
2332 std::unique_ptr<Mapping<dim, spacedim>>
2335 return std_cxx14::make_unique<
2340 template <
int dim,
int spacedim,
typename VectorType,
typename DoFHandlerType>
2353 dof_cell->get_dof_indices(data.local_dof_indices);
2355 for (
unsigned int i = 0; i < data.local_dof_values.size(); ++i)
2356 data.local_dof_values[i] =
2357 internal::ElementAccess<VectorType>::get(*
euler_vector,
2358 data.local_dof_indices[i]);
2362 #define SPLIT_INSTANTIATIONS_COUNT 2 2363 #ifndef SPLIT_INSTANTIATIONS_INDEX 2364 # define SPLIT_INSTANTIATIONS_INDEX 0 2366 #include "mapping_fe_field.inst" 2369 DEAL_II_NAMESPACE_CLOSE
Transformed quadrature weights.
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
static const unsigned int invalid_unsigned_int
#define AssertDimension(dim1, dim2)
Number determinant(const SymmetricTensor< 2, dim, Number > &)
Contravariant transformation.
const std::vector< Point< dim > > & get_points() const
unsigned int n_shape_functions
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
SmartPointer< const VectorType, MappingFEField< dim, spacedim, VectorType, DoFHandlerType > > euler_vector
const std::vector< double > & get_weights() const
virtual void compute_shapes_virtual(const std::vector< Point< dim >> &unit_points, typename MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data) const
Outer normal vector, not normalized.
std::vector< unsigned int > fe_to_real
Determinant of the Jacobian.
Transformed quadrature points.
#define AssertThrow(cond, exc)
numbers::NumberTraits< Number >::real_type norm() const
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static DataSetDescriptor cell()
Point< spacedim > do_transform_unit_to_real_cell(const InternalData &mdata) const
unsigned int get_degree() const
FEValues< dim, spacedim > fe_values
SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
std::vector< Tensor< 1, dim > > shape_derivatives
static ::ExceptionBase & ExcMessage(std::string arg1)
const double & shape(const unsigned int qpoint, const unsigned int shape_nr) const
#define Assert(cond, exc)
Threads::Mutex fe_values_mutex
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
virtual void transform(const ArrayView< const Tensor< 1, dim >> &input, const MappingType type, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim >> &output) const override
numbers::NumberTraits< Number >::real_type norm_square() const
unsigned int size() const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
Gradient of volume element.
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
std::vector< double > local_dof_values
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
const Tensor< 3, dim > & third_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
virtual std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
unsigned int size() const
std::vector< double > shape_values
friend class MappingFEField
static Point< dim > project_to_unit_cell(const Point< dim > &p)
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
const Tensor< 2, dim > & second_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
const Tensor< 4, dim > & fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
SmartPointer< const DoFHandlerType, MappingFEField< dim, spacedim, VectorType, DoFHandlerType > > euler_dof_handler
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
static DataSetDescriptor subface(const unsigned int face_no, const unsigned int subface_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points, const internal::SubfaceCase< dim > ref_case=internal::SubfaceCase< dim >::case_isotropic)
const Tensor< 1, dim > & derivative(const unsigned int qpoint, const unsigned int shape_nr) const
static ::ExceptionBase & ExcNotImplemented()
ComponentMask get_component_mask() const
void update_internal_dofs(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data) const
typename ActiveSelector::cell_iterator cell_iterator
virtual bool preserves_vertex_locations() const override
static ::ExceptionBase & ExcInactiveCell()
Point< dim > do_transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_p_unit, InternalData &mdata) const
static double subface_ratio(const internal::SubfaceCase< dim > &subface_case, const unsigned int subface_no)
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
Covariant transformation.
virtual std::size_t memory_consumption() const override
InternalData(const FiniteElement< dim, spacedim > &fe, const ComponentMask &mask)
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
static ::ExceptionBase & ExcInternalError()
static DataSetDescriptor face(const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)