Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
fe_values.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2020 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 
24 #include <deal.II/base/point.h>
29 
32 
33 #include <deal.II/fe/fe.h>
36 #include <deal.II/fe/mapping.h>
37 
38 #include <deal.II/grid/tria.h>
40 
41 #include <deal.II/hp/dof_handler.h>
42 
43 #include <algorithm>
44 #include <memory>
45 #include <type_traits>
46 
47 
48 // dummy include in order to have the
49 // definition of PetscScalar available
50 // without including other PETSc stuff
51 #ifdef DEAL_II_WITH_PETSC
52 # include <petsc.h>
53 #endif
54 
56 
57 // Forward declaration
58 #ifndef DOXYGEN
59 template <int dim, int spacedim = dim>
60 class FEValuesBase;
61 #endif
62 
63 namespace internal
64 {
69  template <int dim, class NumberType = double>
70  struct CurlType;
71 
78  template <class NumberType>
79  struct CurlType<1, NumberType>
80  {
82  };
83 
90  template <class NumberType>
91  struct CurlType<2, NumberType>
92  {
94  };
95 
102  template <class NumberType>
103  struct CurlType<3, NumberType>
104  {
106  };
107 } // namespace internal
108 
109 
110 
132 namespace FEValuesViews
133 {
145  template <int dim, int spacedim = dim>
146  class Scalar
147  {
148  public:
155 
162 
169 
176 
181  template <typename Number>
182  struct OutputType
183  {
188  using value_type =
189  typename ProductType<Number,
191 
196  using gradient_type = typename ProductType<
197  Number,
199 
204  using laplacian_type =
205  typename ProductType<Number,
207 
212  using hessian_type = typename ProductType<
213  Number,
215 
220  using third_derivative_type = typename ProductType<
221  Number,
223  };
224 
230  {
240 
249  unsigned int row_index;
250  };
251 
255  Scalar();
256 
262  Scalar(const FEValuesBase<dim, spacedim> &fe_values_base,
263  const unsigned int component);
264 
269  Scalar &
270  operator=(const Scalar<dim, spacedim> &) = delete;
271 
285  value_type
286  value(const unsigned int shape_function, const unsigned int q_point) const;
287 
299  gradient(const unsigned int shape_function,
300  const unsigned int q_point) const;
301 
313  hessian(const unsigned int shape_function,
314  const unsigned int q_point) const;
315 
327  third_derivative(const unsigned int shape_function,
328  const unsigned int q_point) const;
329 
347  template <class InputVector>
348  void
350  const InputVector &fe_function,
351  std::vector<typename ProductType<value_type,
352  typename InputVector::value_type>::type>
353  &values) const;
354 
377  template <class InputVector>
378  void
380  const InputVector &dof_values,
381  std::vector<
383  &values) const;
384 
402  template <class InputVector>
403  void
405  const InputVector &fe_function,
406  std::vector<typename ProductType<gradient_type,
407  typename InputVector::value_type>::type>
408  &gradients) const;
409 
413  template <class InputVector>
414  void
416  const InputVector &dof_values,
417  std::vector<
419  &gradients) const;
420 
438  template <class InputVector>
439  void
441  const InputVector &fe_function,
442  std::vector<typename ProductType<hessian_type,
443  typename InputVector::value_type>::type>
444  &hessians) const;
445 
449  template <class InputVector>
450  void
452  const InputVector &dof_values,
453  std::vector<
455  &hessians) const;
456 
457 
476  template <class InputVector>
477  void
479  const InputVector &fe_function,
480  std::vector<typename ProductType<value_type,
481  typename InputVector::value_type>::type>
482  &laplacians) const;
483 
487  template <class InputVector>
488  void
490  const InputVector &dof_values,
491  std::vector<
493  &laplacians) const;
494 
495 
514  template <class InputVector>
515  void
517  const InputVector &fe_function,
518  std::vector<typename ProductType<third_derivative_type,
519  typename InputVector::value_type>::type>
520  &third_derivatives) const;
521 
525  template <class InputVector>
526  void
528  const InputVector & dof_values,
529  std::vector<typename OutputType<typename InputVector::value_type>::
530  third_derivative_type> &third_derivatives) const;
531 
532 
533  private:
538 
543  const unsigned int component;
544 
548  std::vector<ShapeFunctionData> shape_function_data;
549  };
550 
551 
552 
582  template <int dim, int spacedim = dim>
583  class Vector
584  {
585  public:
592 
602 
614 
621 
628  using curl_type = typename ::internal::CurlType<spacedim>::type;
629 
636 
643 
648  template <typename Number>
649  struct OutputType
650  {
655  using value_type =
656  typename ProductType<Number,
658 
663  using gradient_type = typename ProductType<
664  Number,
666 
671  using symmetric_gradient_type = typename ProductType<
672  Number,
674 
679  using divergence_type = typename ProductType<
680  Number,
682 
687  using laplacian_type =
688  typename ProductType<Number,
690 
695  using curl_type =
696  typename ProductType<Number,
698 
703  using hessian_type = typename ProductType<
704  Number,
706 
711  using third_derivative_type = typename ProductType<
712  Number,
714  };
715 
721  {
731 
741  unsigned int row_index[spacedim];
742 
753  };
754 
758  Vector();
759 
768  Vector(const FEValuesBase<dim, spacedim> &fe_values_base,
769  const unsigned int first_vector_component);
770 
775  Vector &
776  operator=(const Vector<dim, spacedim> &) = delete;
777 
794  value_type
795  value(const unsigned int shape_function, const unsigned int q_point) const;
796 
811  gradient(const unsigned int shape_function,
812  const unsigned int q_point) const;
813 
830  symmetric_gradient(const unsigned int shape_function,
831  const unsigned int q_point) const;
832 
844  divergence(const unsigned int shape_function,
845  const unsigned int q_point) const;
846 
867  curl_type
868  curl(const unsigned int shape_function, const unsigned int q_point) const;
869 
881  hessian(const unsigned int shape_function,
882  const unsigned int q_point) const;
883 
895  third_derivative(const unsigned int shape_function,
896  const unsigned int q_point) const;
897 
915  template <class InputVector>
916  void
918  const InputVector &fe_function,
919  std::vector<typename ProductType<value_type,
920  typename InputVector::value_type>::type>
921  &values) const;
922 
945  template <class InputVector>
946  void
948  const InputVector &dof_values,
949  std::vector<
951  &values) const;
952 
970  template <class InputVector>
971  void
973  const InputVector &fe_function,
974  std::vector<typename ProductType<gradient_type,
975  typename InputVector::value_type>::type>
976  &gradients) const;
977 
981  template <class InputVector>
982  void
984  const InputVector &dof_values,
985  std::vector<
987  &gradients) const;
988 
1012  template <class InputVector>
1013  void
1015  const InputVector &fe_function,
1016  std::vector<typename ProductType<symmetric_gradient_type,
1017  typename InputVector::value_type>::type>
1018  &symmetric_gradients) const;
1019 
1023  template <class InputVector>
1024  void
1026  const InputVector & dof_values,
1027  std::vector<typename OutputType<typename InputVector::value_type>::
1028  symmetric_gradient_type> &symmetric_gradients) const;
1029 
1048  template <class InputVector>
1049  void
1051  const InputVector &fe_function,
1052  std::vector<typename ProductType<divergence_type,
1053  typename InputVector::value_type>::type>
1054  &divergences) const;
1055 
1059  template <class InputVector>
1060  void
1062  const InputVector &dof_values,
1063  std::vector<
1065  &divergences) const;
1066 
1085  template <class InputVector>
1086  void
1088  const InputVector &fe_function,
1089  std::vector<
1091  &curls) const;
1092 
1096  template <class InputVector>
1097  void
1099  const InputVector &dof_values,
1100  std::vector<
1102  &curls) const;
1103 
1121  template <class InputVector>
1122  void
1124  const InputVector &fe_function,
1125  std::vector<typename ProductType<hessian_type,
1126  typename InputVector::value_type>::type>
1127  &hessians) const;
1128 
1132  template <class InputVector>
1133  void
1135  const InputVector &dof_values,
1136  std::vector<
1138  &hessians) const;
1139 
1158  template <class InputVector>
1159  void
1161  const InputVector &fe_function,
1162  std::vector<typename ProductType<value_type,
1163  typename InputVector::value_type>::type>
1164  &laplacians) const;
1165 
1169  template <class InputVector>
1170  void
1172  const InputVector &dof_values,
1173  std::vector<
1175  &laplacians) const;
1176 
1195  template <class InputVector>
1196  void
1198  const InputVector &fe_function,
1199  std::vector<typename ProductType<third_derivative_type,
1200  typename InputVector::value_type>::type>
1201  &third_derivatives) const;
1202 
1206  template <class InputVector>
1207  void
1209  const InputVector & dof_values,
1210  std::vector<typename OutputType<typename InputVector::value_type>::
1211  third_derivative_type> &third_derivatives) const;
1212 
1213  private:
1218 
1223  const unsigned int first_vector_component;
1224 
1228  std::vector<ShapeFunctionData> shape_function_data;
1229  };
1230 
1231 
1232  template <int rank, int dim, int spacedim = dim>
1234 
1259  template <int dim, int spacedim>
1260  class SymmetricTensor<2, dim, spacedim>
1261  {
1262  public:
1270 
1281 
1286  template <typename Number>
1287  struct OutputType
1288  {
1293  using value_type = typename ProductType<
1294  Number,
1296 
1301  using divergence_type = typename ProductType<
1302  Number,
1304  };
1305 
1310  struct ShapeFunctionData
1311  {
1320  bool is_nonzero_shape_function_component
1321  [value_type::n_independent_components];
1322 
1332  unsigned int row_index[value_type::n_independent_components];
1333 
1343 
1348  };
1349 
1353  SymmetricTensor();
1354 
1364  SymmetricTensor(const FEValuesBase<dim, spacedim> &fe_values_base,
1365  const unsigned int first_tensor_component);
1366 
1371  SymmetricTensor &
1372  operator=(const SymmetricTensor<2, dim, spacedim> &) = delete;
1373 
1391  value_type
1392  value(const unsigned int shape_function, const unsigned int q_point) const;
1393 
1408  divergence(const unsigned int shape_function,
1409  const unsigned int q_point) const;
1410 
1428  template <class InputVector>
1429  void
1430  get_function_values(
1431  const InputVector &fe_function,
1432  std::vector<typename ProductType<value_type,
1433  typename InputVector::value_type>::type>
1434  &values) const;
1435 
1458  template <class InputVector>
1459  void
1460  get_function_values_from_local_dof_values(
1461  const InputVector &dof_values,
1462  std::vector<
1464  &values) const;
1465 
1487  template <class InputVector>
1488  void
1489  get_function_divergences(
1490  const InputVector &fe_function,
1491  std::vector<typename ProductType<divergence_type,
1492  typename InputVector::value_type>::type>
1493  &divergences) const;
1494 
1498  template <class InputVector>
1499  void
1500  get_function_divergences_from_local_dof_values(
1501  const InputVector &dof_values,
1502  std::vector<
1504  &divergences) const;
1505 
1506  private:
1511 
1516  const unsigned int first_tensor_component;
1517 
1521  std::vector<ShapeFunctionData> shape_function_data;
1522  };
1523 
1524 
1525  template <int rank, int dim, int spacedim = dim>
1526  class Tensor;
1527 
1548  template <int dim, int spacedim>
1549  class Tensor<2, dim, spacedim>
1550  {
1551  public:
1557 
1562 
1568 
1573  template <typename Number>
1574  struct OutputType
1575  {
1580  using value_type = typename ProductType<
1581  Number,
1583 
1588  using divergence_type = typename ProductType<
1589  Number,
1591 
1596  using gradient_type = typename ProductType<
1597  Number,
1599  };
1600 
1605  struct ShapeFunctionData
1606  {
1615  bool is_nonzero_shape_function_component
1616  [value_type::n_independent_components];
1617 
1627  unsigned int row_index[value_type::n_independent_components];
1628 
1638 
1643  };
1644 
1648  Tensor();
1649 
1659  Tensor(const FEValuesBase<dim, spacedim> &fe_values_base,
1660  const unsigned int first_tensor_component);
1661 
1662 
1667  Tensor &
1668  operator=(const Tensor<2, dim, spacedim> &) = delete;
1669 
1686  value_type
1687  value(const unsigned int shape_function, const unsigned int q_point) const;
1688 
1703  divergence(const unsigned int shape_function,
1704  const unsigned int q_point) const;
1705 
1720  gradient(const unsigned int shape_function,
1721  const unsigned int q_point) const;
1722 
1740  template <class InputVector>
1741  void
1742  get_function_values(
1743  const InputVector &fe_function,
1744  std::vector<typename ProductType<value_type,
1745  typename InputVector::value_type>::type>
1746  &values) const;
1747 
1770  template <class InputVector>
1771  void
1772  get_function_values_from_local_dof_values(
1773  const InputVector &dof_values,
1774  std::vector<
1776  &values) const;
1777 
1799  template <class InputVector>
1800  void
1801  get_function_divergences(
1802  const InputVector &fe_function,
1803  std::vector<typename ProductType<divergence_type,
1804  typename InputVector::value_type>::type>
1805  &divergences) const;
1806 
1810  template <class InputVector>
1811  void
1812  get_function_divergences_from_local_dof_values(
1813  const InputVector &dof_values,
1814  std::vector<
1816  &divergences) const;
1817 
1834  template <class InputVector>
1835  void
1836  get_function_gradients(
1837  const InputVector &fe_function,
1838  std::vector<typename ProductType<gradient_type,
1839  typename InputVector::value_type>::type>
1840  &gradients) const;
1841 
1845  template <class InputVector>
1846  void
1847  get_function_gradients_from_local_dof_values(
1848  const InputVector &dof_values,
1849  std::vector<
1851  &gradients) const;
1852 
1853  private:
1858 
1863  const unsigned int first_tensor_component;
1864 
1868  std::vector<ShapeFunctionData> shape_function_data;
1869  };
1870 
1871 } // namespace FEValuesViews
1872 
1873 
1874 namespace internal
1875 {
1876  namespace FEValuesViews
1877  {
1884  template <int dim, int spacedim, typename Extractor>
1885  struct ViewType
1886  {};
1887 
1895  template <int dim, int spacedim>
1896  struct ViewType<dim, spacedim, FEValuesExtractors::Scalar>
1897  {
1898  using type = typename ::FEValuesViews::Scalar<dim, spacedim>;
1899  };
1900 
1908  template <int dim, int spacedim>
1909  struct ViewType<dim, spacedim, FEValuesExtractors::Vector>
1910  {
1911  using type = typename ::FEValuesViews::Vector<dim, spacedim>;
1912  };
1913 
1921  template <int dim, int spacedim, int rank>
1922  struct ViewType<dim, spacedim, FEValuesExtractors::Tensor<rank>>
1923  {
1924  using type = typename ::FEValuesViews::Tensor<rank, dim, spacedim>;
1925  };
1926 
1934  template <int dim, int spacedim, int rank>
1935  struct ViewType<dim, spacedim, FEValuesExtractors::SymmetricTensor<rank>>
1936  {
1937  using type =
1938  typename ::FEValuesViews::SymmetricTensor<rank, dim, spacedim>;
1939  };
1940 
1948  template <int dim, int spacedim>
1949  struct Cache
1950  {
1955  std::vector<::FEValuesViews::Scalar<dim, spacedim>> scalars;
1956  std::vector<::FEValuesViews::Vector<dim, spacedim>> vectors;
1957  std::vector<::FEValuesViews::SymmetricTensor<2, dim, spacedim>>
1959  std::vector<::FEValuesViews::Tensor<2, dim, spacedim>>
1961 
1965  Cache(const FEValuesBase<dim, spacedim> &fe_values);
1966  };
1967  } // namespace FEValuesViews
1968 } // namespace internal
1969 
1970 namespace FEValuesViews
1971 {
1978  template <int dim, int spacedim, typename Extractor>
1979  using View = typename ::internal::FEValuesViews::
1980  ViewType<dim, spacedim, Extractor>::type;
1981 } // namespace FEValuesViews
1982 
1983 
2084 template <int dim, int spacedim>
2085 class FEValuesBase : public Subscriptor
2086 {
2087 public:
2091  static const unsigned int dimension = dim;
2092 
2096  static const unsigned int space_dimension = spacedim;
2097 
2101  const unsigned int n_quadrature_points;
2102 
2108  const unsigned int dofs_per_cell;
2109 
2110 
2118  FEValuesBase(const unsigned int n_q_points,
2119  const unsigned int dofs_per_cell,
2120  const UpdateFlags update_flags,
2123 
2128  FEValuesBase &
2129  operator=(const FEValuesBase &) = delete;
2130 
2135  FEValuesBase(const FEValuesBase &) = delete;
2136 
2140  virtual ~FEValuesBase() override;
2141 
2142 
2146 
2147 
2168  const double &
2169  shape_value(const unsigned int function_no,
2170  const unsigned int point_no) const;
2171 
2192  double
2193  shape_value_component(const unsigned int function_no,
2194  const unsigned int point_no,
2195  const unsigned int component) const;
2196 
2222  const Tensor<1, spacedim> &
2223  shape_grad(const unsigned int function_no,
2224  const unsigned int quadrature_point) const;
2225 
2243  shape_grad_component(const unsigned int function_no,
2244  const unsigned int point_no,
2245  const unsigned int component) const;
2246 
2266  const Tensor<2, spacedim> &
2267  shape_hessian(const unsigned int function_no,
2268  const unsigned int point_no) const;
2269 
2287  shape_hessian_component(const unsigned int function_no,
2288  const unsigned int point_no,
2289  const unsigned int component) const;
2290 
2310  const Tensor<3, spacedim> &
2311  shape_3rd_derivative(const unsigned int function_no,
2312  const unsigned int point_no) const;
2313 
2331  shape_3rd_derivative_component(const unsigned int function_no,
2332  const unsigned int point_no,
2333  const unsigned int component) const;
2334 
2336 
2338 
2375  template <class InputVector>
2376  void
2378  const InputVector & fe_function,
2379  std::vector<typename InputVector::value_type> &values) const;
2380 
2394  template <class InputVector>
2395  void
2397  const InputVector & fe_function,
2398  std::vector<Vector<typename InputVector::value_type>> &values) const;
2399 
2418  template <class InputVector>
2419  void
2421  const InputVector & fe_function,
2423  std::vector<typename InputVector::value_type> & values) const;
2424 
2446  template <class InputVector>
2447  void
2449  const InputVector & fe_function,
2451  std::vector<Vector<typename InputVector::value_type>> &values) const;
2452 
2453 
2484  template <class InputVector>
2485  void
2487  const InputVector & fe_function,
2489  ArrayView<std::vector<typename InputVector::value_type>> values,
2490  const bool quadrature_points_fastest) const;
2491 
2493 
2495 
2532  template <class InputVector>
2533  void
2535  const InputVector &fe_function,
2537  &gradients) const;
2538 
2555  template <class InputVector>
2556  void
2558  const InputVector &fe_function,
2559  std::vector<
2561  &gradients) const;
2562 
2569  template <class InputVector>
2570  void
2572  const InputVector & fe_function,
2575  &gradients) const;
2576 
2583  template <class InputVector>
2584  void
2586  const InputVector & fe_function,
2588  ArrayView<
2590  gradients,
2591  bool quadrature_points_fastest = false) const;
2592 
2594 
2598 
2636  template <class InputVector>
2637  void
2639  const InputVector &fe_function,
2641  &hessians) const;
2642 
2660  template <class InputVector>
2661  void
2663  const InputVector &fe_function,
2664  std::vector<
2666  & hessians,
2667  bool quadrature_points_fastest = false) const;
2668 
2673  template <class InputVector>
2674  void
2676  const InputVector & fe_function,
2679  &hessians) const;
2680 
2687  template <class InputVector>
2688  void
2690  const InputVector & fe_function,
2692  ArrayView<
2694  hessians,
2695  bool quadrature_points_fastest = false) const;
2696 
2737  template <class InputVector>
2738  void
2740  const InputVector & fe_function,
2741  std::vector<typename InputVector::value_type> &laplacians) const;
2742 
2762  template <class InputVector>
2763  void
2765  const InputVector & fe_function,
2766  std::vector<Vector<typename InputVector::value_type>> &laplacians) const;
2767 
2774  template <class InputVector>
2775  void
2777  const InputVector & fe_function,
2779  std::vector<typename InputVector::value_type> & laplacians) const;
2780 
2787  template <class InputVector>
2788  void
2790  const InputVector & fe_function,
2792  std::vector<Vector<typename InputVector::value_type>> &laplacians) const;
2793 
2800  template <class InputVector>
2801  void
2803  const InputVector & fe_function,
2805  std::vector<std::vector<typename InputVector::value_type>> &laplacians,
2806  bool quadrature_points_fastest = false) const;
2807 
2809 
2811 
2850  template <class InputVector>
2851  void
2853  const InputVector &fe_function,
2855  &third_derivatives) const;
2856 
2875  template <class InputVector>
2876  void
2878  const InputVector &fe_function,
2879  std::vector<
2881  & third_derivatives,
2882  bool quadrature_points_fastest = false) const;
2883 
2888  template <class InputVector>
2889  void
2891  const InputVector & fe_function,
2894  &third_derivatives) const;
2895 
2902  template <class InputVector>
2903  void
2905  const InputVector & fe_function,
2907  ArrayView<
2909  third_derivatives,
2910  bool quadrature_points_fastest = false) const;
2912 
2914 
2915 
2940  dof_indices() const;
2941 
2975  dof_indices_starting_at(const unsigned int start_dof_index) const;
2976 
3008  dof_indices_ending_at(const unsigned int end_dof_index) const;
3009 
3011 
3013 
3014 
3036  quadrature_point_indices() const;
3037 
3043  const Point<spacedim> &
3044  quadrature_point(const unsigned int q) const;
3045 
3051  const std::vector<Point<spacedim>> &
3052  get_quadrature_points() const;
3053 
3069  double
3070  JxW(const unsigned int quadrature_point) const;
3071 
3075  const std::vector<double> &
3076  get_JxW_values() const;
3077 
3085  jacobian(const unsigned int quadrature_point) const;
3086 
3093  const std::vector<DerivativeForm<1, dim, spacedim>> &
3094  get_jacobians() const;
3095 
3104  jacobian_grad(const unsigned int quadrature_point) const;
3105 
3112  const std::vector<DerivativeForm<2, dim, spacedim>> &
3113  get_jacobian_grads() const;
3114 
3123  const Tensor<3, spacedim> &
3124  jacobian_pushed_forward_grad(const unsigned int quadrature_point) const;
3125 
3132  const std::vector<Tensor<3, spacedim>> &
3134 
3143  jacobian_2nd_derivative(const unsigned int quadrature_point) const;
3144 
3151  const std::vector<DerivativeForm<3, dim, spacedim>> &
3153 
3163  const Tensor<4, spacedim> &
3165  const unsigned int quadrature_point) const;
3166 
3173  const std::vector<Tensor<4, spacedim>> &
3175 
3185  jacobian_3rd_derivative(const unsigned int quadrature_point) const;
3186 
3193  const std::vector<DerivativeForm<4, dim, spacedim>> &
3195 
3205  const Tensor<5, spacedim> &
3207  const unsigned int quadrature_point) const;
3208 
3215  const std::vector<Tensor<5, spacedim>> &
3217 
3225  inverse_jacobian(const unsigned int quadrature_point) const;
3226 
3233  const std::vector<DerivativeForm<1, spacedim, dim>> &
3234  get_inverse_jacobians() const;
3235 
3255  const Tensor<1, spacedim> &
3256  normal_vector(const unsigned int i) const;
3257 
3265  const std::vector<Tensor<1, spacedim>> &
3266  get_normal_vectors() const;
3267 
3269 
3271 
3272 
3282  operator[](const FEValuesExtractors::Scalar &scalar) const;
3283 
3293  operator[](const FEValuesExtractors::Vector &vector) const;
3294 
3306 
3307 
3317  operator[](const FEValuesExtractors::Tensor<2> &tensor) const;
3318 
3320 
3322 
3323 
3327  const Mapping<dim, spacedim> &
3328  get_mapping() const;
3329 
3334  get_fe() const;
3335 
3339  UpdateFlags
3340  get_update_flags() const;
3341 
3346  get_cell() const;
3347 
3354  get_cell_similarity() const;
3355 
3360  std::size_t
3361  memory_consumption() const;
3363 
3364 
3373  std::string,
3374  << "You are requesting information from an FEValues/FEFaceValues/FESubfaceValues "
3375  << "object for which this kind of information has not been computed. What "
3376  << "information these objects compute is determined by the update_* flags you "
3377  << "pass to the constructor. Here, the operation you are attempting requires "
3378  << "the <" << arg1
3379  << "> flag to be set, but it was apparently not specified "
3380  << "upon construction.");
3381 
3390  "The FiniteElement you provided to FEValues and the FiniteElement that belongs "
3391  "to the DoFHandler that provided the cell iterator do not match.");
3398  int,
3399  << "The shape function with index " << arg1
3400  << " is not primitive, i.e. it is vector-valued and "
3401  << "has more than one non-zero vector component. This "
3402  << "function cannot be called for these shape functions. "
3403  << "Maybe you want to use the same function with the "
3404  << "_component suffix?");
3405 
3414  "The given FiniteElement is not a primitive element but the requested operation "
3415  "only works for those. See FiniteElement::is_primitive() for more information.");
3416 
3417 protected:
3450  class CellIteratorBase;
3451 
3456  template <typename CI>
3458  class TriaCellIterator;
3459 
3465  std::unique_ptr<const CellIteratorBase> present_cell;
3466 
3474  boost::signals2::connection tria_listener_refinement;
3475 
3483  boost::signals2::connection tria_listener_mesh_transform;
3484 
3490  void
3492 
3502  void
3504  const typename Triangulation<dim, spacedim>::cell_iterator &cell);
3505 
3511 
3517  std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
3519 
3526 
3527 
3535 
3541  std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
3543 
3549  spacedim>
3551 
3552 
3557 
3566  UpdateFlags
3568 
3575 
3581  void
3583  const typename Triangulation<dim, spacedim>::cell_iterator &cell);
3584 
3585 private:
3590 
3591  // Make the view classes friends of this class, since they access internal
3592  // data.
3593  template <int, int>
3595  template <int, int>
3597  template <int, int, int>
3599  template <int, int, int>
3601 };
3602 
3603 
3604 
3615 template <int dim, int spacedim = dim>
3616 class FEValues : public FEValuesBase<dim, spacedim>
3617 {
3618 public:
3623  static const unsigned int integral_dimension = dim;
3624 
3629  FEValues(const Mapping<dim, spacedim> & mapping,
3630  const FiniteElement<dim, spacedim> &fe,
3631  const Quadrature<dim> & quadrature,
3632  const UpdateFlags update_flags);
3633 
3640  const Quadrature<dim> & quadrature,
3641  const UpdateFlags update_flags);
3642 
3649  template <template <int, int> class DoFHandlerType, bool level_dof_access>
3650  void
3651  reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
3652  level_dof_access>> &cell);
3653 
3667  void
3669 
3674  const Quadrature<dim> &
3675  get_quadrature() const;
3676 
3681  std::size_t
3682  memory_consumption() const;
3683 
3698  const FEValues<dim, spacedim> &
3699  get_present_fe_values() const;
3700 
3701 private:
3706 
3710  void
3711  initialize(const UpdateFlags update_flags);
3712 
3719  void
3720  do_reinit();
3721 };
3722 
3723 
3734 template <int dim, int spacedim = dim>
3735 class FEFaceValuesBase : public FEValuesBase<dim, spacedim>
3736 {
3737 public:
3742  static const unsigned int integral_dimension = dim - 1;
3743 
3755  FEFaceValuesBase(const unsigned int n_q_points,
3756  const unsigned int dofs_per_cell,
3757  const UpdateFlags update_flags,
3761 
3769  const Tensor<1, spacedim> &
3770  boundary_form(const unsigned int i) const;
3771 
3778  const std::vector<Tensor<1, spacedim>> &
3779  get_boundary_forms() const;
3780 
3785  unsigned int
3786  get_face_index() const;
3787 
3792  const Quadrature<dim - 1> &
3793  get_quadrature() const;
3794 
3799  std::size_t
3800  memory_consumption() const;
3801 
3802 protected:
3807  unsigned int present_face_index;
3808 
3812  const Quadrature<dim - 1> quadrature;
3813 };
3814 
3815 
3816 
3831 template <int dim, int spacedim = dim>
3832 class FEFaceValues : public FEFaceValuesBase<dim, spacedim>
3833 {
3834 public:
3839  static const unsigned int dimension = dim;
3840 
3841  static const unsigned int space_dimension = spacedim;
3842 
3847  static const unsigned int integral_dimension = dim - 1;
3848 
3853  FEFaceValues(const Mapping<dim, spacedim> & mapping,
3854  const FiniteElement<dim, spacedim> &fe,
3855  const Quadrature<dim - 1> & quadrature,
3856  const UpdateFlags update_flags);
3857 
3864  const Quadrature<dim - 1> & quadrature,
3865  const UpdateFlags update_flags);
3866 
3871  template <template <int, int> class DoFHandlerType, bool level_dof_access>
3872  void
3873  reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
3874  level_dof_access>> &cell,
3875  const unsigned int face_no);
3876 
3883  template <template <int, int> class DoFHandlerType, bool level_dof_access>
3884  void
3885  reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
3886  level_dof_access>> & cell,
3887  const typename Triangulation<dim, spacedim>::face_iterator &face);
3888 
3902  void
3904  const unsigned int face_no);
3905 
3906  /*
3907  * Reinitialize the gradients, Jacobi determinants, etc for the given face
3908  * on a given cell of type "iterator into a Triangulation object", and the
3909  * given finite element. Since iterators into a triangulation alone only
3910  * convey information about the geometry of a cell, but not about degrees of
3911  * freedom possibly associated with this cell, you will not be able to call
3912  * some functions of this class if they need information about degrees of
3913  * freedom. These functions are, above all, the
3914  * <tt>get_function_value/gradients/hessians/third_derivatives</tt>
3915  * functions. If you want to call these functions, you have to call the @p
3916  * reinit variants that take iterators into DoFHandler or other DoF handler
3917  * type objects.
3918  *
3919  * @note @p face must be one of @p cell's face iterators.
3920  */
3921  void
3923  const typename Triangulation<dim, spacedim>::face_iterator &face);
3924 
3940  get_present_fe_values() const;
3941 
3942 private:
3946  void
3947  initialize(const UpdateFlags update_flags);
3948 
3955  void
3956  do_reinit(const unsigned int face_no);
3957 };
3958 
3959 
3977 template <int dim, int spacedim = dim>
3978 class FESubfaceValues : public FEFaceValuesBase<dim, spacedim>
3979 {
3980 public:
3984  static const unsigned int dimension = dim;
3985 
3989  static const unsigned int space_dimension = spacedim;
3990 
3995  static const unsigned int integral_dimension = dim - 1;
3996 
4001  FESubfaceValues(const Mapping<dim, spacedim> & mapping,
4002  const FiniteElement<dim, spacedim> &fe,
4003  const Quadrature<dim - 1> & face_quadrature,
4004  const UpdateFlags update_flags);
4005 
4012  const Quadrature<dim - 1> & face_quadrature,
4013  const UpdateFlags update_flags);
4014 
4021  template <template <int, int> class DoFHandlerType, bool level_dof_access>
4022  void
4023  reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
4024  level_dof_access>> &cell,
4025  const unsigned int face_no,
4026  const unsigned int subface_no);
4027 
4032  template <template <int, int> class DoFHandlerType, bool level_dof_access>
4033  void
4034  reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
4035  level_dof_access>> & cell,
4036  const typename Triangulation<dim, spacedim>::face_iterator &face,
4037  const typename Triangulation<dim, spacedim>::face_iterator &subface);
4038 
4052  void
4054  const unsigned int face_no,
4055  const unsigned int subface_no);
4056 
4076  void
4078  const typename Triangulation<dim, spacedim>::face_iterator &face,
4079  const typename Triangulation<dim, spacedim>::face_iterator &subface);
4080 
4096  get_present_fe_values() const;
4097 
4104 
4111 
4112 private:
4116  void
4117  initialize(const UpdateFlags update_flags);
4118 
4125  void
4126  do_reinit(const unsigned int face_no, const unsigned int subface_no);
4127 };
4128 
4129 
4130 #ifndef DOXYGEN
4131 
4132 
4133 /*------------------------ Inline functions: namespace FEValuesViews --------*/
4134 
4135 namespace FEValuesViews
4136 {
4137  template <int dim, int spacedim>
4138  inline typename Scalar<dim, spacedim>::value_type
4139  Scalar<dim, spacedim>::value(const unsigned int shape_function,
4140  const unsigned int q_point) const
4141  {
4142  AssertIndexRange(shape_function, fe_values->fe->dofs_per_cell);
4143  Assert(
4144  fe_values->update_flags & update_values,
4146  "update_values"))));
4147 
4148  // an adaptation of the FEValuesBase::shape_value_component function
4149  // except that here we know the component as fixed and we have
4150  // pre-computed and cached a bunch of information. See the comments there.
4151  if (shape_function_data[shape_function].is_nonzero_shape_function_component)
4152  return fe_values->finite_element_output.shape_values(
4153  shape_function_data[shape_function].row_index, q_point);
4154  else
4155  return 0;
4156  }
4157 
4158 
4159 
4160  template <int dim, int spacedim>
4161  inline typename Scalar<dim, spacedim>::gradient_type
4162  Scalar<dim, spacedim>::gradient(const unsigned int shape_function,
4163  const unsigned int q_point) const
4164  {
4165  AssertIndexRange(shape_function, fe_values->fe->dofs_per_cell);
4166  Assert(fe_values->update_flags & update_gradients,
4168  "update_gradients")));
4169 
4170  // an adaptation of the FEValuesBase::shape_grad_component
4171  // function except that here we know the component as fixed and we have
4172  // pre-computed and cached a bunch of information. See the comments there.
4173  if (shape_function_data[shape_function].is_nonzero_shape_function_component)
4174  return fe_values->finite_element_output
4175  .shape_gradients[shape_function_data[shape_function].row_index]
4176  [q_point];
4177  else
4178  return gradient_type();
4179  }
4180 
4181 
4182 
4183  template <int dim, int spacedim>
4184  inline typename Scalar<dim, spacedim>::hessian_type
4185  Scalar<dim, spacedim>::hessian(const unsigned int shape_function,
4186  const unsigned int q_point) const
4187  {
4188  AssertIndexRange(shape_function, fe_values->fe->dofs_per_cell);
4189  Assert(fe_values->update_flags & update_hessians,
4191  "update_hessians")));
4192 
4193  // an adaptation of the FEValuesBase::shape_hessian_component
4194  // function except that here we know the component as fixed and we have
4195  // pre-computed and cached a bunch of information. See the comments there.
4196  if (shape_function_data[shape_function].is_nonzero_shape_function_component)
4197  return fe_values->finite_element_output
4198  .shape_hessians[shape_function_data[shape_function].row_index][q_point];
4199  else
4200  return hessian_type();
4201  }
4202 
4203 
4204 
4205  template <int dim, int spacedim>
4207  Scalar<dim, spacedim>::third_derivative(const unsigned int shape_function,
4208  const unsigned int q_point) const
4209  {
4210  AssertIndexRange(shape_function, fe_values->fe->dofs_per_cell);
4211  Assert(fe_values->update_flags & update_3rd_derivatives,
4213  "update_3rd_derivatives")));
4214 
4215  // an adaptation of the FEValuesBase::shape_3rdderivative_component
4216  // function except that here we know the component as fixed and we have
4217  // pre-computed and cached a bunch of information. See the comments there.
4218  if (shape_function_data[shape_function].is_nonzero_shape_function_component)
4219  return fe_values->finite_element_output
4220  .shape_3rd_derivatives[shape_function_data[shape_function].row_index]
4221  [q_point];
4222  else
4223  return third_derivative_type();
4224  }
4225 
4226 
4227 
4228  template <int dim, int spacedim>
4229  inline typename Vector<dim, spacedim>::value_type
4230  Vector<dim, spacedim>::value(const unsigned int shape_function,
4231  const unsigned int q_point) const
4232  {
4233  AssertIndexRange(shape_function, fe_values->fe->dofs_per_cell);
4234  Assert(fe_values->update_flags & update_values,
4236  "update_values")));
4237 
4238  // same as for the scalar case except that we have one more index
4239  const int snc =
4240  shape_function_data[shape_function].single_nonzero_component;
4241  if (snc == -2)
4242  return value_type();
4243  else if (snc != -1)
4244  {
4245  value_type return_value;
4246  return_value[shape_function_data[shape_function]
4247  .single_nonzero_component_index] =
4248  fe_values->finite_element_output.shape_values(snc, q_point);
4249  return return_value;
4250  }
4251  else
4252  {
4253  value_type return_value;
4254  for (unsigned int d = 0; d < dim; ++d)
4255  if (shape_function_data[shape_function]
4256  .is_nonzero_shape_function_component[d])
4257  return_value[d] = fe_values->finite_element_output.shape_values(
4258  shape_function_data[shape_function].row_index[d], q_point);
4259 
4260  return return_value;
4261  }
4262  }
4263 
4264 
4265 
4266  template <int dim, int spacedim>
4267  inline typename Vector<dim, spacedim>::gradient_type
4268  Vector<dim, spacedim>::gradient(const unsigned int shape_function,
4269  const unsigned int q_point) const
4270  {
4271  AssertIndexRange(shape_function, fe_values->fe->dofs_per_cell);
4272  Assert(fe_values->update_flags & update_gradients,
4274  "update_gradients")));
4275 
4276  // same as for the scalar case except that we have one more index
4277  const int snc =
4278  shape_function_data[shape_function].single_nonzero_component;
4279  if (snc == -2)
4280  return gradient_type();
4281  else if (snc != -1)
4282  {
4283  gradient_type return_value;
4284  return_value[shape_function_data[shape_function]
4285  .single_nonzero_component_index] =
4286  fe_values->finite_element_output.shape_gradients[snc][q_point];
4287  return return_value;
4288  }
4289  else
4290  {
4291  gradient_type return_value;
4292  for (unsigned int d = 0; d < dim; ++d)
4293  if (shape_function_data[shape_function]
4294  .is_nonzero_shape_function_component[d])
4295  return_value[d] =
4296  fe_values->finite_element_output.shape_gradients
4297  [shape_function_data[shape_function].row_index[d]][q_point];
4298 
4299  return return_value;
4300  }
4301  }
4302 
4303 
4304 
4305  template <int dim, int spacedim>
4307  Vector<dim, spacedim>::divergence(const unsigned int shape_function,
4308  const unsigned int q_point) const
4309  {
4310  // this function works like in the case above
4311  AssertIndexRange(shape_function, fe_values->fe->dofs_per_cell);
4312  Assert(fe_values->update_flags & update_gradients,
4314  "update_gradients")));
4315 
4316  // same as for the scalar case except that we have one more index
4317  const int snc =
4318  shape_function_data[shape_function].single_nonzero_component;
4319  if (snc == -2)
4320  return divergence_type();
4321  else if (snc != -1)
4322  return fe_values->finite_element_output
4323  .shape_gradients[snc][q_point][shape_function_data[shape_function]
4324  .single_nonzero_component_index];
4325  else
4326  {
4327  divergence_type return_value = 0;
4328  for (unsigned int d = 0; d < dim; ++d)
4329  if (shape_function_data[shape_function]
4330  .is_nonzero_shape_function_component[d])
4331  return_value +=
4332  fe_values->finite_element_output.shape_gradients
4333  [shape_function_data[shape_function].row_index[d]][q_point][d];
4334 
4335  return return_value;
4336  }
4337  }
4338 
4339 
4340 
4341  template <int dim, int spacedim>
4342  inline typename Vector<dim, spacedim>::curl_type
4343  Vector<dim, spacedim>::curl(const unsigned int shape_function,
4344  const unsigned int q_point) const
4345  {
4346  // this function works like in the case above
4347 
4348  AssertIndexRange(shape_function, fe_values->fe->dofs_per_cell);
4349  Assert(fe_values->update_flags & update_gradients,
4351  "update_gradients")));
4352  // same as for the scalar case except that we have one more index
4353  const int snc =
4354  shape_function_data[shape_function].single_nonzero_component;
4355 
4356  if (snc == -2)
4357  return curl_type();
4358 
4359  else
4360  switch (dim)
4361  {
4362  case 1:
4363  {
4364  Assert(false,
4365  ExcMessage(
4366  "Computing the curl in 1d is not a useful operation"));
4367  return curl_type();
4368  }
4369 
4370  case 2:
4371  {
4372  if (snc != -1)
4373  {
4374  curl_type return_value;
4375 
4376  // the single nonzero component can only be zero or one in 2d
4377  if (shape_function_data[shape_function]
4378  .single_nonzero_component_index == 0)
4379  return_value[0] =
4380  -1.0 * fe_values->finite_element_output
4381  .shape_gradients[snc][q_point][1];
4382  else
4383  return_value[0] = fe_values->finite_element_output
4384  .shape_gradients[snc][q_point][0];
4385 
4386  return return_value;
4387  }
4388 
4389  else
4390  {
4391  curl_type return_value;
4392 
4393  return_value[0] = 0.0;
4394 
4395  if (shape_function_data[shape_function]
4396  .is_nonzero_shape_function_component[0])
4397  return_value[0] -=
4398  fe_values->finite_element_output
4399  .shape_gradients[shape_function_data[shape_function]
4400  .row_index[0]][q_point][1];
4401 
4402  if (shape_function_data[shape_function]
4403  .is_nonzero_shape_function_component[1])
4404  return_value[0] +=
4405  fe_values->finite_element_output
4406  .shape_gradients[shape_function_data[shape_function]
4407  .row_index[1]][q_point][0];
4408 
4409  return return_value;
4410  }
4411  }
4412 
4413  case 3:
4414  {
4415  if (snc != -1)
4416  {
4417  curl_type return_value;
4418 
4419  switch (shape_function_data[shape_function]
4420  .single_nonzero_component_index)
4421  {
4422  case 0:
4423  {
4424  return_value[0] = 0;
4425  return_value[1] = fe_values->finite_element_output
4426  .shape_gradients[snc][q_point][2];
4427  return_value[2] =
4428  -1.0 * fe_values->finite_element_output
4429  .shape_gradients[snc][q_point][1];
4430  return return_value;
4431  }
4432 
4433  case 1:
4434  {
4435  return_value[0] =
4436  -1.0 * fe_values->finite_element_output
4437  .shape_gradients[snc][q_point][2];
4438  return_value[1] = 0;
4439  return_value[2] = fe_values->finite_element_output
4440  .shape_gradients[snc][q_point][0];
4441  return return_value;
4442  }
4443 
4444  default:
4445  {
4446  return_value[0] = fe_values->finite_element_output
4447  .shape_gradients[snc][q_point][1];
4448  return_value[1] =
4449  -1.0 * fe_values->finite_element_output
4450  .shape_gradients[snc][q_point][0];
4451  return_value[2] = 0;
4452  return return_value;
4453  }
4454  }
4455  }
4456 
4457  else
4458  {
4459  curl_type return_value;
4460 
4461  for (unsigned int i = 0; i < dim; ++i)
4462  return_value[i] = 0.0;
4463 
4464  if (shape_function_data[shape_function]
4465  .is_nonzero_shape_function_component[0])
4466  {
4467  return_value[1] +=
4468  fe_values->finite_element_output
4469  .shape_gradients[shape_function_data[shape_function]
4470  .row_index[0]][q_point][2];
4471  return_value[2] -=
4472  fe_values->finite_element_output
4473  .shape_gradients[shape_function_data[shape_function]
4474  .row_index[0]][q_point][1];
4475  }
4476 
4477  if (shape_function_data[shape_function]
4478  .is_nonzero_shape_function_component[1])
4479  {
4480  return_value[0] -=
4481  fe_values->finite_element_output
4482  .shape_gradients[shape_function_data[shape_function]
4483  .row_index[1]][q_point][2];
4484  return_value[2] +=
4485  fe_values->finite_element_output
4486  .shape_gradients[shape_function_data[shape_function]
4487  .row_index[1]][q_point][0];
4488  }
4489 
4490  if (shape_function_data[shape_function]
4491  .is_nonzero_shape_function_component[2])
4492  {
4493  return_value[0] +=
4494  fe_values->finite_element_output
4495  .shape_gradients[shape_function_data[shape_function]
4496  .row_index[2]][q_point][1];
4497  return_value[1] -=
4498  fe_values->finite_element_output
4499  .shape_gradients[shape_function_data[shape_function]
4500  .row_index[2]][q_point][0];
4501  }
4502 
4503  return return_value;
4504  }
4505  }
4506  }
4507  // should not end up here
4508  Assert(false, ExcInternalError());
4509  return curl_type();
4510  }
4511 
4512 
4513 
4514  template <int dim, int spacedim>
4515  inline typename Vector<dim, spacedim>::hessian_type
4516  Vector<dim, spacedim>::hessian(const unsigned int shape_function,
4517  const unsigned int q_point) const
4518  {
4519  // this function works like in the case above
4520  AssertIndexRange(shape_function, fe_values->fe->dofs_per_cell);
4521  Assert(fe_values->update_flags & update_hessians,
4523  "update_hessians")));
4524 
4525  // same as for the scalar case except that we have one more index
4526  const int snc =
4527  shape_function_data[shape_function].single_nonzero_component;
4528  if (snc == -2)
4529  return hessian_type();
4530  else if (snc != -1)
4531  {
4532  hessian_type return_value;
4533  return_value[shape_function_data[shape_function]
4534  .single_nonzero_component_index] =
4535  fe_values->finite_element_output.shape_hessians[snc][q_point];
4536  return return_value;
4537  }
4538  else
4539  {
4540  hessian_type return_value;
4541  for (unsigned int d = 0; d < dim; ++d)
4542  if (shape_function_data[shape_function]
4543  .is_nonzero_shape_function_component[d])
4544  return_value[d] =
4545  fe_values->finite_element_output.shape_hessians
4546  [shape_function_data[shape_function].row_index[d]][q_point];
4547 
4548  return return_value;
4549  }
4550  }
4551 
4552 
4553 
4554  template <int dim, int spacedim>
4556  Vector<dim, spacedim>::third_derivative(const unsigned int shape_function,
4557  const unsigned int q_point) const
4558  {
4559  // this function works like in the case above
4560  AssertIndexRange(shape_function, fe_values->fe->dofs_per_cell);
4561  Assert(fe_values->update_flags & update_3rd_derivatives,
4563  "update_3rd_derivatives")));
4564 
4565  // same as for the scalar case except that we have one more index
4566  const int snc =
4567  shape_function_data[shape_function].single_nonzero_component;
4568  if (snc == -2)
4569  return third_derivative_type();
4570  else if (snc != -1)
4571  {
4572  third_derivative_type return_value;
4573  return_value[shape_function_data[shape_function]
4574  .single_nonzero_component_index] =
4575  fe_values->finite_element_output.shape_3rd_derivatives[snc][q_point];
4576  return return_value;
4577  }
4578  else
4579  {
4580  third_derivative_type return_value;
4581  for (unsigned int d = 0; d < dim; ++d)
4582  if (shape_function_data[shape_function]
4583  .is_nonzero_shape_function_component[d])
4584  return_value[d] =
4585  fe_values->finite_element_output.shape_3rd_derivatives
4586  [shape_function_data[shape_function].row_index[d]][q_point];
4587 
4588  return return_value;
4589  }
4590  }
4591 
4592 
4593 
4594  namespace internal
4595  {
4600  inline ::SymmetricTensor<2, 1>
4601  symmetrize_single_row(const unsigned int n, const ::Tensor<1, 1> &t)
4602  {
4603  AssertIndexRange(n, 1);
4604  (void)n;
4605 
4606  return {{t[0]}};
4607  }
4608 
4609 
4610 
4611  inline ::SymmetricTensor<2, 2>
4612  symmetrize_single_row(const unsigned int n, const ::Tensor<1, 2> &t)
4613  {
4614  switch (n)
4615  {
4616  case 0:
4617  {
4618  return {{t[0], 0, t[1] / 2}};
4619  }
4620  case 1:
4621  {
4622  return {{0, t[1], t[0] / 2}};
4623  }
4624  default:
4625  {
4626  AssertIndexRange(n, 2);
4627  return {};
4628  }
4629  }
4630  }
4631 
4632 
4633 
4634  inline ::SymmetricTensor<2, 3>
4635  symmetrize_single_row(const unsigned int n, const ::Tensor<1, 3> &t)
4636  {
4637  switch (n)
4638  {
4639  case 0:
4640  {
4641  return {{t[0], 0, 0, t[1] / 2, t[2] / 2, 0}};
4642  }
4643  case 1:
4644  {
4645  return {{0, t[1], 0, t[0] / 2, 0, t[2] / 2}};
4646  }
4647  case 2:
4648  {
4649  return {{0, 0, t[2], 0, t[0] / 2, t[1] / 2}};
4650  }
4651  default:
4652  {
4653  AssertIndexRange(n, 3);
4654  return {};
4655  }
4656  }
4657  }
4658  } // namespace internal
4659 
4660 
4661 
4662  template <int dim, int spacedim>
4664  Vector<dim, spacedim>::symmetric_gradient(const unsigned int shape_function,
4665  const unsigned int q_point) const
4666  {
4667  AssertIndexRange(shape_function, fe_values->fe->dofs_per_cell);
4668  Assert(fe_values->update_flags & update_gradients,
4670  "update_gradients")));
4671 
4672  // same as for the scalar case except that we have one more index
4673  const int snc =
4674  shape_function_data[shape_function].single_nonzero_component;
4675  if (snc == -2)
4676  return symmetric_gradient_type();
4677  else if (snc != -1)
4678  return internal::symmetrize_single_row(
4679  shape_function_data[shape_function].single_nonzero_component_index,
4680  fe_values->finite_element_output.shape_gradients[snc][q_point]);
4681  else
4682  {
4683  gradient_type return_value;
4684  for (unsigned int d = 0; d < dim; ++d)
4685  if (shape_function_data[shape_function]
4686  .is_nonzero_shape_function_component[d])
4687  return_value[d] =
4688  fe_values->finite_element_output.shape_gradients
4689  [shape_function_data[shape_function].row_index[d]][q_point];
4690 
4691  return symmetrize(return_value);
4692  }
4693  }
4694 
4695 
4696 
4697  template <int dim, int spacedim>
4699  SymmetricTensor<2, dim, spacedim>::value(const unsigned int shape_function,
4700  const unsigned int q_point) const
4701  {
4702  AssertIndexRange(shape_function, fe_values->fe->dofs_per_cell);
4703  Assert(fe_values->update_flags & update_values,
4705  "update_values")));
4706 
4707  // similar to the vector case where we have more then one index and we need
4708  // to convert between unrolled and component indexing for tensors
4709  const int snc =
4710  shape_function_data[shape_function].single_nonzero_component;
4711 
4712  if (snc == -2)
4713  {
4714  // shape function is zero for the selected components
4715  return value_type();
4716  }
4717  else if (snc != -1)
4718  {
4719  value_type return_value;
4720  const unsigned int comp =
4721  shape_function_data[shape_function].single_nonzero_component_index;
4722  return_value[value_type::unrolled_to_component_indices(comp)] =
4723  fe_values->finite_element_output.shape_values(snc, q_point);
4724  return return_value;
4725  }
4726  else
4727  {
4728  value_type return_value;
4729  for (unsigned int d = 0; d < value_type::n_independent_components; ++d)
4730  if (shape_function_data[shape_function]
4731  .is_nonzero_shape_function_component[d])
4732  return_value[value_type::unrolled_to_component_indices(d)] =
4733  fe_values->finite_element_output.shape_values(
4734  shape_function_data[shape_function].row_index[d], q_point);
4735  return return_value;
4736  }
4737  }
4738 
4739 
4740 
4741  template <int dim, int spacedim>
4744  const unsigned int shape_function,
4745  const unsigned int q_point) const
4746  {
4747  AssertIndexRange(shape_function, fe_values->fe->dofs_per_cell);
4748  Assert(fe_values->update_flags & update_gradients,
4750  "update_gradients")));
4751 
4752  const int snc =
4753  shape_function_data[shape_function].single_nonzero_component;
4754 
4755  if (snc == -2)
4756  {
4757  // shape function is zero for the selected components
4758  return divergence_type();
4759  }
4760  else if (snc != -1)
4761  {
4762  // we have a single non-zero component when the symmetric tensor is
4763  // represented in unrolled form. this implies we potentially have
4764  // two non-zero components when represented in component form! we
4765  // will only have one non-zero entry if the non-zero component lies on
4766  // the diagonal of the tensor.
4767  //
4768  // the divergence of a second-order tensor is a first order tensor.
4769  //
4770  // assume the second-order tensor is A with components A_{ij}. then
4771  // A_{ij} = A_{ji} and there is only one (if diagonal) or two non-zero
4772  // entries in the tensorial representation. define the
4773  // divergence as:
4774  // b_i \dealcoloneq \dfrac{\partial phi_{ij}}{\partial x_j}.
4775  // (which is incidentally also
4776  // b_j \dealcoloneq \dfrac{\partial phi_{ij}}{\partial x_i}).
4777  // In both cases, a sum is implied.
4778  //
4779  // Now, we know the nonzero component in unrolled form: it is indicated
4780  // by 'snc'. we can figure out which tensor components belong to this:
4781  const unsigned int comp =
4782  shape_function_data[shape_function].single_nonzero_component_index;
4783  const unsigned int ii =
4784  value_type::unrolled_to_component_indices(comp)[0];
4785  const unsigned int jj =
4786  value_type::unrolled_to_component_indices(comp)[1];
4787 
4788  // given the form of the divergence above, if ii=jj there is only a
4789  // single nonzero component of the full tensor and the gradient
4790  // equals
4791  // b_ii \dealcoloneq \dfrac{\partial phi_{ii,ii}}{\partial x_ii}.
4792  // all other entries of 'b' are zero
4793  //
4794  // on the other hand, if ii!=jj, then there are two nonzero entries in
4795  // the full tensor and
4796  // b_ii \dealcoloneq \dfrac{\partial phi_{ii,jj}}{\partial x_ii}.
4797  // b_jj \dealcoloneq \dfrac{\partial phi_{ii,jj}}{\partial x_jj}.
4798  // again, all other entries of 'b' are zero
4799  const ::Tensor<1, spacedim> &phi_grad =
4800  fe_values->finite_element_output.shape_gradients[snc][q_point];
4801 
4802  divergence_type return_value;
4803  return_value[ii] = phi_grad[jj];
4804 
4805  if (ii != jj)
4806  return_value[jj] = phi_grad[ii];
4807 
4808  return return_value;
4809  }
4810  else
4811  {
4812  Assert(false, ExcNotImplemented());
4813  divergence_type return_value;
4814  return return_value;
4815  }
4816  }
4817 
4818 
4819 
4820  template <int dim, int spacedim>
4821  inline typename Tensor<2, dim, spacedim>::value_type
4822  Tensor<2, dim, spacedim>::value(const unsigned int shape_function,
4823  const unsigned int q_point) const
4824  {
4825  AssertIndexRange(shape_function, fe_values->fe->dofs_per_cell);
4826  Assert(fe_values->update_flags & update_values,
4828  "update_values")));
4829 
4830  // similar to the vector case where we have more then one index and we need
4831  // to convert between unrolled and component indexing for tensors
4832  const int snc =
4833  shape_function_data[shape_function].single_nonzero_component;
4834 
4835  if (snc == -2)
4836  {
4837  // shape function is zero for the selected components
4838  return value_type();
4839  }
4840  else if (snc != -1)
4841  {
4842  value_type return_value;
4843  const unsigned int comp =
4844  shape_function_data[shape_function].single_nonzero_component_index;
4845  const TableIndices<2> indices =
4847  return_value[indices] =
4848  fe_values->finite_element_output.shape_values(snc, q_point);
4849  return return_value;
4850  }
4851  else
4852  {
4853  value_type return_value;
4854  for (unsigned int d = 0; d < dim * dim; ++d)
4855  if (shape_function_data[shape_function]
4856  .is_nonzero_shape_function_component[d])
4857  {
4858  const TableIndices<2> indices =
4860  return_value[indices] =
4861  fe_values->finite_element_output.shape_values(
4862  shape_function_data[shape_function].row_index[d], q_point);
4863  }
4864  return return_value;
4865  }
4866  }
4867 
4868 
4869 
4870  template <int dim, int spacedim>
4872  Tensor<2, dim, spacedim>::divergence(const unsigned int shape_function,
4873  const unsigned int q_point) const
4874  {
4875  AssertIndexRange(shape_function, fe_values->fe->dofs_per_cell);
4876  Assert(fe_values->update_flags & update_gradients,
4878  "update_gradients")));
4879 
4880  const int snc =
4881  shape_function_data[shape_function].single_nonzero_component;
4882 
4883  if (snc == -2)
4884  {
4885  // shape function is zero for the selected components
4886  return divergence_type();
4887  }
4888  else if (snc != -1)
4889  {
4890  // we have a single non-zero component when the tensor is
4891  // represented in unrolled form.
4892  //
4893  // the divergence of a second-order tensor is a first order tensor.
4894  //
4895  // assume the second-order tensor is A with components A_{ij},
4896  // then divergence is d_i := \frac{\partial A_{ij}}{\partial x_j}
4897  //
4898  // Now, we know the nonzero component in unrolled form: it is indicated
4899  // by 'snc'. we can figure out which tensor components belong to this:
4900  const unsigned int comp =
4901  shape_function_data[shape_function].single_nonzero_component_index;
4902  const TableIndices<2> indices =
4904  const unsigned int ii = indices[0];
4905  const unsigned int jj = indices[1];
4906 
4907  const ::Tensor<1, spacedim> &phi_grad =
4908  fe_values->finite_element_output.shape_gradients[snc][q_point];
4909 
4910  divergence_type return_value;
4911  // note that we contract \nabla from the right
4912  return_value[ii] = phi_grad[jj];
4913 
4914  return return_value;
4915  }
4916  else
4917  {
4918  Assert(false, ExcNotImplemented());
4919  divergence_type return_value;
4920  return return_value;
4921  }
4922  }
4923 
4924 
4925 
4926  template <int dim, int spacedim>
4928  Tensor<2, dim, spacedim>::gradient(const unsigned int shape_function,
4929  const unsigned int q_point) const
4930  {
4931  AssertIndexRange(shape_function, fe_values->fe->dofs_per_cell);
4932  Assert(fe_values->update_flags & update_gradients,
4934  "update_gradients")));
4935 
4936  const int snc =
4937  shape_function_data[shape_function].single_nonzero_component;
4938 
4939  if (snc == -2)
4940  {
4941  // shape function is zero for the selected components
4942  return gradient_type();
4943  }
4944  else if (snc != -1)
4945  {
4946  // we have a single non-zero component when the tensor is
4947  // represented in unrolled form.
4948  //
4949  // the gradient of a second-order tensor is a third order tensor.
4950  //
4951  // assume the second-order tensor is A with components A_{ij},
4952  // then gradient is B_{ijk} := \frac{\partial A_{ij}}{\partial x_k}
4953  //
4954  // Now, we know the nonzero component in unrolled form: it is indicated
4955  // by 'snc'. we can figure out which tensor components belong to this:
4956  const unsigned int comp =
4957  shape_function_data[shape_function].single_nonzero_component_index;
4958  const TableIndices<2> indices =
4960  const unsigned int ii = indices[0];
4961  const unsigned int jj = indices[1];
4962 
4963  const ::Tensor<1, spacedim> &phi_grad =
4964  fe_values->finite_element_output.shape_gradients[snc][q_point];
4965 
4966  gradient_type return_value;
4967  return_value[ii][jj] = phi_grad;
4968 
4969  return return_value;
4970  }
4971  else
4972  {
4973  Assert(false, ExcNotImplemented());
4974  gradient_type return_value;
4975  return return_value;
4976  }
4977  }
4978 
4979 } // namespace FEValuesViews
4980 
4981 
4982 
4983 /*---------------------- Inline functions: FEValuesBase ---------------------*/
4984 
4985 
4986 
4987 template <int dim, int spacedim>
4989  operator[](const FEValuesExtractors::Scalar &scalar) const
4990 {
4991  AssertIndexRange(scalar.component, fe_values_views_cache.scalars.size());
4992 
4993  return fe_values_views_cache.scalars[scalar.component];
4994 }
4995 
4996 
4997 
4998 template <int dim, int spacedim>
5000  operator[](const FEValuesExtractors::Vector &vector) const
5001 {
5003  fe_values_views_cache.vectors.size());
5004 
5005  return fe_values_views_cache.vectors[vector.first_vector_component];
5006 }
5007 
5008 
5009 
5010 template <int dim, int spacedim>
5014 {
5015  Assert(
5016  tensor.first_tensor_component <
5017  fe_values_views_cache.symmetric_second_order_tensors.size(),
5019  0,
5020  fe_values_views_cache.symmetric_second_order_tensors.size()));
5021 
5022  return fe_values_views_cache
5023  .symmetric_second_order_tensors[tensor.first_tensor_component];
5024 }
5025 
5026 
5027 
5028 template <int dim, int spacedim>
5031  operator[](const FEValuesExtractors::Tensor<2> &tensor) const
5032 {
5034  fe_values_views_cache.second_order_tensors.size());
5035 
5036  return fe_values_views_cache
5037  .second_order_tensors[tensor.first_tensor_component];
5038 }
5039 
5040 
5041 
5042 template <int dim, int spacedim>
5043 inline const double &
5044 FEValuesBase<dim, spacedim>::shape_value(const unsigned int i,
5045  const unsigned int j) const
5046 {
5047  AssertIndexRange(i, fe->dofs_per_cell);
5048  Assert(this->update_flags & update_values,
5049  ExcAccessToUninitializedField("update_values"));
5050  Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
5051  Assert(present_cell.get() != nullptr,
5052  ExcMessage("FEValues object is not reinit'ed to any cell"));
5053  // if the entire FE is primitive,
5054  // then we can take a short-cut:
5055  if (fe->is_primitive())
5056  return this->finite_element_output.shape_values(i, j);
5057  else
5058  {
5059  // otherwise, use the mapping
5060  // between shape function
5061  // numbers and rows. note that
5062  // by the assertions above, we
5063  // know that this particular
5064  // shape function is primitive,
5065  // so we can call
5066  // system_to_component_index
5067  const unsigned int row =
5068  this->finite_element_output
5069  .shape_function_to_row_table[i * fe->n_components() +
5070  fe->system_to_component_index(i).first];
5071  return this->finite_element_output.shape_values(row, j);
5072  }
5073 }
5074 
5075 
5076 
5077 template <int dim, int spacedim>
5078 inline double
5080  const unsigned int i,
5081  const unsigned int j,
5082  const unsigned int component) const
5083 {
5084  AssertIndexRange(i, fe->dofs_per_cell);
5085  Assert(this->update_flags & update_values,
5086  ExcAccessToUninitializedField("update_values"));
5087  AssertIndexRange(component, fe->n_components());
5088  Assert(present_cell.get() != nullptr,
5089  ExcMessage("FEValues object is not reinit'ed to any cell"));
5090 
5091  // check whether the shape function
5092  // is non-zero at all within
5093  // this component:
5094  if (fe->get_nonzero_components(i)[component] == false)
5095  return 0;
5096 
5097  // look up the right row in the
5098  // table and take the data from
5099  // there
5100  const unsigned int row =
5101  this->finite_element_output
5102  .shape_function_to_row_table[i * fe->n_components() + component];
5103  return this->finite_element_output.shape_values(row, j);
5104 }
5105 
5106 
5107 
5108 template <int dim, int spacedim>
5109 inline const Tensor<1, spacedim> &
5110 FEValuesBase<dim, spacedim>::shape_grad(const unsigned int i,
5111  const unsigned int j) const
5112 {
5113  AssertIndexRange(i, fe->dofs_per_cell);
5114  Assert(this->update_flags & update_gradients,
5115  ExcAccessToUninitializedField("update_gradients"));
5116  Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
5117  Assert(present_cell.get() != nullptr,
5118  ExcMessage("FEValues object is not reinit'ed to any cell"));
5119  // if the entire FE is primitive,
5120  // then we can take a short-cut:
5121  if (fe->is_primitive())
5122  return this->finite_element_output.shape_gradients[i][j];
5123  else
5124  {
5125  // otherwise, use the mapping
5126  // between shape function
5127  // numbers and rows. note that
5128  // by the assertions above, we
5129  // know that this particular
5130  // shape function is primitive,
5131  // so we can call
5132  // system_to_component_index
5133  const unsigned int row =
5134  this->finite_element_output
5135  .shape_function_to_row_table[i * fe->n_components() +
5136  fe->system_to_component_index(i).first];
5137  return this->finite_element_output.shape_gradients[row][j];
5138  }
5139 }
5140 
5141 
5142 
5143 template <int dim, int spacedim>
5144 inline Tensor<1, spacedim>
5146  const unsigned int i,
5147  const unsigned int j,
5148  const unsigned int component) const
5149 {
5150  AssertIndexRange(i, fe->dofs_per_cell);
5151  Assert(this->update_flags & update_gradients,
5152  ExcAccessToUninitializedField("update_gradients"));
5153  AssertIndexRange(component, fe->n_components());
5154  Assert(present_cell.get() != nullptr,
5155  ExcMessage("FEValues object is not reinit'ed to any cell"));
5156  // check whether the shape function
5157  // is non-zero at all within
5158  // this component:
5159  if (fe->get_nonzero_components(i)[component] == false)
5160  return Tensor<1, spacedim>();
5161 
5162  // look up the right row in the
5163  // table and take the data from
5164  // there
5165  const unsigned int row =
5166  this->finite_element_output
5167  .shape_function_to_row_table[i * fe->n_components() + component];
5168  return this->finite_element_output.shape_gradients[row][j];
5169 }
5170 
5171 
5172 
5173 template <int dim, int spacedim>
5174 inline const Tensor<2, spacedim> &
5175 FEValuesBase<dim, spacedim>::shape_hessian(const unsigned int i,
5176  const unsigned int j) const
5177 {
5178  AssertIndexRange(i, fe->dofs_per_cell);
5179  Assert(this->update_flags & update_hessians,
5180  ExcAccessToUninitializedField("update_hessians"));
5181  Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
5182  Assert(present_cell.get() != nullptr,
5183  ExcMessage("FEValues object is not reinit'ed to any cell"));
5184  // if the entire FE is primitive,
5185  // then we can take a short-cut:
5186  if (fe->is_primitive())
5187  return this->finite_element_output.shape_hessians[i][j];
5188  else
5189  {
5190  // otherwise, use the mapping
5191  // between shape function
5192  // numbers and rows. note that
5193  // by the assertions above, we
5194  // know that this particular
5195  // shape function is primitive,
5196  // so we can call
5197  // system_to_component_index
5198  const unsigned int row =
5199  this->finite_element_output
5200  .shape_function_to_row_table[i * fe->n_components() +
5201  fe->system_to_component_index(i).first];
5202  return this->finite_element_output.shape_hessians[row][j];
5203  }
5204 }
5205 
5206 
5207 
5208 template <int dim, int spacedim>
5209 inline Tensor<2, spacedim>
5211  const unsigned int i,
5212  const unsigned int j,
5213  const unsigned int component) const
5214 {
5215  AssertIndexRange(i, fe->dofs_per_cell);
5216  Assert(this->update_flags & update_hessians,
5217  ExcAccessToUninitializedField("update_hessians"));
5218  AssertIndexRange(component, fe->n_components());
5219  Assert(present_cell.get() != nullptr,
5220  ExcMessage("FEValues object is not reinit'ed to any cell"));
5221  // check whether the shape function
5222  // is non-zero at all within
5223  // this component:
5224  if (fe->get_nonzero_components(i)[component] == false)
5225  return Tensor<2, spacedim>();
5226 
5227  // look up the right row in the
5228  // table and take the data from
5229  // there
5230  const unsigned int row =
5231  this->finite_element_output
5232  .shape_function_to_row_table[i * fe->n_components() + component];
5233  return this->finite_element_output.shape_hessians[row][j];
5234 }
5235 
5236 
5237 
5238 template <int dim, int spacedim>
5239 inline const Tensor<3, spacedim> &
5241  const unsigned int j) const
5242 {
5243  AssertIndexRange(i, fe->dofs_per_cell);
5244  Assert(this->update_flags & update_3rd_derivatives,
5245  ExcAccessToUninitializedField("update_3rd_derivatives"));
5246  Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
5247  Assert(present_cell.get() != nullptr,
5248  ExcMessage("FEValues object is not reinit'ed to any cell"));
5249  // if the entire FE is primitive,
5250  // then we can take a short-cut:
5251  if (fe->is_primitive())
5252  return this->finite_element_output.shape_3rd_derivatives[i][j];
5253  else
5254  {
5255  // otherwise, use the mapping
5256  // between shape function
5257  // numbers and rows. note that
5258  // by the assertions above, we
5259  // know that this particular
5260  // shape function is primitive,
5261  // so we can call
5262  // system_to_component_index
5263  const unsigned int row =
5264  this->finite_element_output
5265  .shape_function_to_row_table[i * fe->n_components() +
5266  fe->system_to_component_index(i).first];
5267  return this->finite_element_output.shape_3rd_derivatives[row][j];
5268  }
5269 }
5270 
5271 
5272 
5273 template <int dim, int spacedim>
5274 inline Tensor<3, spacedim>
5276  const unsigned int i,
5277  const unsigned int j,
5278  const unsigned int component) const
5279 {
5280  AssertIndexRange(i, fe->dofs_per_cell);
5281  Assert(this->update_flags & update_3rd_derivatives,
5282  ExcAccessToUninitializedField("update_3rd_derivatives"));
5283  AssertIndexRange(component, fe->n_components());
5284  Assert(present_cell.get() != nullptr,
5285  ExcMessage("FEValues object is not reinit'ed to any cell"));
5286  // check whether the shape function
5287  // is non-zero at all within
5288  // this component:
5289  if (fe->get_nonzero_components(i)[component] == false)
5290  return Tensor<3, spacedim>();
5291 
5292  // look up the right row in the
5293  // table and take the data from
5294  // there
5295  const unsigned int row =
5296  this->finite_element_output
5297  .shape_function_to_row_table[i * fe->n_components() + component];
5298  return this->finite_element_output.shape_3rd_derivatives[row][j];
5299 }
5300 
5301 
5302 
5303 template <int dim, int spacedim>
5304 inline const FiniteElement<dim, spacedim> &
5306 {
5307  return *fe;
5308 }
5309 
5310 
5311 
5312 template <int dim, int spacedim>
5313 inline const Mapping<dim, spacedim> &
5315 {
5316  return *mapping;
5317 }
5318 
5319 
5320 
5321 template <int dim, int spacedim>
5322 inline UpdateFlags
5324 {
5325  return this->update_flags;
5326 }
5327 
5328 
5329 
5330 template <int dim, int spacedim>
5331 inline const std::vector<Point<spacedim>> &
5333 {
5334  Assert(this->update_flags & update_quadrature_points,
5335  ExcAccessToUninitializedField("update_quadrature_points"));
5336  Assert(present_cell.get() != nullptr,
5337  ExcMessage("FEValues object is not reinit'ed to any cell"));
5338  return this->mapping_output.quadrature_points;
5339 }
5340 
5341 
5342 
5343 template <int dim, int spacedim>
5344 inline const std::vector<double> &
5346 {
5347  Assert(this->update_flags & update_JxW_values,
5348  ExcAccessToUninitializedField("update_JxW_values"));
5349  Assert(present_cell.get() != nullptr,
5350  ExcMessage("FEValues object is not reinit'ed to any cell"));
5351  return this->mapping_output.JxW_values;
5352 }
5353 
5354 
5355 
5356 template <int dim, int spacedim>
5357 inline const std::vector<DerivativeForm<1, dim, spacedim>> &
5359 {
5360  Assert(this->update_flags & update_jacobians,
5361  ExcAccessToUninitializedField("update_jacobians"));
5362  Assert(present_cell.get() != nullptr,
5363  ExcMessage("FEValues object is not reinit'ed to any cell"));
5364  return this->mapping_output.jacobians;
5365 }
5366 
5367 
5368 
5369 template <int dim, int spacedim>
5370 inline const std::vector<DerivativeForm<2, dim, spacedim>> &
5372 {
5373  Assert(this->update_flags & update_jacobian_grads,
5374  ExcAccessToUninitializedField("update_jacobians_grads"));
5375  Assert(present_cell.get() != nullptr,
5376  ExcMessage("FEValues object is not reinit'ed to any cell"));
5377  return this->mapping_output.jacobian_grads;
5378 }
5379 
5380 
5381 
5382 template <int dim, int spacedim>
5383 inline const Tensor<3, spacedim> &
5385  const unsigned int i) const
5386 {
5387  Assert(this->update_flags & update_jacobian_pushed_forward_grads,
5388  ExcAccessToUninitializedField("update_jacobian_pushed_forward_grads"));
5389  Assert(present_cell.get() != nullptr,
5390  ExcMessage("FEValues object is not reinit'ed to any cell"));
5391  return this->mapping_output.jacobian_pushed_forward_grads[i];
5392 }
5393 
5394 
5395 
5396 template <int dim, int spacedim>
5397 inline const std::vector<Tensor<3, spacedim>> &
5399 {
5400  Assert(this->update_flags & update_jacobian_pushed_forward_grads,
5401  ExcAccessToUninitializedField("update_jacobian_pushed_forward_grads"));
5402  Assert(present_cell.get() != nullptr,
5403  ExcMessage("FEValues object is not reinit'ed to any cell"));
5404  return this->mapping_output.jacobian_pushed_forward_grads;
5405 }
5406 
5407 
5408 
5409 template <int dim, int spacedim>
5410 inline const DerivativeForm<3, dim, spacedim> &
5411 FEValuesBase<dim, spacedim>::jacobian_2nd_derivative(const unsigned int i) const
5412 {
5413  Assert(this->update_flags & update_jacobian_2nd_derivatives,
5414  ExcAccessToUninitializedField("update_jacobian_2nd_derivatives"));
5415  Assert(present_cell.get() != nullptr,
5416  ExcMessage("FEValues object is not reinit'ed to any cell"));
5417  return this->mapping_output.jacobian_2nd_derivatives[i];
5418 }
5419 
5420 
5421 
5422 template <int dim, int spacedim>
5423 inline const std::vector<DerivativeForm<3, dim, spacedim>> &
5425 {
5426  Assert(this->update_flags & update_jacobian_2nd_derivatives,
5427  ExcAccessToUninitializedField("update_jacobian_2nd_derivatives"));
5428  Assert(present_cell.get() != nullptr,
5429  ExcMessage("FEValues object is not reinit'ed to any cell"));
5430  return this->mapping_output.jacobian_2nd_derivatives;
5431 }
5432 
5433 
5434 
5435 template <int dim, int spacedim>
5436 inline const Tensor<4, spacedim> &
5438  const unsigned int i) const
5439 {
5442  "update_jacobian_pushed_forward_2nd_derivatives"));
5443  Assert(present_cell.get() != nullptr,
5444  ExcMessage("FEValues object is not reinit'ed to any cell"));
5445  return this->mapping_output.jacobian_pushed_forward_2nd_derivatives[i];
5446 }
5447 
5448 
5449 
5450 template <int dim, int spacedim>
5451 inline const std::vector<Tensor<4, spacedim>> &
5453 {
5456  "update_jacobian_pushed_forward_2nd_derivatives"));
5457  Assert(present_cell.get() != nullptr,
5458  ExcMessage("FEValues object is not reinit'ed to any cell"));
5459  return this->mapping_output.jacobian_pushed_forward_2nd_derivatives;
5460 }
5461 
5462 
5463 
5464 template <int dim, int spacedim>
5465 inline const DerivativeForm<4, dim, spacedim> &
5466 FEValuesBase<dim, spacedim>::jacobian_3rd_derivative(const unsigned int i) const
5467 {
5468  Assert(this->update_flags & update_jacobian_3rd_derivatives,
5469  ExcAccessToUninitializedField("update_jacobian_3rd_derivatives"));
5470  Assert(present_cell.get() != nullptr,
5471  ExcMessage("FEValues object is not reinit'ed to any cell"));
5472  return this->mapping_output.jacobian_3rd_derivatives[i];
5473 }
5474 
5475 
5476 
5477 template <int dim, int spacedim>
5478 inline const std::vector<DerivativeForm<4, dim, spacedim>> &
5480 {
5481  Assert(this->update_flags & update_jacobian_3rd_derivatives,
5482  ExcAccessToUninitializedField("update_jacobian_3rd_derivatives"));
5483  Assert(present_cell.get() != nullptr,
5484  ExcMessage("FEValues object is not reinit'ed to any cell"));
5485  return this->mapping_output.jacobian_3rd_derivatives;
5486 }
5487 
5488 
5489 
5490 template <int dim, int spacedim>
5491 inline const Tensor<5, spacedim> &
5493  const unsigned int i) const
5494 {
5497  "update_jacobian_pushed_forward_3rd_derivatives"));
5498  Assert(present_cell.get() != nullptr,
5499  ExcMessage("FEValues object is not reinit'ed to any cell"));
5500  return this->mapping_output.jacobian_pushed_forward_3rd_derivatives[i];
5501 }
5502 
5503 
5504 
5505 template <int dim, int spacedim>
5506 inline const std::vector<Tensor<5, spacedim>> &
5508 {
5511  "update_jacobian_pushed_forward_3rd_derivatives"));
5512  Assert(present_cell.get() != nullptr,
5513  ExcMessage("FEValues object is not reinit'ed to any cell"));
5514  return this->mapping_output.jacobian_pushed_forward_3rd_derivatives;
5515 }
5516 
5517 
5518 
5519 template <int dim, int spacedim>
5520 inline const std::vector<DerivativeForm<1, spacedim, dim>> &
5522 {
5523  Assert(this->update_flags & update_inverse_jacobians,
5524  ExcAccessToUninitializedField("update_inverse_jacobians"));
5525  Assert(present_cell.get() != nullptr,
5526  ExcMessage("FEValues object is not reinit'ed to any cell"));
5527  return this->mapping_output.inverse_jacobians;
5528 }
5529 
5530 
5531 
5532 template <int dim, int spacedim>
5535 {
5536  return {0U, dofs_per_cell};
5537 }
5538 
5539 
5540 
5541 template <int dim, int spacedim>
5544  const unsigned int start_dof_index) const
5545 {
5546  Assert(start_dof_index <= dofs_per_cell,
5547  ExcIndexRange(start_dof_index, 0, dofs_per_cell + 1));
5548  return {start_dof_index, dofs_per_cell};
5549 }
5550 
5551 
5552 
5553 template <int dim, int spacedim>
5556  const unsigned int end_dof_index) const
5557 {
5558  Assert(end_dof_index < dofs_per_cell,
5559  ExcIndexRange(end_dof_index, 0, dofs_per_cell));
5560  return {0U, end_dof_index + 1};
5561 }
5562 
5563 
5564 
5565 template <int dim, int spacedim>
5568 {
5569  return {0U, n_quadrature_points};
5570 }
5571 
5572 
5573 
5574 template <int dim, int spacedim>
5575 inline const Point<spacedim> &
5576 FEValuesBase<dim, spacedim>::quadrature_point(const unsigned int i) const
5577 {
5578  Assert(this->update_flags & update_quadrature_points,
5579  ExcAccessToUninitializedField("update_quadrature_points"));
5580  AssertIndexRange(i, this->mapping_output.quadrature_points.size());
5581  Assert(present_cell.get() != nullptr,
5582  ExcMessage("FEValues object is not reinit'ed to any cell"));
5583 
5584  return this->mapping_output.quadrature_points[i];
5585 }
5586 
5587 
5588 
5589 template <int dim, int spacedim>
5590 inline double
5591 FEValuesBase<dim, spacedim>::JxW(const unsigned int i) const
5592 {
5593  Assert(this->update_flags & update_JxW_values,
5594  ExcAccessToUninitializedField("update_JxW_values"));
5595  AssertIndexRange(i, this->mapping_output.JxW_values.size());
5596  Assert(present_cell.get() != nullptr,
5597  ExcMessage("FEValues object is not reinit'ed to any cell"));
5598 
5599  return this->mapping_output.JxW_values[i];
5600 }
5601 
5602 
5603 
5604 template <int dim, int spacedim>
5605 inline const DerivativeForm<1, dim, spacedim> &
5606 FEValuesBase<dim, spacedim>::jacobian(const unsigned int i) const
5607 {
5608  Assert(this->update_flags & update_jacobians,
5609  ExcAccessToUninitializedField("update_jacobians"));
5610  AssertIndexRange(i, this->mapping_output.jacobians.size());
5611  Assert(present_cell.get() != nullptr,
5612  ExcMessage("FEValues object is not reinit'ed to any cell"));
5613 
5614  return this->mapping_output.jacobians[i];
5615 }
5616 
5617 
5618 
5619 template <int dim, int spacedim>
5620 inline const DerivativeForm<2, dim, spacedim> &
5621 FEValuesBase<dim, spacedim>::jacobian_grad(const unsigned int i) const
5622 {
5623  Assert(this->update_flags & update_jacobian_grads,
5624  ExcAccessToUninitializedField("update_jacobians_grads"));
5625  AssertIndexRange(i, this->mapping_output.jacobian_grads.size());
5626  Assert(present_cell.get() != nullptr,
5627  ExcMessage("FEValues object is not reinit'ed to any cell"));
5628 
5629  return this->mapping_output.jacobian_grads[i];
5630 }
5631 
5632 
5633 
5634 template <int dim, int spacedim>
5635 inline const DerivativeForm<1, spacedim, dim> &
5636 FEValuesBase<dim, spacedim>::inverse_jacobian(const unsigned int i) const
5637 {
5638  Assert(this->update_flags & update_inverse_jacobians,
5639  ExcAccessToUninitializedField("update_inverse_jacobians"));
5640  AssertIndexRange(i, this->mapping_output.inverse_jacobians.size());
5641  Assert(present_cell.get() != nullptr,
5642  ExcMessage("FEValues object is not reinit'ed to any cell"));
5643 
5644  return this->mapping_output.inverse_jacobians[i];
5645 }
5646 
5647 
5648 
5649 template <int dim, int spacedim>
5650 inline const Tensor<1, spacedim> &
5651 FEValuesBase<dim, spacedim>::normal_vector(const unsigned int i) const
5652 {
5653  Assert(this->update_flags & update_normal_vectors,
5655  "update_normal_vectors")));
5656  AssertIndexRange(i, this->mapping_output.normal_vectors.size());
5657  Assert(present_cell.get() != nullptr,
5658  ExcMessage("FEValues object is not reinit'ed to any cell"));
5659 
5660  return this->mapping_output.normal_vectors[i];
5661 }
5662 
5663 
5664 
5665 /*--------------------- Inline functions: FEValues --------------------------*/
5666 
5667 
5668 template <int dim, int spacedim>
5669 inline const Quadrature<dim> &
5671 {
5672  return quadrature;
5673 }
5674 
5675 
5676 
5677 template <int dim, int spacedim>
5678 inline const FEValues<dim, spacedim> &
5680 {
5681  return *this;
5682 }
5683 
5684 
5685 /*---------------------- Inline functions: FEFaceValuesBase -----------------*/
5686 
5687 
5688 template <int dim, int spacedim>
5689 inline unsigned int
5691 {
5692  return present_face_index;
5693 }
5694 
5695 
5696 /*----------------------- Inline functions: FE*FaceValues -------------------*/
5697 
5698 template <int dim, int spacedim>
5699 inline const Quadrature<dim - 1> &
5701 {
5702  return quadrature;
5703 }
5704 
5705 
5706 
5707 template <int dim, int spacedim>
5708 inline const FEFaceValues<dim, spacedim> &
5710 {
5711  return *this;
5712 }
5713 
5714 
5715 
5716 template <int dim, int spacedim>
5717 inline const FESubfaceValues<dim, spacedim> &
5719 {
5720  return *this;
5721 }
5722 
5723 
5724 
5725 template <int dim, int spacedim>
5726 inline const Tensor<1, spacedim> &
5727 FEFaceValuesBase<dim, spacedim>::boundary_form(const unsigned int i) const
5728 {
5729  AssertIndexRange(i, this->mapping_output.boundary_forms.size());
5730  Assert(this->update_flags & update_boundary_forms,
5732  "update_boundary_forms")));
5733 
5734  return this->mapping_output.boundary_forms[i];
5735 }
5736 
5737 #endif // DOXYGEN
5738 
5740 
5741 #endif
FEValuesBase::ExcFEDontMatch
static ::ExceptionBase & ExcFEDontMatch()
FEValuesExtractors::SymmetricTensor::first_tensor_component
unsigned int first_tensor_component
Definition: fe_values_extractors.h:204
FEValuesViews::Vector::fe_values
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
Definition: fe_values.h:1217
FEValuesBase::shape_value
const double & shape_value(const unsigned int function_no, const unsigned int point_no) const
internal::FEValuesViews::Cache::Cache
Cache(const FEValuesBase< dim, spacedim > &fe_values)
Definition: fe_values.cc:2591
FEValuesViews::Tensor< 2, dim, spacedim >::ShapeFunctionData::single_nonzero_component
int single_nonzero_component
Definition: fe_values.h:1637
FEValuesViews::Vector::OutputType::divergence_type
typename ProductType< Number, typename Vector< dim, spacedim >::divergence_type >::type divergence_type
Definition: fe_values.h:681
DeclExceptionMsg
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:496
FEValuesViews::Scalar::third_derivative_type
::Tensor< 3, spacedim > third_derivative_type
Definition: fe_values.h:175
FEValuesBase::shape_3rd_derivative_component
Tensor< 3, spacedim > shape_3rd_derivative_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
FEValues::memory_consumption
std::size_t memory_consumption() const
Definition: fe_values.cc:4594
derivative_form.h
iota_view.h
FEValuesExtractors
Definition: fe_values_extractors.h:82
FEValuesViews::Vector::OutputType::laplacian_type
typename ProductType< Number, typename Vector< dim, spacedim >::value_type >::type laplacian_type
Definition: fe_values.h:689
FEValuesViews::SymmetricTensor< 2, dim, spacedim >::ShapeFunctionData::single_nonzero_component_index
unsigned int single_nonzero_component_index
Definition: fe_values.h:1347
internal::FEValuesViews::Cache::vectors
std::vector<::FEValuesViews::Vector< dim, spacedim > > vectors
Definition: fe_values.h:1956
FEValuesViews::Vector::third_derivative
third_derivative_type third_derivative(const unsigned int shape_function, const unsigned int q_point) const
FEValuesViews::Scalar::OutputType
Definition: fe_values.h:182
internal::CurlType
Definition: fe_values.h:70
LinearAlgebraDealII::Vector
Vector< double > Vector
Definition: generic_linear_algebra.h:43
update_quadrature_points
@ update_quadrature_points
Transformed quadrature points.
Definition: fe_update_flags.h:122
update_jacobian_pushed_forward_3rd_derivatives
@ update_jacobian_pushed_forward_3rd_derivatives
Definition: fe_update_flags.h:215
FEValuesBase::get_jacobian_pushed_forward_3rd_derivatives
const std::vector< Tensor< 5, spacedim > > & get_jacobian_pushed_forward_3rd_derivatives() const
internal::ExcAccessToUninitializedField
static ::ExceptionBase & ExcAccessToUninitializedField()
FEValuesViews
Definition: fe_values.h:132
fe_update_flags.h
FEValuesBase::present_cell
std::unique_ptr< const CellIteratorBase > present_cell
Definition: fe_values.h:3458
FEValuesBase::tria_listener_refinement
boost::signals2::connection tria_listener_refinement
Definition: fe_values.h:3474
FEValuesViews::Scalar::hessian_type
::Tensor< 2, spacedim > hessian_type
Definition: fe_values.h:168
FEValuesBase::jacobian_3rd_derivative
const DerivativeForm< 4, dim, spacedim > & jacobian_3rd_derivative(const unsigned int quadrature_point) const
internal::FEValuesImplementation::FiniteElementRelatedData
Definition: fe_update_flags.h:524
TableIndices
Definition: table_indices.h:45
FEValuesBase::get_function_hessians
void get_function_hessians(const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type >> &hessians) const
Definition: fe_values.cc:3826
FEFaceValues::initialize
void initialize(const UpdateFlags update_flags)
Definition: fe_values.cc:4693
FEValuesViews::Tensor< 2, dim, spacedim >::OutputType::value_type
typename ProductType< Number, typename Tensor< 2, dim, spacedim >::value_type >::type value_type
Definition: fe_values.h:1582
FEValuesViews::Scalar
Definition: fe_values.h:146
FEValuesViews::Vector::ShapeFunctionData
Definition: fe_values.h:720
SymmetricTensor< 2, spacedim >
FEValuesViews::Scalar::component
const unsigned int component
Definition: fe_values.h:543
FEValuesViews::Scalar::ShapeFunctionData::is_nonzero_shape_function_component
bool is_nonzero_shape_function_component
Definition: fe_values.h:239
FEValuesViews::Scalar::gradient
gradient_type gradient(const unsigned int shape_function, const unsigned int q_point) const
FEValuesBase::quadrature_point_indices
std_cxx20::ranges::iota_view< unsigned int, unsigned int > quadrature_point_indices() const
FEValuesViews::Scalar::shape_function_data
std::vector< ShapeFunctionData > shape_function_data
Definition: fe_values.h:548
FEValuesViews::Vector
Definition: fe_values.h:583
FEValuesBase::get_normal_vectors
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors() const
Definition: fe_values.cc:4199
FEValuesViews::SymmetricTensor< 2, dim, spacedim >::OutputType::divergence_type
typename ProductType< Number, typename SymmetricTensor< 2, dim, spacedim >::divergence_type >::type divergence_type
Definition: fe_values.h:1303
FEValuesViews::SymmetricTensor< 2, dim, spacedim >::OutputType::value_type
typename ProductType< Number, typename SymmetricTensor< 2, dim, spacedim >::value_type >::type value_type
Definition: fe_values.h:1295
internal::FEValuesImplementation::MappingRelatedData
Definition: fe_update_flags.h:418
update_3rd_derivatives
@ update_3rd_derivatives
Third derivatives of shape functions.
Definition: fe_update_flags.h:96
FEValuesViews::Vector::ShapeFunctionData::single_nonzero_component
int single_nonzero_component
Definition: fe_values.h:751
FEValuesBase::ExcAccessToUninitializedField
static ::ExceptionBase & ExcAccessToUninitializedField(std::string arg1)
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
FEValuesViews::Vector::ShapeFunctionData::is_nonzero_shape_function_component
bool is_nonzero_shape_function_component[spacedim]
Definition: fe_values.h:730
FEFaceValues::integral_dimension
static const unsigned int integral_dimension
Definition: fe_values.h:3847
tria.h
FEValuesBase::~FEValuesBase
virtual ~FEValuesBase() override
Definition: fe_values.cc:3099
FEValuesBase::shape_value_component
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
ArrayView
Definition: array_view.h:77
FEValuesBase::get_jacobian_pushed_forward_grads
const std::vector< Tensor< 3, spacedim > > & get_jacobian_pushed_forward_grads() const
FEValuesViews::Scalar::fe_values
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
Definition: fe_values.h:537
FEValuesViews::Vector::third_derivative_type
::Tensor< 4, spacedim > third_derivative_type
Definition: fe_values.h:642
FEValuesViews::SymmetricTensor
Definition: fe_values.h:1233
FEValuesViews::Tensor
Definition: fe_values.h:1526
FEValuesBase::update_flags
UpdateFlags update_flags
Definition: fe_values.h:3556
FEValuesExtractors::Scalar
Definition: fe_values_extractors.h:95
FEValuesBase::tria_listener_mesh_transform
boost::signals2::connection tria_listener_mesh_transform
Definition: fe_values.h:3483
tria_iterator.h
FEValuesBase::dof_indices_starting_at
std_cxx20::ranges::iota_view< unsigned int, unsigned int > dof_indices_starting_at(const unsigned int start_dof_index) const
FEValuesBase::normal_vector
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const
FEFaceValuesBase::get_face_index
unsigned int get_face_index() const
FEValuesViews::Vector::get_function_curls_from_local_dof_values
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:2104
FEValues::FEValues
FEValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags update_flags)
Definition: fe_values.cc:4390
FEValuesViews::Vector::get_function_symmetric_gradients_from_local_dof_values
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:1992
FEValuesViews::Vector::curl_type
typename ::internal::CurlType< spacedim >::type curl_type
Definition: fe_values.h:628
FEValuesBase::cell_similarity
CellSimilarity::Similarity cell_similarity
Definition: fe_values.h:3574
LAPACKSupport::U
static const char U
Definition: lapack_support.h:167
FEValuesViews::Vector::get_function_hessians_from_local_dof_values
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:2160
FEValuesViews::Tensor< 2, dim, spacedim >
Definition: fe_values.h:1549
FEValuesViews::Tensor< 2, dim, spacedim >::OutputType::divergence_type
typename ProductType< Number, typename Tensor< 2, dim, spacedim >::divergence_type >::type divergence_type
Definition: fe_values.h:1590
FEValuesViews::Scalar::OutputType::gradient_type
typename ProductType< Number, typename Scalar< dim, spacedim >::gradient_type >::type gradient_type
Definition: fe_values.h:198
FEValuesViews::Vector::OutputType
Definition: fe_values.h:649
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
ProductType::type
typename internal::ProductTypeImpl< typename std::decay< T >::type, typename std::decay< U >::type >::type type
Definition: template_constraints.h:426
FEValuesViews::Vector::get_function_third_derivatives
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:2249
FEValuesBase::jacobian_pushed_forward_2nd_derivative
const Tensor< 4, spacedim > & jacobian_pushed_forward_2nd_derivative(const unsigned int quadrature_point) const
FEValuesViews::SymmetricTensor< 2, dim, spacedim >::shape_function_data
std::vector< ShapeFunctionData > shape_function_data
Definition: fe_values.h:1521
FEValuesBase::get_function_laplacians
void get_function_laplacians(const InputVector &fe_function, std::vector< typename InputVector::value_type > &laplacians) const
Definition: fe_values.cc:3939
FEValuesViews::Scalar::get_function_gradients
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:1626
FEValuesBase::ExcShapeFunctionNotPrimitive
static ::ExceptionBase & ExcShapeFunctionNotPrimitive(int arg1)
FEValuesViews::Vector::get_function_hessians
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:2129
FEValues::quadrature
const Quadrature< dim > quadrature
Definition: fe_values.h:3705
CellSimilarity::Similarity
Similarity
Definition: fe_update_flags.h:378
FEValues::get_present_fe_values
const FEValues< dim, spacedim > & get_present_fe_values() const
FESubfaceValues::reinit
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access >> &cell, const unsigned int face_no, const unsigned int subface_no)
FEValuesViews::Tensor< 2, dim, spacedim >::ShapeFunctionData::single_nonzero_component_index
unsigned int single_nonzero_component_index
Definition: fe_values.h:1642
internal::FEValuesViews::ViewType< dim, spacedim, FEValuesExtractors::SymmetricTensor< rank > >::type
typename ::FEValuesViews::SymmetricTensor< rank, dim, spacedim > type
Definition: fe_values.h:1938
petsc.h
FEValuesViews::Vector::get_function_gradients_from_local_dof_values
void get_function_gradients_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::gradient_type > &gradients) const
Definition: fe_values.cc:1936
FESubfaceValues::initialize
void initialize(const UpdateFlags update_flags)
Definition: fe_values.cc:4892
FEValuesViews::Vector::gradient_type
::Tensor< 2, spacedim > gradient_type
Definition: fe_values.h:601
FEValuesBase::ExcFENotPrimitive
static ::ExceptionBase & ExcFENotPrimitive()
FEValuesViews::Vector::get_function_divergences
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:2016
FEValuesViews::Scalar::OutputType::laplacian_type
typename ProductType< Number, typename Scalar< dim, spacedim >::value_type >::type laplacian_type
Definition: fe_values.h:206
update_values
@ update_values
Shape function values.
Definition: fe_update_flags.h:78
FEValuesBase::maybe_invalidate_previous_present_cell
void maybe_invalidate_previous_present_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
Definition: fe_values.cc:4269
FEValuesViews::Tensor< 2, dim, spacedim >::OutputType::gradient_type
typename ProductType< Number, typename Tensor< 2, dim, spacedim >::gradient_type >::type gradient_type
Definition: fe_values.h:1598
FEValuesViews::Scalar::get_function_values
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:1569
Physics::Elasticity::Kinematics::d
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
FEValuesViews::SymmetricTensor< 2, dim, spacedim >::first_tensor_component
const unsigned int first_tensor_component
Definition: fe_values.h:1516
FEValuesViews::Scalar::get_function_laplacians_from_local_dof_values
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:1769
Vector::value_type
Number value_type
Definition: vector.h:124
FEFaceValuesBase::quadrature
const Quadrature< dim - 1 > quadrature
Definition: fe_values.h:3812
FESubfaceValues::get_present_fe_values
const FESubfaceValues< dim, spacedim > & get_present_fe_values() const
FEValuesBase::get_function_gradients
void get_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type >> &gradients) const
Definition: fe_values.cc:3713
FEValuesBase::invalidate_present_cell
void invalidate_present_cell()
Definition: fe_values.cc:4251
FEValuesBase::mapping_output
::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > mapping_output
Definition: fe_values.h:3525
FEFaceValuesBase::integral_dimension
static const unsigned int integral_dimension
Definition: fe_values.h:3742
FEValuesBase::get_jacobian_grads
const std::vector< DerivativeForm< 2, dim, spacedim > > & get_jacobian_grads() const
update_jacobian_pushed_forward_grads
@ update_jacobian_pushed_forward_grads
Definition: fe_update_flags.h:197
FEFaceValues::dimension
static const unsigned int dimension
Definition: fe_values.h:3839
FEValuesBase::finite_element_output
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > finite_element_output
Definition: fe_values.h:3550
ProductType
Definition: template_constraints.h:422
FEValues
Definition: fe.h:38
FEValuesViews::Scalar::gradient_type
::Tensor< 1, spacedim > gradient_type
Definition: fe_values.h:161
FEValuesBase::dofs_per_cell
const unsigned int dofs_per_cell
Definition: fe_values.h:2108
FEValuesViews::Vector::first_vector_component
const unsigned int first_vector_component
Definition: fe_values.h:1223
FEValuesViews::Vector::shape_function_data
std::vector< ShapeFunctionData > shape_function_data
Definition: fe_values.h:1228
FEValuesViews::Scalar::OutputType::third_derivative_type
typename ProductType< Number, typename Scalar< dim, spacedim >::third_derivative_type >::type third_derivative_type
Definition: fe_values.h:222
FEValuesViews::Vector::ShapeFunctionData::row_index
unsigned int row_index[spacedim]
Definition: fe_values.h:741
Subscriptor
Definition: subscriptor.h:62
internal::FEValuesViews::Cache
Definition: fe_values.h:1949
DerivativeForm< 1, dim, spacedim >
FEValuesViews::Vector::ShapeFunctionData::single_nonzero_component_index
unsigned int single_nonzero_component_index
Definition: fe_values.h:752
FEValuesViews::Vector::Vector
Vector()
Definition: fe_values.cc:258
FEValuesBase::CellIterator
Definition: fe_values.h:3457
FEValuesViews::Vector::get_function_divergences_from_local_dof_values
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:2048
fe_values_extractors.h
FESubfaceValues::do_reinit
void do_reinit(const unsigned int face_no, const unsigned int subface_no)
Definition: fe_values.cc:5049
FEValuesViews::Vector::get_function_laplacians
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:2185
FEValuesViews::Scalar::third_derivative
third_derivative_type third_derivative(const unsigned int shape_function, const unsigned int q_point) const
FEValuesExtractors::SymmetricTensor
Definition: fe_values_extractors.h:199
FEValuesBase::memory_consumption
std::size_t memory_consumption() const
Definition: fe_values.cc:4212
internal::FEValuesViews::Cache::symmetric_second_order_tensors
std::vector<::FEValuesViews::SymmetricTensor< 2, dim, spacedim > > symmetric_second_order_tensors
Definition: fe_values.h:1958
Tensor::unrolled_to_component_indices
static constexpr TableIndices< rank_ > unrolled_to_component_indices(const unsigned int i)
FEValuesExtractors::Vector::first_vector_component
unsigned int first_vector_component
Definition: fe_values_extractors.h:155
Mapping
Abstract base class for mapping classes.
Definition: mapping.h:302
FEValuesBase::get_cell
const Triangulation< dim, spacedim >::cell_iterator get_cell() const
Definition: fe_values.cc:4190
FEFaceValuesBase::get_boundary_forms
const std::vector< Tensor< 1, spacedim > > & get_boundary_forms() const
Definition: fe_values.cc:4625
FEValuesBase::jacobian_grad
const DerivativeForm< 2, dim, spacedim > & jacobian_grad(const unsigned int quadrature_point) const
FEValuesBase::shape_grad
const Tensor< 1, spacedim > & shape_grad(const unsigned int function_no, const unsigned int quadrature_point) const
FEValuesViews::Vector::symmetric_gradient
symmetric_gradient_type symmetric_gradient(const unsigned int shape_function, const unsigned int q_point) const
FEValuesViews::Scalar::OutputType::hessian_type
typename ProductType< Number, typename Scalar< dim, spacedim >::hessian_type >::type hessian_type
Definition: fe_values.h:214
FEValuesExtractors::Vector
Definition: fe_values_extractors.h:150
StandardExceptions::ExcIndexRange
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
FEValuesViews::Scalar::get_function_hessians_from_local_dof_values
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:1713
FEValuesBase::space_dimension
static const unsigned int space_dimension
Definition: fe_values.h:2096
FEValuesBase::mapping
const SmartPointer< const Mapping< dim, spacedim >, FEValuesBase< dim, spacedim > > mapping
Definition: fe_values.h:3510
FEValuesBase::get_inverse_jacobians
const std::vector< DerivativeForm< 1, spacedim, dim > > & get_inverse_jacobians() const
FEValuesViews::Vector::get_function_symmetric_gradients
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:1961
FEValuesViews::Scalar::operator=
Scalar & operator=(const Scalar< dim, spacedim > &)=delete
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
double
FEFaceValues::do_reinit
void do_reinit(const unsigned int face_no)
Definition: fe_values.cc:4815
FEValuesBase::jacobian_pushed_forward_3rd_derivative
const Tensor< 5, spacedim > & jacobian_pushed_forward_3rd_derivative(const unsigned int quadrature_point) const
FEValuesBase::fe_data
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > fe_data
Definition: fe_values.h:3542
subscriptor.h
FEValuesBase::check_cell_similarity
void check_cell_similarity(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
Definition: fe_values.cc:4312
FEFaceValuesBase::boundary_form
const Tensor< 1, spacedim > & boundary_form(const unsigned int i) const
Tensor
Definition: tensor.h:450
FEValuesBase::jacobian_2nd_derivative
const DerivativeForm< 3, dim, spacedim > & jacobian_2nd_derivative(const unsigned int quadrature_point) const
FEValuesViews::Scalar::get_function_hessians
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:1682
FEValuesViews::Vector::OutputType::symmetric_gradient_type
typename ProductType< Number, typename Vector< dim, spacedim >::symmetric_gradient_type >::type symmetric_gradient_type
Definition: fe_values.h:673
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
FEValuesViews::Scalar::get_function_laplacians
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:1738
update_jacobians
@ update_jacobians
Volume element.
Definition: fe_update_flags.h:150
update_gradients
@ update_gradients
Shape function gradients.
Definition: fe_update_flags.h:84
FEValuesBase::fe_values_views_cache
::internal::FEValuesViews::Cache< dim, spacedim > fe_values_views_cache
Definition: fe_values.h:3589
FEValuesBase::shape_grad_component
Tensor< 1, spacedim > shape_grad_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
FEValues::reinit
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access >> &cell)
FEValuesBase::get_quadrature_points
const std::vector< Point< spacedim > > & get_quadrature_points() const
FESubfaceValues::integral_dimension
static const unsigned int integral_dimension
Definition: fe_values.h:3995
fe.h
FEValuesViews::Vector::get_function_curls
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:2073
FEFaceValues::space_dimension
static const unsigned int space_dimension
Definition: fe_values.h:3841
DoFCellAccessor
Definition: dof_accessor.h:1321
DeclException1
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:518
FEValuesViews::Scalar::ShapeFunctionData::row_index
unsigned int row_index
Definition: fe_values.h:249
FEValuesViews::SymmetricTensor< 2, dim, spacedim >::fe_values
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
Definition: fe_values.h:1510
update_jacobian_pushed_forward_2nd_derivatives
@ update_jacobian_pushed_forward_2nd_derivatives
Definition: fe_update_flags.h:206
FEFaceValuesBase::FEFaceValuesBase
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:4605
FEValuesViews::Vector::OutputType::third_derivative_type
typename ProductType< Number, typename Vector< dim, spacedim >::third_derivative_type >::type third_derivative_type
Definition: fe_values.h:713
FEValuesViews::Vector::OutputType::value_type
typename ProductType< Number, typename Vector< dim, spacedim >::value_type >::type value_type
Definition: fe_values.h:657
FEValuesBase::JxW
double JxW(const unsigned int quadrature_point) const
FEValuesViews::Vector::OutputType::gradient_type
typename ProductType< Number, typename Vector< dim, spacedim >::gradient_type >::type gradient_type
Definition: fe_values.h:665
FEValuesBase::compute_update_flags
UpdateFlags compute_update_flags(const UpdateFlags update_flags) const
Definition: fe_values.cc:4232
FESubfaceValues::dimension
static const unsigned int dimension
Definition: fe_values.h:3984
FEValues::integral_dimension
static const unsigned int integral_dimension
Definition: fe_values.h:3623
update_jacobian_3rd_derivatives
@ update_jacobian_3rd_derivatives
Definition: fe_update_flags.h:210
FiniteElement
Definition: fe.h:648
FEValues::do_reinit
void do_reinit()
Definition: fe_values.cc:4560
FEFaceValuesBase::memory_consumption
std::size_t memory_consumption() const
Definition: fe_values.cc:4637
FEValuesBase
Definition: fe.h:36
exceptions.h
UpdateFlags
UpdateFlags
Definition: fe_update_flags.h:66
FEValuesBase::get_jacobian_2nd_derivatives
const std::vector< DerivativeForm< 3, dim, spacedim > > & get_jacobian_2nd_derivatives() const
FEFaceValues::reinit
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access >> &cell, const unsigned int face_no)
FEValuesBase::shape_hessian
const Tensor< 2, spacedim > & shape_hessian(const unsigned int function_no, const unsigned int point_no) const
FEValuesViews::Vector::hessian
hessian_type hessian(const unsigned int shape_function, const unsigned int q_point) const
update_jacobian_grads
@ update_jacobian_grads
Gradient of volume element.
Definition: fe_update_flags.h:155
SmartPointer
Definition: smartpointer.h:68
FEValuesBase::inverse_jacobian
const DerivativeForm< 1, spacedim, dim > & inverse_jacobian(const unsigned int quadrature_point) const
value
static const bool value
Definition: dof_tools_constraints.cc:433
FESubfaceValues
Definition: fe.h:42
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
update_JxW_values
@ update_JxW_values
Transformed quadrature weights.
Definition: fe_update_flags.h:129
update_normal_vectors
@ update_normal_vectors
Normal vectors.
Definition: fe_update_flags.h:136
mapping.h
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
update_hessians
@ update_hessians
Second derivatives of shape functions.
Definition: fe_update_flags.h:90
FEValuesBase::get_fe
const FiniteElement< dim, spacedim > & get_fe() const
std_cxx20::ranges::iota_view
boost::integer_range< IncrementableType > iota_view
Definition: iota_view.h:48
symmetric_tensor.h
FESubfaceValues::ExcReinitCalledWithBoundaryFace
static ::ExceptionBase & ExcReinitCalledWithBoundaryFace()
FEValuesBase::operator=
FEValuesBase & operator=(const FEValuesBase &)=delete
internal::FEValuesViews::ViewType< dim, spacedim, FEValuesExtractors::Scalar >::type
typename ::FEValuesViews::Scalar< dim, spacedim > type
Definition: fe_values.h:1898
FEValuesBase::get_mapping
const Mapping< dim, spacedim > & get_mapping() const
FEValuesBase::get_function_values
void get_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
Definition: fe_values.cc:3572
internal::NumberType
Definition: numbers.h:700
FEValuesViews::Vector::hessian_type
::Tensor< 3, spacedim > hessian_type
Definition: fe_values.h:635
FEFaceValuesBase::present_face_index
unsigned int present_face_index
Definition: fe_values.h:3807
FEValuesViews::Scalar::value_type
double value_type
Definition: fe_values.h:154
dof_handler.h
dof_accessor.h
DeclException0
#define DeclException0(Exception0)
Definition: exceptions.h:473
internal::FEValuesViews::ViewType
Definition: fe_values.h:1885
FEValuesBase::mapping_data
std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > mapping_data
Definition: fe_values.h:3518
internal::FEValuesViews::Cache::scalars
std::vector<::FEValuesViews::Scalar< dim, spacedim > > scalars
Definition: fe_values.h:1955
FEValuesViews::Scalar::get_function_third_derivatives_from_local_dof_values
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:1825
FEValuesViews::Scalar::get_function_gradients_from_local_dof_values
void get_function_gradients_from_local_dof_values(const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::gradient_type > &gradients) const
Definition: fe_values.cc:1657
FEValuesBase::jacobian
const DerivativeForm< 1, dim, spacedim > & jacobian(const unsigned int quadrature_point) const
FEValuesBase::get_jacobian_3rd_derivatives
const std::vector< DerivativeForm< 4, dim, spacedim > > & get_jacobian_3rd_derivatives() const
Tensor::value_type
typename Tensor< rank_ - 1, dim, Number >::tensor_type value_type
Definition: tensor.h:484
FEValuesExtractors::Tensor::first_tensor_component
unsigned int first_tensor_component
Definition: fe_values_extractors.h:253
FEValuesViews::Vector::divergence
divergence_type divergence(const unsigned int shape_function, const unsigned int q_point) const
FEValuesViews::Scalar::ShapeFunctionData
Definition: fe_values.h:229
Point< spacedim >
FEFaceValuesBase
Definition: fe_values.h:3735
FEValuesViews::Vector::get_function_values
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:1849
quadrature.h
update_inverse_jacobians
@ update_inverse_jacobians
Volume element.
Definition: fe_update_flags.h:161
FEValuesViews::Tensor< 2, dim, spacedim >::shape_function_data
std::vector< ShapeFunctionData > shape_function_data
Definition: fe_values.h:1868
FEValuesExtractors::Scalar::component
unsigned int component
Definition: fe_values_extractors.h:100
FEValuesBase::get_update_flags
UpdateFlags get_update_flags() const
FEValuesViews::SymmetricTensor< 2, dim, spacedim >::ShapeFunctionData::single_nonzero_component
int single_nonzero_component
Definition: fe_values.h:1342
FEValuesViews::Scalar::Scalar
Scalar()
Definition: fe_values.cc:180
FEValuesBase::n_quadrature_points
const unsigned int n_quadrature_points
Definition: fe_values.h:2101
FEValuesViews::Vector::get_function_third_derivatives_from_local_dof_values
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:2280
FEValuesViews::Tensor< 2, dim, spacedim >::fe_values
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
Definition: fe_values.h:1857
config.h
FEValuesBase::quadrature_point
const Point< spacedim > & quadrature_point(const unsigned int q) const
FEValuesBase::get_jacobian_pushed_forward_2nd_derivatives
const std::vector< Tensor< 4, spacedim > > & get_jacobian_pushed_forward_2nd_derivatives() const
FEValuesBase::shape_3rd_derivative
const Tensor< 3, spacedim > & shape_3rd_derivative(const unsigned int function_no, const unsigned int point_no) const
FEValuesViews::Vector::gradient
gradient_type gradient(const unsigned int shape_function, const unsigned int q_point) const
FEFaceValues::get_present_fe_values
const FEFaceValues< dim, spacedim > & get_present_fe_values() const
FEValuesBase::dof_indices
std_cxx20::ranges::iota_view< unsigned int, unsigned int > dof_indices() const
FEValuesBase::fe
const SmartPointer< const FiniteElement< dim, spacedim >, FEValuesBase< dim, spacedim > > fe
Definition: fe_values.h:3534
FEValuesBase::shape_hessian_component
Tensor< 2, spacedim > shape_hessian_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
FEValuesBase::jacobian_pushed_forward_grad
const Tensor< 3, spacedim > & jacobian_pushed_forward_grad(const unsigned int quadrature_point) const
MeshWorker::internal::OutputType
typename FEValuesViews::View< dim, spacedim, Extractor >::template OutputType< NumberType > OutputType
Definition: scratch_data.h:46
internal::FEValuesViews::ViewType< dim, spacedim, FEValuesExtractors::Tensor< rank > >::type
typename ::FEValuesViews::Tensor< rank, dim, spacedim > type
Definition: fe_values.h:1924
FEValuesViews::Scalar::get_function_third_derivatives
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:1794
FESubfaceValues::FESubfaceValues
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:4856
FEValuesViews::Vector::operator=
Vector & operator=(const Vector< dim, spacedim > &)=delete
FEValuesExtractors::Tensor
Definition: fe_values_extractors.h:248
FEFaceValues::FEFaceValues
FEFaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature, const UpdateFlags update_flags)
Definition: fe_values.cc:4657
internal
Definition: aligned_vector.h:369
FEValuesViews::Scalar::hessian
hessian_type hessian(const unsigned int shape_function, const unsigned int q_point) const
FEValuesBase::dimension
static const unsigned int dimension
Definition: fe_values.h:2091
FEValues::initialize
void initialize(const UpdateFlags update_flags)
Definition: fe_values.cc:4424
FEFaceValues
Definition: fe.h:40
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
FEValuesViews::Scalar::get_function_values_from_local_dof_values
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:1601
FEValuesViews::Vector::curl
curl_type curl(const unsigned int shape_function, const unsigned int q_point) const
FEValuesViews::Vector::symmetric_gradient_type
::SymmetricTensor< 2, spacedim > symmetric_gradient_type
Definition: fe_values.h:613
FEValuesBase::operator[]
const FEValuesViews::Scalar< dim, spacedim > & operator[](const FEValuesExtractors::Scalar &scalar) const
FEValuesViews::Tensor< 2, dim, spacedim >::first_tensor_component
const unsigned int first_tensor_component
Definition: fe_values.h:1863
FEValuesViews::Vector::get_function_laplacians_from_local_dof_values
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:2221
value_type
Quadrature
Definition: quadrature.h:85
FEValuesViews::Vector::OutputType::curl_type
typename ProductType< Number, typename Vector< dim, spacedim >::curl_type >::type curl_type
Definition: fe_values.h:697
FEValuesBase::get_function_third_derivatives
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:4076
FEFaceValuesBase::get_quadrature
const Quadrature< dim - 1 > & get_quadrature() const
FEValuesViews::Vector::get_function_values_from_local_dof_values
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:1880
FEValuesViews::Vector::OutputType::hessian_type
typename ProductType< Number, typename Vector< dim, spacedim >::hessian_type >::type hessian_type
Definition: fe_values.h:705
FEValuesViews::Scalar::OutputType::value_type
typename ProductType< Number, typename Scalar< dim, spacedim >::value_type >::type value_type
Definition: fe_values.h:190
TriaIterator
Definition: tria_iterator.h:578
FEValuesViews::Vector::divergence_type
double divergence_type
Definition: fe_values.h:620
FESubfaceValues::space_dimension
static const unsigned int space_dimension
Definition: fe_values.h:3989
FEValuesBase::get_jacobians
const std::vector< DerivativeForm< 1, dim, spacedim > > & get_jacobians() const
FEValuesViews::View
typename ::internal::FEValuesViews::ViewType< dim, spacedim, Extractor >::type View
Definition: fe_values.h:1980
FEValues::get_quadrature
const Quadrature< dim > & get_quadrature() const
FEValuesBase::get_JxW_values
const std::vector< double > & get_JxW_values() const
FEValuesBase::FEValuesBase
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:3076
FEValuesViews::Vector::get_function_gradients
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:1905
internal::FEValuesViews::Cache::second_order_tensors
std::vector<::FEValuesViews::Tensor< 2, dim, spacedim > > second_order_tensors
Definition: fe_values.h:1960
FEValuesViews::Scalar::value
value_type value(const unsigned int shape_function, const unsigned int q_point) const
FEValuesBase::dof_indices_ending_at
std_cxx20::ranges::iota_view< unsigned int, unsigned int > dof_indices_ending_at(const unsigned int end_dof_index) const
FEValuesBase::get_cell_similarity
CellSimilarity::Similarity get_cell_similarity() const
Definition: fe_values.cc:4367
update_boundary_forms
@ update_boundary_forms
Outer normal vector, not normalized.
Definition: fe_update_flags.h:103
FEValuesViews::SymmetricTensor< 2, dim, spacedim >
Definition: fe_values.h:1260
FEValuesViews::Vector::value
value_type value(const unsigned int shape_function, const unsigned int q_point) const
FESubfaceValues::ExcFaceHasNoSubfaces
static ::ExceptionBase & ExcFaceHasNoSubfaces()
dof_handler.h
point.h
update_jacobian_2nd_derivatives
@ update_jacobian_2nd_derivatives
Definition: fe_update_flags.h:201
symmetrize
constexpr SymmetricTensor< 2, dim, Number > symmetrize(const Tensor< 2, dim, Number > &t)
Definition: symmetric_tensor.h:3547
internal::FEValuesViews::ViewType< dim, spacedim, FEValuesExtractors::Vector >::type
typename ::FEValuesViews::Vector< dim, spacedim > type
Definition: fe_values.h:1911