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_evaluation_data.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2020 - 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
17#ifndef dealii_matrix_free_fe_evaluation_data_h
18#define dealii_matrix_free_fe_evaluation_data_h
19
20
21#include <deal.II/base/config.h>
22
30#include <deal.II/base/tensor.h>
32
37
38
40
41
42
43namespace internal
44{
46
49 std::string,
50 << "You are requesting information from an FEEvaluation/FEFaceEvaluation "
51 << "object for which this kind of information has not been computed. What "
52 << "information these objects compute is determined by the update_* flags you "
53 << "pass to MatrixFree::reinit() via MatrixFree::AdditionalData. Here, "
54 << "the operation you are attempting requires the <" << arg1
55 << "> flag to be set, but it was apparently not specified "
56 << "upon initialization.");
57} // namespace internal
58
59// forward declarations
60template <int dim,
61 int fe_degree,
62 int n_q_points_1d = fe_degree + 1,
63 int n_components_ = 1,
64 typename Number = double,
65 typename VectorizedArrayType = VectorizedArray<Number>>
66class FEEvaluation;
67
68template <int dim,
69 int n_components_,
70 typename Number,
71 bool = false,
72 typename = VectorizedArray<Number>>
74
75namespace internal
76{
77 namespace MatrixFreeFunctions
78 {
79 template <int, typename>
80 class MappingDataOnTheFly;
81 }
82} // namespace internal
83
84
113template <int dim, typename Number, bool is_face>
115{
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;
125
128
129 static constexpr unsigned int n_lanes = sizeof(Number) / sizeof(ScalarNumber);
130
139 const bool is_interior_face = true);
140
144 FEEvaluationData(const FEEvaluationData &other) = default;
145
151
157 void
159 const unsigned int n_components);
160
165 void
168
173
182
187 Number
188 JxW(const unsigned int q_point) const;
189
195 quadrature_point(const unsigned int q) const;
196
209 inverse_jacobian(const unsigned int q_point) const;
210
224 get_normal_vector(const unsigned int q_point) const;
225
227
240 const Number *
242
251 Number *
253
264 const Number *
266
277 Number *
279
291 const Number *
293
305 Number *
307
320 const Number *
322
335 Number *
337
339
344
353 unsigned int
355
364
368 const ShapeInfoType &
370
376
382 const std::vector<unsigned int> &
384
388 unsigned int
390
394 unsigned int
396
401 unsigned int
403
408 unsigned int
410
414 unsigned int
416
429 std::uint8_t
430 get_face_no(const unsigned int v = 0) const;
431
441 unsigned int
443
452 std::uint8_t
453 get_face_orientation(const unsigned int v = 0) const;
454
464
473 bool
475
480 const std::array<unsigned int, n_lanes> &
482 {
483// implemented inline to avoid compilation problems on Windows
484#ifdef DEBUG
486#endif
487 return cell_ids;
488 }
489
494 const std::array<unsigned int, n_lanes> &
496 {
497// implemented inline to avoid compilation problems on Windows
498#ifdef DEBUG
500#endif
501 return face_ids;
502 }
503
508 unsigned int
510 {
511// implemented inline to avoid compilation problems on Windows
512#ifdef DEBUG
514#endif
515
516 return cell;
517 }
518
523 const std::array<unsigned int, n_lanes> &
525 {
526// implemented inline to avoid compilation problems on Windows
527#ifdef DEBUG
529#endif
530
531 if (!is_face || dof_access_index ==
533 return cell_ids;
534 else
535 return face_ids;
536 }
537
545
547
552
559 Number
561
568 void
569 set_cell_data(AlignedVector<Number> &array, const Number &value) const;
570
575 template <typename T>
576 std::array<T, Number::size()>
578 const AlignedVector<std::array<T, Number::size()>> &array) const;
579
584 template <typename T>
585 void
586 set_cell_data(AlignedVector<std::array<T, Number::size()>> &array,
587 const std::array<T, Number::size()> & value) const;
588
597 Number
599
608 void
609 set_face_data(AlignedVector<Number> &array, const Number &value) const;
610
615 template <typename T>
616 std::array<T, Number::size()>
618 const AlignedVector<std::array<T, Number::size()>> &array) const;
619
624 template <typename T>
625 void
626 set_face_data(AlignedVector<std::array<T, Number::size()>> &array,
627 const std::array<T, Number::size()> & value) const;
628
630
636 {
640 unsigned int active_fe_index;
641 unsigned int active_quad_index;
643 };
644
645protected:
652 FEEvaluationData(const InitializationData &initialization_data,
653 const bool is_interior_face,
654 const unsigned int quad_no,
655 const unsigned int first_selected_component);
656
662 const std::shared_ptr<
664 & mapping_data,
665 const unsigned int n_fe_components,
666 const unsigned int first_selected_component);
667
675
682
690
696 const unsigned int quad_no;
697
702 const unsigned int n_fe_components;
703
708 const unsigned int first_selected_component;
709
713 const unsigned int active_fe_index;
714
719 const unsigned int active_quad_index;
720
727
731 const unsigned int n_quadrature_points;
732
738
752
757 const Tensor<1, dim *(dim + 1) / 2, Tensor<1, dim, Number>>
759
766 const Number *J_value;
767
772
777
782
790
801 Number *values_dofs;
802
812 Number *values_quad;
813
826
839
850
856
863
870
877
884
891
898
905
910 unsigned int cell;
911
919
925
934 std::array<std::uint8_t, n_lanes> face_numbers;
935
944 std::array<std::uint8_t, n_lanes> face_orientations;
945
953 unsigned int subface_index;
954
962
967 std::array<unsigned int, n_lanes> cell_ids;
968
973 std::array<unsigned int, n_lanes> face_ids;
974
980 std::shared_ptr<
983
984 // Make FEEvaluation and FEEvaluationBase objects friends for access to
985 // protected member mapped_geometry.
986 template <int, int, typename, bool, typename>
987 friend class FEEvaluationBase;
988
989 template <int, int, int, int, typename, typename>
990 friend class FEEvaluation;
991};
992
993
998template <int dim,
999 typename Number,
1000 bool is_face,
1001 typename VectorizedArrayType = VectorizedArray<Number>>
1004
1005
1006/*----------------------- Inline functions ----------------------------------*/
1007
1008#ifndef DOXYGEN
1009
1010template <int dim, typename Number, bool is_face>
1012 const ShapeInfoType &shape_info,
1013 const bool is_interior_face)
1015 InitializationData{&shape_info, nullptr, nullptr, 0, 0, nullptr},
1017 0,
1018 0)
1019{}
1020
1021
1022template <int dim, typename Number, bool is_face>
1024 const InitializationData &initialization_data,
1025 const bool is_interior_face,
1026 const unsigned int quad_no,
1027 const unsigned int first_selected_component)
1028 : data(initialization_data.shape_info)
1029 , dof_info(initialization_data.dof_info)
1030 , mapping_data(initialization_data.mapping_data)
1031 , quad_no(quad_no)
1032 , n_fe_components(dof_info != nullptr ? dof_info->start_components.back() : 0)
1034 , active_fe_index(initialization_data.active_fe_index)
1035 , active_quad_index(initialization_data.active_quad_index)
1036 , descriptor(initialization_data.descriptor)
1037 , n_quadrature_points(descriptor == nullptr ?
1038 (is_face ? data->n_q_points_face : data->n_q_points) :
1040 , quadrature_points(nullptr)
1041 , jacobian(nullptr)
1042 , jacobian_gradients(nullptr)
1043 , J_value(nullptr)
1044 , normal_vectors(nullptr)
1045 , normal_x_jacobian(nullptr)
1047 descriptor != nullptr ? descriptor->quadrature_weights.begin() : nullptr)
1048# ifdef DEBUG
1049 , is_reinitialized(false)
1050 , dof_values_initialized(false)
1054 , values_quad_submitted(false)
1056# endif
1057 , cell(numbers::invalid_unsigned_int)
1060 is_face ?
1062 internal::MatrixFreeFunctions::DoFInfo::dof_access_face_interior :
1063 internal::MatrixFreeFunctions::DoFInfo::dof_access_face_exterior) :
1064 internal::MatrixFreeFunctions::DoFInfo::dof_access_cell)
1065 , subface_index(0)
1066 , cell_type(internal::MatrixFreeFunctions::general)
1067{}
1068
1069
1070
1071template <int dim, typename Number, bool is_face>
1073 const std::shared_ptr<
1076 const unsigned int n_fe_components,
1077 const unsigned int first_selected_component)
1078 : data(nullptr)
1079 , dof_info(nullptr)
1080 , mapping_data(nullptr)
1081 , quad_no(numbers::invalid_unsigned_int)
1084 , active_fe_index(numbers::invalid_unsigned_int)
1085 , active_quad_index(numbers::invalid_unsigned_int)
1086 , descriptor(nullptr)
1088 mapped_geometry->get_data_storage().descriptor[0].n_q_points)
1089 , quadrature_points(nullptr)
1090 , jacobian(nullptr)
1091 , jacobian_gradients(nullptr)
1092 , J_value(nullptr)
1093 , normal_vectors(nullptr)
1094 , normal_x_jacobian(nullptr)
1095 , quadrature_weights(nullptr)
1096 , cell(0)
1097 , cell_type(internal::MatrixFreeFunctions::general)
1098 , interior_face(true)
1099 , dof_access_index(internal::MatrixFreeFunctions::DoFInfo::dof_access_cell)
1101 , is_reinitialized(false)
1102{
1103 mapping_data = &mapped_geometry->get_data_storage();
1104 jacobian = mapped_geometry->get_data_storage().jacobians[0].begin();
1105 J_value = mapped_geometry->get_data_storage().JxW_values.begin();
1107 mapped_geometry->get_data_storage().jacobian_gradients[0].begin();
1109 mapped_geometry->get_data_storage().quadrature_points.begin();
1110}
1111
1112
1113
1114template <int dim, typename Number, bool is_face>
1117{
1124
1125 data = other.data;
1126 dof_info = other.dof_info;
1127 mapping_data = other.mapping_data;
1128 descriptor = other.descriptor;
1129 jacobian = nullptr;
1130 J_value = nullptr;
1131 normal_vectors = nullptr;
1132 normal_x_jacobian = nullptr;
1133 jacobian_gradients = nullptr;
1134 quadrature_points = nullptr;
1136
1137# ifdef DEBUG
1138 is_reinitialized = false;
1139 dof_values_initialized = false;
1143 values_quad_submitted = false;
1145# endif
1146
1150 is_face ?
1151 (is_interior_face() ?
1155 face_numbers[0] = 0;
1156 face_orientations[0] = 0;
1157 subface_index = 0;
1159
1160 return *this;
1161}
1162
1163
1164
1165template <int dim, typename Number, bool is_face>
1166inline void
1169 const unsigned int n_components)
1170{
1172
1173 const unsigned int tensor_dofs_per_component =
1174 Utilities::fixed_power<dim>(data->data.front().fe_degree + 1);
1176
1177 const unsigned int size_scratch_data =
1178 std::max(tensor_dofs_per_component + 1, dofs_per_component) * n_components *
1179 3 +
1181 const unsigned int size_data_arrays =
1183 (n_components * ((dim * (dim + 1)) / 2 + 2 * dim + 1) *
1185
1186 const unsigned int allocated_size = size_scratch_data + size_data_arrays;
1187 scratch_data_array->resize_fast(allocated_size);
1188 scratch_data.reinit(scratch_data_array->begin() + size_data_arrays,
1189 size_scratch_data);
1190
1191 // set the pointers to the correct position in the data array
1202}
1203
1204
1205
1206template <int dim, typename Number, bool is_face>
1207inline void
1210{
1211 Assert(is_face == true,
1212 ExcMessage("Faces can only be set if the is_face template parameter "
1213 "is true"));
1214 face_numbers[0] =
1216 subface_index = is_interior_face() == true ?
1218 face.subface_index;
1219
1220 // First check if interior or exterior cell has non-standard orientation
1221 // (i.e. the third bit is one or not). Then set zero if this cell has
1222 // standard-orientation else copy the first three bits
1223 // (which is equivalent to modulo 8). See also the documentation of
1224 // internal::MatrixFreeFunctions::FaceToCellTopology::face_orientation.
1225 face_orientations[0] = (is_interior_face() == (face.face_orientation >= 8)) ?
1226 (face.face_orientation % 8) :
1227 0;
1228
1229 if (is_interior_face())
1230 cell_ids = face.cells_interior;
1231 else
1232 cell_ids = face.cells_exterior;
1233}
1234
1235
1236
1237template <int dim, typename Number, bool is_face>
1240 const unsigned int q_point) const
1241{
1243 Assert(normal_vectors != nullptr,
1245 "update_normal_vectors"));
1247 return normal_vectors[0];
1248 else
1249 return normal_vectors[q_point];
1250}
1251
1252
1253
1254template <int dim, typename Number, bool is_face>
1255inline DEAL_II_ALWAYS_INLINE Number
1256FEEvaluationData<dim, Number, is_face>::JxW(const unsigned int q_point) const
1257{
1259 Assert(J_value != nullptr,
1261 "update_values|update_gradients"));
1263 {
1265 return J_value[0] * quadrature_weights[q_point];
1266 }
1267 else
1268 return J_value[q_point];
1269}
1270
1271
1272
1273template <int dim, typename Number, bool is_face>
1276 const unsigned int q) const
1277{
1279 Assert(this->quadrature_points != nullptr,
1281 "update_quadrature_points"));
1282
1283 // Cartesian/affine mesh: only first vertex of cell is stored, we must
1284 // compute it through the Jacobian (which is stored in non-inverted and
1285 // non-transposed form as index '1' in the jacobian field)
1286 if (is_face == false &&
1288 {
1289 Assert(this->jacobian != nullptr, ExcNotInitialized());
1291
1292 const Tensor<2, dim, Number> &jac = this->jacobian[1];
1294 for (unsigned int d = 0; d < dim; ++d)
1295 point[d] += jac[d][d] * static_cast<typename Number::value_type>(
1296 this->descriptor->quadrature.point(q)[d]);
1297 else
1298 for (unsigned int d = 0; d < dim; ++d)
1299 for (unsigned int e = 0; e < dim; ++e)
1300 point[d] += jac[d][e] * static_cast<typename Number::value_type>(
1301 this->descriptor->quadrature.point(q)[e]);
1302 return point;
1303 }
1304 else
1305 return this->quadrature_points[q];
1306}
1307
1308
1309
1310template <int dim, typename Number, bool is_face>
1313 const unsigned int q_point) const
1314{
1316 Assert(jacobian != nullptr,
1318 "update_gradients"));
1320 return jacobian[0];
1321 else
1322 return jacobian[q_point];
1323}
1324
1325
1326
1327template <int dim, typename Number, bool is_face>
1328inline const Number *
1330{
1331 return values_dofs;
1332}
1333
1334
1335
1336template <int dim, typename Number, bool is_face>
1337inline Number *
1339{
1340# ifdef DEBUG
1342# endif
1343 return values_dofs;
1344}
1345
1346
1347
1348template <int dim, typename Number, bool is_face>
1349inline const Number *
1351{
1352# ifdef DEBUG
1354# endif
1355 return values_quad;
1356}
1357
1358
1359
1360template <int dim, typename Number, bool is_face>
1361inline Number *
1363{
1364# ifdef DEBUG
1366 values_quad_submitted = true;
1367# endif
1368 return values_quad;
1369}
1370
1371
1372
1373template <int dim, typename Number, bool is_face>
1374inline const Number *
1376{
1377# ifdef DEBUG
1380# endif
1381 return gradients_quad;
1382}
1383
1384
1385
1386template <int dim, typename Number, bool is_face>
1387inline Number *
1389{
1390# ifdef DEBUG
1393# endif
1394 return gradients_quad;
1395}
1396
1397
1398
1399template <int dim, typename Number, bool is_face>
1400inline const Number *
1402{
1403# ifdef DEBUG
1405# endif
1406 return hessians_quad;
1407}
1408
1409
1410
1411template <int dim, typename Number, bool is_face>
1412inline Number *
1414{
1415# ifdef DEBUG
1417# endif
1418 return hessians_quad;
1419}
1420
1421
1422
1423template <int dim, typename Number, bool is_face>
1424inline unsigned int
1426{
1427 Assert(mapping_data != nullptr, ExcInternalError());
1428
1429 if (dof_info == nullptr)
1430 return 0;
1431 else
1432 {
1435 }
1436}
1437
1438
1439
1440template <int dim, typename Number, bool is_face>
1443{
1444# ifdef DEBUG
1446# endif
1447 return cell_type;
1448}
1449
1450
1451
1452template <int dim, typename Number, bool is_face>
1455{
1456 Assert(data != nullptr, ExcInternalError());
1457 return *data;
1458}
1459
1460
1461
1462template <int dim, typename Number, bool is_face>
1465{
1466 Assert(dof_info != nullptr,
1467 ExcMessage(
1468 "FEEvaluation was not initialized with a MatrixFree object!"));
1469 return *dof_info;
1470}
1471
1472
1473
1474template <int dim, typename Number, bool is_face>
1475inline const std::vector<unsigned int> &
1477{
1479}
1480
1481
1482
1483template <int dim, typename Number, bool is_face>
1484inline unsigned int
1486{
1487 return quad_no;
1488}
1489
1490
1491
1492template <int dim, typename Number, bool is_face>
1493inline unsigned int
1495{
1496 if (is_face && dof_access_index ==
1499 else
1500 return cell;
1501}
1502
1503
1504
1505template <int dim, typename Number, bool is_face>
1506inline unsigned int
1508{
1509 return active_fe_index;
1510}
1511
1512
1513
1514template <int dim, typename Number, bool is_face>
1515inline unsigned int
1517{
1518 return active_quad_index;
1519}
1520
1521
1522
1523template <int dim, typename Number, bool is_face>
1524inline unsigned int
1526{
1528}
1529
1530
1531
1532template <int dim, typename Number, bool is_face>
1533inline ArrayView<Number>
1535{
1536 return scratch_data;
1537}
1538
1539
1540
1541template <int dim, typename Number, bool is_face>
1542std::uint8_t
1544{
1545 Assert(is_face, ExcNotInitialized());
1549 is_interior_face() == false),
1550 ExcMessage("All face numbers can only be queried for ECL at exterior "
1551 "faces. Use get_face_no() in other cases."));
1552
1553 return face_numbers[v];
1554}
1555
1556
1557
1558template <int dim, typename Number, bool is_face>
1559inline unsigned int
1561{
1562 return subface_index;
1563}
1564
1565
1566
1567template <int dim, typename Number, bool is_face>
1568std::uint8_t
1570 const unsigned int v) const
1571{
1572 Assert(is_face, ExcNotInitialized());
1576 is_interior_face() == false),
1577 ExcMessage("All face numbers can only be queried for ECL at exterior "
1578 "faces. Use get_face_no() in other cases."));
1579
1580 return face_orientations[v];
1581}
1582
1583
1584
1585template <int dim, typename Number, bool is_face>
1588{
1589 return dof_access_index;
1590}
1591
1592
1593
1594template <int dim, typename Number, bool is_face>
1595inline bool
1597{
1598 return interior_face;
1599}
1600
1601
1602
1603template <int dim, typename Number, bool is_face>
1606{
1607 return {0U, n_quadrature_points};
1608}
1609
1610
1611
1612namespace internal
1613{
1614 template <std::size_t N,
1615 typename VectorizedArrayType2,
1616 typename GlobalVectorType,
1617 typename FU>
1618 inline void
1619 process_cell_or_face_data(const std::array<unsigned int, N> indices,
1620 GlobalVectorType & array,
1621 VectorizedArrayType2 & out,
1622 const FU & fu)
1623 {
1624 for (unsigned int i = 0; i < N; ++i)
1625 if (indices[i] != numbers::invalid_unsigned_int)
1626 {
1627 AssertIndexRange(indices[i] / N, array.size());
1628 fu(out[i], array[indices[i] / N][indices[i] % N]);
1629 }
1630 }
1631} // namespace internal
1632
1633
1634
1635template <int dim, typename Number, bool is_face>
1636inline Number
1638 const AlignedVector<Number> &array) const
1639{
1640 Number out = Number(1.);
1641 internal::process_cell_or_face_data(this->get_cell_ids(),
1642 array,
1643 out,
1644 [](auto &local, const auto &global) {
1645 local = global;
1646 });
1647 return out;
1648}
1649
1650
1651
1652template <int dim, typename Number, bool is_face>
1653inline void
1655 AlignedVector<Number> &array,
1656 const Number & in) const
1657{
1658 internal::process_cell_or_face_data(this->get_cell_ids(),
1659 array,
1660 in,
1661 [](const auto &local, auto &global) {
1662 global = local;
1663 });
1664}
1665
1666
1667
1668template <int dim, typename Number, bool is_face>
1669template <typename T>
1670inline std::array<T, Number::size()>
1672 const AlignedVector<std::array<T, Number::size()>> &array) const
1673{
1674 std::array<T, Number::size()> out;
1675 internal::process_cell_or_face_data(this->get_cell_ids(),
1676 array,
1677 out,
1678 [](auto &local, const auto &global) {
1679 local = global;
1680 });
1681 return out;
1682}
1683
1684
1685
1686template <int dim, typename Number, bool is_face>
1687template <typename T>
1688inline void
1690 AlignedVector<std::array<T, Number::size()>> &array,
1691 const std::array<T, Number::size()> & in) const
1692{
1693 internal::process_cell_or_face_data(this->get_cell_ids(),
1694 array,
1695 in,
1696 [](const auto &local, auto &global) {
1697 global = local;
1698 });
1699}
1700
1701
1702
1703template <int dim, typename Number, bool is_face>
1704inline Number
1706 const AlignedVector<Number> &array) const
1707{
1708 Number out = Number(1.);
1709 internal::process_cell_or_face_data(this->get_face_ids(),
1710 array,
1711 out,
1712 [](auto &local, const auto &global) {
1713 local = global;
1714 });
1715 return out;
1716}
1717
1718
1719
1720template <int dim, typename Number, bool is_face>
1721inline void
1723 AlignedVector<Number> &array,
1724 const Number & in) const
1725{
1726 internal::process_cell_or_face_data(this->get_face_ids(),
1727 array,
1728 in,
1729 [](const auto &local, auto &global) {
1730 global = local;
1731 });
1732}
1733
1734
1735
1736template <int dim, typename Number, bool is_face>
1737template <typename T>
1738inline std::array<T, Number::size()>
1740 const AlignedVector<std::array<T, Number::size()>> &array) const
1741{
1742 std::array<T, Number::size()> out;
1743 internal::process_cell_or_face_data(this->get_face_ids(),
1744 array,
1745 out,
1746 [](auto &local, const auto &global) {
1747 local = global;
1748 });
1749 return out;
1750}
1751
1752
1753
1754template <int dim, typename Number, bool is_face>
1755template <typename T>
1756inline void
1758 AlignedVector<std::array<T, Number::size()>> &array,
1759 const std::array<T, Number::size()> & in) const
1760{
1761 internal::process_cell_or_face_data(this->get_face_ids(),
1762 array,
1763 in,
1764 [](const auto &local, auto &global) {
1765 global = local;
1766 });
1767}
1768
1769
1770
1771#endif // ifndef DOXYGEN
1772
1773
1775
1776#endif
void resize_fast(const size_type new_size)
iterator begin()
size_type size() const
void reinit(value_type *starting_element, const std::size_t n_elements)
Definition: array_view.h:414
AlignedVector< VectorizedArrayType > * scratch_data_array
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
std::array< T, Number::size()> read_cell_data(const AlignedVector< std::array< T, Number::size()> > &array) const
Number * begin_gradients()
Number * begin_hessians()
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
Number * begin_values()
ArrayView< Number > scratch_data
internal::MatrixFreeFunctions::DoFInfo::DoFAccessIndex get_dof_access_index() const
const Point< dim, Number > * quadrature_points
unsigned int subface_index
FEEvaluationData(const FEEvaluationData &other)=default
const ShapeInfoType & get_shape_info() const
void reinit_face(const internal::MatrixFreeFunctions::FaceToCellTopology< n_lanes > &face)
Number read_face_data(const AlignedVector< Number > &array) const
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
void set_cell_data(AlignedVector< Number > &array, const Number &value) const
const ShapeInfoType * data
Number * gradients_from_hessians_quad
std::shared_ptr< internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number > > mapped_geometry
internal::MatrixFreeFunctions::DoFInfo DoFInfo
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
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
Number read_cell_data(const AlignedVector< Number > &array) const
const unsigned int active_quad_index
unsigned int get_active_fe_index() 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)
ArrayView< Number > get_scratch_data() const
Number * begin_dof_values()
void set_cell_data(AlignedVector< std::array< T, Number::size()> > &array, const std::array< T, Number::size()> &value) const
static constexpr unsigned int dimension
internal::MatrixFreeFunctions::ShapeInfo< VectorizedArrayType > ShapeInfoType
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
std::array< T, Number::size()> read_face_data(const AlignedVector< std::array< T, Number::size()> > &array) const
typename internal::VectorizedArrayTrait< Number >::value_type ScalarNumber
std::array< unsigned int, n_lanes > face_ids
void set_face_data(AlignedVector< std::array< T, Number::size()> > &array, const std::array< T, Number::size()> &value) const
FEEvaluationData(const InitializationData &initialization_data, const bool is_interior_face, const unsigned int quad_no, const unsigned int first_selected_component)
void set_face_data(AlignedVector< Number > &array, const Number &value) const
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
Definition: tensor.h:503
#define DEAL_II_ALWAYS_INLINE
Definition: config.h:102
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
#define DeclException0(Exception0)
Definition: exceptions.h:464
#define Assert(cond, exc)
Definition: exceptions.h:1473
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1667
#define AssertIndexRange(index, range)
Definition: exceptions.h:1732
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcAccessToUninitializedField()
static ::ExceptionBase & ExcNotInitialized()
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:509
static ::ExceptionBase & ExcMatrixFreeAccessToUninitializedMappingField(std::string arg1)
static ::ExceptionBase & ExcMessage(std::string arg1)
@ general
No special properties.
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition: utilities.cc:190
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:201
boost::integer_range< IncrementableType > iota_view
Definition: iota_view.h:46
::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:61
std::array< unsigned int, vectorization_width > cells_exterior
Definition: face_info.h:76
std::vector< UnivariateShapeData< Number > > data
Definition: shape_info.h:423
std::vector< unsigned int > lexicographic_numbering
Definition: shape_info.h:417