Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
fe_values.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2019 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_fe_values_h
17 #define dealii_fe_values_h
18 
19 
20 #include <deal.II/base/config.h>
21 
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>
28 
29 #include <deal.II/dofs/dof_accessor.h>
30 #include <deal.II/dofs/dof_handler.h>
31 
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>
36 
37 #include <deal.II/grid/tria.h>
38 #include <deal.II/grid/tria_iterator.h>
39 
40 #include <deal.II/hp/dof_handler.h>
41 
42 #include <algorithm>
43 #include <memory>
44 #include <type_traits>
45 
46 
47 // dummy include in order to have the
48 // definition of PetscScalar available
49 // without including other PETSc stuff
50 #ifdef DEAL_II_WITH_PETSC
51 # include <petsc.h>
52 #endif
53 
54 DEAL_II_NAMESPACE_OPEN
55 
56 template <int dim, int spacedim = dim>
57 class FEValuesBase;
58 
59 namespace internal
60 {
65  template <int dim, class NumberType = double>
66  struct CurlType;
67 
74  template <class NumberType>
75  struct CurlType<1, NumberType>
76  {
78  };
79 
86  template <class NumberType>
87  struct CurlType<2, NumberType>
88  {
90  };
91 
98  template <class NumberType>
99  struct CurlType<3, NumberType>
100  {
102  };
103 } // namespace internal
104 
105 
106 
128 namespace FEValuesViews
129 {
141  template <int dim, int spacedim = dim>
142  class Scalar
143  {
144  public:
150  using value_type = double;
151 
158 
165 
172 
177  template <typename Number>
178  struct OutputType
179  {
184  using value_type =
185  typename ProductType<Number,
187 
192  using gradient_type = typename ProductType<
193  Number,
195 
200  using laplacian_type =
201  typename ProductType<Number,
203 
208  using hessian_type = typename ProductType<
209  Number,
211 
216  using third_derivative_type = typename ProductType<
217  Number,
219  };
220 
226  {
236 
245  unsigned int row_index;
246  };
247 
251  Scalar();
252 
258  Scalar(const FEValuesBase<dim, spacedim> &fe_values_base,
259  const unsigned int component);
260 
265  Scalar &
266  operator=(const Scalar<dim, spacedim> &) = delete;
267 
281  value_type
282  value(const unsigned int shape_function, const unsigned int q_point) const;
283 
295  gradient(const unsigned int shape_function,
296  const unsigned int q_point) const;
297 
309  hessian(const unsigned int shape_function,
310  const unsigned int q_point) const;
311 
323  third_derivative(const unsigned int shape_function,
324  const unsigned int q_point) const;
325 
343  template <class InputVector>
344  void
346  const InputVector &fe_function,
347  std::vector<typename ProductType<value_type,
348  typename InputVector::value_type>::type>
349  &values) const;
350 
373  template <class InputVector>
374  void
376  const InputVector &dof_values,
377  std::vector<
379  &values) const;
380 
398  template <class InputVector>
399  void
401  const InputVector &fe_function,
402  std::vector<typename ProductType<gradient_type,
403  typename InputVector::value_type>::type>
404  &gradients) const;
405 
409  template <class InputVector>
410  void
412  const InputVector &dof_values,
413  std::vector<
415  &gradients) const;
416 
434  template <class InputVector>
435  void
437  const InputVector &fe_function,
438  std::vector<typename ProductType<hessian_type,
439  typename InputVector::value_type>::type>
440  &hessians) const;
441 
445  template <class InputVector>
446  void
448  const InputVector &dof_values,
449  std::vector<
451  &hessians) const;
452 
453 
472  template <class InputVector>
473  void
475  const InputVector &fe_function,
476  std::vector<typename ProductType<value_type,
477  typename InputVector::value_type>::type>
478  &laplacians) const;
479 
483  template <class InputVector>
484  void
486  const InputVector &dof_values,
487  std::vector<
489  &laplacians) const;
490 
491 
510  template <class InputVector>
511  void
513  const InputVector &fe_function,
514  std::vector<typename ProductType<third_derivative_type,
515  typename InputVector::value_type>::type>
516  &third_derivatives) const;
517 
521  template <class InputVector>
522  void
524  const InputVector & dof_values,
525  std::vector<typename OutputType<typename InputVector::value_type>::
526  third_derivative_type> &third_derivatives) const;
527 
528 
529  private:
534 
539  const unsigned int component;
540 
544  std::vector<ShapeFunctionData> shape_function_data;
545  };
546 
547 
548 
578  template <int dim, int spacedim = dim>
579  class Vector
580  {
581  public:
588 
598 
610 
616  using divergence_type = double;
617 
624  using curl_type = typename ::internal::CurlType<spacedim>::type;
625 
632 
639 
644  template <typename Number>
645  struct OutputType
646  {
651  using value_type =
652  typename ProductType<Number,
654 
659  using gradient_type = typename ProductType<
660  Number,
662 
667  using symmetric_gradient_type = typename ProductType<
668  Number,
670 
675  using divergence_type = typename ProductType<
676  Number,
678 
683  using laplacian_type =
684  typename ProductType<Number,
686 
691  using curl_type =
692  typename ProductType<Number,
694 
699  using hessian_type = typename ProductType<
700  Number,
702 
707  using third_derivative_type = typename ProductType<
708  Number,
710  };
711 
717  {
727 
737  unsigned int row_index[spacedim];
738 
748  unsigned int single_nonzero_component_index;
749  };
750 
754  Vector();
755 
764  Vector(const FEValuesBase<dim, spacedim> &fe_values_base,
765  const unsigned int first_vector_component);
766 
771  Vector &
772  operator=(const Vector<dim, spacedim> &) = delete;
773 
790  value_type
791  value(const unsigned int shape_function, const unsigned int q_point) const;
792 
807  gradient(const unsigned int shape_function,
808  const unsigned int q_point) const;
809 
826  symmetric_gradient(const unsigned int shape_function,
827  const unsigned int q_point) const;
828 
840  divergence(const unsigned int shape_function,
841  const unsigned int q_point) const;
842 
863  curl_type
864  curl(const unsigned int shape_function, const unsigned int q_point) const;
865 
877  hessian(const unsigned int shape_function,
878  const unsigned int q_point) const;
879 
891  third_derivative(const unsigned int shape_function,
892  const unsigned int q_point) const;
893 
911  template <class InputVector>
912  void
914  const InputVector &fe_function,
915  std::vector<typename ProductType<value_type,
916  typename InputVector::value_type>::type>
917  &values) const;
918 
941  template <class InputVector>
942  void
944  const InputVector &dof_values,
945  std::vector<
947  &values) const;
948 
966  template <class InputVector>
967  void
969  const InputVector &fe_function,
970  std::vector<typename ProductType<gradient_type,
971  typename InputVector::value_type>::type>
972  &gradients) const;
973 
977  template <class InputVector>
978  void
980  const InputVector &dof_values,
981  std::vector<
983  &gradients) const;
984 
1008  template <class InputVector>
1009  void
1011  const InputVector &fe_function,
1012  std::vector<typename ProductType<symmetric_gradient_type,
1013  typename InputVector::value_type>::type>
1014  &symmetric_gradients) const;
1015 
1019  template <class InputVector>
1020  void
1022  const InputVector & dof_values,
1023  std::vector<typename OutputType<typename InputVector::value_type>::
1024  symmetric_gradient_type> &symmetric_gradients) const;
1025 
1044  template <class InputVector>
1045  void
1047  const InputVector &fe_function,
1048  std::vector<typename ProductType<divergence_type,
1049  typename InputVector::value_type>::type>
1050  &divergences) const;
1051 
1055  template <class InputVector>
1056  void
1058  const InputVector &dof_values,
1059  std::vector<
1061  &divergences) const;
1062 
1081  template <class InputVector>
1082  void
1084  const InputVector &fe_function,
1085  std::vector<
1086  typename ProductType<curl_type, typename InputVector::value_type>::type>
1087  &curls) const;
1088 
1092  template <class InputVector>
1093  void
1095  const InputVector &dof_values,
1096  std::vector<
1098  &curls) const;
1099 
1117  template <class InputVector>
1118  void
1120  const InputVector &fe_function,
1121  std::vector<typename ProductType<hessian_type,
1122  typename InputVector::value_type>::type>
1123  &hessians) const;
1124 
1128  template <class InputVector>
1129  void
1131  const InputVector &dof_values,
1132  std::vector<
1134  &hessians) const;
1135 
1154  template <class InputVector>
1155  void
1157  const InputVector &fe_function,
1158  std::vector<typename ProductType<value_type,
1159  typename InputVector::value_type>::type>
1160  &laplacians) const;
1161 
1165  template <class InputVector>
1166  void
1168  const InputVector &dof_values,
1169  std::vector<
1171  &laplacians) const;
1172 
1191  template <class InputVector>
1192  void
1194  const InputVector &fe_function,
1195  std::vector<typename ProductType<third_derivative_type,
1196  typename InputVector::value_type>::type>
1197  &third_derivatives) const;
1198 
1202  template <class InputVector>
1203  void
1205  const InputVector & dof_values,
1206  std::vector<typename OutputType<typename InputVector::value_type>::
1207  third_derivative_type> &third_derivatives) const;
1208 
1209  private:
1214 
1219  const unsigned int first_vector_component;
1220 
1224  std::vector<ShapeFunctionData> shape_function_data;
1225  };
1226 
1227 
1228  template <int rank, int dim, int spacedim = dim>
1229  class SymmetricTensor;
1230 
1255  template <int dim, int spacedim>
1256  class SymmetricTensor<2, dim, spacedim>
1257  {
1258  public:
1266 
1277 
1282  template <typename Number>
1283  struct OutputType
1284  {
1289  using value_type = typename ProductType<
1290  Number,
1292 
1297  using divergence_type = typename ProductType<
1298  Number,
1300  };
1301 
1306  struct ShapeFunctionData
1307  {
1316  bool is_nonzero_shape_function_component
1317  [value_type::n_independent_components];
1318 
1328  unsigned int row_index[value_type::n_independent_components];
1329 
1339 
1344  };
1345 
1349  SymmetricTensor();
1350 
1360  SymmetricTensor(const FEValuesBase<dim, spacedim> &fe_values_base,
1361  const unsigned int first_tensor_component);
1362 
1367  SymmetricTensor &
1368  operator=(const SymmetricTensor<2, dim, spacedim> &) = delete;
1369 
1387  value_type
1388  value(const unsigned int shape_function, const unsigned int q_point) const;
1389 
1404  divergence(const unsigned int shape_function,
1405  const unsigned int q_point) const;
1406 
1424  template <class InputVector>
1425  void
1426  get_function_values(
1427  const InputVector &fe_function,
1428  std::vector<typename ProductType<value_type,
1429  typename InputVector::value_type>::type>
1430  &values) const;
1431 
1454  template <class InputVector>
1455  void
1456  get_function_values_from_local_dof_values(
1457  const InputVector &dof_values,
1458  std::vector<
1460  &values) const;
1461 
1483  template <class InputVector>
1484  void
1485  get_function_divergences(
1486  const InputVector &fe_function,
1487  std::vector<typename ProductType<divergence_type,
1488  typename InputVector::value_type>::type>
1489  &divergences) const;
1490 
1494  template <class InputVector>
1495  void
1496  get_function_divergences_from_local_dof_values(
1497  const InputVector &dof_values,
1498  std::vector<
1500  &divergences) const;
1501 
1502  private:
1507 
1512  const unsigned int first_tensor_component;
1513 
1517  std::vector<ShapeFunctionData> shape_function_data;
1518  };
1519 
1520 
1521  template <int rank, int dim, int spacedim = dim>
1522  class Tensor;
1523 
1544  template <int dim, int spacedim>
1545  class Tensor<2, dim, spacedim>
1546  {
1547  public:
1553 
1558 
1564 
1569  template <typename Number>
1570  struct OutputType
1571  {
1576  using value_type = typename ProductType<
1577  Number,
1579 
1584  using divergence_type = typename ProductType<
1585  Number,
1587 
1592  using gradient_type = typename ProductType<
1593  Number,
1595  };
1596 
1601  struct ShapeFunctionData
1602  {
1611  bool is_nonzero_shape_function_component
1612  [value_type::n_independent_components];
1613 
1623  unsigned int row_index[value_type::n_independent_components];
1624 
1634 
1639  };
1640 
1644  Tensor();
1645 
1655  Tensor(const FEValuesBase<dim, spacedim> &fe_values_base,
1656  const unsigned int first_tensor_component);
1657 
1658 
1663  Tensor &
1664  operator=(const Tensor<2, dim, spacedim> &) = delete;
1665 
1682  value_type
1683  value(const unsigned int shape_function, const unsigned int q_point) const;
1684 
1699  divergence(const unsigned int shape_function,
1700  const unsigned int q_point) const;
1701 
1716  gradient(const unsigned int shape_function,
1717  const unsigned int q_point) const;
1718 
1736  template <class InputVector>
1737  void
1738  get_function_values(
1739  const InputVector &fe_function,
1740  std::vector<typename ProductType<value_type,
1741  typename InputVector::value_type>::type>
1742  &values) const;
1743 
1766  template <class InputVector>
1767  void
1768  get_function_values_from_local_dof_values(
1769  const InputVector &dof_values,
1770  std::vector<
1772  &values) const;
1773 
1795  template <class InputVector>
1796  void
1797  get_function_divergences(
1798  const InputVector &fe_function,
1799  std::vector<typename ProductType<divergence_type,
1800  typename InputVector::value_type>::type>
1801  &divergences) const;
1802 
1806  template <class InputVector>
1807  void
1808  get_function_divergences_from_local_dof_values(
1809  const InputVector &dof_values,
1810  std::vector<
1812  &divergences) const;
1813 
1830  template <class InputVector>
1831  void
1832  get_function_gradients(
1833  const InputVector &fe_function,
1834  std::vector<typename ProductType<gradient_type,
1835  typename InputVector::value_type>::type>
1836  &gradients) const;
1837 
1841  template <class InputVector>
1842  void
1843  get_function_gradients_from_local_dof_values(
1844  const InputVector &dof_values,
1845  std::vector<
1847  &gradients) const;
1848 
1849  private:
1854 
1859  const unsigned int first_tensor_component;
1860 
1864  std::vector<ShapeFunctionData> shape_function_data;
1865  };
1866 
1867 } // namespace FEValuesViews
1868 
1869 
1870 namespace internal
1871 {
1872  namespace FEValuesViews
1873  {
1880  template <int dim, int spacedim, typename Extractor>
1881  struct ViewType
1882  {};
1883 
1891  template <int dim, int spacedim>
1892  struct ViewType<dim, spacedim, FEValuesExtractors::Scalar>
1893  {
1894  using type = typename ::FEValuesViews::Scalar<dim, spacedim>;
1895  };
1896 
1904  template <int dim, int spacedim>
1905  struct ViewType<dim, spacedim, FEValuesExtractors::Vector>
1906  {
1907  using type = typename ::FEValuesViews::Vector<dim, spacedim>;
1908  };
1909 
1917  template <int dim, int spacedim, int rank>
1918  struct ViewType<dim, spacedim, FEValuesExtractors::Tensor<rank>>
1919  {
1920  using type = typename ::FEValuesViews::Tensor<rank, dim, spacedim>;
1921  };
1922 
1930  template <int dim, int spacedim, int rank>
1931  struct ViewType<dim, spacedim, FEValuesExtractors::SymmetricTensor<rank>>
1932  {
1933  using type =
1934  typename ::FEValuesViews::SymmetricTensor<rank, dim, spacedim>;
1935  };
1936 
1944  template <int dim, int spacedim>
1945  struct Cache
1946  {
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;
1957 
1961  Cache(const FEValuesBase<dim, spacedim> &fe_values);
1962  };
1963  } // namespace FEValuesViews
1964 } // namespace internal
1965 
1966 namespace FEValuesViews
1967 {
1974  template <int dim, int spacedim, typename Extractor>
1975  using View =
1977 } // namespace FEValuesViews
1978 
1979 
2082 template <int dim, int spacedim>
2083 class FEValuesBase : public Subscriptor
2084 {
2085 public:
2089  static const unsigned int dimension = dim;
2090 
2094  static const unsigned int space_dimension = spacedim;
2095 
2099  const unsigned int n_quadrature_points;
2100 
2106  const unsigned int dofs_per_cell;
2107 
2108 
2116  FEValuesBase(const unsigned int n_q_points,
2117  const unsigned int dofs_per_cell,
2118  const UpdateFlags update_flags,
2121 
2126  FEValuesBase &
2127  operator=(const FEValuesBase &) = delete;
2128 
2133  FEValuesBase(const FEValuesBase &) = delete;
2134 
2138  virtual ~FEValuesBase() override;
2139 
2140 
2144 
2145 
2166  const double &
2167  shape_value(const unsigned int function_no,
2168  const unsigned int point_no) const;
2169 
2190  double
2191  shape_value_component(const unsigned int function_no,
2192  const unsigned int point_no,
2193  const unsigned int component) const;
2194 
2220  const Tensor<1, spacedim> &
2221  shape_grad(const unsigned int function_no,
2222  const unsigned int quadrature_point) const;
2223 
2241  shape_grad_component(const unsigned int function_no,
2242  const unsigned int point_no,
2243  const unsigned int component) const;
2244 
2264  const Tensor<2, spacedim> &
2265  shape_hessian(const unsigned int function_no,
2266  const unsigned int point_no) const;
2267 
2285  shape_hessian_component(const unsigned int function_no,
2286  const unsigned int point_no,
2287  const unsigned int component) const;
2288 
2308  const Tensor<3, spacedim> &
2309  shape_3rd_derivative(const unsigned int function_no,
2310  const unsigned int point_no) const;
2311 
2329  shape_3rd_derivative_component(const unsigned int function_no,
2330  const unsigned int point_no,
2331  const unsigned int component) const;
2332 
2334 
2336 
2373  template <class InputVector>
2374  void
2376  const InputVector & fe_function,
2377  std::vector<typename InputVector::value_type> &values) const;
2378 
2392  template <class InputVector>
2393  void
2395  const InputVector & fe_function,
2396  std::vector<Vector<typename InputVector::value_type>> &values) const;
2397 
2416  template <class InputVector>
2417  void
2419  const InputVector & fe_function,
2421  std::vector<typename InputVector::value_type> & values) const;
2422 
2444  template <class InputVector>
2445  void
2447  const InputVector & fe_function,
2449  std::vector<Vector<typename InputVector::value_type>> &values) const;
2450 
2451 
2482  template <class InputVector>
2483  void
2485  const InputVector & fe_function,
2487  ArrayView<std::vector<typename InputVector::value_type>> values,
2488  const bool quadrature_points_fastest) const;
2489 
2491 
2493 
2530  template <class InputVector>
2531  void
2533  const InputVector &fe_function,
2535  &gradients) const;
2536 
2553  template <class InputVector>
2554  void
2556  const InputVector &fe_function,
2557  std::vector<
2559  &gradients) const;
2560 
2567  template <class InputVector>
2568  void
2570  const InputVector & fe_function,
2573  &gradients) const;
2574 
2581  template <class InputVector>
2582  void
2584  const InputVector & fe_function,
2586  ArrayView<
2588  gradients,
2589  bool quadrature_points_fastest = false) const;
2590 
2592 
2596 
2634  template <class InputVector>
2635  void
2637  const InputVector &fe_function,
2639  &hessians) const;
2640 
2658  template <class InputVector>
2659  void
2661  const InputVector &fe_function,
2662  std::vector<
2664  & hessians,
2665  bool quadrature_points_fastest = false) const;
2666 
2671  template <class InputVector>
2672  void
2674  const InputVector & fe_function,
2677  &hessians) const;
2678 
2685  template <class InputVector>
2686  void
2688  const InputVector & fe_function,
2690  ArrayView<
2692  hessians,
2693  bool quadrature_points_fastest = false) const;
2694 
2735  template <class InputVector>
2736  void
2738  const InputVector & fe_function,
2739  std::vector<typename InputVector::value_type> &laplacians) const;
2740 
2760  template <class InputVector>
2761  void
2763  const InputVector & fe_function,
2764  std::vector<Vector<typename InputVector::value_type>> &laplacians) const;
2765 
2772  template <class InputVector>
2773  void
2775  const InputVector & fe_function,
2777  std::vector<typename InputVector::value_type> & laplacians) const;
2778 
2785  template <class InputVector>
2786  void
2788  const InputVector & fe_function,
2790  std::vector<Vector<typename InputVector::value_type>> &laplacians) const;
2791 
2798  template <class InputVector>
2799  void
2801  const InputVector & fe_function,
2803  std::vector<std::vector<typename InputVector::value_type>> &laplacians,
2804  bool quadrature_points_fastest = false) const;
2805 
2807 
2809 
2848  template <class InputVector>
2849  void
2851  const InputVector &fe_function,
2853  &third_derivatives) const;
2854 
2873  template <class InputVector>
2874  void
2876  const InputVector &fe_function,
2877  std::vector<
2879  & third_derivatives,
2880  bool quadrature_points_fastest = false) const;
2881 
2886  template <class InputVector>
2887  void
2889  const InputVector & fe_function,
2892  &third_derivatives) const;
2893 
2900  template <class InputVector>
2901  void
2903  const InputVector & fe_function,
2905  ArrayView<
2907  third_derivatives,
2908  bool quadrature_points_fastest = false) const;
2910 
2912 
2913 
2919  const Point<spacedim> &
2920  quadrature_point(const unsigned int q) const;
2921 
2927  const std::vector<Point<spacedim>> &
2928  get_quadrature_points() const;
2929 
2945  double
2946  JxW(const unsigned int quadrature_point) const;
2947 
2951  const std::vector<double> &
2952  get_JxW_values() const;
2953 
2961  jacobian(const unsigned int quadrature_point) const;
2962 
2969  const std::vector<DerivativeForm<1, dim, spacedim>> &
2970  get_jacobians() const;
2971 
2980  jacobian_grad(const unsigned int quadrature_point) const;
2981 
2988  const std::vector<DerivativeForm<2, dim, spacedim>> &
2989  get_jacobian_grads() const;
2990 
2999  const Tensor<3, spacedim> &
3000  jacobian_pushed_forward_grad(const unsigned int quadrature_point) const;
3001 
3008  const std::vector<Tensor<3, spacedim>> &
3010 
3019  jacobian_2nd_derivative(const unsigned int quadrature_point) const;
3020 
3027  const std::vector<DerivativeForm<3, dim, spacedim>> &
3029 
3039  const Tensor<4, spacedim> &
3041  const unsigned int quadrature_point) const;
3042 
3049  const std::vector<Tensor<4, spacedim>> &
3051 
3061  jacobian_3rd_derivative(const unsigned int quadrature_point) const;
3062 
3069  const std::vector<DerivativeForm<4, dim, spacedim>> &
3071 
3081  const Tensor<5, spacedim> &
3083  const unsigned int quadrature_point) const;
3084 
3091  const std::vector<Tensor<5, spacedim>> &
3093 
3101  inverse_jacobian(const unsigned int quadrature_point) const;
3102 
3109  const std::vector<DerivativeForm<1, spacedim, dim>> &
3110  get_inverse_jacobians() const;
3111 
3125  const Tensor<1, spacedim> &
3126  normal_vector(const unsigned int i) const;
3127 
3138  DEAL_II_DEPRECATED
3139  const std::vector<Tensor<1, spacedim>> &
3140  get_all_normal_vectors() const;
3141 
3149  const std::vector<Tensor<1, spacedim>> &
3150  get_normal_vectors() const;
3151 
3153 
3155 
3156 
3166  operator[](const FEValuesExtractors::Scalar &scalar) const;
3167 
3177  operator[](const FEValuesExtractors::Vector &vector) const;
3178 
3190 
3191 
3201  operator[](const FEValuesExtractors::Tensor<2> &tensor) const;
3202 
3204 
3206 
3207 
3211  const Mapping<dim, spacedim> &
3212  get_mapping() const;
3213 
3218  get_fe() const;
3219 
3223  UpdateFlags
3224  get_update_flags() const;
3225 
3230  get_cell() const;
3231 
3238  get_cell_similarity() const;
3239 
3244  std::size_t
3245  memory_consumption() const;
3247 
3248 
3257  std::string,
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 "
3262  << "the <" << arg1
3263  << "> flag to be set, but it was apparently not specified "
3264  << "upon construction.");
3265 
3274  "The FiniteElement you provided to FEValues and the FiniteElement that belongs "
3275  "to the DoFHandler that provided the cell iterator do not match.");
3282  int,
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?");
3289 
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.");
3300 
3301 protected:
3330  class CellIteratorBase;
3331 
3336  template <typename CI>
3339 
3345  std::unique_ptr<const CellIteratorBase> present_cell;
3346 
3354  boost::signals2::connection tria_listener_refinement;
3355 
3363  boost::signals2::connection tria_listener_mesh_transform;
3364 
3370  void
3372 
3382  void
3384  const typename Triangulation<dim, spacedim>::cell_iterator &cell);
3385 
3391 
3397  std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
3399 
3406 
3407 
3415 
3421  std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
3423 
3429  spacedim>
3431 
3432 
3437 
3446  UpdateFlags
3448 
3455 
3461  void
3463  const typename Triangulation<dim, spacedim>::cell_iterator &cell);
3464 
3465 private:
3470 
3475  template <int, int>
3477  template <int, int>
3478  friend class FEValuesViews::Vector;
3479  template <int, int, int>
3480  friend class FEValuesViews::SymmetricTensor;
3481  template <int, int, int>
3482  friend class FEValuesViews::Tensor;
3483 };
3484 
3485 
3486 
3497 template <int dim, int spacedim = dim>
3498 class FEValues : public FEValuesBase<dim, spacedim>
3499 {
3500 public:
3505  static const unsigned int integral_dimension = dim;
3506 
3513  const Quadrature<dim> & quadrature,
3514  const UpdateFlags update_flags);
3515 
3522  const Quadrature<dim> & quadrature,
3523  const UpdateFlags update_flags);
3524 
3531  template <template <int, int> class DoFHandlerType, bool level_dof_access>
3532  void
3533  reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
3534  level_dof_access>> &cell);
3535 
3549  void
3551 
3556  const Quadrature<dim> &
3557  get_quadrature() const;
3558 
3563  std::size_t
3564  memory_consumption() const;
3565 
3580  const FEValues<dim, spacedim> &
3581  get_present_fe_values() const;
3582 
3583 private:
3588 
3592  void
3594 
3601  void
3602  do_reinit();
3603 };
3604 
3605 
3616 template <int dim, int spacedim = dim>
3617 class FEFaceValuesBase : public FEValuesBase<dim, spacedim>
3618 {
3619 public:
3624  static const unsigned int integral_dimension = dim - 1;
3625 
3637  FEFaceValuesBase(const unsigned int n_q_points,
3638  const unsigned int dofs_per_cell,
3639  const UpdateFlags update_flags,
3643 
3651  const Tensor<1, spacedim> &
3652  boundary_form(const unsigned int i) const;
3653 
3660  const std::vector<Tensor<1, spacedim>> &
3661  get_boundary_forms() const;
3662 
3667  unsigned int
3668  get_face_index() const;
3669 
3674  const Quadrature<dim - 1> &
3675  get_quadrature() const;
3676 
3681  std::size_t
3682  memory_consumption() const;
3683 
3684 protected:
3689  unsigned int present_face_index;
3690 
3694  const Quadrature<dim - 1> quadrature;
3695 };
3696 
3697 
3698 
3713 template <int dim, int spacedim = dim>
3714 class FEFaceValues : public FEFaceValuesBase<dim, spacedim>
3715 {
3716 public:
3721  static const unsigned int dimension = dim;
3722 
3723  static const unsigned int space_dimension = spacedim;
3724 
3729  static const unsigned int integral_dimension = dim - 1;
3730 
3738  const UpdateFlags update_flags);
3739 
3747  const UpdateFlags update_flags);
3748 
3753  template <template <int, int> class DoFHandlerType, bool level_dof_access>
3754  void
3755  reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
3756  level_dof_access>> &cell,
3757  const unsigned int face_no);
3758 
3772  void
3774  const unsigned int face_no);
3775 
3791  get_present_fe_values() const;
3792 
3793 private:
3797  void
3799 
3806  void
3807  do_reinit(const unsigned int face_no);
3808 };
3809 
3810 
3828 template <int dim, int spacedim = dim>
3829 class FESubfaceValues : public FEFaceValuesBase<dim, spacedim>
3830 {
3831 public:
3835  static const unsigned int dimension = dim;
3836 
3840  static const unsigned int space_dimension = spacedim;
3841 
3846  static const unsigned int integral_dimension = dim - 1;
3847 
3854  const Quadrature<dim - 1> & face_quadrature,
3855  const UpdateFlags update_flags);
3856 
3863  const Quadrature<dim - 1> & face_quadrature,
3864  const UpdateFlags update_flags);
3865 
3872  template <template <int, int> class DoFHandlerType, bool level_dof_access>
3873  void
3874  reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
3875  level_dof_access>> &cell,
3876  const unsigned int face_no,
3877  const unsigned int subface_no);
3878 
3892  void
3894  const unsigned int face_no,
3895  const unsigned int subface_no);
3896 
3912  get_present_fe_values() const;
3913 
3920 
3927 
3928 private:
3932  void
3934 
3941  void
3942  do_reinit(const unsigned int face_no, const unsigned int subface_no);
3943 };
3944 
3945 
3946 #ifndef DOXYGEN
3947 
3948 
3949 /*------------------------ Inline functions: namespace FEValuesViews --------*/
3950 
3951 namespace FEValuesViews
3952 {
3953  template <int dim, int spacedim>
3954  inline typename Scalar<dim, spacedim>::value_type
3955  Scalar<dim, spacedim>::value(const unsigned int shape_function,
3956  const unsigned int q_point) const
3957  {
3958  Assert(shape_function < fe_values->fe->dofs_per_cell,
3959  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
3960  Assert(
3961  fe_values->update_flags & update_values,
3963  "update_values"))));
3964 
3965  // an adaptation of the FEValuesBase::shape_value_component function
3966  // except that here we know the component as fixed and we have
3967  // pre-computed and cached a bunch of information. See the comments there.
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);
3971  else
3972  return 0;
3973  }
3974 
3975 
3976 
3977  template <int dim, int spacedim>
3978  inline typename Scalar<dim, spacedim>::gradient_type
3979  Scalar<dim, spacedim>::gradient(const unsigned int shape_function,
3980  const unsigned int q_point) const
3981  {
3982  Assert(shape_function < fe_values->fe->dofs_per_cell,
3983  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
3984  Assert(fe_values->update_flags & update_gradients,
3986  "update_gradients")));
3987 
3988  // an adaptation of the FEValuesBase::shape_grad_component
3989  // function except that here we know the component as fixed and we have
3990  // pre-computed and cached a bunch of information. See the comments there.
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]
3994  [q_point];
3995  else
3996  return gradient_type();
3997  }
3998 
3999 
4000 
4001  template <int dim, int spacedim>
4002  inline typename Scalar<dim, spacedim>::hessian_type
4003  Scalar<dim, spacedim>::hessian(const unsigned int shape_function,
4004  const unsigned int q_point) const
4005  {
4006  Assert(shape_function < fe_values->fe->dofs_per_cell,
4007  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4008  Assert(fe_values->update_flags & update_hessians,
4010  "update_hessians")));
4011 
4012  // an adaptation of the FEValuesBase::shape_hessian_component
4013  // function except that here we know the component as fixed and we have
4014  // pre-computed and cached a bunch of information. See the comments there.
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];
4018  else
4019  return hessian_type();
4020  }
4021 
4022 
4023 
4024  template <int dim, int spacedim>
4026  Scalar<dim, spacedim>::third_derivative(const unsigned int shape_function,
4027  const unsigned int q_point) const
4028  {
4029  Assert(shape_function < fe_values->fe->dofs_per_cell,
4030  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4031  Assert(fe_values->update_flags & update_3rd_derivatives,
4033  "update_3rd_derivatives")));
4034 
4035  // an adaptation of the FEValuesBase::shape_3rdderivative_component
4036  // function except that here we know the component as fixed and we have
4037  // pre-computed and cached a bunch of information. See the comments there.
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]
4041  [q_point];
4042  else
4043  return third_derivative_type();
4044  }
4045 
4046 
4047 
4048  template <int dim, int spacedim>
4049  inline typename Vector<dim, spacedim>::value_type
4050  Vector<dim, spacedim>::value(const unsigned int shape_function,
4051  const unsigned int q_point) const
4052  {
4053  Assert(shape_function < fe_values->fe->dofs_per_cell,
4054  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4055  Assert(fe_values->update_flags & update_values,
4057  "update_values")));
4058 
4059  // same as for the scalar case except that we have one more index
4060  const int snc =
4061  shape_function_data[shape_function].single_nonzero_component;
4062  if (snc == -2)
4063  return value_type();
4064  else if (snc != -1)
4065  {
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;
4071  }
4072  else
4073  {
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);
4080 
4081  return return_value;
4082  }
4083  }
4084 
4085 
4086 
4087  template <int dim, int spacedim>
4088  inline typename Vector<dim, spacedim>::gradient_type
4089  Vector<dim, spacedim>::gradient(const unsigned int shape_function,
4090  const unsigned int q_point) const
4091  {
4092  Assert(shape_function < fe_values->fe->dofs_per_cell,
4093  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4094  Assert(fe_values->update_flags & update_gradients,
4096  "update_gradients")));
4097 
4098  // same as for the scalar case except that we have one more index
4099  const int snc =
4100  shape_function_data[shape_function].single_nonzero_component;
4101  if (snc == -2)
4102  return gradient_type();
4103  else if (snc != -1)
4104  {
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;
4110  }
4111  else
4112  {
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])
4117  return_value[d] =
4118  fe_values->finite_element_output.shape_gradients
4119  [shape_function_data[shape_function].row_index[d]][q_point];
4120 
4121  return return_value;
4122  }
4123  }
4124 
4125 
4126 
4127  template <int dim, int spacedim>
4129  Vector<dim, spacedim>::divergence(const unsigned int shape_function,
4130  const unsigned int q_point) const
4131  {
4132  // this function works like in the case above
4133  Assert(shape_function < fe_values->fe->dofs_per_cell,
4134  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4135  Assert(fe_values->update_flags & update_gradients,
4137  "update_gradients")));
4138 
4139  // same as for the scalar case except that we have one more index
4140  const int snc =
4141  shape_function_data[shape_function].single_nonzero_component;
4142  if (snc == -2)
4143  return divergence_type();
4144  else if (snc != -1)
4145  return fe_values->finite_element_output
4146  .shape_gradients[snc][q_point][shape_function_data[shape_function]
4147  .single_nonzero_component_index];
4148  else
4149  {
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])
4154  return_value +=
4155  fe_values->finite_element_output.shape_gradients
4156  [shape_function_data[shape_function].row_index[d]][q_point][d];
4157 
4158  return return_value;
4159  }
4160  }
4161 
4162 
4163 
4164  template <int dim, int spacedim>
4165  inline typename Vector<dim, spacedim>::curl_type
4166  Vector<dim, spacedim>::curl(const unsigned int shape_function,
4167  const unsigned int q_point) const
4168  {
4169  // this function works like in the case above
4170 
4171  Assert(shape_function < fe_values->fe->dofs_per_cell,
4172  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4173  Assert(fe_values->update_flags & update_gradients,
4175  "update_gradients")));
4176  // same as for the scalar case except that we have one more index
4177  const int snc =
4178  shape_function_data[shape_function].single_nonzero_component;
4179 
4180  if (snc == -2)
4181  return curl_type();
4182 
4183  else
4184  switch (dim)
4185  {
4186  case 1:
4187  {
4188  Assert(false,
4189  ExcMessage(
4190  "Computing the curl in 1d is not a useful operation"));
4191  return curl_type();
4192  }
4193 
4194  case 2:
4195  {
4196  if (snc != -1)
4197  {
4198  curl_type return_value;
4199 
4200  // the single nonzero component can only be zero or one in 2d
4201  if (shape_function_data[shape_function]
4202  .single_nonzero_component_index == 0)
4203  return_value[0] =
4204  -1.0 * fe_values->finite_element_output
4205  .shape_gradients[snc][q_point][1];
4206  else
4207  return_value[0] = fe_values->finite_element_output
4208  .shape_gradients[snc][q_point][0];
4209 
4210  return return_value;
4211  }
4212 
4213  else
4214  {
4215  curl_type return_value;
4216 
4217  return_value[0] = 0.0;
4218 
4219  if (shape_function_data[shape_function]
4220  .is_nonzero_shape_function_component[0])
4221  return_value[0] -=
4222  fe_values->finite_element_output
4223  .shape_gradients[shape_function_data[shape_function]
4224  .row_index[0]][q_point][1];
4225 
4226  if (shape_function_data[shape_function]
4227  .is_nonzero_shape_function_component[1])
4228  return_value[0] +=
4229  fe_values->finite_element_output
4230  .shape_gradients[shape_function_data[shape_function]
4231  .row_index[1]][q_point][0];
4232 
4233  return return_value;
4234  }
4235  }
4236 
4237  case 3:
4238  {
4239  if (snc != -1)
4240  {
4241  curl_type return_value;
4242 
4243  switch (shape_function_data[shape_function]
4244  .single_nonzero_component_index)
4245  {
4246  case 0:
4247  {
4248  return_value[0] = 0;
4249  return_value[1] = fe_values->finite_element_output
4250  .shape_gradients[snc][q_point][2];
4251  return_value[2] =
4252  -1.0 * fe_values->finite_element_output
4253  .shape_gradients[snc][q_point][1];
4254  return return_value;
4255  }
4256 
4257  case 1:
4258  {
4259  return_value[0] =
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;
4266  }
4267 
4268  default:
4269  {
4270  return_value[0] = fe_values->finite_element_output
4271  .shape_gradients[snc][q_point][1];
4272  return_value[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;
4277  }
4278  }
4279  }
4280 
4281  else
4282  {
4283  curl_type return_value;
4284 
4285  for (unsigned int i = 0; i < dim; ++i)
4286  return_value[i] = 0.0;
4287 
4288  if (shape_function_data[shape_function]
4289  .is_nonzero_shape_function_component[0])
4290  {
4291  return_value[1] +=
4292  fe_values->finite_element_output
4293  .shape_gradients[shape_function_data[shape_function]
4294  .row_index[0]][q_point][2];
4295  return_value[2] -=
4296  fe_values->finite_element_output
4297  .shape_gradients[shape_function_data[shape_function]
4298  .row_index[0]][q_point][1];
4299  }
4300 
4301  if (shape_function_data[shape_function]
4302  .is_nonzero_shape_function_component[1])
4303  {
4304  return_value[0] -=
4305  fe_values->finite_element_output
4306  .shape_gradients[shape_function_data[shape_function]
4307  .row_index[1]][q_point][2];
4308  return_value[2] +=
4309  fe_values->finite_element_output
4310  .shape_gradients[shape_function_data[shape_function]
4311  .row_index[1]][q_point][0];
4312  }
4313 
4314  if (shape_function_data[shape_function]
4315  .is_nonzero_shape_function_component[2])
4316  {
4317  return_value[0] +=
4318  fe_values->finite_element_output
4319  .shape_gradients[shape_function_data[shape_function]
4320  .row_index[2]][q_point][1];
4321  return_value[1] -=
4322  fe_values->finite_element_output
4323  .shape_gradients[shape_function_data[shape_function]
4324  .row_index[2]][q_point][0];
4325  }
4326 
4327  return return_value;
4328  }
4329  }
4330  }
4331  // should not end up here
4332  Assert(false, ExcInternalError());
4333  return curl_type();
4334  }
4335 
4336 
4337 
4338  template <int dim, int spacedim>
4339  inline typename Vector<dim, spacedim>::hessian_type
4340  Vector<dim, spacedim>::hessian(const unsigned int shape_function,
4341  const unsigned int q_point) const
4342  {
4343  // this function works like in the case above
4344  Assert(shape_function < fe_values->fe->dofs_per_cell,
4345  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4346  Assert(fe_values->update_flags & update_hessians,
4348  "update_hessians")));
4349 
4350  // same as for the scalar case except that we have one more index
4351  const int snc =
4352  shape_function_data[shape_function].single_nonzero_component;
4353  if (snc == -2)
4354  return hessian_type();
4355  else if (snc != -1)
4356  {
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;
4362  }
4363  else
4364  {
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])
4369  return_value[d] =
4370  fe_values->finite_element_output.shape_hessians
4371  [shape_function_data[shape_function].row_index[d]][q_point];
4372 
4373  return return_value;
4374  }
4375  }
4376 
4377 
4378 
4379  template <int dim, int spacedim>
4381  Vector<dim, spacedim>::third_derivative(const unsigned int shape_function,
4382  const unsigned int q_point) const
4383  {
4384  // this function works like in the case above
4385  Assert(shape_function < fe_values->fe->dofs_per_cell,
4386  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4387  Assert(fe_values->update_flags & update_3rd_derivatives,
4389  "update_3rd_derivatives")));
4390 
4391  // same as for the scalar case except that we have one more index
4392  const int snc =
4393  shape_function_data[shape_function].single_nonzero_component;
4394  if (snc == -2)
4395  return third_derivative_type();
4396  else if (snc != -1)
4397  {
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;
4403  }
4404  else
4405  {
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])
4410  return_value[d] =
4411  fe_values->finite_element_output.shape_3rd_derivatives
4412  [shape_function_data[shape_function].row_index[d]][q_point];
4413 
4414  return return_value;
4415  }
4416  }
4417 
4418 
4419 
4420  namespace internal
4421  {
4426  inline ::SymmetricTensor<2, 1>
4427  symmetrize_single_row(const unsigned int n, const ::Tensor<1, 1> &t)
4428  {
4429  Assert(n < 1, ExcIndexRange(n, 0, 1));
4430  (void)n;
4431 
4432  return {{t[0]}};
4433  }
4434 
4435 
4436 
4437  inline ::SymmetricTensor<2, 2>
4438  symmetrize_single_row(const unsigned int n, const ::Tensor<1, 2> &t)
4439  {
4440  switch (n)
4441  {
4442  case 0:
4443  {
4444  return {{t[0], 0, t[1] / 2}};
4445  }
4446  case 1:
4447  {
4448  return {{0, t[1], t[0] / 2}};
4449  }
4450  default:
4451  {
4452  Assert(false, ExcIndexRange(n, 0, 2));
4453  return {};
4454  }
4455  }
4456  }
4457 
4458 
4459 
4460  inline ::SymmetricTensor<2, 3>
4461  symmetrize_single_row(const unsigned int n, const ::Tensor<1, 3> &t)
4462  {
4463  switch (n)
4464  {
4465  case 0:
4466  {
4467  return {{t[0], 0, 0, t[1] / 2, t[2] / 2, 0}};
4468  }
4469  case 1:
4470  {
4471  return {{0, t[1], 0, t[0] / 2, 0, t[2] / 2}};
4472  }
4473  case 2:
4474  {
4475  return {{0, 0, t[2], 0, t[0] / 2, t[1] / 2}};
4476  }
4477  default:
4478  {
4479  Assert(false, ExcIndexRange(n, 0, 3));
4480  return {};
4481  }
4482  }
4483  }
4484  } // namespace internal
4485 
4486 
4487 
4488  template <int dim, int spacedim>
4490  Vector<dim, spacedim>::symmetric_gradient(const unsigned int shape_function,
4491  const unsigned int q_point) const
4492  {
4493  Assert(shape_function < fe_values->fe->dofs_per_cell,
4494  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4495  Assert(fe_values->update_flags & update_gradients,
4497  "update_gradients")));
4498 
4499  // same as for the scalar case except that we have one more index
4500  const int snc =
4501  shape_function_data[shape_function].single_nonzero_component;
4502  if (snc == -2)
4503  return symmetric_gradient_type();
4504  else if (snc != -1)
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]);
4508  else
4509  {
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])
4514  return_value[d] =
4515  fe_values->finite_element_output.shape_gradients
4516  [shape_function_data[shape_function].row_index[d]][q_point];
4517 
4518  return symmetrize(return_value);
4519  }
4520  }
4521 
4522 
4523 
4524  template <int dim, int spacedim>
4526  SymmetricTensor<2, dim, spacedim>::value(const unsigned int shape_function,
4527  const unsigned int q_point) const
4528  {
4529  Assert(shape_function < fe_values->fe->dofs_per_cell,
4530  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4531  Assert(fe_values->update_flags & update_values,
4533  "update_values")));
4534 
4535  // similar to the vector case where we have more then one index and we need
4536  // to convert between unrolled and component indexing for tensors
4537  const int snc =
4538  shape_function_data[shape_function].single_nonzero_component;
4539 
4540  if (snc == -2)
4541  {
4542  // shape function is zero for the selected components
4543  return value_type();
4544  }
4545  else if (snc != -1)
4546  {
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;
4553  }
4554  else
4555  {
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;
4564  }
4565  }
4566 
4567 
4568 
4569  template <int dim, int spacedim>
4572  const unsigned int shape_function,
4573  const unsigned int q_point) const
4574  {
4575  Assert(shape_function < fe_values->fe->dofs_per_cell,
4576  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4577  Assert(fe_values->update_flags & update_gradients,
4579  "update_gradients")));
4580 
4581  const int snc =
4582  shape_function_data[shape_function].single_nonzero_component;
4583 
4584  if (snc == -2)
4585  {
4586  // shape function is zero for the selected components
4587  return divergence_type();
4588  }
4589  else if (snc != -1)
4590  {
4591  // we have a single non-zero component when the symmetric tensor is
4592  // represented in unrolled form. this implies we potentially have
4593  // two non-zero components when represented in component form! we
4594  // will only have one non-zero entry if the non-zero component lies on
4595  // the diagonal of the tensor.
4596  //
4597  // the divergence of a second-order tensor is a first order tensor.
4598  //
4599  // assume the second-order tensor is A with components A_{ij}. then
4600  // A_{ij} = A_{ji} and there is only one (if diagonal) or two non-zero
4601  // entries in the tensorial representation. define the
4602  // divergence as:
4603  // b_i \dealcoloneq \dfrac{\partial phi_{ij}}{\partial x_j}.
4604  // (which is incidentally also
4605  // b_j \dealcoloneq \dfrac{\partial phi_{ij}}{\partial x_i}).
4606  // In both cases, a sum is implied.
4607  //
4608  // Now, we know the nonzero component in unrolled form: it is indicated
4609  // by 'snc'. we can figure out which tensor components belong to this:
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];
4616 
4617  // given the form of the divergence above, if ii=jj there is only a
4618  // single nonzero component of the full tensor and the gradient
4619  // equals
4620  // b_ii \dealcoloneq \dfrac{\partial phi_{ii,ii}}{\partial x_ii}.
4621  // all other entries of 'b' are zero
4622  //
4623  // on the other hand, if ii!=jj, then there are two nonzero entries in
4624  // the full tensor and
4625  // b_ii \dealcoloneq \dfrac{\partial phi_{ii,jj}}{\partial x_ii}.
4626  // b_jj \dealcoloneq \dfrac{\partial phi_{ii,jj}}{\partial x_jj}.
4627  // again, all other entries of 'b' are zero
4628  const ::Tensor<1, spacedim> &phi_grad =
4629  fe_values->finite_element_output.shape_gradients[snc][q_point];
4630 
4631  divergence_type return_value;
4632  return_value[ii] = phi_grad[jj];
4633 
4634  if (ii != jj)
4635  return_value[jj] = phi_grad[ii];
4636 
4637  return return_value;
4638  }
4639  else
4640  {
4641  Assert(false, ExcNotImplemented());
4642  divergence_type return_value;
4643  return return_value;
4644  }
4645  }
4646 
4647 
4648 
4649  template <int dim, int spacedim>
4650  inline typename Tensor<2, dim, spacedim>::value_type
4651  Tensor<2, dim, spacedim>::value(const unsigned int shape_function,
4652  const unsigned int q_point) const
4653  {
4654  Assert(shape_function < fe_values->fe->dofs_per_cell,
4655  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4656  Assert(fe_values->update_flags & update_values,
4658  "update_values")));
4659 
4660  // similar to the vector case where we have more then one index and we need
4661  // to convert between unrolled and component indexing for tensors
4662  const int snc =
4663  shape_function_data[shape_function].single_nonzero_component;
4664 
4665  if (snc == -2)
4666  {
4667  // shape function is zero for the selected components
4668  return value_type();
4669  }
4670  else if (snc != -1)
4671  {
4672  value_type return_value;
4673  const unsigned int comp =
4674  shape_function_data[shape_function].single_nonzero_component_index;
4675  const TableIndices<2> indices =
4677  return_value[indices] =
4678  fe_values->finite_element_output.shape_values(snc, q_point);
4679  return return_value;
4680  }
4681  else
4682  {
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])
4687  {
4688  const TableIndices<2> indices =
4690  return_value[indices] =
4691  fe_values->finite_element_output.shape_values(
4692  shape_function_data[shape_function].row_index[d], q_point);
4693  }
4694  return return_value;
4695  }
4696  }
4697 
4698 
4699 
4700  template <int dim, int spacedim>
4702  Tensor<2, dim, spacedim>::divergence(const unsigned int shape_function,
4703  const unsigned int q_point) const
4704  {
4705  Assert(shape_function < fe_values->fe->dofs_per_cell,
4706  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4707  Assert(fe_values->update_flags & update_gradients,
4709  "update_gradients")));
4710 
4711  const int snc =
4712  shape_function_data[shape_function].single_nonzero_component;
4713 
4714  if (snc == -2)
4715  {
4716  // shape function is zero for the selected components
4717  return divergence_type();
4718  }
4719  else if (snc != -1)
4720  {
4721  // we have a single non-zero component when the tensor is
4722  // represented in unrolled form.
4723  //
4724  // the divergence of a second-order tensor is a first order tensor.
4725  //
4726  // assume the second-order tensor is A with components A_{ij},
4727  // then divergence is d_i := \frac{\partial A_{ij}}{\partial x_j}
4728  //
4729  // Now, we know the nonzero component in unrolled form: it is indicated
4730  // by 'snc'. we can figure out which tensor components belong to this:
4731  const unsigned int comp =
4732  shape_function_data[shape_function].single_nonzero_component_index;
4733  const TableIndices<2> indices =
4735  const unsigned int ii = indices[0];
4736  const unsigned int jj = indices[1];
4737 
4738  const ::Tensor<1, spacedim> &phi_grad =
4739  fe_values->finite_element_output.shape_gradients[snc][q_point];
4740 
4741  divergence_type return_value;
4742  // note that we contract \nabla from the right
4743  return_value[ii] = phi_grad[jj];
4744 
4745  return return_value;
4746  }
4747  else
4748  {
4749  Assert(false, ExcNotImplemented());
4750  divergence_type return_value;
4751  return return_value;
4752  }
4753  }
4754 
4755 
4756 
4757  template <int dim, int spacedim>
4759  Tensor<2, dim, spacedim>::gradient(const unsigned int shape_function,
4760  const unsigned int q_point) const
4761  {
4762  Assert(shape_function < fe_values->fe->dofs_per_cell,
4763  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4764  Assert(fe_values->update_flags & update_gradients,
4766  "update_gradients")));
4767 
4768  const int snc =
4769  shape_function_data[shape_function].single_nonzero_component;
4770 
4771  if (snc == -2)
4772  {
4773  // shape function is zero for the selected components
4774  return gradient_type();
4775  }
4776  else if (snc != -1)
4777  {
4778  // we have a single non-zero component when the tensor is
4779  // represented in unrolled form.
4780  //
4781  // the gradient of a second-order tensor is a third order tensor.
4782  //
4783  // assume the second-order tensor is A with components A_{ij},
4784  // then gradient is B_{ijk} := \frac{\partial A_{ij}}{\partial x_k}
4785  //
4786  // Now, we know the nonzero component in unrolled form: it is indicated
4787  // by 'snc'. we can figure out which tensor components belong to this:
4788  const unsigned int comp =
4789  shape_function_data[shape_function].single_nonzero_component_index;
4790  const TableIndices<2> indices =
4792  const unsigned int ii = indices[0];
4793  const unsigned int jj = indices[1];
4794 
4795  const ::Tensor<1, spacedim> &phi_grad =
4796  fe_values->finite_element_output.shape_gradients[snc][q_point];
4797 
4798  gradient_type return_value;
4799  return_value[ii][jj] = phi_grad;
4800 
4801  return return_value;
4802  }
4803  else
4804  {
4805  Assert(false, ExcNotImplemented());
4806  gradient_type return_value;
4807  return return_value;
4808  }
4809  }
4810 
4811 } // namespace FEValuesViews
4812 
4813 
4814 
4815 /*---------------------- Inline functions: FEValuesBase ---------------------*/
4816 
4817 
4818 
4819 template <int dim, int spacedim>
4821  operator[](const FEValuesExtractors::Scalar &scalar) const
4822 {
4823  Assert(scalar.component < fe_values_views_cache.scalars.size(),
4824  ExcIndexRange(scalar.component,
4825  0,
4826  fe_values_views_cache.scalars.size()));
4827 
4828  return fe_values_views_cache.scalars[scalar.component];
4829 }
4830 
4831 
4832 
4833 template <int dim, int spacedim>
4835  operator[](const FEValuesExtractors::Vector &vector) const
4836 {
4837  Assert(vector.first_vector_component < fe_values_views_cache.vectors.size(),
4839  0,
4840  fe_values_views_cache.vectors.size()));
4841 
4842  return fe_values_views_cache.vectors[vector.first_vector_component];
4843 }
4844 
4845 
4846 
4847 template <int dim, int spacedim>
4851 {
4852  Assert(
4853  tensor.first_tensor_component <
4854  fe_values_views_cache.symmetric_second_order_tensors.size(),
4856  0,
4857  fe_values_views_cache.symmetric_second_order_tensors.size()));
4858 
4859  return fe_values_views_cache
4860  .symmetric_second_order_tensors[tensor.first_tensor_component];
4861 }
4862 
4863 
4864 
4865 template <int dim, int spacedim>
4868  operator[](const FEValuesExtractors::Tensor<2> &tensor) const
4869 {
4871  fe_values_views_cache.second_order_tensors.size(),
4873  0,
4874  fe_values_views_cache.second_order_tensors.size()));
4875 
4876  return fe_values_views_cache
4877  .second_order_tensors[tensor.first_tensor_component];
4878 }
4879 
4880 
4881 
4882 template <int dim, int spacedim>
4883 inline const double &
4884 FEValuesBase<dim, spacedim>::shape_value(const unsigned int i,
4885  const unsigned int j) const
4886 {
4887  Assert(i < fe->dofs_per_cell, ExcIndexRange(i, 0, fe->dofs_per_cell));
4888  Assert(this->update_flags & update_values,
4889  ExcAccessToUninitializedField("update_values"));
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"));
4893  // if the entire FE is primitive,
4894  // then we can take a short-cut:
4895  if (fe->is_primitive())
4896  return this->finite_element_output.shape_values(i, j);
4897  else
4898  {
4899  // otherwise, use the mapping
4900  // between shape function
4901  // numbers and rows. note that
4902  // by the assertions above, we
4903  // know that this particular
4904  // shape function is primitive,
4905  // so we can call
4906  // system_to_component_index
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);
4912  }
4913 }
4914 
4915 
4916 
4917 template <int dim, int spacedim>
4918 inline double
4920  const unsigned int i,
4921  const unsigned int j,
4922  const unsigned int component) const
4923 {
4924  Assert(i < fe->dofs_per_cell, ExcIndexRange(i, 0, fe->dofs_per_cell));
4925  Assert(this->update_flags & update_values,
4926  ExcAccessToUninitializedField("update_values"));
4927  Assert(component < fe->n_components(),
4928  ExcIndexRange(component, 0, fe->n_components()));
4929  Assert(present_cell.get() != nullptr,
4930  ExcMessage("FEValues object is not reinit'ed to any cell"));
4931 
4932  // check whether the shape function
4933  // is non-zero at all within
4934  // this component:
4935  if (fe->get_nonzero_components(i)[component] == false)
4936  return 0;
4937 
4938  // look up the right row in the
4939  // table and take the data from
4940  // there
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);
4945 }
4946 
4947 
4948 
4949 template <int dim, int spacedim>
4950 inline const Tensor<1, spacedim> &
4951 FEValuesBase<dim, spacedim>::shape_grad(const unsigned int i,
4952  const unsigned int j) const
4953 {
4954  Assert(i < fe->dofs_per_cell, ExcIndexRange(i, 0, fe->dofs_per_cell));
4955  Assert(this->update_flags & update_gradients,
4956  ExcAccessToUninitializedField("update_gradients"));
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"));
4960  // if the entire FE is primitive,
4961  // then we can take a short-cut:
4962  if (fe->is_primitive())
4963  return this->finite_element_output.shape_gradients[i][j];
4964  else
4965  {
4966  // otherwise, use the mapping
4967  // between shape function
4968  // numbers and rows. note that
4969  // by the assertions above, we
4970  // know that this particular
4971  // shape function is primitive,
4972  // so we can call
4973  // system_to_component_index
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];
4979  }
4980 }
4981 
4982 
4983 
4984 template <int dim, int spacedim>
4985 inline Tensor<1, spacedim>
4987  const unsigned int i,
4988  const unsigned int j,
4989  const unsigned int component) const
4990 {
4991  Assert(i < fe->dofs_per_cell, ExcIndexRange(i, 0, fe->dofs_per_cell));
4992  Assert(this->update_flags & update_gradients,
4993  ExcAccessToUninitializedField("update_gradients"));
4994  Assert(component < fe->n_components(),
4995  ExcIndexRange(component, 0, fe->n_components()));
4996  Assert(present_cell.get() != nullptr,
4997  ExcMessage("FEValues object is not reinit'ed to any cell"));
4998  // check whether the shape function
4999  // is non-zero at all within
5000  // this component:
5001  if (fe->get_nonzero_components(i)[component] == false)
5002  return Tensor<1, spacedim>();
5003 
5004  // look up the right row in the
5005  // table and take the data from
5006  // there
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];
5011 }
5012 
5013 
5014 
5015 template <int dim, int spacedim>
5016 inline const Tensor<2, spacedim> &
5017 FEValuesBase<dim, spacedim>::shape_hessian(const unsigned int i,
5018  const unsigned int j) const
5019 {
5020  Assert(i < fe->dofs_per_cell, ExcIndexRange(i, 0, fe->dofs_per_cell));
5021  Assert(this->update_flags & update_hessians,
5022  ExcAccessToUninitializedField("update_hessians"));
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"));
5026  // if the entire FE is primitive,
5027  // then we can take a short-cut:
5028  if (fe->is_primitive())
5029  return this->finite_element_output.shape_hessians[i][j];
5030  else
5031  {
5032  // otherwise, use the mapping
5033  // between shape function
5034  // numbers and rows. note that
5035  // by the assertions above, we
5036  // know that this particular
5037  // shape function is primitive,
5038  // so we can call
5039  // system_to_component_index
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];
5045  }
5046 }
5047 
5048 
5049 
5050 template <int dim, int spacedim>
5051 inline Tensor<2, spacedim>
5053  const unsigned int i,
5054  const unsigned int j,
5055  const unsigned int component) const
5056 {
5057  Assert(i < fe->dofs_per_cell, ExcIndexRange(i, 0, fe->dofs_per_cell));
5058  Assert(this->update_flags & update_hessians,
5059  ExcAccessToUninitializedField("update_hessians"));
5060  Assert(component < fe->n_components(),
5061  ExcIndexRange(component, 0, fe->n_components()));
5062  Assert(present_cell.get() != nullptr,
5063  ExcMessage("FEValues object is not reinit'ed to any cell"));
5064  // check whether the shape function
5065  // is non-zero at all within
5066  // this component:
5067  if (fe->get_nonzero_components(i)[component] == false)
5068  return Tensor<2, spacedim>();
5069 
5070  // look up the right row in the
5071  // table and take the data from
5072  // there
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];
5077 }
5078 
5079 
5080 
5081 template <int dim, int spacedim>
5082 inline const Tensor<3, spacedim> &
5084  const unsigned int j) const
5085 {
5086  Assert(i < fe->dofs_per_cell, ExcIndexRange(i, 0, fe->dofs_per_cell));
5087  Assert(this->update_flags & update_hessians,
5088  ExcAccessToUninitializedField("update_3rd_derivatives"));
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"));
5092  // if the entire FE is primitive,
5093  // then we can take a short-cut:
5094  if (fe->is_primitive())
5095  return this->finite_element_output.shape_3rd_derivatives[i][j];
5096  else
5097  {
5098  // otherwise, use the mapping
5099  // between shape function
5100  // numbers and rows. note that
5101  // by the assertions above, we
5102  // know that this particular
5103  // shape function is primitive,
5104  // so we can call
5105  // system_to_component_index
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];
5111  }
5112 }
5113 
5114 
5115 
5116 template <int dim, int spacedim>
5117 inline Tensor<3, spacedim>
5119  const unsigned int i,
5120  const unsigned int j,
5121  const unsigned int component) const
5122 {
5123  Assert(i < fe->dofs_per_cell, ExcIndexRange(i, 0, fe->dofs_per_cell));
5124  Assert(this->update_flags & update_hessians,
5125  ExcAccessToUninitializedField("update_3rd_derivatives"));
5126  Assert(component < fe->n_components(),
5127  ExcIndexRange(component, 0, fe->n_components()));
5128  Assert(present_cell.get() != nullptr,
5129  ExcMessage("FEValues object is not reinit'ed to any cell"));
5130  // check whether the shape function
5131  // is non-zero at all within
5132  // this component:
5133  if (fe->get_nonzero_components(i)[component] == false)
5134  return Tensor<3, spacedim>();
5135 
5136  // look up the right row in the
5137  // table and take the data from
5138  // there
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];
5143 }
5144 
5145 
5146 
5147 template <int dim, int spacedim>
5148 inline const FiniteElement<dim, spacedim> &
5150 {
5151  return *fe;
5152 }
5153 
5154 
5155 
5156 template <int dim, int spacedim>
5157 inline const Mapping<dim, spacedim> &
5159 {
5160  return *mapping;
5161 }
5162 
5163 
5164 
5165 template <int dim, int spacedim>
5166 inline UpdateFlags
5168 {
5169  return this->update_flags;
5170 }
5171 
5172 
5173 
5174 template <int dim, int spacedim>
5175 inline const std::vector<Point<spacedim>> &
5177 {
5178  Assert(this->update_flags & update_quadrature_points,
5179  ExcAccessToUninitializedField("update_quadrature_points"));
5180  Assert(present_cell.get() != nullptr,
5181  ExcMessage("FEValues object is not reinit'ed to any cell"));
5182  return this->mapping_output.quadrature_points;
5183 }
5184 
5185 
5186 
5187 template <int dim, int spacedim>
5188 inline const std::vector<double> &
5190 {
5191  Assert(this->update_flags & update_JxW_values,
5192  ExcAccessToUninitializedField("update_JxW_values"));
5193  Assert(present_cell.get() != nullptr,
5194  ExcMessage("FEValues object is not reinit'ed to any cell"));
5195  return this->mapping_output.JxW_values;
5196 }
5197 
5198 
5199 
5200 template <int dim, int spacedim>
5201 inline const std::vector<DerivativeForm<1, dim, spacedim>> &
5203 {
5204  Assert(this->update_flags & update_jacobians,
5205  ExcAccessToUninitializedField("update_jacobians"));
5206  Assert(present_cell.get() != nullptr,
5207  ExcMessage("FEValues object is not reinit'ed to any cell"));
5208  return this->mapping_output.jacobians;
5209 }
5210 
5211 
5212 
5213 template <int dim, int spacedim>
5214 inline const std::vector<DerivativeForm<2, dim, spacedim>> &
5216 {
5217  Assert(this->update_flags & update_jacobian_grads,
5218  ExcAccessToUninitializedField("update_jacobians_grads"));
5219  Assert(present_cell.get() != nullptr,
5220  ExcMessage("FEValues object is not reinit'ed to any cell"));
5221  return this->mapping_output.jacobian_grads;
5222 }
5223 
5224 
5225 
5226 template <int dim, int spacedim>
5227 inline const Tensor<3, spacedim> &
5229  const unsigned int i) const
5230 {
5231  Assert(this->update_flags & update_jacobian_pushed_forward_grads,
5232  ExcAccessToUninitializedField("update_jacobian_pushed_forward_grads"));
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];
5236 }
5237 
5238 
5239 
5240 template <int dim, int spacedim>
5241 inline const std::vector<Tensor<3, spacedim>> &
5243 {
5244  Assert(this->update_flags & update_jacobian_pushed_forward_grads,
5245  ExcAccessToUninitializedField("update_jacobian_pushed_forward_grads"));
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;
5249 }
5250 
5251 
5252 
5253 template <int dim, int spacedim>
5254 inline const DerivativeForm<3, dim, spacedim> &
5255 FEValuesBase<dim, spacedim>::jacobian_2nd_derivative(const unsigned int i) const
5256 {
5257  Assert(this->update_flags & update_jacobian_2nd_derivatives,
5258  ExcAccessToUninitializedField("update_jacobian_2nd_derivatives"));
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];
5262 }
5263 
5264 
5265 
5266 template <int dim, int spacedim>
5267 inline const std::vector<DerivativeForm<3, dim, spacedim>> &
5269 {
5270  Assert(this->update_flags & update_jacobian_2nd_derivatives,
5271  ExcAccessToUninitializedField("update_jacobian_2nd_derivatives"));
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;
5275 }
5276 
5277 
5278 
5279 template <int dim, int spacedim>
5280 inline const Tensor<4, spacedim> &
5282  const unsigned int i) const
5283 {
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];
5290 }
5291 
5292 
5293 
5294 template <int dim, int spacedim>
5295 inline const std::vector<Tensor<4, spacedim>> &
5297 {
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;
5304 }
5305 
5306 
5307 
5308 template <int dim, int spacedim>
5309 inline const DerivativeForm<4, dim, spacedim> &
5310 FEValuesBase<dim, spacedim>::jacobian_3rd_derivative(const unsigned int i) const
5311 {
5312  Assert(this->update_flags & update_jacobian_3rd_derivatives,
5313  ExcAccessToUninitializedField("update_jacobian_3rd_derivatives"));
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];
5317 }
5318 
5319 
5320 
5321 template <int dim, int spacedim>
5322 inline const std::vector<DerivativeForm<4, dim, spacedim>> &
5324 {
5325  Assert(this->update_flags & update_jacobian_3rd_derivatives,
5326  ExcAccessToUninitializedField("update_jacobian_3rd_derivatives"));
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;
5330 }
5331 
5332 
5333 
5334 template <int dim, int spacedim>
5335 inline const Tensor<5, spacedim> &
5337  const unsigned int i) const
5338 {
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];
5345 }
5346 
5347 
5348 
5349 template <int dim, int spacedim>
5350 inline const std::vector<Tensor<5, spacedim>> &
5352 {
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;
5359 }
5360 
5361 
5362 
5363 template <int dim, int spacedim>
5364 inline const std::vector<DerivativeForm<1, spacedim, dim>> &
5366 {
5367  Assert(this->update_flags & update_inverse_jacobians,
5368  ExcAccessToUninitializedField("update_inverse_jacobians"));
5369  Assert(present_cell.get() != nullptr,
5370  ExcMessage("FEValues object is not reinit'ed to any cell"));
5371  return this->mapping_output.inverse_jacobians;
5372 }
5373 
5374 
5375 
5376 template <int dim, int spacedim>
5377 inline const Point<spacedim> &
5378 FEValuesBase<dim, spacedim>::quadrature_point(const unsigned int i) const
5379 {
5380  Assert(this->update_flags & update_quadrature_points,
5381  ExcAccessToUninitializedField("update_quadrature_points"));
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"));
5386 
5387  return this->mapping_output.quadrature_points[i];
5388 }
5389 
5390 
5391 
5392 template <int dim, int spacedim>
5393 inline double
5394 FEValuesBase<dim, spacedim>::JxW(const unsigned int i) const
5395 {
5396  Assert(this->update_flags & update_JxW_values,
5397  ExcAccessToUninitializedField("update_JxW_values"));
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"));
5402 
5403  return this->mapping_output.JxW_values[i];
5404 }
5405 
5406 
5407 
5408 template <int dim, int spacedim>
5409 inline const DerivativeForm<1, dim, spacedim> &
5410 FEValuesBase<dim, spacedim>::jacobian(const unsigned int i) const
5411 {
5412  Assert(this->update_flags & update_jacobians,
5413  ExcAccessToUninitializedField("update_jacobians"));
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"));
5418 
5419  return this->mapping_output.jacobians[i];
5420 }
5421 
5422 
5423 
5424 template <int dim, int spacedim>
5425 inline const DerivativeForm<2, dim, spacedim> &
5426 FEValuesBase<dim, spacedim>::jacobian_grad(const unsigned int i) const
5427 {
5428  Assert(this->update_flags & update_jacobian_grads,
5429  ExcAccessToUninitializedField("update_jacobians_grads"));
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"));
5434 
5435  return this->mapping_output.jacobian_grads[i];
5436 }
5437 
5438 
5439 
5440 template <int dim, int spacedim>
5441 inline const DerivativeForm<1, spacedim, dim> &
5442 FEValuesBase<dim, spacedim>::inverse_jacobian(const unsigned int i) const
5443 {
5444  Assert(this->update_flags & update_inverse_jacobians,
5445  ExcAccessToUninitializedField("update_inverse_jacobians"));
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"));
5450 
5451  return this->mapping_output.inverse_jacobians[i];
5452 }
5453 
5454 
5455 
5456 template <int dim, int spacedim>
5457 inline const Tensor<1, spacedim> &
5458 FEValuesBase<dim, spacedim>::normal_vector(const unsigned int i) const
5459 {
5460  Assert(this->update_flags & update_normal_vectors,
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"));
5467 
5468  return this->mapping_output.normal_vectors[i];
5469 }
5470 
5471 
5472 
5473 /*--------------------- Inline functions: FEValues --------------------------*/
5474 
5475 
5476 template <int dim, int spacedim>
5477 inline const Quadrature<dim> &
5479 {
5480  return quadrature;
5481 }
5482 
5483 
5484 
5485 template <int dim, int spacedim>
5486 inline const FEValues<dim, spacedim> &
5488 {
5489  return *this;
5490 }
5491 
5492 
5493 /*---------------------- Inline functions: FEFaceValuesBase -----------------*/
5494 
5495 
5496 template <int dim, int spacedim>
5497 inline unsigned int
5499 {
5500  return present_face_index;
5501 }
5502 
5503 
5504 /*----------------------- Inline functions: FE*FaceValues -------------------*/
5505 
5506 template <int dim, int spacedim>
5507 inline const Quadrature<dim - 1> &
5509 {
5510  return quadrature;
5511 }
5512 
5513 
5514 
5515 template <int dim, int spacedim>
5516 inline const FEFaceValues<dim, spacedim> &
5518 {
5519  return *this;
5520 }
5521 
5522 
5523 
5524 template <int dim, int spacedim>
5525 inline const FESubfaceValues<dim, spacedim> &
5527 {
5528  return *this;
5529 }
5530 
5531 
5532 
5533 template <int dim, int spacedim>
5534 inline const Tensor<1, spacedim> &
5535 FEFaceValuesBase<dim, spacedim>::boundary_form(const unsigned int i) const
5536 {
5537  Assert(i < this->mapping_output.boundary_forms.size(),
5538  ExcIndexRange(i, 0, this->mapping_output.boundary_forms.size()));
5539  Assert(this->update_flags & update_boundary_forms,
5541  "update_boundary_forms")));
5542 
5543  return this->mapping_output.boundary_forms[i];
5544 }
5545 
5546 #endif // DOXYGEN
5547 
5548 DEAL_II_NAMESPACE_CLOSE
5549 
5550 #endif
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
Definition: fe_values.cc:1800
Transformed quadrature weights.
constexpr Tensor()=default
virtual ~FEValuesBase() override
Definition: fe_values.cc:3105
Shape function values.
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
Definition: fe_values.h:693
::SymmetricTensor< 2, spacedim > symmetric_gradient_type
Definition: fe_values.h:609
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)
Definition: fe_values.cc:4940
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
Definition: fe_values.h:3422
CellSimilarity::Similarity cell_similarity
Definition: fe_values.h:3454
typename ProductType< Number, typename Tensor< 2, dim, spacedim >::gradient_type >::type gradient_type
Definition: fe_values.h:1594
::Tensor< 1, spacedim > divergence_type
Definition: fe_values.h:1557
static const unsigned int dimension
Definition: fe_values.h:3835
Cache(const FEValuesBase< dim, spacedim > &fe_values)
Definition: fe_values.cc:2597
typename ProductType< Number, typename SymmetricTensor< 2, dim, spacedim >::value_type >::type value_type
Definition: fe_values.h:1291
unsigned int present_face_index
Definition: fe_values.h:3689
std::vector< ShapeFunctionData > shape_function_data
Definition: fe_values.h:1517
::Tensor< 2, spacedim > gradient_type
Definition: fe_values.h:597
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
Definition: fe_values.h:533
::Tensor< 2, spacedim > hessian_type
Definition: fe_values.h:164
static ::ExceptionBase & ExcFaceHasNoSubfaces()
const unsigned int dofs_per_cell
Definition: fe_values.h:2106
typename ::internal::CurlType< spacedim >::type curl_type
Definition: fe_values.h:624
const unsigned int component
Definition: fe_values.h:539
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
Definition: fe_values.cc:1998
void get_function_values(const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &values) const
Definition: fe_values.cc:1575
gradient_type gradient(const unsigned int shape_function, const unsigned int q_point) const
const Quadrature< dim - 1 > quadrature
Definition: fe_values.h:3694
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)
Definition: fe_values.cc:3082
Volume element.
void get_function_third_derivatives(const InputVector &fe_function, std::vector< typename ProductType< third_derivative_type, typename InputVector::value_type >::type > &third_derivatives) const
Definition: fe_values.cc:2255
static const unsigned int dimension
Definition: fe_values.h:3721
static ::ExceptionBase & ExcAccessToUninitializedField(std::string arg1)
const Mapping< dim, spacedim > & get_mapping() const
void get_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
Definition: fe_values.cc:3578
static ::ExceptionBase & ExcReinitCalledWithBoundaryFace()
Outer normal vector, not normalized.
typename ProductType< Number, typename Scalar< dim, spacedim >::hessian_type >::type hessian_type
Definition: fe_values.h:210
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)
Definition: fe_values.cc:4627
const FiniteElement< dim, spacedim > & get_fe() const
std::unique_ptr< const CellIteratorBase > present_cell
Definition: fe_values.h:3338
::Tensor< 4, spacedim > third_derivative_type
Definition: fe_values.h:638
::Tensor< 3, spacedim > gradient_type
Definition: fe_values.h:1563
static ::ExceptionBase & ExcFEDontMatch()
void get_function_divergences(const InputVector &fe_function, std::vector< typename ProductType< divergence_type, typename InputVector::value_type >::type > &divergences) const
Definition: fe_values.cc:2022
void get_function_laplacians(const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &laplacians) const
Definition: fe_values.cc:2191
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)
Definition: fe_values.cc:4814
std::size_t memory_consumption() const
Definition: fe_values.cc:4659
typename ProductType< Number, typename Scalar< dim, spacedim >::value_type >::type value_type
Definition: fe_values.h:186
unsigned int row_index[spacedim]
Definition: fe_values.h:737
void check_cell_similarity(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
Definition: fe_values.cc:4334
static const unsigned int space_dimension
Definition: fe_values.h:2094
const SmartPointer< const Mapping< dim, spacedim >, FEValuesBase< dim, spacedim > > mapping
Definition: fe_values.h:3390
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
Definition: fe_values.h:3469
void do_reinit(const unsigned int face_no, const unsigned int subface_no)
Definition: fe_values.cc:5013
static const unsigned int integral_dimension
Definition: fe_values.h:3846
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
Definition: fe_values.cc:2054
void get_function_curls_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::curl_type > &curls) const
Definition: fe_values.cc:2110
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
Definition: fe_values.h:661
static const unsigned int integral_dimension
Definition: fe_values.h:3729
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
Definition: fe_values.cc:1744
void get_function_hessians_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::hessian_type > &hessians) const
Definition: fe_values.cc:2166
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
Definition: fe_values.cc:1632
typename ProductType< Number, typename Scalar< dim, spacedim >::gradient_type >::type gradient_type
Definition: fe_values.h:194
typename internal::FEValuesViews::ViewType< dim, spacedim, Extractor >::type View
Definition: fe_values.h:1976
void get_function_hessians(const InputVector &fe_function, std::vector< typename ProductType< hessian_type, typename InputVector::value_type >::type > &hessians) const
Definition: fe_values.cc:2135
void get_function_laplacians(const InputVector &fe_function, std::vector< typename InputVector::value_type > &laplacians) const
Definition: fe_values.cc:3945
typename ProductType< Number, typename Tensor< 2, dim, spacedim >::value_type >::type value_type
Definition: fe_values.h:1578
static ::ExceptionBase & ExcMessage(std::string arg1)
::Tensor< 1, spacedim > gradient_type
Definition: fe_values.h:157
void get_function_gradients_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::gradient_type > &gradients) const
Definition: fe_values.cc:1942
Number value_type
Definition: vector.h:123
typename ProductType< Number, typename Vector< dim, spacedim >::third_derivative_type >::type third_derivative_type
Definition: fe_values.h:709
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:518
const double & shape_value(const unsigned int function_no, const unsigned int point_no) const
CellSimilarity::Similarity get_cell_similarity() const
Definition: fe_values.cc:4389
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
Definition: fe_values.cc:4250
std::vector< ShapeFunctionData > shape_function_data
Definition: fe_values.h:1864
void get_function_laplacians_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::laplacian_type > &laplacians) const
Definition: fe_values.cc:1775
unsigned int get_face_index() const
#define Assert(cond, exc)
Definition: exceptions.h:1407
UpdateFlags
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)
Definition: fe_values.cc:4855
Abstract base class for mapping classes.
Definition: dof_tools.h:57
Scalar & operator=(const Scalar< dim, spacedim > &)=delete
std::vector< ShapeFunctionData > shape_function_data
Definition: fe_values.h:1224
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:496
const Quadrature< dim > quadrature
Definition: fe_values.h:3587
const unsigned int first_vector_component
Definition: fe_values.h:1219
const std::vector< DerivativeForm< 1, dim, spacedim > > & get_jacobians() const
#define DeclException0(Exception0)
Definition: exceptions.h:473
void invalidate_present_cell()
Definition: fe_values.cc:4269
std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > mapping_data
Definition: fe_values.h:3398
const DerivativeForm< 4, dim, spacedim > & jacobian_3rd_derivative(const unsigned int quadrature_point) const
typename Tensor< rank_ - 1, dim, Number >::tensor_type value_type
Definition: tensor.h:427
static const unsigned int integral_dimension
Definition: fe_values.h:3505
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
Definition: fe_values.cc:1719
bool is_nonzero_shape_function_component[spacedim]
Definition: fe_values.h:726
void get_function_symmetric_gradients(const InputVector &fe_function, std::vector< typename ProductType< symmetric_gradient_type, typename InputVector::value_type >::type > &symmetric_gradients) const
Definition: fe_values.cc:1967
static const unsigned int integral_dimension
Definition: fe_values.h:3624
Second derivatives of shape functions.
static ::ExceptionBase & ExcFENotPrimitive()
Gradient of volume element.
const std::vector< Tensor< 1, spacedim > > & get_all_normal_vectors() const
Definition: fe_values.cc:4205
typename ProductType< Number, typename SymmetricTensor< 2, dim, spacedim >::divergence_type >::type divergence_type
Definition: fe_values.h:1299
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)
Definition: tensor.h:1416
divergence_type divergence(const unsigned int shape_function, const unsigned int q_point) const
std::vector<::FEValuesViews::Scalar< dim, spacedim > > scalars
Definition: fe_values.h:1951
std::size_t memory_consumption() const
Definition: fe_values.cc:4616
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access >> &cell)
Definition: fe_values.cc:4555
::SymmetricTensor< 2, spacedim > value_type
Definition: fe_values.h:1265
void get_function_hessians(const InputVector &fe_function, std::vector< typename ProductType< hessian_type, typename InputVector::value_type >::type > &hessians) const
Definition: fe_values.cc:1688
const Tensor< 1, spacedim > & boundary_form(const unsigned int i) const
const Quadrature< dim > & get_quadrature() const
const unsigned int n_quadrature_points
Definition: fe_values.h:2099
void get_function_hessians(const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type >> &hessians) const
Definition: fe_values.cc:3832
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)
Definition: fe_values.cc:4679
void get_function_curls(const InputVector &fe_function, std::vector< typename ProductType< curl_type, typename InputVector::value_type >::type > &curls) const
Definition: fe_values.cc:2079
void get_function_third_derivatives(const InputVector &fe_function, std::vector< Tensor< 3, spacedim, typename InputVector::value_type >> &third_derivatives) const
Definition: fe_values.cc:4082
static const unsigned int dimension
Definition: fe_values.h:2089
boost::signals2::connection tria_listener_mesh_transform
Definition: fe_values.h:3363
typename ProductType< Number, typename Vector< dim, spacedim >::value_type >::type value_type
Definition: fe_values.h:653
static const unsigned int space_dimension
Definition: fe_values.h:3840
void get_function_laplacians_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::laplacian_type > &laplacians) const
Definition: fe_values.cc:2227
const std::vector< Tensor< 3, spacedim > > & get_jacobian_pushed_forward_grads() const
Definition: mpi.h:90
const std::vector< double > & get_JxW_values() const
double JxW(const unsigned int quadrature_point) const
typename ProductType< Number, typename Vector< dim, spacedim >::symmetric_gradient_type >::type symmetric_gradient_type
Definition: fe_values.h:669
typename ProductType< Number, typename Vector< dim, spacedim >::divergence_type >::type divergence_type
Definition: fe_values.h:677
void initialize(const UpdateFlags update_flags)
Definition: fe_values.cc:4446
typename ProductType< Number, typename Vector< dim, spacedim >::hessian_type >::type hessian_type
Definition: fe_values.h:701
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
Definition: fe_values.cc:2286
Normal vectors.
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)
Definition: fe_values.cc:4287
typename ProductType< Number, typename Tensor< 2, dim, spacedim >::divergence_type >::type divergence_type
Definition: fe_values.h:1586
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
Definition: fe_values.h:1853
void get_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type >> &gradients) const
Definition: fe_values.cc:3719
Definition: fe.h:38
void initialize(const UpdateFlags update_flags)
Definition: fe_values.cc:4891
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
Definition: fe_values.h:1506
void get_function_gradients_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::gradient_type > &gradients) const
Definition: fe_values.cc:1663
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)
Definition: fe_values.cc:4762
const Triangulation< dim, spacedim >::cell_iterator get_cell() const
Definition: fe_values.cc:4196
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
Definition: fe_values.h:1213
void get_function_values(const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &values) const
Definition: fe_values.cc:1855
const Tensor< 1, spacedim > & shape_grad(const unsigned int function_no, const unsigned int quadrature_point) const
boost::signals2::connection tria_listener_refinement
Definition: fe_values.h:3354
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > finite_element_output
Definition: fe_values.h:3430
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
Definition: fe_values.h:685
const std::vector< Tensor< 1, spacedim > > & get_boundary_forms() const
Definition: fe_values.cc:4647
::Tensor< 3, spacedim > hessian_type
Definition: fe_values.h:631
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
void initialize(const UpdateFlags update_flags)
Definition: fe_values.cc:4715
third_derivative_type third_derivative(const unsigned int shape_function, const unsigned int q_point) const
::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > mapping_output
Definition: fe_values.h:3405
FEValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags update_flags)
Definition: fe_values.cc:4412
void get_function_gradients(const InputVector &fe_function, std::vector< typename ProductType< gradient_type, typename InputVector::value_type >::type > &gradients) const
Definition: fe_values.cc:1911
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
Definition: fe_values.cc:1831
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
Definition: fe_values.cc:4217
void do_reinit()
Definition: fe_values.cc:4582
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
Definition: fe_values.cc:1607
typename ProductType< Number, typename Scalar< dim, spacedim >::value_type >::type laplacian_type
Definition: fe_values.h:202
::Tensor< 3, spacedim > third_derivative_type
Definition: fe_values.h:171
typename ProductType< Number, typename Scalar< dim, spacedim >::third_derivative_type >::type third_derivative_type
Definition: fe_values.h:218
std::size_t memory_consumption() const
Definition: fe_values.cc:4230
void get_function_values_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::value_type > &values) const
Definition: fe_values.cc:1886
std::vector< ShapeFunctionData > shape_function_data
Definition: fe_values.h:544
UpdateFlags update_flags
Definition: fe_values.h:3436
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
Definition: fe_values.h:3414
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