Reference documentation for deal.II version 9.1.0-pre
fe_values.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_fe_values_h
17 #define dealii_fe_values_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #include <deal.II/base/derivative_form.h>
23 #include <deal.II/base/exceptions.h>
24 #include <deal.II/base/point.h>
25 #include <deal.II/base/quadrature.h>
26 #include <deal.II/base/subscriptor.h>
28 #include <deal.II/base/table.h>
29 #include <deal.II/base/vector_slice.h>
30 
31 #include <deal.II/dofs/dof_accessor.h>
32 #include <deal.II/dofs/dof_handler.h>
33 
34 #include <deal.II/fe/fe.h>
35 #include <deal.II/fe/fe_update_flags.h>
36 #include <deal.II/fe/fe_values_extractors.h>
37 #include <deal.II/fe/mapping.h>
38 
39 #include <deal.II/grid/tria.h>
40 #include <deal.II/grid/tria_iterator.h>
41 
42 #include <deal.II/hp/dof_handler.h>
43 
44 #include <algorithm>
45 #include <memory>
46 #include <type_traits>
47 
48 
49 // dummy include in order to have the
50 // definition of PetscScalar available
51 // without including other PETSc stuff
52 #ifdef DEAL_II_WITH_PETSC
53 # include <petsc.h>
54 #endif
55 
56 DEAL_II_NAMESPACE_OPEN
57 
58 template <int dim, int spacedim = dim>
59 class FEValuesBase;
60 
61 namespace internal
62 {
67  template <int dim, class NumberType = double>
68  struct CurlType;
69 
76  template <class NumberType>
77  struct CurlType<1, NumberType>
78  {
80  };
81 
88  template <class NumberType>
89  struct CurlType<2, NumberType>
90  {
92  };
93 
100  template <class NumberType>
101  struct CurlType<3, NumberType>
102  {
104  };
105 } // namespace internal
106 
107 
108 
130 namespace FEValuesViews
131 {
143  template <int dim, int spacedim = dim>
144  class Scalar
145  {
146  public:
152  using value_type = double;
153 
160 
167 
174 
179  template <typename Number>
180  struct OutputType
181  {
186  using value_type =
187  typename ProductType<Number,
189 
194  using gradient_type = typename ProductType<
195  Number,
197 
202  using laplacian_type =
203  typename ProductType<Number,
205 
210  using hessian_type = typename ProductType<
211  Number,
213 
218  using third_derivative_type = typename ProductType<
219  Number,
221  };
222 
228  {
238 
247  unsigned int row_index;
248  };
249 
253  Scalar();
254 
260  Scalar(const FEValuesBase<dim, spacedim> &fe_values_base,
261  const unsigned int component);
262 
267  Scalar &
268  operator=(const Scalar<dim, spacedim> &);
269 
283  value_type
284  value(const unsigned int shape_function, const unsigned int q_point) const;
285 
297  gradient(const unsigned int shape_function,
298  const unsigned int q_point) const;
299 
311  hessian(const unsigned int shape_function,
312  const unsigned int q_point) const;
313 
325  third_derivative(const unsigned int shape_function,
326  const unsigned int q_point) const;
327 
345  template <class InputVector>
346  void
347  get_function_values(
348  const InputVector &fe_function,
349  std::vector<typename ProductType<value_type,
350  typename InputVector::value_type>::type>
351  &values) const;
352 
375  template <class InputVector>
376  void
377  get_function_values_from_local_dof_values(
378  const InputVector &dof_values,
379  std::vector<
381  &values) const;
382 
400  template <class InputVector>
401  void
402  get_function_gradients(
403  const InputVector &fe_function,
404  std::vector<typename ProductType<gradient_type,
405  typename InputVector::value_type>::type>
406  &gradients) const;
407 
411  template <class InputVector>
412  void
413  get_function_gradients_from_local_dof_values(
414  const InputVector &dof_values,
415  std::vector<
417  &gradients) const;
418 
436  template <class InputVector>
437  void
438  get_function_hessians(
439  const InputVector &fe_function,
440  std::vector<typename ProductType<hessian_type,
441  typename InputVector::value_type>::type>
442  &hessians) const;
443 
447  template <class InputVector>
448  void
449  get_function_hessians_from_local_dof_values(
450  const InputVector &dof_values,
451  std::vector<
453  &hessians) const;
454 
455 
474  template <class InputVector>
475  void
476  get_function_laplacians(
477  const InputVector &fe_function,
478  std::vector<typename ProductType<value_type,
479  typename InputVector::value_type>::type>
480  &laplacians) const;
481 
485  template <class InputVector>
486  void
487  get_function_laplacians_from_local_dof_values(
488  const InputVector &dof_values,
489  std::vector<
491  &laplacians) const;
492 
493 
512  template <class InputVector>
513  void
514  get_function_third_derivatives(
515  const InputVector &fe_function,
516  std::vector<typename ProductType<third_derivative_type,
517  typename InputVector::value_type>::type>
518  &third_derivatives) const;
519 
523  template <class InputVector>
524  void
525  get_function_third_derivatives_from_local_dof_values(
526  const InputVector & dof_values,
527  std::vector<typename OutputType<typename InputVector::value_type>::
528  third_derivative_type> &third_derivatives) const;
529 
530 
531  private:
536 
541  const unsigned int component;
542 
546  std::vector<ShapeFunctionData> shape_function_data;
547  };
548 
549 
550 
580  template <int dim, int spacedim = dim>
581  class Vector
582  {
583  public:
590 
600 
612 
618  using divergence_type = double;
619 
626  using curl_type = typename ::internal::CurlType<spacedim>::type;
627 
634 
641 
646  template <typename Number>
647  struct OutputType
648  {
653  using value_type =
654  typename ProductType<Number,
656 
661  using gradient_type = typename ProductType<
662  Number,
664 
669  using symmetric_gradient_type = typename ProductType<
670  Number,
672 
677  using divergence_type = typename ProductType<
678  Number,
680 
685  using laplacian_type =
686  typename ProductType<Number,
688 
693  using curl_type =
694  typename ProductType<Number,
696 
701  using hessian_type = typename ProductType<
702  Number,
704 
709  using third_derivative_type = typename ProductType<
710  Number,
712  };
713 
719  {
728  bool is_nonzero_shape_function_component[spacedim];
729 
739  unsigned int row_index[spacedim];
740 
750  unsigned int single_nonzero_component_index;
751  };
752 
756  Vector();
757 
766  Vector(const FEValuesBase<dim, spacedim> &fe_values_base,
767  const unsigned int first_vector_component);
768 
773  Vector &
774  operator=(const Vector<dim, spacedim> &);
775 
792  value_type
793  value(const unsigned int shape_function, const unsigned int q_point) const;
794 
809  gradient(const unsigned int shape_function,
810  const unsigned int q_point) const;
811 
828  symmetric_gradient(const unsigned int shape_function,
829  const unsigned int q_point) const;
830 
842  divergence(const unsigned int shape_function,
843  const unsigned int q_point) const;
844 
865  curl_type
866  curl(const unsigned int shape_function, const unsigned int q_point) const;
867 
879  hessian(const unsigned int shape_function,
880  const unsigned int q_point) const;
881 
893  third_derivative(const unsigned int shape_function,
894  const unsigned int q_point) const;
895 
913  template <class InputVector>
914  void
915  get_function_values(
916  const InputVector &fe_function,
917  std::vector<typename ProductType<value_type,
918  typename InputVector::value_type>::type>
919  &values) const;
920 
943  template <class InputVector>
944  void
945  get_function_values_from_local_dof_values(
946  const InputVector &dof_values,
947  std::vector<
949  &values) const;
950 
968  template <class InputVector>
969  void
970  get_function_gradients(
971  const InputVector &fe_function,
972  std::vector<typename ProductType<gradient_type,
973  typename InputVector::value_type>::type>
974  &gradients) const;
975 
979  template <class InputVector>
980  void
981  get_function_gradients_from_local_dof_values(
982  const InputVector &dof_values,
983  std::vector<
985  &gradients) const;
986 
1010  template <class InputVector>
1011  void
1012  get_function_symmetric_gradients(
1013  const InputVector &fe_function,
1014  std::vector<typename ProductType<symmetric_gradient_type,
1015  typename InputVector::value_type>::type>
1016  &symmetric_gradients) const;
1017 
1021  template <class InputVector>
1022  void
1023  get_function_symmetric_gradients_from_local_dof_values(
1024  const InputVector & dof_values,
1025  std::vector<typename OutputType<typename InputVector::value_type>::
1026  symmetric_gradient_type> &symmetric_gradients) const;
1027 
1046  template <class InputVector>
1047  void
1048  get_function_divergences(
1049  const InputVector &fe_function,
1050  std::vector<typename ProductType<divergence_type,
1051  typename InputVector::value_type>::type>
1052  &divergences) const;
1053 
1057  template <class InputVector>
1058  void
1059  get_function_divergences_from_local_dof_values(
1060  const InputVector &dof_values,
1061  std::vector<
1063  &divergences) const;
1064 
1083  template <class InputVector>
1084  void
1085  get_function_curls(
1086  const InputVector &fe_function,
1087  std::vector<
1088  typename ProductType<curl_type, typename InputVector::value_type>::type>
1089  &curls) const;
1090 
1094  template <class InputVector>
1095  void
1096  get_function_curls_from_local_dof_values(
1097  const InputVector &dof_values,
1098  std::vector<
1100  &curls) const;
1101 
1119  template <class InputVector>
1120  void
1121  get_function_hessians(
1122  const InputVector &fe_function,
1123  std::vector<typename ProductType<hessian_type,
1124  typename InputVector::value_type>::type>
1125  &hessians) const;
1126 
1130  template <class InputVector>
1131  void
1132  get_function_hessians_from_local_dof_values(
1133  const InputVector &dof_values,
1134  std::vector<
1136  &hessians) const;
1137 
1156  template <class InputVector>
1157  void
1158  get_function_laplacians(
1159  const InputVector &fe_function,
1160  std::vector<typename ProductType<value_type,
1161  typename InputVector::value_type>::type>
1162  &laplacians) const;
1163 
1167  template <class InputVector>
1168  void
1169  get_function_laplacians_from_local_dof_values(
1170  const InputVector &dof_values,
1171  std::vector<
1173  &laplacians) const;
1174 
1193  template <class InputVector>
1194  void
1195  get_function_third_derivatives(
1196  const InputVector &fe_function,
1197  std::vector<typename ProductType<third_derivative_type,
1198  typename InputVector::value_type>::type>
1199  &third_derivatives) const;
1200 
1204  template <class InputVector>
1205  void
1206  get_function_third_derivatives_from_local_dof_values(
1207  const InputVector & dof_values,
1208  std::vector<typename OutputType<typename InputVector::value_type>::
1209  third_derivative_type> &third_derivatives) const;
1210 
1211  private:
1216 
1221  const unsigned int first_vector_component;
1222 
1226  std::vector<ShapeFunctionData> shape_function_data;
1227  };
1228 
1229 
1230  template <int rank, int dim, int spacedim = dim>
1231  class SymmetricTensor;
1232 
1257  template <int dim, int spacedim>
1258  class SymmetricTensor<2, dim, spacedim>
1259  {
1260  public:
1268 
1279 
1284  template <typename Number>
1285  struct OutputType
1286  {
1291  using value_type = typename ProductType<
1292  Number,
1293  typename SymmetricTensor<2, dim, spacedim>::value_type>::type;
1294 
1299  using divergence_type = typename ProductType<
1300  Number,
1301  typename SymmetricTensor<2, dim, spacedim>::divergence_type>::type;
1302  };
1303 
1308  struct ShapeFunctionData
1309  {
1318  bool is_nonzero_shape_function_component
1319  [value_type::n_independent_components];
1320 
1330  unsigned int row_index[value_type::n_independent_components];
1331 
1341 
1346  };
1347 
1351  SymmetricTensor();
1352 
1362  SymmetricTensor(const FEValuesBase<dim, spacedim> &fe_values_base,
1363  const unsigned int first_tensor_component);
1364 
1369  SymmetricTensor &
1370  operator=(const SymmetricTensor<2, dim, spacedim> &);
1371 
1389  value_type
1390  value(const unsigned int shape_function, const unsigned int q_point) const;
1391 
1406  divergence(const unsigned int shape_function,
1407  const unsigned int q_point) const;
1408 
1426  template <class InputVector>
1427  void
1428  get_function_values(
1429  const InputVector &fe_function,
1430  std::vector<typename ProductType<value_type,
1431  typename InputVector::value_type>::type>
1432  &values) const;
1433 
1456  template <class InputVector>
1457  void
1458  get_function_values_from_local_dof_values(
1459  const InputVector &dof_values,
1460  std::vector<
1462  &values) const;
1463 
1485  template <class InputVector>
1486  void
1487  get_function_divergences(
1488  const InputVector &fe_function,
1489  std::vector<typename ProductType<divergence_type,
1490  typename InputVector::value_type>::type>
1491  &divergences) const;
1492 
1496  template <class InputVector>
1497  void
1498  get_function_divergences_from_local_dof_values(
1499  const InputVector &dof_values,
1500  std::vector<
1502  &divergences) const;
1503 
1504  private:
1509 
1514  const unsigned int first_tensor_component;
1515 
1519  std::vector<ShapeFunctionData> shape_function_data;
1520  };
1521 
1522 
1523  template <int rank, int dim, int spacedim = dim>
1524  class Tensor;
1525 
1546  template <int dim, int spacedim>
1547  class Tensor<2, dim, spacedim>
1548  {
1549  public:
1555 
1560 
1566 
1571  template <typename Number>
1572  struct OutputType
1573  {
1578  using value_type = typename ProductType<
1579  Number,
1581 
1586  using divergence_type = typename ProductType<
1587  Number,
1588  typename Tensor<2, dim, spacedim>::divergence_type>::type;
1589 
1594  using gradient_type = typename ProductType<
1595  Number,
1596  typename Tensor<2, dim, spacedim>::gradient_type>::type;
1597  };
1598 
1603  struct ShapeFunctionData
1604  {
1613  bool is_nonzero_shape_function_component
1614  [value_type::n_independent_components];
1615 
1625  unsigned int row_index[value_type::n_independent_components];
1626 
1636 
1641  };
1642 
1646  Tensor();
1647 
1657  Tensor(const FEValuesBase<dim, spacedim> &fe_values_base,
1658  const unsigned int first_tensor_component);
1659 
1660 
1665  Tensor &
1666  operator=(const Tensor<2, dim, spacedim> &);
1667 
1684  value_type
1685  value(const unsigned int shape_function, const unsigned int q_point) const;
1686 
1701  divergence(const unsigned int shape_function,
1702  const unsigned int q_point) const;
1703 
1718  gradient(const unsigned int shape_function,
1719  const unsigned int q_point) const;
1720 
1738  template <class InputVector>
1739  void
1740  get_function_values(
1741  const InputVector &fe_function,
1742  std::vector<typename ProductType<value_type,
1743  typename InputVector::value_type>::type>
1744  &values) const;
1745 
1768  template <class InputVector>
1769  void
1770  get_function_values_from_local_dof_values(
1771  const InputVector &dof_values,
1772  std::vector<
1774  &values) const;
1775 
1797  template <class InputVector>
1798  void
1799  get_function_divergences(
1800  const InputVector &fe_function,
1801  std::vector<typename ProductType<divergence_type,
1802  typename InputVector::value_type>::type>
1803  &divergences) const;
1804 
1808  template <class InputVector>
1809  void
1810  get_function_divergences_from_local_dof_values(
1811  const InputVector &dof_values,
1812  std::vector<
1814  &divergences) const;
1815 
1832  template <class InputVector>
1833  void
1834  get_function_gradients(
1835  const InputVector &fe_function,
1836  std::vector<typename ProductType<gradient_type,
1837  typename InputVector::value_type>::type>
1838  &gradients) const;
1839 
1843  template <class InputVector>
1844  void
1845  get_function_gradients_from_local_dof_values(
1846  const InputVector &dof_values,
1847  std::vector<
1849  &gradients) const;
1850 
1851  private:
1856 
1861  const unsigned int first_tensor_component;
1862 
1866  std::vector<ShapeFunctionData> shape_function_data;
1867  };
1868 
1869 } // namespace FEValuesViews
1870 
1871 
1872 namespace internal
1873 {
1874  namespace FEValuesViews
1875  {
1883  template <int dim, int spacedim>
1884  struct Cache
1885  {
1890  std::vector<::FEValuesViews::Scalar<dim, spacedim>> scalars;
1891  std::vector<::FEValuesViews::Vector<dim, spacedim>> vectors;
1892  std::vector<::FEValuesViews::SymmetricTensor<2, dim, spacedim>>
1893  symmetric_second_order_tensors;
1894  std::vector<::FEValuesViews::Tensor<2, dim, spacedim>>
1895  second_order_tensors;
1896 
1900  Cache(const FEValuesBase<dim, spacedim> &fe_values);
1901  };
1902  } // namespace FEValuesViews
1903 } // namespace internal
1904 
1905 
1906 
2009 template <int dim, int spacedim>
2010 class FEValuesBase : public Subscriptor
2011 {
2012 public:
2016  static const unsigned int dimension = dim;
2017 
2021  static const unsigned int space_dimension = spacedim;
2022 
2026  const unsigned int n_quadrature_points;
2027 
2033  const unsigned int dofs_per_cell;
2034 
2035 
2043  FEValuesBase(const unsigned int n_q_points,
2044  const unsigned int dofs_per_cell,
2045  const UpdateFlags update_flags,
2046  const Mapping<dim, spacedim> & mapping,
2047  const FiniteElement<dim, spacedim> &fe);
2048 
2049 
2053  virtual ~FEValuesBase() override;
2054 
2055 
2058 
2059 
2080  const double &
2081  shape_value(const unsigned int function_no,
2082  const unsigned int point_no) const;
2083 
2104  double
2105  shape_value_component(const unsigned int function_no,
2106  const unsigned int point_no,
2107  const unsigned int component) const;
2108 
2134  const Tensor<1, spacedim> &
2135  shape_grad(const unsigned int function_no,
2136  const unsigned int quadrature_point) const;
2137 
2155  shape_grad_component(const unsigned int function_no,
2156  const unsigned int point_no,
2157  const unsigned int component) const;
2158 
2178  const Tensor<2, spacedim> &
2179  shape_hessian(const unsigned int function_no,
2180  const unsigned int point_no) const;
2181 
2199  shape_hessian_component(const unsigned int function_no,
2200  const unsigned int point_no,
2201  const unsigned int component) const;
2202 
2222  const Tensor<3, spacedim> &
2223  shape_3rd_derivative(const unsigned int function_no,
2224  const unsigned int point_no) const;
2225 
2243  shape_3rd_derivative_component(const unsigned int function_no,
2244  const unsigned int point_no,
2245  const unsigned int component) const;
2246 
2248 
2250 
2289  template <class InputVector>
2290  void
2291  get_function_values(
2292  const InputVector & fe_function,
2293  std::vector<typename InputVector::value_type> &values) const;
2294 
2308  template <class InputVector>
2309  void
2310  get_function_values(
2311  const InputVector & fe_function,
2312  std::vector<Vector<typename InputVector::value_type>> &values) const;
2313 
2332  template <class InputVector>
2333  void
2334  get_function_values(
2335  const InputVector & fe_function,
2336  const VectorSlice<const std::vector<types::global_dof_index>> &indices,
2337  std::vector<typename InputVector::value_type> &values) const;
2338 
2360  template <class InputVector>
2361  void
2362  get_function_values(
2363  const InputVector & fe_function,
2364  const VectorSlice<const std::vector<types::global_dof_index>> &indices,
2365  std::vector<Vector<typename InputVector::value_type>> &values) const;
2366 
2367 
2398  template <class InputVector>
2399  void
2400  get_function_values(
2401  const InputVector & fe_function,
2402  const VectorSlice<const std::vector<types::global_dof_index>> &indices,
2403  VectorSlice<std::vector<std::vector<typename InputVector::value_type>>>
2404  values,
2405  const bool quadrature_points_fastest) const;
2406 
2408 
2410 
2449  template <class InputVector>
2450  void
2451  get_function_gradients(
2452  const InputVector &fe_function,
2454  &gradients) const;
2455 
2472  template <class InputVector>
2473  void
2474  get_function_gradients(
2475  const InputVector &fe_function,
2476  std::vector<
2478  &gradients) const;
2479 
2486  template <class InputVector>
2487  void
2488  get_function_gradients(
2489  const InputVector & fe_function,
2490  const VectorSlice<const std::vector<types::global_dof_index>> &indices,
2492  &gradients) const;
2493 
2500  template <class InputVector>
2501  void
2502  get_function_gradients(
2503  const InputVector & fe_function,
2504  const VectorSlice<const std::vector<types::global_dof_index>> &indices,
2505  VectorSlice<std::vector<
2507  gradients,
2508  bool quadrature_points_fastest = false) const;
2509 
2511 
2514 
2554  template <class InputVector>
2555  void
2556  get_function_hessians(
2557  const InputVector &fe_function,
2559  &hessians) const;
2560 
2578  template <class InputVector>
2579  void
2580  get_function_hessians(
2581  const InputVector &fe_function,
2582  std::vector<
2584  & hessians,
2585  bool quadrature_points_fastest = false) const;
2586 
2591  template <class InputVector>
2592  void
2593  get_function_hessians(
2594  const InputVector & fe_function,
2595  const VectorSlice<const std::vector<types::global_dof_index>> &indices,
2597  &hessians) const;
2598 
2605  template <class InputVector>
2606  void
2607  get_function_hessians(
2608  const InputVector & fe_function,
2609  const VectorSlice<const std::vector<types::global_dof_index>> &indices,
2610  VectorSlice<std::vector<
2612  hessians,
2613  bool quadrature_points_fastest = false) const;
2614 
2657  template <class InputVector>
2658  void
2659  get_function_laplacians(
2660  const InputVector & fe_function,
2661  std::vector<typename InputVector::value_type> &laplacians) const;
2662 
2682  template <class InputVector>
2683  void
2684  get_function_laplacians(
2685  const InputVector & fe_function,
2686  std::vector<Vector<typename InputVector::value_type>> &laplacians) const;
2687 
2694  template <class InputVector>
2695  void
2696  get_function_laplacians(
2697  const InputVector & fe_function,
2698  const VectorSlice<const std::vector<types::global_dof_index>> &indices,
2699  std::vector<typename InputVector::value_type> &laplacians) const;
2700 
2707  template <class InputVector>
2708  void
2709  get_function_laplacians(
2710  const InputVector & fe_function,
2711  const VectorSlice<const std::vector<types::global_dof_index>> &indices,
2712  std::vector<Vector<typename InputVector::value_type>> &laplacians) const;
2713 
2720  template <class InputVector>
2721  void
2722  get_function_laplacians(
2723  const InputVector & fe_function,
2724  const VectorSlice<const std::vector<types::global_dof_index>> &indices,
2725  std::vector<std::vector<typename InputVector::value_type>> & laplacians,
2726  bool quadrature_points_fastest = false) const;
2727 
2729 
2731 
2772  template <class InputVector>
2773  void
2774  get_function_third_derivatives(
2775  const InputVector &fe_function,
2777  &third_derivatives) const;
2778 
2797  template <class InputVector>
2798  void
2799  get_function_third_derivatives(
2800  const InputVector &fe_function,
2801  std::vector<
2803  & third_derivatives,
2804  bool quadrature_points_fastest = false) const;
2805 
2810  template <class InputVector>
2811  void
2812  get_function_third_derivatives(
2813  const InputVector & fe_function,
2814  const VectorSlice<const std::vector<types::global_dof_index>> &indices,
2816  &third_derivatives) const;
2817 
2824  template <class InputVector>
2825  void
2826  get_function_third_derivatives(
2827  const InputVector & fe_function,
2828  const VectorSlice<const std::vector<types::global_dof_index>> &indices,
2829  VectorSlice<std::vector<
2831  third_derivatives,
2832  bool quadrature_points_fastest = false) const;
2834 
2836 
2837 
2843  const Point<spacedim> &
2844  quadrature_point(const unsigned int q) const;
2845 
2851  const std::vector<Point<spacedim>> &
2852  get_quadrature_points() const;
2853 
2869  double
2870  JxW(const unsigned int quadrature_point) const;
2871 
2875  const std::vector<double> &
2876  get_JxW_values() const;
2877 
2885  jacobian(const unsigned int quadrature_point) const;
2886 
2893  const std::vector<DerivativeForm<1, dim, spacedim>> &
2894  get_jacobians() const;
2895 
2904  jacobian_grad(const unsigned int quadrature_point) const;
2905 
2912  const std::vector<DerivativeForm<2, dim, spacedim>> &
2913  get_jacobian_grads() const;
2914 
2923  const Tensor<3, spacedim> &
2924  jacobian_pushed_forward_grad(const unsigned int quadrature_point) const;
2925 
2932  const std::vector<Tensor<3, spacedim>> &
2933  get_jacobian_pushed_forward_grads() const;
2934 
2943  jacobian_2nd_derivative(const unsigned int quadrature_point) const;
2944 
2951  const std::vector<DerivativeForm<3, dim, spacedim>> &
2952  get_jacobian_2nd_derivatives() const;
2953 
2963  const Tensor<4, spacedim> &
2964  jacobian_pushed_forward_2nd_derivative(
2965  const unsigned int quadrature_point) const;
2966 
2973  const std::vector<Tensor<4, spacedim>> &
2974  get_jacobian_pushed_forward_2nd_derivatives() const;
2975 
2985  jacobian_3rd_derivative(const unsigned int quadrature_point) const;
2986 
2993  const std::vector<DerivativeForm<4, dim, spacedim>> &
2994  get_jacobian_3rd_derivatives() const;
2995 
3005  const Tensor<5, spacedim> &
3006  jacobian_pushed_forward_3rd_derivative(
3007  const unsigned int quadrature_point) const;
3008 
3015  const std::vector<Tensor<5, spacedim>> &
3016  get_jacobian_pushed_forward_3rd_derivatives() const;
3017 
3025  inverse_jacobian(const unsigned int quadrature_point) const;
3026 
3033  const std::vector<DerivativeForm<1, spacedim, dim>> &
3034  get_inverse_jacobians() const;
3035 
3049  const Tensor<1, spacedim> &
3050  normal_vector(const unsigned int i) const;
3051 
3062  DEAL_II_DEPRECATED
3063  const std::vector<Tensor<1, spacedim>> &
3064  get_all_normal_vectors() const;
3065 
3073  const std::vector<Tensor<1, spacedim>> &
3074  get_normal_vectors() const;
3075 
3077 
3079 
3080 
3090  operator[](const FEValuesExtractors::Scalar &scalar) const;
3091 
3101  operator[](const FEValuesExtractors::Vector &vector) const;
3102 
3113  operator[](const FEValuesExtractors::SymmetricTensor<2> &tensor) const;
3114 
3115 
3125  operator[](const FEValuesExtractors::Tensor<2> &tensor) const;
3126 
3128 
3130 
3131 
3135  const Mapping<dim, spacedim> &
3136  get_mapping() const;
3137 
3142  get_fe() const;
3143 
3147  UpdateFlags
3148  get_update_flags() const;
3149 
3154  get_cell() const;
3155 
3162  get_cell_similarity() const;
3163 
3168  std::size_t
3169  memory_consumption() const;
3171 
3172 
3181  std::string,
3182  << "You are requesting information from an FEValues/FEFaceValues/FESubfaceValues "
3183  << "object for which this kind of information has not been computed. What "
3184  << "information these objects compute is determined by the update_* flags you "
3185  << "pass to the constructor. Here, the operation you are attempting requires "
3186  << "the <" << arg1
3187  << "> flag to be set, but it was apparently not specified "
3188  << "upon construction.");
3189 
3197  ExcFEDontMatch,
3198  "The FiniteElement you provided to FEValues and the FiniteElement that belongs "
3199  "to the DoFHandler that provided the cell iterator do not match.");
3205  DeclException1(ExcShapeFunctionNotPrimitive,
3206  int,
3207  << "The shape function with index " << arg1
3208  << " is not primitive, i.e. it is vector-valued and "
3209  << "has more than one non-zero vector component. This "
3210  << "function cannot be called for these shape functions. "
3211  << "Maybe you want to use the same function with the "
3212  << "_component suffix?");
3213 
3222  "The given FiniteElement is not a primitive element but the requested operation "
3223  "only works for those. See FiniteElement::is_primitive() for more information.");
3224 
3225 protected:
3254  class CellIteratorBase;
3255 
3260  template <typename CI>
3263 
3269  std::unique_ptr<const CellIteratorBase> present_cell;
3270 
3278  boost::signals2::connection tria_listener_refinement;
3279 
3287  boost::signals2::connection tria_listener_mesh_transform;
3288 
3294  void
3295  invalidate_present_cell();
3296 
3306  void
3307  maybe_invalidate_previous_present_cell(
3308  const typename Triangulation<dim, spacedim>::cell_iterator &cell);
3309 
3315 
3321  std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
3323 
3330 
3331 
3339 
3345  std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
3347 
3353  spacedim>
3355 
3356 
3361 
3370  UpdateFlags
3371  compute_update_flags(const UpdateFlags update_flags) const;
3372 
3379 
3385  void
3386  check_cell_similarity(
3387  const typename Triangulation<dim, spacedim>::cell_iterator &cell);
3388 
3389 private:
3394  FEValuesBase(const FEValuesBase &);
3395 
3400  FEValuesBase &
3401  operator=(const FEValuesBase &);
3402 
3407 
3412  template <int, int>
3414  template <int, int>
3415  friend class FEValuesViews::Vector;
3416  template <int, int, int>
3417  friend class FEValuesViews::SymmetricTensor;
3418  template <int, int, int>
3419  friend class FEValuesViews::Tensor;
3420 };
3421 
3422 
3423 
3434 template <int dim, int spacedim = dim>
3435 class FEValues : public FEValuesBase<dim, spacedim>
3436 {
3437 public:
3442  static const unsigned int integral_dimension = dim;
3443 
3448  FEValues(const Mapping<dim, spacedim> & mapping,
3449  const FiniteElement<dim, spacedim> &fe,
3450  const Quadrature<dim> & quadrature,
3451  const UpdateFlags update_flags);
3452 
3459  const Quadrature<dim> & quadrature,
3460  const UpdateFlags update_flags);
3461 
3468  template <template <int, int> class DoFHandlerType, bool level_dof_access>
3469  void
3470  reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
3471  level_dof_access>> &cell);
3472 
3486  void
3487  reinit(const typename Triangulation<dim, spacedim>::cell_iterator &cell);
3488 
3493  const Quadrature<dim> &
3494  get_quadrature() const;
3495 
3500  std::size_t
3501  memory_consumption() const;
3502 
3517  const FEValues<dim, spacedim> &
3518  get_present_fe_values() const;
3519 
3520 private:
3525 
3529  void
3530  initialize(const UpdateFlags update_flags);
3531 
3538  void
3539  do_reinit();
3540 };
3541 
3542 
3553 template <int dim, int spacedim = dim>
3554 class FEFaceValuesBase : public FEValuesBase<dim, spacedim>
3555 {
3556 public:
3561  static const unsigned int integral_dimension = dim - 1;
3562 
3574  FEFaceValuesBase(const unsigned int n_q_points,
3575  const unsigned int dofs_per_cell,
3576  const UpdateFlags update_flags,
3577  const Mapping<dim, spacedim> & mapping,
3578  const FiniteElement<dim, spacedim> &fe,
3579  const Quadrature<dim - 1> & quadrature);
3580 
3588  const Tensor<1, spacedim> &
3589  boundary_form(const unsigned int i) const;
3590 
3597  const std::vector<Tensor<1, spacedim>> &
3598  get_boundary_forms() const;
3599 
3604  unsigned int
3605  get_face_index() const;
3606 
3611  const Quadrature<dim - 1> &
3612  get_quadrature() const;
3613 
3618  std::size_t
3619  memory_consumption() const;
3620 
3621 protected:
3626  unsigned int present_face_index;
3627 
3631  const Quadrature<dim - 1> quadrature;
3632 };
3633 
3634 
3635 
3650 template <int dim, int spacedim = dim>
3651 class FEFaceValues : public FEFaceValuesBase<dim, spacedim>
3652 {
3653 public:
3658  static const unsigned int dimension = dim;
3659 
3660  static const unsigned int space_dimension = spacedim;
3661 
3666  static const unsigned int integral_dimension = dim - 1;
3667 
3672  FEFaceValues(const Mapping<dim, spacedim> & mapping,
3673  const FiniteElement<dim, spacedim> &fe,
3674  const Quadrature<dim - 1> & quadrature,
3675  const UpdateFlags update_flags);
3676 
3683  const Quadrature<dim - 1> & quadrature,
3684  const UpdateFlags update_flags);
3685 
3690  template <template <int, int> class DoFHandlerType, bool level_dof_access>
3691  void
3692  reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
3693  level_dof_access>> &cell,
3694  const unsigned int face_no);
3695 
3709  void
3710  reinit(const typename Triangulation<dim, spacedim>::cell_iterator &cell,
3711  const unsigned int face_no);
3712 
3728  get_present_fe_values() const;
3729 
3730 private:
3734  void
3735  initialize(const UpdateFlags update_flags);
3736 
3743  void
3744  do_reinit(const unsigned int face_no);
3745 };
3746 
3747 
3765 template <int dim, int spacedim = dim>
3766 class FESubfaceValues : public FEFaceValuesBase<dim, spacedim>
3767 {
3768 public:
3772  static const unsigned int dimension = dim;
3773 
3777  static const unsigned int space_dimension = spacedim;
3778 
3783  static const unsigned int integral_dimension = dim - 1;
3784 
3789  FESubfaceValues(const Mapping<dim, spacedim> & mapping,
3790  const FiniteElement<dim, spacedim> &fe,
3791  const Quadrature<dim - 1> & face_quadrature,
3792  const UpdateFlags update_flags);
3793 
3800  const Quadrature<dim - 1> & face_quadrature,
3801  const UpdateFlags update_flags);
3802 
3809  template <template <int, int> class DoFHandlerType, bool level_dof_access>
3810  void
3811  reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
3812  level_dof_access>> &cell,
3813  const unsigned int face_no,
3814  const unsigned int subface_no);
3815 
3829  void
3830  reinit(const typename Triangulation<dim, spacedim>::cell_iterator &cell,
3831  const unsigned int face_no,
3832  const unsigned int subface_no);
3833 
3849  get_present_fe_values() const;
3850 
3856  DeclException0(ExcReinitCalledWithBoundaryFace);
3857 
3863  DeclException0(ExcFaceHasNoSubfaces);
3864 
3865 private:
3869  void
3870  initialize(const UpdateFlags update_flags);
3871 
3878  void
3879  do_reinit(const unsigned int face_no, const unsigned int subface_no);
3880 };
3881 
3882 
3883 #ifndef DOXYGEN
3884 
3885 
3886 /*------------------------ Inline functions: namespace FEValuesViews --------*/
3887 
3888 namespace FEValuesViews
3889 {
3890  template <int dim, int spacedim>
3891  inline typename Scalar<dim, spacedim>::value_type
3892  Scalar<dim, spacedim>::value(const unsigned int shape_function,
3893  const unsigned int q_point) const
3894  {
3895  Assert(shape_function < fe_values->fe->dofs_per_cell,
3896  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
3897  Assert(
3898  fe_values->update_flags & update_values,
3900  "update_values"))));
3901 
3902  // an adaptation of the FEValuesBase::shape_value_component function
3903  // except that here we know the component as fixed and we have
3904  // pre-computed and cached a bunch of information. See the comments there.
3905  if (shape_function_data[shape_function].is_nonzero_shape_function_component)
3906  return fe_values->finite_element_output.shape_values(
3907  shape_function_data[shape_function].row_index, q_point);
3908  else
3909  return 0;
3910  }
3911 
3912 
3913 
3914  template <int dim, int spacedim>
3915  inline typename Scalar<dim, spacedim>::gradient_type
3916  Scalar<dim, spacedim>::gradient(const unsigned int shape_function,
3917  const unsigned int q_point) const
3918  {
3919  Assert(shape_function < fe_values->fe->dofs_per_cell,
3920  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
3921  Assert(fe_values->update_flags & update_gradients,
3923  "update_gradients")));
3924 
3925  // an adaptation of the FEValuesBase::shape_grad_component
3926  // function except that here we know the component as fixed and we have
3927  // pre-computed and cached a bunch of information. See the comments there.
3928  if (shape_function_data[shape_function].is_nonzero_shape_function_component)
3929  return fe_values->finite_element_output
3930  .shape_gradients[shape_function_data[shape_function].row_index]
3931  [q_point];
3932  else
3933  return gradient_type();
3934  }
3935 
3936 
3937 
3938  template <int dim, int spacedim>
3939  inline typename Scalar<dim, spacedim>::hessian_type
3940  Scalar<dim, spacedim>::hessian(const unsigned int shape_function,
3941  const unsigned int q_point) const
3942  {
3943  Assert(shape_function < fe_values->fe->dofs_per_cell,
3944  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
3945  Assert(fe_values->update_flags & update_hessians,
3947  "update_hessians")));
3948 
3949  // an adaptation of the FEValuesBase::shape_hessian_component
3950  // function except that here we know the component as fixed and we have
3951  // pre-computed and cached a bunch of information. See the comments there.
3952  if (shape_function_data[shape_function].is_nonzero_shape_function_component)
3953  return fe_values->finite_element_output
3954  .shape_hessians[shape_function_data[shape_function].row_index][q_point];
3955  else
3956  return hessian_type();
3957  }
3958 
3959 
3960 
3961  template <int dim, int spacedim>
3962  inline typename Scalar<dim, spacedim>::third_derivative_type
3963  Scalar<dim, spacedim>::third_derivative(const unsigned int shape_function,
3964  const unsigned int q_point) const
3965  {
3966  Assert(shape_function < fe_values->fe->dofs_per_cell,
3967  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
3968  Assert(fe_values->update_flags & update_3rd_derivatives,
3970  "update_3rd_derivatives")));
3971 
3972  // an adaptation of the FEValuesBase::shape_3rdderivative_component
3973  // function except that here we know the component as fixed and we have
3974  // pre-computed and cached a bunch of information. See the comments there.
3975  if (shape_function_data[shape_function].is_nonzero_shape_function_component)
3976  return fe_values->finite_element_output
3977  .shape_3rd_derivatives[shape_function_data[shape_function].row_index]
3978  [q_point];
3979  else
3980  return third_derivative_type();
3981  }
3982 
3983 
3984 
3985  template <int dim, int spacedim>
3986  inline typename Vector<dim, spacedim>::value_type
3987  Vector<dim, spacedim>::value(const unsigned int shape_function,
3988  const unsigned int q_point) const
3989  {
3990  Assert(shape_function < fe_values->fe->dofs_per_cell,
3991  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
3992  Assert(fe_values->update_flags & update_values,
3994  "update_values")));
3995 
3996  // same as for the scalar case except that we have one more index
3997  const int snc =
3998  shape_function_data[shape_function].single_nonzero_component;
3999  if (snc == -2)
4000  return value_type();
4001  else if (snc != -1)
4002  {
4003  value_type return_value;
4004  return_value[shape_function_data[shape_function]
4005  .single_nonzero_component_index] =
4006  fe_values->finite_element_output.shape_values(snc, q_point);
4007  return return_value;
4008  }
4009  else
4010  {
4011  value_type return_value;
4012  for (unsigned int d = 0; d < dim; ++d)
4013  if (shape_function_data[shape_function]
4014  .is_nonzero_shape_function_component[d])
4015  return_value[d] = fe_values->finite_element_output.shape_values(
4016  shape_function_data[shape_function].row_index[d], q_point);
4017 
4018  return return_value;
4019  }
4020  }
4021 
4022 
4023 
4024  template <int dim, int spacedim>
4025  inline typename Vector<dim, spacedim>::gradient_type
4026  Vector<dim, spacedim>::gradient(const unsigned int shape_function,
4027  const unsigned int q_point) const
4028  {
4029  Assert(shape_function < fe_values->fe->dofs_per_cell,
4030  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4031  Assert(fe_values->update_flags & update_gradients,
4033  "update_gradients")));
4034 
4035  // same as for the scalar case except that we have one more index
4036  const int snc =
4037  shape_function_data[shape_function].single_nonzero_component;
4038  if (snc == -2)
4039  return gradient_type();
4040  else if (snc != -1)
4041  {
4042  gradient_type return_value;
4043  return_value[shape_function_data[shape_function]
4044  .single_nonzero_component_index] =
4045  fe_values->finite_element_output.shape_gradients[snc][q_point];
4046  return return_value;
4047  }
4048  else
4049  {
4050  gradient_type return_value;
4051  for (unsigned int d = 0; d < dim; ++d)
4052  if (shape_function_data[shape_function]
4053  .is_nonzero_shape_function_component[d])
4054  return_value[d] =
4055  fe_values->finite_element_output.shape_gradients
4056  [shape_function_data[shape_function].row_index[d]][q_point];
4057 
4058  return return_value;
4059  }
4060  }
4061 
4062 
4063 
4064  template <int dim, int spacedim>
4065  inline typename Vector<dim, spacedim>::divergence_type
4066  Vector<dim, spacedim>::divergence(const unsigned int shape_function,
4067  const unsigned int q_point) const
4068  {
4069  // this function works like in the case above
4070  Assert(shape_function < fe_values->fe->dofs_per_cell,
4071  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4072  Assert(fe_values->update_flags & update_gradients,
4074  "update_gradients")));
4075 
4076  // same as for the scalar case except that we have one more index
4077  const int snc =
4078  shape_function_data[shape_function].single_nonzero_component;
4079  if (snc == -2)
4080  return divergence_type();
4081  else if (snc != -1)
4082  return fe_values->finite_element_output
4083  .shape_gradients[snc][q_point][shape_function_data[shape_function]
4084  .single_nonzero_component_index];
4085  else
4086  {
4087  divergence_type return_value = 0;
4088  for (unsigned int d = 0; d < dim; ++d)
4089  if (shape_function_data[shape_function]
4090  .is_nonzero_shape_function_component[d])
4091  return_value +=
4092  fe_values->finite_element_output.shape_gradients
4093  [shape_function_data[shape_function].row_index[d]][q_point][d];
4094 
4095  return return_value;
4096  }
4097  }
4098 
4099 
4100 
4101  template <int dim, int spacedim>
4102  inline typename Vector<dim, spacedim>::curl_type
4103  Vector<dim, spacedim>::curl(const unsigned int shape_function,
4104  const unsigned int q_point) const
4105  {
4106  // this function works like in the case above
4107 
4108  Assert(shape_function < fe_values->fe->dofs_per_cell,
4109  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4110  Assert(fe_values->update_flags & update_gradients,
4112  "update_gradients")));
4113  // same as for the scalar case except that we have one more index
4114  const int snc =
4115  shape_function_data[shape_function].single_nonzero_component;
4116 
4117  if (snc == -2)
4118  return curl_type();
4119 
4120  else
4121  switch (dim)
4122  {
4123  case 1:
4124  {
4125  Assert(false,
4126  ExcMessage(
4127  "Computing the curl in 1d is not a useful operation"));
4128  return curl_type();
4129  }
4130 
4131  case 2:
4132  {
4133  if (snc != -1)
4134  {
4135  curl_type return_value;
4136 
4137  // the single nonzero component can only be zero or one in 2d
4138  if (shape_function_data[shape_function]
4139  .single_nonzero_component_index == 0)
4140  return_value[0] =
4141  -1.0 * fe_values->finite_element_output
4142  .shape_gradients[snc][q_point][1];
4143  else
4144  return_value[0] = fe_values->finite_element_output
4145  .shape_gradients[snc][q_point][0];
4146 
4147  return return_value;
4148  }
4149 
4150  else
4151  {
4152  curl_type return_value;
4153 
4154  return_value[0] = 0.0;
4155 
4156  if (shape_function_data[shape_function]
4157  .is_nonzero_shape_function_component[0])
4158  return_value[0] -=
4159  fe_values->finite_element_output
4160  .shape_gradients[shape_function_data[shape_function]
4161  .row_index[0]][q_point][1];
4162 
4163  if (shape_function_data[shape_function]
4164  .is_nonzero_shape_function_component[1])
4165  return_value[0] +=
4166  fe_values->finite_element_output
4167  .shape_gradients[shape_function_data[shape_function]
4168  .row_index[1]][q_point][0];
4169 
4170  return return_value;
4171  }
4172  }
4173 
4174  case 3:
4175  {
4176  if (snc != -1)
4177  {
4178  curl_type return_value;
4179 
4180  switch (shape_function_data[shape_function]
4181  .single_nonzero_component_index)
4182  {
4183  case 0:
4184  {
4185  return_value[0] = 0;
4186  return_value[1] = fe_values->finite_element_output
4187  .shape_gradients[snc][q_point][2];
4188  return_value[2] =
4189  -1.0 * fe_values->finite_element_output
4190  .shape_gradients[snc][q_point][1];
4191  return return_value;
4192  }
4193 
4194  case 1:
4195  {
4196  return_value[0] =
4197  -1.0 * fe_values->finite_element_output
4198  .shape_gradients[snc][q_point][2];
4199  return_value[1] = 0;
4200  return_value[2] = fe_values->finite_element_output
4201  .shape_gradients[snc][q_point][0];
4202  return return_value;
4203  }
4204 
4205  default:
4206  {
4207  return_value[0] = fe_values->finite_element_output
4208  .shape_gradients[snc][q_point][1];
4209  return_value[1] =
4210  -1.0 * fe_values->finite_element_output
4211  .shape_gradients[snc][q_point][0];
4212  return_value[2] = 0;
4213  return return_value;
4214  }
4215  }
4216  }
4217 
4218  else
4219  {
4220  curl_type return_value;
4221 
4222  for (unsigned int i = 0; i < dim; ++i)
4223  return_value[i] = 0.0;
4224 
4225  if (shape_function_data[shape_function]
4226  .is_nonzero_shape_function_component[0])
4227  {
4228  return_value[1] +=
4229  fe_values->finite_element_output
4230  .shape_gradients[shape_function_data[shape_function]
4231  .row_index[0]][q_point][2];
4232  return_value[2] -=
4233  fe_values->finite_element_output
4234  .shape_gradients[shape_function_data[shape_function]
4235  .row_index[0]][q_point][1];
4236  }
4237 
4238  if (shape_function_data[shape_function]
4239  .is_nonzero_shape_function_component[1])
4240  {
4241  return_value[0] -=
4242  fe_values->finite_element_output
4243  .shape_gradients[shape_function_data[shape_function]
4244  .row_index[1]][q_point][2];
4245  return_value[2] +=
4246  fe_values->finite_element_output
4247  .shape_gradients[shape_function_data[shape_function]
4248  .row_index[1]][q_point][0];
4249  }
4250 
4251  if (shape_function_data[shape_function]
4252  .is_nonzero_shape_function_component[2])
4253  {
4254  return_value[0] +=
4255  fe_values->finite_element_output
4256  .shape_gradients[shape_function_data[shape_function]
4257  .row_index[2]][q_point][1];
4258  return_value[1] -=
4259  fe_values->finite_element_output
4260  .shape_gradients[shape_function_data[shape_function]
4261  .row_index[2]][q_point][0];
4262  }
4263 
4264  return return_value;
4265  }
4266  }
4267  }
4268  // should not end up here
4269  Assert(false, ExcInternalError());
4270  return curl_type();
4271  }
4272 
4273 
4274 
4275  template <int dim, int spacedim>
4276  inline typename Vector<dim, spacedim>::hessian_type
4277  Vector<dim, spacedim>::hessian(const unsigned int shape_function,
4278  const unsigned int q_point) const
4279  {
4280  // this function works like in the case above
4281  Assert(shape_function < fe_values->fe->dofs_per_cell,
4282  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4283  Assert(fe_values->update_flags & update_hessians,
4285  "update_hessians")));
4286 
4287  // same as for the scalar case except that we have one more index
4288  const int snc =
4289  shape_function_data[shape_function].single_nonzero_component;
4290  if (snc == -2)
4291  return hessian_type();
4292  else if (snc != -1)
4293  {
4294  hessian_type return_value;
4295  return_value[shape_function_data[shape_function]
4296  .single_nonzero_component_index] =
4297  fe_values->finite_element_output.shape_hessians[snc][q_point];
4298  return return_value;
4299  }
4300  else
4301  {
4302  hessian_type return_value;
4303  for (unsigned int d = 0; d < dim; ++d)
4304  if (shape_function_data[shape_function]
4305  .is_nonzero_shape_function_component[d])
4306  return_value[d] =
4307  fe_values->finite_element_output.shape_hessians
4308  [shape_function_data[shape_function].row_index[d]][q_point];
4309 
4310  return return_value;
4311  }
4312  }
4313 
4314 
4315 
4316  template <int dim, int spacedim>
4317  inline typename Vector<dim, spacedim>::third_derivative_type
4318  Vector<dim, spacedim>::third_derivative(const unsigned int shape_function,
4319  const unsigned int q_point) const
4320  {
4321  // this function works like in the case above
4322  Assert(shape_function < fe_values->fe->dofs_per_cell,
4323  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4324  Assert(fe_values->update_flags & update_3rd_derivatives,
4326  "update_3rd_derivatives")));
4327 
4328  // same as for the scalar case except that we have one more index
4329  const int snc =
4330  shape_function_data[shape_function].single_nonzero_component;
4331  if (snc == -2)
4332  return third_derivative_type();
4333  else if (snc != -1)
4334  {
4335  third_derivative_type return_value;
4336  return_value[shape_function_data[shape_function]
4337  .single_nonzero_component_index] =
4338  fe_values->finite_element_output.shape_3rd_derivatives[snc][q_point];
4339  return return_value;
4340  }
4341  else
4342  {
4343  third_derivative_type return_value;
4344  for (unsigned int d = 0; d < dim; ++d)
4345  if (shape_function_data[shape_function]
4346  .is_nonzero_shape_function_component[d])
4347  return_value[d] =
4348  fe_values->finite_element_output.shape_3rd_derivatives
4349  [shape_function_data[shape_function].row_index[d]][q_point];
4350 
4351  return return_value;
4352  }
4353  }
4354 
4355 
4356 
4357  namespace internal
4358  {
4363  inline ::SymmetricTensor<2, 1>
4364  symmetrize_single_row(const unsigned int n, const ::Tensor<1, 1> &t)
4365  {
4366  Assert(n < 1, ExcIndexRange(n, 0, 1));
4367  (void)n;
4368 
4369  const double array[1] = {t[0]};
4370  return ::SymmetricTensor<2, 1>(array);
4371  }
4372 
4373 
4374 
4375  inline ::SymmetricTensor<2, 2>
4376  symmetrize_single_row(const unsigned int n, const ::Tensor<1, 2> &t)
4377  {
4378  switch (n)
4379  {
4380  case 0:
4381  {
4382  const double array[3] = {t[0], 0, t[1] / 2};
4383  return ::SymmetricTensor<2, 2>(array);
4384  }
4385  case 1:
4386  {
4387  const double array[3] = {0, t[1], t[0] / 2};
4388  return ::SymmetricTensor<2, 2>(array);
4389  }
4390  default:
4391  {
4392  Assert(false, ExcIndexRange(n, 0, 2));
4393  return ::SymmetricTensor<2, 2>();
4394  }
4395  }
4396  }
4397 
4398 
4399 
4400  inline ::SymmetricTensor<2, 3>
4401  symmetrize_single_row(const unsigned int n, const ::Tensor<1, 3> &t)
4402  {
4403  switch (n)
4404  {
4405  case 0:
4406  {
4407  const double array[6] = {t[0], 0, 0, t[1] / 2, t[2] / 2, 0};
4408  return ::SymmetricTensor<2, 3>(array);
4409  }
4410  case 1:
4411  {
4412  const double array[6] = {0, t[1], 0, t[0] / 2, 0, t[2] / 2};
4413  return ::SymmetricTensor<2, 3>(array);
4414  }
4415  case 2:
4416  {
4417  const double array[6] = {0, 0, t[2], 0, t[0] / 2, t[1] / 2};
4418  return ::SymmetricTensor<2, 3>(array);
4419  }
4420  default:
4421  {
4422  Assert(false, ExcIndexRange(n, 0, 3));
4423  return ::SymmetricTensor<2, 3>();
4424  }
4425  }
4426  }
4427  } // namespace internal
4428 
4429 
4430 
4431  template <int dim, int spacedim>
4432  inline typename Vector<dim, spacedim>::symmetric_gradient_type
4433  Vector<dim, spacedim>::symmetric_gradient(const unsigned int shape_function,
4434  const unsigned int q_point) const
4435  {
4436  Assert(shape_function < fe_values->fe->dofs_per_cell,
4437  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4438  Assert(fe_values->update_flags & update_gradients,
4440  "update_gradients")));
4441 
4442  // same as for the scalar case except that we have one more index
4443  const int snc =
4444  shape_function_data[shape_function].single_nonzero_component;
4445  if (snc == -2)
4446  return symmetric_gradient_type();
4447  else if (snc != -1)
4448  return internal::symmetrize_single_row(
4449  shape_function_data[shape_function].single_nonzero_component_index,
4450  fe_values->finite_element_output.shape_gradients[snc][q_point]);
4451  else
4452  {
4453  gradient_type return_value;
4454  for (unsigned int d = 0; d < dim; ++d)
4455  if (shape_function_data[shape_function]
4456  .is_nonzero_shape_function_component[d])
4457  return_value[d] =
4458  fe_values->finite_element_output.shape_gradients
4459  [shape_function_data[shape_function].row_index[d]][q_point];
4460 
4461  return symmetrize(return_value);
4462  }
4463  }
4464 
4465 
4466 
4467  template <int dim, int spacedim>
4469  SymmetricTensor<2, dim, spacedim>::value(const unsigned int shape_function,
4470  const unsigned int q_point) const
4471  {
4472  Assert(shape_function < fe_values->fe->dofs_per_cell,
4473  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4474  Assert(fe_values->update_flags & update_values,
4476  "update_values")));
4477 
4478  // similar to the vector case where we have more then one index and we need
4479  // to convert between unrolled and component indexing for tensors
4480  const int snc =
4481  shape_function_data[shape_function].single_nonzero_component;
4482 
4483  if (snc == -2)
4484  {
4485  // shape function is zero for the selected components
4486  return value_type();
4487  }
4488  else if (snc != -1)
4489  {
4490  value_type return_value;
4491  const unsigned int comp =
4492  shape_function_data[shape_function].single_nonzero_component_index;
4493  return_value[value_type::unrolled_to_component_indices(comp)] =
4494  fe_values->finite_element_output.shape_values(snc, q_point);
4495  return return_value;
4496  }
4497  else
4498  {
4499  value_type return_value;
4500  for (unsigned int d = 0; d < value_type::n_independent_components; ++d)
4501  if (shape_function_data[shape_function]
4502  .is_nonzero_shape_function_component[d])
4503  return_value[value_type::unrolled_to_component_indices(d)] =
4504  fe_values->finite_element_output.shape_values(
4505  shape_function_data[shape_function].row_index[d], q_point);
4506  return return_value;
4507  }
4508  }
4509 
4510 
4511 
4512  template <int dim, int spacedim>
4515  const unsigned int shape_function,
4516  const unsigned int q_point) const
4517  {
4518  Assert(shape_function < fe_values->fe->dofs_per_cell,
4519  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4520  Assert(fe_values->update_flags & update_gradients,
4522  "update_gradients")));
4523 
4524  const int snc =
4525  shape_function_data[shape_function].single_nonzero_component;
4526 
4527  if (snc == -2)
4528  {
4529  // shape function is zero for the selected components
4530  return divergence_type();
4531  }
4532  else if (snc != -1)
4533  {
4534  // we have a single non-zero component when the symmetric tensor is
4535  // represented in unrolled form. this implies we potentially have
4536  // two non-zero components when represented in component form! we
4537  // will only have one non-zero entry if the non-zero component lies on
4538  // the diagonal of the tensor.
4539  //
4540  // the divergence of a second-order tensor is a first order tensor.
4541  //
4542  // assume the second-order tensor is A with components A_{ij}. then
4543  // A_{ij} = A_{ji} and there is only one (if diagonal) or two non-zero
4544  // entries in the tensorial representation. define the
4545  // divergence as:
4546  // b_i := \dfrac{\partial phi_{ij}}{\partial x_j}.
4547  // (which is incidentally also
4548  // b_j := \dfrac{\partial phi_{ij}}{\partial x_i}).
4549  // In both cases, a sum is implied.
4550  //
4551  // Now, we know the nonzero component in unrolled form: it is indicated
4552  // by 'snc'. we can figure out which tensor components belong to this:
4553  const unsigned int comp =
4554  shape_function_data[shape_function].single_nonzero_component_index;
4555  const unsigned int ii =
4556  value_type::unrolled_to_component_indices(comp)[0];
4557  const unsigned int jj =
4558  value_type::unrolled_to_component_indices(comp)[1];
4559 
4560  // given the form of the divergence above, if ii=jj there is only a
4561  // single nonzero component of the full tensor and the gradient
4562  // equals
4563  // b_ii := \dfrac{\partial phi_{ii,ii}}{\partial x_ii}.
4564  // all other entries of 'b' are zero
4565  //
4566  // on the other hand, if ii!=jj, then there are two nonzero entries in
4567  // the full tensor and
4568  // b_ii := \dfrac{\partial phi_{ii,jj}}{\partial x_ii}.
4569  // b_jj := \dfrac{\partial phi_{ii,jj}}{\partial x_jj}.
4570  // again, all other entries of 'b' are zero
4571  const ::Tensor<1, spacedim> &phi_grad =
4572  fe_values->finite_element_output.shape_gradients[snc][q_point];
4573 
4574  divergence_type return_value;
4575  return_value[ii] = phi_grad[jj];
4576 
4577  if (ii != jj)
4578  return_value[jj] = phi_grad[ii];
4579 
4580  return return_value;
4581  }
4582  else
4583  {
4584  Assert(false, ExcNotImplemented());
4585  divergence_type return_value;
4586  return return_value;
4587  }
4588  }
4589 
4590 
4591 
4592  template <int dim, int spacedim>
4593  inline typename Tensor<2, dim, spacedim>::value_type
4594  Tensor<2, dim, spacedim>::value(const unsigned int shape_function,
4595  const unsigned int q_point) const
4596  {
4597  Assert(shape_function < fe_values->fe->dofs_per_cell,
4598  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4599  Assert(fe_values->update_flags & update_values,
4601  "update_values")));
4602 
4603  // similar to the vector case where we have more then one index and we need
4604  // to convert between unrolled and component indexing for tensors
4605  const int snc =
4606  shape_function_data[shape_function].single_nonzero_component;
4607 
4608  if (snc == -2)
4609  {
4610  // shape function is zero for the selected components
4611  return value_type();
4612  }
4613  else if (snc != -1)
4614  {
4615  value_type return_value;
4616  const unsigned int comp =
4617  shape_function_data[shape_function].single_nonzero_component_index;
4618  const TableIndices<2> indices =
4620  return_value[indices] =
4621  fe_values->finite_element_output.shape_values(snc, q_point);
4622  return return_value;
4623  }
4624  else
4625  {
4626  value_type return_value;
4627  for (unsigned int d = 0; d < dim * dim; ++d)
4628  if (shape_function_data[shape_function]
4629  .is_nonzero_shape_function_component[d])
4630  {
4631  const TableIndices<2> indices =
4633  return_value[indices] =
4634  fe_values->finite_element_output.shape_values(
4635  shape_function_data[shape_function].row_index[d], q_point);
4636  }
4637  return return_value;
4638  }
4639  }
4640 
4641 
4642 
4643  template <int dim, int spacedim>
4645  Tensor<2, dim, spacedim>::divergence(const unsigned int shape_function,
4646  const unsigned int q_point) const
4647  {
4648  Assert(shape_function < fe_values->fe->dofs_per_cell,
4649  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4650  Assert(fe_values->update_flags & update_gradients,
4652  "update_gradients")));
4653 
4654  const int snc =
4655  shape_function_data[shape_function].single_nonzero_component;
4656 
4657  if (snc == -2)
4658  {
4659  // shape function is zero for the selected components
4660  return divergence_type();
4661  }
4662  else if (snc != -1)
4663  {
4664  // we have a single non-zero component when the tensor is
4665  // represented in unrolled form.
4666  //
4667  // the divergence of a second-order tensor is a first order tensor.
4668  //
4669  // assume the second-order tensor is A with components A_{ij},
4670  // then divergence is d_i := \frac{\partial A_{ij}}{\partial x_j}
4671  //
4672  // Now, we know the nonzero component in unrolled form: it is indicated
4673  // by 'snc'. we can figure out which tensor components belong to this:
4674  const unsigned int comp =
4675  shape_function_data[shape_function].single_nonzero_component_index;
4676  const TableIndices<2> indices =
4678  const unsigned int ii = indices[0];
4679  const unsigned int jj = indices[1];
4680 
4681  const ::Tensor<1, spacedim> &phi_grad =
4682  fe_values->finite_element_output.shape_gradients[snc][q_point];
4683 
4684  divergence_type return_value;
4685  // note that we contract \nabla from the right
4686  return_value[ii] = phi_grad[jj];
4687 
4688  return return_value;
4689  }
4690  else
4691  {
4692  Assert(false, ExcNotImplemented());
4693  divergence_type return_value;
4694  return return_value;
4695  }
4696  }
4697 
4698 
4699 
4700  template <int dim, int spacedim>
4702  Tensor<2, dim, spacedim>::gradient(const unsigned int shape_function,
4703  const unsigned int q_point) const
4704  {
4705  Assert(shape_function < fe_values->fe->dofs_per_cell,
4706  ExcIndexRange(shape_function, 0, fe_values->fe->dofs_per_cell));
4707  Assert(fe_values->update_flags & update_gradients,
4709  "update_gradients")));
4710 
4711  const int snc =
4712  shape_function_data[shape_function].single_nonzero_component;
4713 
4714  if (snc == -2)
4715  {
4716  // shape function is zero for the selected components
4717  return gradient_type();
4718  }
4719  else if (snc != -1)
4720  {
4721  // we have a single non-zero component when the tensor is
4722  // represented in unrolled form.
4723  //
4724  // the gradient of a second-order tensor is a third order tensor.
4725  //
4726  // assume the second-order tensor is A with components A_{ij},
4727  // then gradient is B_{ijk} := \frac{\partial A_{ij}}{\partial x_k}
4728  //
4729  // Now, we know the nonzero component in unrolled form: it is indicated
4730  // by 'snc'. we can figure out which tensor components belong to this:
4731  const unsigned int comp =
4732  shape_function_data[shape_function].single_nonzero_component_index;
4733  const TableIndices<2> indices =
4735  const unsigned int ii = indices[0];
4736  const unsigned int jj = indices[1];
4737 
4738  const ::Tensor<1, spacedim> &phi_grad =
4739  fe_values->finite_element_output.shape_gradients[snc][q_point];
4740 
4741  gradient_type return_value;
4742  return_value[ii][jj] = phi_grad;
4743 
4744  return return_value;
4745  }
4746  else
4747  {
4748  Assert(false, ExcNotImplemented());
4749  gradient_type return_value;
4750  return return_value;
4751  }
4752  }
4753 
4754 } // namespace FEValuesViews
4755 
4756 
4757 
4758 /*---------------------- Inline functions: FEValuesBase ---------------------*/
4759 
4760 
4761 
4762 template <int dim, int spacedim>
4764  operator[](const FEValuesExtractors::Scalar &scalar) const
4765 {
4766  Assert(scalar.component < fe_values_views_cache.scalars.size(),
4767  ExcIndexRange(scalar.component,
4768  0,
4769  fe_values_views_cache.scalars.size()));
4770 
4771  return fe_values_views_cache.scalars[scalar.component];
4772 }
4773 
4774 
4775 
4776 template <int dim, int spacedim>
4778  operator[](const FEValuesExtractors::Vector &vector) const
4779 {
4780  Assert(vector.first_vector_component < fe_values_views_cache.vectors.size(),
4782  0,
4783  fe_values_views_cache.vectors.size()));
4784 
4785  return fe_values_views_cache.vectors[vector.first_vector_component];
4786 }
4787 
4788 
4789 
4790 template <int dim, int spacedim>
4794 {
4795  Assert(
4796  tensor.first_tensor_component <
4797  fe_values_views_cache.symmetric_second_order_tensors.size(),
4799  0,
4800  fe_values_views_cache.symmetric_second_order_tensors.size()));
4801 
4802  return fe_values_views_cache
4803  .symmetric_second_order_tensors[tensor.first_tensor_component];
4804 }
4805 
4806 
4807 
4808 template <int dim, int spacedim>
4811  operator[](const FEValuesExtractors::Tensor<2> &tensor) const
4812 {
4814  fe_values_views_cache.second_order_tensors.size(),
4816  0,
4817  fe_values_views_cache.second_order_tensors.size()));
4818 
4819  return fe_values_views_cache
4820  .second_order_tensors[tensor.first_tensor_component];
4821 }
4822 
4823 
4824 
4825 template <int dim, int spacedim>
4826 inline const double &
4827 FEValuesBase<dim, spacedim>::shape_value(const unsigned int i,
4828  const unsigned int j) const
4829 {
4830  Assert(i < fe->dofs_per_cell, ExcIndexRange(i, 0, fe->dofs_per_cell));
4831  Assert(this->update_flags & update_values,
4832  ExcAccessToUninitializedField("update_values"));
4833  Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
4834  Assert(present_cell.get() != nullptr,
4835  ExcMessage("FEValues object is not reinit'ed to any cell"));
4836  // if the entire FE is primitive,
4837  // then we can take a short-cut:
4838  if (fe->is_primitive())
4839  return this->finite_element_output.shape_values(i, j);
4840  else
4841  {
4842  // otherwise, use the mapping
4843  // between shape function
4844  // numbers and rows. note that
4845  // by the assertions above, we
4846  // know that this particular
4847  // shape function is primitive,
4848  // so we can call
4849  // system_to_component_index
4850  const unsigned int row =
4851  this->finite_element_output
4852  .shape_function_to_row_table[i * fe->n_components() +
4853  fe->system_to_component_index(i).first];
4854  return this->finite_element_output.shape_values(row, j);
4855  }
4856 }
4857 
4858 
4859 
4860 template <int dim, int spacedim>
4861 inline double
4863  const unsigned int i,
4864  const unsigned int j,
4865  const unsigned int component) const
4866 {
4867  Assert(i < fe->dofs_per_cell, ExcIndexRange(i, 0, fe->dofs_per_cell));
4868  Assert(this->update_flags & update_values,
4869  ExcAccessToUninitializedField("update_values"));
4870  Assert(component < fe->n_components(),
4871  ExcIndexRange(component, 0, fe->n_components()));
4872  Assert(present_cell.get() != nullptr,
4873  ExcMessage("FEValues object is not reinit'ed to any cell"));
4874 
4875  // check whether the shape function
4876  // is non-zero at all within
4877  // this component:
4878  if (fe->get_nonzero_components(i)[component] == false)
4879  return 0;
4880 
4881  // look up the right row in the
4882  // table and take the data from
4883  // there
4884  const unsigned int row =
4885  this->finite_element_output
4886  .shape_function_to_row_table[i * fe->n_components() + component];
4887  return this->finite_element_output.shape_values(row, j);
4888 }
4889 
4890 
4891 
4892 template <int dim, int spacedim>
4893 inline const Tensor<1, spacedim> &
4894 FEValuesBase<dim, spacedim>::shape_grad(const unsigned int i,
4895  const unsigned int j) const
4896 {
4897  Assert(i < fe->dofs_per_cell, ExcIndexRange(i, 0, fe->dofs_per_cell));
4898  Assert(this->update_flags & update_gradients,
4899  ExcAccessToUninitializedField("update_gradients"));
4900  Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
4901  Assert(present_cell.get() != nullptr,
4902  ExcMessage("FEValues object is not reinit'ed to any cell"));
4903  // if the entire FE is primitive,
4904  // then we can take a short-cut:
4905  if (fe->is_primitive())
4906  return this->finite_element_output.shape_gradients[i][j];
4907  else
4908  {
4909  // otherwise, use the mapping
4910  // between shape function
4911  // numbers and rows. note that
4912  // by the assertions above, we
4913  // know that this particular
4914  // shape function is primitive,
4915  // so we can call
4916  // system_to_component_index
4917  const unsigned int row =
4918  this->finite_element_output
4919  .shape_function_to_row_table[i * fe->n_components() +
4920  fe->system_to_component_index(i).first];
4921  return this->finite_element_output.shape_gradients[row][j];
4922  }
4923 }
4924 
4925 
4926 
4927 template <int dim, int spacedim>
4928 inline Tensor<1, spacedim>
4930  const unsigned int i,
4931  const unsigned int j,
4932  const unsigned int component) const
4933 {
4934  Assert(i < fe->dofs_per_cell, ExcIndexRange(i, 0, fe->dofs_per_cell));
4935  Assert(this->update_flags & update_gradients,
4936  ExcAccessToUninitializedField("update_gradients"));
4937  Assert(component < fe->n_components(),
4938  ExcIndexRange(component, 0, fe->n_components()));
4939  Assert(present_cell.get() != nullptr,
4940  ExcMessage("FEValues object is not reinit'ed to any cell"));
4941  // check whether the shape function
4942  // is non-zero at all within
4943  // this component:
4944  if (fe->get_nonzero_components(i)[component] == false)
4945  return Tensor<1, spacedim>();
4946 
4947  // look up the right row in the
4948  // table and take the data from
4949  // there
4950  const unsigned int row =
4951  this->finite_element_output
4952  .shape_function_to_row_table[i * fe->n_components() + component];
4953  return this->finite_element_output.shape_gradients[row][j];
4954 }
4955 
4956 
4957 
4958 template <int dim, int spacedim>
4959 inline const Tensor<2, spacedim> &
4960 FEValuesBase<dim, spacedim>::shape_hessian(const unsigned int i,
4961  const unsigned int j) const
4962 {
4963  Assert(i < fe->dofs_per_cell, ExcIndexRange(i, 0, fe->dofs_per_cell));
4964  Assert(this->update_flags & update_hessians,
4965  ExcAccessToUninitializedField("update_hessians"));
4966  Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
4967  Assert(present_cell.get() != nullptr,
4968  ExcMessage("FEValues object is not reinit'ed to any cell"));
4969  // if the entire FE is primitive,
4970  // then we can take a short-cut:
4971  if (fe->is_primitive())
4972  return this->finite_element_output.shape_hessians[i][j];
4973  else
4974  {
4975  // otherwise, use the mapping
4976  // between shape function
4977  // numbers and rows. note that
4978  // by the assertions above, we
4979  // know that this particular
4980  // shape function is primitive,
4981  // so we can call
4982  // system_to_component_index
4983  const unsigned int row =
4984  this->finite_element_output
4985  .shape_function_to_row_table[i * fe->n_components() +
4986  fe->system_to_component_index(i).first];
4987  return this->finite_element_output.shape_hessians[row][j];
4988  }
4989 }
4990 
4991 
4992 
4993 template <int dim, int spacedim>
4994 inline Tensor<2, spacedim>
4996  const unsigned int i,
4997  const unsigned int j,
4998  const unsigned int component) const
4999 {
5000  Assert(i < fe->dofs_per_cell, ExcIndexRange(i, 0, fe->dofs_per_cell));
5001  Assert(this->update_flags & update_hessians,
5002  ExcAccessToUninitializedField("update_hessians"));
5003  Assert(component < fe->n_components(),
5004  ExcIndexRange(component, 0, fe->n_components()));
5005  Assert(present_cell.get() != nullptr,
5006  ExcMessage("FEValues object is not reinit'ed to any cell"));
5007  // check whether the shape function
5008  // is non-zero at all within
5009  // this component:
5010  if (fe->get_nonzero_components(i)[component] == false)
5011  return Tensor<2, spacedim>();
5012 
5013  // look up the right row in the
5014  // table and take the data from
5015  // there
5016  const unsigned int row =
5017  this->finite_element_output
5018  .shape_function_to_row_table[i * fe->n_components() + component];
5019  return this->finite_element_output.shape_hessians[row][j];
5020 }
5021 
5022 
5023 
5024 template <int dim, int spacedim>
5025 inline const Tensor<3, spacedim> &
5027  const unsigned int j) const
5028 {
5029  Assert(i < fe->dofs_per_cell, ExcIndexRange(i, 0, fe->dofs_per_cell));
5030  Assert(this->update_flags & update_hessians,
5031  ExcAccessToUninitializedField("update_3rd_derivatives"));
5032  Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
5033  Assert(present_cell.get() != nullptr,
5034  ExcMessage("FEValues object is not reinit'ed to any cell"));
5035  // if the entire FE is primitive,
5036  // then we can take a short-cut:
5037  if (fe->is_primitive())
5038  return this->finite_element_output.shape_3rd_derivatives[i][j];
5039  else
5040  {
5041  // otherwise, use the mapping
5042  // between shape function
5043  // numbers and rows. note that
5044  // by the assertions above, we
5045  // know that this particular
5046  // shape function is primitive,
5047  // so we can call
5048  // system_to_component_index
5049  const unsigned int row =
5050  this->finite_element_output
5051  .shape_function_to_row_table[i * fe->n_components() +
5052  fe->system_to_component_index(i).first];
5053  return this->finite_element_output.shape_3rd_derivatives[row][j];
5054  }
5055 }
5056 
5057 
5058 
5059 template <int dim, int spacedim>
5060 inline Tensor<3, spacedim>
5062  const unsigned int i,
5063  const unsigned int j,
5064  const unsigned int component) const
5065 {
5066  Assert(i < fe->dofs_per_cell, ExcIndexRange(i, 0, fe->dofs_per_cell));
5067  Assert(this->update_flags & update_hessians,
5068  ExcAccessToUninitializedField("update_3rd_derivatives"));
5069  Assert(component < fe->n_components(),
5070  ExcIndexRange(component, 0, fe->n_components()));
5071  Assert(present_cell.get() != nullptr,
5072  ExcMessage("FEValues object is not reinit'ed to any cell"));
5073  // check whether the shape function
5074  // is non-zero at all within
5075  // this component:
5076  if (fe->get_nonzero_components(i)[component] == false)
5077  return Tensor<3, spacedim>();
5078 
5079  // look up the right row in the
5080  // table and take the data from
5081  // there
5082  const unsigned int row =
5083  this->finite_element_output
5084  .shape_function_to_row_table[i * fe->n_components() + component];
5085  return this->finite_element_output.shape_3rd_derivatives[row][j];
5086 }
5087 
5088 
5089 
5090 template <int dim, int spacedim>
5091 inline const FiniteElement<dim, spacedim> &
5093 {
5094  return *fe;
5095 }
5096 
5097 
5098 
5099 template <int dim, int spacedim>
5100 inline const Mapping<dim, spacedim> &
5102 {
5103  return *mapping;
5104 }
5105 
5106 
5107 
5108 template <int dim, int spacedim>
5109 inline UpdateFlags
5111 {
5112  return this->update_flags;
5113 }
5114 
5115 
5116 
5117 template <int dim, int spacedim>
5118 inline const std::vector<Point<spacedim>> &
5120 {
5121  Assert(this->update_flags & update_quadrature_points,
5122  ExcAccessToUninitializedField("update_quadrature_points"));
5123  Assert(present_cell.get() != nullptr,
5124  ExcMessage("FEValues object is not reinit'ed to any cell"));
5125  return this->mapping_output.quadrature_points;
5126 }
5127 
5128 
5129 
5130 template <int dim, int spacedim>
5131 inline const std::vector<double> &
5133 {
5134  Assert(this->update_flags & update_JxW_values,
5135  ExcAccessToUninitializedField("update_JxW_values"));
5136  Assert(present_cell.get() != nullptr,
5137  ExcMessage("FEValues object is not reinit'ed to any cell"));
5138  return this->mapping_output.JxW_values;
5139 }
5140 
5141 
5142 
5143 template <int dim, int spacedim>
5144 inline const std::vector<DerivativeForm<1, dim, spacedim>> &
5146 {
5147  Assert(this->update_flags & update_jacobians,
5148  ExcAccessToUninitializedField("update_jacobians"));
5149  Assert(present_cell.get() != nullptr,
5150  ExcMessage("FEValues object is not reinit'ed to any cell"));
5151  return this->mapping_output.jacobians;
5152 }
5153 
5154 
5155 
5156 template <int dim, int spacedim>
5157 inline const std::vector<DerivativeForm<2, dim, spacedim>> &
5159 {
5160  Assert(this->update_flags & update_jacobian_grads,
5161  ExcAccessToUninitializedField("update_jacobians_grads"));
5162  Assert(present_cell.get() != nullptr,
5163  ExcMessage("FEValues object is not reinit'ed to any cell"));
5164  return this->mapping_output.jacobian_grads;
5165 }
5166 
5167 
5168 
5169 template <int dim, int spacedim>
5170 inline const Tensor<3, spacedim> &
5172  const unsigned int i) const
5173 {
5174  Assert(this->update_flags & update_jacobian_pushed_forward_grads,
5175  ExcAccessToUninitializedField("update_jacobian_pushed_forward_grads"));
5176  Assert(present_cell.get() != nullptr,
5177  ExcMessage("FEValues object is not reinit'ed to any cell"));
5178  return this->mapping_output.jacobian_pushed_forward_grads[i];
5179 }
5180 
5181 
5182 
5183 template <int dim, int spacedim>
5184 inline const std::vector<Tensor<3, spacedim>> &
5186 {
5187  Assert(this->update_flags & update_jacobian_pushed_forward_grads,
5188  ExcAccessToUninitializedField("update_jacobian_pushed_forward_grads"));
5189  Assert(present_cell.get() != nullptr,
5190  ExcMessage("FEValues object is not reinit'ed to any cell"));
5191  return this->mapping_output.jacobian_pushed_forward_grads;
5192 }
5193 
5194 
5195 
5196 template <int dim, int spacedim>
5197 inline const DerivativeForm<3, dim, spacedim> &
5198 FEValuesBase<dim, spacedim>::jacobian_2nd_derivative(const unsigned int i) const
5199 {
5200  Assert(this->update_flags & update_jacobian_2nd_derivatives,
5201  ExcAccessToUninitializedField("update_jacobian_2nd_derivatives"));
5202  Assert(present_cell.get() != nullptr,
5203  ExcMessage("FEValues object is not reinit'ed to any cell"));
5204  return this->mapping_output.jacobian_2nd_derivatives[i];
5205 }
5206 
5207 
5208 
5209 template <int dim, int spacedim>
5210 inline const std::vector<DerivativeForm<3, dim, spacedim>> &
5212 {
5213  Assert(this->update_flags & update_jacobian_2nd_derivatives,
5214  ExcAccessToUninitializedField("update_jacobian_2nd_derivatives"));
5215  Assert(present_cell.get() != nullptr,
5216  ExcMessage("FEValues object is not reinit'ed to any cell"));
5217  return this->mapping_output.jacobian_2nd_derivatives;
5218 }
5219 
5220 
5221 
5222 template <int dim, int spacedim>
5223 inline const Tensor<4, spacedim> &
5225  const unsigned int i) const
5226 {
5229  "update_jacobian_pushed_forward_2nd_derivatives"));
5230  Assert(present_cell.get() != nullptr,
5231  ExcMessage("FEValues object is not reinit'ed to any cell"));
5232  return this->mapping_output.jacobian_pushed_forward_2nd_derivatives[i];
5233 }
5234 
5235 
5236 
5237 template <int dim, int spacedim>
5238 inline const std::vector<Tensor<4, spacedim>> &
5240 {
5241  Assert(this->update_flags & update_jacobian_pushed_forward_2nd_derivatives,
5243  "update_jacobian_pushed_forward_2nd_derivatives"));
5244  Assert(present_cell.get() != nullptr,
5245  ExcMessage("FEValues object is not reinit'ed to any cell"));
5246  return this->mapping_output.jacobian_pushed_forward_2nd_derivatives;
5247 }
5248 
5249 
5250 
5251 template <int dim, int spacedim>
5252 inline const DerivativeForm<4, dim, spacedim> &
5253 FEValuesBase<dim, spacedim>::jacobian_3rd_derivative(const unsigned int i) const
5254 {
5255  Assert(this->update_flags & update_jacobian_3rd_derivatives,
5256  ExcAccessToUninitializedField("update_jacobian_3rd_derivatives"));
5257  Assert(present_cell.get() != nullptr,
5258  ExcMessage("FEValues object is not reinit'ed to any cell"));
5259  return this->mapping_output.jacobian_3rd_derivatives[i];
5260 }
5261 
5262 
5263 
5264 template <int dim, int spacedim>
5265 inline const std::vector<DerivativeForm<4, dim, spacedim>> &
5267 {
5268  Assert(this->update_flags & update_jacobian_3rd_derivatives,
5269  ExcAccessToUninitializedField("update_jacobian_3rd_derivatives"));
5270  Assert(present_cell.get() != nullptr,
5271  ExcMessage("FEValues object is not reinit'ed to any cell"));
5272  return this->mapping_output.jacobian_3rd_derivatives;
5273 }
5274 
5275 
5276 
5277 template <int dim, int spacedim>
5278 inline const Tensor<5, spacedim> &
5280  const unsigned int i) const
5281 {
5284  "update_jacobian_pushed_forward_3rd_derivatives"));
5285  Assert(present_cell.get() != nullptr,
5286  ExcMessage("FEValues object is not reinit'ed to any cell"));
5287  return this->mapping_output.jacobian_pushed_forward_3rd_derivatives[i];
5288 }
5289 
5290 
5291 
5292 template <int dim, int spacedim>
5293 inline const std::vector<Tensor<5, spacedim>> &
5295 {
5296  Assert(this->update_flags & update_jacobian_pushed_forward_3rd_derivatives,
5298  "update_jacobian_pushed_forward_3rd_derivatives"));
5299  Assert(present_cell.get() != nullptr,
5300  ExcMessage("FEValues object is not reinit'ed to any cell"));
5301  return this->mapping_output.jacobian_pushed_forward_3rd_derivatives;
5302 }
5303 
5304 
5305 
5306 template <int dim, int spacedim>
5307 inline const std::vector<DerivativeForm<1, spacedim, dim>> &
5309 {
5310  Assert(this->update_flags & update_inverse_jacobians,
5311  ExcAccessToUninitializedField("update_inverse_jacobians"));
5312  Assert(present_cell.get() != nullptr,
5313  ExcMessage("FEValues object is not reinit'ed to any cell"));
5314  return this->mapping_output.inverse_jacobians;
5315 }
5316 
5317 
5318 
5319 template <int dim, int spacedim>
5320 inline const Point<spacedim> &
5321 FEValuesBase<dim, spacedim>::quadrature_point(const unsigned int i) const
5322 {
5323  Assert(this->update_flags & update_quadrature_points,
5324  ExcAccessToUninitializedField("update_quadrature_points"));
5325  Assert(i < this->mapping_output.quadrature_points.size(),
5326  ExcIndexRange(i, 0, this->mapping_output.quadrature_points.size()));
5327  Assert(present_cell.get() != nullptr,
5328  ExcMessage("FEValues object is not reinit'ed to any cell"));
5329 
5330  return this->mapping_output.quadrature_points[i];
5331 }
5332 
5333 
5334 
5335 template <int dim, int spacedim>
5336 inline double
5337 FEValuesBase<dim, spacedim>::JxW(const unsigned int i) const
5338 {
5339  Assert(this->update_flags & update_JxW_values,
5340  ExcAccessToUninitializedField("update_JxW_values"));
5341  Assert(i < this->mapping_output.JxW_values.size(),
5342  ExcIndexRange(i, 0, this->mapping_output.JxW_values.size()));
5343  Assert(present_cell.get() != nullptr,
5344  ExcMessage("FEValues object is not reinit'ed to any cell"));
5345 
5346  return this->mapping_output.JxW_values[i];
5347 }
5348 
5349 
5350 
5351 template <int dim, int spacedim>
5352 inline const DerivativeForm<1, dim, spacedim> &
5353 FEValuesBase<dim, spacedim>::jacobian(const unsigned int i) const
5354 {
5355  Assert(this->update_flags & update_jacobians,
5356  ExcAccessToUninitializedField("update_jacobians"));
5357  Assert(i < this->mapping_output.jacobians.size(),
5358  ExcIndexRange(i, 0, this->mapping_output.jacobians.size()));
5359  Assert(present_cell.get() != nullptr,
5360  ExcMessage("FEValues object is not reinit'ed to any cell"));
5361 
5362  return this->mapping_output.jacobians[i];
5363 }
5364 
5365 
5366 
5367 template <int dim, int spacedim>
5368 inline const DerivativeForm<2, dim, spacedim> &
5369 FEValuesBase<dim, spacedim>::jacobian_grad(const unsigned int i) const
5370 {
5371  Assert(this->update_flags & update_jacobian_grads,
5372  ExcAccessToUninitializedField("update_jacobians_grads"));
5373  Assert(i < this->mapping_output.jacobian_grads.size(),
5374  ExcIndexRange(i, 0, this->mapping_output.jacobian_grads.size()));
5375  Assert(present_cell.get() != nullptr,
5376  ExcMessage("FEValues object is not reinit'ed to any cell"));
5377 
5378  return this->mapping_output.jacobian_grads[i];
5379 }
5380 
5381 
5382 
5383 template <int dim, int spacedim>
5384 inline const DerivativeForm<1, spacedim, dim> &
5385 FEValuesBase<dim, spacedim>::inverse_jacobian(const unsigned int i) const
5386 {
5387  Assert(this->update_flags & update_inverse_jacobians,
5388  ExcAccessToUninitializedField("update_inverse_jacobians"));
5389  Assert(i < this->mapping_output.inverse_jacobians.size(),
5390  ExcIndexRange(i, 0, this->mapping_output.inverse_jacobians.size()));
5391  Assert(present_cell.get() != nullptr,
5392  ExcMessage("FEValues object is not reinit'ed to any cell"));
5393 
5394  return this->mapping_output.inverse_jacobians[i];
5395 }
5396 
5397 
5398 
5399 template <int dim, int spacedim>
5400 inline const Tensor<1, spacedim> &
5401 FEValuesBase<dim, spacedim>::normal_vector(const unsigned int i) const
5402 {
5403  Assert(this->update_flags & update_normal_vectors,
5405  "update_normal_vectors")));
5406  Assert(i < this->mapping_output.normal_vectors.size(),
5407  ExcIndexRange(i, 0, this->mapping_output.normal_vectors.size()));
5408  Assert(present_cell.get() != nullptr,
5409  ExcMessage("FEValues object is not reinit'ed to any cell"));
5410 
5411  return this->mapping_output.normal_vectors[i];
5412 }
5413 
5414 
5415 
5416 /*--------------------- Inline functions: FEValues --------------------------*/
5417 
5418 
5419 template <int dim, int spacedim>
5420 inline const Quadrature<dim> &
5422 {
5423  return quadrature;
5424 }
5425 
5426 
5427 
5428 template <int dim, int spacedim>
5429 inline const FEValues<dim, spacedim> &
5431 {
5432  return *this;
5433 }
5434 
5435 
5436 /*---------------------- Inline functions: FEFaceValuesBase -----------------*/
5437 
5438 
5439 template <int dim, int spacedim>
5440 inline unsigned int
5442 {
5443  return present_face_index;
5444 }
5445 
5446 
5447 /*----------------------- Inline functions: FE*FaceValues -------------------*/
5448 
5449 template <int dim, int spacedim>
5450 inline const Quadrature<dim - 1> &
5452 {
5453  return quadrature;
5454 }
5455 
5456 
5457 
5458 template <int dim, int spacedim>
5459 inline const FEFaceValues<dim, spacedim> &
5461 {
5462  return *this;
5463 }
5464 
5465 
5466 
5467 template <int dim, int spacedim>
5468 inline const FESubfaceValues<dim, spacedim> &
5470 {
5471  return *this;
5472 }
5473 
5474 
5475 
5476 template <int dim, int spacedim>
5477 inline const Tensor<1, spacedim> &
5478 FEFaceValuesBase<dim, spacedim>::boundary_form(const unsigned int i) const
5479 {
5480  Assert(i < this->mapping_output.boundary_forms.size(),
5481  ExcIndexRange(i, 0, this->mapping_output.boundary_forms.size()));
5482  Assert(this->update_flags & update_boundary_forms,
5484  "update_boundary_forms")));
5485 
5486  return this->mapping_output.boundary_forms[i];
5487 }
5488 
5489 #endif // DOXYGEN
5490 
5491 DEAL_II_NAMESPACE_CLOSE
5492 
5493 #endif
Transformed quadrature weights.
const std::vector< Tensor< 4, spacedim > > & get_jacobian_pushed_forward_2nd_derivatives() const
const Tensor< 3, spacedim > & shape_3rd_derivative(const unsigned int function_no, const unsigned int point_no) const
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
Shape function values.
const DerivativeForm< 1, spacedim, dim > & inverse_jacobian(const unsigned int quadrature_point) const
typename ProductType< Number, typename Vector< dim, spacedim >::curl_type >::type curl_type
Definition: fe_values.h:695
const Tensor< 1, spacedim > & shape_grad(const unsigned int function_no, const unsigned int quadrature_point) const
const DerivativeForm< 1, dim, spacedim > & jacobian(const unsigned int quadrature_point) const
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > fe_data
Definition: fe_values.h:3346
CellSimilarity::Similarity cell_similarity
Definition: fe_values.h:3378
typename ProductType< Number, typename Tensor< 2, dim, spacedim >::gradient_type >::type gradient_type
Definition: fe_values.h:1596
const FEFaceValues< dim, spacedim > & get_present_fe_values() const
typename ProductType< Number, typename SymmetricTensor< 2, dim, spacedim >::value_type >::type value_type
Definition: fe_values.h:1293
unsigned int present_face_index
Definition: fe_values.h:3626
std::vector< ShapeFunctionData > shape_function_data
Definition: fe_values.h:1519
static::ExceptionBase & ExcAccessToUninitializedField()
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
Definition: fe_values.h:535
const DerivativeForm< 3, dim, spacedim > & jacobian_2nd_derivative(const unsigned int quadrature_point) const
const unsigned int dofs_per_cell
Definition: fe_values.h:2033
typename::internal::CurlType< spacedim >::type curl_type
Definition: fe_values.h:626
const unsigned int component
Definition: fe_values.h:541
const Quadrature< dim-1 > quadrature
Definition: fe_values.h:3631
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const
Volume element.
const Tensor< 3, spacedim > & jacobian_pushed_forward_grad(const unsigned int quadrature_point) const
unsigned int get_face_index() const
Outer normal vector, not normalized.
typename ProductType< Number, typename Scalar< dim, spacedim >::hessian_type >::type hessian_type
Definition: fe_values.h:212
std::unique_ptr< const CellIteratorBase > present_cell
Definition: fe_values.h:3262
Transformed quadrature points.
const std::vector< DerivativeForm< 4, dim, spacedim > > & get_jacobian_3rd_derivatives() const
typename ProductType< Number, typename Scalar< dim, spacedim >::value_type >::type value_type
Definition: fe_values.h:188
const SmartPointer< const Mapping< dim, spacedim >, FEValuesBase< dim, spacedim > > mapping
Definition: fe_values.h:3314
static::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
const std::vector< double > & get_JxW_values() const
::internal::FEValuesViews::Cache< dim, spacedim > fe_values_views_cache
Definition: fe_values.h:3406
const std::vector< DerivativeForm< 1, dim, spacedim > > & get_jacobians() const
UpdateFlags get_update_flags() const
static::ExceptionBase & ExcFENotPrimitive()
const Quadrature< dim > & get_quadrature() const
typename ProductType< Number, typename Vector< dim, spacedim >::gradient_type >::type gradient_type
Definition: fe_values.h:663
const std::vector< DerivativeForm< 3, dim, spacedim > > & get_jacobian_2nd_derivatives() const
typename ProductType< Number, typename Scalar< dim, spacedim >::gradient_type >::type gradient_type
Definition: fe_values.h:196
typename ProductType< Number, typename Tensor< 2, dim, spacedim >::value_type >::type value_type
Definition: fe_values.h:1580
static::ExceptionBase & ExcMessage(std::string arg1)
const std::vector< Point< spacedim > > & get_quadrature_points() const
const std::vector< Tensor< 5, spacedim > > & get_jacobian_pushed_forward_3rd_derivatives() const
const Quadrature< dim-1 > & get_quadrature() const
Number value_type
Definition: vector.h:123
typename ProductType< Number, typename Vector< dim, spacedim >::third_derivative_type >::type third_derivative_type
Definition: fe_values.h:711
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:408
const Tensor< 2, spacedim > & shape_hessian(const unsigned int function_no, const unsigned int point_no) const
Third derivatives of shape functions.
std::vector< ShapeFunctionData > shape_function_data
Definition: fe_values.h:1866
#define Assert(cond, exc)
Definition: exceptions.h:1227
UpdateFlags
Tensor< 3, spacedim > shape_3rd_derivative_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
Abstract base class for mapping classes.
Definition: dof_tools.h:57
std::vector< ShapeFunctionData > shape_function_data
Definition: fe_values.h:1226
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:397
const Quadrature< dim > quadrature
Definition: fe_values.h:3524
const FEValuesViews::Scalar< dim, spacedim > & operator[](const FEValuesExtractors::Scalar &scalar) const
const unsigned int first_vector_component
Definition: fe_values.h:1221
#define DeclException0(Exception0)
Definition: exceptions.h:385
const FiniteElement< dim, spacedim > & get_fe() const
const std::vector< Tensor< 3, spacedim > > & get_jacobian_pushed_forward_grads() const
std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > mapping_data
Definition: fe_values.h:3322
typename Tensor< rank_-1, dim, Number >::tensor_type value_type
Definition: tensor.h:427
SymmetricTensor< 2, dim, Number > symmetrize(const Tensor< 2, dim, Number > &t)
Tensor()
Definition: tensor.h:1015
const FEValues< dim, spacedim > & get_present_fe_values() const
Second derivatives of shape functions.
Gradient of volume element.
typename ProductType< Number, typename SymmetricTensor< 2, dim, spacedim >::divergence_type >::type divergence_type
Definition: fe_values.h:1301
static TableIndices< rank_ > unrolled_to_component_indices(const unsigned int i)
Definition: tensor.h:1399
std::vector<::FEValuesViews::Scalar< dim, spacedim > > scalars
Definition: fe_values.h:1890
const unsigned int n_quadrature_points
Definition: fe_values.h:2026
const double & shape_value(const unsigned int function_no, const unsigned int point_no) const
const FESubfaceValues< dim, spacedim > & get_present_fe_values() const
boost::signals2::connection tria_listener_mesh_transform
Definition: fe_values.h:3287
typename ProductType< Number, typename Vector< dim, spacedim >::value_type >::type value_type
Definition: fe_values.h:655
const Tensor< 5, spacedim > & jacobian_pushed_forward_3rd_derivative(const unsigned int quadrature_point) const
const Tensor< 4, spacedim > & jacobian_pushed_forward_2nd_derivative(const unsigned int quadrature_point) const
Definition: mpi.h:55
typename ProductType< Number, typename Vector< dim, spacedim >::symmetric_gradient_type >::type symmetric_gradient_type
Definition: fe_values.h:671
typename ProductType< Number, typename Vector< dim, spacedim >::divergence_type >::type divergence_type
Definition: fe_values.h:679
typename ProductType< Number, typename Vector< dim, spacedim >::hessian_type >::type hessian_type
Definition: fe_values.h:703
Shape function gradients.
Normal vectors.
typename ProductType< Number, typename Tensor< 2, dim, spacedim >::divergence_type >::type divergence_type
Definition: fe_values.h:1588
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
Definition: fe_values.h:1855
Definition: fe.h:36
const std::vector< DerivativeForm< 1, spacedim, dim > > & get_inverse_jacobians() const
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
Definition: fe_values.h:1508
static::ExceptionBase & ExcNotImplemented()
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
Definition: fe_values.h:1215
boost::signals2::connection tria_listener_refinement
Definition: fe_values.h:3278
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > finite_element_output
Definition: fe_values.h:3354
const DerivativeForm< 4, dim, spacedim > & jacobian_3rd_derivative(const unsigned int quadrature_point) const
typename ProductType< Number, typename Vector< dim, spacedim >::value_type >::type laplacian_type
Definition: fe_values.h:687
const std::vector< DerivativeForm< 2, dim, spacedim > > & get_jacobian_grads() const
const Point< spacedim > & quadrature_point(const unsigned int q) const
::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > mapping_output
Definition: fe_values.h:3329
const Mapping< dim, spacedim > & get_mapping() const
Tensor< 1, spacedim > shape_grad_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
typename ProductType< Number, typename Scalar< dim, spacedim >::value_type >::type laplacian_type
Definition: fe_values.h:204
double JxW(const unsigned int quadrature_point) const
typename ProductType< Number, typename Scalar< dim, spacedim >::third_derivative_type >::type third_derivative_type
Definition: fe_values.h:220
std::vector< ShapeFunctionData > shape_function_data
Definition: fe_values.h:546
const DerivativeForm< 2, dim, spacedim > & jacobian_grad(const unsigned int quadrature_point) const
UpdateFlags update_flags
Definition: fe_values.h:3360
static::ExceptionBase & ExcInternalError()
const SmartPointer< const FiniteElement< dim, spacedim >, FEValuesBase< dim, spacedim > > fe
Definition: fe_values.h:3338
const Tensor< 1, spacedim > & boundary_form(const unsigned int i) const
Tensor< 2, spacedim > shape_hessian_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const