deal.II version GIT relicensing-2017-g7ae82fad8b 2024-10-22 19:20:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
fe_values_base.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_fe_values_base_h
16#define dealii_fe_values_base_h
17
18
19#include <deal.II/base/config.h>
20
24#include <deal.II/base/point.h>
28
31
32#include <deal.II/fe/fe.h>
36#include <deal.II/fe/mapping.h>
38
39#include <deal.II/grid/tria.h>
41
43
45
46#include <algorithm>
47#include <memory>
48#include <optional>
49#include <type_traits>
50#include <variant>
51
53
153template <int dim, int spacedim>
155{
156public:
160 static constexpr unsigned int dimension = dim;
161
165 static constexpr unsigned int space_dimension = spacedim;
166
174 const unsigned int n_quadrature_points;
175
185 const unsigned int max_n_quadrature_points;
186
192 const unsigned int dofs_per_cell;
193
194
202 FEValuesBase(const unsigned int n_q_points,
203 const unsigned int dofs_per_cell,
207
213 operator=(const FEValuesBase &) = delete;
214
219 FEValuesBase(const FEValuesBase &) = delete;
220
224 virtual ~FEValuesBase() override;
225
246 void
248
252
274 const double &
275 shape_value(const unsigned int i, const unsigned int q_point) const;
276
297 double
298 shape_value_component(const unsigned int i,
299 const unsigned int q_point,
300 const unsigned int component) const;
301
327 const Tensor<1, spacedim> &
328 shape_grad(const unsigned int i, const unsigned int q_point) const;
329
347 shape_grad_component(const unsigned int i,
348 const unsigned int q_point,
349 const unsigned int component) const;
350
370 const Tensor<2, spacedim> &
371 shape_hessian(const unsigned int i, const unsigned int q_point) const;
372
390 shape_hessian_component(const unsigned int i,
391 const unsigned int q_point,
392 const unsigned int component) const;
393
413 const Tensor<3, spacedim> &
414 shape_3rd_derivative(const unsigned int i, const unsigned int q_point) const;
415
433 shape_3rd_derivative_component(const unsigned int i,
434 const unsigned int q_point,
435 const unsigned int component) const;
436
439
485 template <typename Number>
486 void
487 get_function_values(const ReadVector<Number> &fe_function,
488 std::vector<Number> &values) const;
489
503 template <typename Number>
504 void
505 get_function_values(const ReadVector<Number> &fe_function,
506 std::vector<Vector<Number>> &values) const;
507
564 template <typename Number>
565 void
566 get_function_values(const ReadVector<Number> &fe_function,
568 std::vector<Number> &values) const;
569
578 template <typename Number>
579 void
580 get_function_values(const ReadVector<Number> &fe_function,
582 std::vector<Vector<Number>> &values) const;
583
584
606 template <typename Number>
607 void
608 get_function_values(const ReadVector<Number> &fe_function,
610 ArrayView<std::vector<Number>> values,
611 const bool quadrature_points_fastest) const;
612
615
659 template <typename Number>
660 void
662 const ReadVector<Number> &fe_function,
663 std::vector<Tensor<1, spacedim, Number>> &gradients) const;
664
681 template <typename Number>
682 void
684 const ReadVector<Number> &fe_function,
685 std::vector<std::vector<Tensor<1, spacedim, Number>>> &gradients) const;
686
695 template <typename Number>
696 void
698 const ReadVector<Number> &fe_function,
700 std::vector<Tensor<1, spacedim, Number>> &gradients) const;
701
710 template <typename Number>
711 void
713 const ReadVector<Number> &fe_function,
715 ArrayView<std::vector<Tensor<1, spacedim, Number>>> gradients,
716 const bool quadrature_points_fastest = false) const;
717
722
761 template <typename Number>
762 void
764 const ReadVector<Number> &fe_function,
765 std::vector<Tensor<2, spacedim, Number>> &hessians) const;
766
784 template <typename Number>
785 void
787 const ReadVector<Number> &fe_function,
788 std::vector<std::vector<Tensor<2, spacedim, Number>>> &hessians,
789 const bool quadrature_points_fastest = false) const;
790
799 template <typename Number>
800 void
802 const ReadVector<Number> &fe_function,
804 std::vector<Tensor<2, spacedim, Number>> &hessians) const;
805
814 template <typename Number>
815 void
817 const ReadVector<Number> &fe_function,
819 ArrayView<std::vector<Tensor<2, spacedim, Number>>> hessians,
820 const bool quadrature_points_fastest = false) const;
821
862 template <typename Number>
863 void
865 std::vector<Number> &laplacians) const;
866
886 template <typename Number>
887 void
889 std::vector<Vector<Number>> &laplacians) const;
890
899 template <typename Number>
900 void
902 const ReadVector<Number> &fe_function,
904 std::vector<Number> &laplacians) const;
905
914 template <typename Number>
915 void
917 const ReadVector<Number> &fe_function,
919 std::vector<Vector<Number>> &laplacians) const;
920
929 template <typename Number>
930 void
932 const ReadVector<Number> &fe_function,
934 std::vector<std::vector<Number>> &laplacians,
935 const bool quadrature_points_fastest = false) const;
936
939
979 template <typename Number>
980 void
982 const ReadVector<Number> &fe_function,
983 std::vector<Tensor<3, spacedim, Number>> &third_derivatives) const;
984
1004 template <typename Number>
1005 void
1007 const ReadVector<Number> &fe_function,
1008 std::vector<std::vector<Tensor<3, spacedim, Number>>> &third_derivatives,
1009 const bool quadrature_points_fastest = false) const;
1010
1019 template <typename Number>
1020 void
1022 const ReadVector<Number> &fe_function,
1024 std::vector<Tensor<3, spacedim, Number>> &third_derivatives) const;
1025
1034 template <typename Number>
1035 void
1037 const ReadVector<Number> &fe_function,
1039 ArrayView<std::vector<Tensor<3, spacedim, Number>>> third_derivatives,
1040 const bool quadrature_points_fastest = false) const;
1044
1071
1105 dof_indices_starting_at(const unsigned int start_dof_index) const;
1106
1138 dof_indices_ending_at(const unsigned int end_dof_index) const;
1139
1143
1167
1174 const Point<spacedim> &
1175 quadrature_point(const unsigned int q_point) const;
1176
1182 const std::vector<Point<spacedim>> &
1184
1201 double
1202 JxW(const unsigned int q_point) const;
1203
1207 const std::vector<double> &
1209
1217 jacobian(const unsigned int q_point) const;
1218
1225 const std::vector<DerivativeForm<1, dim, spacedim>> &
1227
1236 jacobian_grad(const unsigned int q_point) const;
1237
1244 const std::vector<DerivativeForm<2, dim, spacedim>> &
1246
1255 const Tensor<3, spacedim> &
1256 jacobian_pushed_forward_grad(const unsigned int q_point) const;
1257
1264 const std::vector<Tensor<3, spacedim>> &
1266
1275 jacobian_2nd_derivative(const unsigned int q_point) const;
1276
1283 const std::vector<DerivativeForm<3, dim, spacedim>> &
1285
1295 const Tensor<4, spacedim> &
1296 jacobian_pushed_forward_2nd_derivative(const unsigned int q_point) const;
1297
1304 const std::vector<Tensor<4, spacedim>> &
1306
1316 jacobian_3rd_derivative(const unsigned int q_point) const;
1317
1324 const std::vector<DerivativeForm<4, dim, spacedim>> &
1326
1336 const Tensor<5, spacedim> &
1337 jacobian_pushed_forward_3rd_derivative(const unsigned int q_point) const;
1338
1345 const std::vector<Tensor<5, spacedim>> &
1347
1355 inverse_jacobian(const unsigned int q_point) const;
1356
1363 const std::vector<DerivativeForm<1, spacedim, dim>> &
1365
1385 const Tensor<1, spacedim> &
1386 normal_vector(const unsigned int q_point) const;
1387
1395 const std::vector<Tensor<1, spacedim>> &
1396 get_normal_vectors() const;
1397
1401
1413
1424
1436
1437
1448
1452
1459
1464 get_fe() const;
1465
1471
1476 get_cell() const;
1477
1484 get_cell_similarity() const;
1485
1490 std::size_t
1491 memory_consumption() const;
1503 std::string,
1504 << "You are requesting information from an FEValues/FEFaceValues/FESubfaceValues "
1505 << "object for which this kind of information has not been computed. What "
1506 << "information these objects compute is determined by the update_* flags you "
1507 << "pass to the constructor. Here, the operation you are attempting requires "
1508 << "the <" << arg1
1509 << "> flag to be set, but it was apparently not specified "
1510 << "upon construction.");
1511
1518 "FEValues object is not reinit'ed to any cell");
1519
1528 "The FiniteElement you provided to FEValues and the FiniteElement that belongs "
1529 "to the DoFHandler that provided the cell iterator do not match.");
1536 int,
1537 << "The shape function with index " << arg1
1538 << " is not primitive, i.e. it is vector-valued and "
1539 << "has more than one non-zero vector component. This "
1540 << "function cannot be called for these shape functions. "
1541 << "Maybe you want to use the same function with the "
1542 << "_component suffix?");
1543
1552 "The given FiniteElement is not a primitive element but the requested operation "
1553 "only works for those. See FiniteElement::is_primitive() for more information.");
1554
1555protected:
1573 {
1574 public:
1577 "You have previously called the FEValues::reinit() function with a "
1578 "cell iterator of type Triangulation<dim,spacedim>::cell_iterator. However, "
1579 "when you do this, you cannot call some functions in the FEValues "
1580 "class, such as the get_function_values/gradients/hessians/third_derivatives "
1581 "functions. If you need these functions, then you need to call "
1582 "FEValues::reinit() with an iterator type that allows to extract "
1583 "degrees of freedom, such as DoFHandler<dim,spacedim>::cell_iterator.");
1584
1590
1596
1602
1608
1612 bool
1613 is_initialized() const;
1614
1621 operator typename Triangulation<dim, spacedim>::cell_iterator() const;
1622
1628 n_dofs_for_dof_handler() const;
1629
1634 template <typename Number>
1635 void
1637 Vector<Number> &out) const;
1638
1639 private:
1645 std::optional<
1646 std::variant<typename Triangulation<dim, spacedim>::cell_iterator,
1650 };
1651
1658
1666 boost::signals2::connection tria_listener_refinement;
1667
1675 boost::signals2::connection tria_listener_mesh_transform;
1676
1682 void
1684
1694 void
1696 const typename Triangulation<dim, spacedim>::cell_iterator &cell);
1697
1704
1710 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
1712
1719
1727
1733 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1735
1741 spacedim>
1743
1744
1749
1760
1767
1773 void
1775 const typename Triangulation<dim, spacedim>::cell_iterator &cell);
1776
1777private:
1782
1787
1788 // Make the view classes friends of this class, since they access internal
1789 // data.
1790 template <int, int>
1792 template <int, int>
1794 template <int, int, int>
1796 template <int, int, int>
1798};
1799
1800#ifndef DOXYGEN
1801
1802/*---------------------- Inline functions: FEValuesBase ---------------------*/
1803
1804template <int dim, int spacedim>
1807 const FEValuesExtractors::Scalar &scalar) const
1808{
1809 AssertIndexRange(scalar.component, fe_values_views_cache.scalars.size());
1810
1811 return fe_values_views_cache.scalars[scalar.component].value_or_initialize(
1812 [scalar, this]() {
1813 return FEValuesViews::Scalar<dim, spacedim>(*this, scalar.component);
1814 });
1815}
1816
1817
1818
1819template <int dim, int spacedim>
1822 const FEValuesExtractors::Vector &vector) const
1823{
1825 fe_values_views_cache.vectors.size());
1826
1827 return fe_values_views_cache.vectors[vector.first_vector_component]
1828 .value_or_initialize([vector, this]() {
1830 *this, vector.first_vector_component);
1831 });
1832}
1833
1834
1835
1836template <int dim, int spacedim>
1839 const FEValuesExtractors::SymmetricTensor<2> &tensor) const
1840{
1841 Assert(
1842 tensor.first_tensor_component <
1843 fe_values_views_cache.symmetric_second_order_tensors.size(),
1845 0,
1846 fe_values_views_cache.symmetric_second_order_tensors.size()));
1847
1848 return fe_values_views_cache
1850 .value_or_initialize([tensor, this]() {
1852 *this, tensor.first_tensor_component);
1853 });
1854}
1855
1856
1857
1858template <int dim, int spacedim>
1861 const FEValuesExtractors::Tensor<2> &tensor) const
1862{
1864 fe_values_views_cache.second_order_tensors.size());
1865
1866 return fe_values_views_cache
1868 .value_or_initialize([tensor, this]() {
1870 *this, tensor.first_tensor_component);
1871 });
1872}
1873
1874
1875
1876template <int dim, int spacedim>
1877inline const double &
1878FEValuesBase<dim, spacedim>::shape_value(const unsigned int i,
1879 const unsigned int q_point) const
1880{
1881 AssertIndexRange(i, fe->n_dofs_per_cell());
1882 Assert(this->update_flags & update_values,
1883 ExcAccessToUninitializedField("update_values"));
1884 Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
1885 Assert(present_cell.is_initialized(), ExcNotReinited());
1886 // if the entire FE is primitive,
1887 // then we can take a short-cut:
1888 if (fe->is_primitive())
1889 return this->finite_element_output.shape_values(i, q_point);
1890 else
1891 {
1892 // otherwise, use the mapping
1893 // between shape function
1894 // numbers and rows. note that
1895 // by the assertions above, we
1896 // know that this particular
1897 // shape function is primitive,
1898 // so we can call
1899 // system_to_component_index
1900 const unsigned int row =
1901 this->finite_element_output
1902 .shape_function_to_row_table[i * fe->n_components() +
1903 fe->system_to_component_index(i).first];
1904 return this->finite_element_output.shape_values(row, q_point);
1905 }
1906}
1907
1908
1909
1910template <int dim, int spacedim>
1911inline double
1913 const unsigned int i,
1914 const unsigned int q_point,
1915 const unsigned int component) const
1916{
1917 AssertIndexRange(i, fe->n_dofs_per_cell());
1918 Assert(this->update_flags & update_values,
1919 ExcAccessToUninitializedField("update_values"));
1920 AssertIndexRange(component, fe->n_components());
1921 Assert(present_cell.is_initialized(), ExcNotReinited());
1922
1923 // check whether the shape function
1924 // is non-zero at all within
1925 // this component:
1926 if (fe->get_nonzero_components(i)[component] == false)
1927 return 0;
1928
1929 // look up the right row in the
1930 // table and take the data from
1931 // there
1932 const unsigned int row =
1933 this->finite_element_output
1934 .shape_function_to_row_table[i * fe->n_components() + component];
1935 return this->finite_element_output.shape_values(row, q_point);
1936}
1937
1938
1939
1940template <int dim, int spacedim>
1941inline const Tensor<1, spacedim> &
1942FEValuesBase<dim, spacedim>::shape_grad(const unsigned int i,
1943 const unsigned int q_point) const
1944{
1945 AssertIndexRange(i, fe->n_dofs_per_cell());
1946 Assert(this->update_flags & update_gradients,
1947 ExcAccessToUninitializedField("update_gradients"));
1948 Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
1949 Assert(present_cell.is_initialized(), ExcNotReinited());
1950 // if the entire FE is primitive,
1951 // then we can take a short-cut:
1952 if (fe->is_primitive())
1953 return this->finite_element_output.shape_gradients[i][q_point];
1954 else
1955 {
1956 // otherwise, use the mapping
1957 // between shape function
1958 // numbers and rows. note that
1959 // by the assertions above, we
1960 // know that this particular
1961 // shape function is primitive,
1962 // so we can call
1963 // system_to_component_index
1964 const unsigned int row =
1965 this->finite_element_output
1966 .shape_function_to_row_table[i * fe->n_components() +
1967 fe->system_to_component_index(i).first];
1968 return this->finite_element_output.shape_gradients[row][q_point];
1969 }
1970}
1971
1972
1973
1974template <int dim, int spacedim>
1977 const unsigned int i,
1978 const unsigned int q_point,
1979 const unsigned int component) const
1980{
1981 AssertIndexRange(i, fe->n_dofs_per_cell());
1982 Assert(this->update_flags & update_gradients,
1983 ExcAccessToUninitializedField("update_gradients"));
1984 AssertIndexRange(component, fe->n_components());
1985 Assert(present_cell.is_initialized(), ExcNotReinited());
1986 // check whether the shape function
1987 // is non-zero at all within
1988 // this component:
1989 if (fe->get_nonzero_components(i)[component] == false)
1990 return Tensor<1, spacedim>();
1991
1992 // look up the right row in the
1993 // table and take the data from
1994 // there
1995 const unsigned int row =
1996 this->finite_element_output
1997 .shape_function_to_row_table[i * fe->n_components() + component];
1998 return this->finite_element_output.shape_gradients[row][q_point];
1999}
2000
2001
2002
2003template <int dim, int spacedim>
2004inline const Tensor<2, spacedim> &
2006 const unsigned int q_point) const
2007{
2008 AssertIndexRange(i, fe->n_dofs_per_cell());
2009 Assert(this->update_flags & update_hessians,
2010 ExcAccessToUninitializedField("update_hessians"));
2011 Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
2012 Assert(present_cell.is_initialized(), ExcNotReinited());
2013 // if the entire FE is primitive,
2014 // then we can take a short-cut:
2015 if (fe->is_primitive())
2016 return this->finite_element_output.shape_hessians[i][q_point];
2017 else
2018 {
2019 // otherwise, use the mapping
2020 // between shape function
2021 // numbers and rows. note that
2022 // by the assertions above, we
2023 // know that this particular
2024 // shape function is primitive,
2025 // so we can call
2026 // system_to_component_index
2027 const unsigned int row =
2028 this->finite_element_output
2029 .shape_function_to_row_table[i * fe->n_components() +
2030 fe->system_to_component_index(i).first];
2031 return this->finite_element_output.shape_hessians[row][q_point];
2032 }
2033}
2034
2035
2036
2037template <int dim, int spacedim>
2040 const unsigned int i,
2041 const unsigned int q_point,
2042 const unsigned int component) const
2043{
2044 AssertIndexRange(i, fe->n_dofs_per_cell());
2045 Assert(this->update_flags & update_hessians,
2046 ExcAccessToUninitializedField("update_hessians"));
2047 AssertIndexRange(component, fe->n_components());
2048 Assert(present_cell.is_initialized(), ExcNotReinited());
2049 // check whether the shape function
2050 // is non-zero at all within
2051 // this component:
2052 if (fe->get_nonzero_components(i)[component] == false)
2053 return Tensor<2, spacedim>();
2054
2055 // look up the right row in the
2056 // table and take the data from
2057 // there
2058 const unsigned int row =
2059 this->finite_element_output
2060 .shape_function_to_row_table[i * fe->n_components() + component];
2061 return this->finite_element_output.shape_hessians[row][q_point];
2062}
2063
2064
2065
2066template <int dim, int spacedim>
2067inline const Tensor<3, spacedim> &
2069 const unsigned int i,
2070 const unsigned int q_point) const
2071{
2072 AssertIndexRange(i, fe->n_dofs_per_cell());
2073 Assert(this->update_flags & update_3rd_derivatives,
2074 ExcAccessToUninitializedField("update_3rd_derivatives"));
2075 Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
2076 Assert(present_cell.is_initialized(), ExcNotReinited());
2077 // if the entire FE is primitive,
2078 // then we can take a short-cut:
2079 if (fe->is_primitive())
2080 return this->finite_element_output.shape_3rd_derivatives[i][q_point];
2081 else
2082 {
2083 // otherwise, use the mapping
2084 // between shape function
2085 // numbers and rows. note that
2086 // by the assertions above, we
2087 // know that this particular
2088 // shape function is primitive,
2089 // so we can call
2090 // system_to_component_index
2091 const unsigned int row =
2092 this->finite_element_output
2093 .shape_function_to_row_table[i * fe->n_components() +
2094 fe->system_to_component_index(i).first];
2095 return this->finite_element_output.shape_3rd_derivatives[row][q_point];
2096 }
2097}
2098
2099
2100
2101template <int dim, int spacedim>
2104 const unsigned int i,
2105 const unsigned int q_point,
2106 const unsigned int component) const
2107{
2108 AssertIndexRange(i, fe->n_dofs_per_cell());
2109 Assert(this->update_flags & update_3rd_derivatives,
2110 ExcAccessToUninitializedField("update_3rd_derivatives"));
2111 AssertIndexRange(component, fe->n_components());
2112 Assert(present_cell.is_initialized(), ExcNotReinited());
2113 // check whether the shape function
2114 // is non-zero at all within
2115 // this component:
2116 if (fe->get_nonzero_components(i)[component] == false)
2117 return Tensor<3, spacedim>();
2118
2119 // look up the right row in the
2120 // table and take the data from
2121 // there
2122 const unsigned int row =
2123 this->finite_element_output
2124 .shape_function_to_row_table[i * fe->n_components() + component];
2125 return this->finite_element_output.shape_3rd_derivatives[row][q_point];
2126}
2127
2128
2129
2130template <int dim, int spacedim>
2131inline const FiniteElement<dim, spacedim> &
2133{
2134 return *fe;
2135}
2136
2137
2138
2139template <int dim, int spacedim>
2140inline const Mapping<dim, spacedim> &
2142{
2143 return *mapping;
2144}
2145
2146
2147
2148template <int dim, int spacedim>
2149inline UpdateFlags
2151{
2152 return this->update_flags;
2153}
2154
2155
2156
2157template <int dim, int spacedim>
2158inline const std::vector<Point<spacedim>> &
2160{
2161 Assert(this->update_flags & update_quadrature_points,
2162 ExcAccessToUninitializedField("update_quadrature_points"));
2163 Assert(present_cell.is_initialized(), ExcNotReinited());
2164 return this->mapping_output.quadrature_points;
2165}
2166
2167
2168
2169template <int dim, int spacedim>
2170inline const std::vector<double> &
2172{
2173 Assert(this->update_flags & update_JxW_values,
2174 ExcAccessToUninitializedField("update_JxW_values"));
2175 Assert(present_cell.is_initialized(), ExcNotReinited());
2176 return this->mapping_output.JxW_values;
2177}
2178
2179
2180
2181template <int dim, int spacedim>
2182inline const std::vector<DerivativeForm<1, dim, spacedim>> &
2184{
2185 Assert(this->update_flags & update_jacobians,
2186 ExcAccessToUninitializedField("update_jacobians"));
2187 Assert(present_cell.is_initialized(), ExcNotReinited());
2188 return this->mapping_output.jacobians;
2189}
2190
2191
2192
2193template <int dim, int spacedim>
2194inline const std::vector<DerivativeForm<2, dim, spacedim>> &
2196{
2197 Assert(this->update_flags & update_jacobian_grads,
2198 ExcAccessToUninitializedField("update_jacobians_grads"));
2199 Assert(present_cell.is_initialized(), ExcNotReinited());
2200 return this->mapping_output.jacobian_grads;
2201}
2202
2203
2204
2205template <int dim, int spacedim>
2206inline const Tensor<3, spacedim> &
2208 const unsigned int q_point) const
2209{
2210 Assert(this->update_flags & update_jacobian_pushed_forward_grads,
2211 ExcAccessToUninitializedField("update_jacobian_pushed_forward_grads"));
2212 Assert(present_cell.is_initialized(), ExcNotReinited());
2213 return this->mapping_output.jacobian_pushed_forward_grads[q_point];
2214}
2215
2216
2217
2218template <int dim, int spacedim>
2219inline const std::vector<Tensor<3, spacedim>> &
2221{
2222 Assert(this->update_flags & update_jacobian_pushed_forward_grads,
2223 ExcAccessToUninitializedField("update_jacobian_pushed_forward_grads"));
2224 Assert(present_cell.is_initialized(), ExcNotReinited());
2225 return this->mapping_output.jacobian_pushed_forward_grads;
2226}
2227
2228
2229
2230template <int dim, int spacedim>
2233 const unsigned int q_point) const
2234{
2235 Assert(this->update_flags & update_jacobian_2nd_derivatives,
2236 ExcAccessToUninitializedField("update_jacobian_2nd_derivatives"));
2237 Assert(present_cell.is_initialized(), ExcNotReinited());
2238 return this->mapping_output.jacobian_2nd_derivatives[q_point];
2239}
2240
2241
2242
2243template <int dim, int spacedim>
2244inline const std::vector<DerivativeForm<3, dim, spacedim>> &
2246{
2247 Assert(this->update_flags & update_jacobian_2nd_derivatives,
2248 ExcAccessToUninitializedField("update_jacobian_2nd_derivatives"));
2249 Assert(present_cell.is_initialized(), ExcNotReinited());
2250 return this->mapping_output.jacobian_2nd_derivatives;
2251}
2252
2253
2254
2255template <int dim, int spacedim>
2256inline const Tensor<4, spacedim> &
2258 const unsigned int q_point) const
2259{
2261 ExcAccessToUninitializedField(
2262 "update_jacobian_pushed_forward_2nd_derivatives"));
2263 Assert(present_cell.is_initialized(), ExcNotReinited());
2264 return this->mapping_output.jacobian_pushed_forward_2nd_derivatives[q_point];
2265}
2266
2267
2268
2269template <int dim, int spacedim>
2270inline const std::vector<Tensor<4, spacedim>> &
2272{
2274 ExcAccessToUninitializedField(
2275 "update_jacobian_pushed_forward_2nd_derivatives"));
2276 Assert(present_cell.is_initialized(), ExcNotReinited());
2277 return this->mapping_output.jacobian_pushed_forward_2nd_derivatives;
2278}
2279
2280
2281
2282template <int dim, int spacedim>
2285 const unsigned int q_point) const
2286{
2287 Assert(this->update_flags & update_jacobian_3rd_derivatives,
2288 ExcAccessToUninitializedField("update_jacobian_3rd_derivatives"));
2289 Assert(present_cell.is_initialized(), ExcNotReinited());
2290 return this->mapping_output.jacobian_3rd_derivatives[q_point];
2291}
2292
2293
2294
2295template <int dim, int spacedim>
2296inline const std::vector<DerivativeForm<4, dim, spacedim>> &
2298{
2299 Assert(this->update_flags & update_jacobian_3rd_derivatives,
2300 ExcAccessToUninitializedField("update_jacobian_3rd_derivatives"));
2301 Assert(present_cell.is_initialized(), ExcNotReinited());
2302 return this->mapping_output.jacobian_3rd_derivatives;
2303}
2304
2305
2306
2307template <int dim, int spacedim>
2308inline const Tensor<5, spacedim> &
2310 const unsigned int q_point) const
2311{
2313 ExcAccessToUninitializedField(
2314 "update_jacobian_pushed_forward_3rd_derivatives"));
2315 Assert(present_cell.is_initialized(), ExcNotReinited());
2316 return this->mapping_output.jacobian_pushed_forward_3rd_derivatives[q_point];
2317}
2318
2319
2320
2321template <int dim, int spacedim>
2322inline const std::vector<Tensor<5, spacedim>> &
2324{
2326 ExcAccessToUninitializedField(
2327 "update_jacobian_pushed_forward_3rd_derivatives"));
2328 Assert(present_cell.is_initialized(), ExcNotReinited());
2329 return this->mapping_output.jacobian_pushed_forward_3rd_derivatives;
2330}
2331
2332
2333
2334template <int dim, int spacedim>
2335inline const std::vector<DerivativeForm<1, spacedim, dim>> &
2337{
2338 Assert(this->update_flags & update_inverse_jacobians,
2339 ExcAccessToUninitializedField("update_inverse_jacobians"));
2340 Assert(present_cell.is_initialized(), ExcNotReinited());
2341 return this->mapping_output.inverse_jacobians;
2342}
2343
2344
2345
2346template <int dim, int spacedim>
2349{
2350 return {0U, dofs_per_cell};
2351}
2352
2353
2354
2355template <int dim, int spacedim>
2358 const unsigned int start_dof_index) const
2359{
2360 Assert(start_dof_index <= dofs_per_cell,
2361 ExcIndexRange(start_dof_index, 0, dofs_per_cell + 1));
2362 return {start_dof_index, dofs_per_cell};
2363}
2364
2365
2366
2367template <int dim, int spacedim>
2370 const unsigned int end_dof_index) const
2371{
2372 Assert(end_dof_index < dofs_per_cell,
2373 ExcIndexRange(end_dof_index, 0, dofs_per_cell));
2374 return {0U, end_dof_index + 1};
2375}
2376
2377
2378
2379template <int dim, int spacedim>
2382{
2383 return {0U, n_quadrature_points};
2384}
2385
2386
2387
2388template <int dim, int spacedim>
2389inline const Point<spacedim> &
2390FEValuesBase<dim, spacedim>::quadrature_point(const unsigned int q_point) const
2391{
2392 Assert(this->update_flags & update_quadrature_points,
2393 ExcAccessToUninitializedField("update_quadrature_points"));
2394 AssertIndexRange(q_point, this->mapping_output.quadrature_points.size());
2395 Assert(present_cell.is_initialized(), ExcNotReinited());
2396
2397 return this->mapping_output.quadrature_points[q_point];
2398}
2399
2400
2401
2402template <int dim, int spacedim>
2403inline double
2404FEValuesBase<dim, spacedim>::JxW(const unsigned int q_point) const
2405{
2406 Assert(this->update_flags & update_JxW_values,
2407 ExcAccessToUninitializedField("update_JxW_values"));
2408 AssertIndexRange(q_point, this->mapping_output.JxW_values.size());
2409 Assert(present_cell.is_initialized(), ExcNotReinited());
2410
2411 return this->mapping_output.JxW_values[q_point];
2412}
2413
2414
2415
2416template <int dim, int spacedim>
2418FEValuesBase<dim, spacedim>::jacobian(const unsigned int q_point) const
2419{
2420 Assert(this->update_flags & update_jacobians,
2421 ExcAccessToUninitializedField("update_jacobians"));
2422 AssertIndexRange(q_point, this->mapping_output.jacobians.size());
2423 Assert(present_cell.is_initialized(), ExcNotReinited());
2424
2425 return this->mapping_output.jacobians[q_point];
2426}
2427
2428
2429
2430template <int dim, int spacedim>
2432FEValuesBase<dim, spacedim>::jacobian_grad(const unsigned int q_point) const
2433{
2434 Assert(this->update_flags & update_jacobian_grads,
2435 ExcAccessToUninitializedField("update_jacobians_grads"));
2436 AssertIndexRange(q_point, this->mapping_output.jacobian_grads.size());
2437 Assert(present_cell.is_initialized(), ExcNotReinited());
2438
2439 return this->mapping_output.jacobian_grads[q_point];
2440}
2441
2442
2443
2444template <int dim, int spacedim>
2446FEValuesBase<dim, spacedim>::inverse_jacobian(const unsigned int q_point) const
2447{
2448 Assert(this->update_flags & update_inverse_jacobians,
2449 ExcAccessToUninitializedField("update_inverse_jacobians"));
2450 AssertIndexRange(q_point, this->mapping_output.inverse_jacobians.size());
2451 Assert(present_cell.is_initialized(), ExcNotReinited());
2452
2453 return this->mapping_output.inverse_jacobians[q_point];
2454}
2455
2456
2457
2458template <int dim, int spacedim>
2459inline const Tensor<1, spacedim> &
2460FEValuesBase<dim, spacedim>::normal_vector(const unsigned int q_point) const
2461{
2462 Assert(this->update_flags & update_normal_vectors,
2464 "update_normal_vectors")));
2465 AssertIndexRange(q_point, this->mapping_output.normal_vectors.size());
2466 Assert(present_cell.is_initialized(), ExcNotReinited());
2467
2468 return this->mapping_output.normal_vectors[q_point];
2469}
2470
2471#endif
2472
2474
2475#endif
typename LevelSelector::cell_iterator level_cell_iterator
types::global_dof_index n_dofs_for_dof_handler() const
std::optional< std::variant< typename Triangulation< dim, spacedim >::cell_iterator, typename DoFHandler< dim, spacedim >::cell_iterator, typename DoFHandler< dim, spacedim >::level_cell_iterator > > cell
void get_interpolated_dof_values(const ReadVector< Number > &in, Vector< Number > &out) const
CellSimilarity::Similarity cell_similarity
const Tensor< 5, spacedim > & jacobian_pushed_forward_3rd_derivative(const unsigned int q_point) const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > dof_indices_ending_at(const unsigned int end_dof_index) const
const FEValuesViews::Vector< dim, spacedim > & operator[](const FEValuesExtractors::Vector &vector) const
const std::vector< double > & get_JxW_values() const
const DerivativeForm< 1, dim, spacedim > & jacobian(const unsigned int q_point) const
::internal::FEValuesViews::Cache< dim, spacedim > fe_values_views_cache
FEValuesBase(const FEValuesBase &)=delete
const FEValuesViews::Scalar< dim, spacedim > & operator[](const FEValuesExtractors::Scalar &scalar) const
Triangulation< dim, spacedim >::cell_iterator get_cell() const
boost::signals2::connection tria_listener_mesh_transform
void get_function_values(const ReadVector< Number > &fe_function, std::vector< Number > &values) const
const std::vector< Point< spacedim > > & get_quadrature_points() const
const std::vector< DerivativeForm< 1, dim, spacedim > > & get_jacobians() const
const DerivativeForm< 2, dim, spacedim > & jacobian_grad(const unsigned int q_point) const
internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > mapping_output
virtual ~FEValuesBase() override
const std::vector< DerivativeForm< 3, dim, spacedim > > & get_jacobian_2nd_derivatives() const
UpdateFlags get_update_flags() const
const FEValuesViews::Tensor< 2, dim, spacedim > & operator[](const FEValuesExtractors::Tensor< 2 > &tensor) const
const unsigned int dofs_per_cell
void check_cell_similarity(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
Tensor< 3, spacedim > shape_3rd_derivative_component(const unsigned int i, const unsigned int q_point, const unsigned int component) const
CellIteratorWrapper present_cell
const Tensor< 1, spacedim > & normal_vector(const unsigned int q_point) const
UpdateFlags update_flags
static constexpr unsigned int space_dimension
Tensor< 2, spacedim > shape_hessian_component(const unsigned int i, const unsigned int q_point, const unsigned int component) const
void get_function_hessians(const ReadVector< Number > &fe_function, std::vector< Tensor< 2, spacedim, Number > > &hessians) const
void always_allow_check_for_cell_similarity(const bool allow)
void get_function_laplacians(const ReadVector< Number > &fe_function, std::vector< Number > &laplacians) const
const Point< spacedim > & quadrature_point(const unsigned int q_point) const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > dof_indices_starting_at(const unsigned int start_dof_index) const
double shape_value_component(const unsigned int i, const unsigned int q_point, const unsigned int component) const
const Mapping< dim, spacedim > & get_mapping() const
const unsigned int n_quadrature_points
CellSimilarity::Similarity get_cell_similarity() const
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors() const
std::size_t memory_consumption() const
const Tensor< 3, spacedim > & jacobian_pushed_forward_grad(const unsigned int q_point) const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > dof_indices() const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > quadrature_point_indices() const
void get_function_third_derivatives(const ReadVector< Number > &fe_function, std::vector< Tensor< 3, spacedim, Number > > &third_derivatives) const
const Tensor< 2, spacedim > & shape_hessian(const unsigned int i, const unsigned int q_point) const
const std::vector< Tensor< 5, spacedim > > & get_jacobian_pushed_forward_3rd_derivatives() const
const Tensor< 4, spacedim > & jacobian_pushed_forward_2nd_derivative(const unsigned int q_point) const
const Tensor< 3, spacedim > & shape_3rd_derivative(const unsigned int i, const unsigned int q_point) const
UpdateFlags compute_update_flags(const UpdateFlags update_flags) const
std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > mapping_data
void get_function_gradients(const ReadVector< Number > &fe_function, std::vector< Tensor< 1, spacedim, Number > > &gradients) const
const std::vector< DerivativeForm< 2, dim, spacedim > > & get_jacobian_grads() const
boost::signals2::connection tria_listener_refinement
const std::vector< Tensor< 3, spacedim > > & get_jacobian_pushed_forward_grads() const
const DerivativeForm< 4, dim, spacedim > & jacobian_3rd_derivative(const unsigned int q_point) const
const DerivativeForm< 1, spacedim, dim > & inverse_jacobian(const unsigned int q_point) const
void invalidate_present_cell()
const std::vector< Tensor< 4, spacedim > > & get_jacobian_pushed_forward_2nd_derivatives() const
const DerivativeForm< 3, dim, spacedim > & jacobian_2nd_derivative(const unsigned int q_point) const
const ObserverPointer< const Mapping< dim, spacedim >, FEValuesBase< dim, spacedim > > mapping
Tensor< 1, spacedim > shape_grad_component(const unsigned int i, const unsigned int q_point, const unsigned int component) const
bool check_for_cell_similarity_allowed
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > finite_element_output
const FEValuesViews::SymmetricTensor< 2, dim, spacedim > & operator[](const FEValuesExtractors::SymmetricTensor< 2 > &tensor) const
const ObserverPointer< const FiniteElement< dim, spacedim >, FEValuesBase< dim, spacedim > > fe
const Tensor< 1, spacedim > & shape_grad(const unsigned int i, const unsigned int q_point) const
const std::vector< DerivativeForm< 4, dim, spacedim > > & get_jacobian_3rd_derivatives() const
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > fe_data
const FiniteElement< dim, spacedim > & get_fe() const
double JxW(const unsigned int q_point) const
const double & shape_value(const unsigned int i, const unsigned int q_point) const
void maybe_invalidate_previous_present_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
FEValuesBase & operator=(const FEValuesBase &)=delete
const std::vector< DerivativeForm< 1, spacedim, dim > > & get_inverse_jacobians() const
const unsigned int max_n_quadrature_points
static constexpr unsigned int dimension
Abstract base class for mapping classes.
Definition mapping.h:318
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
static ::ExceptionBase & ExcAccessToUninitializedField(std::string arg1)
#define Assert(cond, exc)
static ::ExceptionBase & ExcNotReinited()
static ::ExceptionBase & ExcNeedsDoFHandler()
static ::ExceptionBase & ExcShapeFunctionNotPrimitive(int arg1)
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:489
static ::ExceptionBase & ExcIndexRange(std::size_t arg1, std::size_t arg2, std::size_t arg3)
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:511
static ::ExceptionBase & ExcFEDontMatch()
static ::ExceptionBase & ExcFENotPrimitive()
typename ActiveSelector::cell_iterator cell_iterator
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition tria.h:1556
UpdateFlags
@ update_jacobian_pushed_forward_2nd_derivatives
@ update_jacobian_pushed_forward_grads
@ update_hessians
Second derivatives of shape functions.
@ update_jacobian_3rd_derivatives
@ update_values
Shape function values.
@ update_jacobian_grads
Gradient of volume element.
@ update_normal_vectors
Normal vectors.
@ update_3rd_derivatives
Third derivatives of shape functions.
@ update_JxW_values
Transformed quadrature weights.
@ update_jacobians
Volume element.
@ update_inverse_jacobians
Volume element.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
@ update_jacobian_pushed_forward_3rd_derivatives
@ update_jacobian_2nd_derivatives
boost::integer_range< IncrementableType > iota_view
Definition iota_view.h:45
std::vector< Lazy<::FEValuesViews::Vector< dim, spacedim > > > vectors
std::vector< Lazy<::FEValuesViews::SymmetricTensor< 2, dim, spacedim > > > symmetric_second_order_tensors
std::vector< Lazy<::FEValuesViews::Tensor< 2, dim, spacedim > > > second_order_tensors
std::vector< Lazy<::FEValuesViews::Scalar< dim, spacedim > > > scalars