16 #ifndef dealii_fe_values_h 17 #define dealii_fe_values_h 20 #include <deal.II/base/config.h> 22 #include <deal.II/base/derivative_form.h> 23 #include <deal.II/base/exceptions.h> 24 #include <deal.II/base/point.h> 25 #include <deal.II/base/quadrature.h> 26 #include <deal.II/base/subscriptor.h> 29 #include <deal.II/dofs/dof_accessor.h> 30 #include <deal.II/dofs/dof_handler.h> 32 #include <deal.II/fe/fe.h> 33 #include <deal.II/fe/fe_update_flags.h> 34 #include <deal.II/fe/fe_values_extractors.h> 35 #include <deal.II/fe/mapping.h> 37 #include <deal.II/grid/tria.h> 38 #include <deal.II/grid/tria_iterator.h> 40 #include <deal.II/hp/dof_handler.h> 44 #include <type_traits> 50 #ifdef DEAL_II_WITH_PETSC 54 DEAL_II_NAMESPACE_OPEN
56 template <
int dim,
int spacedim = dim>
65 template <
int dim,
class NumberType =
double>
74 template <
class NumberType>
86 template <
class NumberType>
98 template <
class NumberType>
141 template <
int dim,
int spacedim = dim>
177 template <
typename Number>
282 value(
const unsigned int shape_function,
const unsigned int q_point)
const;
295 gradient(
const unsigned int shape_function,
296 const unsigned int q_point)
const;
309 hessian(
const unsigned int shape_function,
310 const unsigned int q_point)
const;
324 const unsigned int q_point)
const;
343 template <
class InputVector>
346 const InputVector &fe_function,
348 typename InputVector::value_type>::type>
373 template <
class InputVector>
376 const InputVector &dof_values,
398 template <
class InputVector>
401 const InputVector &fe_function,
403 typename InputVector::value_type>::type>
409 template <
class InputVector>
412 const InputVector &dof_values,
434 template <
class InputVector>
437 const InputVector &fe_function,
439 typename InputVector::value_type>::type>
445 template <
class InputVector>
448 const InputVector &dof_values,
472 template <
class InputVector>
475 const InputVector &fe_function,
477 typename InputVector::value_type>::type>
483 template <
class InputVector>
486 const InputVector &dof_values,
510 template <
class InputVector>
513 const InputVector &fe_function,
515 typename InputVector::value_type>::type>
516 &third_derivatives)
const;
521 template <
class InputVector>
524 const InputVector & dof_values,
578 template <
int dim,
int spacedim = dim>
624 using curl_type = typename ::internal::CurlType<spacedim>::type;
644 template <
typename Number>
748 unsigned int single_nonzero_component_index;
791 value(
const unsigned int shape_function,
const unsigned int q_point)
const;
807 gradient(
const unsigned int shape_function,
808 const unsigned int q_point)
const;
827 const unsigned int q_point)
const;
841 const unsigned int q_point)
const;
864 curl(
const unsigned int shape_function,
const unsigned int q_point)
const;
877 hessian(
const unsigned int shape_function,
878 const unsigned int q_point)
const;
892 const unsigned int q_point)
const;
911 template <
class InputVector>
914 const InputVector &fe_function,
916 typename InputVector::value_type>::type>
941 template <
class InputVector>
944 const InputVector &dof_values,
966 template <
class InputVector>
969 const InputVector &fe_function,
971 typename InputVector::value_type>::type>
977 template <
class InputVector>
980 const InputVector &dof_values,
1008 template <
class InputVector>
1011 const InputVector &fe_function,
1013 typename InputVector::value_type>::type>
1014 &symmetric_gradients)
const;
1019 template <
class InputVector>
1022 const InputVector & dof_values,
1044 template <
class InputVector>
1047 const InputVector &fe_function,
1049 typename InputVector::value_type>::type>
1050 &divergences)
const;
1055 template <
class InputVector>
1058 const InputVector &dof_values,
1061 &divergences)
const;
1081 template <
class InputVector>
1084 const InputVector &fe_function,
1086 typename ProductType<curl_type, typename InputVector::value_type>::type>
1092 template <
class InputVector>
1095 const InputVector &dof_values,
1117 template <
class InputVector>
1120 const InputVector &fe_function,
1122 typename InputVector::value_type>::type>
1128 template <
class InputVector>
1131 const InputVector &dof_values,
1154 template <
class InputVector>
1157 const InputVector &fe_function,
1159 typename InputVector::value_type>::type>
1165 template <
class InputVector>
1168 const InputVector &dof_values,
1191 template <
class InputVector>
1194 const InputVector &fe_function,
1196 typename InputVector::value_type>::type>
1197 &third_derivatives)
const;
1202 template <
class InputVector>
1205 const InputVector & dof_values,
1228 template <
int rank,
int dim,
int spacedim = dim>
1255 template <
int dim,
int spacedim>
1282 template <
typename Number>
1306 struct ShapeFunctionData
1316 bool is_nonzero_shape_function_component
1317 [value_type::n_independent_components];
1328 unsigned int row_index[value_type::n_independent_components];
1361 const unsigned int first_tensor_component);
1388 value(
const unsigned int shape_function,
const unsigned int q_point)
const;
1404 divergence(
const unsigned int shape_function,
1405 const unsigned int q_point)
const;
1424 template <
class InputVector>
1426 get_function_values(
1427 const InputVector &fe_function,
1429 typename InputVector::value_type>::type>
1454 template <
class InputVector>
1456 get_function_values_from_local_dof_values(
1457 const InputVector &dof_values,
1483 template <
class InputVector>
1485 get_function_divergences(
1486 const InputVector &fe_function,
1488 typename InputVector::value_type>::type>
1489 &divergences)
const;
1494 template <
class InputVector>
1496 get_function_divergences_from_local_dof_values(
1497 const InputVector &dof_values,
1500 &divergences)
const;
1521 template <
int rank,
int dim,
int spacedim = dim>
1544 template <
int dim,
int spacedim>
1569 template <
typename Number>
1601 struct ShapeFunctionData
1611 bool is_nonzero_shape_function_component
1612 [value_type::n_independent_components];
1623 unsigned int row_index[value_type::n_independent_components];
1656 const unsigned int first_tensor_component);
1683 value(
const unsigned int shape_function,
const unsigned int q_point)
const;
1699 divergence(
const unsigned int shape_function,
1700 const unsigned int q_point)
const;
1716 gradient(
const unsigned int shape_function,
1717 const unsigned int q_point)
const;
1736 template <
class InputVector>
1738 get_function_values(
1739 const InputVector &fe_function,
1741 typename InputVector::value_type>::type>
1766 template <
class InputVector>
1768 get_function_values_from_local_dof_values(
1769 const InputVector &dof_values,
1795 template <
class InputVector>
1797 get_function_divergences(
1798 const InputVector &fe_function,
1800 typename InputVector::value_type>::type>
1801 &divergences)
const;
1806 template <
class InputVector>
1808 get_function_divergences_from_local_dof_values(
1809 const InputVector &dof_values,
1812 &divergences)
const;
1830 template <
class InputVector>
1832 get_function_gradients(
1833 const InputVector &fe_function,
1835 typename InputVector::value_type>::type>
1841 template <
class InputVector>
1843 get_function_gradients_from_local_dof_values(
1844 const InputVector &dof_values,
1880 template <
int dim,
int spacedim,
typename Extractor>
1891 template <
int dim,
int spacedim>
1894 using type = typename ::FEValuesViews::Scalar<dim, spacedim>;
1904 template <
int dim,
int spacedim>
1907 using type = typename ::FEValuesViews::Vector<dim, spacedim>;
1917 template <
int dim,
int spacedim,
int rank>
1920 using type = typename ::FEValuesViews::Tensor<rank, dim, spacedim>;
1930 template <
int dim,
int spacedim,
int rank>
1934 typename ::FEValuesViews::SymmetricTensor<rank, dim, spacedim>;
1944 template <
int dim,
int spacedim>
1951 std::vector<::FEValuesViews::Scalar<dim, spacedim>>
scalars;
1952 std::vector<::FEValuesViews::Vector<dim, spacedim>> vectors;
1953 std::vector<::FEValuesViews::SymmetricTensor<2, dim, spacedim>>
1954 symmetric_second_order_tensors;
1955 std::vector<::FEValuesViews::Tensor<2, dim, spacedim>>
1956 second_order_tensors;
1974 template <
int dim,
int spacedim,
typename Extractor>
2082 template <
int dim,
int spacedim>
2168 const unsigned int point_no)
const;
2192 const unsigned int point_no,
2193 const unsigned int component)
const;
2242 const unsigned int point_no,
2243 const unsigned int component)
const;
2266 const unsigned int point_no)
const;
2286 const unsigned int point_no,
2287 const unsigned int component)
const;
2310 const unsigned int point_no)
const;
2330 const unsigned int point_no,
2331 const unsigned int component)
const;
2373 template <
class InputVector>
2376 const InputVector & fe_function,
2377 std::vector<typename InputVector::value_type> &values)
const;
2392 template <
class InputVector>
2395 const InputVector & fe_function,
2396 std::vector<Vector<typename InputVector::value_type>> &values)
const;
2416 template <
class InputVector>
2419 const InputVector & fe_function,
2421 std::vector<typename InputVector::value_type> & values)
const;
2444 template <
class InputVector>
2447 const InputVector & fe_function,
2449 std::vector<Vector<typename InputVector::value_type>> &values)
const;
2482 template <
class InputVector>
2485 const InputVector & fe_function,
2487 ArrayView<std::vector<typename InputVector::value_type>> values,
2488 const bool quadrature_points_fastest)
const;
2530 template <
class InputVector>
2533 const InputVector &fe_function,
2553 template <
class InputVector>
2556 const InputVector &fe_function,
2567 template <
class InputVector>
2570 const InputVector & fe_function,
2581 template <
class InputVector>
2584 const InputVector & fe_function,
2589 bool quadrature_points_fastest =
false)
const;
2634 template <
class InputVector>
2637 const InputVector &fe_function,
2658 template <
class InputVector>
2661 const InputVector &fe_function,
2665 bool quadrature_points_fastest =
false)
const;
2671 template <
class InputVector>
2674 const InputVector & fe_function,
2685 template <
class InputVector>
2688 const InputVector & fe_function,
2693 bool quadrature_points_fastest =
false)
const;
2735 template <
class InputVector>
2738 const InputVector & fe_function,
2739 std::vector<typename InputVector::value_type> &laplacians)
const;
2760 template <
class InputVector>
2763 const InputVector & fe_function,
2764 std::vector<Vector<typename InputVector::value_type>> &laplacians)
const;
2772 template <
class InputVector>
2775 const InputVector & fe_function,
2777 std::vector<typename InputVector::value_type> & laplacians)
const;
2785 template <
class InputVector>
2788 const InputVector & fe_function,
2790 std::vector<Vector<typename InputVector::value_type>> &laplacians)
const;
2798 template <
class InputVector>
2801 const InputVector & fe_function,
2803 std::vector<std::vector<typename InputVector::value_type>> &laplacians,
2804 bool quadrature_points_fastest =
false)
const;
2848 template <
class InputVector>
2851 const InputVector &fe_function,
2853 &third_derivatives)
const;
2873 template <
class InputVector>
2876 const InputVector &fe_function,
2879 & third_derivatives,
2880 bool quadrature_points_fastest =
false)
const;
2886 template <
class InputVector>
2889 const InputVector & fe_function,
2892 &third_derivatives)
const;
2900 template <
class InputVector>
2903 const InputVector & fe_function,
2908 bool quadrature_points_fastest =
false)
const;
2927 const std::vector<Point<spacedim>> &
2951 const std::vector<double> &
2969 const std::vector<DerivativeForm<1, dim, spacedim>> &
2988 const std::vector<DerivativeForm<2, dim, spacedim>> &
3008 const std::vector<Tensor<3, spacedim>> &
3027 const std::vector<DerivativeForm<3, dim, spacedim>> &
3049 const std::vector<Tensor<4, spacedim>> &
3069 const std::vector<DerivativeForm<4, dim, spacedim>> &
3091 const std::vector<Tensor<5, spacedim>> &
3109 const std::vector<DerivativeForm<1, spacedim, dim>> &
3139 const std::vector<Tensor<1, spacedim>> &
3149 const std::vector<Tensor<1, spacedim>> &
3258 <<
"You are requesting information from an FEValues/FEFaceValues/FESubfaceValues " 3259 <<
"object for which this kind of information has not been computed. What " 3260 <<
"information these objects compute is determined by the update_* flags you " 3261 <<
"pass to the constructor. Here, the operation you are attempting requires " 3263 <<
"> flag to be set, but it was apparently not specified " 3264 <<
"upon construction.");
3274 "The FiniteElement you provided to FEValues and the FiniteElement that belongs " 3275 "to the DoFHandler that provided the cell iterator do not match.");
3283 <<
"The shape function with index " << arg1
3284 <<
" is not primitive, i.e. it is vector-valued and " 3285 <<
"has more than one non-zero vector component. This " 3286 <<
"function cannot be called for these shape functions. " 3287 <<
"Maybe you want to use the same function with the " 3288 <<
"_component suffix?");
3298 "The given FiniteElement is not a primitive element but the requested operation " 3299 "only works for those. See FiniteElement::is_primitive() for more information.");
3330 class CellIteratorBase;
3336 template <
typename CI>
3397 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
3421 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
3479 template <
int,
int,
int>
3480 friend class FEValuesViews::SymmetricTensor;
3481 template <
int,
int,
int>
3482 friend class FEValuesViews::Tensor;
3497 template <
int dim,
int spacedim = dim>
3531 template <
template <
int,
int>
class DoFHandlerType,
bool level_dof_access>
3534 level_dof_access>> &cell);
3616 template <
int dim,
int spacedim = dim>
3660 const std::vector<Tensor<1, spacedim>> &
3713 template <
int dim,
int spacedim = dim>
3723 static const unsigned int space_dimension = spacedim;
3753 template <
template <
int,
int>
class DoFHandlerType,
bool level_dof_access>
3756 level_dof_access>> &cell,
3757 const unsigned int face_no);
3774 const unsigned int face_no);
3828 template <
int dim,
int spacedim = dim>
3872 template <
template <
int,
int>
class DoFHandlerType,
bool level_dof_access>
3875 level_dof_access>> &cell,
3876 const unsigned int face_no,
3877 const unsigned int subface_no);
3894 const unsigned int face_no,
3895 const unsigned int subface_no);
3942 do_reinit(
const unsigned int face_no,
const unsigned int subface_no);
3953 template <
int dim,
int spacedim>
3956 const unsigned int q_point)
const 3958 Assert(shape_function < fe_values->fe->dofs_per_cell,
3959 ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
3963 "update_values"))));
3968 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
3969 return fe_values->finite_element_output.shape_values(
3970 shape_function_data[shape_function].row_index, q_point);
3977 template <
int dim,
int spacedim>
3980 const unsigned int q_point)
const 3982 Assert(shape_function < fe_values->fe->dofs_per_cell,
3983 ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
3986 "update_gradients")));
3991 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
3992 return fe_values->finite_element_output
3993 .shape_gradients[shape_function_data[shape_function].row_index]
3996 return gradient_type();
4001 template <
int dim,
int spacedim>
4004 const unsigned int q_point)
const 4006 Assert(shape_function < fe_values->fe->dofs_per_cell,
4007 ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4010 "update_hessians")));
4015 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
4016 return fe_values->finite_element_output
4017 .shape_hessians[shape_function_data[shape_function].row_index][q_point];
4019 return hessian_type();
4024 template <
int dim,
int spacedim>
4027 const unsigned int q_point)
const 4029 Assert(shape_function < fe_values->fe->dofs_per_cell,
4030 ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4033 "update_3rd_derivatives")));
4038 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
4039 return fe_values->finite_element_output
4040 .shape_3rd_derivatives[shape_function_data[shape_function].row_index]
4043 return third_derivative_type();
4048 template <
int dim,
int spacedim>
4051 const unsigned int q_point)
const 4053 Assert(shape_function < fe_values->fe->dofs_per_cell,
4054 ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4061 shape_function_data[shape_function].single_nonzero_component;
4063 return value_type();
4066 value_type return_value;
4067 return_value[shape_function_data[shape_function]
4068 .single_nonzero_component_index] =
4069 fe_values->finite_element_output.shape_values(snc, q_point);
4070 return return_value;
4074 value_type return_value;
4075 for (
unsigned int d = 0;
d < dim; ++
d)
4076 if (shape_function_data[shape_function]
4077 .is_nonzero_shape_function_component[d])
4078 return_value[
d] = fe_values->finite_element_output.shape_values(
4079 shape_function_data[shape_function].row_index[d], q_point);
4081 return return_value;
4087 template <
int dim,
int spacedim>
4090 const unsigned int q_point)
const 4092 Assert(shape_function < fe_values->fe->dofs_per_cell,
4093 ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4096 "update_gradients")));
4100 shape_function_data[shape_function].single_nonzero_component;
4102 return gradient_type();
4105 gradient_type return_value;
4106 return_value[shape_function_data[shape_function]
4107 .single_nonzero_component_index] =
4108 fe_values->finite_element_output.shape_gradients[snc][q_point];
4109 return return_value;
4113 gradient_type return_value;
4114 for (
unsigned int d = 0;
d < dim; ++
d)
4115 if (shape_function_data[shape_function]
4116 .is_nonzero_shape_function_component[d])
4118 fe_values->finite_element_output.shape_gradients
4119 [shape_function_data[shape_function].row_index[
d]][q_point];
4121 return return_value;
4127 template <
int dim,
int spacedim>
4130 const unsigned int q_point)
const 4133 Assert(shape_function < fe_values->fe->dofs_per_cell,
4134 ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4137 "update_gradients")));
4141 shape_function_data[shape_function].single_nonzero_component;
4143 return divergence_type();
4145 return fe_values->finite_element_output
4146 .shape_gradients[snc][q_point][shape_function_data[shape_function]
4147 .single_nonzero_component_index];
4150 divergence_type return_value = 0;
4151 for (
unsigned int d = 0;
d < dim; ++
d)
4152 if (shape_function_data[shape_function]
4153 .is_nonzero_shape_function_component[d])
4155 fe_values->finite_element_output.shape_gradients
4156 [shape_function_data[shape_function].row_index[
d]][q_point][
d];
4158 return return_value;
4164 template <
int dim,
int spacedim>
4167 const unsigned int q_point)
const 4171 Assert(shape_function < fe_values->fe->dofs_per_cell,
4172 ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4175 "update_gradients")));
4178 shape_function_data[shape_function].single_nonzero_component;
4190 "Computing the curl in 1d is not a useful operation"));
4198 curl_type return_value;
4201 if (shape_function_data[shape_function]
4202 .single_nonzero_component_index == 0)
4204 -1.0 * fe_values->finite_element_output
4205 .shape_gradients[snc][q_point][1];
4207 return_value[0] = fe_values->finite_element_output
4208 .shape_gradients[snc][q_point][0];
4210 return return_value;
4215 curl_type return_value;
4217 return_value[0] = 0.0;
4219 if (shape_function_data[shape_function]
4220 .is_nonzero_shape_function_component[0])
4222 fe_values->finite_element_output
4223 .shape_gradients[shape_function_data[shape_function]
4224 .row_index[0]][q_point][1];
4226 if (shape_function_data[shape_function]
4227 .is_nonzero_shape_function_component[1])
4229 fe_values->finite_element_output
4230 .shape_gradients[shape_function_data[shape_function]
4231 .row_index[1]][q_point][0];
4233 return return_value;
4241 curl_type return_value;
4243 switch (shape_function_data[shape_function]
4244 .single_nonzero_component_index)
4248 return_value[0] = 0;
4249 return_value[1] = fe_values->finite_element_output
4250 .shape_gradients[snc][q_point][2];
4252 -1.0 * fe_values->finite_element_output
4253 .shape_gradients[snc][q_point][1];
4254 return return_value;
4260 -1.0 * fe_values->finite_element_output
4261 .shape_gradients[snc][q_point][2];
4262 return_value[1] = 0;
4263 return_value[2] = fe_values->finite_element_output
4264 .shape_gradients[snc][q_point][0];
4265 return return_value;
4270 return_value[0] = fe_values->finite_element_output
4271 .shape_gradients[snc][q_point][1];
4273 -1.0 * fe_values->finite_element_output
4274 .shape_gradients[snc][q_point][0];
4275 return_value[2] = 0;
4276 return return_value;
4283 curl_type return_value;
4285 for (
unsigned int i = 0; i < dim; ++i)
4286 return_value[i] = 0.0;
4288 if (shape_function_data[shape_function]
4289 .is_nonzero_shape_function_component[0])
4292 fe_values->finite_element_output
4293 .shape_gradients[shape_function_data[shape_function]
4294 .row_index[0]][q_point][2];
4296 fe_values->finite_element_output
4297 .shape_gradients[shape_function_data[shape_function]
4298 .row_index[0]][q_point][1];
4301 if (shape_function_data[shape_function]
4302 .is_nonzero_shape_function_component[1])
4305 fe_values->finite_element_output
4306 .shape_gradients[shape_function_data[shape_function]
4307 .row_index[1]][q_point][2];
4309 fe_values->finite_element_output
4310 .shape_gradients[shape_function_data[shape_function]
4311 .row_index[1]][q_point][0];
4314 if (shape_function_data[shape_function]
4315 .is_nonzero_shape_function_component[2])
4318 fe_values->finite_element_output
4319 .shape_gradients[shape_function_data[shape_function]
4320 .row_index[2]][q_point][1];
4322 fe_values->finite_element_output
4323 .shape_gradients[shape_function_data[shape_function]
4324 .row_index[2]][q_point][0];
4327 return return_value;
4338 template <
int dim,
int spacedim>
4341 const unsigned int q_point)
const 4344 Assert(shape_function < fe_values->fe->dofs_per_cell,
4345 ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4348 "update_hessians")));
4352 shape_function_data[shape_function].single_nonzero_component;
4354 return hessian_type();
4357 hessian_type return_value;
4358 return_value[shape_function_data[shape_function]
4359 .single_nonzero_component_index] =
4360 fe_values->finite_element_output.shape_hessians[snc][q_point];
4361 return return_value;
4365 hessian_type return_value;
4366 for (
unsigned int d = 0;
d < dim; ++
d)
4367 if (shape_function_data[shape_function]
4368 .is_nonzero_shape_function_component[d])
4370 fe_values->finite_element_output.shape_hessians
4371 [shape_function_data[shape_function].row_index[
d]][q_point];
4373 return return_value;
4379 template <
int dim,
int spacedim>
4382 const unsigned int q_point)
const 4385 Assert(shape_function < fe_values->fe->dofs_per_cell,
4386 ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4389 "update_3rd_derivatives")));
4393 shape_function_data[shape_function].single_nonzero_component;
4395 return third_derivative_type();
4398 third_derivative_type return_value;
4399 return_value[shape_function_data[shape_function]
4400 .single_nonzero_component_index] =
4401 fe_values->finite_element_output.shape_3rd_derivatives[snc][q_point];
4402 return return_value;
4406 third_derivative_type return_value;
4407 for (
unsigned int d = 0;
d < dim; ++
d)
4408 if (shape_function_data[shape_function]
4409 .is_nonzero_shape_function_component[d])
4411 fe_values->finite_element_output.shape_3rd_derivatives
4412 [shape_function_data[shape_function].row_index[
d]][q_point];
4414 return return_value;
4426 inline ::SymmetricTensor<2, 1>
4427 symmetrize_single_row(
const unsigned int n, const ::Tensor<1, 1> &t)
4437 inline ::SymmetricTensor<2, 2>
4438 symmetrize_single_row(
const unsigned int n, const ::Tensor<1, 2> &t)
4444 return {{t[0], 0, t[1] / 2}};
4448 return {{0, t[1], t[0] / 2}};
4460 inline ::SymmetricTensor<2, 3>
4461 symmetrize_single_row(
const unsigned int n, const ::Tensor<1, 3> &t)
4467 return {{t[0], 0, 0, t[1] / 2, t[2] / 2, 0}};
4471 return {{0, t[1], 0, t[0] / 2, 0, t[2] / 2}};
4475 return {{0, 0, t[2], 0, t[0] / 2, t[1] / 2}};
4488 template <
int dim,
int spacedim>
4491 const unsigned int q_point)
const 4493 Assert(shape_function < fe_values->fe->dofs_per_cell,
4494 ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4497 "update_gradients")));
4501 shape_function_data[shape_function].single_nonzero_component;
4503 return symmetric_gradient_type();
4505 return internal::symmetrize_single_row(
4506 shape_function_data[shape_function].single_nonzero_component_index,
4507 fe_values->finite_element_output.shape_gradients[snc][q_point]);
4510 gradient_type return_value;
4511 for (
unsigned int d = 0;
d < dim; ++
d)
4512 if (shape_function_data[shape_function]
4513 .is_nonzero_shape_function_component[d])
4515 fe_values->finite_element_output.shape_gradients
4516 [shape_function_data[shape_function].row_index[
d]][q_point];
4524 template <
int dim,
int spacedim>
4527 const unsigned int q_point)
const 4529 Assert(shape_function < fe_values->fe->dofs_per_cell,
4530 ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4538 shape_function_data[shape_function].single_nonzero_component;
4543 return value_type();
4547 value_type return_value;
4548 const unsigned int comp =
4549 shape_function_data[shape_function].single_nonzero_component_index;
4550 return_value[value_type::unrolled_to_component_indices(comp)] =
4551 fe_values->finite_element_output.shape_values(snc, q_point);
4552 return return_value;
4556 value_type return_value;
4557 for (
unsigned int d = 0;
d < value_type::n_independent_components; ++
d)
4558 if (shape_function_data[shape_function]
4559 .is_nonzero_shape_function_component[d])
4560 return_value[value_type::unrolled_to_component_indices(d)] =
4561 fe_values->finite_element_output.shape_values(
4562 shape_function_data[shape_function].row_index[d], q_point);
4563 return return_value;
4569 template <
int dim,
int spacedim>
4572 const unsigned int shape_function,
4573 const unsigned int q_point)
const 4575 Assert(shape_function < fe_values->fe->dofs_per_cell,
4576 ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4579 "update_gradients")));
4582 shape_function_data[shape_function].single_nonzero_component;
4587 return divergence_type();
4610 const unsigned int comp =
4611 shape_function_data[shape_function].single_nonzero_component_index;
4612 const unsigned int ii =
4613 value_type::unrolled_to_component_indices(comp)[0];
4614 const unsigned int jj =
4615 value_type::unrolled_to_component_indices(comp)[1];
4628 const ::Tensor<1, spacedim> &phi_grad =
4629 fe_values->finite_element_output.shape_gradients[snc][q_point];
4631 divergence_type return_value;
4632 return_value[ii] = phi_grad[jj];
4635 return_value[jj] = phi_grad[ii];
4637 return return_value;
4642 divergence_type return_value;
4643 return return_value;
4649 template <
int dim,
int spacedim>
4652 const unsigned int q_point)
const 4654 Assert(shape_function < fe_values->fe->dofs_per_cell,
4655 ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4663 shape_function_data[shape_function].single_nonzero_component;
4668 return value_type();
4672 value_type return_value;
4673 const unsigned int comp =
4674 shape_function_data[shape_function].single_nonzero_component_index;
4677 return_value[indices] =
4678 fe_values->finite_element_output.shape_values(snc, q_point);
4679 return return_value;
4683 value_type return_value;
4684 for (
unsigned int d = 0;
d < dim * dim; ++
d)
4685 if (shape_function_data[shape_function]
4686 .is_nonzero_shape_function_component[d])
4690 return_value[indices] =
4691 fe_values->finite_element_output.shape_values(
4692 shape_function_data[shape_function].row_index[d], q_point);
4694 return return_value;
4700 template <
int dim,
int spacedim>
4703 const unsigned int q_point)
const 4705 Assert(shape_function < fe_values->fe->dofs_per_cell,
4706 ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4709 "update_gradients")));
4712 shape_function_data[shape_function].single_nonzero_component;
4717 return divergence_type();
4731 const unsigned int comp =
4732 shape_function_data[shape_function].single_nonzero_component_index;
4735 const unsigned int ii = indices[0];
4736 const unsigned int jj = indices[1];
4738 const ::Tensor<1, spacedim> &phi_grad =
4739 fe_values->finite_element_output.shape_gradients[snc][q_point];
4741 divergence_type return_value;
4743 return_value[ii] = phi_grad[jj];
4745 return return_value;
4750 divergence_type return_value;
4751 return return_value;
4757 template <
int dim,
int spacedim>
4760 const unsigned int q_point)
const 4762 Assert(shape_function < fe_values->fe->dofs_per_cell,
4763 ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4766 "update_gradients")));
4769 shape_function_data[shape_function].single_nonzero_component;
4774 return gradient_type();
4788 const unsigned int comp =
4789 shape_function_data[shape_function].single_nonzero_component_index;
4792 const unsigned int ii = indices[0];
4793 const unsigned int jj = indices[1];
4795 const ::Tensor<1, spacedim> &phi_grad =
4796 fe_values->finite_element_output.shape_gradients[snc][q_point];
4798 gradient_type return_value;
4799 return_value[ii][jj] = phi_grad;
4801 return return_value;
4806 gradient_type return_value;
4807 return return_value;
4819 template <
int dim,
int spacedim>
4826 fe_values_views_cache.scalars.size()));
4828 return fe_values_views_cache.scalars[scalar.
component];
4833 template <
int dim,
int spacedim>
4840 fe_values_views_cache.vectors.size()));
4847 template <
int dim,
int spacedim>
4854 fe_values_views_cache.symmetric_second_order_tensors.size(),
4857 fe_values_views_cache.symmetric_second_order_tensors.size()));
4859 return fe_values_views_cache
4865 template <
int dim,
int spacedim>
4871 fe_values_views_cache.second_order_tensors.size(),
4874 fe_values_views_cache.second_order_tensors.size()));
4876 return fe_values_views_cache
4882 template <
int dim,
int spacedim>
4883 inline const double &
4885 const unsigned int j)
const 4890 Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
4891 Assert(present_cell.get() !=
nullptr,
4892 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
4895 if (fe->is_primitive())
4896 return this->finite_element_output.shape_values(i, j);
4907 const unsigned int row =
4908 this->finite_element_output
4909 .shape_function_to_row_table[i * fe->n_components() +
4910 fe->system_to_component_index(i).first];
4911 return this->finite_element_output.shape_values(row, j);
4917 template <
int dim,
int spacedim>
4920 const unsigned int i,
4921 const unsigned int j,
4922 const unsigned int component)
const 4929 Assert(present_cell.get() !=
nullptr,
4930 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
4935 if (fe->get_nonzero_components(i)[component] ==
false)
4941 const unsigned int row =
4942 this->finite_element_output
4943 .shape_function_to_row_table[i * fe->n_components() + component];
4944 return this->finite_element_output.shape_values(row, j);
4949 template <
int dim,
int spacedim>
4952 const unsigned int j)
const 4957 Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
4958 Assert(present_cell.get() !=
nullptr,
4959 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
4962 if (fe->is_primitive())
4963 return this->finite_element_output.shape_gradients[i][j];
4974 const unsigned int row =
4975 this->finite_element_output
4976 .shape_function_to_row_table[i * fe->n_components() +
4977 fe->system_to_component_index(i).first];
4978 return this->finite_element_output.shape_gradients[row][j];
4984 template <
int dim,
int spacedim>
4987 const unsigned int i,
4988 const unsigned int j,
4989 const unsigned int component)
const 4996 Assert(present_cell.get() !=
nullptr,
4997 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5001 if (fe->get_nonzero_components(i)[component] ==
false)
5007 const unsigned int row =
5008 this->finite_element_output
5009 .shape_function_to_row_table[i * fe->n_components() + component];
5010 return this->finite_element_output.shape_gradients[row][j];
5015 template <
int dim,
int spacedim>
5018 const unsigned int j)
const 5023 Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
5024 Assert(present_cell.get() !=
nullptr,
5025 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5028 if (fe->is_primitive())
5029 return this->finite_element_output.shape_hessians[i][j];
5040 const unsigned int row =
5041 this->finite_element_output
5042 .shape_function_to_row_table[i * fe->n_components() +
5043 fe->system_to_component_index(i).first];
5044 return this->finite_element_output.shape_hessians[row][j];
5050 template <
int dim,
int spacedim>
5053 const unsigned int i,
5054 const unsigned int j,
5055 const unsigned int component)
const 5062 Assert(present_cell.get() !=
nullptr,
5063 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5067 if (fe->get_nonzero_components(i)[component] ==
false)
5073 const unsigned int row =
5074 this->finite_element_output
5075 .shape_function_to_row_table[i * fe->n_components() + component];
5076 return this->finite_element_output.shape_hessians[row][j];
5081 template <
int dim,
int spacedim>
5084 const unsigned int j)
const 5089 Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
5090 Assert(present_cell.get() !=
nullptr,
5091 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5094 if (fe->is_primitive())
5095 return this->finite_element_output.shape_3rd_derivatives[i][j];
5106 const unsigned int row =
5107 this->finite_element_output
5108 .shape_function_to_row_table[i * fe->n_components() +
5109 fe->system_to_component_index(i).first];
5110 return this->finite_element_output.shape_3rd_derivatives[row][j];
5116 template <
int dim,
int spacedim>
5119 const unsigned int i,
5120 const unsigned int j,
5121 const unsigned int component)
const 5128 Assert(present_cell.get() !=
nullptr,
5129 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5133 if (fe->get_nonzero_components(i)[component] ==
false)
5139 const unsigned int row =
5140 this->finite_element_output
5141 .shape_function_to_row_table[i * fe->n_components() + component];
5142 return this->finite_element_output.shape_3rd_derivatives[row][j];
5147 template <
int dim,
int spacedim>
5156 template <
int dim,
int spacedim>
5165 template <
int dim,
int spacedim>
5169 return this->update_flags;
5174 template <
int dim,
int spacedim>
5175 inline const std::vector<Point<spacedim>> &
5180 Assert(present_cell.get() !=
nullptr,
5181 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5182 return this->mapping_output.quadrature_points;
5187 template <
int dim,
int spacedim>
5188 inline const std::vector<double> &
5193 Assert(present_cell.get() !=
nullptr,
5194 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5195 return this->mapping_output.JxW_values;
5200 template <
int dim,
int spacedim>
5201 inline const std::vector<DerivativeForm<1, dim, spacedim>> &
5206 Assert(present_cell.get() !=
nullptr,
5207 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5208 return this->mapping_output.jacobians;
5213 template <
int dim,
int spacedim>
5214 inline const std::vector<DerivativeForm<2, dim, spacedim>> &
5219 Assert(present_cell.get() !=
nullptr,
5220 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5221 return this->mapping_output.jacobian_grads;
5226 template <
int dim,
int spacedim>
5229 const unsigned int i)
const 5233 Assert(present_cell.get() !=
nullptr,
5234 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5235 return this->mapping_output.jacobian_pushed_forward_grads[i];
5240 template <
int dim,
int spacedim>
5241 inline const std::vector<Tensor<3, spacedim>> &
5246 Assert(present_cell.get() !=
nullptr,
5247 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5248 return this->mapping_output.jacobian_pushed_forward_grads;
5253 template <
int dim,
int spacedim>
5259 Assert(present_cell.get() !=
nullptr,
5260 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5261 return this->mapping_output.jacobian_2nd_derivatives[i];
5266 template <
int dim,
int spacedim>
5267 inline const std::vector<DerivativeForm<3, dim, spacedim>> &
5272 Assert(present_cell.get() !=
nullptr,
5273 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5274 return this->mapping_output.jacobian_2nd_derivatives;
5279 template <
int dim,
int spacedim>
5282 const unsigned int i)
const 5286 "update_jacobian_pushed_forward_2nd_derivatives"));
5287 Assert(present_cell.get() !=
nullptr,
5288 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5289 return this->mapping_output.jacobian_pushed_forward_2nd_derivatives[i];
5294 template <
int dim,
int spacedim>
5295 inline const std::vector<Tensor<4, spacedim>> &
5300 "update_jacobian_pushed_forward_2nd_derivatives"));
5301 Assert(present_cell.get() !=
nullptr,
5302 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5303 return this->mapping_output.jacobian_pushed_forward_2nd_derivatives;
5308 template <
int dim,
int spacedim>
5314 Assert(present_cell.get() !=
nullptr,
5315 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5316 return this->mapping_output.jacobian_3rd_derivatives[i];
5321 template <
int dim,
int spacedim>
5322 inline const std::vector<DerivativeForm<4, dim, spacedim>> &
5327 Assert(present_cell.get() !=
nullptr,
5328 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5329 return this->mapping_output.jacobian_3rd_derivatives;
5334 template <
int dim,
int spacedim>
5337 const unsigned int i)
const 5341 "update_jacobian_pushed_forward_3rd_derivatives"));
5342 Assert(present_cell.get() !=
nullptr,
5343 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5344 return this->mapping_output.jacobian_pushed_forward_3rd_derivatives[i];
5349 template <
int dim,
int spacedim>
5350 inline const std::vector<Tensor<5, spacedim>> &
5355 "update_jacobian_pushed_forward_3rd_derivatives"));
5356 Assert(present_cell.get() !=
nullptr,
5357 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5358 return this->mapping_output.jacobian_pushed_forward_3rd_derivatives;
5363 template <
int dim,
int spacedim>
5364 inline const std::vector<DerivativeForm<1, spacedim, dim>> &
5369 Assert(present_cell.get() !=
nullptr,
5370 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5371 return this->mapping_output.inverse_jacobians;
5376 template <
int dim,
int spacedim>
5382 Assert(i < this->mapping_output.quadrature_points.size(),
5383 ExcIndexRange(i, 0, this->mapping_output.quadrature_points.size()));
5384 Assert(present_cell.get() !=
nullptr,
5385 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5387 return this->mapping_output.quadrature_points[i];
5392 template <
int dim,
int spacedim>
5398 Assert(i < this->mapping_output.JxW_values.size(),
5399 ExcIndexRange(i, 0, this->mapping_output.JxW_values.size()));
5400 Assert(present_cell.get() !=
nullptr,
5401 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5403 return this->mapping_output.JxW_values[i];
5408 template <
int dim,
int spacedim>
5414 Assert(i < this->mapping_output.jacobians.size(),
5415 ExcIndexRange(i, 0, this->mapping_output.jacobians.size()));
5416 Assert(present_cell.get() !=
nullptr,
5417 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5419 return this->mapping_output.jacobians[i];
5424 template <
int dim,
int spacedim>
5430 Assert(i < this->mapping_output.jacobian_grads.size(),
5431 ExcIndexRange(i, 0, this->mapping_output.jacobian_grads.size()));
5432 Assert(present_cell.get() !=
nullptr,
5433 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5435 return this->mapping_output.jacobian_grads[i];
5440 template <
int dim,
int spacedim>
5446 Assert(i < this->mapping_output.inverse_jacobians.size(),
5447 ExcIndexRange(i, 0, this->mapping_output.inverse_jacobians.size()));
5448 Assert(present_cell.get() !=
nullptr,
5449 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5451 return this->mapping_output.inverse_jacobians[i];
5456 template <
int dim,
int spacedim>
5462 "update_normal_vectors")));
5463 Assert(i < this->mapping_output.normal_vectors.size(),
5464 ExcIndexRange(i, 0, this->mapping_output.normal_vectors.size()));
5465 Assert(present_cell.get() !=
nullptr,
5466 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
5468 return this->mapping_output.normal_vectors[i];
5476 template <
int dim,
int spacedim>
5485 template <
int dim,
int spacedim>
5496 template <
int dim,
int spacedim>
5500 return present_face_index;
5506 template <
int dim,
int spacedim>
5515 template <
int dim,
int spacedim>
5524 template <
int dim,
int spacedim>
5533 template <
int dim,
int spacedim>
5537 Assert(i < this->mapping_output.boundary_forms.size(),
5538 ExcIndexRange(i, 0, this->mapping_output.boundary_forms.size()));
5541 "update_boundary_forms")));
5543 return this->mapping_output.boundary_forms[i];
5548 DEAL_II_NAMESPACE_CLOSE
Vector & operator=(const Vector< dim, spacedim > &)=delete
void get_function_third_derivatives(const InputVector &fe_function, std::vector< typename ProductType< third_derivative_type, typename InputVector::value_type >::type > &third_derivatives) const
Transformed quadrature weights.
constexpr Tensor()=default
virtual ~FEValuesBase() override
value_type value(const unsigned int shape_function, const unsigned int q_point) const
typename ProductType< Number, typename Vector< dim, spacedim >::curl_type >::type curl_type
::SymmetricTensor< 2, spacedim > symmetric_gradient_type
const FEFaceValues< dim, spacedim > & get_present_fe_values() const
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access >> &cell, const unsigned int face_no, const unsigned int subface_no)
const FEValuesViews::Scalar< dim, spacedim > & operator[](const FEValuesExtractors::Scalar &scalar) const
const Tensor< 3, spacedim > & jacobian_pushed_forward_grad(const unsigned int quadrature_point) const
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > fe_data
CellSimilarity::Similarity cell_similarity
typename ProductType< Number, typename Tensor< 2, dim, spacedim >::gradient_type >::type gradient_type
::Tensor< 1, spacedim > divergence_type
static const unsigned int dimension
Cache(const FEValuesBase< dim, spacedim > &fe_values)
typename ProductType< Number, typename SymmetricTensor< 2, dim, spacedim >::value_type >::type value_type
unsigned int present_face_index
std::vector< ShapeFunctionData > shape_function_data
::Tensor< 2, spacedim > gradient_type
static ::ExceptionBase & ExcAccessToUninitializedField()
value_type value(const unsigned int shape_function, const unsigned int q_point) const
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
::Tensor< 2, spacedim > hessian_type
static ::ExceptionBase & ExcFaceHasNoSubfaces()
::Tensor< 1, spacedim > divergence_type
const unsigned int dofs_per_cell
typename ::internal::CurlType< spacedim >::type curl_type
const unsigned int component
void get_function_symmetric_gradients_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::symmetric_gradient_type > &symmetric_gradients) const
void get_function_values(const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &values) const
gradient_type gradient(const unsigned int shape_function, const unsigned int q_point) const
const Quadrature< dim - 1 > quadrature
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)
void get_function_third_derivatives(const InputVector &fe_function, std::vector< typename ProductType< third_derivative_type, typename InputVector::value_type >::type > &third_derivatives) const
static const unsigned int dimension
static ::ExceptionBase & ExcAccessToUninitializedField(std::string arg1)
int single_nonzero_component
const Mapping< dim, spacedim > & get_mapping() const
void get_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
static ::ExceptionBase & ExcReinitCalledWithBoundaryFace()
Outer normal vector, not normalized.
typename ProductType< Number, typename Scalar< dim, spacedim >::hessian_type >::type hessian_type
const DerivativeForm< 2, dim, spacedim > & jacobian_grad(const unsigned int quadrature_point) const
FEFaceValuesBase(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, const Quadrature< dim - 1 > &quadrature)
const FiniteElement< dim, spacedim > & get_fe() const
std::unique_ptr< const CellIteratorBase > present_cell
::Tensor< 4, spacedim > third_derivative_type
::Tensor< 3, spacedim > gradient_type
static ::ExceptionBase & ExcFEDontMatch()
void get_function_divergences(const InputVector &fe_function, std::vector< typename ProductType< divergence_type, typename InputVector::value_type >::type > &divergences) const
const unsigned int first_tensor_component
void get_function_laplacians(const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &laplacians) const
UpdateFlags get_update_flags() const
third_derivative_type third_derivative(const unsigned int shape_function, const unsigned int q_point) const
Transformed quadrature points.
void do_reinit(const unsigned int face_no)
std::size_t memory_consumption() const
typename ProductType< Number, typename Scalar< dim, spacedim >::value_type >::type value_type
unsigned int row_index[spacedim]
void check_cell_similarity(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
static const unsigned int space_dimension
const SmartPointer< const Mapping< dim, spacedim >, FEValuesBase< dim, spacedim > > mapping
const DerivativeForm< 1, dim, spacedim > & jacobian(const unsigned int quadrature_point) const
const DerivativeForm< 3, dim, spacedim > & jacobian_2nd_derivative(const unsigned int quadrature_point) const
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcShapeFunctionNotPrimitive(int arg1)
::internal::FEValuesViews::Cache< dim, spacedim > fe_values_views_cache
void do_reinit(const unsigned int face_no, const unsigned int subface_no)
static const unsigned int integral_dimension
LinearAlgebra::distributed::Vector< Number > Vector
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
void get_function_divergences_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::divergence_type > &divergences) const
void get_function_curls_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::curl_type > &curls) const
Tensor< 2, spacedim > shape_hessian_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
typename ProductType< Number, typename Vector< dim, spacedim >::gradient_type >::type gradient_type
static const unsigned int integral_dimension
const Quadrature< dim - 1 > & get_quadrature() const
void get_function_laplacians(const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &laplacians) const
void get_function_hessians_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::hessian_type > &hessians) const
const std::vector< Point< spacedim > > & get_quadrature_points() const
const Point< spacedim > & quadrature_point(const unsigned int q) const
const std::vector< Tensor< 4, spacedim > > & get_jacobian_pushed_forward_2nd_derivatives() const
void get_function_gradients(const InputVector &fe_function, std::vector< typename ProductType< gradient_type, typename InputVector::value_type >::type > &gradients) const
typename ProductType< Number, typename Scalar< dim, spacedim >::gradient_type >::type gradient_type
typename internal::FEValuesViews::ViewType< dim, spacedim, Extractor >::type View
void get_function_hessians(const InputVector &fe_function, std::vector< typename ProductType< hessian_type, typename InputVector::value_type >::type > &hessians) const
void get_function_laplacians(const InputVector &fe_function, std::vector< typename InputVector::value_type > &laplacians) const
typename ProductType< Number, typename Tensor< 2, dim, spacedim >::value_type >::type value_type
static ::ExceptionBase & ExcMessage(std::string arg1)
::Tensor< 1, spacedim > gradient_type
const unsigned int first_tensor_component
void get_function_gradients_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::gradient_type > &gradients) const
typename ProductType< Number, typename Vector< dim, spacedim >::third_derivative_type >::type third_derivative_type
#define DeclException1(Exception1, type1, outsequence)
const double & shape_value(const unsigned int function_no, const unsigned int point_no) const
CellSimilarity::Similarity get_cell_similarity() const
Tensor< 3, spacedim > shape_3rd_derivative_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
Third derivatives of shape functions.
const FEValues< dim, spacedim > & get_present_fe_values() const
UpdateFlags compute_update_flags(const UpdateFlags update_flags) const
std::vector< ShapeFunctionData > shape_function_data
void get_function_laplacians_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::laplacian_type > &laplacians) const
unsigned int get_face_index() const
#define Assert(cond, exc)
const Tensor< 2, spacedim > & shape_hessian(const unsigned int function_no, const unsigned int point_no) const
FESubfaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &face_quadrature, const UpdateFlags update_flags)
Abstract base class for mapping classes.
Scalar & operator=(const Scalar< dim, spacedim > &)=delete
std::vector< ShapeFunctionData > shape_function_data
#define DeclExceptionMsg(Exception, defaulttext)
const Quadrature< dim > quadrature
const unsigned int first_vector_component
const std::vector< DerivativeForm< 1, dim, spacedim > > & get_jacobians() const
#define DeclException0(Exception0)
void invalidate_present_cell()
unsigned int single_nonzero_component_index
std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > mapping_data
const DerivativeForm< 4, dim, spacedim > & jacobian_3rd_derivative(const unsigned int quadrature_point) const
typename Tensor< rank_ - 1, dim, Number >::tensor_type value_type
static const unsigned int integral_dimension
const std::vector< DerivativeForm< 4, dim, spacedim > > & get_jacobian_3rd_derivatives() const
SymmetricTensor< 2, dim, Number > symmetrize(const Tensor< 2, dim, Number > &t)
symmetric_gradient_type symmetric_gradient(const unsigned int shape_function, const unsigned int q_point) const
void get_function_hessians_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::hessian_type > &hessians) const
bool is_nonzero_shape_function_component[spacedim]
void get_function_symmetric_gradients(const InputVector &fe_function, std::vector< typename ProductType< symmetric_gradient_type, typename InputVector::value_type >::type > &symmetric_gradients) const
static const unsigned int integral_dimension
Second derivatives of shape functions.
static ::ExceptionBase & ExcFENotPrimitive()
Gradient of volume element.
const std::vector< Tensor< 1, spacedim > > & get_all_normal_vectors() const
typename ProductType< Number, typename SymmetricTensor< 2, dim, spacedim >::divergence_type >::type divergence_type
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
static TableIndices< rank_ > unrolled_to_component_indices(const unsigned int i)
divergence_type divergence(const unsigned int shape_function, const unsigned int q_point) const
std::vector<::FEValuesViews::Scalar< dim, spacedim > > scalars
std::size_t memory_consumption() const
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access >> &cell)
::SymmetricTensor< 2, spacedim > value_type
void get_function_hessians(const InputVector &fe_function, std::vector< typename ProductType< hessian_type, typename InputVector::value_type >::type > &hessians) const
const Tensor< 1, spacedim > & boundary_form(const unsigned int i) const
const Quadrature< dim > & get_quadrature() const
const unsigned int n_quadrature_points
void get_function_hessians(const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type >> &hessians) const
Tensor< 1, spacedim > shape_grad_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
FEFaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature, const UpdateFlags update_flags)
void get_function_curls(const InputVector &fe_function, std::vector< typename ProductType< curl_type, typename InputVector::value_type >::type > &curls) const
void get_function_third_derivatives(const InputVector &fe_function, std::vector< Tensor< 3, spacedim, typename InputVector::value_type >> &third_derivatives) const
static const unsigned int dimension
boost::signals2::connection tria_listener_mesh_transform
typename ProductType< Number, typename Vector< dim, spacedim >::value_type >::type value_type
static const unsigned int space_dimension
int single_nonzero_component
void get_function_laplacians_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::laplacian_type > &laplacians) const
const std::vector< Tensor< 3, spacedim > > & get_jacobian_pushed_forward_grads() const
const std::vector< double > & get_JxW_values() const
unsigned int single_nonzero_component_index
double JxW(const unsigned int quadrature_point) const
typename ProductType< Number, typename Vector< dim, spacedim >::symmetric_gradient_type >::type symmetric_gradient_type
typename ProductType< Number, typename Vector< dim, spacedim >::divergence_type >::type divergence_type
void initialize(const UpdateFlags update_flags)
typename ProductType< Number, typename Vector< dim, spacedim >::hessian_type >::type hessian_type
Shape function gradients.
void get_function_third_derivatives_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::third_derivative_type > &third_derivatives) const
const std::vector< DerivativeForm< 2, dim, spacedim > > & get_jacobian_grads() const
void maybe_invalidate_previous_present_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
typename ProductType< Number, typename Tensor< 2, dim, spacedim >::divergence_type >::type divergence_type
hessian_type hessian(const unsigned int shape_function, const unsigned int q_point) const
const Tensor< 4, spacedim > & jacobian_pushed_forward_2nd_derivative(const unsigned int quadrature_point) const
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
void get_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type >> &gradients) const
void initialize(const UpdateFlags update_flags)
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
void get_function_gradients_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::gradient_type > &gradients) const
const DerivativeForm< 1, spacedim, dim > & inverse_jacobian(const unsigned int quadrature_point) const
curl_type curl(const unsigned int shape_function, const unsigned int q_point) const
static ::ExceptionBase & ExcNotImplemented()
const std::vector< DerivativeForm< 3, dim, spacedim > > & get_jacobian_2nd_derivatives() const
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access >> &cell, const unsigned int face_no)
bool is_nonzero_shape_function_component
const Triangulation< dim, spacedim >::cell_iterator get_cell() const
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
void get_function_values(const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &values) const
const Tensor< 1, spacedim > & shape_grad(const unsigned int function_no, const unsigned int quadrature_point) const
boost::signals2::connection tria_listener_refinement
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > finite_element_output
gradient_type gradient(const unsigned int shape_function, const unsigned int q_point) const
typename ProductType< Number, typename Vector< dim, spacedim >::value_type >::type laplacian_type
const std::vector< Tensor< 1, spacedim > > & get_boundary_forms() const
::Tensor< 3, spacedim > hessian_type
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
void initialize(const UpdateFlags update_flags)
third_derivative_type third_derivative(const unsigned int shape_function, const unsigned int q_point) const
::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > mapping_output
FEValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags update_flags)
void get_function_gradients(const InputVector &fe_function, std::vector< typename ProductType< gradient_type, typename InputVector::value_type >::type > &gradients) const
hessian_type hessian(const unsigned int shape_function, const unsigned int q_point) const
void get_function_third_derivatives_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::third_derivative_type > &third_derivatives) const
const std::vector< Tensor< 5, spacedim > > & get_jacobian_pushed_forward_3rd_derivatives() const
const FESubfaceValues< dim, spacedim > & get_present_fe_values() const
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors() const
const std::vector< DerivativeForm< 1, spacedim, dim > > & get_inverse_jacobians() const
FEValuesBase & operator=(const FEValuesBase &)=delete
void get_function_values_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::value_type > &values) const
int single_nonzero_component
typename ProductType< Number, typename Scalar< dim, spacedim >::value_type >::type laplacian_type
::Tensor< 3, spacedim > third_derivative_type
typename ProductType< Number, typename Scalar< dim, spacedim >::third_derivative_type >::type third_derivative_type
std::size_t memory_consumption() const
void get_function_values_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::value_type > &values) const
std::vector< ShapeFunctionData > shape_function_data
const Tensor< 3, spacedim > & shape_3rd_derivative(const unsigned int function_no, const unsigned int point_no) const
static ::ExceptionBase & ExcInternalError()
const SmartPointer< const FiniteElement< dim, spacedim >, FEValuesBase< dim, spacedim > > fe
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const
const Tensor< 5, spacedim > & jacobian_pushed_forward_3rd_derivative(const unsigned int quadrature_point) const