Reference documentation for deal.II version 9.6.0
\(\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_coupling_values.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2023 - 2024 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_coupling_values_h
17#define dealii_fe_coupling_values_h
18
19#include <deal.II/base/config.h>
20
22
27
29
31
32// Forward declaration
33#ifndef DOXYGEN
34template <int dim, int spacedim>
35class FEValuesBase;
36#endif
37
39{
51 {
92 const unsigned int n_inner_quadrature_points =
94 const std::vector<unsigned int> &dof_renumbering = {},
95 const std::vector<unsigned int> &quadrature_renumbering = {});
96
100 RenumberingData(const RenumberingData &other) = delete;
101
105 const unsigned int n_inner_dofs;
106
110 const unsigned int n_dofs;
111
116 const unsigned int n_inner_quadrature_points;
117
121 const unsigned int n_quadrature_points;
122
126 const std::vector<unsigned int> dof_renumbering;
127
131 const std::vector<unsigned int> quadrature_renumbering;
132
145 };
146
168 template <typename ViewType>
170 {
171 public:
175 using value_type = typename ViewType::value_type;
176
180 using gradient_type = typename ViewType::gradient_type;
181
186 template <typename Number>
188 typename ViewType::template solution_value_type<Number>;
189
194 template <typename Number>
196 typename ViewType::template solution_gradient_type<Number>;
197
211 RenumberedView(const ViewType &view, const RenumberingData &data);
212
230 value(const unsigned int shape_function, const unsigned int q_point) const;
231
249 gradient(const unsigned int shape_function,
250 const unsigned int q_point) const;
251
266 template <typename Number>
267 void
269 std::vector<solution_value_type<Number>> &values) const;
270
291 template <class InputVector>
292 void
294 const InputVector &dof_values,
296 &values) const;
297
313 template <typename Number>
314 void
316 const ReadVector<Number> &fe_function,
317 std::vector<solution_gradient_type<Number>> &gradients) const;
318
339 template <class InputVector>
340 void
342 const InputVector &dof_values,
344 &gradients) const;
345
346 private:
352
353
370 template <typename Number>
371 std::string
372 get_unique_container_name(const std::string &prefix,
373 const unsigned int size,
374 const Number &exemplar_number) const;
375
380 template <typename InputVector>
381 const InputVector &
382 outer_to_inner_dofs(const InputVector &outer_vector) const;
383
388 template <typename ValueType>
389 std::vector<ValueType> &
390 outer_to_inner_values(std::vector<ValueType> &outer_values) const;
391
397 template <typename ValueType>
398 void
399 inner_to_outer_values(const std::vector<ValueType> &inner_values,
400 std::vector<ValueType> &outer_values) const;
401
405 const ViewType &view;
406 };
407} // namespace FEValuesViews
408
409
435{
444
455 unrolled,
456
467 matching,
468
477 reorder,
478
489};
490
535{
543
557};
558
730template <int dim1, int dim2 = dim1, int spacedim = dim1>
732{
733public:
739
776 const FEValuesBase<dim1, spacedim> &fe_values_1,
777 const FEValuesBase<dim2, spacedim> &fe_values_2,
781
816 void
818 const FEValuesBase<dim1, spacedim> &fe_values_1,
819 const FEValuesBase<dim2, spacedim> &fe_values_2,
823
824
829 template <typename Extractor>
831 get_first_extractor(const Extractor &extractor) const;
832
837 template <typename Extractor>
839 get_second_extractor(const Extractor &extractor) const;
840
846 double
847 JxW(const unsigned int quadrature_point) const;
848
856 std::pair<Point<spacedim>, Point<spacedim>>
857 quadrature_point(const unsigned int quadrature_point) const;
858
868
876
884
899
940 std::vector<types::global_dof_index>
942 const std::vector<types::global_dof_index> &dof_indices_1,
943 const std::vector<types::global_dof_index> &dof_indices_2,
944 const types::global_dof_index dofs_offset_1 = 0,
945 const types::global_dof_index dofs_offset_2 = 0) const;
946
952 std::pair<unsigned int, unsigned int>
953 coupling_dof_to_dof_indices(const unsigned int coupling_dof_index) const;
954
960 std::pair<unsigned int, unsigned int>
962 const unsigned int quadrature_point) const;
963
974 template <typename Extractor>
978 const FEValuesExtractors::FirstCoupling<Extractor> &extractor) const;
979
985 template <typename Extractor>
989 const FEValuesExtractors::SecondCoupling<Extractor> &extractor) const;
990
1000 unsigned int
1002
1007 unsigned int
1009
1014 unsigned int
1016
1020 unsigned int
1022
1023private:
1028
1033
1038
1043
1047 std::unique_ptr<const FEValuesViews::RenumberingData> first_renumbering_data;
1048
1052 std::unique_ptr<const FEValuesViews::RenumberingData> second_renumbering_data;
1053
1058
1064 unsigned int n_coupling_dofs_;
1065};
1066
1067
1068#ifndef DOXYGEN
1069
1070
1071/*------------------------ Inline functions: namespace FEValuesViews --------*/
1072
1073namespace FEValuesViews
1074{
1076 const unsigned int n_inner_dofs,
1077 const unsigned int n_inner_quadrature_points,
1078 const std::vector<unsigned int> &dof_renumbering,
1079 const std::vector<unsigned int> &quadrature_renumbering)
1080 : n_inner_dofs(n_inner_dofs)
1081 , n_dofs(dof_renumbering.empty() ? n_inner_dofs : dof_renumbering.size())
1082 , n_inner_quadrature_points(n_inner_quadrature_points)
1083 , n_quadrature_points(quadrature_renumbering.empty() ?
1084 n_inner_quadrature_points :
1085 quadrature_renumbering.size())
1086 , dof_renumbering(dof_renumbering)
1087 , quadrature_renumbering(quadrature_renumbering)
1088 {
1089// Check that the renumbering vectors are valid.
1090# ifdef DEBUG
1091 // While for dofs we admit invalid values, this is not the case for
1092 // quadrature points.
1093 for (const auto i : dof_renumbering)
1094 Assert(i < n_inner_dofs || i == numbers::invalid_unsigned_int,
1095 ExcIndexRange(i, 0, n_inner_dofs));
1096
1097 for (const auto q : quadrature_renumbering)
1098 AssertIndexRange(q, n_inner_quadrature_points);
1099# endif
1100 }
1101
1102
1103
1104 template <typename ViewType>
1105 RenumberedView<ViewType>::RenumberedView(const ViewType &view,
1106 const RenumberingData &data)
1107 : view(view)
1108 , data(data)
1109 {}
1110
1111
1112
1113 template <typename ViewType>
1114 template <typename Number>
1115 inline std::string
1116 RenumberedView<ViewType>::get_unique_container_name(
1117 const std::string &prefix,
1118 const unsigned int size,
1119 const Number &exemplar_number) const
1120 {
1121 return prefix + "_" + Utilities::int_to_string(size) + "_" +
1122 Utilities::type_to_string(exemplar_number);
1123 }
1124
1125
1126
1127 template <typename ViewType>
1128 typename RenumberedView<ViewType>::value_type
1129 RenumberedView<ViewType>::value(const unsigned int shape_function,
1130 const unsigned int q_point) const
1131 {
1132 AssertIndexRange(shape_function, data.n_dofs);
1133 AssertIndexRange(q_point, data.n_quadrature_points);
1134
1135 const auto inner_shape_function = data.dof_renumbering.empty() ?
1136 shape_function :
1137 data.dof_renumbering[shape_function];
1138 const auto inner_q_point = data.quadrature_renumbering.empty() ?
1139 q_point :
1140 data.quadrature_renumbering[q_point];
1141 if (inner_shape_function == numbers::invalid_unsigned_int)
1142 return value_type(0);
1143 else
1144 {
1145 AssertIndexRange(inner_shape_function, data.n_inner_dofs);
1146 AssertIndexRange(inner_q_point, data.n_inner_quadrature_points);
1147 return view.value(inner_shape_function, inner_q_point);
1148 }
1149 }
1150
1151
1152
1153 template <typename ViewType>
1154 typename RenumberedView<ViewType>::gradient_type
1155 RenumberedView<ViewType>::gradient(const unsigned int shape_function,
1156 const unsigned int q_point) const
1157 {
1158 AssertIndexRange(shape_function, data.n_dofs);
1159 AssertIndexRange(q_point, data.n_quadrature_points);
1160
1161 const auto inner_shape_function = data.dof_renumbering.empty() ?
1162 shape_function :
1163 data.dof_renumbering[shape_function];
1164 const auto inner_q_point = data.quadrature_renumbering.empty() ?
1165 q_point :
1166 data.quadrature_renumbering[q_point];
1167 if (inner_shape_function == numbers::invalid_unsigned_int)
1168 return gradient_type();
1169 else
1170 return view.gradient(inner_shape_function, inner_q_point);
1171 }
1172
1173
1174
1175 template <typename ViewType>
1176 template <typename ValueType>
1177 std::vector<ValueType> &
1178 RenumberedView<ViewType>::outer_to_inner_values(
1179 std::vector<ValueType> &outer_values) const
1180 {
1181 AssertDimension(outer_values.size(), data.n_quadrature_points);
1182 if (data.quadrature_renumbering.empty())
1183 {
1184 return outer_values;
1185 }
1186 else
1187 {
1188 const auto name =
1189 get_unique_container_name("RenumberedView::outer_to_inner_values",
1190 data.n_inner_quadrature_points,
1191 outer_values[0]);
1192 auto &inner_values =
1193 data.data_storage.get()
1194 .template get_or_add_object_with_name<std::vector<ValueType>>(
1195 name, data.n_inner_quadrature_points);
1196 return inner_values;
1197 }
1198 }
1199
1200
1201
1202 template <typename ViewType>
1203 template <typename VectorType>
1204 const VectorType &
1205 RenumberedView<ViewType>::outer_to_inner_dofs(
1206 const VectorType &outer_dofs) const
1207 {
1208 AssertDimension(outer_dofs.size(), data.n_dofs);
1209 if (data.dof_renumbering.empty())
1210 {
1211 return outer_dofs;
1212 }
1213 else
1214 {
1215 const auto name =
1216 get_unique_container_name("RenumberedView::outer_to_inner_dofs",
1217 data.n_inner_dofs,
1218 outer_dofs[0]);
1219
1220 auto &inner_dofs = data.data_storage.get()
1221 .template get_or_add_object_with_name<VectorType>(
1222 name, data.n_inner_dofs);
1223 for (unsigned int i = 0; i < data.n_dofs; ++i)
1224 {
1225 const auto inner_i = data.dof_renumbering[i];
1226 if (inner_i != numbers::invalid_unsigned_int)
1227 {
1228 AssertIndexRange(inner_i, data.n_inner_dofs);
1229 inner_dofs[inner_i] = outer_dofs[i];
1230 }
1231 }
1232 return inner_dofs;
1233 }
1234 }
1235
1236
1237
1238 template <typename ViewType>
1239 template <typename ValueType>
1240 void
1241 RenumberedView<ViewType>::inner_to_outer_values(
1242 const std::vector<ValueType> &inner_values,
1243 std::vector<ValueType> &outer_values) const
1244 {
1245 AssertDimension(outer_values.size(), data.n_quadrature_points);
1246 AssertDimension(inner_values.size(), data.n_inner_quadrature_points);
1247 if (data.quadrature_renumbering.empty())
1248 {
1249 Assert(&inner_values == &outer_values, ExcInternalError());
1250 return;
1251 }
1252 for (unsigned int i = 0; i < data.quadrature_renumbering.size(); ++i)
1253 {
1254 outer_values[i] = inner_values[data.quadrature_renumbering[i]];
1255 }
1256 }
1257
1258
1259
1260 template <typename ViewType>
1261 template <typename Number>
1262 void
1263 RenumberedView<ViewType>::get_function_values(
1264 const ReadVector<Number> &fe_function,
1265 std::vector<solution_value_type<Number>> &values) const
1266 {
1267 auto &inner_values = outer_to_inner_values(values);
1268 view.get_function_values(fe_function, inner_values);
1269 inner_to_outer_values(inner_values, values);
1270 }
1271
1272
1273
1274 template <typename ViewType>
1275 template <typename InputVector>
1276 void
1277 RenumberedView<ViewType>::get_function_values_from_local_dof_values(
1278 const InputVector &dof_values,
1279 std::vector<solution_value_type<typename InputVector::value_type>> &values)
1280 const
1281 {
1282 const auto &inner_dof_values = outer_to_inner_dofs(dof_values);
1283 auto &inner_values = outer_to_inner_values(values);
1284
1285 view.get_function_values_from_local_dof_values(inner_dof_values,
1286 inner_values);
1287 inner_to_outer_values(inner_values, values);
1288 }
1289
1290
1291 template <typename ViewType>
1292 template <typename Number>
1293 void
1294 RenumberedView<ViewType>::get_function_gradients(
1295 const ReadVector<Number> &fe_function,
1296 std::vector<solution_gradient_type<Number>> &gradients) const
1297 {
1298 auto &inner_gradients = outer_to_inner_values(gradients);
1299 view.get_function_gradients(fe_function, inner_gradients);
1300 inner_to_outer_values(inner_gradients, gradients);
1301 }
1302
1303
1304
1305 template <typename ViewType>
1306 template <typename InputVector>
1307 void
1308 RenumberedView<ViewType>::get_function_gradients_from_local_dof_values(
1309 const InputVector &dof_values,
1310 std::vector<solution_gradient_type<typename InputVector::value_type>>
1311 &gradients) const
1312 {
1313 const auto &inner_dof_values = outer_to_inner_dofs(dof_values);
1314 auto &inner_gradients = outer_to_inner_values(gradients);
1315
1316 view.get_function_gradients_from_local_dof_values(inner_dof_values,
1317 inner_gradients);
1318 inner_to_outer_values(inner_gradients, gradients);
1319 }
1320} // namespace FEValuesViews
1321
1322/*-------------- Inline functions FECouplingValues ---------------------*/
1323
1324template <int dim1, int dim2, int spacedim>
1325std::vector<types::global_dof_index>
1327 const std::vector<types::global_dof_index> &dof_indices_1,
1328 const std::vector<types::global_dof_index> &dof_indices_2,
1329 const types::global_dof_index dofs_offset_1,
1330 const types::global_dof_index dofs_offset_2) const
1331{
1333 n_coupling_dofs_ != numbers::invalid_unsigned_int,
1334 ExcMessage(
1335 "Dofs are independent. You cannot ask for coupling dof indices."));
1336 AssertDimension(dof_indices_1.size(), first_fe_values->dofs_per_cell);
1337 AssertDimension(dof_indices_2.size(), second_fe_values->dofs_per_cell);
1338
1339 std::vector<types::global_dof_index> coupling_dof_indices(
1340 dof_indices_1.size() + dof_indices_2.size());
1341 unsigned int idx = 0;
1342 for (const auto &i : dof_indices_1)
1343 coupling_dof_indices[idx++] = i + dofs_offset_1;
1344 for (const auto &i : dof_indices_2)
1345 coupling_dof_indices[idx++] = i + dofs_offset_2;
1346 return coupling_dof_indices;
1347}
1348
1349
1350
1351template <int dim1, int dim2, int spacedim>
1353 : first_fe_values(nullptr)
1354 , second_fe_values(nullptr)
1355 , quadrature_coupling_type(QuadratureCouplingType::unrolled)
1356 , n_quadrature_points_(0)
1357 , n_coupling_dofs_(numbers::invalid_unsigned_int)
1358{}
1359
1360
1361
1362template <int dim1, int dim2, int spacedim>
1364 const FEValuesBase<dim1, spacedim> &fe_values_1,
1365 const FEValuesBase<dim2, spacedim> &fe_values_2,
1366 const DoFCouplingType &dof_coupling_type,
1367 const QuadratureCouplingType &quadrature_coupling_type)
1368 : FECouplingValues() // delegate to other constructor
1369{
1370 reinit(fe_values_1, fe_values_2, dof_coupling_type, quadrature_coupling_type);
1371}
1372
1373
1374
1375template <int dim1, int dim2, int spacedim>
1376void
1378 const FEValuesBase<dim1, spacedim> &fe_values_1,
1379 const FEValuesBase<dim2, spacedim> &fe_values_2,
1380 const DoFCouplingType &dof_coupling_type,
1381 const QuadratureCouplingType &quadrature_coupling_type)
1382{
1383 first_fe_values = &fe_values_1;
1384 second_fe_values = &fe_values_2;
1385 this->dof_coupling_type = dof_coupling_type;
1386 this->quadrature_coupling_type = quadrature_coupling_type;
1387
1388 // Gather information about the inner objects
1389 unsigned int first_n_inner_dofs = fe_values_1.dofs_per_cell;
1390 unsigned int second_n_inner_dofs = fe_values_2.dofs_per_cell;
1391
1392 unsigned int first_n_inner_quadrature_points =
1393 fe_values_1.n_quadrature_points;
1394 unsigned int second_n_inner_quadrature_points =
1395 fe_values_2.n_quadrature_points;
1396
1397 // Initialize counters and renumbering vectors to zero
1398 unsigned int first_n_dofs = 0;
1399 unsigned int second_n_dofs = 0;
1400
1401 unsigned int first_n_quadrature_points = 0;
1402 unsigned int second_n_quadrature_points = 0;
1403
1404 std::vector<unsigned int> first_dofs_map;
1405 std::vector<unsigned int> second_dofs_map;
1406
1407 std::vector<unsigned int> first_quad_map;
1408 std::vector<unsigned int> second_quad_map;
1409
1410 // Local relative tolerance for comparing quadrature points
1411 const double tol = std::min(fe_values_1.get_cell()->diameter(),
1412 fe_values_2.get_cell()->diameter()) *
1413 1e-6;
1414
1415 // Now fill the renumbering vectors
1416 switch (dof_coupling_type)
1417 {
1419 {
1420 n_coupling_dofs_ = numbers::invalid_unsigned_int;
1421 first_dofs_map.clear();
1422 second_dofs_map.clear();
1423 first_n_dofs = first_fe_values->dofs_per_cell;
1424 second_n_dofs = second_fe_values->dofs_per_cell;
1425 break;
1426 }
1428 {
1429 n_coupling_dofs_ =
1430 fe_values_1.dofs_per_cell + fe_values_2.dofs_per_cell;
1431
1432 first_n_dofs = n_coupling_dofs_;
1433 second_n_dofs = n_coupling_dofs_;
1434
1435 first_dofs_map.resize(n_coupling_dofs_);
1436 second_dofs_map.resize(n_coupling_dofs_);
1437
1438 unsigned int idx = 0;
1439 for (const unsigned int &i : fe_values_1.dof_indices())
1440 {
1441 first_dofs_map[idx] = i;
1442 second_dofs_map[idx++] = numbers::invalid_unsigned_int;
1443 }
1444 for (const unsigned int &i : fe_values_2.dof_indices())
1445 {
1446 first_dofs_map[idx] = numbers::invalid_unsigned_int;
1447 second_dofs_map[idx++] = i;
1448 }
1449 // Make sure we have the right number of dofs
1450 AssertDimension(idx, n_coupling_dofs_);
1451 break;
1452 }
1453 default:
1455 }
1456 // Compute the quadrature maps
1457 switch (quadrature_coupling_type)
1458 {
1460 {
1461 const auto &quadrature_points_1 = fe_values_1.get_quadrature_points();
1462 const auto &quadrature_points_2 = fe_values_2.get_quadrature_points();
1463
1464 n_quadrature_points_ =
1465 quadrature_points_1.size() * quadrature_points_2.size();
1466
1467 first_n_quadrature_points = n_quadrature_points_;
1468 second_n_quadrature_points = n_quadrature_points_;
1469
1470 first_quad_map.resize(n_quadrature_points_);
1471 second_quad_map.resize(n_quadrature_points_);
1472
1473 unsigned int idx = 0;
1474 for (const unsigned int &i : fe_values_1.quadrature_point_indices())
1475 for (const unsigned int &j : fe_values_2.quadrature_point_indices())
1476 {
1477 first_quad_map[idx] = i;
1478 second_quad_map[idx] = j;
1479 ++idx;
1480 }
1481 break;
1482 }
1484 {
1485 Assert(fe_values_1.get_quadrature_points().size() ==
1486 fe_values_2.get_quadrature_points().size(),
1487 ExcMessage("The two FEValuesBase objects must have the same "
1488 "number of quadrature points"));
1489
1490 n_quadrature_points_ = fe_values_1.get_quadrature_points().size();
1491
1492 first_n_quadrature_points = n_quadrature_points_;
1493 second_n_quadrature_points = n_quadrature_points_;
1494
1495 first_quad_map.clear();
1496 second_quad_map.clear();
1497
1498 break;
1499 }
1501 {
1502 const auto &quadrature_points_1 = fe_values_1.get_quadrature_points();
1503 const auto &quadrature_points_2 = fe_values_2.get_quadrature_points();
1504
1505 Assert(quadrature_points_1.size() == quadrature_points_2.size(),
1506 ExcMessage("The two FEValuesBase objects must have the same "
1507 "number of quadrature points"));
1508
1509 for (const unsigned int &i : fe_values_1.quadrature_point_indices())
1510 {
1511 Assert(quadrature_points_1[i].distance(quadrature_points_2[i]) <
1512 tol,
1513 ExcMessage(
1514 "The two FEValuesBase objects must have the same "
1515 "quadrature points"));
1516 }
1517
1518 n_quadrature_points_ = fe_values_1.get_quadrature_points().size();
1519
1520 first_n_quadrature_points = n_quadrature_points_;
1521 second_n_quadrature_points = n_quadrature_points_;
1522
1523 first_quad_map.clear();
1524 second_quad_map.clear();
1525
1526 break;
1527 }
1529 {
1530 const auto &quadrature_points_1 = fe_values_1.get_quadrature_points();
1531 const auto &quadrature_points_2 = fe_values_2.get_quadrature_points();
1532
1533 Assert(quadrature_points_1.size() == quadrature_points_2.size(),
1534 ExcMessage("The two FEValuesBase objects must have the same "
1535 "number of quadrature points"));
1536
1537 n_quadrature_points_ = fe_values_1.get_quadrature_points().size();
1538
1539 first_n_quadrature_points = n_quadrature_points_;
1540 second_n_quadrature_points = n_quadrature_points_;
1541
1542 // The first is the id. The second is renumbered.
1543 first_quad_map.clear();
1544 second_quad_map.resize(fe_values_1.get_quadrature_points().size());
1545
1546 // [TODO]: Avoid quadratic complexity here
1547 for (const unsigned int &i : fe_values_1.quadrature_point_indices())
1548 {
1550 for (const unsigned int &j :
1551 fe_values_2.quadrature_point_indices())
1552 if (quadrature_points_1[i].distance(quadrature_points_2[j]) <
1553 tol)
1554 {
1555 id = i;
1556 second_quad_map[i] = j;
1557 break;
1558 }
1560 ExcMessage(
1561 "The two FEValuesBase objects must have the same "
1562 "quadrature points, even if not in the same order."));
1563 }
1564 break;
1565 }
1567 {
1568 const auto &quadrature_points_1 = fe_values_1.get_quadrature_points();
1569 const auto &quadrature_points_2 = fe_values_2.get_quadrature_points();
1570
1571 // [TODO]: Avoid quadratic complexity here
1572 for (const unsigned int &i : fe_values_1.quadrature_point_indices())
1573 {
1574 for (const unsigned int &j :
1575 fe_values_2.quadrature_point_indices())
1576 if (quadrature_points_1[i].distance(quadrature_points_2[j]) <
1577 tol)
1578 {
1579 first_quad_map.emplace_back(i);
1580 second_quad_map.emplace_back(j);
1581 break;
1582 }
1583 }
1584 n_quadrature_points_ = first_quad_map.size();
1585
1586 first_n_quadrature_points = n_quadrature_points_;
1587 second_n_quadrature_points = n_quadrature_points_;
1588 break;
1589 }
1590 default:
1591 Assert(false, ExcInternalError());
1592 }
1593
1594 first_renumbering_data = std::make_unique<FEValuesViews::RenumberingData>(
1595 first_n_inner_dofs,
1596 first_n_inner_quadrature_points,
1597 first_dofs_map,
1598 first_quad_map);
1599
1600 second_renumbering_data = std::make_unique<FEValuesViews::RenumberingData>(
1601 second_n_inner_dofs,
1602 second_n_inner_quadrature_points,
1603 second_dofs_map,
1604 second_quad_map);
1605}
1606
1607
1608template <int dim1, int dim2, int spacedim>
1609inline double
1610FECouplingValues<dim1, dim2, spacedim>::JxW(const unsigned int q) const
1611{
1612 AssertIndexRange(q, n_quadrature_points_);
1613
1614 const auto first_q = first_renumbering_data->quadrature_renumbering.empty() ?
1615 q :
1616 first_renumbering_data->quadrature_renumbering[q];
1617
1618 if (quadrature_coupling_type == QuadratureCouplingType::tensor_product ||
1619 quadrature_coupling_type == QuadratureCouplingType::unrolled)
1620 {
1621 const auto second_q =
1622 second_renumbering_data->quadrature_renumbering.empty() ?
1623 q :
1624 second_renumbering_data->quadrature_renumbering[q];
1625 return first_fe_values->JxW(first_q) * second_fe_values->JxW(second_q);
1626 }
1627 else
1628 return first_fe_values->JxW(first_q);
1629}
1630
1631
1632
1633template <int dim1, int dim2, int spacedim>
1636{
1637 return {0U, n_quadrature_points_};
1638}
1639
1640
1641
1642template <int dim1, int dim2, int spacedim>
1645{
1646 AssertThrow(n_coupling_dofs_ != numbers::invalid_unsigned_int,
1647 ExcMessage(
1648 "Dofs are independent. You cannot ask for coupling dofs."));
1649 return {0U, n_coupling_dofs_};
1650}
1651
1652
1653
1654template <int dim1, int dim2, int spacedim>
1657{
1658 return {0U, n_first_dofs()};
1659}
1660
1661
1662
1663template <int dim1, int dim2, int spacedim>
1666{
1667 return {0U, n_second_dofs()};
1668}
1669
1670
1671
1672template <int dim1, int dim2, int spacedim>
1673inline unsigned int
1675{
1676 AssertThrow(n_coupling_dofs_ != numbers::invalid_unsigned_int,
1677 ExcMessage(
1678 "Dofs are independent. You cannot ask for coupling dofs."));
1679 return n_coupling_dofs_;
1680}
1681
1682
1683
1684template <int dim1, int dim2, int spacedim>
1685inline unsigned int
1687{
1688 return first_renumbering_data->n_dofs;
1689}
1690
1691
1692
1693template <int dim1, int dim2, int spacedim>
1694inline unsigned int
1696{
1697 return second_renumbering_data->n_dofs;
1698}
1699
1700
1701
1702template <int dim1, int dim2, int spacedim>
1703inline unsigned int
1705{
1706 return n_quadrature_points_;
1707}
1708
1709
1710
1711template <int dim1, int dim2, int spacedim>
1712std::pair<Point<spacedim>, Point<spacedim>>
1714 const unsigned int quadrature_point) const
1715{
1716 AssertIndexRange(quadrature_point, n_quadrature_points_);
1717 const auto first_q = first_fe_values->quadrature_point(
1718 first_renumbering_data->quadrature_renumbering.empty() ?
1719 quadrature_point :
1720 first_renumbering_data->quadrature_renumbering[quadrature_point]);
1721
1722 const auto second_q = second_fe_values->quadrature_point(
1723 second_renumbering_data->quadrature_renumbering.empty() ?
1724 quadrature_point :
1725 second_renumbering_data->quadrature_renumbering[quadrature_point]);
1726
1727 return {first_q, second_q};
1728}
1729
1730
1731
1732template <int dim1, int dim2, int spacedim>
1733std::pair<unsigned int, unsigned int>
1736 const unsigned int quadrature_point) const
1737{
1738 AssertIndexRange(quadrature_point, n_quadrature_points_);
1739 const auto first_id =
1740 first_renumbering_data->quadrature_renumbering.empty() ?
1741 quadrature_point :
1742 first_renumbering_data->quadrature_renumbering[quadrature_point];
1743
1744 const auto second_id =
1745 second_renumbering_data->quadrature_renumbering.empty() ?
1746 quadrature_point :
1747 second_renumbering_data->quadrature_renumbering[quadrature_point];
1748 return std::make_pair(first_id, second_id);
1749}
1750
1751
1752
1753template <int dim1, int dim2, int spacedim>
1754std::pair<unsigned int, unsigned int>
1756 const unsigned int coupling_dof_index) const
1757{
1758 AssertIndexRange(coupling_dof_index, n_coupling_dofs_);
1759 const auto first_id =
1760 first_renumbering_data->dof_renumbering.empty() ?
1761 coupling_dof_index :
1762 first_renumbering_data->dof_renumbering[coupling_dof_index];
1763
1764 const auto second_id =
1765 second_renumbering_data->dof_renumbering.empty() ?
1766 coupling_dof_index :
1767 second_renumbering_data->dof_renumbering[coupling_dof_index];
1768 return std::make_pair(first_id, second_id);
1769}
1770
1771
1772
1773template <int dim1, int dim2, int spacedim>
1774template <typename Extractor>
1777 const Extractor &extractor) const
1778{
1780}
1781
1782
1783
1784template <int dim1, int dim2, int spacedim>
1785template <typename Extractor>
1788 const Extractor &extractor) const
1789{
1791}
1792
1793
1794
1795template <int dim1, int dim2, int spacedim>
1796template <typename Extractor>
1800 const FEValuesExtractors::FirstCoupling<Extractor> &extractor) const
1801{
1804 (*first_fe_values)[extractor.extractor], *first_renumbering_data);
1805}
1806
1807
1808
1809template <int dim1, int dim2, int spacedim>
1810template <typename Extractor>
1814 const FEValuesExtractors::SecondCoupling<Extractor> &extractor) const
1815{
1818 (*second_fe_values)[extractor.extractor], *second_renumbering_data);
1819}
1820
1821#endif // DOXYGEN
1822
1824
1825#endif
unsigned int n_quadrature_points() const
std::pair< Point< spacedim >, Point< spacedim > > quadrature_point(const unsigned int quadrature_point) const
const FEValuesViews::RenumberedView< FEValuesViews::View< dim2, spacedim, Extractor > > operator[](const FEValuesExtractors::SecondCoupling< Extractor > &extractor) const
std::vector< types::global_dof_index > get_coupling_dof_indices(const std::vector< types::global_dof_index > &dof_indices_1, const std::vector< types::global_dof_index > &dof_indices_2, const types::global_dof_index dofs_offset_1=0, const types::global_dof_index dofs_offset_2=0) const
FEValuesExtractors::FirstCoupling< Extractor > get_first_extractor(const Extractor &extractor) const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > first_dof_indices() const
FECouplingValues(const FEValuesBase< dim1, spacedim > &fe_values_1, const FEValuesBase< dim2, spacedim > &fe_values_2, const DoFCouplingType &dof_coupling_type=DoFCouplingType::independent, const QuadratureCouplingType &quadrature_coupling_type=QuadratureCouplingType::tensor_product)
unsigned int n_quadrature_points_
std_cxx20::ranges::iota_view< unsigned int, unsigned int > quadrature_point_indices() const
unsigned int n_coupling_dofs_
SmartPointer< const FEValuesBase< dim1, spacedim > > first_fe_values
unsigned int n_coupling_dofs() const
QuadratureCouplingType quadrature_coupling_type
const FEValuesViews::RenumberedView< FEValuesViews::View< dim1, spacedim, Extractor > > operator[](const FEValuesExtractors::FirstCoupling< Extractor > &extractor) const
std::pair< unsigned int, unsigned int > coupling_quadrature_to_quadrature_indices(const unsigned int quadrature_point) const
double JxW(const unsigned int quadrature_point) const
std::unique_ptr< const FEValuesViews::RenumberingData > first_renumbering_data
std_cxx20::ranges::iota_view< unsigned int, unsigned int > coupling_dof_indices() const
unsigned int n_second_dofs() const
FEValuesExtractors::SecondCoupling< Extractor > get_second_extractor(const Extractor &extractor) const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > second_dof_indices() const
DoFCouplingType dof_coupling_type
unsigned int n_first_dofs() const
void reinit(const FEValuesBase< dim1, spacedim > &fe_values_1, const FEValuesBase< dim2, spacedim > &fe_values_2, const DoFCouplingType &dof_coupling_type=DoFCouplingType::independent, const QuadratureCouplingType &quadrature_coupling_type=QuadratureCouplingType::tensor_product)
std::pair< unsigned int, unsigned int > coupling_dof_to_dof_indices(const unsigned int coupling_dof_index) const
SmartPointer< const FEValuesBase< dim2, spacedim > > second_fe_values
std::unique_ptr< const FEValuesViews::RenumberingData > second_renumbering_data
Triangulation< dim, spacedim >::cell_iterator get_cell() const
const std::vector< Point< spacedim > > & get_quadrature_points() const
const unsigned int dofs_per_cell
const unsigned int n_quadrature_points
void get_function_values_from_local_dof_values(const InputVector &dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
void get_function_gradients_from_local_dof_values(const InputVector &dof_values, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
const InputVector & outer_to_inner_dofs(const InputVector &outer_vector) const
gradient_type gradient(const unsigned int shape_function, const unsigned int q_point) const
typename ViewType::template solution_gradient_type< Number > solution_gradient_type
typename ViewType::value_type value_type
typename ViewType::template solution_value_type< Number > solution_value_type
void get_function_gradients(const ReadVector< Number > &fe_function, std::vector< solution_gradient_type< Number > > &gradients) const
void inner_to_outer_values(const std::vector< ValueType > &inner_values, std::vector< ValueType > &outer_values) const
std::vector< ValueType > & outer_to_inner_values(std::vector< ValueType > &outer_values) const
typename ViewType::gradient_type gradient_type
void get_function_values(const ReadVector< Number > &fe_function, std::vector< solution_value_type< Number > > &values) const
RenumberedView(const ViewType &view, const RenumberingData &data)
std::string get_unique_container_name(const std::string &prefix, const unsigned int size, const Number &exemplar_number) const
value_type value(const unsigned int shape_function, const unsigned int q_point) const
Definition point.h:111
A class that provides a separate storage location on each thread that accesses the object.
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
QuadratureCouplingType
DoFCouplingType
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcIndexRange(std::size_t arg1, std::size_t arg2, std::size_t arg3)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
typename ::internal::FEValuesViews:: ViewType< dim, spacedim, Extractor >::type View
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
std::string type_to_string(const T &t)
Definition utilities.h:1023
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition utilities.cc:470
static const unsigned int invalid_unsigned_int
Definition types.h:220
boost::integer_range< IncrementableType > iota_view
Definition iota_view.h:45
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
const std::vector< unsigned int > dof_renumbering
RenumberingData(const unsigned int n_inner_dofs=numbers::invalid_unsigned_int, const unsigned int n_inner_quadrature_points=numbers::invalid_unsigned_int, const std::vector< unsigned int > &dof_renumbering={}, const std::vector< unsigned int > &quadrature_renumbering={})
const unsigned int n_inner_quadrature_points
const std::vector< unsigned int > quadrature_renumbering
RenumberingData(const RenumberingData &other)=delete
Threads::ThreadLocalStorage< GeneralDataStorage > data_storage