48#include <boost/container/small_vector.hpp>
60 template <
class VectorType>
62 const VectorType & vector,
79 template <
int dim,
int spacedim>
80 inline std::vector<unsigned int>
83 std::vector<unsigned int> shape_function_to_row_table(
92 unsigned int nth_nonzero_component = 0;
97 row + nth_nonzero_component;
98 ++nth_nonzero_component;
103 return shape_function_to_row_table;
110 template <
typename Number,
typename T =
void>
125 template <
typename Number>
128 typename
std::enable_if<
129 Differentiation::AD::is_ad_number<Number>::value>::type>
132 value(
const Number & )
144 template <
int dim,
int spacedim>
146 const unsigned int component)
147 : fe_values(&fe_values)
148 , component(component)
149 , shape_function_data(this->fe_values->fe->n_dofs_per_cell())
157 const std::vector<unsigned int> shape_function_to_row_table =
164 if (is_primitive ==
true)
181 template <
int dim,
int spacedim>
184 , component(
numbers::invalid_unsigned_int)
189 template <
int dim,
int spacedim>
191 const unsigned int first_vector_component)
192 : fe_values(&fe_values)
193 , first_vector_component(first_vector_component)
194 , shape_function_data(this->fe_values->fe->n_dofs_per_cell())
202 const std::vector<unsigned int> shape_function_to_row_table =
205 for (
unsigned int d = 0; d < spacedim; ++d)
213 if (is_primitive ==
true)
223 shape_function_to_row_table[i * fe.
n_components() + component];
232 unsigned int n_nonzero_components = 0;
233 for (
unsigned int d = 0; d < spacedim; ++d)
236 ++n_nonzero_components;
238 if (n_nonzero_components == 0)
240 else if (n_nonzero_components > 1)
244 for (
unsigned int d = 0; d < spacedim; ++d)
246 .is_nonzero_shape_function_component[d] ==
true)
259 template <
int dim,
int spacedim>
262 , first_vector_component(
numbers::invalid_unsigned_int)
267 template <
int dim,
int spacedim>
270 const unsigned int first_tensor_component)
271 : fe_values(&fe_values)
272 , first_tensor_component(first_tensor_component)
273 , shape_function_data(this->fe_values->fe->n_dofs_per_cell())
276 Assert(first_tensor_component + (dim * dim + dim) / 2 - 1 <
279 first_tensor_component +
286 const std::vector<unsigned int> shape_function_to_row_table =
289 for (
unsigned int d = 0;
290 d < ::SymmetricTensor<2, dim>::n_independent_components;
293 const unsigned int component = first_tensor_component + d;
299 if (is_primitive ==
true)
300 shape_function_data[i].is_nonzero_shape_function_component[d] =
303 shape_function_data[i].is_nonzero_shape_function_component[d] =
306 if (shape_function_data[i].is_nonzero_shape_function_component[d] ==
308 shape_function_data[i].row_index[d] =
309 shape_function_to_row_table[i * fe.
n_components() + component];
311 shape_function_data[i].row_index[d] =
318 unsigned int n_nonzero_components = 0;
319 for (
unsigned int d = 0;
320 d < ::SymmetricTensor<2, dim>::n_independent_components;
322 if (shape_function_data[i].is_nonzero_shape_function_component[d] ==
324 ++n_nonzero_components;
326 if (n_nonzero_components == 0)
327 shape_function_data[i].single_nonzero_component = -2;
328 else if (n_nonzero_components > 1)
329 shape_function_data[i].single_nonzero_component = -1;
332 for (
unsigned int d = 0;
333 d < ::SymmetricTensor<2, dim>::n_independent_components;
335 if (shape_function_data[i]
336 .is_nonzero_shape_function_component[d] ==
true)
338 shape_function_data[i].single_nonzero_component =
339 shape_function_data[i].row_index[d];
340 shape_function_data[i].single_nonzero_component_index = d;
349 template <
int dim,
int spacedim>
352 , first_tensor_component(
numbers::invalid_unsigned_int)
357 template <
int dim,
int spacedim>
359 const unsigned int first_tensor_component)
360 : fe_values(&fe_values)
361 , first_tensor_component(first_tensor_component)
362 , shape_function_data(this->fe_values->fe->n_dofs_per_cell())
369 const std::vector<unsigned int> shape_function_to_row_table =
372 for (
unsigned int d = 0; d < dim * dim; ++d)
374 const unsigned int component = first_tensor_component + d;
380 if (is_primitive ==
true)
381 shape_function_data[i].is_nonzero_shape_function_component[d] =
384 shape_function_data[i].is_nonzero_shape_function_component[d] =
387 if (shape_function_data[i].is_nonzero_shape_function_component[d] ==
389 shape_function_data[i].row_index[d] =
390 shape_function_to_row_table[i * fe.
n_components() + component];
392 shape_function_data[i].row_index[d] =
399 unsigned int n_nonzero_components = 0;
400 for (
unsigned int d = 0; d < dim * dim; ++d)
401 if (shape_function_data[i].is_nonzero_shape_function_component[d] ==
403 ++n_nonzero_components;
405 if (n_nonzero_components == 0)
406 shape_function_data[i].single_nonzero_component = -2;
407 else if (n_nonzero_components > 1)
408 shape_function_data[i].single_nonzero_component = -1;
411 for (
unsigned int d = 0; d < dim * dim; ++d)
412 if (shape_function_data[i]
413 .is_nonzero_shape_function_component[d] ==
true)
415 shape_function_data[i].single_nonzero_component =
416 shape_function_data[i].row_index[d];
417 shape_function_data[i].single_nonzero_component_index = d;
426 template <
int dim,
int spacedim>
429 , first_tensor_component(
numbers::invalid_unsigned_int)
440 template <
int dim,
int spacedim,
typename Number>
446 &shape_function_data,
449 const unsigned int dofs_per_cell = dof_values.
size();
450 const unsigned int n_quadrature_points = values.size();
452 std::fill(values.begin(),
456 for (
unsigned int shape_function = 0; shape_function < dofs_per_cell;
458 if (shape_function_data[shape_function]
459 .is_nonzero_shape_function_component)
461 const Number &value = dof_values[shape_function];
465 if (::internal::CheckForZero<Number>::value(value) ==
true)
468 const double *shape_value_ptr =
469 &shape_values(shape_function_data[shape_function].row_index, 0);
470 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
472 values[q_point] += value * (*shape_value_ptr++);
480 template <
int order,
int dim,
int spacedim,
typename Number>
486 &shape_function_data,
491 const unsigned int dofs_per_cell = dof_values.
size();
492 const unsigned int n_quadrature_points = derivatives.size();
499 for (
unsigned int shape_function = 0; shape_function < dofs_per_cell;
501 if (shape_function_data[shape_function]
502 .is_nonzero_shape_function_component)
504 const Number &value = dof_values[shape_function];
508 if (::internal::CheckForZero<Number>::value(value) ==
true)
511 const ::Tensor<order, spacedim> *shape_derivative_ptr =
512 &shape_derivatives[shape_function_data[shape_function].row_index]
514 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
516 derivatives[q_point] += value * (*shape_derivative_ptr++);
522 template <
int dim,
int spacedim,
typename Number>
528 &shape_function_data,
532 const unsigned int dofs_per_cell = dof_values.
size();
533 const unsigned int n_quadrature_points = laplacians.size();
541 for (
unsigned int shape_function = 0; shape_function < dofs_per_cell;
543 if (shape_function_data[shape_function]
544 .is_nonzero_shape_function_component)
546 const Number &value = dof_values[shape_function];
550 if (::internal::CheckForZero<Number>::value(value) ==
true)
553 const ::Tensor<2, spacedim> *shape_hessian_ptr =
554 &shape_hessians[shape_function_data[shape_function].row_index][0];
555 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
557 laplacians[q_point] += value *
trace(*shape_hessian_ptr++);
565 template <
int dim,
int spacedim,
typename Number>
571 &shape_function_data,
576 const unsigned int dofs_per_cell = dof_values.
size();
577 const unsigned int n_quadrature_points = values.size();
584 for (
unsigned int shape_function = 0; shape_function < dofs_per_cell;
588 shape_function_data[shape_function].single_nonzero_component;
594 const Number &value = dof_values[shape_function];
598 if (::internal::CheckForZero<Number>::value(value) ==
true)
603 const unsigned int comp = shape_function_data[shape_function]
604 .single_nonzero_component_index;
605 const double *shape_value_ptr = &shape_values(snc, 0);
606 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
608 values[q_point][comp] += value * (*shape_value_ptr++);
611 for (
unsigned int d = 0; d < spacedim; ++d)
612 if (shape_function_data[shape_function]
613 .is_nonzero_shape_function_component[d])
615 const double *shape_value_ptr = &shape_values(
616 shape_function_data[shape_function].row_index[d], 0);
617 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
619 values[q_point][d] += value * (*shape_value_ptr++);
626 template <
int order,
int dim,
int spacedim,
typename Number>
632 &shape_function_data,
637 const unsigned int dofs_per_cell = dof_values.
size();
638 const unsigned int n_quadrature_points = derivatives.size();
646 for (
unsigned int shape_function = 0; shape_function < dofs_per_cell;
650 shape_function_data[shape_function].single_nonzero_component;
656 const Number &value = dof_values[shape_function];
660 if (::internal::CheckForZero<Number>::value(value) ==
true)
665 const unsigned int comp = shape_function_data[shape_function]
666 .single_nonzero_component_index;
667 const ::Tensor<order, spacedim> *shape_derivative_ptr =
668 &shape_derivatives[snc][0];
669 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
671 derivatives[q_point][comp] += value * (*shape_derivative_ptr++);
674 for (
unsigned int d = 0; d < spacedim; ++d)
675 if (shape_function_data[shape_function]
676 .is_nonzero_shape_function_component[d])
678 const ::Tensor<order, spacedim> *shape_derivative_ptr =
679 &shape_derivatives[shape_function_data[shape_function]
681 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
683 derivatives[q_point][d] +=
684 value * (*shape_derivative_ptr++);
691 template <
int dim,
int spacedim,
typename Number>
697 &shape_function_data,
701 &symmetric_gradients)
703 const unsigned int dofs_per_cell = dof_values.
size();
704 const unsigned int n_quadrature_points = symmetric_gradients.size();
707 symmetric_gradients.begin(),
708 symmetric_gradients.end(),
712 for (
unsigned int shape_function = 0; shape_function < dofs_per_cell;
716 shape_function_data[shape_function].single_nonzero_component;
722 const Number &value = dof_values[shape_function];
726 if (::internal::CheckForZero<Number>::value(value) ==
true)
731 const unsigned int comp = shape_function_data[shape_function]
732 .single_nonzero_component_index;
733 const ::Tensor<1, spacedim> *shape_gradient_ptr =
734 &shape_gradients[snc][0];
735 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
737 symmetric_gradients[q_point] +=
739 symmetrize_single_row(comp, *shape_gradient_ptr++));
742 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
747 for (
unsigned int d = 0; d < spacedim; ++d)
748 if (shape_function_data[shape_function]
749 .is_nonzero_shape_function_component[d])
752 shape_gradients[shape_function_data[shape_function]
753 .row_index[d]][q_point];
754 symmetric_gradients[q_point] +=
symmetrize(grad);
761 template <
int dim,
int spacedim,
typename Number>
767 &shape_function_data,
769 template solution_divergence_type<Number>> &divergences)
771 const unsigned int dofs_per_cell = dof_values.
size();
772 const unsigned int n_quadrature_points = divergences.size();
778 spacedim>::template solution_divergence_type<Number>());
780 for (
unsigned int shape_function = 0; shape_function < dofs_per_cell;
784 shape_function_data[shape_function].single_nonzero_component;
790 const Number &value = dof_values[shape_function];
794 if (::internal::CheckForZero<Number>::value(value) ==
true)
799 const unsigned int comp = shape_function_data[shape_function]
800 .single_nonzero_component_index;
801 const ::Tensor<1, spacedim> *shape_gradient_ptr =
802 &shape_gradients[snc][0];
803 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
805 divergences[q_point] += value * (*shape_gradient_ptr++)[comp];
808 for (
unsigned int d = 0; d < spacedim; ++d)
809 if (shape_function_data[shape_function]
810 .is_nonzero_shape_function_component[d])
812 const ::Tensor<1, spacedim> *shape_gradient_ptr =
813 &shape_gradients[shape_function_data[shape_function]
815 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
817 divergences[q_point] += value * (*shape_gradient_ptr++)[d];
824 template <
int dim,
int spacedim,
typename Number>
830 &shape_function_data,
833 typename ::internal::CurlType<spacedim>::type>::type> &curls)
835 const unsigned int dofs_per_cell = dof_values.
size();
836 const unsigned int n_quadrature_points = curls.size();
838 std::fill(curls.begin(),
842 typename ::internal::CurlType<spacedim>::type>::type());
850 "Computing the curl in 1d is not a useful operation"));
856 for (
unsigned int shape_function = 0;
857 shape_function < dofs_per_cell;
860 const int snc = shape_function_data[shape_function]
861 .single_nonzero_component;
867 const Number &value = dof_values[shape_function];
871 if (::internal::CheckForZero<Number>::value(value) ==
877 const ::Tensor<1, spacedim> *shape_gradient_ptr =
878 &shape_gradients[snc][0];
880 Assert(shape_function_data[shape_function]
881 .single_nonzero_component >= 0,
884 if (shape_function_data[shape_function]
885 .single_nonzero_component_index == 0)
886 for (
unsigned int q_point = 0;
887 q_point < n_quadrature_points;
890 value * (*shape_gradient_ptr++)[1];
892 for (
unsigned int q_point = 0;
893 q_point < n_quadrature_points;
896 value * (*shape_gradient_ptr++)[0];
904 if (shape_function_data[shape_function]
905 .is_nonzero_shape_function_component[0])
908 spacedim> *shape_gradient_ptr =
909 &shape_gradients[shape_function_data[shape_function]
912 for (
unsigned int q_point = 0;
913 q_point < n_quadrature_points;
916 value * (*shape_gradient_ptr++)[1];
919 if (shape_function_data[shape_function]
920 .is_nonzero_shape_function_component[1])
923 spacedim> *shape_gradient_ptr =
924 &shape_gradients[shape_function_data[shape_function]
927 for (
unsigned int q_point = 0;
928 q_point < n_quadrature_points;
931 value * (*shape_gradient_ptr++)[0];
940 for (
unsigned int shape_function = 0;
941 shape_function < dofs_per_cell;
944 const int snc = shape_function_data[shape_function]
945 .single_nonzero_component;
951 const Number &value = dof_values[shape_function];
955 if (::internal::CheckForZero<Number>::value(value) ==
961 const ::Tensor<1, spacedim> *shape_gradient_ptr =
962 &shape_gradients[snc][0];
964 switch (shape_function_data[shape_function]
965 .single_nonzero_component_index)
969 for (
unsigned int q_point = 0;
970 q_point < n_quadrature_points;
974 value * (*shape_gradient_ptr)[2];
976 value * (*shape_gradient_ptr++)[1];
984 for (
unsigned int q_point = 0;
985 q_point < n_quadrature_points;
989 value * (*shape_gradient_ptr)[2];
991 value * (*shape_gradient_ptr++)[0];
999 for (
unsigned int q_point = 0;
1000 q_point < n_quadrature_points;
1003 curls[q_point][0] +=
1004 value * (*shape_gradient_ptr)[1];
1005 curls[q_point][1] -=
1006 value * (*shape_gradient_ptr++)[0];
1022 if (shape_function_data[shape_function]
1023 .is_nonzero_shape_function_component[0])
1026 spacedim> *shape_gradient_ptr =
1027 &shape_gradients[shape_function_data[shape_function]
1030 for (
unsigned int q_point = 0;
1031 q_point < n_quadrature_points;
1034 curls[q_point][1] +=
1035 value * (*shape_gradient_ptr)[2];
1036 curls[q_point][2] -=
1037 value * (*shape_gradient_ptr++)[1];
1041 if (shape_function_data[shape_function]
1042 .is_nonzero_shape_function_component[1])
1045 spacedim> *shape_gradient_ptr =
1046 &shape_gradients[shape_function_data[shape_function]
1049 for (
unsigned int q_point = 0;
1050 q_point < n_quadrature_points;
1053 curls[q_point][0] -=
1054 value * (*shape_gradient_ptr)[2];
1055 curls[q_point][2] +=
1056 value * (*shape_gradient_ptr++)[0];
1060 if (shape_function_data[shape_function]
1061 .is_nonzero_shape_function_component[2])
1064 spacedim> *shape_gradient_ptr =
1065 &shape_gradients[shape_function_data[shape_function]
1068 for (
unsigned int q_point = 0;
1069 q_point < n_quadrature_points;
1072 curls[q_point][0] +=
1073 value * (*shape_gradient_ptr)[1];
1074 curls[q_point][1] -=
1075 value * (*shape_gradient_ptr++)[0];
1086 template <
int dim,
int spacedim,
typename Number>
1092 &shape_function_data,
1094 template solution_laplacian_type<Number>> &laplacians)
1096 const unsigned int dofs_per_cell = dof_values.
size();
1097 const unsigned int n_quadrature_points = laplacians.size();
1103 spacedim>::template solution_laplacian_type<Number>());
1105 for (
unsigned int shape_function = 0; shape_function < dofs_per_cell;
1109 shape_function_data[shape_function].single_nonzero_component;
1115 const Number &value = dof_values[shape_function];
1119 if (::internal::CheckForZero<Number>::value(value) ==
true)
1124 const unsigned int comp = shape_function_data[shape_function]
1125 .single_nonzero_component_index;
1126 const ::Tensor<2, spacedim> *shape_hessian_ptr =
1127 &shape_hessians[snc][0];
1128 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
1130 laplacians[q_point][comp] +=
1131 value *
trace(*shape_hessian_ptr++);
1134 for (
unsigned int d = 0; d < spacedim; ++d)
1135 if (shape_function_data[shape_function]
1136 .is_nonzero_shape_function_component[d])
1138 const ::Tensor<2, spacedim> *shape_hessian_ptr =
1139 &shape_hessians[shape_function_data[shape_function]
1141 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
1143 laplacians[q_point][d] +=
1144 value *
trace(*shape_hessian_ptr++);
1153 template <
int dim,
int spacedim,
typename Number>
1157 const ::Table<2, double> &shape_values,
1160 &shape_function_data,
1166 const unsigned int dofs_per_cell = dof_values.
size();
1167 const unsigned int n_quadrature_points = values.size();
1175 for (
unsigned int shape_function = 0; shape_function < dofs_per_cell;
1179 shape_function_data[shape_function].single_nonzero_component;
1185 const Number &value = dof_values[shape_function];
1189 if (::internal::CheckForZero<Number>::value(value) ==
true)
1196 shape_function_data[shape_function]
1197 .single_nonzero_component_index);
1198 const double *shape_value_ptr = &shape_values(snc, 0);
1199 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
1201 values[q_point][comp] += value * (*shape_value_ptr++);
1204 for (
unsigned int d = 0;
1208 if (shape_function_data[shape_function]
1209 .is_nonzero_shape_function_component[d])
1214 const double *shape_value_ptr = &shape_values(
1215 shape_function_data[shape_function].row_index[d], 0);
1216 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
1218 values[q_point][comp] += value * (*shape_value_ptr++);
1225 template <
int dim,
int spacedim,
typename Number>
1232 &shape_function_data,
1236 const unsigned int dofs_per_cell = dof_values.
size();
1237 const unsigned int n_quadrature_points = divergences.size();
1239 std::fill(divergences.begin(),
1244 for (
unsigned int shape_function = 0; shape_function < dofs_per_cell;
1248 shape_function_data[shape_function].single_nonzero_component;
1254 const Number &value = dof_values[shape_function];
1258 if (::internal::CheckForZero<Number>::value(value) ==
true)
1263 const unsigned int comp = shape_function_data[shape_function]
1264 .single_nonzero_component_index;
1266 const ::Tensor<1, spacedim> *shape_gradient_ptr =
1267 &shape_gradients[snc][0];
1274 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
1275 ++q_point, ++shape_gradient_ptr)
1277 divergences[q_point][ii] += value * (*shape_gradient_ptr)[jj];
1280 divergences[q_point][jj] +=
1281 value * (*shape_gradient_ptr)[ii];
1286 for (
unsigned int d = 0;
1289 spacedim>::n_independent_components;
1291 if (shape_function_data[shape_function]
1292 .is_nonzero_shape_function_component[d])
1304 const unsigned int comp =
1305 shape_function_data[shape_function]
1306 .single_nonzero_component_index;
1308 const ::Tensor<1, spacedim> *shape_gradient_ptr =
1309 &shape_gradients[shape_function_data[shape_function]
1311 for (
unsigned int q_point = 0;
1312 q_point < n_quadrature_points;
1313 ++q_point, ++shape_gradient_ptr)
1315 for (
unsigned int j = 0; j < spacedim; ++j)
1317 const unsigned int vector_component =
1321 divergences[q_point][vector_component] +=
1322 value * (*shape_gradient_ptr++)[j];
1332 template <
int dim,
int spacedim,
typename Number>
1336 const ::Table<2, double> &shape_values,
1338 &shape_function_data,
1343 const unsigned int dofs_per_cell = dof_values.
size();
1344 const unsigned int n_quadrature_points = values.size();
1351 for (
unsigned int shape_function = 0; shape_function < dofs_per_cell;
1355 shape_function_data[shape_function].single_nonzero_component;
1361 const Number &value = dof_values[shape_function];
1365 if (::internal::CheckForZero<Number>::value(value) ==
true)
1370 const unsigned int comp = shape_function_data[shape_function]
1371 .single_nonzero_component_index;
1377 const double *shape_value_ptr = &shape_values(snc, 0);
1378 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
1380 values[q_point][indices] += value * (*shape_value_ptr++);
1383 for (
unsigned int d = 0; d < dim * dim; ++d)
1384 if (shape_function_data[shape_function]
1385 .is_nonzero_shape_function_component[d])
1391 const double *shape_value_ptr = &shape_values(
1392 shape_function_data[shape_function].row_index[d], 0);
1393 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
1395 values[q_point][indices] += value * (*shape_value_ptr++);
1402 template <
int dim,
int spacedim,
typename Number>
1408 &shape_function_data,
1410 template solution_divergence_type<Number>> &divergences)
1412 const unsigned int dofs_per_cell = dof_values.
size();
1413 const unsigned int n_quadrature_points = divergences.size();
1416 divergences.begin(),
1421 for (
unsigned int shape_function = 0; shape_function < dofs_per_cell;
1425 shape_function_data[shape_function].single_nonzero_component;
1431 const Number &value = dof_values[shape_function];
1435 if (::internal::CheckForZero<Number>::value(value) ==
true)
1440 const unsigned int comp = shape_function_data[shape_function]
1441 .single_nonzero_component_index;
1443 const ::Tensor<1, spacedim> *shape_gradient_ptr =
1444 &shape_gradients[snc][0];
1449 const unsigned int ii = indices[0];
1450 const unsigned int jj = indices[1];
1452 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
1453 ++q_point, ++shape_gradient_ptr)
1455 divergences[q_point][ii] += value * (*shape_gradient_ptr)[jj];
1460 for (
unsigned int d = 0; d < dim * dim; ++d)
1461 if (shape_function_data[shape_function]
1462 .is_nonzero_shape_function_component[d])
1472 template <
int dim,
int spacedim,
typename Number>
1478 &shape_function_data,
1480 template solution_gradient_type<Number>> &gradients)
1482 const unsigned int dofs_per_cell = dof_values.
size();
1483 const unsigned int n_quadrature_points = gradients.size();
1491 for (
unsigned int shape_function = 0; shape_function < dofs_per_cell;
1495 shape_function_data[shape_function].single_nonzero_component;
1501 const Number &value = dof_values[shape_function];
1505 if (::internal::CheckForZero<Number>::value(value) ==
true)
1510 const unsigned int comp = shape_function_data[shape_function]
1511 .single_nonzero_component_index;
1513 const ::Tensor<1, spacedim> *shape_gradient_ptr =
1514 &shape_gradients[snc][0];
1519 const unsigned int ii = indices[0];
1520 const unsigned int jj = indices[1];
1522 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
1523 ++q_point, ++shape_gradient_ptr)
1525 gradients[q_point][ii][jj] += value * (*shape_gradient_ptr);
1530 for (
unsigned int d = 0; d < dim * dim; ++d)
1531 if (shape_function_data[shape_function]
1532 .is_nonzero_shape_function_component[d])
1544 template <
int dim,
int spacedim>
1545 template <
class InputVector>
1548 const InputVector &fe_function,
1566 internal::do_function_values<dim, spacedim>(
1569 shape_function_data,
1575 template <
int dim,
int spacedim>
1576 template <
class InputVector>
1579 const InputVector &dof_values,
1590 internal::do_function_values<dim, spacedim>(
1593 shape_function_data,
1599 template <
int dim,
int spacedim>
1600 template <
class InputVector>
1603 const InputVector &fe_function,
1609 "update_gradients")));
1620 internal::do_function_derivatives<1, dim, spacedim>(
1623 shape_function_data,
1629 template <
int dim,
int spacedim>
1630 template <
class InputVector>
1633 const InputVector &dof_values,
1639 "update_gradients")));
1644 internal::do_function_derivatives<1, dim, spacedim>(
1647 shape_function_data,
1653 template <
int dim,
int spacedim>
1654 template <
class InputVector>
1657 const InputVector &fe_function,
1663 "update_hessians")));
1674 internal::do_function_derivatives<2, dim, spacedim>(
1677 shape_function_data,
1683 template <
int dim,
int spacedim>
1684 template <
class InputVector>
1687 const InputVector &dof_values,
1693 "update_hessians")));
1698 internal::do_function_derivatives<2, dim, spacedim>(
1701 shape_function_data,
1707 template <
int dim,
int spacedim>
1708 template <
class InputVector>
1711 const InputVector &fe_function,
1717 "update_hessians")));
1728 internal::do_function_laplacians<dim, spacedim>(
1731 shape_function_data,
1737 template <
int dim,
int spacedim>
1738 template <
class InputVector>
1741 const InputVector &dof_values,
1747 "update_hessians")));
1752 internal::do_function_laplacians<dim, spacedim>(
1755 shape_function_data,
1761 template <
int dim,
int spacedim>
1762 template <
class InputVector>
1765 const InputVector &fe_function,
1768 &third_derivatives)
const
1772 "update_3rd_derivatives")));
1783 internal::do_function_derivatives<3, dim, spacedim>(
1786 shape_function_data,
1792 template <
int dim,
int spacedim>
1793 template <
class InputVector>
1796 const InputVector &dof_values,
1799 &third_derivatives)
const
1803 "update_3rd_derivatives")));
1808 internal::do_function_derivatives<3, dim, spacedim>(
1811 shape_function_data,
1817 template <
int dim,
int spacedim>
1818 template <
class InputVector>
1821 const InputVector &fe_function,
1838 internal::do_function_values<dim, spacedim>(
1841 shape_function_data,
1847 template <
int dim,
int spacedim>
1848 template <
class InputVector>
1851 const InputVector &dof_values,
1862 internal::do_function_values<dim, spacedim>(
1865 shape_function_data,
1871 template <
int dim,
int spacedim>
1872 template <
class InputVector>
1875 const InputVector &fe_function,
1881 "update_gradients")));
1892 internal::do_function_derivatives<1, dim, spacedim>(
1895 shape_function_data,
1901 template <
int dim,
int spacedim>
1902 template <
class InputVector>
1905 const InputVector &dof_values,
1911 "update_gradients")));
1916 internal::do_function_derivatives<1, dim, spacedim>(
1919 shape_function_data,
1925 template <
int dim,
int spacedim>
1926 template <
class InputVector>
1929 const InputVector &fe_function,
1932 &symmetric_gradients)
const
1936 "update_gradients")));
1947 internal::do_function_symmetric_gradients<dim, spacedim>(
1950 shape_function_data,
1951 symmetric_gradients);
1956 template <
int dim,
int spacedim>
1957 template <
class InputVector>
1960 const InputVector &dof_values,
1963 &symmetric_gradients)
const
1967 "update_gradients")));
1972 internal::do_function_symmetric_gradients<dim, spacedim>(
1975 shape_function_data,
1976 symmetric_gradients);
1981 template <
int dim,
int spacedim>
1982 template <
class InputVector>
1985 const InputVector &fe_function,
1991 "update_gradients")));
2003 internal::do_function_divergences<dim, spacedim>(
2006 shape_function_data,
2012 template <
int dim,
int spacedim>
2013 template <
class InputVector>
2016 const InputVector &dof_values,
2022 "update_gradients")));
2027 internal::do_function_divergences<dim, spacedim>(
2030 shape_function_data,
2036 template <
int dim,
int spacedim>
2037 template <
class InputVector>
2040 const InputVector &fe_function,
2046 "update_gradients")));
2048 ExcMessage(
"FEValues object is not reinited to any cell"));
2057 internal::do_function_curls<dim, spacedim>(
2060 shape_function_data,
2066 template <
int dim,
int spacedim>
2067 template <
class InputVector>
2070 const InputVector &dof_values,
2076 "update_gradients")));
2078 ExcMessage(
"FEValues object is not reinited to any cell"));
2081 internal::do_function_curls<dim, spacedim>(
2084 shape_function_data,
2090 template <
int dim,
int spacedim>
2091 template <
class InputVector>
2094 const InputVector &fe_function,
2100 "update_hessians")));
2111 internal::do_function_derivatives<2, dim, spacedim>(
2114 shape_function_data,
2120 template <
int dim,
int spacedim>
2121 template <
class InputVector>
2124 const InputVector &dof_values,
2130 "update_hessians")));
2135 internal::do_function_derivatives<2, dim, spacedim>(
2138 shape_function_data,
2144 template <
int dim,
int spacedim>
2145 template <
class InputVector>
2148 const InputVector &fe_function,
2154 "update_hessians")));
2170 internal::do_function_laplacians<dim, spacedim>(
2173 shape_function_data,
2179 template <
int dim,
int spacedim>
2180 template <
class InputVector>
2183 const InputVector &dof_values,
2189 "update_hessians")));
2197 internal::do_function_laplacians<dim, spacedim>(
2200 shape_function_data,
2206 template <
int dim,
int spacedim>
2207 template <
class InputVector>
2210 const InputVector &fe_function,
2213 &third_derivatives)
const
2217 "update_3rd_derivatives")));
2228 internal::do_function_derivatives<3, dim, spacedim>(
2231 shape_function_data,
2237 template <
int dim,
int spacedim>
2238 template <
class InputVector>
2241 const InputVector &dof_values,
2244 &third_derivatives)
const
2248 "update_3rd_derivatives")));
2253 internal::do_function_derivatives<3, dim, spacedim>(
2256 shape_function_data,
2262 template <
int dim,
int spacedim>
2263 template <
class InputVector>
2266 const InputVector &fe_function,
2283 internal::do_function_values<dim, spacedim>(
2286 shape_function_data,
2292 template <
int dim,
int spacedim>
2293 template <
class InputVector>
2296 const InputVector &dof_values,
2307 internal::do_function_values<dim, spacedim>(
2310 shape_function_data,
2316 template <
int dim,
int spacedim>
2317 template <
class InputVector>
2320 const InputVector &fe_function,
2326 "update_gradients")));
2338 internal::do_function_divergences<dim, spacedim>(
2341 shape_function_data,
2347 template <
int dim,
int spacedim>
2348 template <
class InputVector>
2352 const InputVector &dof_values,
2358 "update_gradients")));
2363 internal::do_function_divergences<dim, spacedim>(
2366 shape_function_data,
2372 template <
int dim,
int spacedim>
2373 template <
class InputVector>
2376 const InputVector &fe_function,
2393 internal::do_function_values<dim, spacedim>(
2396 shape_function_data,
2402 template <
int dim,
int spacedim>
2403 template <
class InputVector>
2406 const InputVector &dof_values,
2417 internal::do_function_values<dim, spacedim>(
2420 shape_function_data,
2426 template <
int dim,
int spacedim>
2427 template <
class InputVector>
2430 const InputVector &fe_function,
2436 "update_gradients")));
2448 internal::do_function_divergences<dim, spacedim>(
2451 shape_function_data,
2457 template <
int dim,
int spacedim>
2458 template <
class InputVector>
2461 const InputVector &dof_values,
2462 std::vector<solution_divergence_type<typename InputVector::value_type>>
2467 "update_gradients")));
2472 internal::do_function_divergences<dim, spacedim>(
2475 shape_function_data,
2481 template <
int dim,
int spacedim>
2482 template <
class InputVector>
2485 const InputVector &fe_function,
2491 "update_gradients")));
2503 internal::do_function_gradients<dim, spacedim>(
2506 shape_function_data,
2512 template <
int dim,
int spacedim>
2513 template <
class InputVector>
2516 const InputVector &dof_values,
2522 "update_gradients")));
2527 internal::do_function_gradients<dim, spacedim>(
2530 shape_function_data,
2541 template <
int dim,
int spacedim>
2547 scalars.reserve(n_scalars);
2548 for (
unsigned int component = 0; component < n_scalars; ++component)
2549 scalars.emplace_back(fe_values, component);
2554 const unsigned int n_vectors =
2559 vectors.reserve(n_vectors);
2560 for (
unsigned int component = 0; component < n_vectors; ++component)
2561 vectors.emplace_back(fe_values, component);
2564 const unsigned int n_symmetric_second_order_tensors =
2570 symmetric_second_order_tensors.reserve(n_symmetric_second_order_tensors);
2571 for (
unsigned int component = 0;
2572 component < n_symmetric_second_order_tensors;
2574 symmetric_second_order_tensors.emplace_back(fe_values, component);
2578 const unsigned int n_second_order_tensors =
2583 second_order_tensors.reserve(n_second_order_tensors);
2584 for (
unsigned int component = 0; component < n_second_order_tensors;
2586 second_order_tensors.emplace_back(fe_values, component);
2595template <
int dim,
int spacedim>
2597 : initialized(false)
2598 , cell(typename
Triangulation<dim, spacedim>::cell_iterator(nullptr, -1, -1))
2599 , dof_handler(nullptr)
2600 , level_dof_access(false)
2605template <
int dim,
int spacedim>
2610 , dof_handler(nullptr)
2611 , level_dof_access(false)
2616template <
int dim,
int spacedim>
2625template <
int dim,
int spacedim>
2636template <
int dim,
int spacedim>
2642 Assert(dof_handler !=
nullptr, ExcNeedsDoFHandler());
2644 return dof_handler->n_dofs();
2649template <
int dim,
int spacedim>
2650template <
typename VectorType>
2653 const VectorType & in,
2657 Assert(dof_handler !=
nullptr, ExcNeedsDoFHandler());
2659 if (level_dof_access)
2675template <
int dim,
int spacedim>
2682 Assert(dof_handler !=
nullptr, ExcNeedsDoFHandler());
2686 &cell->get_triangulation(), cell->level(), cell->index(), dof_handler);
2689 cell_dofs.
get_fe().n_dofs_per_cell());
2692 for (
unsigned int i = 0; i < cell_dofs.
get_fe().n_dofs_per_cell(); ++i)
2700 namespace FEValuesImplementation
2702 template <
int dim,
int spacedim>
2705 const unsigned int n_quadrature_points,
2709 this->quadrature_points.resize(
2710 n_quadrature_points,
2714 this->JxW_values.resize(n_quadrature_points,
2715 numbers::signaling_nan<double>());
2718 this->jacobians.resize(
2719 n_quadrature_points,
2723 this->jacobian_grads.resize(
2724 n_quadrature_points,
2728 this->jacobian_pushed_forward_grads.resize(
2732 this->jacobian_2nd_derivatives.resize(
2733 n_quadrature_points,
2737 this->jacobian_pushed_forward_2nd_derivatives.resize(
2741 this->jacobian_3rd_derivatives.resize(n_quadrature_points);
2744 this->jacobian_pushed_forward_3rd_derivatives.resize(
2748 this->inverse_jacobians.resize(
2749 n_quadrature_points,
2753 this->boundary_forms.resize(
2757 this->normal_vectors.resize(
2763 template <
int dim,
int spacedim>
2774 jacobian_pushed_forward_2nd_derivatives) +
2777 jacobian_pushed_forward_3rd_derivatives) +
2786 template <
int dim,
int spacedim>
2789 const unsigned int n_quadrature_points,
2796 this->shape_function_to_row_table =
2801 unsigned int n_nonzero_shape_components = 0;
2811 this->shape_values.reinit(n_nonzero_shape_components,
2812 n_quadrature_points);
2813 this->shape_values.fill(numbers::signaling_nan<double>());
2818 this->shape_gradients.reinit(n_nonzero_shape_components,
2819 n_quadrature_points);
2820 this->shape_gradients.fill(
2826 this->shape_hessians.reinit(n_nonzero_shape_components,
2827 n_quadrature_points);
2828 this->shape_hessians.fill(
2834 this->shape_3rd_derivatives.reinit(n_nonzero_shape_components,
2835 n_quadrature_points);
2836 this->shape_3rd_derivatives.fill(
2843 template <
int dim,
int spacedim>
2862template <
int dim,
int spacedim>
2864 const unsigned int n_q_points,
2873 ,
fe(&
fe, typeid(*this).name())
2878 ExcMessage(
"There is nothing useful you can do with an FEValues "
2879 "object when using a quadrature formula with zero "
2880 "quadrature points!"));
2886template <
int dim,
int spacedim>
2889 tria_listener_refinement.disconnect();
2890 tria_listener_mesh_transform.disconnect();
2904 template <
typename Number,
typename Number2>
2907 const ::Table<2, double> &shape_values,
2908 std::vector<Number> & values)
2911 const unsigned int dofs_per_cell = shape_values.n_rows();
2912 const unsigned int n_quadrature_points = values.size();
2915 std::fill_n(values.begin(),
2916 n_quadrature_points,
2926 for (
unsigned int shape_func = 0; shape_func < dofs_per_cell; ++shape_func)
2928 const Number2
value = dof_values_ptr[shape_func];
2936 const double *shape_value_ptr = &shape_values(shape_func, 0);
2937 for (
unsigned int point = 0; point < n_quadrature_points; ++point)
2938 values[point] +=
value * (*shape_value_ptr++);
2944 template <
int dim,
int spacedim,
typename VectorType>
2947 const typename VectorType::value_type *dof_values_ptr,
2948 const ::Table<2, double> & shape_values,
2950 const std::vector<unsigned int> & shape_function_to_row_table,
2952 const bool quadrature_points_fastest =
false,
2953 const unsigned int component_multiple = 1)
2955 using Number =
typename VectorType::value_type;
2957 for (
unsigned int i = 0; i < values.size(); ++i)
2958 std::fill_n(values[i].begin(),
2960 typename VectorType::value_type());
2965 if (dofs_per_cell == 0)
2968 const unsigned int n_quadrature_points =
2969 quadrature_points_fastest ? values[0].size() : values.size();
2973 const unsigned result_components = n_components * component_multiple;
2974 (void)result_components;
2975 if (quadrature_points_fastest)
2978 for (
unsigned int i = 0; i < values.size(); ++i)
2984 for (
unsigned int i = 0; i < values.size(); ++i)
2991 for (
unsigned int mc = 0; mc < component_multiple; ++mc)
2992 for (
unsigned int shape_func = 0; shape_func < dofs_per_cell;
2995 const Number &
value = dof_values_ptr[shape_func + mc * dofs_per_cell];
2999 if (::internal::CheckForZero<Number>::value(
value) ==
true)
3004 const unsigned int comp =
3007 const unsigned int row =
3008 shape_function_to_row_table[shape_func * n_components + comp];
3010 const double *shape_value_ptr = &shape_values(row, 0);
3012 if (quadrature_points_fastest)
3014 VectorType &values_comp = values[comp];
3015 for (
unsigned int point = 0; point < n_quadrature_points;
3017 values_comp[point] +=
value * (*shape_value_ptr++);
3020 for (
unsigned int point = 0; point < n_quadrature_points;
3022 values[point][comp] +=
value * (*shape_value_ptr++);
3025 for (
unsigned int c = 0; c < n_components; ++c)
3030 const unsigned int row =
3031 shape_function_to_row_table[shape_func * n_components + c];
3033 const double * shape_value_ptr = &shape_values(row, 0);
3034 const unsigned int comp = c + mc * n_components;
3036 if (quadrature_points_fastest)
3038 VectorType &values_comp = values[comp];
3039 for (
unsigned int point = 0; point < n_quadrature_points;
3041 values_comp[point] +=
value * (*shape_value_ptr++);
3044 for (
unsigned int point = 0; point < n_quadrature_points;
3046 values[point][comp] +=
value * (*shape_value_ptr++);
3055 template <
int order,
int spacedim,
typename Number>
3058 const Number * dof_values_ptr,
3062 const unsigned int dofs_per_cell = shape_derivatives.size()[0];
3063 const unsigned int n_quadrature_points = derivatives.size();
3066 std::fill_n(derivatives.begin(),
3067 n_quadrature_points,
3077 for (
unsigned int shape_func = 0; shape_func < dofs_per_cell; ++shape_func)
3079 const Number &
value = dof_values_ptr[shape_func];
3083 if (::internal::CheckForZero<Number>::value(
value) ==
true)
3087 &shape_derivatives[shape_func][0];
3088 for (
unsigned int point = 0; point < n_quadrature_points; ++point)
3089 derivatives[point] +=
value * (*shape_derivative_ptr++);
3095 template <
int order,
int dim,
int spacedim,
typename Number>
3098 const Number * dof_values_ptr,
3101 const std::vector<unsigned int> &shape_function_to_row_table,
3103 const bool quadrature_points_fastest =
false,
3104 const unsigned int component_multiple = 1)
3107 for (
unsigned int i = 0; i < derivatives.size(); ++i)
3108 std::fill_n(derivatives[i].begin(),
3109 derivatives[i].size(),
3115 if (dofs_per_cell == 0)
3119 const unsigned int n_quadrature_points =
3120 quadrature_points_fastest ? derivatives[0].size() : derivatives.size();
3124 const unsigned result_components = n_components * component_multiple;
3125 (void)result_components;
3126 if (quadrature_points_fastest)
3129 for (
unsigned int i = 0; i < derivatives.size(); ++i)
3135 for (
unsigned int i = 0; i < derivatives.size(); ++i)
3142 for (
unsigned int mc = 0; mc < component_multiple; ++mc)
3143 for (
unsigned int shape_func = 0; shape_func < dofs_per_cell;
3146 const Number &
value = dof_values_ptr[shape_func + mc * dofs_per_cell];
3150 if (::internal::CheckForZero<Number>::value(
value) ==
true)
3155 const unsigned int comp =
3158 const unsigned int row =
3159 shape_function_to_row_table[shape_func * n_components + comp];
3162 &shape_derivatives[row][0];
3164 if (quadrature_points_fastest)
3165 for (
unsigned int point = 0; point < n_quadrature_points;
3167 derivatives[comp][point] +=
value * (*shape_derivative_ptr++);
3169 for (
unsigned int point = 0; point < n_quadrature_points;
3171 derivatives[point][comp] +=
value * (*shape_derivative_ptr++);
3174 for (
unsigned int c = 0; c < n_components; ++c)
3179 const unsigned int row =
3180 shape_function_to_row_table[shape_func * n_components + c];
3183 &shape_derivatives[row][0];
3184 const unsigned int comp = c + mc * n_components;
3186 if (quadrature_points_fastest)
3187 for (
unsigned int point = 0; point < n_quadrature_points;
3189 derivatives[comp][point] +=
3190 value * (*shape_derivative_ptr++);
3192 for (
unsigned int point = 0; point < n_quadrature_points;
3194 derivatives[point][comp] +=
3195 value * (*shape_derivative_ptr++);
3202 template <
int spacedim,
typename Number,
typename Number2>
3205 const Number2 * dof_values_ptr,
3207 std::vector<Number> & laplacians)
3209 const unsigned int dofs_per_cell = shape_hessians.size()[0];
3210 const unsigned int n_quadrature_points = laplacians.size();
3213 std::fill_n(laplacians.begin(),
3214 n_quadrature_points,
3220 for (
unsigned int shape_func = 0; shape_func < dofs_per_cell; ++shape_func)
3222 const Number2
value = dof_values_ptr[shape_func];
3231 &shape_hessians[shape_func][0];
3232 for (
unsigned int point = 0; point < n_quadrature_points; ++point)
3233 laplacians[point] +=
value *
trace(*shape_hessian_ptr++);
3239 template <
int dim,
int spacedim,
typename VectorType,
typename Number>
3242 const Number * dof_values_ptr,
3245 const std::vector<unsigned int> & shape_function_to_row_table,
3246 std::vector<VectorType> & laplacians,
3247 const bool quadrature_points_fastest =
false,
3248 const unsigned int component_multiple = 1)
3251 for (
unsigned int i = 0; i < laplacians.size(); ++i)
3252 std::fill_n(laplacians[i].begin(),
3253 laplacians[i].size(),
3254 typename VectorType::value_type());
3259 if (dofs_per_cell == 0)
3263 const unsigned int n_quadrature_points = laplacians.size();
3267 const unsigned result_components = n_components * component_multiple;
3268 (void)result_components;
3269 if (quadrature_points_fastest)
3272 for (
unsigned int i = 0; i < laplacians.size(); ++i)
3278 for (
unsigned int i = 0; i < laplacians.size(); ++i)
3285 for (
unsigned int mc = 0; mc < component_multiple; ++mc)
3286 for (
unsigned int shape_func = 0; shape_func < dofs_per_cell;
3289 const Number &
value = dof_values_ptr[shape_func + mc * dofs_per_cell];
3293 if (::internal::CheckForZero<Number>::value(
value) ==
true)
3298 const unsigned int comp =
3301 const unsigned int row =
3302 shape_function_to_row_table[shape_func * n_components + comp];
3305 &shape_hessians[row][0];
3306 if (quadrature_points_fastest)
3308 VectorType &laplacians_comp = laplacians[comp];
3309 for (
unsigned int point = 0; point < n_quadrature_points;
3311 laplacians_comp[point] +=
3315 for (
unsigned int point = 0; point < n_quadrature_points;
3317 laplacians[point][comp] +=
3321 for (
unsigned int c = 0; c < n_components; ++c)
3326 const unsigned int row =
3327 shape_function_to_row_table[shape_func * n_components + c];
3330 &shape_hessians[row][0];
3331 const unsigned int comp = c + mc * n_components;
3333 if (quadrature_points_fastest)
3335 VectorType &laplacians_comp = laplacians[comp];
3336 for (
unsigned int point = 0; point < n_quadrature_points;
3338 laplacians_comp[point] +=
3342 for (
unsigned int point = 0; point < n_quadrature_points;
3344 laplacians[point][comp] +=
3353template <
int dim,
int spacedim>
3354template <
class InputVector>
3357 const InputVector & fe_function,
3358 std::vector<typename InputVector::value_type> &values)
const
3360 using Number =
typename InputVector::value_type;
3362 ExcAccessToUninitializedField(
"update_values"));
3364 Assert(present_cell.is_initialized(), ExcNotReinited());
3365 AssertDimension(fe_function.size(), present_cell.n_dofs_for_dof_handler());
3369 present_cell.get_interpolated_dof_values(fe_function, dof_values);
3371 this->finite_element_output.shape_values,
3377template <
int dim,
int spacedim>
3378template <
class InputVector>
3381 const InputVector & fe_function,
3383 std::vector<typename InputVector::value_type> & values)
const
3385 using Number =
typename InputVector::value_type;
3387 ExcAccessToUninitializedField(
"update_values"));
3391 boost::container::small_vector<Number, 200> dof_values(dofs_per_cell);
3392 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
3395 this->finite_element_output.shape_values,
3401template <
int dim,
int spacedim>
3402template <
class InputVector>
3405 const InputVector & fe_function,
3408 using Number =
typename InputVector::value_type;
3409 Assert(present_cell.is_initialized(), ExcNotReinited());
3412 ExcAccessToUninitializedField(
"update_values"));
3413 AssertDimension(fe_function.size(), present_cell.n_dofs_for_dof_handler());
3417 present_cell.get_interpolated_dof_values(fe_function, dof_values);
3420 this->finite_element_output.shape_values,
3422 this->finite_element_output.shape_function_to_row_table,
3428template <
int dim,
int spacedim>
3429template <
class InputVector>
3432 const InputVector & fe_function,
3436 using Number =
typename InputVector::value_type;
3442 ExcAccessToUninitializedField(
"update_values"));
3444 boost::container::small_vector<Number, 200> dof_values(dofs_per_cell);
3445 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
3449 this->finite_element_output.shape_values,
3451 this->finite_element_output.shape_function_to_row_table,
3454 indices.
size() / dofs_per_cell);
3459template <
int dim,
int spacedim>
3460template <
class InputVector>
3463 const InputVector & fe_function,
3465 ArrayView<std::vector<typename InputVector::value_type>> values,
3466 const bool quadrature_points_fastest)
const
3468 using Number =
typename InputVector::value_type;
3470 ExcAccessToUninitializedField(
"update_values"));
3477 boost::container::small_vector<Number, 200> dof_values(indices.
size());
3478 for (
unsigned int i = 0; i < indices.
size(); ++i)
3482 this->finite_element_output.shape_values,
3484 this->finite_element_output.shape_function_to_row_table,
3486 quadrature_points_fastest,
3487 indices.
size() / dofs_per_cell);
3492template <
int dim,
int spacedim>
3493template <
class InputVector>
3496 const InputVector &fe_function,
3500 using Number =
typename InputVector::value_type;
3502 ExcAccessToUninitializedField(
"update_gradients"));
3504 Assert(present_cell.is_initialized(), ExcNotReinited());
3505 AssertDimension(fe_function.size(), present_cell.n_dofs_for_dof_handler());
3509 present_cell.get_interpolated_dof_values(fe_function, dof_values);
3511 this->finite_element_output.shape_gradients,
3517template <
int dim,
int spacedim>
3518template <
class InputVector>
3521 const InputVector & fe_function,
3526 using Number =
typename InputVector::value_type;
3528 ExcAccessToUninitializedField(
"update_gradients"));
3532 boost::container::small_vector<Number, 200> dof_values(dofs_per_cell);
3533 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
3536 this->finite_element_output.shape_gradients,
3542template <
int dim,
int spacedim>
3543template <
class InputVector>
3546 const InputVector &fe_function,
3551 using Number =
typename InputVector::value_type;
3553 ExcAccessToUninitializedField(
"update_gradients"));
3554 Assert(present_cell.is_initialized(), ExcNotReinited());
3555 AssertDimension(fe_function.size(), present_cell.n_dofs_for_dof_handler());
3559 present_cell.get_interpolated_dof_values(fe_function, dof_values);
3562 this->finite_element_output.shape_gradients,
3564 this->finite_element_output.shape_function_to_row_table,
3570template <
int dim,
int spacedim>
3571template <
class InputVector>
3574 const InputVector & fe_function,
3578 const bool quadrature_points_fastest)
const
3580 using Number =
typename InputVector::value_type;
3586 ExcAccessToUninitializedField(
"update_gradients"));
3588 boost::container::small_vector<Number, 200> dof_values(indices.
size());
3589 for (
unsigned int i = 0; i < indices.
size(); ++i)
3593 this->finite_element_output.shape_gradients,
3595 this->finite_element_output.shape_function_to_row_table,
3597 quadrature_points_fastest,
3598 indices.
size() / dofs_per_cell);
3603template <
int dim,
int spacedim>
3604template <
class InputVector>
3607 const InputVector &fe_function,
3611 using Number =
typename InputVector::value_type;
3614 ExcAccessToUninitializedField(
"update_hessians"));
3615 Assert(present_cell.is_initialized(), ExcNotReinited());
3616 AssertDimension(fe_function.size(), present_cell.n_dofs_for_dof_handler());
3620 present_cell.get_interpolated_dof_values(fe_function, dof_values);
3622 this->finite_element_output.shape_hessians,
3628template <
int dim,
int spacedim>
3629template <
class InputVector>
3632 const InputVector & fe_function,
3637 using Number =
typename InputVector::value_type;
3639 ExcAccessToUninitializedField(
"update_hessians"));
3640 AssertDimension(fe_function.size(), present_cell.n_dofs_for_dof_handler());
3643 boost::container::small_vector<Number, 200> dof_values(dofs_per_cell);
3644 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
3647 this->finite_element_output.shape_hessians,
3653template <
int dim,
int spacedim>
3654template <
class InputVector>
3657 const InputVector &fe_function,
3661 const bool quadrature_points_fastest)
const
3663 using Number =
typename InputVector::value_type;
3665 ExcAccessToUninitializedField(
"update_hessians"));
3666 Assert(present_cell.is_initialized(), ExcNotReinited());
3667 AssertDimension(fe_function.size(), present_cell.n_dofs_for_dof_handler());
3671 present_cell.get_interpolated_dof_values(fe_function, dof_values);
3674 this->finite_element_output.shape_hessians,
3676 this->finite_element_output.shape_function_to_row_table,
3678 quadrature_points_fastest);
3683template <
int dim,
int spacedim>
3684template <
class InputVector>
3687 const InputVector & fe_function,
3691 const bool quadrature_points_fastest)
const
3693 using Number =
typename InputVector::value_type;
3695 ExcAccessToUninitializedField(
"update_hessians"));
3699 boost::container::small_vector<Number, 200> dof_values(indices.
size());
3700 for (
unsigned int i = 0; i < indices.
size(); ++i)
3704 this->finite_element_output.shape_hessians,
3706 this->finite_element_output.shape_function_to_row_table,
3708 quadrature_points_fastest,
3709 indices.
size() / dofs_per_cell);
3714template <
int dim,
int spacedim>
3715template <
class InputVector>
3718 const InputVector & fe_function,
3719 std::vector<typename InputVector::value_type> &laplacians)
const
3721 using Number =
typename InputVector::value_type;
3723 ExcAccessToUninitializedField(
"update_hessians"));
3725 Assert(present_cell.is_initialized(), ExcNotReinited());
3730 present_cell.get_interpolated_dof_values(fe_function, dof_values);
3732 this->finite_element_output.shape_hessians,
3738template <
int dim,
int spacedim>
3739template <
class InputVector>
3742 const InputVector & fe_function,
3744 std::vector<typename InputVector::value_type> & laplacians)
const
3746 using Number =
typename InputVector::value_type;
3748 ExcAccessToUninitializedField(
"update_hessians"));
3752 boost::container::small_vector<Number, 200> dof_values(dofs_per_cell);
3753 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
3756 this->finite_element_output.shape_hessians,
3762template <
int dim,
int spacedim>
3763template <
class InputVector>
3766 const InputVector & fe_function,
3769 using Number =
typename InputVector::value_type;
3770 Assert(present_cell.is_initialized(), ExcNotReinited());
3772 ExcAccessToUninitializedField(
"update_hessians"));
3773 AssertDimension(fe_function.size(), present_cell.n_dofs_for_dof_handler());
3777 present_cell.get_interpolated_dof_values(fe_function, dof_values);
3780 this->finite_element_output.shape_hessians,
3782 this->finite_element_output.shape_function_to_row_table,
3788template <
int dim,
int spacedim>
3789template <
class InputVector>
3792 const InputVector & fe_function,
3796 using Number =
typename InputVector::value_type;
3802 ExcAccessToUninitializedField(
"update_hessians"));
3804 boost::container::small_vector<Number, 200> dof_values(indices.
size());
3805 for (
unsigned int i = 0; i < indices.
size(); ++i)
3809 this->finite_element_output.shape_hessians,
3811 this->finite_element_output.shape_function_to_row_table,
3814 indices.
size() / dofs_per_cell);
3819template <
int dim,
int spacedim>
3820template <
class InputVector>
3823 const InputVector & fe_function,
3825 std::vector<std::vector<typename InputVector::value_type>> &laplacians,
3826 const bool quadrature_points_fastest)
const
3828 using Number =
typename InputVector::value_type;
3832 ExcAccessToUninitializedField(
"update_hessians"));
3834 boost::container::small_vector<Number, 200> dof_values(indices.
size());
3835 for (
unsigned int i = 0; i < indices.
size(); ++i)
3839 this->finite_element_output.shape_hessians,
3841 this->finite_element_output.shape_function_to_row_table,
3843 quadrature_points_fastest,
3844 indices.
size() / dofs_per_cell);
3849template <
int dim,
int spacedim>
3850template <
class InputVector>
3853 const InputVector &fe_function,
3855 &third_derivatives)
const
3857 using Number =
typename InputVector::value_type;
3860 ExcAccessToUninitializedField(
"update_3rd_derivatives"));
3861 Assert(present_cell.is_initialized(), ExcNotReinited());
3866 present_cell.get_interpolated_dof_values(fe_function, dof_values);
3869 this->finite_element_output.shape_3rd_derivatives,
3875template <
int dim,
int spacedim>
3876template <
class InputVector>
3879 const InputVector & fe_function,
3882 &third_derivatives)
const
3884 using Number =
typename InputVector::value_type;
3886 ExcAccessToUninitializedField(
"update_3rd_derivatives"));
3887 AssertDimension(fe_function.size(), present_cell.n_dofs_for_dof_handler());
3890 boost::container::small_vector<Number, 200> dof_values(dofs_per_cell);
3891 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
3895 this->finite_element_output.shape_3rd_derivatives,
3901template <
int dim,
int spacedim>
3902template <
class InputVector>
3905 const InputVector &fe_function,
3908 & third_derivatives,
3909 const bool quadrature_points_fastest)
const
3911 using Number =
typename InputVector::value_type;
3913 ExcAccessToUninitializedField(
"update_3rd_derivatives"));
3914 Assert(present_cell.is_initialized(), ExcNotReinited());
3915 AssertDimension(fe_function.size(), present_cell.n_dofs_for_dof_handler());
3919 present_cell.get_interpolated_dof_values(fe_function, dof_values);
3922 this->finite_element_output.shape_3rd_derivatives,
3924 this->finite_element_output.shape_function_to_row_table,
3926 quadrature_points_fastest);
3931template <
int dim,
int spacedim>
3932template <
class InputVector>
3935 const InputVector & fe_function,
3939 const bool quadrature_points_fastest)
const
3941 using Number =
typename InputVector::value_type;
3943 ExcAccessToUninitializedField(
"update_3rd_derivatives"));
3947 boost::container::small_vector<Number, 200> dof_values(indices.
size());
3948 for (
unsigned int i = 0; i < indices.
size(); ++i)
3952 this->finite_element_output.shape_3rd_derivatives,
3954 this->finite_element_output.shape_function_to_row_table,
3956 quadrature_points_fastest,
3957 indices.
size() / dofs_per_cell);
3962template <
int dim,
int spacedim>
3966 return present_cell;
3971template <
int dim,
int spacedim>
3972const std::vector<Tensor<1, spacedim>> &
3977 "update_normal_vectors")));
3979 return this->mapping_output.normal_vectors;
3984template <
int dim,
int spacedim>
3988 return (
sizeof(this->update_flags) +
3991 sizeof(cell_similarity) +
4005template <
int dim,
int spacedim>
4017 flags |= mapping->requires_update_flags(flags);
4024template <
int dim,
int spacedim>
4035 tria_listener_refinement.disconnect();
4036 tria_listener_mesh_transform.disconnect();
4042template <
int dim,
int spacedim>
4047 if (present_cell.is_initialized())
4049 if (&cell->get_triangulation() !=
4053 ->get_triangulation())
4060 invalidate_present_cell();
4061 tria_listener_refinement =
4062 cell->get_triangulation().signals.any_change.connect(
4063 [
this]() { this->invalidate_present_cell(); });
4064 tria_listener_mesh_transform =
4065 cell->get_triangulation().signals.mesh_movement.connect(
4066 [
this]() { this->invalidate_present_cell(); });
4074 tria_listener_refinement =
4075 cell->get_triangulation().signals.post_refinement.connect(
4076 [
this]() { this->invalidate_present_cell(); });
4077 tria_listener_mesh_transform =
4078 cell->get_triangulation().signals.mesh_movement.connect(
4079 [
this]() { this->invalidate_present_cell(); });
4085template <
int dim,
int spacedim>
4112 if (this->present_cell.is_initialized() ==
false)
4121 (cell->is_translation_of(
4123 &
>(this->present_cell)) ?
4130 &
>(this->present_cell)
4131 ->direction_flag() != cell->direction_flag())
4140template <
int dim,
int spacedim>
4144 return cell_similarity;
4149template <
int dim,
int spacedim>
4154template <
int dim,
int spacedim>
4159template <
int dim,
int spacedim>
4164template <
int dim,
int spacedim>
4170 fe.n_dofs_per_cell(),
4181template <
int dim,
int spacedim>
4186 :
FEValues(mapping, fe, q[0], update_flags)
4193template <
int dim,
int spacedim>
4199 fe.n_dofs_per_cell(),
4210template <
int dim,
int spacedim>
4221template <
int dim,
int spacedim>
4227 if (dim != spacedim - 1)
4229 ExcMessage(
"You can only pass the 'update_normal_vectors' "
4230 "flag to FEFaceValues or FESubfaceValues objects, "
4231 "but not to an FEValues object unless the "
4232 "triangulation it refers to is embedded in a higher "
4233 "dimensional space."));
4235 const UpdateFlags flags = this->compute_update_flags(update_flags);
4239 this->mapping_output.initialize(this->max_n_quadrature_points, flags);
4240 this->finite_element_output.initialize(this->max_n_quadrature_points,
4247 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>>
4252 this->finite_element_output);
4256 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>>
4260 [&]() {
return this->mapping->get_data(flags, quadrature); });
4262 this->update_flags = flags;
4265 this->fe_data = std::move(fe_get_data.return_value());
4267 this->mapping_data = std::move(mapping_get_data.return_value());
4269 this->mapping_data =
4270 std::make_unique<typename Mapping<dim, spacedim>::InternalDataBase>();
4275template <
int dim,
int spacedim>
4281 Assert(this->get_mapping().is_compatible_with(cell->reference_cell()),
4283 "You are trying to call FEValues::reinit() with a cell of type " +
4284 cell->reference_cell().to_string() +
4285 " with a Mapping that is not compatible with it."));
4289 this->maybe_invalidate_previous_present_cell(cell);
4290 this->check_cell_similarity(cell);
4292 this->present_cell = {cell};
4302template <
int dim,
int spacedim>
4315 Assert(this->get_mapping().is_compatible_with(cell->reference_cell()),
4317 "You are trying to call FEValues::reinit() with a cell of type " +
4318 cell->reference_cell().to_string() +
4319 " with a Mapping that is not compatible with it."));
4321 this->maybe_invalidate_previous_present_cell(cell);
4322 this->check_cell_similarity(cell);
4324 this->present_cell = {cell};
4334template <
int dim,
int spacedim>
4344 this->cell_similarity =
4345 this->get_mapping().fill_fe_values(this->present_cell,
4346 this->cell_similarity,
4348 *this->mapping_data,
4349 this->mapping_output);
4356 this->get_fe().fill_fe_values(this->present_cell,
4357 this->cell_similarity,
4359 this->get_mapping(),
4360 *this->mapping_data,
4361 this->mapping_output,
4363 this->finite_element_output);
4368template <
int dim,
int spacedim>
4380template <
int dim,
int spacedim>
4382 const unsigned int dofs_per_cell,
4391 hp::QCollection<dim - 1>(quadrature))
4396template <
int dim,
int spacedim>
4398 const unsigned int dofs_per_cell,
4403 :
FEValuesBase<dim, spacedim>(quadrature.max_n_quadrature_points(),
4408 , present_face_index(
numbers::invalid_unsigned_int)
4409 , quadrature(quadrature)
4412 quadrature.
size() ==
fe.reference_cell().n_faces(),
4418template <
int dim,
int spacedim>
4419const std::vector<Tensor<1, spacedim>> &
4424 "update_boundary_forms")));
4425 return this->mapping_output.boundary_forms;
4430template <
int dim,
int spacedim>
4441template <
int dim,
int spacedim>
4446template <
int dim,
int spacedim>
4451template <
int dim,
int spacedim>
4459 hp::QCollection<dim - 1>(quadrature),
4465template <
int dim,
int spacedim>
4482template <
int dim,
int spacedim>
4488 hp::QCollection<dim - 1>(quadrature),
4494template <
int dim,
int spacedim>
4500 fe.n_dofs_per_cell(),
4511template <
int dim,
int spacedim>
4515 const UpdateFlags flags = this->compute_update_flags(update_flags);
4519 this->mapping_output.initialize(this->max_n_quadrature_points, flags);
4520 this->finite_element_output.initialize(this->max_n_quadrature_points,
4527 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase> (
4536 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> (
4543 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>>
4549 this->finite_element_output);
4551 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>>
4559 this->update_flags = flags;
4562 this->fe_data = std::move(fe_get_data.return_value());
4564 this->mapping_data = std::move(mapping_get_data.return_value());
4566 this->mapping_data =
4567 std::make_unique<typename Mapping<dim, spacedim>::InternalDataBase>();
4572template <
int dim,
int spacedim>
4577 const unsigned int face_no)
4583 cell->get_dof_handler().get_fe(cell->active_fe_index())),
4588 this->maybe_invalidate_previous_present_cell(cell);
4589 this->present_cell = {cell};
4599template <
int dim,
int spacedim>
4606 const auto face_n = cell->face_iterator_to_index(face);
4612template <
int dim,
int spacedim>
4616 const unsigned int face_no)
4620 this->maybe_invalidate_previous_present_cell(cell);
4621 this->present_cell = {cell};
4631template <
int dim,
int spacedim>
4637 const auto face_n = cell->face_iterator_to_index(face);
4643template <
int dim,
int spacedim>
4647 this->present_face_no = face_no;
4652 this->present_face_index = cell->face_index(face_no);
4656 this->get_mapping().fill_fe_face_values(this->present_cell,
4659 *this->mapping_data,
4660 this->mapping_output);
4663 this->get_fe().fill_fe_face_values(this->present_cell,
4666 this->get_mapping(),
4667 *this->mapping_data,
4668 this->mapping_output,
4670 this->finite_element_output);
4672 const_cast<unsigned int &
>(this->n_quadrature_points) =
4673 this->quadrature[this->quadrature.
size() == 1 ? 0 : face_no].
size();
4680template <
int dim,
int spacedim>
4685template <
int dim,
int spacedim>
4690template <
int dim,
int spacedim>
4707template <
int dim,
int spacedim>
4720template <
int dim,
int spacedim>
4726 fe.n_dofs_per_cell(),
4737template <
int dim,
int spacedim>
4749template <
int dim,
int spacedim>
4753 const UpdateFlags flags = this->compute_update_flags(update_flags);
4757 this->mapping_output.initialize(this->max_n_quadrature_points, flags);
4758 this->finite_element_output.initialize(this->max_n_quadrature_points,
4766 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>>
4772 this->quadrature[0],
4773 this->finite_element_output);
4775 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>>
4782 this->quadrature[0]);
4784 this->update_flags = flags;
4787 this->fe_data = std::move(fe_get_data.return_value());
4789 this->mapping_data = std::move(mapping_get_data.return_value());
4791 this->mapping_data =
4792 std::make_unique<typename Mapping<dim, spacedim>::InternalDataBase>();
4797template <
int dim,
int spacedim>
4802 const unsigned int face_no,
4803 const unsigned int subface_no)
4809 cell->get_dof_handler().get_fe(cell->active_fe_index())),
4816 Assert(cell->face(face_no)->has_children() ||
4821 Assert(!cell->face(face_no)->has_children() ||
4822 subface_no < cell->face(face_no)->n_active_descendants(),
4825 cell->face(face_no)->n_active_descendants()));
4826 Assert(cell->has_children() ==
false,
4827 ExcMessage(
"You can't use subface data for cells that are "
4828 "already refined. Iterate over their children "
4829 "instead in these cases."));
4831 this->maybe_invalidate_previous_present_cell(cell);
4832 this->present_cell = {cell};
4837 do_reinit(face_no, subface_no);
4842template <
int dim,
int spacedim>
4851 cell->face_iterator_to_index(face),
4852 face->child_iterator_to_index(subface));
4857template <
int dim,
int spacedim>
4861 const unsigned int face_no,
4862 const unsigned int subface_no)
4870 (cell->has_periodic_neighbor(face_no) ?
4871 cell->periodic_neighbor(face_no)
4872 ->face(cell->periodic_neighbor_face_no(face_no))
4874 cell->face(face_no)->n_children()));
4876 this->maybe_invalidate_previous_present_cell(cell);
4877 this->present_cell = {cell};
4882 do_reinit(face_no, subface_no);
4887template <
int dim,
int spacedim>
4895 cell->face_iterator_to_index(face),
4896 face->child_iterator_to_index(subface));
4901template <
int dim,
int spacedim>
4904 const unsigned int subface_no)
4906 this->present_face_no = face_no;
4912 if (!cell->face(face_no)->has_children())
4915 this->present_face_index = cell->face_index(face_no);
4917 this->present_face_index = cell->face(face_no)->child_index(subface_no);
4923 switch (cell->subface_case(face_no))
4928 subface_index = cell->face(face_no)->child_index(subface_no);
4932 subface_index = cell->face(face_no)
4933 ->child(subface_no / 2)
4934 ->child_index(subface_no % 2);
4943 cell->face(face_no)->child(0)->child_index(subface_no);
4946 subface_index = cell->face(face_no)->child_index(1);
4957 subface_index = cell->face(face_no)->child_index(0);
4962 cell->face(face_no)->child(1)->child_index(subface_no - 1);
4974 this->present_face_index = subface_index;
4980 this->get_mapping().fill_fe_subface_values(this->present_cell,
4983 this->quadrature[0],
4984 *this->mapping_data,
4985 this->mapping_output);
4988 this->get_fe().fill_fe_subface_values(this->present_cell,
4991 this->quadrature[0],
4992 this->get_mapping(),
4993 *this->mapping_data,
4994 this->mapping_output,
4996 this->finite_element_output);
5001#define SPLIT_INSTANTIATIONS_COUNT 6
5002#ifndef SPLIT_INSTANTIATIONS_INDEX
5003# define SPLIT_INSTANTIATIONS_INDEX 0
5005#include "fe_values.inst"
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
value_type * data() const noexcept
const FiniteElement< dimension_, space_dimension_ > & get_fe() const
void get_interpolated_dof_values(const InputVector &values, Vector< number > &interpolated_values, const unsigned int fe_index=DoFHandler< dimension_, space_dimension_ >::invalid_fe_index) const
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
FEFaceValuesBase(const unsigned int dofs_per_cell, const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature)
std::size_t memory_consumption() const
const std::vector< Tensor< 1, spacedim > > & get_boundary_forms() const
void initialize(const UpdateFlags update_flags)
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell, const unsigned int face_no)
void do_reinit(const unsigned int face_no)
FEFaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature, const UpdateFlags update_flags)
FESubfaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &face_quadrature, const UpdateFlags update_flags)
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell, const unsigned int face_no, const unsigned int subface_no)
void initialize(const UpdateFlags update_flags)
void do_reinit(const unsigned int face_no, const unsigned int subface_no)
bool is_initialized() const
void get_interpolated_dof_values(const VectorType &in, Vector< typename VectorType::value_type > &out) const
types::global_dof_index n_dofs_for_dof_handler() const
CellSimilarity::Similarity cell_similarity
CellIteratorContainer present_cell
::internal::FEValuesViews::Cache< dim, spacedim > fe_values_views_cache
void get_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
FEValuesBase(const unsigned int n_q_points, const unsigned int dofs_per_cell, const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe)
virtual ~FEValuesBase() override
const SmartPointer< const Mapping< dim, spacedim >, FEValuesBase< dim, spacedim > > mapping
const unsigned int dofs_per_cell
void check_cell_similarity(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
const SmartPointer< const FiniteElement< dim, spacedim >, FEValuesBase< dim, spacedim > > fe
const unsigned int n_quadrature_points
CellSimilarity::Similarity get_cell_similarity() const
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors() const
std::size_t memory_consumption() const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > dof_indices() const
const Triangulation< dim, spacedim >::cell_iterator get_cell() const
void get_function_laplacians(const InputVector &fe_function, std::vector< typename InputVector::value_type > &laplacians) const
UpdateFlags compute_update_flags(const UpdateFlags update_flags) const
void get_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &gradients) const
void invalidate_present_cell()
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > finite_element_output
void get_function_hessians(const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type > > &hessians) const
const FiniteElement< dim, spacedim > & get_fe() const
void maybe_invalidate_previous_present_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
const unsigned int max_n_quadrature_points
void get_function_third_derivatives(const InputVector &fe_function, std::vector< Tensor< 3, spacedim, typename InputVector::value_type > > &third_derivatives) const
void get_function_third_derivatives_from_local_dof_values(const InputVector &dof_values, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &third_derivatives) const
typename ProductType< Number, hessian_type >::type solution_hessian_type
const unsigned int component
void get_function_hessians_from_local_dof_values(const InputVector &dof_values, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
void get_function_values_from_local_dof_values(const InputVector &dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
std::vector< ShapeFunctionData > shape_function_data
void get_function_laplacians(const InputVector &fe_function, std::vector< solution_laplacian_type< typename InputVector::value_type > > &laplacians) const
typename ProductType< Number, value_type >::type solution_laplacian_type
void get_function_third_derivatives(const InputVector &fe_function, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &third_derivatives) const
void get_function_hessians(const InputVector &fe_function, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
typename ProductType< Number, third_derivative_type >::type solution_third_derivative_type
typename ProductType< Number, value_type >::type solution_value_type
typename ProductType< Number, gradient_type >::type solution_gradient_type
void get_function_laplacians_from_local_dof_values(const InputVector &dof_values, std::vector< solution_laplacian_type< typename InputVector::value_type > > &laplacians) const
void get_function_gradients(const InputVector &fe_function, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
void get_function_gradients_from_local_dof_values(const InputVector &dof_values, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
void get_function_values(const InputVector &fe_function, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
typename ProductType< Number, value_type >::type solution_value_type
typename ProductType< Number, divergence_type >::type solution_divergence_type
typename ProductType< Number, value_type >::type solution_value_type
typename ProductType< Number, divergence_type >::type solution_divergence_type
typename ProductType< Number, gradient_type >::type solution_gradient_type
void get_function_hessians(const InputVector &fe_function, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
typename ProductType< Number, third_derivative_type >::type solution_third_derivative_type
void get_function_symmetric_gradients_from_local_dof_values(const InputVector &dof_values, std::vector< solution_symmetric_gradient_type< typename InputVector::value_type > > &symmetric_gradients) const
void get_function_gradients_from_local_dof_values(const InputVector &dof_values, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
typename ProductType< Number, divergence_type >::type solution_divergence_type
typename ProductType< Number, hessian_type >::type solution_hessian_type
typename ProductType< Number, symmetric_gradient_type >::type solution_symmetric_gradient_type
void get_function_laplacians_from_local_dof_values(const InputVector &dof_values, std::vector< solution_laplacian_type< typename InputVector::value_type > > &laplacians) const
typename ProductType< Number, gradient_type >::type solution_gradient_type
void get_function_third_derivatives(const InputVector &fe_function, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &third_derivatives) const
typename ProductType< Number, value_type >::type solution_value_type
void get_function_curls_from_local_dof_values(const InputVector &dof_values, std::vector< solution_curl_type< typename InputVector::value_type > > &curls) const
void get_function_values(const InputVector &fe_function, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
const unsigned int first_vector_component
typename ProductType< Number, curl_type >::type solution_curl_type
std::vector< ShapeFunctionData > shape_function_data
void get_function_gradients(const InputVector &fe_function, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
void get_function_symmetric_gradients(const InputVector &fe_function, std::vector< solution_symmetric_gradient_type< typename InputVector::value_type > > &symmetric_gradients) const
void get_function_hessians_from_local_dof_values(const InputVector &dof_values, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
void get_function_laplacians(const InputVector &fe_function, std::vector< solution_laplacian_type< typename InputVector::value_type > > &laplacians) const
typename ProductType< Number, value_type >::type solution_laplacian_type
void get_function_values_from_local_dof_values(const InputVector &dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
void get_function_curls(const InputVector &fe_function, std::vector< solution_curl_type< typename InputVector::value_type > > &curls) const
void get_function_divergences_from_local_dof_values(const InputVector &dof_values, std::vector< solution_divergence_type< typename InputVector::value_type > > &divergences) const
void get_function_third_derivatives_from_local_dof_values(const InputVector &dof_values, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &third_derivatives) const
void get_function_divergences(const InputVector &fe_function, std::vector< solution_divergence_type< typename InputVector::value_type > > &divergences) const
FEValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags update_flags)
void initialize(const UpdateFlags update_flags)
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
std::size_t memory_consumption() const
unsigned int n_dofs_per_cell() const
unsigned int n_components() const
virtual std::unique_ptr< InternalDataBase > get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
const ComponentMask & get_nonzero_components(const unsigned int i) const
bool is_primitive() const
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const =0
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
unsigned int n_nonzero_components(const unsigned int i) const
bool is_element(const size_type index) const
Abstract base class for mapping classes.
static unsigned int n_threads()
static constexpr unsigned int component_to_unrolled_index(const TableIndices< rank_ > &indices)
static constexpr TableIndices< rank_ > unrolled_to_component_indices(const unsigned int i)
constexpr SymmetricTensor()=default
static constexpr TableIndices< rank_ > unrolled_to_component_indices(const unsigned int i)
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
@ update_jacobian_pushed_forward_2nd_derivatives
@ update_jacobian_pushed_forward_grads
@ update_hessians
Second derivatives of shape functions.
@ update_jacobian_3rd_derivatives
@ update_values
Shape function values.
@ update_jacobian_grads
Gradient of volume element.
@ update_normal_vectors
Normal vectors.
@ update_3rd_derivatives
Third derivatives of shape functions.
@ update_JxW_values
Transformed quadrature weights.
@ update_jacobians
Volume element.
@ update_inverse_jacobians
Volume element.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
@ update_default
No update.
@ update_jacobian_pushed_forward_3rd_derivatives
@ update_boundary_forms
Outer normal vector, not normalized.
@ update_jacobian_2nd_derivatives
static ::ExceptionBase & ExcNotMultiple(int arg1, int arg2)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcNotReinited()
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcIndexRange(std::size_t arg1, std::size_t arg2, std::size_t arg3)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Task< RT > new_task(const std::function< RT()> &function)
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
void do_function_laplacians(const ArrayView< Number > &dof_values, const Table< 2, ::Tensor< 2, spacedim > > &shape_hessians, const std::vector< typename Scalar< dim, spacedim >::ShapeFunctionData > &shape_function_data, std::vector< typename Scalar< dim, spacedim >::template solution_laplacian_type< Number > > &laplacians)
void do_function_values(const ArrayView< Number > &dof_values, const Table< 2, double > &shape_values, const std::vector< typename Scalar< dim, spacedim >::ShapeFunctionData > &shape_function_data, std::vector< typename ProductType< Number, double >::type > &values)
void do_function_gradients(const ArrayView< Number > &dof_values, const Table< 2, ::Tensor< 1, spacedim > > &shape_gradients, const std::vector< typename Tensor< 2, dim, spacedim >::ShapeFunctionData > &shape_function_data, std::vector< typename Tensor< 2, dim, spacedim >::template solution_gradient_type< Number > > &gradients)
void do_function_divergences(const ArrayView< Number > &dof_values, const Table< 2, ::Tensor< 1, spacedim > > &shape_gradients, const std::vector< typename Vector< dim, spacedim >::ShapeFunctionData > &shape_function_data, std::vector< typename Vector< dim, spacedim >::template solution_divergence_type< Number > > &divergences)
void do_function_curls(const ArrayView< Number > &dof_values, const Table< 2, ::Tensor< 1, spacedim > > &shape_gradients, const std::vector< typename Vector< dim, spacedim >::ShapeFunctionData > &shape_function_data, std::vector< typename ProductType< Number, typename ::internal::CurlType< spacedim >::type >::type > &curls)
void do_function_symmetric_gradients(const ArrayView< Number > &dof_values, const Table< 2, ::Tensor< 1, spacedim > > &shape_gradients, const std::vector< typename Vector< dim, spacedim >::ShapeFunctionData > &shape_function_data, std::vector< typename ProductType< Number, ::SymmetricTensor< 2, spacedim > >::type > &symmetric_gradients)
void do_function_derivatives(const ArrayView< Number > &dof_values, const Table< 2, ::Tensor< order, spacedim > > &shape_derivatives, const std::vector< typename Scalar< dim, spacedim >::ShapeFunctionData > &shape_function_data, std::vector< typename ProductType< Number, ::Tensor< order, spacedim > >::type > &derivatives)
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
void do_function_derivatives(const Number *dof_values_ptr, const ::Table< 2, Tensor< order, spacedim > > &shape_derivatives, std::vector< Tensor< order, spacedim, Number > > &derivatives)
void do_function_values(const Number2 *dof_values_ptr, const ::Table< 2, double > &shape_values, std::vector< Number > &values)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
void do_function_laplacians(const Number2 *dof_values_ptr, const ::Table< 2, Tensor< 2, spacedim > > &shape_hessians, std::vector< Number > &laplacians)
VectorType::value_type get_vector_element(const VectorType &vector, const types::global_dof_index cell_number)
std::vector< unsigned int > make_shape_function_to_row_table(const FiniteElement< dim, spacedim > &fe)
static const unsigned int invalid_unsigned_int
typename internal::ProductTypeImpl< typename std::decay< T >::type, typename std::decay< U >::type >::type type
static VectorType::value_type get(const VectorType &V, const types::global_dof_index i)
Cache(const FEValuesBase< dim, spacedim > &fe_values)
static constexpr const T & value(const T &t)
constexpr SymmetricTensor< 2, dim, Number > symmetrize(const Tensor< 2, dim, Number > &t)
constexpr Number trace(const SymmetricTensor< 2, dim2, Number > &)