deal.II version GIT relicensing-2206-gaa53ff9447 2024-12-02 09:10:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
fe_evaluation_data.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2021 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
16#ifndef dealii_matrix_free_fe_evaluation_data_h
17#define dealii_matrix_free_fe_evaluation_data_h
18
19
20#include <deal.II/base/config.h>
21
30#include <deal.II/base/tensor.h>
32
36
37
39
40
41
42namespace internal
43{
45
48 std::string,
49 << "You are requesting information from an FEEvaluation/FEFaceEvaluation "
50 << "object for which this kind of information has not been computed. What "
51 << "information these objects compute is determined by the update_* flags "
52 << "you pass to MatrixFree::reinit() via MatrixFree::AdditionalData. "
53 << "Here, the operation you are attempting requires the <" << arg1
54 << "> flag to be set, but it was apparently not specified "
55 << "upon initialization.");
56} // namespace internal
57
58// forward declarations
59template <int dim,
60 int fe_degree,
61 int n_q_points_1d = fe_degree + 1,
62 int n_components_ = 1,
63 typename Number = double,
64 typename VectorizedArrayType = VectorizedArray<Number>>
65class FEEvaluation;
66
67template <int dim,
68 int n_components_,
69 typename Number,
70 bool = false,
71 typename = VectorizedArray<Number>>
73
74namespace internal
75{
76 namespace MatrixFreeFunctions
77 {
78 template <int, typename>
79 class MappingDataOnTheFly;
80 }
81} // namespace internal
82
83
112template <int dim, typename Number, bool is_face>
114{
117 using MappingInfoStorageType = internal::MatrixFreeFunctions::
118 MappingInfoStorage<(is_face ? dim - 1 : dim), dim, Number>;
120
121public:
122 static constexpr unsigned int dimension = dim;
123
124 using NumberType = Number;
128
129 static constexpr unsigned int n_lanes =
131
140 const bool is_interior_face = true);
141
145 FEEvaluationData(const FEEvaluationData &other) = default;
146
152
156 virtual ~FEEvaluationData() = default;
157
163 void
165 const unsigned int n_components);
166
171 void
174
188
193 Number
194 JxW(const unsigned int q_point) const;
195
201 quadrature_point(const unsigned int q) const;
202
215 inverse_jacobian(const unsigned int q_point) const;
216
230 normal_vector(const unsigned int q_point) const;
231
237 DEAL_II_DEPRECATED_WITH_COMMENT("Use normal_vector() instead.")
238 Tensor<1, dim, Number>
239 get_normal_vector(const unsigned int q_point) const;
240
255 const Number *
257
266 Number *
268
279 const Number *
281
292 Number *
294
306 const Number *
308
320 Number *
322
335 const Number *
337
350 Number *
352
368 unsigned int
370
377 internal::MatrixFreeFunctions::GeometryType
379
383 const ShapeInfoType &
385
389 const internal::MatrixFreeFunctions::DoFInfo &
391
397 const std::vector<unsigned int> &
399
403 unsigned int
405
409 unsigned int
411
416 unsigned int
418
423 unsigned int
425
429 unsigned int
431
444 std::uint8_t
445 get_face_no(const unsigned int v = 0) const;
446
456 unsigned int
458
467 std::uint8_t
468 get_face_orientation(const unsigned int v = 0) const;
469
477 internal::MatrixFreeFunctions::DoFInfo::DoFAccessIndex
479
488 bool
490
495 const std::array<unsigned int, n_lanes> &
497 {
498// implemented inline to avoid compilation problems on Windows
499#ifdef DEBUG
501#endif
502 return cell_ids;
503 }
504
509 const std::array<unsigned int, n_lanes> &
511 {
512// implemented inline to avoid compilation problems on Windows
513#ifdef DEBUG
515#endif
516 return face_ids;
517 }
518
523 unsigned int
525 {
526// implemented inline to avoid compilation problems on Windows
527#ifdef DEBUG
529#endif
530
531 return cell;
532 }
533
538 const std::array<unsigned int, n_lanes> &
540 {
541// implemented inline to avoid compilation problems on Windows
542#ifdef DEBUG
544#endif
545
546 if (!is_face || dof_access_index ==
548 return cell_ids;
549 else
550 return face_ids;
551 }
552
560
578 template <typename T>
579 T
580 read_cell_data(const AlignedVector<T> &array) const;
581
588 template <typename T>
589 void
590 set_cell_data(AlignedVector<T> &array, const T &value) const;
591
600 template <typename T>
601 T
602 read_face_data(const AlignedVector<T> &array) const;
603
612 template <typename T>
613 void
614 set_face_data(AlignedVector<T> &array, const T &value) const;
615
631
632protected:
639 FEEvaluationData(const InitializationData &initialization_data,
640 const bool is_interior_face,
641 const unsigned int quad_no,
642 const unsigned int first_selected_component);
643
649 const std::shared_ptr<
652 const unsigned int n_fe_components,
653 const unsigned int first_selected_component);
654
662
669
677
683 const unsigned int quad_no;
684
689 const unsigned int n_fe_components;
690
695 const unsigned int first_selected_component;
696
700 const unsigned int active_fe_index;
701
706 const unsigned int active_quad_index;
707
714
718 const unsigned int n_quadrature_points;
719
725
739
744 const Tensor<1, dim *(dim + 1) / 2, Tensor<1, dim, Number>>
746
751 const Tensor<1, dim *(dim + 1) / 2, Tensor<1, dim, Number>>
753
760 const Number *J_value;
761
766
771
776
784
795 Number *values_dofs;
796
806 Number *values_quad;
807
821
834
847
858
864
871
878
885
892
899
906
913
918 unsigned int cell;
919
927
933
942 std::array<std::uint8_t, n_lanes> face_numbers;
943
952 std::array<std::uint8_t, n_lanes> face_orientations;
953
961 unsigned int subface_index;
962
970
975 std::array<unsigned int, n_lanes> cell_ids;
976
981 std::array<unsigned int, n_lanes> face_ids;
982
988 std::shared_ptr<
991
997
998 // Make FEEvaluation and FEEvaluationBase objects friends for access to
999 // protected member mapped_geometry.
1000 template <int, int, typename, bool, typename>
1001 friend class FEEvaluationBase;
1002
1003 template <int, int, int, int, typename, typename>
1004 friend class FEEvaluation;
1005};
1006
1007
1012template <int dim,
1013 typename Number,
1014 bool is_face,
1015 typename VectorizedArrayType = VectorizedArray<Number>>
1018
1019
1020/*----------------------- Inline functions ----------------------------------*/
1021
1022#ifndef DOXYGEN
1023
1024template <int dim, typename Number, bool is_face>
1026 const ShapeInfoType &shape_info,
1027 const bool is_interior_face)
1029 InitializationData{&shape_info, nullptr, nullptr, 0, 0, nullptr},
1031 0,
1032 0)
1033{}
1034
1035
1036template <int dim, typename Number, bool is_face>
1038 const InitializationData &initialization_data,
1039 const bool is_interior_face,
1040 const unsigned int quad_no,
1041 const unsigned int first_selected_component)
1042 : data(initialization_data.shape_info)
1043 , dof_info(initialization_data.dof_info)
1044 , mapping_data(initialization_data.mapping_data)
1045 , quad_no(quad_no)
1046 , n_fe_components(dof_info != nullptr ? dof_info->start_components.back() : 0)
1048 , active_fe_index(initialization_data.active_fe_index)
1049 , active_quad_index(initialization_data.active_quad_index)
1050 , descriptor(initialization_data.descriptor)
1051 , n_quadrature_points(descriptor == nullptr ?
1052 (is_face ? data->n_q_points_face : data->n_q_points) :
1054 , quadrature_points(nullptr)
1055 , jacobian(nullptr)
1056 , jacobian_gradients(nullptr)
1058 , J_value(nullptr)
1059 , normal_vectors(nullptr)
1060 , normal_x_jacobian(nullptr)
1062 descriptor != nullptr ? descriptor->quadrature_weights.begin() : nullptr)
1063# ifdef DEBUG
1064 , is_reinitialized(false)
1065 , dof_values_initialized(false)
1069 , values_quad_submitted(false)
1071# endif
1072 , cell(numbers::invalid_unsigned_int)
1075 is_face ?
1077 internal::MatrixFreeFunctions::DoFInfo::dof_access_face_interior :
1078 internal::MatrixFreeFunctions::DoFInfo::dof_access_face_exterior) :
1079 internal::MatrixFreeFunctions::DoFInfo::dof_access_cell)
1080 , subface_index(0)
1081 , cell_type(internal::MatrixFreeFunctions::general)
1083{}
1084
1085
1086
1087template <int dim, typename Number, bool is_face>
1089 const std::shared_ptr<
1092 const unsigned int n_fe_components,
1093 const unsigned int first_selected_component)
1094 : data(nullptr)
1095 , dof_info(nullptr)
1096 , mapping_data(nullptr)
1097 , quad_no(numbers::invalid_unsigned_int)
1100 , active_fe_index(numbers::invalid_unsigned_int)
1101 , active_quad_index(numbers::invalid_unsigned_int)
1102 , descriptor(nullptr)
1104 mapped_geometry->get_data_storage().descriptor[0].n_q_points)
1105 , quadrature_points(nullptr)
1106 , jacobian(nullptr)
1107 , jacobian_gradients(nullptr)
1109 , J_value(nullptr)
1110 , normal_vectors(nullptr)
1111 , normal_x_jacobian(nullptr)
1112 , quadrature_weights(nullptr)
1113 , cell(0)
1114 , cell_type(internal::MatrixFreeFunctions::general)
1115 , interior_face(true)
1116 , dof_access_index(internal::MatrixFreeFunctions::DoFInfo::dof_access_cell)
1118 , is_reinitialized(false)
1120{
1121 mapping_data = &mapped_geometry->get_data_storage();
1122 jacobian = mapped_geometry->get_data_storage().jacobians[0].begin();
1123 J_value = mapped_geometry->get_data_storage().JxW_values.begin();
1125 mapped_geometry->get_data_storage().jacobian_gradients[0].begin();
1127 .jacobian_gradients_non_inverse[0]
1128 .begin();
1130 mapped_geometry->get_data_storage().quadrature_points.begin();
1131}
1132
1133
1134
1135template <int dim, typename Number, bool is_face>
1138{
1145
1146 data = other.data;
1147 dof_info = other.dof_info;
1148 mapping_data = other.mapping_data;
1149 descriptor = other.descriptor;
1150 jacobian = nullptr;
1151 J_value = nullptr;
1152 normal_vectors = nullptr;
1153 normal_x_jacobian = nullptr;
1154 jacobian_gradients = nullptr;
1156 quadrature_points = nullptr;
1158
1159# ifdef DEBUG
1160 is_reinitialized = false;
1161 dof_values_initialized = false;
1165 values_quad_submitted = false;
1167# endif
1168
1172 is_face ?
1173 (is_interior_face() ?
1176 internal::MatrixFreeFunctions::DoFInfo::dof_access_cell;
1177 face_numbers[0] = 0;
1178 face_orientations[0] = 0;
1179 subface_index = 0;
1182
1183 return *this;
1184}
1185
1186
1187
1188template <int dim, typename Number, bool is_face>
1189inline void
1192 const unsigned int n_components)
1193{
1195
1196 const unsigned int tensor_dofs_per_component =
1197 Utilities::fixed_power<dim>(data->data.front().fe_degree + 1);
1198 const unsigned int dofs_per_component = data->dofs_per_component_on_cell;
1199
1200 const unsigned int size_scratch_data =
1201 std::max(tensor_dofs_per_component + 1, dofs_per_component) * n_components *
1202 3 +
1204 const unsigned int size_data_arrays =
1206 (n_components * ((dim * (dim + 1)) / 2 + 2 * dim + 2) *
1208
1209 // include 12 extra fields to insert some padding between values, gradients
1210 // and hessians, which helps to reduce the probability of cache conflicts
1211 const unsigned int allocated_size = size_scratch_data + size_data_arrays + 12;
1212# ifdef DEBUG
1214 scratch_data_array->resize(allocated_size,
1215 Number(numbers::signaling_nan<ScalarNumber>()));
1216# else
1217 scratch_data_array->resize_fast(allocated_size);
1218# endif
1219 scratch_data.reinit(scratch_data_array->begin() + size_data_arrays + 12,
1220 size_scratch_data);
1221
1222 // set the pointers to the correct position in the data array
1224 values_quad =
1227 scratch_data_array->begin() + 6 +
1230 scratch_data_array->begin() + 8 +
1233 scratch_data_array->begin() + 12 +
1236 scratch_data_array->begin() + 12 +
1238}
1239
1240
1241
1242template <int dim, typename Number, bool is_face>
1243inline void
1246{
1247 Assert(is_face == true,
1248 ExcMessage("Faces can only be set if the is_face template parameter "
1249 "is true"));
1250 face_numbers[0] =
1252 subface_index = is_interior_face() == true ?
1254 face.subface_index;
1255
1256 // First check if interior or exterior cell has non-standard orientation
1257 // (i.e. the third bit is one or not). Then set zero if this cell has
1258 // standard-orientation else copy the first three bits
1259 // (which is equivalent to modulo 8). See also the documentation of
1260 // internal::MatrixFreeFunctions::FaceToCellTopology::face_orientation.
1261 face_orientations[0] = (is_interior_face() == (face.face_orientation >= 8)) ?
1262 (face.face_orientation % 8) :
1263 0;
1264
1265 if (is_interior_face())
1266 cell_ids = face.cells_interior;
1267 else
1268 cell_ids = face.cells_exterior;
1269}
1270
1271
1272
1273template <int dim, typename Number, bool is_face>
1276 const unsigned int q_point) const
1277{
1279 Assert(normal_vectors != nullptr,
1281 "update_normal_vectors"));
1283 return normal_vectors[0];
1284 else
1285 return normal_vectors[q_point];
1286}
1287
1288
1289
1290// This function is deprecated.
1291template <int dim, typename Number, bool is_face>
1294 const unsigned int q_point) const
1295{
1296 return normal_vector(q_point);
1297}
1298
1299
1300
1301template <int dim, typename Number, bool is_face>
1302inline DEAL_II_ALWAYS_INLINE Number
1303FEEvaluationData<dim, Number, is_face>::JxW(const unsigned int q_point) const
1304{
1306 Assert(J_value != nullptr,
1308 "update_values|update_gradients"));
1310 {
1312 return J_value[0] * quadrature_weights[q_point];
1313 }
1314 else
1315 return J_value[q_point];
1316}
1317
1318
1319
1320template <int dim, typename Number, bool is_face>
1323 const unsigned int q) const
1324{
1326 Assert(this->quadrature_points != nullptr,
1328 "update_quadrature_points"));
1329
1330 // Cartesian/affine mesh: only first vertex of cell is stored, we must
1331 // compute it through the Jacobian (which is stored in non-inverted and
1332 // non-transposed form as index '1' in the jacobian field)
1333 if (is_face == false &&
1335 {
1336 Assert(this->jacobian != nullptr, ExcNotInitialized());
1338
1339 const Tensor<2, dim, Number> &jac = this->jacobian[1];
1341 for (unsigned int d = 0; d < dim; ++d)
1342 point[d] += jac[d][d] * static_cast<typename Number::value_type>(
1343 this->descriptor->quadrature.point(q)[d]);
1344 else
1345 for (unsigned int d = 0; d < dim; ++d)
1346 for (unsigned int e = 0; e < dim; ++e)
1347 point[d] += jac[d][e] * static_cast<typename Number::value_type>(
1348 this->descriptor->quadrature.point(q)[e]);
1349 return point;
1350 }
1351 else
1352 return this->quadrature_points[q];
1353}
1354
1355
1356
1357template <int dim, typename Number, bool is_face>
1360 const unsigned int q_point) const
1361{
1363 Assert(jacobian != nullptr,
1365 "update_gradients"));
1367 return jacobian[0];
1368 else
1369 return jacobian[q_point];
1370}
1371
1372
1373
1374template <int dim, typename Number, bool is_face>
1375inline const Number *
1377{
1378 return values_dofs;
1379}
1380
1381
1382
1383template <int dim, typename Number, bool is_face>
1384inline Number *
1386{
1387# ifdef DEBUG
1389# endif
1390 return values_dofs;
1391}
1392
1393
1394
1395template <int dim, typename Number, bool is_face>
1396inline const Number *
1398{
1399# ifdef DEBUG
1401# endif
1402 return values_quad;
1403}
1404
1405
1406
1407template <int dim, typename Number, bool is_face>
1408inline Number *
1410{
1411# ifdef DEBUG
1413 values_quad_submitted = true;
1414# endif
1415 return values_quad;
1416}
1417
1418
1419
1420template <int dim, typename Number, bool is_face>
1421inline const Number *
1423{
1424# ifdef DEBUG
1427# endif
1428 return gradients_quad;
1429}
1430
1431
1432
1433template <int dim, typename Number, bool is_face>
1434inline Number *
1436{
1437# ifdef DEBUG
1440# endif
1441 return gradients_quad;
1442}
1443
1444
1445
1446template <int dim, typename Number, bool is_face>
1447inline const Number *
1449{
1450# ifdef DEBUG
1452# endif
1453 return hessians_quad;
1454}
1455
1456
1457
1458template <int dim, typename Number, bool is_face>
1459inline Number *
1461{
1462# ifdef DEBUG
1464# endif
1465 return hessians_quad;
1466}
1467
1468
1469
1470template <int dim, typename Number, bool is_face>
1471inline unsigned int
1473{
1474 Assert(mapping_data != nullptr, ExcInternalError());
1475
1476 if (dof_info == nullptr)
1477 return 0;
1478 else
1479 {
1482 }
1483}
1484
1485
1486
1487template <int dim, typename Number, bool is_face>
1490{
1491# ifdef DEBUG
1493# endif
1494 return cell_type;
1495}
1496
1497
1498
1499template <int dim, typename Number, bool is_face>
1502{
1503 Assert(data != nullptr, ExcInternalError());
1504 return *data;
1505}
1506
1507
1508
1509template <int dim, typename Number, bool is_face>
1512{
1513 Assert(dof_info != nullptr,
1514 ExcMessage(
1515 "FEEvaluation was not initialized with a MatrixFree object!"));
1516 return *dof_info;
1517}
1518
1519
1520
1521template <int dim, typename Number, bool is_face>
1522inline const std::vector<unsigned int> &
1524{
1525 return data->lexicographic_numbering;
1526}
1527
1528
1529
1530template <int dim, typename Number, bool is_face>
1531inline unsigned int
1533{
1534 return quad_no;
1535}
1536
1537
1538
1539template <int dim, typename Number, bool is_face>
1540inline unsigned int
1542{
1543 if (is_face && dof_access_index ==
1546 else
1547 return cell;
1548}
1549
1550
1551
1552template <int dim, typename Number, bool is_face>
1553inline unsigned int
1555{
1556 return active_fe_index;
1557}
1558
1559
1560
1561template <int dim, typename Number, bool is_face>
1562inline unsigned int
1564{
1565 return active_quad_index;
1566}
1567
1568
1569
1570template <int dim, typename Number, bool is_face>
1571inline unsigned int
1573{
1575}
1576
1577
1578
1579template <int dim, typename Number, bool is_face>
1580inline ArrayView<Number>
1582{
1583 return scratch_data;
1584}
1585
1586
1587
1588template <int dim, typename Number, bool is_face>
1589inline std::uint8_t
1591{
1592 Assert(is_face, ExcNotInitialized());
1596 is_interior_face() == false),
1597 ExcMessage("All face numbers can only be queried for ECL at exterior "
1598 "faces. Use get_face_no() in other cases."));
1599
1600 return face_numbers[v];
1601}
1602
1603
1604
1605template <int dim, typename Number, bool is_face>
1606inline unsigned int
1608{
1609 return subface_index;
1610}
1611
1612
1613
1614template <int dim, typename Number, bool is_face>
1615inline std::uint8_t
1617 const unsigned int v) const
1618{
1619 Assert(is_face, ExcNotInitialized());
1623 is_interior_face() == false),
1624 ExcMessage("All face numbers can only be queried for ECL at exterior "
1625 "faces. Use get_face_no() in other cases."));
1626
1627 return face_orientations[v];
1628}
1629
1630
1631
1632template <int dim, typename Number, bool is_face>
1635{
1636 return dof_access_index;
1637}
1638
1639
1640
1641template <int dim, typename Number, bool is_face>
1642inline bool
1644{
1645 return interior_face;
1646}
1647
1648
1649
1650template <int dim, typename Number, bool is_face>
1653{
1654 return {0U, n_quadrature_points};
1655}
1656
1657
1658
1659namespace internal
1660{
1661 template <std::size_t N,
1662 typename VectorOfArrayType,
1663 typename ArrayType,
1664 typename FU>
1665 inline void
1666 process_cell_or_face_data(const std::array<unsigned int, N> indices,
1667 VectorOfArrayType &array,
1668 ArrayType &out,
1669 const FU &fu)
1670 {
1671 for (unsigned int i = 0; i < N; ++i)
1672 if (indices[i] != numbers::invalid_unsigned_int)
1673 {
1674 AssertIndexRange(indices[i] / N, array.size());
1675 fu(out[i], array[indices[i] / N][indices[i] % N]);
1676 }
1677 }
1678
1679 template <std::size_t N, typename VectorOfArrayType, typename ArrayType>
1680 inline void
1681 set_valid_element_to_array(const std::array<unsigned int, N> indices,
1682 const VectorOfArrayType &array,
1683 ArrayType &out)
1684 {
1685 AssertDimension(indices.size(), out.size());
1686 AssertDimension(indices.size(), array[0].size());
1687 // set all entries in array to a valid element
1688 int index = 0;
1689 for (; index < N; ++index)
1690 if (indices[index] != numbers::invalid_unsigned_int)
1691 break;
1692 for (unsigned int i = 0; i < N; ++i)
1693 out[i] = array[indices[index] / N][indices[index] % N];
1694 }
1695} // namespace internal
1696
1697
1698
1699template <int dim, typename Number, bool is_face>
1700template <typename T>
1701inline T
1703 const AlignedVector<T> &array) const
1704{
1705 T out;
1706 internal::set_valid_element_to_array(this->get_cell_ids(), array, out);
1707 internal::process_cell_or_face_data(this->get_cell_ids(),
1708 array,
1709 out,
1710 [](auto &local, const auto &global) {
1711 local = global;
1712 });
1713 return out;
1714}
1715
1716
1717
1718template <int dim, typename Number, bool is_face>
1719template <typename T>
1720inline void
1722 const T &in) const
1723{
1724 internal::process_cell_or_face_data(this->get_cell_ids(),
1725 array,
1726 in,
1727 [](const auto &local, auto &global) {
1728 global = local;
1729 });
1730}
1731
1732
1733
1734template <int dim, typename Number, bool is_face>
1735template <typename T>
1736inline T
1738 const AlignedVector<T> &array) const
1739{
1740 T out;
1741 internal::set_valid_element_to_array(this->get_cell_ids(), array, out);
1742 internal::process_cell_or_face_data(this->get_face_ids(),
1743 array,
1744 out,
1745 [](auto &local, const auto &global) {
1746 local = global;
1747 });
1748 return out;
1749}
1750
1751
1752
1753template <int dim, typename Number, bool is_face>
1754template <typename T>
1755inline void
1757 const T &in) const
1758{
1759 internal::process_cell_or_face_data(this->get_face_ids(),
1760 array,
1761 in,
1762 [](const auto &local, auto &global) {
1763 global = local;
1764 });
1765}
1766
1767
1768
1769#endif // ifndef DOXYGEN
1770
1771
1773
1774#endif
void resize_fast(const size_type new_size)
iterator begin()
size_type size() const
void resize(const size_type new_size)
void reinit(value_type *starting_element, const std::size_t n_elements)
Definition array_view.h:521
const unsigned int quad_no
const Number * J_value
internal::MatrixFreeFunctions::GeometryType get_cell_type() const
const Tensor< 2, dim, Number > * jacobian
const MappingInfoStorageType::QuadratureDescriptor * descriptor
const unsigned int n_quadrature_points
const MappingInfoStorageType * mapping_data
internal::MatrixFreeFunctions::DoFInfo::DoFAccessIndex dof_access_index
std::uint8_t get_face_no(const unsigned int v=0) const
ArrayView< Number > scratch_data
internal::MatrixFreeFunctions::DoFInfo::DoFAccessIndex get_dof_access_index() const
const Point< dim, Number > * quadrature_points
unsigned int subface_index
const Tensor< 1, dim *(dim+1)/2, Tensor< 1, dim, Number > > * jacobian_gradients_non_inverse
ScalarNumber shape_info_number_type
FEEvaluationData(const FEEvaluationData &other)=default
const ShapeInfoType & get_shape_info() const
void set_face_data(AlignedVector< T > &array, const T &value) const
internal::MatrixFreeFunctions::ShapeInfo< typename internal::VectorizedArrayTrait< VectorizedArrayType >::value_type > ShapeInfoType
void reinit_face(const internal::MatrixFreeFunctions::FaceToCellTopology< n_lanes > &face)
Number JxW(const unsigned int q_point) const
const std::array< unsigned int, n_lanes > & get_face_ids() const
unsigned int get_first_selected_component() const
const unsigned int n_fe_components
const ShapeInfoType * data
T read_face_data(const AlignedVector< T > &array) const
Number * gradients_from_hessians_quad
std::shared_ptr< internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number > > mapped_geometry
unsigned int get_active_quadrature_index() const
const internal::MatrixFreeFunctions::DoFInfo & get_dof_info() const
const std::array< unsigned int, n_lanes > & get_cell_ids() const
unsigned int get_mapping_data_index_offset() const
void set_data_pointers(AlignedVector< Number > *scratch_data, const unsigned int n_components)
Tensor< 2, dim, Number > inverse_jacobian(const unsigned int q_point) const
std::array< std::uint8_t, n_lanes > face_numbers
const unsigned int active_fe_index
static constexpr unsigned int n_lanes
const Tensor< 1, dim, Number > * normal_vectors
internal::MatrixFreeFunctions::GeometryType cell_type
const Number * begin_gradients() const
const std::array< unsigned int, n_lanes > & get_cell_or_face_ids() const
Tensor< 1, dim, Number > normal_vector(const unsigned int q_point) const
std::array< unsigned int, n_lanes > cell_ids
const DoFInfo * dof_info
unsigned int get_current_cell_index() const
unsigned int get_subface_index() const
Tensor< 1, dim, Number > get_normal_vector(const unsigned int q_point) const
const unsigned int active_quad_index
unsigned int get_active_fe_index() const
Number * values_from_gradients_quad
void set_cell_data(AlignedVector< T > &array, const T &value) const
unsigned int get_quadrature_index() const
std::array< std::uint8_t, n_lanes > face_orientations
const ScalarNumber * quadrature_weights
const Tensor< 1, dim *(dim+1)/2, Tensor< 1, dim, Number > > * jacobian_gradients
const Tensor< 1, dim, Number > * normal_x_jacobian
bool is_interior_face() const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > quadrature_point_indices() const
const std::vector< unsigned int > & get_internal_dof_numbering() const
FEEvaluationData & operator=(const FEEvaluationData &other)
virtual ~FEEvaluationData()=default
ArrayView< Number > get_scratch_data() const
T read_cell_data(const AlignedVector< T > &array) const
static constexpr unsigned int dimension
FEEvaluationData(const std::shared_ptr< internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number > > &mapping_data, const unsigned int n_fe_components, const unsigned int first_selected_component)
FEEvaluationData(const ShapeInfoType &shape_info, const bool is_interior_face=true)
Point< dim, Number > quadrature_point(const unsigned int q) const
const Number * begin_values() const
typename internal::VectorizedArrayTrait< Number >::value_type ScalarNumber
std::array< unsigned int, n_lanes > face_ids
FEEvaluationData(const InitializationData &initialization_data, const bool is_interior_face, const unsigned int quad_no, const unsigned int first_selected_component)
const Number * begin_dof_values() const
std::uint8_t get_face_orientation(const unsigned int v=0) const
const Number * begin_hessians() const
unsigned int get_cell_or_face_batch_id() const
const unsigned int first_selected_component
const unsigned int dofs_per_component
const unsigned int n_q_points
static constexpr unsigned int n_components
Definition point.h:111
const Point< dim > & point(const unsigned int i) const
#define DEAL_II_ALWAYS_INLINE
Definition config.h:109
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_DEPRECATED_WITH_COMMENT(comment)
Definition config.h:206
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
#define DeclException0(Exception0)
Definition exceptions.h:466
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcAccessToUninitializedField()
static ::ExceptionBase & ExcNotInitialized()
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:511
static ::ExceptionBase & ExcMatrixFreeAccessToUninitializedMappingField(std::string arg1)
static ::ExceptionBase & ExcMessage(std::string arg1)
std::vector< index_type > data
Definition mpi.cc:735
@ general
No special properties.
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:191
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
VectorType::value_type * begin(VectorType &V)
static const unsigned int invalid_unsigned_int
Definition types.h:220
boost::integer_range< IncrementableType > iota_view
Definition iota_view.h:45
STL namespace.
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
const MappingInfoStorageType * mapping_data
const MappingInfoStorageType::QuadratureDescriptor * descriptor
std::array< unsigned int, vectorization_width > cells_interior
Definition face_info.h:60
std::array< unsigned int, vectorization_width > cells_exterior
Definition face_info.h:75
static constexpr std::size_t width()