Reference documentation for deal.II version 9.4.1
\(\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.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1998 - 2022 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_fe_values_h
17#define dealii_fe_values_h
18
19
20#include <deal.II/base/config.h>
21
24#include <deal.II/base/point.h>
29
32
33#include <deal.II/fe/fe.h>
36#include <deal.II/fe/mapping.h>
37
38#include <deal.II/grid/tria.h>
40
42
43#include <algorithm>
44#include <memory>
45#include <type_traits>
46
47
48// dummy include in order to have the
49// definition of PetscScalar available
50// without including other PETSc stuff
51#ifdef DEAL_II_WITH_PETSC
52# include <petsc.h>
53#endif
54
56
57// Forward declaration
58#ifndef DOXYGEN
59template <int dim, int spacedim = dim>
60class FEValuesBase;
61#endif
62
63namespace internal
64{
69 template <int dim, class NumberType = double>
70 struct CurlType;
71
78 template <class NumberType>
80 {
82 };
83
90 template <class NumberType>
92 {
94 };
95
102 template <class NumberType>
104 {
106 };
107} // namespace internal
108
109
110
133{
145 template <int dim, int spacedim = dim>
146 class Scalar
147 {
148 public:
154 using value_type = double;
155
162
169
176
183 template <typename Number>
185
192 template <typename Number>
195
202 template <typename Number>
205
212 template <typename Number>
215
222 template <typename Number>
225
232 template <typename Number>
234 {
240 typename ProductType<Number,
241 typename Scalar<dim, spacedim>::value_type>::type;
242
247 using gradient_type = typename ProductType<
248 Number,
250
256 typename ProductType<Number,
257 typename Scalar<dim, spacedim>::value_type>::type;
258
263 using hessian_type = typename ProductType<
264 Number,
266
272 Number,
274 };
275
281 {
291
300 unsigned int row_index;
301 };
302
306 Scalar();
307
313 Scalar(const FEValuesBase<dim, spacedim> &fe_values_base,
314 const unsigned int component);
315
320 Scalar(const Scalar<dim, spacedim> &) = delete;
321
325 // NOLINTNEXTLINE OSX does not compile with noexcept
327
331 ~Scalar() = default;
332
337 Scalar &
339
343 Scalar &
344 operator=(Scalar<dim, spacedim> &&) noexcept = default;
345
360 value(const unsigned int shape_function, const unsigned int q_point) const;
361
373 gradient(const unsigned int shape_function,
374 const unsigned int q_point) const;
375
387 hessian(const unsigned int shape_function,
388 const unsigned int q_point) const;
389
401 third_derivative(const unsigned int shape_function,
402 const unsigned int q_point) const;
403
421 template <class InputVector>
422 void
424 const InputVector &fe_function,
425 std::vector<solution_value_type<typename InputVector::value_type>>
426 &values) const;
427
462 template <class InputVector>
463 void
465 const InputVector &dof_values,
466 std::vector<solution_value_type<typename InputVector::value_type>>
467 &values) const;
468
486 template <class InputVector>
487 void
489 const InputVector &fe_function,
490 std::vector<solution_gradient_type<typename InputVector::value_type>>
491 &gradients) const;
492
499 template <class InputVector>
500 void
502 const InputVector &dof_values,
503 std::vector<solution_gradient_type<typename InputVector::value_type>>
504 &gradients) const;
505
523 template <class InputVector>
524 void
526 const InputVector &fe_function,
527 std::vector<solution_hessian_type<typename InputVector::value_type>>
528 &hessians) const;
529
536 template <class InputVector>
537 void
539 const InputVector &dof_values,
540 std::vector<solution_hessian_type<typename InputVector::value_type>>
541 &hessians) const;
542
543
562 template <class InputVector>
563 void
565 const InputVector &fe_function,
566 std::vector<solution_laplacian_type<typename InputVector::value_type>>
567 &laplacians) const;
568
575 template <class InputVector>
576 void
578 const InputVector &dof_values,
579 std::vector<solution_laplacian_type<typename InputVector::value_type>>
580 &laplacians) const;
581
582
601 template <class InputVector>
602 void
604 const InputVector &fe_function,
605 std::vector<
606 solution_third_derivative_type<typename InputVector::value_type>>
607 &third_derivatives) const;
608
615 template <class InputVector>
616 void
618 const InputVector &dof_values,
619 std::vector<
620 solution_third_derivative_type<typename InputVector::value_type>>
621 &third_derivatives) const;
622
623
624 private:
628 const SmartPointer<const FEValuesBase<dim, spacedim>> fe_values;
629
634 const unsigned int component;
635
640 };
641
642
643
673 template <int dim, int spacedim = dim>
674 class Vector
675 {
676 public:
683
693
705
711 using divergence_type = double;
712
719 using curl_type = typename ::internal::CurlType<spacedim>::type;
720
727
734
741 template <typename Number>
743
750 template <typename Number>
753
760 template <typename Number>
763
770 template <typename Number>
773
780 template <typename Number>
783
790 template <typename Number>
792
799 template <typename Number>
802
809 template <typename Number>
812
819 template <typename Number>
821 {
827 typename ProductType<Number,
828 typename Vector<dim, spacedim>::value_type>::type;
829
834 using gradient_type = typename ProductType<
835 Number,
837
843 Number,
845
850 using divergence_type = typename ProductType<
851 Number,
853
859 typename ProductType<Number,
860 typename Vector<dim, spacedim>::value_type>::type;
861
866 using curl_type =
867 typename ProductType<Number,
868 typename Vector<dim, spacedim>::curl_type>::type;
869
874 using hessian_type = typename ProductType<
875 Number,
877
883 Number,
885 };
886
892 {
902
912 unsigned int row_index[spacedim];
913
924 };
925
929 Vector();
930
939 Vector(const FEValuesBase<dim, spacedim> &fe_values_base,
940 const unsigned int first_vector_component);
941
946 Vector(const Vector<dim, spacedim> &) = delete;
947
951 // NOLINTNEXTLINE OSX does not compile with noexcept
953
957 ~Vector() = default;
958
963 Vector &
965
969 // NOLINTNEXTLINE OSX does not compile with noexcept
970 Vector &
971 operator=(Vector<dim, spacedim> &&) = default; // NOLINT
972
990 value(const unsigned int shape_function, const unsigned int q_point) const;
991
1006 gradient(const unsigned int shape_function,
1007 const unsigned int q_point) const;
1008
1025 symmetric_gradient(const unsigned int shape_function,
1026 const unsigned int q_point) const;
1027
1039 divergence(const unsigned int shape_function,
1040 const unsigned int q_point) const;
1041
1062 curl_type
1063 curl(const unsigned int shape_function, const unsigned int q_point) const;
1064
1076 hessian(const unsigned int shape_function,
1077 const unsigned int q_point) const;
1078
1090 third_derivative(const unsigned int shape_function,
1091 const unsigned int q_point) const;
1092
1110 template <class InputVector>
1111 void
1113 const InputVector &fe_function,
1115 &values) const;
1116
1151 template <class InputVector>
1152 void
1154 const InputVector &dof_values,
1156 &values) const;
1157
1175 template <class InputVector>
1176 void
1178 const InputVector &fe_function,
1180 &gradients) const;
1181
1188 template <class InputVector>
1189 void
1191 const InputVector &dof_values,
1193 &gradients) const;
1194
1218 template <class InputVector>
1219 void
1220 get_function_symmetric_gradients(
1221 const InputVector &fe_function,
1222 std::vector<
1224 &symmetric_gradients) const;
1225
1232 template <class InputVector>
1233 void
1234 get_function_symmetric_gradients_from_local_dof_values(
1235 const InputVector &dof_values,
1236 std::vector<
1238 &symmetric_gradients) const;
1239
1258 template <class InputVector>
1259 void
1260 get_function_divergences(
1261 const InputVector &fe_function,
1263 &divergences) const;
1264
1271 template <class InputVector>
1272 void
1273 get_function_divergences_from_local_dof_values(
1274 const InputVector &dof_values,
1276 &divergences) const;
1277
1296 template <class InputVector>
1297 void
1298 get_function_curls(
1299 const InputVector &fe_function,
1301 const;
1302
1309 template <class InputVector>
1310 void
1311 get_function_curls_from_local_dof_values(
1312 const InputVector &dof_values,
1314 const;
1315
1333 template <class InputVector>
1334 void
1336 const InputVector &fe_function,
1338 &hessians) const;
1339
1346 template <class InputVector>
1347 void
1349 const InputVector &dof_values,
1351 &hessians) const;
1352
1371 template <class InputVector>
1372 void
1374 const InputVector &fe_function,
1376 &laplacians) const;
1377
1384 template <class InputVector>
1385 void
1387 const InputVector &dof_values,
1389 &laplacians) const;
1390
1409 template <class InputVector>
1410 void
1412 const InputVector &fe_function,
1413 std::vector<
1415 &third_derivatives) const;
1416
1423 template <class InputVector>
1424 void
1426 const InputVector &dof_values,
1427 std::vector<
1429 &third_derivatives) const;
1430
1431 private:
1436
1441 const unsigned int first_vector_component;
1442
1446 std::vector<ShapeFunctionData> shape_function_data;
1447 };
1448
1449
1450 template <int rank, int dim, int spacedim = dim>
1452
1475 template <int dim, int spacedim>
1476 class SymmetricTensor<2, dim, spacedim>
1477 {
1478 public:
1486
1497
1504 template <typename Number>
1506
1513 template <typename Number>
1516
1517
1524 template <typename Number>
1525 struct DEAL_II_DEPRECATED OutputType
1526 {
1531 using value_type = typename ProductType<
1532 Number,
1534
1540 Number,
1542 };
1543
1548 struct ShapeFunctionData
1549 {
1558 bool is_nonzero_shape_function_component
1559 [value_type::n_independent_components];
1560
1570 unsigned int row_index[value_type::n_independent_components];
1571
1581
1586 };
1587
1592
1602 SymmetricTensor(const FEValuesBase<dim, spacedim> &fe_values_base,
1603 const unsigned int first_tensor_component);
1604
1610
1614 // NOLINTNEXTLINE OSX does not compile with noexcept
1616
1623
1629
1648 value(const unsigned int shape_function, const unsigned int q_point) const;
1649
1664 divergence(const unsigned int shape_function,
1665 const unsigned int q_point) const;
1666
1684 template <class InputVector>
1685 void
1686 get_function_values(
1687 const InputVector &fe_function,
1688 std::vector<solution_value_type<typename InputVector::value_type>>
1689 &values) const;
1690
1725 template <class InputVector>
1726 void
1727 get_function_values_from_local_dof_values(
1728 const InputVector &dof_values,
1729 std::vector<solution_value_type<typename InputVector::value_type>>
1730 &values) const;
1731
1753 template <class InputVector>
1754 void
1755 get_function_divergences(
1756 const InputVector &fe_function,
1757 std::vector<solution_divergence_type<typename InputVector::value_type>>
1758 &divergences) const;
1759
1766 template <class InputVector>
1767 void
1768 get_function_divergences_from_local_dof_values(
1769 const InputVector &dof_values,
1770 std::vector<solution_divergence_type<typename InputVector::value_type>>
1771 &divergences) const;
1772
1773 private:
1777 const SmartPointer<const FEValuesBase<dim, spacedim>> fe_values;
1778
1783 const unsigned int first_tensor_component;
1784
1788 std::vector<ShapeFunctionData> shape_function_data;
1789 };
1790
1791
1792 template <int rank, int dim, int spacedim = dim>
1793 class Tensor;
1794
1813 template <int dim, int spacedim>
1814 class Tensor<2, dim, spacedim>
1815 {
1816 public:
1822
1827
1833
1840 template <typename Number>
1842
1849 template <typename Number>
1852
1859 template <typename Number>
1862
1863
1870 template <typename Number>
1871 struct DEAL_II_DEPRECATED OutputType
1872 {
1877 using value_type = typename ProductType<
1878 Number,
1880
1886 Number,
1888
1893 using gradient_type = typename ProductType<
1894 Number,
1896 };
1897
1902 struct ShapeFunctionData
1903 {
1912 bool is_nonzero_shape_function_component
1913 [value_type::n_independent_components];
1914
1924 unsigned int row_index[value_type::n_independent_components];
1925
1935
1940 };
1941
1945 Tensor();
1946
1952
1956 // NOLINTNEXTLINE OSX does not compile with noexcept
1958
1962 ~Tensor() = default;
1963
1973 Tensor(const FEValuesBase<dim, spacedim> &fe_values_base,
1974 const unsigned int first_tensor_component);
1975
1976
1981 Tensor &
1983
1987 Tensor &
1988 operator=(Tensor<2, dim, spacedim> &&) = default; // NOLINT
1989
2007 value(const unsigned int shape_function, const unsigned int q_point) const;
2008
2023 divergence(const unsigned int shape_function,
2024 const unsigned int q_point) const;
2025
2040 gradient(const unsigned int shape_function,
2041 const unsigned int q_point) const;
2042
2060 template <class InputVector>
2061 void
2062 get_function_values(
2063 const InputVector &fe_function,
2065 &values) const;
2066
2101 template <class InputVector>
2102 void
2103 get_function_values_from_local_dof_values(
2104 const InputVector &dof_values,
2106 &values) const;
2107
2129 template <class InputVector>
2130 void
2131 get_function_divergences(
2132 const InputVector &fe_function,
2134 &divergences) const;
2135
2142 template <class InputVector>
2143 void
2144 get_function_divergences_from_local_dof_values(
2145 const InputVector &dof_values,
2147 &divergences) const;
2148
2165 template <class InputVector>
2166 void
2167 get_function_gradients(
2168 const InputVector &fe_function,
2170 &gradients) const;
2171
2178 template <class InputVector>
2179 void
2180 get_function_gradients_from_local_dof_values(
2181 const InputVector &dof_values,
2183 &gradients) const;
2184
2185 private:
2190
2195 const unsigned int first_tensor_component;
2196
2200 std::vector<ShapeFunctionData> shape_function_data;
2201 };
2202
2203} // namespace FEValuesViews
2204
2205
2206namespace internal
2207{
2209 {
2214 template <int dim, int spacedim, typename Extractor>
2216 {};
2217
2225 template <int dim, int spacedim>
2226 struct ViewType<dim, spacedim, FEValuesExtractors::Scalar>
2227 {
2228 using type = typename ::FEValuesViews::Scalar<dim, spacedim>;
2229 };
2230
2238 template <int dim, int spacedim>
2239 struct ViewType<dim, spacedim, FEValuesExtractors::Vector>
2240 {
2241 using type = typename ::FEValuesViews::Vector<dim, spacedim>;
2242 };
2243
2251 template <int dim, int spacedim, int rank>
2252 struct ViewType<dim, spacedim, FEValuesExtractors::Tensor<rank>>
2253 {
2254 using type = typename ::FEValuesViews::Tensor<rank, dim, spacedim>;
2255 };
2256
2264 template <int dim, int spacedim, int rank>
2265 struct ViewType<dim, spacedim, FEValuesExtractors::SymmetricTensor<rank>>
2266 {
2267 using type =
2268 typename ::FEValuesViews::SymmetricTensor<rank, dim, spacedim>;
2269 };
2270
2278 template <int dim, int spacedim>
2279 struct Cache
2280 {
2285 std::vector<::FEValuesViews::Scalar<dim, spacedim>> scalars;
2286 std::vector<::FEValuesViews::Vector<dim, spacedim>> vectors;
2287 std::vector<::FEValuesViews::SymmetricTensor<2, dim, spacedim>>
2289 std::vector<::FEValuesViews::Tensor<2, dim, spacedim>>
2291
2295 Cache(const FEValuesBase<dim, spacedim> &fe_values);
2296 };
2297 } // namespace FEValuesViews
2298} // namespace internal
2299
2300namespace FEValuesViews
2301{
2306 template <int dim, int spacedim, typename Extractor>
2307 using View = typename ::internal::FEValuesViews::
2308 ViewType<dim, spacedim, Extractor>::type;
2309} // namespace FEValuesViews
2310
2311
2411template <int dim, int spacedim>
2413{
2414public:
2418 static constexpr unsigned int dimension = dim;
2419
2423 static constexpr unsigned int space_dimension = spacedim;
2424
2432 const unsigned int n_quadrature_points;
2433
2443 const unsigned int max_n_quadrature_points;
2444
2450 const unsigned int dofs_per_cell;
2451
2452
2460 FEValuesBase(const unsigned int n_q_points,
2461 const unsigned int dofs_per_cell,
2462 const UpdateFlags update_flags,
2463 const Mapping<dim, spacedim> & mapping,
2465
2470 FEValuesBase &
2471 operator=(const FEValuesBase &) = delete;
2472
2477 FEValuesBase(const FEValuesBase &) = delete;
2478
2482 virtual ~FEValuesBase() override;
2483
2484
2488
2489
2510 const double &
2511 shape_value(const unsigned int function_no,
2512 const unsigned int point_no) const;
2513
2534 double
2535 shape_value_component(const unsigned int function_no,
2536 const unsigned int point_no,
2537 const unsigned int component) const;
2538
2564 const Tensor<1, spacedim> &
2565 shape_grad(const unsigned int function_no,
2566 const unsigned int quadrature_point) const;
2567
2585 shape_grad_component(const unsigned int function_no,
2586 const unsigned int point_no,
2587 const unsigned int component) const;
2588
2608 const Tensor<2, spacedim> &
2609 shape_hessian(const unsigned int function_no,
2610 const unsigned int point_no) const;
2611
2629 shape_hessian_component(const unsigned int function_no,
2630 const unsigned int point_no,
2631 const unsigned int component) const;
2632
2652 const Tensor<3, spacedim> &
2653 shape_3rd_derivative(const unsigned int function_no,
2654 const unsigned int point_no) const;
2655
2673 shape_3rd_derivative_component(const unsigned int function_no,
2674 const unsigned int point_no,
2675 const unsigned int component) const;
2676
2678
2680
2717 template <class InputVector>
2718 void
2719 get_function_values(
2720 const InputVector & fe_function,
2721 std::vector<typename InputVector::value_type> &values) const;
2722
2736 template <class InputVector>
2737 void
2738 get_function_values(
2739 const InputVector & fe_function,
2740 std::vector<Vector<typename InputVector::value_type>> &values) const;
2741
2798 template <class InputVector>
2799 void
2800 get_function_values(
2801 const InputVector & fe_function,
2803 std::vector<typename InputVector::value_type> & values) const;
2804
2813 template <class InputVector>
2814 void
2815 get_function_values(
2816 const InputVector & fe_function,
2818 std::vector<Vector<typename InputVector::value_type>> &values) const;
2819
2820
2842 template <class InputVector>
2843 void
2844 get_function_values(
2845 const InputVector & fe_function,
2847 ArrayView<std::vector<typename InputVector::value_type>> values,
2848 const bool quadrature_points_fastest) const;
2849
2851
2853
2890 template <class InputVector>
2891 void
2892 get_function_gradients(
2893 const InputVector &fe_function,
2895 &gradients) const;
2896
2913 template <class InputVector>
2914 void
2915 get_function_gradients(
2916 const InputVector &fe_function,
2917 std::vector<
2919 &gradients) const;
2920
2929 template <class InputVector>
2930 void
2931 get_function_gradients(
2932 const InputVector & fe_function,
2935 &gradients) const;
2936
2945 template <class InputVector>
2946 void
2947 get_function_gradients(
2948 const InputVector & fe_function,
2950 ArrayView<
2952 gradients,
2953 const bool quadrature_points_fastest = false) const;
2954
2956
2959
2960
2998 template <class InputVector>
2999 void
3000 get_function_hessians(
3001 const InputVector &fe_function,
3003 &hessians) const;
3004
3022 template <class InputVector>
3023 void
3024 get_function_hessians(
3025 const InputVector &fe_function,
3026 std::vector<
3028 & hessians,
3029 const bool quadrature_points_fastest = false) const;
3030
3039 template <class InputVector>
3040 void
3041 get_function_hessians(
3042 const InputVector & fe_function,
3045 &hessians) const;
3046
3055 template <class InputVector>
3056 void
3057 get_function_hessians(
3058 const InputVector & fe_function,
3060 ArrayView<
3062 hessians,
3063 const bool quadrature_points_fastest = false) const;
3064
3105 template <class InputVector>
3106 void
3107 get_function_laplacians(
3108 const InputVector & fe_function,
3109 std::vector<typename InputVector::value_type> &laplacians) const;
3110
3130 template <class InputVector>
3131 void
3132 get_function_laplacians(
3133 const InputVector & fe_function,
3134 std::vector<Vector<typename InputVector::value_type>> &laplacians) const;
3135
3144 template <class InputVector>
3145 void
3146 get_function_laplacians(
3147 const InputVector & fe_function,
3149 std::vector<typename InputVector::value_type> & laplacians) const;
3150
3159 template <class InputVector>
3160 void
3161 get_function_laplacians(
3162 const InputVector & fe_function,
3164 std::vector<Vector<typename InputVector::value_type>> &laplacians) const;
3165
3174 template <class InputVector>
3175 void
3176 get_function_laplacians(
3177 const InputVector & fe_function,
3179 std::vector<std::vector<typename InputVector::value_type>> &laplacians,
3180 const bool quadrature_points_fastest = false) const;
3181
3183
3185
3224 template <class InputVector>
3225 void
3226 get_function_third_derivatives(
3227 const InputVector &fe_function,
3229 &third_derivatives) const;
3230
3249 template <class InputVector>
3250 void
3251 get_function_third_derivatives(
3252 const InputVector &fe_function,
3253 std::vector<
3255 & third_derivatives,
3256 const bool quadrature_points_fastest = false) const;
3257
3266 template <class InputVector>
3267 void
3268 get_function_third_derivatives(
3269 const InputVector & fe_function,
3272 &third_derivatives) const;
3273
3282 template <class InputVector>
3283 void
3284 get_function_third_derivatives(
3285 const InputVector & fe_function,
3287 ArrayView<
3289 third_derivatives,
3290 const bool quadrature_points_fastest = false) const;
3292
3294
3295
3321
3355 dof_indices_starting_at(const unsigned int start_dof_index) const;
3356
3388 dof_indices_ending_at(const unsigned int end_dof_index) const;
3389
3391
3393
3394
3417
3423 const Point<spacedim> &
3424 quadrature_point(const unsigned int q) const;
3425
3431 const std::vector<Point<spacedim>> &
3433
3449 double
3450 JxW(const unsigned int quadrature_point) const;
3451
3455 const std::vector<double> &
3457
3465 jacobian(const unsigned int quadrature_point) const;
3466
3473 const std::vector<DerivativeForm<1, dim, spacedim>> &
3475
3484 jacobian_grad(const unsigned int quadrature_point) const;
3485
3492 const std::vector<DerivativeForm<2, dim, spacedim>> &
3494
3503 const Tensor<3, spacedim> &
3504 jacobian_pushed_forward_grad(const unsigned int quadrature_point) const;
3505
3512 const std::vector<Tensor<3, spacedim>> &
3514
3523 jacobian_2nd_derivative(const unsigned int quadrature_point) const;
3524
3531 const std::vector<DerivativeForm<3, dim, spacedim>> &
3533
3543 const Tensor<4, spacedim> &
3545 const unsigned int quadrature_point) const;
3546
3553 const std::vector<Tensor<4, spacedim>> &
3555
3565 jacobian_3rd_derivative(const unsigned int quadrature_point) const;
3566
3573 const std::vector<DerivativeForm<4, dim, spacedim>> &
3575
3585 const Tensor<5, spacedim> &
3587 const unsigned int quadrature_point) const;
3588
3595 const std::vector<Tensor<5, spacedim>> &
3597
3605 inverse_jacobian(const unsigned int quadrature_point) const;
3606
3613 const std::vector<DerivativeForm<1, spacedim, dim>> &
3615
3635 const Tensor<1, spacedim> &
3636 normal_vector(const unsigned int i) const;
3637
3645 const std::vector<Tensor<1, spacedim>> &
3646 get_normal_vectors() const;
3647
3649
3651
3652
3663
3674
3686
3687
3698
3700
3702
3703
3709
3714 get_fe() const;
3715
3721
3726 get_cell() const;
3727
3734 get_cell_similarity() const;
3735
3740 std::size_t
3741 memory_consumption() const;
3743
3744
3752 ExcAccessToUninitializedField,
3753 std::string,
3754 << "You are requesting information from an FEValues/FEFaceValues/FESubfaceValues "
3755 << "object for which this kind of information has not been computed. What "
3756 << "information these objects compute is determined by the update_* flags you "
3757 << "pass to the constructor. Here, the operation you are attempting requires "
3758 << "the <" << arg1
3759 << "> flag to be set, but it was apparently not specified "
3760 << "upon construction.");
3761
3767 DeclExceptionMsg(ExcNotReinited,
3768 "FEValues object is not reinit'ed to any cell");
3769
3777 ExcFEDontMatch,
3778 "The FiniteElement you provided to FEValues and the FiniteElement that belongs "
3779 "to the DoFHandler that provided the cell iterator do not match.");
3785 DeclException1(ExcShapeFunctionNotPrimitive,
3786 int,
3787 << "The shape function with index " << arg1
3788 << " is not primitive, i.e. it is vector-valued and "
3789 << "has more than one non-zero vector component. This "
3790 << "function cannot be called for these shape functions. "
3791 << "Maybe you want to use the same function with the "
3792 << "_component suffix?");
3793
3801 ExcFENotPrimitive,
3802 "The given FiniteElement is not a primitive element but the requested operation "
3803 "only works for those. See FiniteElement::is_primitive() for more information.");
3804
3805protected:
3813 {
3814 public:
3816 ExcNeedsDoFHandler,
3817 "You have previously called the FEValues::reinit() function with a "
3818 "cell iterator of type Triangulation<dim,spacedim>::cell_iterator. However, "
3819 "when you do this, you cannot call some functions in the FEValues "
3820 "class, such as the get_function_values/gradients/hessians/third_derivatives "
3821 "functions. If you need these functions, then you need to call "
3822 "FEValues::reinit() with an iterator type that allows to extract "
3823 "degrees of freedom, such as DoFHandler<dim,spacedim>::cell_iterator.");
3824
3829
3833 template <bool lda>
3836
3841 const typename Triangulation<dim, spacedim>::cell_iterator &cell);
3842
3846 bool
3847 is_initialized() const;
3848
3855 operator typename Triangulation<dim, spacedim>::cell_iterator() const;
3856
3862 n_dofs_for_dof_handler() const;
3863
3868 template <typename VectorType>
3869 void
3870 get_interpolated_dof_values(
3871 const VectorType & in,
3873
3878 void
3879 get_interpolated_dof_values(const IndexSet & in,
3880 Vector<IndexSet::value_type> &out) const;
3881
3882 private:
3887 };
3888
3895
3903 boost::signals2::connection tria_listener_refinement;
3904
3912 boost::signals2::connection tria_listener_mesh_transform;
3913
3919 void
3920 invalidate_present_cell();
3921
3931 void
3932 maybe_invalidate_previous_present_cell(
3933 const typename Triangulation<dim, spacedim>::cell_iterator &cell);
3934
3940
3946 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
3948
3955
3956
3964
3970 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
3972
3978 spacedim>
3980
3981
3986
3996 compute_update_flags(const UpdateFlags update_flags) const;
3997
4004
4010 void
4011 check_cell_similarity(
4012 const typename Triangulation<dim, spacedim>::cell_iterator &cell);
4013
4014private:
4019
4020 // Make the view classes friends of this class, since they access internal
4021 // data.
4022 template <int, int>
4024 template <int, int>
4026 template <int, int, int>
4028 template <int, int, int>
4030};
4031
4032
4033
4043template <int dim, int spacedim = dim>
4044class FEValues : public FEValuesBase<dim, spacedim>
4045{
4046public:
4051 static constexpr unsigned int integral_dimension = dim;
4052
4059 const Quadrature<dim> & quadrature,
4060 const UpdateFlags update_flags);
4061
4070 const hp::QCollection<dim> & quadrature,
4071 const UpdateFlags update_flags);
4072
4079 const Quadrature<dim> & quadrature,
4080 const UpdateFlags update_flags);
4081
4089 const hp::QCollection<dim> & quadrature,
4090 const UpdateFlags update_flags);
4091
4098 template <bool level_dof_access>
4099 void
4102
4116 void
4118
4123 const Quadrature<dim> &
4125
4130 std::size_t
4132
4149
4150private:
4155
4159 void
4160 initialize(const UpdateFlags update_flags);
4161
4168 void
4170};
4171
4172
4182template <int dim, int spacedim = dim>
4183class FEFaceValuesBase : public FEValuesBase<dim, spacedim>
4184{
4185public:
4190 static constexpr unsigned int integral_dimension = dim - 1;
4191
4203 FEFaceValuesBase(const unsigned int dofs_per_cell,
4204 const UpdateFlags update_flags,
4205 const Mapping<dim, spacedim> & mapping,
4207 const Quadrature<dim - 1> & quadrature);
4208
4215 FEFaceValuesBase(const unsigned int dofs_per_cell,
4216 const UpdateFlags update_flags,
4217 const Mapping<dim, spacedim> & mapping,
4219 const hp::QCollection<dim - 1> & quadrature);
4220
4228 const Tensor<1, spacedim> &
4229 boundary_form(const unsigned int i) const;
4230
4237 const std::vector<Tensor<1, spacedim>> &
4239
4244 unsigned int
4246
4251 unsigned int
4253
4258 const Quadrature<dim - 1> &
4260
4265 std::size_t
4267
4268protected:
4273 unsigned int present_face_no;
4274
4280
4285};
4286
4287
4288
4302template <int dim, int spacedim = dim>
4303class FEFaceValues : public FEFaceValuesBase<dim, spacedim>
4304{
4305public:
4310 static constexpr unsigned int dimension = dim;
4311
4312 static constexpr unsigned int space_dimension = spacedim;
4313
4318 static constexpr unsigned int integral_dimension = dim - 1;
4319
4326 const Quadrature<dim - 1> & quadrature,
4327 const UpdateFlags update_flags);
4328
4337 const hp::QCollection<dim - 1> & quadrature,
4338 const UpdateFlags update_flags);
4339
4346 const Quadrature<dim - 1> & quadrature,
4347 const UpdateFlags update_flags);
4348
4356 const hp::QCollection<dim - 1> & quadrature,
4357 const UpdateFlags update_flags);
4358
4363 template <bool level_dof_access>
4364 void
4367 const unsigned int face_no);
4368
4375 template <bool level_dof_access>
4376 void
4379 const typename Triangulation<dim, spacedim>::face_iterator & face);
4380
4394 void
4396 const unsigned int face_no);
4397
4398 /*
4399 * Reinitialize the gradients, Jacobi determinants, etc for the given face
4400 * on a given cell of type "iterator into a Triangulation object", and the
4401 * given finite element. Since iterators into a triangulation alone only
4402 * convey information about the geometry of a cell, but not about degrees of
4403 * freedom possibly associated with this cell, you will not be able to call
4404 * some functions of this class if they need information about degrees of
4405 * freedom. These functions are, above all, the
4406 * <tt>get_function_value/gradients/hessians/third_derivatives</tt>
4407 * functions. If you want to call these functions, you have to call the @p
4408 * reinit variants that take iterators into DoFHandler or other DoF handler
4409 * type objects.
4410 *
4411 * @note @p face must be one of @p cell's face iterators.
4412 */
4413 void
4415 const typename Triangulation<dim, spacedim>::face_iterator &face);
4416
4433
4434private:
4438 void
4439 initialize(const UpdateFlags update_flags);
4440
4447 void
4448 do_reinit(const unsigned int face_no);
4449};
4450
4451
4468template <int dim, int spacedim = dim>
4469class FESubfaceValues : public FEFaceValuesBase<dim, spacedim>
4470{
4471public:
4475 static constexpr unsigned int dimension = dim;
4476
4480 static constexpr unsigned int space_dimension = spacedim;
4481
4486 static constexpr unsigned int integral_dimension = dim - 1;
4487
4494 const Quadrature<dim - 1> & face_quadrature,
4495 const UpdateFlags update_flags);
4496
4505 const hp::QCollection<dim - 1> & face_quadrature,
4506 const UpdateFlags update_flags);
4507
4514 const Quadrature<dim - 1> & face_quadrature,
4515 const UpdateFlags update_flags);
4516
4524 const hp::QCollection<dim - 1> & face_quadrature,
4525 const UpdateFlags update_flags);
4526
4533 template <bool level_dof_access>
4534 void
4537 const unsigned int face_no,
4538 const unsigned int subface_no);
4539
4544 template <bool level_dof_access>
4545 void
4548 const typename Triangulation<dim, spacedim>::face_iterator & face,
4549 const typename Triangulation<dim, spacedim>::face_iterator &subface);
4550
4564 void
4566 const unsigned int face_no,
4567 const unsigned int subface_no);
4568
4588 void
4591 const typename Triangulation<dim, spacedim>::face_iterator &subface);
4592
4609
4615 DeclException0(ExcReinitCalledWithBoundaryFace);
4616
4622 DeclException0(ExcFaceHasNoSubfaces);
4623
4624private:
4628 void
4629 initialize(const UpdateFlags update_flags);
4630
4637 void
4638 do_reinit(const unsigned int face_no, const unsigned int subface_no);
4639};
4640
4641
4642#ifndef DOXYGEN
4643
4644
4645/*------------------------ Inline functions: namespace FEValuesViews --------*/
4646
4647namespace FEValuesViews
4648{
4649 template <int dim, int spacedim>
4650 inline typename Scalar<dim, spacedim>::value_type
4651 Scalar<dim, spacedim>::value(const unsigned int shape_function,
4652 const unsigned int q_point) const
4653 {
4654 AssertIndexRange(shape_function, fe_values->fe->n_dofs_per_cell());
4655 Assert(
4656 fe_values->update_flags & update_values,
4658 "update_values"))));
4659
4660 // an adaptation of the FEValuesBase::shape_value_component function
4661 // except that here we know the component as fixed and we have
4662 // pre-computed and cached a bunch of information. See the comments there.
4663 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
4664 return fe_values->finite_element_output.shape_values(
4665 shape_function_data[shape_function].row_index, q_point);
4666 else
4667 return 0;
4668 }
4669
4670
4671
4672 template <int dim, int spacedim>
4673 inline typename Scalar<dim, spacedim>::gradient_type
4674 Scalar<dim, spacedim>::gradient(const unsigned int shape_function,
4675 const unsigned int q_point) const
4676 {
4677 AssertIndexRange(shape_function, fe_values->fe->n_dofs_per_cell());
4678 Assert(fe_values->update_flags & update_gradients,
4680 "update_gradients")));
4681
4682 // an adaptation of the FEValuesBase::shape_grad_component
4683 // function except that here we know the component as fixed and we have
4684 // pre-computed and cached a bunch of information. See the comments there.
4685 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
4686 return fe_values->finite_element_output
4687 .shape_gradients[shape_function_data[shape_function].row_index]
4688 [q_point];
4689 else
4690 return gradient_type();
4691 }
4692
4693
4694
4695 template <int dim, int spacedim>
4696 inline typename Scalar<dim, spacedim>::hessian_type
4697 Scalar<dim, spacedim>::hessian(const unsigned int shape_function,
4698 const unsigned int q_point) const
4699 {
4700 AssertIndexRange(shape_function, fe_values->fe->n_dofs_per_cell());
4701 Assert(fe_values->update_flags & update_hessians,
4703 "update_hessians")));
4704
4705 // an adaptation of the FEValuesBase::shape_hessian_component
4706 // function except that here we know the component as fixed and we have
4707 // pre-computed and cached a bunch of information. See the comments there.
4708 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
4709 return fe_values->finite_element_output
4710 .shape_hessians[shape_function_data[shape_function].row_index][q_point];
4711 else
4712 return hessian_type();
4713 }
4714
4715
4716
4717 template <int dim, int spacedim>
4718 inline typename Scalar<dim, spacedim>::third_derivative_type
4719 Scalar<dim, spacedim>::third_derivative(const unsigned int shape_function,
4720 const unsigned int q_point) const
4721 {
4722 AssertIndexRange(shape_function, fe_values->fe->n_dofs_per_cell());
4723 Assert(fe_values->update_flags & update_3rd_derivatives,
4725 "update_3rd_derivatives")));
4726
4727 // an adaptation of the FEValuesBase::shape_3rdderivative_component
4728 // function except that here we know the component as fixed and we have
4729 // pre-computed and cached a bunch of information. See the comments there.
4730 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
4731 return fe_values->finite_element_output
4732 .shape_3rd_derivatives[shape_function_data[shape_function].row_index]
4733 [q_point];
4734 else
4735 return third_derivative_type();
4736 }
4737
4738
4739
4740 template <int dim, int spacedim>
4741 inline typename Vector<dim, spacedim>::value_type
4742 Vector<dim, spacedim>::value(const unsigned int shape_function,
4743 const unsigned int q_point) const
4744 {
4745 AssertIndexRange(shape_function, fe_values->fe->n_dofs_per_cell());
4746 Assert(fe_values->update_flags & update_values,
4748 "update_values")));
4749
4750 // same as for the scalar case except that we have one more index
4751 const int snc =
4752 shape_function_data[shape_function].single_nonzero_component;
4753 if (snc == -2)
4754 return value_type();
4755 else if (snc != -1)
4756 {
4757 value_type return_value;
4758 return_value[shape_function_data[shape_function]
4759 .single_nonzero_component_index] =
4760 fe_values->finite_element_output.shape_values(snc, q_point);
4761 return return_value;
4762 }
4763 else
4764 {
4765 value_type return_value;
4766 for (unsigned int d = 0; d < dim; ++d)
4767 if (shape_function_data[shape_function]
4768 .is_nonzero_shape_function_component[d])
4769 return_value[d] = fe_values->finite_element_output.shape_values(
4770 shape_function_data[shape_function].row_index[d], q_point);
4771
4772 return return_value;
4773 }
4774 }
4775
4776
4777
4778 template <int dim, int spacedim>
4780 Vector<dim, spacedim>::gradient(const unsigned int shape_function,
4781 const unsigned int q_point) const
4782 {
4783 AssertIndexRange(shape_function, fe_values->fe->n_dofs_per_cell());
4784 Assert(fe_values->update_flags & update_gradients,
4786 "update_gradients")));
4787
4788 // same as for the scalar case except that we have one more index
4789 const int snc =
4790 shape_function_data[shape_function].single_nonzero_component;
4791 if (snc == -2)
4792 return gradient_type();
4793 else if (snc != -1)
4794 {
4795 gradient_type return_value;
4796 return_value[shape_function_data[shape_function]
4797 .single_nonzero_component_index] =
4798 fe_values->finite_element_output.shape_gradients[snc][q_point];
4799 return return_value;
4800 }
4801 else
4802 {
4803 gradient_type return_value;
4804 for (unsigned int d = 0; d < dim; ++d)
4805 if (shape_function_data[shape_function]
4806 .is_nonzero_shape_function_component[d])
4807 return_value[d] =
4808 fe_values->finite_element_output.shape_gradients
4809 [shape_function_data[shape_function].row_index[d]][q_point];
4810
4811 return return_value;
4812 }
4813 }
4814
4815
4816
4817 template <int dim, int spacedim>
4819 Vector<dim, spacedim>::divergence(const unsigned int shape_function,
4820 const unsigned int q_point) const
4821 {
4822 // this function works like in the case above
4823 AssertIndexRange(shape_function, fe_values->fe->n_dofs_per_cell());
4824 Assert(fe_values->update_flags & update_gradients,
4826 "update_gradients")));
4827
4828 // same as for the scalar case except that we have one more index
4829 const int snc =
4830 shape_function_data[shape_function].single_nonzero_component;
4831 if (snc == -2)
4832 return divergence_type();
4833 else if (snc != -1)
4834 return fe_values->finite_element_output
4835 .shape_gradients[snc][q_point][shape_function_data[shape_function]
4836 .single_nonzero_component_index];
4837 else
4838 {
4839 divergence_type return_value = 0;
4840 for (unsigned int d = 0; d < dim; ++d)
4841 if (shape_function_data[shape_function]
4842 .is_nonzero_shape_function_component[d])
4843 return_value +=
4844 fe_values->finite_element_output.shape_gradients
4845 [shape_function_data[shape_function].row_index[d]][q_point][d];
4846
4847 return return_value;
4848 }
4849 }
4850
4851
4852
4853 template <int dim, int spacedim>
4854 inline typename Vector<dim, spacedim>::curl_type
4855 Vector<dim, spacedim>::curl(const unsigned int shape_function,
4856 const unsigned int q_point) const
4857 {
4858 // this function works like in the case above
4859
4860 AssertIndexRange(shape_function, fe_values->fe->n_dofs_per_cell());
4861 Assert(fe_values->update_flags & update_gradients,
4863 "update_gradients")));
4864 // same as for the scalar case except that we have one more index
4865 const int snc =
4866 shape_function_data[shape_function].single_nonzero_component;
4867
4868 if (snc == -2)
4869 return curl_type();
4870
4871 else
4872 switch (dim)
4873 {
4874 case 1:
4875 {
4876 Assert(false,
4877 ExcMessage(
4878 "Computing the curl in 1d is not a useful operation"));
4879 return curl_type();
4880 }
4881
4882 case 2:
4883 {
4884 if (snc != -1)
4885 {
4886 curl_type return_value;
4887
4888 // the single nonzero component can only be zero or one in 2d
4889 if (shape_function_data[shape_function]
4890 .single_nonzero_component_index == 0)
4891 return_value[0] =
4892 -1.0 * fe_values->finite_element_output
4893 .shape_gradients[snc][q_point][1];
4894 else
4895 return_value[0] = fe_values->finite_element_output
4896 .shape_gradients[snc][q_point][0];
4897
4898 return return_value;
4899 }
4900
4901 else
4902 {
4903 curl_type return_value;
4904
4905 return_value[0] = 0.0;
4906
4907 if (shape_function_data[shape_function]
4908 .is_nonzero_shape_function_component[0])
4909 return_value[0] -=
4910 fe_values->finite_element_output
4911 .shape_gradients[shape_function_data[shape_function]
4912 .row_index[0]][q_point][1];
4913
4914 if (shape_function_data[shape_function]
4915 .is_nonzero_shape_function_component[1])
4916 return_value[0] +=
4917 fe_values->finite_element_output
4918 .shape_gradients[shape_function_data[shape_function]
4919 .row_index[1]][q_point][0];
4920
4921 return return_value;
4922 }
4923 }
4924
4925 case 3:
4926 {
4927 if (snc != -1)
4928 {
4929 curl_type return_value;
4930
4931 switch (shape_function_data[shape_function]
4932 .single_nonzero_component_index)
4933 {
4934 case 0:
4935 {
4936 return_value[0] = 0;
4937 return_value[1] = fe_values->finite_element_output
4938 .shape_gradients[snc][q_point][2];
4939 return_value[2] =
4940 -1.0 * fe_values->finite_element_output
4941 .shape_gradients[snc][q_point][1];
4942 return return_value;
4943 }
4944
4945 case 1:
4946 {
4947 return_value[0] =
4948 -1.0 * fe_values->finite_element_output
4949 .shape_gradients[snc][q_point][2];
4950 return_value[1] = 0;
4951 return_value[2] = fe_values->finite_element_output
4952 .shape_gradients[snc][q_point][0];
4953 return return_value;
4954 }
4955
4956 default:
4957 {
4958 return_value[0] = fe_values->finite_element_output
4959 .shape_gradients[snc][q_point][1];
4960 return_value[1] =
4961 -1.0 * fe_values->finite_element_output
4962 .shape_gradients[snc][q_point][0];
4963 return_value[2] = 0;
4964 return return_value;
4965 }
4966 }
4967 }
4968
4969 else
4970 {
4971 curl_type return_value;
4972
4973 for (unsigned int i = 0; i < dim; ++i)
4974 return_value[i] = 0.0;
4975
4976 if (shape_function_data[shape_function]
4977 .is_nonzero_shape_function_component[0])
4978 {
4979 return_value[1] +=
4980 fe_values->finite_element_output
4981 .shape_gradients[shape_function_data[shape_function]
4982 .row_index[0]][q_point][2];
4983 return_value[2] -=
4984 fe_values->finite_element_output
4985 .shape_gradients[shape_function_data[shape_function]
4986 .row_index[0]][q_point][1];
4987 }
4988
4989 if (shape_function_data[shape_function]
4990 .is_nonzero_shape_function_component[1])
4991 {
4992 return_value[0] -=
4993 fe_values->finite_element_output
4994 .shape_gradients[shape_function_data[shape_function]
4995 .row_index[1]][q_point][2];
4996 return_value[2] +=
4997 fe_values->finite_element_output
4998 .shape_gradients[shape_function_data[shape_function]
4999 .row_index[1]][q_point][0];
5000 }
5001
5002 if (shape_function_data[shape_function]
5003 .is_nonzero_shape_function_component[2])
5004 {
5005 return_value[0] +=
5006 fe_values->finite_element_output
5007 .shape_gradients[shape_function_data[shape_function]
5008 .row_index[2]][q_point][1];
5009 return_value[1] -=
5010 fe_values->finite_element_output
5011 .shape_gradients[shape_function_data[shape_function]
5012 .row_index[2]][q_point][0];
5013 }
5014
5015 return return_value;
5016 }
5017 }
5018 }
5019 // should not end up here
5020 Assert(false, ExcInternalError());
5021 return curl_type();
5022 }
5023
5024
5025
5026 template <int dim, int spacedim>
5028 Vector<dim, spacedim>::hessian(const unsigned int shape_function,
5029 const unsigned int q_point) const
5030 {
5031 // this function works like in the case above
5032 AssertIndexRange(shape_function, fe_values->fe->n_dofs_per_cell());
5033 Assert(fe_values->update_flags & update_hessians,
5035 "update_hessians")));
5036
5037 // same as for the scalar case except that we have one more index
5038 const int snc =
5039 shape_function_data[shape_function].single_nonzero_component;
5040 if (snc == -2)
5041 return hessian_type();
5042 else if (snc != -1)
5043 {
5044 hessian_type return_value;
5045 return_value[shape_function_data[shape_function]
5046 .single_nonzero_component_index] =
5047 fe_values->finite_element_output.shape_hessians[snc][q_point];
5048 return return_value;
5049 }
5050 else
5051 {
5052 hessian_type return_value;
5053 for (unsigned int d = 0; d < dim; ++d)
5054 if (shape_function_data[shape_function]
5055 .is_nonzero_shape_function_component[d])
5056 return_value[d] =
5057 fe_values->finite_element_output.shape_hessians
5058 [shape_function_data[shape_function].row_index[d]][q_point];
5059
5060 return return_value;
5061 }
5062 }
5063
5064
5065
5066 template <int dim, int spacedim>
5068 Vector<dim, spacedim>::third_derivative(const unsigned int shape_function,
5069 const unsigned int q_point) const
5070 {
5071 // this function works like in the case above
5072 AssertIndexRange(shape_function, fe_values->fe->n_dofs_per_cell());
5073 Assert(fe_values->update_flags & update_3rd_derivatives,
5075 "update_3rd_derivatives")));
5076
5077 // same as for the scalar case except that we have one more index
5078 const int snc =
5079 shape_function_data[shape_function].single_nonzero_component;
5080 if (snc == -2)
5081 return third_derivative_type();
5082 else if (snc != -1)
5083 {
5084 third_derivative_type return_value;
5085 return_value[shape_function_data[shape_function]
5086 .single_nonzero_component_index] =
5087 fe_values->finite_element_output.shape_3rd_derivatives[snc][q_point];
5088 return return_value;
5089 }
5090 else
5091 {
5092 third_derivative_type return_value;
5093 for (unsigned int d = 0; d < dim; ++d)
5094 if (shape_function_data[shape_function]
5095 .is_nonzero_shape_function_component[d])
5096 return_value[d] =
5097 fe_values->finite_element_output.shape_3rd_derivatives
5098 [shape_function_data[shape_function].row_index[d]][q_point];
5099
5100 return return_value;
5101 }
5102 }
5103
5104
5105
5106 namespace internal
5107 {
5112 inline ::SymmetricTensor<2, 1>
5113 symmetrize_single_row(const unsigned int n, const ::Tensor<1, 1> &t)
5114 {
5115 AssertIndexRange(n, 1);
5116 (void)n;
5117
5118 return {{t[0]}};
5119 }
5120
5121
5122
5123 inline ::SymmetricTensor<2, 2>
5124 symmetrize_single_row(const unsigned int n, const ::Tensor<1, 2> &t)
5125 {
5126 switch (n)
5127 {
5128 case 0:
5129 {
5130 return {{t[0], 0, t[1] / 2}};
5131 }
5132 case 1:
5133 {
5134 return {{0, t[1], t[0] / 2}};
5135 }
5136 default:
5137 {
5138 AssertIndexRange(n, 2);
5139 return {};
5140 }
5141 }
5142 }
5143
5144
5145
5146 inline ::SymmetricTensor<2, 3>
5147 symmetrize_single_row(const unsigned int n, const ::Tensor<1, 3> &t)
5148 {
5149 switch (n)
5150 {
5151 case 0:
5152 {
5153 return {{t[0], 0, 0, t[1] / 2, t[2] / 2, 0}};
5154 }
5155 case 1:
5156 {
5157 return {{0, t[1], 0, t[0] / 2, 0, t[2] / 2}};
5158 }
5159 case 2:
5160 {
5161 return {{0, 0, t[2], 0, t[0] / 2, t[1] / 2}};
5162 }
5163 default:
5164 {
5165 AssertIndexRange(n, 3);
5166 return {};
5167 }
5168 }
5169 }
5170 } // namespace internal
5171
5172
5173
5174 template <int dim, int spacedim>
5176 Vector<dim, spacedim>::symmetric_gradient(const unsigned int shape_function,
5177 const unsigned int q_point) const
5178 {
5179 AssertIndexRange(shape_function, fe_values->fe->n_dofs_per_cell());
5180 Assert(fe_values->update_flags & update_gradients,
5182 "update_gradients")));
5183
5184 // same as for the scalar case except that we have one more index
5185 const int snc =
5186 shape_function_data[shape_function].single_nonzero_component;
5187 if (snc == -2)
5188 return symmetric_gradient_type();
5189 else if (snc != -1)
5190 return internal::symmetrize_single_row(
5191 shape_function_data[shape_function].single_nonzero_component_index,
5192 fe_values->finite_element_output.shape_gradients[snc][q_point]);
5193 else
5194 {
5195 gradient_type return_value;
5196 for (unsigned int d = 0; d < dim; ++d)
5197 if (shape_function_data[shape_function]
5198 .is_nonzero_shape_function_component[d])
5199 return_value[d] =
5200 fe_values->finite_element_output.shape_gradients
5201 [shape_function_data[shape_function].row_index[d]][q_point];
5202
5203 return symmetrize(return_value);
5204 }
5205 }
5206
5207
5208
5209 template <int dim, int spacedim>
5211 SymmetricTensor<2, dim, spacedim>::value(const unsigned int shape_function,
5212 const unsigned int q_point) const
5213 {
5214 AssertIndexRange(shape_function, fe_values->fe->n_dofs_per_cell());
5215 Assert(fe_values->update_flags & update_values,
5217 "update_values")));
5218
5219 // similar to the vector case where we have more then one index and we need
5220 // to convert between unrolled and component indexing for tensors
5221 const int snc =
5222 shape_function_data[shape_function].single_nonzero_component;
5223
5224 if (snc == -2)
5225 {
5226 // shape function is zero for the selected components
5227 return value_type();
5228 }
5229 else if (snc != -1)
5230 {
5231 value_type return_value;
5232 const unsigned int comp =
5233 shape_function_data[shape_function].single_nonzero_component_index;
5234 return_value[value_type::unrolled_to_component_indices(comp)] =
5235 fe_values->finite_element_output.shape_values(snc, q_point);
5236 return return_value;
5237 }
5238 else
5239 {
5240 value_type return_value;
5241 for (unsigned int d = 0; d < value_type::n_independent_components; ++d)
5242 if (shape_function_data[shape_function]
5243 .is_nonzero_shape_function_component[d])
5244 return_value[value_type::unrolled_to_component_indices(d)] =
5245 fe_values->finite_element_output.shape_values(
5246 shape_function_data[shape_function].row_index[d], q_point);
5247 return return_value;
5248 }
5249 }
5250
5251
5252
5253 template <int dim, int spacedim>
5256 const unsigned int shape_function,
5257 const unsigned int q_point) const
5258 {
5259 AssertIndexRange(shape_function, fe_values->fe->n_dofs_per_cell());
5260 Assert(fe_values->update_flags & update_gradients,
5262 "update_gradients")));
5263
5264 const int snc =
5265 shape_function_data[shape_function].single_nonzero_component;
5266
5267 if (snc == -2)
5268 {
5269 // shape function is zero for the selected components
5270 return divergence_type();
5271 }
5272 else if (snc != -1)
5273 {
5274 // we have a single non-zero component when the symmetric tensor is
5275 // represented in unrolled form. this implies we potentially have
5276 // two non-zero components when represented in component form! we
5277 // will only have one non-zero entry if the non-zero component lies on
5278 // the diagonal of the tensor.
5279 //
5280 // the divergence of a second-order tensor is a first order tensor.
5281 //
5282 // assume the second-order tensor is A with components A_{ij}. then
5283 // A_{ij} = A_{ji} and there is only one (if diagonal) or two non-zero
5284 // entries in the tensorial representation. define the
5285 // divergence as:
5286 // b_i \dealcoloneq \dfrac{\partial phi_{ij}}{\partial x_j}.
5287 // (which is incidentally also
5288 // b_j \dealcoloneq \dfrac{\partial phi_{ij}}{\partial x_i}).
5289 // In both cases, a sum is implied.
5290 //
5291 // Now, we know the nonzero component in unrolled form: it is indicated
5292 // by 'snc'. we can figure out which tensor components belong to this:
5293 const unsigned int comp =
5294 shape_function_data[shape_function].single_nonzero_component_index;
5295 const unsigned int ii =
5296 value_type::unrolled_to_component_indices(comp)[0];
5297 const unsigned int jj =
5298 value_type::unrolled_to_component_indices(comp)[1];
5299
5300 // given the form of the divergence above, if ii=jj there is only a
5301 // single nonzero component of the full tensor and the gradient
5302 // equals
5303 // b_ii \dealcoloneq \dfrac{\partial phi_{ii,ii}}{\partial x_ii}.
5304 // all other entries of 'b' are zero
5305 //
5306 // on the other hand, if ii!=jj, then there are two nonzero entries in
5307 // the full tensor and
5308 // b_ii \dealcoloneq \dfrac{\partial phi_{ii,jj}}{\partial x_ii}.
5309 // b_jj \dealcoloneq \dfrac{\partial phi_{ii,jj}}{\partial x_jj}.
5310 // again, all other entries of 'b' are zero
5311 const ::Tensor<1, spacedim> &phi_grad =
5312 fe_values->finite_element_output.shape_gradients[snc][q_point];
5313
5314 divergence_type return_value;
5315 return_value[ii] = phi_grad[jj];
5316
5317 if (ii != jj)
5318 return_value[jj] = phi_grad[ii];
5319
5320 return return_value;
5321 }
5322 else
5323 {
5324 Assert(false, ExcNotImplemented());
5325 divergence_type return_value;
5326 return return_value;
5327 }
5328 }
5329
5330
5331
5332 template <int dim, int spacedim>
5334 Tensor<2, dim, spacedim>::value(const unsigned int shape_function,
5335 const unsigned int q_point) const
5336 {
5337 AssertIndexRange(shape_function, fe_values->fe->n_dofs_per_cell());
5338 Assert(fe_values->update_flags & update_values,
5340 "update_values")));
5341
5342 // similar to the vector case where we have more then one index and we need
5343 // to convert between unrolled and component indexing for tensors
5344 const int snc =
5345 shape_function_data[shape_function].single_nonzero_component;
5346
5347 if (snc == -2)
5348 {
5349 // shape function is zero for the selected components
5350 return value_type();
5351 }
5352 else if (snc != -1)
5353 {
5354 value_type return_value;
5355 const unsigned int comp =
5356 shape_function_data[shape_function].single_nonzero_component_index;
5357 const TableIndices<2> indices =
5359 return_value[indices] =
5360 fe_values->finite_element_output.shape_values(snc, q_point);
5361 return return_value;
5362 }
5363 else
5364 {
5365 value_type return_value;
5366 for (unsigned int d = 0; d < dim * dim; ++d)
5367 if (shape_function_data[shape_function]
5368 .is_nonzero_shape_function_component[d])
5369 {
5370 const TableIndices<2> indices =
5372 return_value[indices] =
5373 fe_values->finite_element_output.shape_values(
5374 shape_function_data[shape_function].row_index[d], q_point);
5375 }
5376 return return_value;
5377 }
5378 }
5379
5380
5381
5382 template <int dim, int spacedim>
5384 Tensor<2, dim, spacedim>::divergence(const unsigned int shape_function,
5385 const unsigned int q_point) const
5386 {
5387 AssertIndexRange(shape_function, fe_values->fe->n_dofs_per_cell());
5388 Assert(fe_values->update_flags & update_gradients,
5390 "update_gradients")));
5391
5392 const int snc =
5393 shape_function_data[shape_function].single_nonzero_component;
5394
5395 if (snc == -2)
5396 {
5397 // shape function is zero for the selected components
5398 return divergence_type();
5399 }
5400 else if (snc != -1)
5401 {
5402 // we have a single non-zero component when the tensor is
5403 // represented in unrolled form.
5404 //
5405 // the divergence of a second-order tensor is a first order tensor.
5406 //
5407 // assume the second-order tensor is A with components A_{ij},
5408 // then divergence is d_i := \frac{\partial A_{ij}}{\partial x_j}
5409 //
5410 // Now, we know the nonzero component in unrolled form: it is indicated
5411 // by 'snc'. we can figure out which tensor components belong to this:
5412 const unsigned int comp =
5413 shape_function_data[shape_function].single_nonzero_component_index;
5414 const TableIndices<2> indices =
5416 const unsigned int ii = indices[0];
5417 const unsigned int jj = indices[1];
5418
5419 const ::Tensor<1, spacedim> &phi_grad =
5420 fe_values->finite_element_output.shape_gradients[snc][q_point];
5421
5422 divergence_type return_value;
5423 // note that we contract \nabla from the right
5424 return_value[ii] = phi_grad[jj];
5425
5426 return return_value;
5427 }
5428 else
5429 {
5430 Assert(false, ExcNotImplemented());
5431 divergence_type return_value;
5432 return return_value;
5433 }
5434 }
5435
5436
5437
5438 template <int dim, int spacedim>
5440 Tensor<2, dim, spacedim>::gradient(const unsigned int shape_function,
5441 const unsigned int q_point) const
5442 {
5443 AssertIndexRange(shape_function, fe_values->fe->n_dofs_per_cell());
5444 Assert(fe_values->update_flags & update_gradients,
5446 "update_gradients")));
5447
5448 const int snc =
5449 shape_function_data[shape_function].single_nonzero_component;
5450
5451 if (snc == -2)
5452 {
5453 // shape function is zero for the selected components
5454 return gradient_type();
5455 }
5456 else if (snc != -1)
5457 {
5458 // we have a single non-zero component when the tensor is
5459 // represented in unrolled form.
5460 //
5461 // the gradient of a second-order tensor is a third order tensor.
5462 //
5463 // assume the second-order tensor is A with components A_{ij},
5464 // then gradient is B_{ijk} := \frac{\partial A_{ij}}{\partial x_k}
5465 //
5466 // Now, we know the nonzero component in unrolled form: it is indicated
5467 // by 'snc'. we can figure out which tensor components belong to this:
5468 const unsigned int comp =
5469 shape_function_data[shape_function].single_nonzero_component_index;
5470 const TableIndices<2> indices =
5472 const unsigned int ii = indices[0];
5473 const unsigned int jj = indices[1];
5474
5475 const ::Tensor<1, spacedim> &phi_grad =
5476 fe_values->finite_element_output.shape_gradients[snc][q_point];
5477
5478 gradient_type return_value;
5479 return_value[ii][jj] = phi_grad;
5480
5481 return return_value;
5482 }
5483 else
5484 {
5485 Assert(false, ExcNotImplemented());
5486 gradient_type return_value;
5487 return return_value;
5488 }
5489 }
5490
5491} // namespace FEValuesViews
5492
5493
5494
5495/*---------------------- Inline functions: FEValuesBase ---------------------*/
5496
5497
5498
5499template <int dim, int spacedim>
5500template <bool lda>
5504 : initialized(true)
5505 , cell(cell)
5506 , dof_handler(&cell->get_dof_handler())
5507 , level_dof_access(lda)
5508{}
5509
5510
5511
5512template <int dim, int spacedim>
5515 const FEValuesExtractors::Scalar &scalar) const
5516{
5517 AssertIndexRange(scalar.component, fe_values_views_cache.scalars.size());
5518
5519 return fe_values_views_cache.scalars[scalar.component];
5520}
5521
5522
5523
5524template <int dim, int spacedim>
5527 const FEValuesExtractors::Vector &vector) const
5528{
5530 fe_values_views_cache.vectors.size());
5531
5532 return fe_values_views_cache.vectors[vector.first_vector_component];
5533}
5534
5535
5536
5537template <int dim, int spacedim>
5540 const FEValuesExtractors::SymmetricTensor<2> &tensor) const
5541{
5542 Assert(
5543 tensor.first_tensor_component <
5544 fe_values_views_cache.symmetric_second_order_tensors.size(),
5546 0,
5547 fe_values_views_cache.symmetric_second_order_tensors.size()));
5548
5549 return fe_values_views_cache
5550 .symmetric_second_order_tensors[tensor.first_tensor_component];
5551}
5552
5553
5554
5555template <int dim, int spacedim>
5558 const FEValuesExtractors::Tensor<2> &tensor) const
5559{
5561 fe_values_views_cache.second_order_tensors.size());
5562
5563 return fe_values_views_cache
5564 .second_order_tensors[tensor.first_tensor_component];
5565}
5566
5567
5568
5569template <int dim, int spacedim>
5570inline const double &
5571FEValuesBase<dim, spacedim>::shape_value(const unsigned int i,
5572 const unsigned int j) const
5573{
5574 AssertIndexRange(i, fe->n_dofs_per_cell());
5575 Assert(this->update_flags & update_values,
5576 ExcAccessToUninitializedField("update_values"));
5577 Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
5578 Assert(present_cell.is_initialized(), ExcNotReinited());
5579 // if the entire FE is primitive,
5580 // then we can take a short-cut:
5581 if (fe->is_primitive())
5582 return this->finite_element_output.shape_values(i, j);
5583 else
5584 {
5585 // otherwise, use the mapping
5586 // between shape function
5587 // numbers and rows. note that
5588 // by the assertions above, we
5589 // know that this particular
5590 // shape function is primitive,
5591 // so we can call
5592 // system_to_component_index
5593 const unsigned int row =
5594 this->finite_element_output
5595 .shape_function_to_row_table[i * fe->n_components() +
5596 fe->system_to_component_index(i).first];
5597 return this->finite_element_output.shape_values(row, j);
5598 }
5599}
5600
5601
5602
5603template <int dim, int spacedim>
5604inline double
5606 const unsigned int i,
5607 const unsigned int j,
5608 const unsigned int component) const
5609{
5610 AssertIndexRange(i, fe->n_dofs_per_cell());
5611 Assert(this->update_flags & update_values,
5612 ExcAccessToUninitializedField("update_values"));
5613 AssertIndexRange(component, fe->n_components());
5614 Assert(present_cell.is_initialized(), ExcNotReinited());
5615
5616 // check whether the shape function
5617 // is non-zero at all within
5618 // this component:
5619 if (fe->get_nonzero_components(i)[component] == false)
5620 return 0;
5621
5622 // look up the right row in the
5623 // table and take the data from
5624 // there
5625 const unsigned int row =
5626 this->finite_element_output
5627 .shape_function_to_row_table[i * fe->n_components() + component];
5628 return this->finite_element_output.shape_values(row, j);
5629}
5630
5631
5632
5633template <int dim, int spacedim>
5634inline const Tensor<1, spacedim> &
5635FEValuesBase<dim, spacedim>::shape_grad(const unsigned int i,
5636 const unsigned int j) const
5637{
5638 AssertIndexRange(i, fe->n_dofs_per_cell());
5639 Assert(this->update_flags & update_gradients,
5640 ExcAccessToUninitializedField("update_gradients"));
5641 Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
5642 Assert(present_cell.is_initialized(), ExcNotReinited());
5643 // if the entire FE is primitive,
5644 // then we can take a short-cut:
5645 if (fe->is_primitive())
5646 return this->finite_element_output.shape_gradients[i][j];
5647 else
5648 {
5649 // otherwise, use the mapping
5650 // between shape function
5651 // numbers and rows. note that
5652 // by the assertions above, we
5653 // know that this particular
5654 // shape function is primitive,
5655 // so we can call
5656 // system_to_component_index
5657 const unsigned int row =
5658 this->finite_element_output
5659 .shape_function_to_row_table[i * fe->n_components() +
5660 fe->system_to_component_index(i).first];
5661 return this->finite_element_output.shape_gradients[row][j];
5662 }
5663}
5664
5665
5666
5667template <int dim, int spacedim>
5670 const unsigned int i,
5671 const unsigned int j,
5672 const unsigned int component) const
5673{
5674 AssertIndexRange(i, fe->n_dofs_per_cell());
5675 Assert(this->update_flags & update_gradients,
5676 ExcAccessToUninitializedField("update_gradients"));
5677 AssertIndexRange(component, fe->n_components());
5678 Assert(present_cell.is_initialized(), ExcNotReinited());
5679 // check whether the shape function
5680 // is non-zero at all within
5681 // this component:
5682 if (fe->get_nonzero_components(i)[component] == false)
5683 return Tensor<1, spacedim>();
5684
5685 // look up the right row in the
5686 // table and take the data from
5687 // there
5688 const unsigned int row =
5689 this->finite_element_output
5690 .shape_function_to_row_table[i * fe->n_components() + component];
5691 return this->finite_element_output.shape_gradients[row][j];
5692}
5693
5694
5695
5696template <int dim, int spacedim>
5697inline const Tensor<2, spacedim> &
5699 const unsigned int j) const
5700{
5701 AssertIndexRange(i, fe->n_dofs_per_cell());
5702 Assert(this->update_flags & update_hessians,
5703 ExcAccessToUninitializedField("update_hessians"));
5704 Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
5705 Assert(present_cell.is_initialized(), ExcNotReinited());
5706 // if the entire FE is primitive,
5707 // then we can take a short-cut:
5708 if (fe->is_primitive())
5709 return this->finite_element_output.shape_hessians[i][j];
5710 else
5711 {
5712 // otherwise, use the mapping
5713 // between shape function
5714 // numbers and rows. note that
5715 // by the assertions above, we
5716 // know that this particular
5717 // shape function is primitive,
5718 // so we can call
5719 // system_to_component_index
5720 const unsigned int row =
5721 this->finite_element_output
5722 .shape_function_to_row_table[i * fe->n_components() +
5723 fe->system_to_component_index(i).first];
5724 return this->finite_element_output.shape_hessians[row][j];
5725 }
5726}
5727
5728
5729
5730template <int dim, int spacedim>
5733 const unsigned int i,
5734 const unsigned int j,
5735 const unsigned int component) const
5736{
5737 AssertIndexRange(i, fe->n_dofs_per_cell());
5738 Assert(this->update_flags & update_hessians,
5739 ExcAccessToUninitializedField("update_hessians"));
5740 AssertIndexRange(component, fe->n_components());
5741 Assert(present_cell.is_initialized(), ExcNotReinited());
5742 // check whether the shape function
5743 // is non-zero at all within
5744 // this component:
5745 if (fe->get_nonzero_components(i)[component] == false)
5746 return Tensor<2, spacedim>();
5747
5748 // look up the right row in the
5749 // table and take the data from
5750 // there
5751 const unsigned int row =
5752 this->finite_element_output
5753 .shape_function_to_row_table[i * fe->n_components() + component];
5754 return this->finite_element_output.shape_hessians[row][j];
5755}
5756
5757
5758
5759template <int dim, int spacedim>
5760inline const Tensor<3, spacedim> &
5762 const unsigned int j) const
5763{
5764 AssertIndexRange(i, fe->n_dofs_per_cell());
5765 Assert(this->update_flags & update_3rd_derivatives,
5766 ExcAccessToUninitializedField("update_3rd_derivatives"));
5767 Assert(fe->is_primitive(i), ExcShapeFunctionNotPrimitive(i));
5768 Assert(present_cell.is_initialized(), ExcNotReinited());
5769 // if the entire FE is primitive,
5770 // then we can take a short-cut:
5771 if (fe->is_primitive())
5772 return this->finite_element_output.shape_3rd_derivatives[i][j];
5773 else
5774 {
5775 // otherwise, use the mapping
5776 // between shape function
5777 // numbers and rows. note that
5778 // by the assertions above, we
5779 // know that this particular
5780 // shape function is primitive,
5781 // so we can call
5782 // system_to_component_index
5783 const unsigned int row =
5784 this->finite_element_output
5785 .shape_function_to_row_table[i * fe->n_components() +
5786 fe->system_to_component_index(i).first];
5787 return this->finite_element_output.shape_3rd_derivatives[row][j];
5788 }
5789}
5790
5791
5792
5793template <int dim, int spacedim>
5796 const unsigned int i,
5797 const unsigned int j,
5798 const unsigned int component) const
5799{
5800 AssertIndexRange(i, fe->n_dofs_per_cell());
5801 Assert(this->update_flags & update_3rd_derivatives,
5802 ExcAccessToUninitializedField("update_3rd_derivatives"));
5803 AssertIndexRange(component, fe->n_components());
5804 Assert(present_cell.is_initialized(), ExcNotReinited());
5805 // check whether the shape function
5806 // is non-zero at all within
5807 // this component:
5808 if (fe->get_nonzero_components(i)[component] == false)
5809 return Tensor<3, spacedim>();
5810
5811 // look up the right row in the
5812 // table and take the data from
5813 // there
5814 const unsigned int row =
5815 this->finite_element_output
5816 .shape_function_to_row_table[i * fe->n_components() + component];
5817 return this->finite_element_output.shape_3rd_derivatives[row][j];
5818}
5819
5820
5821
5822template <int dim, int spacedim>
5823inline const FiniteElement<dim, spacedim> &
5825{
5826 return *fe;
5827}
5828
5829
5830
5831template <int dim, int spacedim>
5832inline const Mapping<dim, spacedim> &
5834{
5835 return *mapping;
5836}
5837
5838
5839
5840template <int dim, int spacedim>
5841inline UpdateFlags
5843{
5844 return this->update_flags;
5845}
5846
5847
5848
5849template <int dim, int spacedim>
5850inline const std::vector<Point<spacedim>> &
5852{
5853 Assert(this->update_flags & update_quadrature_points,
5854 ExcAccessToUninitializedField("update_quadrature_points"));
5855 Assert(present_cell.is_initialized(), ExcNotReinited());
5856 return this->mapping_output.quadrature_points;
5857}
5858
5859
5860
5861template <int dim, int spacedim>
5862inline const std::vector<double> &
5864{
5865 Assert(this->update_flags & update_JxW_values,
5866 ExcAccessToUninitializedField("update_JxW_values"));
5867 Assert(present_cell.is_initialized(), ExcNotReinited());
5868 return this->mapping_output.JxW_values;
5869}
5870
5871
5872
5873template <int dim, int spacedim>
5874inline const std::vector<DerivativeForm<1, dim, spacedim>> &
5876{
5877 Assert(this->update_flags & update_jacobians,
5878 ExcAccessToUninitializedField("update_jacobians"));
5879 Assert(present_cell.is_initialized(), ExcNotReinited());
5880 return this->mapping_output.jacobians;
5881}
5882
5883
5884
5885template <int dim, int spacedim>
5886inline const std::vector<DerivativeForm<2, dim, spacedim>> &
5888{
5889 Assert(this->update_flags & update_jacobian_grads,
5890 ExcAccessToUninitializedField("update_jacobians_grads"));
5891 Assert(present_cell.is_initialized(), ExcNotReinited());
5892 return this->mapping_output.jacobian_grads;
5893}
5894
5895
5896
5897template <int dim, int spacedim>
5898inline const Tensor<3, spacedim> &
5900 const unsigned int i) const
5901{
5902 Assert(this->update_flags & update_jacobian_pushed_forward_grads,
5903 ExcAccessToUninitializedField("update_jacobian_pushed_forward_grads"));
5904 Assert(present_cell.is_initialized(), ExcNotReinited());
5905 return this->mapping_output.jacobian_pushed_forward_grads[i];
5906}
5907
5908
5909
5910template <int dim, int spacedim>
5911inline const std::vector<Tensor<3, spacedim>> &
5913{
5914 Assert(this->update_flags & update_jacobian_pushed_forward_grads,
5915 ExcAccessToUninitializedField("update_jacobian_pushed_forward_grads"));
5916 Assert(present_cell.is_initialized(), ExcNotReinited());
5917 return this->mapping_output.jacobian_pushed_forward_grads;
5918}
5919
5920
5921
5922template <int dim, int spacedim>
5925{
5926 Assert(this->update_flags & update_jacobian_2nd_derivatives,
5927 ExcAccessToUninitializedField("update_jacobian_2nd_derivatives"));
5928 Assert(present_cell.is_initialized(), ExcNotReinited());
5929 return this->mapping_output.jacobian_2nd_derivatives[i];
5930}
5931
5932
5933
5934template <int dim, int spacedim>
5935inline const std::vector<DerivativeForm<3, dim, spacedim>> &
5937{
5938 Assert(this->update_flags & update_jacobian_2nd_derivatives,
5939 ExcAccessToUninitializedField("update_jacobian_2nd_derivatives"));
5940 Assert(present_cell.is_initialized(), ExcNotReinited());
5941 return this->mapping_output.jacobian_2nd_derivatives;
5942}
5943
5944
5945
5946template <int dim, int spacedim>
5947inline const Tensor<4, spacedim> &
5949 const unsigned int i) const
5950{
5952 ExcAccessToUninitializedField(
5953 "update_jacobian_pushed_forward_2nd_derivatives"));
5954 Assert(present_cell.is_initialized(), ExcNotReinited());
5955 return this->mapping_output.jacobian_pushed_forward_2nd_derivatives[i];
5956}
5957
5958
5959
5960template <int dim, int spacedim>
5961inline const std::vector<Tensor<4, spacedim>> &
5963{
5965 ExcAccessToUninitializedField(
5966 "update_jacobian_pushed_forward_2nd_derivatives"));
5967 Assert(present_cell.is_initialized(), ExcNotReinited());
5968 return this->mapping_output.jacobian_pushed_forward_2nd_derivatives;
5969}
5970
5971
5972
5973template <int dim, int spacedim>
5976{
5977 Assert(this->update_flags & update_jacobian_3rd_derivatives,
5978 ExcAccessToUninitializedField("update_jacobian_3rd_derivatives"));
5979 Assert(present_cell.is_initialized(), ExcNotReinited());
5980 return this->mapping_output.jacobian_3rd_derivatives[i];
5981}
5982
5983
5984
5985template <int dim, int spacedim>
5986inline const std::vector<DerivativeForm<4, dim, spacedim>> &
5988{
5989 Assert(this->update_flags & update_jacobian_3rd_derivatives,
5990 ExcAccessToUninitializedField("update_jacobian_3rd_derivatives"));
5991 Assert(present_cell.is_initialized(), ExcNotReinited());
5992 return this->mapping_output.jacobian_3rd_derivatives;
5993}
5994
5995
5996
5997template <int dim, int spacedim>
5998inline const Tensor<5, spacedim> &
6000 const unsigned int i) const
6001{
6003 ExcAccessToUninitializedField(
6004 "update_jacobian_pushed_forward_3rd_derivatives"));
6005 Assert(present_cell.is_initialized(), ExcNotReinited());
6006 return this->mapping_output.jacobian_pushed_forward_3rd_derivatives[i];
6007}
6008
6009
6010
6011template <int dim, int spacedim>
6012inline const std::vector<Tensor<5, spacedim>> &
6014{
6016 ExcAccessToUninitializedField(
6017 "update_jacobian_pushed_forward_3rd_derivatives"));
6018 Assert(present_cell.is_initialized(), ExcNotReinited());
6019 return this->mapping_output.jacobian_pushed_forward_3rd_derivatives;
6020}
6021
6022
6023
6024template <int dim, int spacedim>
6025inline const std::vector<DerivativeForm<1, spacedim, dim>> &
6027{
6028 Assert(this->update_flags & update_inverse_jacobians,
6029 ExcAccessToUninitializedField("update_inverse_jacobians"));
6030 Assert(present_cell.is_initialized(), ExcNotReinited());
6031 return this->mapping_output.inverse_jacobians;
6032}
6033
6034
6035
6036template <int dim, int spacedim>
6039{
6040 return {0U, dofs_per_cell};
6041}
6042
6043
6044
6045template <int dim, int spacedim>
6048 const unsigned int start_dof_index) const
6049{
6050 Assert(start_dof_index <= dofs_per_cell,
6051 ExcIndexRange(start_dof_index, 0, dofs_per_cell + 1));
6052 return {start_dof_index, dofs_per_cell};
6053}
6054
6055
6056
6057template <int dim, int spacedim>
6060 const unsigned int end_dof_index) const
6061{
6062 Assert(end_dof_index < dofs_per_cell,
6063 ExcIndexRange(end_dof_index, 0, dofs_per_cell));
6064 return {0U, end_dof_index + 1};
6065}
6066
6067
6068
6069template <int dim, int spacedim>
6072{
6073 return {0U, n_quadrature_points};
6074}
6075
6076
6077
6078template <int dim, int spacedim>
6079inline const Point<spacedim> &
6080FEValuesBase<dim, spacedim>::quadrature_point(const unsigned int i) const
6081{
6082 Assert(this->update_flags & update_quadrature_points,
6083 ExcAccessToUninitializedField("update_quadrature_points"));
6084 AssertIndexRange(i, this->mapping_output.quadrature_points.size());
6085 Assert(present_cell.is_initialized(), ExcNotReinited());
6086
6087 return this->mapping_output.quadrature_points[i];
6088}
6089
6090
6091
6092template <int dim, int spacedim>
6093inline double
6094FEValuesBase<dim, spacedim>::JxW(const unsigned int i) const
6095{
6096 Assert(this->update_flags & update_JxW_values,
6097 ExcAccessToUninitializedField("update_JxW_values"));
6098 AssertIndexRange(i, this->mapping_output.JxW_values.size());
6099 Assert(present_cell.is_initialized(), ExcNotReinited());
6100
6101 return this->mapping_output.JxW_values[i];
6102}
6103
6104
6105
6106template <int dim, int spacedim>
6108FEValuesBase<dim, spacedim>::jacobian(const unsigned int i) const
6109{
6110 Assert(this->update_flags & update_jacobians,
6111 ExcAccessToUninitializedField("update_jacobians"));
6112 AssertIndexRange(i, this->mapping_output.jacobians.size());
6113 Assert(present_cell.is_initialized(), ExcNotReinited());
6114
6115 return this->mapping_output.jacobians[i];
6116}
6117
6118
6119
6120template <int dim, int spacedim>
6122FEValuesBase<dim, spacedim>::jacobian_grad(const unsigned int i) const
6123{
6124 Assert(this->update_flags & update_jacobian_grads,
6125 ExcAccessToUninitializedField("update_jacobians_grads"));
6126 AssertIndexRange(i, this->mapping_output.jacobian_grads.size());
6127 Assert(present_cell.is_initialized(), ExcNotReinited());
6128
6129 return this->mapping_output.jacobian_grads[i];
6130}
6131
6132
6133
6134template <int dim, int spacedim>
6136FEValuesBase<dim, spacedim>::inverse_jacobian(const unsigned int i) const
6137{
6138 Assert(this->update_flags & update_inverse_jacobians,
6139 ExcAccessToUninitializedField("update_inverse_jacobians"));
6140 AssertIndexRange(i, this->mapping_output.inverse_jacobians.size());
6141 Assert(present_cell.is_initialized(), ExcNotReinited());
6142
6143 return this->mapping_output.inverse_jacobians[i];
6144}
6145
6146
6147
6148template <int dim, int spacedim>
6149inline const Tensor<1, spacedim> &
6150FEValuesBase<dim, spacedim>::normal_vector(const unsigned int i) const
6151{
6152 Assert(this->update_flags & update_normal_vectors,
6154 "update_normal_vectors")));
6155 AssertIndexRange(i, this->mapping_output.normal_vectors.size());
6156 Assert(present_cell.is_initialized(), ExcNotReinited());
6157
6158 return this->mapping_output.normal_vectors[i];
6159}
6160
6161
6162
6163/*--------------------- Inline functions: FEValues --------------------------*/
6164
6165
6166template <int dim, int spacedim>
6167inline const Quadrature<dim> &
6169{
6170 return quadrature;
6171}
6172
6173
6174
6175template <int dim, int spacedim>
6176inline const FEValues<dim, spacedim> &
6178{
6179 return *this;
6180}
6181
6182
6183/*---------------------- Inline functions: FEFaceValuesBase -----------------*/
6184
6185
6186template <int dim, int spacedim>
6187inline unsigned int
6189{
6190 return present_face_no;
6191}
6192
6193
6194template <int dim, int spacedim>
6195inline unsigned int
6197{
6198 return present_face_index;
6199}
6200
6201
6202/*----------------------- Inline functions: FE*FaceValues -------------------*/
6203
6204template <int dim, int spacedim>
6205inline const Quadrature<dim - 1> &
6207{
6208 return quadrature[quadrature.size() == 1 ? 0 : present_face_no];
6209}
6210
6211
6212
6213template <int dim, int spacedim>
6214inline const FEFaceValues<dim, spacedim> &
6216{
6217 return *this;
6218}
6219
6220
6221
6222template <int dim, int spacedim>
6223inline const FESubfaceValues<dim, spacedim> &
6225{
6226 return *this;
6227}
6228
6229
6230
6231template <int dim, int spacedim>
6232inline const Tensor<1, spacedim> &
6233FEFaceValuesBase<dim, spacedim>::boundary_form(const unsigned int i) const
6234{
6235 AssertIndexRange(i, this->mapping_output.boundary_forms.size());
6236 Assert(this->update_flags & update_boundary_forms,
6238 "update_boundary_forms")));
6239
6240 return this->mapping_output.boundary_forms[i];
6241}
6242
6243#endif // DOXYGEN
6244
6246
6247#endif
const Tensor< 1, spacedim > & boundary_form(const unsigned int i) const
const Quadrature< dim - 1 > & get_quadrature() const
unsigned int get_face_index() const
unsigned int present_face_no
Definition: fe_values.h:4273
FEFaceValuesBase(const unsigned int dofs_per_cell, const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const hp::QCollection< dim - 1 > &quadrature)
unsigned int present_face_index
Definition: fe_values.h:4279
FEFaceValuesBase(const unsigned int dofs_per_cell, const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature)
std::size_t memory_consumption() const
const std::vector< Tensor< 1, spacedim > > & get_boundary_forms() const
unsigned int get_face_number() const
const hp::QCollection< dim - 1 > quadrature
Definition: fe_values.h:4284
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell, const typename Triangulation< dim, spacedim >::face_iterator &face)
const FEFaceValues< dim, spacedim > & get_present_fe_values() const
FEFaceValues(const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature, const UpdateFlags update_flags)
void initialize(const UpdateFlags update_flags)
void reinit(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no)
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell, const unsigned int face_no)
void do_reinit(const unsigned int face_no)
FEFaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature, const UpdateFlags update_flags)
FEFaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const hp::QCollection< dim - 1 > &quadrature, const UpdateFlags update_flags)
FEFaceValues(const FiniteElement< dim, spacedim > &fe, const hp::QCollection< dim - 1 > &quadrature, const UpdateFlags update_flags)
void reinit(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::face_iterator &face)
FESubfaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &face_quadrature, const UpdateFlags update_flags)
void reinit(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no)
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell, const unsigned int face_no, const unsigned int subface_no)
void initialize(const UpdateFlags update_flags)
const FESubfaceValues< dim, spacedim > & get_present_fe_values() const
FESubfaceValues(const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &face_quadrature, const UpdateFlags update_flags)
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell, const typename Triangulation< dim, spacedim >::face_iterator &face, const typename Triangulation< dim, spacedim >::face_iterator &subface)
void do_reinit(const unsigned int face_no, const unsigned int subface_no)
FESubfaceValues(const FiniteElement< dim, spacedim > &fe, const hp::QCollection< dim - 1 > &face_quadrature, const UpdateFlags update_flags)
FESubfaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const hp::QCollection< dim - 1 > &face_quadrature, const UpdateFlags update_flags)
void reinit(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::face_iterator &face, const typename Triangulation< dim, spacedim >::face_iterator &subface)
const DoFHandler< dim, spacedim > * dof_handler
Definition: fe_values.h:3885
CellIteratorContainer(const TriaIterator< DoFCellAccessor< dim, spacedim, lda > > &cell)
Triangulation< dim, spacedim >::cell_iterator cell
Definition: fe_values.h:3884
CellSimilarity::Similarity cell_similarity
Definition: fe_values.h:4003
std_cxx20::ranges::iota_view< unsigned int, unsigned int > dof_indices_ending_at(const unsigned int end_dof_index) const
const DerivativeForm< 2, dim, spacedim > & jacobian_grad(const unsigned int quadrature_point) const
const FEValuesViews::Vector< dim, spacedim > & operator[](const FEValuesExtractors::Vector &vector) const
const std::vector< double > & get_JxW_values() const
CellIteratorContainer present_cell
Definition: fe_values.h:3894
::internal::FEValuesViews::Cache< dim, spacedim > fe_values_views_cache
Definition: fe_values.h:4018
FEValuesBase(const FEValuesBase &)=delete
const FEValuesViews::Scalar< dim, spacedim > & operator[](const FEValuesExtractors::Scalar &scalar) const
const Tensor< 4, spacedim > & jacobian_pushed_forward_2nd_derivative(const unsigned int quadrature_point) const
boost::signals2::connection tria_listener_mesh_transform
Definition: fe_values.h:3912
const std::vector< Point< spacedim > > & get_quadrature_points() const
const std::vector< DerivativeForm< 1, dim, spacedim > > & get_jacobians() const
const SmartPointer< const Mapping< dim, spacedim >, FEValuesBase< dim, spacedim > > mapping
Definition: fe_values.h:3939
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
Definition: fe_values.h:2450
UpdateFlags update_flags
Definition: fe_values.h:3985
const SmartPointer< const FiniteElement< dim, spacedim >, FEValuesBase< dim, spacedim > > fe
Definition: fe_values.h:3963
std_cxx20::ranges::iota_view< unsigned int, unsigned int > dof_indices_starting_at(const unsigned int start_dof_index) const
const Mapping< dim, spacedim > & get_mapping() const
const unsigned int n_quadrature_points
Definition: fe_values.h:2432
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
Tensor< 2, spacedim > shape_hessian_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
const std::vector< Tensor< 5, spacedim > > & get_jacobian_pushed_forward_3rd_derivatives() const
const DerivativeForm< 1, dim, spacedim > & jacobian(const unsigned int quadrature_point) const
::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > mapping_output
Definition: fe_values.h:3954
std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > mapping_data
Definition: fe_values.h:3947
const std::vector< DerivativeForm< 2, dim, spacedim > > & get_jacobian_grads() const
Tensor< 1, spacedim > shape_grad_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
const Tensor< 2, spacedim > & shape_hessian(const unsigned int function_no, const unsigned int point_no) const
boost::signals2::connection tria_listener_refinement
Definition: fe_values.h:3903
const std::vector< Tensor< 3, spacedim > > & get_jacobian_pushed_forward_grads() const
const Tensor< 5, spacedim > & jacobian_pushed_forward_3rd_derivative(const unsigned int quadrature_point) const
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const
const std::vector< Tensor< 4, spacedim > > & get_jacobian_pushed_forward_2nd_derivatives() const
double JxW(const unsigned int quadrature_point) const
const Point< spacedim > & quadrature_point(const unsigned int q) const
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > finite_element_output
Definition: fe_values.h:3979
const Tensor< 3, spacedim > & jacobian_pushed_forward_grad(const unsigned int quadrature_point) const
const FEValuesViews::SymmetricTensor< 2, dim, spacedim > & operator[](const FEValuesExtractors::SymmetricTensor< 2 > &tensor) const
const std::vector< DerivativeForm< 4, dim, spacedim > > & get_jacobian_3rd_derivatives() const
const DerivativeForm< 4, dim, spacedim > & jacobian_3rd_derivative(const unsigned int quadrature_point) const
const Tensor< 3, spacedim > & shape_3rd_derivative(const unsigned int function_no, const unsigned int point_no) const
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > fe_data
Definition: fe_values.h:3971
const FiniteElement< dim, spacedim > & get_fe() const
const Tensor< 1, spacedim > & shape_grad(const unsigned int function_no, const unsigned int quadrature_point) const
FEValuesBase & operator=(const FEValuesBase &)=delete
const DerivativeForm< 3, dim, spacedim > & jacobian_2nd_derivative(const unsigned int quadrature_point) const
const DerivativeForm< 1, spacedim, dim > & inverse_jacobian(const unsigned int quadrature_point) const
const std::vector< DerivativeForm< 1, spacedim, dim > > & get_inverse_jacobians() const
const unsigned int max_n_quadrature_points
Definition: fe_values.h:2443
const double & shape_value(const unsigned int function_no, const unsigned int point_no) const
Tensor< 3, spacedim > shape_3rd_derivative_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
void get_function_third_derivatives_from_local_dof_values(const InputVector &dof_values, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &third_derivatives) const
Definition: fe_values.cc:1795
Scalar & operator=(const Scalar< dim, spacedim > &)=delete
typename ProductType< Number, hessian_type >::type solution_hessian_type
Definition: fe_values.h:214
value_type value(const unsigned int shape_function, const unsigned int q_point) const
const unsigned int component
Definition: fe_values.h:634
void get_function_hessians_from_local_dof_values(const InputVector &dof_values, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
Definition: fe_values.cc:1686
void get_function_values_from_local_dof_values(const InputVector &dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
Definition: fe_values.cc:1578
gradient_type gradient(const unsigned int shape_function, const unsigned int q_point) const
std::vector< ShapeFunctionData > shape_function_data
Definition: fe_values.h:639
void get_function_laplacians(const InputVector &fe_function, std::vector< solution_laplacian_type< typename InputVector::value_type > > &laplacians) const
Definition: fe_values.cc:1710
typename ProductType< Number, value_type >::type solution_laplacian_type
Definition: fe_values.h:204
void get_function_third_derivatives(const InputVector &fe_function, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &third_derivatives) const
Definition: fe_values.cc:1764
Scalar & operator=(Scalar< dim, spacedim > &&) noexcept=default
third_derivative_type third_derivative(const unsigned int shape_function, const unsigned int q_point) const
void get_function_hessians(const InputVector &fe_function, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
Definition: fe_values.cc:1656
typename ProductType< Number, third_derivative_type >::type solution_third_derivative_type
Definition: fe_values.h:224
typename ProductType< Number, value_type >::type solution_value_type
Definition: fe_values.h:184
typename ProductType< Number, gradient_type >::type solution_gradient_type
Definition: fe_values.h:194
void get_function_laplacians_from_local_dof_values(const InputVector &dof_values, std::vector< solution_laplacian_type< typename InputVector::value_type > > &laplacians) const
Definition: fe_values.cc:1740
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
Definition: fe_values.h:628
void get_function_gradients(const InputVector &fe_function, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
Definition: fe_values.cc:1602
void get_function_gradients_from_local_dof_values(const InputVector &dof_values, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
Definition: fe_values.cc:1632
Scalar(Scalar< dim, spacedim > &&)=default
Scalar(const Scalar< dim, spacedim > &)=delete
hessian_type hessian(const unsigned int shape_function, const unsigned int q_point) const
void get_function_values(const InputVector &fe_function, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
Definition: fe_values.cc:1547
SymmetricTensor(const SymmetricTensor< 2, dim, spacedim > &)=delete
SymmetricTensor(SymmetricTensor< 2, dim, spacedim > &&)=default
SymmetricTensor & operator=(SymmetricTensor< 2, dim, spacedim > &&) noexcept=default
typename ProductType< Number, value_type >::type solution_value_type
Definition: fe_values.h:1505
SymmetricTensor & operator=(const SymmetricTensor< 2, dim, spacedim > &)=delete
typename ProductType< Number, divergence_type >::type solution_divergence_type
Definition: fe_values.h:1515
typename ProductType< Number, value_type >::type solution_value_type
Definition: fe_values.h:1841
typename ProductType< Number, divergence_type >::type solution_divergence_type
Definition: fe_values.h:1851
typename ProductType< Number, gradient_type >::type solution_gradient_type
Definition: fe_values.h:1861
divergence_type divergence(const unsigned int shape_function, const unsigned int q_point) const
Tensor(const Tensor< 2, dim, spacedim > &)=delete
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
Definition: fe_values.h:2189
Tensor & operator=(const Tensor< 2, dim, spacedim > &)=delete
Tensor & operator=(Tensor< 2, dim, spacedim > &&)=default
value_type value(const unsigned int shape_function, const unsigned int q_point) const
gradient_type gradient(const unsigned int shape_function, const unsigned int q_point) const
Tensor(Tensor< 2, dim, spacedim > &&)=default
std::vector< ShapeFunctionData > shape_function_data
Definition: fe_values.h:2200
typename ProductType< Number, third_derivative_type >::type solution_third_derivative_type
Definition: fe_values.h:811
typename ProductType< Number, divergence_type >::type solution_divergence_type
Definition: fe_values.h:772
Vector(const Vector< dim, spacedim > &)=delete
hessian_type hessian(const unsigned int shape_function, const unsigned int q_point) const
typename ProductType< Number, hessian_type >::type solution_hessian_type
Definition: fe_values.h:801
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
Definition: fe_values.h:1435
Vector & operator=(Vector< dim, spacedim > &&)=default
typename ProductType< Number, symmetric_gradient_type >::type solution_symmetric_gradient_type
Definition: fe_values.h:762
gradient_type gradient(const unsigned int shape_function, const unsigned int q_point) const
divergence_type divergence(const unsigned int shape_function, const unsigned int q_point) const
Vector & operator=(const Vector< dim, spacedim > &)=delete
third_derivative_type third_derivative(const unsigned int shape_function, const unsigned int q_point) const
typename ProductType< Number, gradient_type >::type solution_gradient_type
Definition: fe_values.h:752
typename ProductType< Number, value_type >::type solution_value_type
Definition: fe_values.h:742
symmetric_gradient_type symmetric_gradient(const unsigned int shape_function, const unsigned int q_point) const
Vector(Vector< dim, spacedim > &&)=default
typename ::internal::CurlType< spacedim >::type curl_type
Definition: fe_values.h:719
const unsigned int first_vector_component
Definition: fe_values.h:1441
typename ProductType< Number, curl_type >::type solution_curl_type
Definition: fe_values.h:791
std::vector< ShapeFunctionData > shape_function_data
Definition: fe_values.h:1446
value_type value(const unsigned int shape_function, const unsigned int q_point) const
curl_type curl(const unsigned int shape_function, const unsigned int q_point) const
typename ProductType< Number, value_type >::type solution_laplacian_type
Definition: fe_values.h:782
FEValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags update_flags)
void do_reinit()
const Quadrature< dim > & get_quadrature() const
const Quadrature< dim > quadrature
Definition: fe_values.h:4154
FEValues(const FiniteElement< dim, spacedim > &fe, const hp::QCollection< dim > &quadrature, const UpdateFlags update_flags)
const FEValues< dim, spacedim > & get_present_fe_values() const
FEValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const hp::QCollection< dim > &quadrature, const UpdateFlags update_flags)
FEValues(const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags update_flags)
void initialize(const UpdateFlags update_flags)
void reinit(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
std::size_t memory_consumption() const
Abstract base class for mapping classes.
Definition: mapping.h:311
Definition: point.h:111
Definition: tensor.h:503
typename Tensor< rank_ - 1, dim, Number >::tensor_type value_type
Definition: tensor.h:536
static constexpr TableIndices< rank_ > unrolled_to_component_indices(const unsigned int i)
Definition: vector.h:109
#define DEAL_II_DEPRECATED
Definition: config.h:164
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
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_boundary_forms
Outer normal vector, not normalized.
@ update_jacobian_2nd_derivatives
#define DeclException0(Exception0)
Definition: exceptions.h:464
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
Definition: exceptions.h:1473
#define AssertIndexRange(index, range)
Definition: exceptions.h:1732
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:487
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcIndexRange(std::size_t arg1, std::size_t arg2, std::size_t arg3)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:509
static ::ExceptionBase & ExcMessage(std::string arg1)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition: tria.h:1355
Number value_type
Definition: vector.h:121
typename ::internal::FEValuesViews::ViewType< dim, spacedim, Extractor >::type View
Definition: fe_values.h:2308
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
boost::integer_range< IncrementableType > iota_view
Definition: iota_view.h:46
STL namespace.
typename ProductType< Number, typename Scalar< dim, spacedim >::value_type >::type laplacian_type
Definition: fe_values.h:257
typename ProductType< Number, typename Scalar< dim, spacedim >::value_type >::type value_type
Definition: fe_values.h:241
typename ProductType< Number, typename Scalar< dim, spacedim >::gradient_type >::type gradient_type
Definition: fe_values.h:249
typename ProductType< Number, typename Scalar< dim, spacedim >::hessian_type >::type hessian_type
Definition: fe_values.h:265
typename ProductType< Number, typename Scalar< dim, spacedim >::third_derivative_type >::type third_derivative_type
Definition: fe_values.h:273
typename ProductType< Number, typename SymmetricTensor< 2, dim, spacedim >::divergence_type >::type divergence_type
Definition: fe_values.h:1541
typename ProductType< Number, typename SymmetricTensor< 2, dim, spacedim >::value_type >::type value_type
Definition: fe_values.h:1533
typename ProductType< Number, typename Tensor< 2, dim, spacedim >::gradient_type >::type gradient_type
Definition: fe_values.h:1895
typename ProductType< Number, typename Tensor< 2, dim, spacedim >::value_type >::type value_type
Definition: fe_values.h:1879
typename ProductType< Number, typename Tensor< 2, dim, spacedim >::divergence_type >::type divergence_type
Definition: fe_values.h:1887
typename ProductType< Number, typename Vector< dim, spacedim >::third_derivative_type >::type third_derivative_type
Definition: fe_values.h:884
typename ProductType< Number, typename Vector< dim, spacedim >::value_type >::type laplacian_type
Definition: fe_values.h:860
typename ProductType< Number, typename Vector< dim, spacedim >::hessian_type >::type hessian_type
Definition: fe_values.h:876
typename ProductType< Number, typename Vector< dim, spacedim >::symmetric_gradient_type >::type symmetric_gradient_type
Definition: fe_values.h:844
typename ProductType< Number, typename Vector< dim, spacedim >::gradient_type >::type gradient_type
Definition: fe_values.h:836
typename ProductType< Number, typename Vector< dim, spacedim >::curl_type >::type curl_type
Definition: fe_values.h:868
typename ProductType< Number, typename Vector< dim, spacedim >::divergence_type >::type divergence_type
Definition: fe_values.h:852
typename ProductType< Number, typename Vector< dim, spacedim >::value_type >::type value_type
Definition: fe_values.h:828
typename internal::ProductTypeImpl< typename std::decay< T >::type, typename std::decay< U >::type >::type type
std::vector<::FEValuesViews::Scalar< dim, spacedim > > scalars
Definition: fe_values.h:2285
std::vector<::FEValuesViews::Vector< dim, spacedim > > vectors
Definition: fe_values.h:2286
std::vector<::FEValuesViews::SymmetricTensor< 2, dim, spacedim > > symmetric_second_order_tensors
Definition: fe_values.h:2288
std::vector<::FEValuesViews::Tensor< 2, dim, spacedim > > second_order_tensors
Definition: fe_values.h:2290
typename ::FEValuesViews::SymmetricTensor< rank, dim, spacedim > type
Definition: fe_values.h:2268
typename ::FEValuesViews::Tensor< rank, dim, spacedim > type
Definition: fe_values.h:2254
constexpr SymmetricTensor< 2, dim, Number > symmetrize(const Tensor< 2, dim, Number > &t)