15#ifndef dealii_fe_coupling_values_h
16#define dealii_fe_coupling_values_h
33template <
int dim,
int spacedim>
167 template <
typename ViewType>
185 template <
typename Number>
193 template <
typename Number>
229 value(
const unsigned int shape_function,
const unsigned int q_point)
const;
249 const unsigned int q_point)
const;
265 template <
typename Number>
290 template <
class InputVector>
293 const InputVector &dof_values,
312 template <
typename Number>
338 template <
class InputVector>
341 const InputVector &dof_values,
369 template <
typename Number>
372 const unsigned int size,
373 const Number &exemplar_number)
const;
379 template <
typename InputVector>
387 template <
typename ValueType>
388 std::vector<ValueType> &
396 template <
typename ValueType>
399 std::vector<ValueType> &outer_values)
const;
729template <
int dim1,
int dim2 = dim1,
int spacedim = dim1>
828 template <
typename Extractor>
836 template <
typename Extractor>
939 std::vector<types::global_dof_index>
941 const std::vector<types::global_dof_index> &dof_indices_1,
942 const std::vector<types::global_dof_index> &dof_indices_2,
951 std::pair<unsigned int, unsigned int>
959 std::pair<unsigned int, unsigned int>
973 template <
typename Extractor>
984 template <
typename Extractor>
1075 const unsigned int n_inner_dofs,
1076 const unsigned int n_inner_quadrature_points,
1077 const std::vector<unsigned int> &dof_renumbering,
1078 const std::vector<unsigned int> &quadrature_renumbering)
1079 : n_inner_dofs(n_inner_dofs)
1080 , n_dofs(dof_renumbering.empty() ? n_inner_dofs : dof_renumbering.size())
1081 , n_inner_quadrature_points(n_inner_quadrature_points)
1082 , n_quadrature_points(quadrature_renumbering.empty() ?
1083 n_inner_quadrature_points :
1084 quadrature_renumbering.size())
1085 , dof_renumbering(dof_renumbering)
1086 , quadrature_renumbering(quadrature_renumbering)
1092 for (
const auto i : dof_renumbering)
1096 for (
const auto q : quadrature_renumbering)
1103 template <
typename ViewType>
1105 const RenumberingData &
data)
1112 template <
typename ViewType>
1113 template <
typename Number>
1115 RenumberedView<ViewType>::get_unique_container_name(
1116 const std::string &prefix,
1117 const unsigned int size,
1118 const Number &exemplar_number)
const
1126 template <
typename ViewType>
1127 typename RenumberedView<ViewType>::value_type
1128 RenumberedView<ViewType>::value(
const unsigned int shape_function,
1129 const unsigned int q_point)
const
1134 const auto inner_shape_function =
data.dof_renumbering.empty() ?
1136 data.dof_renumbering[shape_function];
1137 const auto inner_q_point =
data.quadrature_renumbering.empty() ?
1139 data.quadrature_renumbering[q_point];
1141 return value_type(0);
1146 return view.value(inner_shape_function, inner_q_point);
1152 template <
typename ViewType>
1153 typename RenumberedView<ViewType>::gradient_type
1154 RenumberedView<ViewType>::gradient(
const unsigned int shape_function,
1155 const unsigned int q_point)
const
1160 const auto inner_shape_function =
data.dof_renumbering.empty() ?
1162 data.dof_renumbering[shape_function];
1163 const auto inner_q_point =
data.quadrature_renumbering.empty() ?
1165 data.quadrature_renumbering[q_point];
1167 return gradient_type();
1169 return view.gradient(inner_shape_function, inner_q_point);
1174 template <
typename ViewType>
1175 template <
typename ValueType>
1176 std::vector<ValueType> &
1177 RenumberedView<ViewType>::outer_to_inner_values(
1178 std::vector<ValueType> &outer_values)
const
1181 if (
data.quadrature_renumbering.empty())
1183 return outer_values;
1188 get_unique_container_name(
"RenumberedView::outer_to_inner_values",
1189 data.n_inner_quadrature_points,
1191 auto &inner_values =
1192 data.data_storage.get()
1193 .template get_or_add_object_with_name<std::vector<ValueType>>(
1194 name,
data.n_inner_quadrature_points);
1195 return inner_values;
1201 template <
typename ViewType>
1202 template <
typename VectorType>
1204 RenumberedView<ViewType>::outer_to_inner_dofs(
1205 const VectorType &outer_dofs)
const
1208 if (
data.dof_renumbering.empty())
1215 get_unique_container_name(
"RenumberedView::outer_to_inner_dofs",
1219 auto &inner_dofs =
data.data_storage.get()
1220 .template get_or_add_object_with_name<VectorType>(
1221 name,
data.n_inner_dofs);
1222 for (
unsigned int i = 0; i <
data.n_dofs; ++i)
1224 const auto inner_i =
data.dof_renumbering[i];
1228 inner_dofs[inner_i] = outer_dofs[i];
1237 template <
typename ViewType>
1238 template <
typename ValueType>
1240 RenumberedView<ViewType>::inner_to_outer_values(
1241 const std::vector<ValueType> &inner_values,
1242 std::vector<ValueType> &outer_values)
const
1246 if (
data.quadrature_renumbering.empty())
1251 for (
unsigned int i = 0; i <
data.quadrature_renumbering.size(); ++i)
1253 outer_values[i] = inner_values[
data.quadrature_renumbering[i]];
1259 template <
typename ViewType>
1260 template <
typename Number>
1262 RenumberedView<ViewType>::get_function_values(
1264 std::vector<solution_value_type<Number>> &values)
const
1266 auto &inner_values = outer_to_inner_values(values);
1267 view.get_function_values(fe_function, inner_values);
1268 inner_to_outer_values(inner_values, values);
1273 template <
typename ViewType>
1274 template <
typename InputVector>
1276 RenumberedView<ViewType>::get_function_values_from_local_dof_values(
1277 const InputVector &dof_values,
1278 std::vector<solution_value_type<typename InputVector::value_type>> &values)
1281 const auto &inner_dof_values = outer_to_inner_dofs(dof_values);
1282 auto &inner_values = outer_to_inner_values(values);
1284 view.get_function_values_from_local_dof_values(inner_dof_values,
1286 inner_to_outer_values(inner_values, values);
1290 template <
typename ViewType>
1291 template <
typename Number>
1293 RenumberedView<ViewType>::get_function_gradients(
1295 std::vector<solution_gradient_type<Number>> &gradients)
const
1297 auto &inner_gradients = outer_to_inner_values(gradients);
1298 view.get_function_gradients(fe_function, inner_gradients);
1299 inner_to_outer_values(inner_gradients, gradients);
1304 template <
typename ViewType>
1305 template <
typename InputVector>
1307 RenumberedView<ViewType>::get_function_gradients_from_local_dof_values(
1308 const InputVector &dof_values,
1309 std::vector<solution_gradient_type<typename InputVector::value_type>>
1312 const auto &inner_dof_values = outer_to_inner_dofs(dof_values);
1313 auto &inner_gradients = outer_to_inner_values(gradients);
1315 view.get_function_gradients_from_local_dof_values(inner_dof_values,
1317 inner_to_outer_values(inner_gradients, gradients);
1323template <
int dim1,
int dim2,
int spacedim>
1324std::vector<types::global_dof_index>
1326 const std::vector<types::global_dof_index> &dof_indices_1,
1327 const std::vector<types::global_dof_index> &dof_indices_2,
1334 "Dofs are independent. You cannot ask for coupling dof indices."));
1335 AssertDimension(dof_indices_1.size(), first_fe_values->dofs_per_cell);
1336 AssertDimension(dof_indices_2.size(), second_fe_values->dofs_per_cell);
1338 std::vector<types::global_dof_index> coupling_dof_indices(
1339 dof_indices_1.size() + dof_indices_2.size());
1340 unsigned int idx = 0;
1341 for (
const auto &i : dof_indices_1)
1342 coupling_dof_indices[idx++] = i + dofs_offset_1;
1343 for (
const auto &i : dof_indices_2)
1344 coupling_dof_indices[idx++] = i + dofs_offset_2;
1345 return coupling_dof_indices;
1350template <
int dim1,
int dim2,
int spacedim>
1352 : first_fe_values(nullptr)
1353 , second_fe_values(nullptr)
1355 , n_quadrature_points_(0)
1356 , n_coupling_dofs_(
numbers::invalid_unsigned_int)
1361template <
int dim1,
int dim2,
int spacedim>
1374template <
int dim1,
int dim2,
int spacedim>
1382 first_fe_values = &fe_values_1;
1383 second_fe_values = &fe_values_2;
1384 this->dof_coupling_type = dof_coupling_type;
1385 this->quadrature_coupling_type = quadrature_coupling_type;
1388 unsigned int first_n_inner_dofs = fe_values_1.
dofs_per_cell;
1389 unsigned int second_n_inner_dofs = fe_values_2.
dofs_per_cell;
1391 unsigned int first_n_inner_quadrature_points =
1393 unsigned int second_n_inner_quadrature_points =
1397 unsigned int first_n_dofs = 0;
1398 unsigned int second_n_dofs = 0;
1400 unsigned int first_n_quadrature_points = 0;
1401 unsigned int second_n_quadrature_points = 0;
1403 std::vector<unsigned int> first_dofs_map;
1404 std::vector<unsigned int> second_dofs_map;
1406 std::vector<unsigned int> first_quad_map;
1407 std::vector<unsigned int> second_quad_map;
1411 fe_values_2.
get_cell()->diameter()) *
1415 switch (dof_coupling_type)
1420 first_dofs_map.clear();
1421 second_dofs_map.clear();
1422 first_n_dofs = first_fe_values->dofs_per_cell;
1423 second_n_dofs = second_fe_values->dofs_per_cell;
1431 first_n_dofs = n_coupling_dofs_;
1432 second_n_dofs = n_coupling_dofs_;
1434 first_dofs_map.resize(n_coupling_dofs_);
1435 second_dofs_map.resize(n_coupling_dofs_);
1437 unsigned int idx = 0;
1438 for (
const unsigned int &i : fe_values_1.dof_indices())
1440 first_dofs_map[idx] = i;
1443 for (
const unsigned int &i : fe_values_2.dof_indices())
1446 second_dofs_map[idx++] = i;
1456 switch (quadrature_coupling_type)
1463 n_quadrature_points_ =
1464 quadrature_points_1.size() * quadrature_points_2.size();
1466 first_n_quadrature_points = n_quadrature_points_;
1467 second_n_quadrature_points = n_quadrature_points_;
1469 first_quad_map.resize(n_quadrature_points_);
1470 second_quad_map.resize(n_quadrature_points_);
1472 unsigned int idx = 0;
1473 for (
const unsigned int &i : fe_values_1.quadrature_point_indices())
1474 for (const unsigned
int &j : fe_values_2.quadrature_point_indices())
1476 first_quad_map[idx] = i;
1477 second_quad_map[idx] = j;
1486 ExcMessage(
"The two FEValuesBase objects must have the same "
1487 "number of quadrature points"));
1491 first_n_quadrature_points = n_quadrature_points_;
1492 second_n_quadrature_points = n_quadrature_points_;
1494 first_quad_map.clear();
1495 second_quad_map.clear();
1504 Assert(quadrature_points_1.size() == quadrature_points_2.size(),
1505 ExcMessage(
"The two FEValuesBase objects must have the same "
1506 "number of quadrature points"));
1508 for (
const unsigned int &i : fe_values_1.quadrature_point_indices())
1510 Assert(quadrature_points_1[i].distance(quadrature_points_2[i]) <
1513 "The two FEValuesBase objects must have the same "
1514 "quadrature points"));
1519 first_n_quadrature_points = n_quadrature_points_;
1520 second_n_quadrature_points = n_quadrature_points_;
1522 first_quad_map.clear();
1523 second_quad_map.clear();
1532 Assert(quadrature_points_1.size() == quadrature_points_2.size(),
1533 ExcMessage(
"The two FEValuesBase objects must have the same "
1534 "number of quadrature points"));
1538 first_n_quadrature_points = n_quadrature_points_;
1539 second_n_quadrature_points = n_quadrature_points_;
1542 first_quad_map.clear();
1546 for (
const unsigned int &i : fe_values_1.quadrature_point_indices())
1549 for (
const unsigned int &j :
1550 fe_values_2.quadrature_point_indices())
1551 if (quadrature_points_1[i].distance(quadrature_points_2[j]) <
1555 second_quad_map[i] = j;
1560 "The two FEValuesBase objects must have the same "
1561 "quadrature points, even if not in the same order."));
1571 for (
const unsigned int &i : fe_values_1.quadrature_point_indices())
1573 for (
const unsigned int &j :
1574 fe_values_2.quadrature_point_indices())
1575 if (quadrature_points_1[i].distance(quadrature_points_2[j]) <
1578 first_quad_map.emplace_back(i);
1579 second_quad_map.emplace_back(j);
1583 n_quadrature_points_ = first_quad_map.size();
1585 first_n_quadrature_points = n_quadrature_points_;
1586 second_n_quadrature_points = n_quadrature_points_;
1593 first_renumbering_data = std::make_unique<FEValuesViews::RenumberingData>(
1595 first_n_inner_quadrature_points,
1599 second_renumbering_data = std::make_unique<FEValuesViews::RenumberingData>(
1600 second_n_inner_dofs,
1601 second_n_inner_quadrature_points,
1607template <
int dim1,
int dim2,
int spacedim>
1613 const auto first_q = first_renumbering_data->quadrature_renumbering.empty() ?
1615 first_renumbering_data->quadrature_renumbering[q];
1620 const auto second_q =
1621 second_renumbering_data->quadrature_renumbering.empty() ?
1623 second_renumbering_data->quadrature_renumbering[q];
1624 return first_fe_values->JxW(first_q) * second_fe_values->JxW(second_q);
1627 return first_fe_values->JxW(first_q);
1632template <
int dim1,
int dim2,
int spacedim>
1636 return {0U, n_quadrature_points_};
1641template <
int dim1,
int dim2,
int spacedim>
1647 "Dofs are independent. You cannot ask for coupling dofs."));
1648 return {0U, n_coupling_dofs_};
1653template <
int dim1,
int dim2,
int spacedim>
1657 return {0U, n_first_dofs()};
1662template <
int dim1,
int dim2,
int spacedim>
1666 return {0U, n_second_dofs()};
1671template <
int dim1,
int dim2,
int spacedim>
1677 "Dofs are independent. You cannot ask for coupling dofs."));
1678 return n_coupling_dofs_;
1683template <
int dim1,
int dim2,
int spacedim>
1687 return first_renumbering_data->n_dofs;
1692template <
int dim1,
int dim2,
int spacedim>
1696 return second_renumbering_data->n_dofs;
1701template <
int dim1,
int dim2,
int spacedim>
1705 return n_quadrature_points_;
1710template <
int dim1,
int dim2,
int spacedim>
1713 const unsigned int quadrature_point)
const
1716 const auto first_q = first_fe_values->quadrature_point(
1717 first_renumbering_data->quadrature_renumbering.empty() ?
1719 first_renumbering_data->quadrature_renumbering[quadrature_point]);
1721 const auto second_q = second_fe_values->quadrature_point(
1722 second_renumbering_data->quadrature_renumbering.empty() ?
1724 second_renumbering_data->quadrature_renumbering[quadrature_point]);
1726 return {first_q, second_q};
1731template <
int dim1,
int dim2,
int spacedim>
1732std::pair<unsigned int, unsigned int>
1735 const unsigned int quadrature_point)
const
1738 const auto first_id =
1739 first_renumbering_data->quadrature_renumbering.empty() ?
1741 first_renumbering_data->quadrature_renumbering[quadrature_point];
1743 const auto second_id =
1744 second_renumbering_data->quadrature_renumbering.empty() ?
1746 second_renumbering_data->quadrature_renumbering[quadrature_point];
1747 return std::make_pair(first_id, second_id);
1752template <
int dim1,
int dim2,
int spacedim>
1753std::pair<unsigned int, unsigned int>
1755 const unsigned int coupling_dof_index)
const
1758 const auto first_id =
1759 first_renumbering_data->dof_renumbering.empty() ?
1760 coupling_dof_index :
1761 first_renumbering_data->dof_renumbering[coupling_dof_index];
1763 const auto second_id =
1764 second_renumbering_data->dof_renumbering.empty() ?
1765 coupling_dof_index :
1766 second_renumbering_data->dof_renumbering[coupling_dof_index];
1767 return std::make_pair(first_id, second_id);
1772template <
int dim1,
int dim2,
int spacedim>
1773template <
typename Extractor>
1776 const Extractor &extractor)
const
1783template <
int dim1,
int dim2,
int spacedim>
1784template <
typename Extractor>
1787 const Extractor &extractor)
const
1794template <
int dim1,
int dim2,
int spacedim>
1795template <
typename Extractor>
1803 (*first_fe_values)[extractor.
extractor], *first_renumbering_data);
1808template <
int dim1,
int dim2,
int spacedim>
1809template <
typename Extractor>
1817 (*second_fe_values)[extractor.
extractor], *second_renumbering_data);
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_
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
ObserverPointer< const FEValuesBase< dim1, spacedim > > first_fe_values
double JxW(const unsigned int quadrature_point) const
ObserverPointer< const FEValuesBase< dim2, spacedim > > second_fe_values
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
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::value_type 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)
typename ViewType::template solution_value_type< Number > solution_value_type
std::string get_unique_container_name(const std::string &prefix, const unsigned int size, const Number &exemplar_number) const
typename ViewType::template solution_gradient_type< Number > solution_gradient_type
value_type value(const unsigned int shape_function, const unsigned int q_point) const
const RenumberingData & data
A class that provides a separate storage location on each thread that accesses the object.
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
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)
std::vector< index_type > data
typename ::internal::FEValuesViews::ViewType< dim, spacedim, Extractor >::type View
Tpetra::Vector< Number, LO, GO, NodeType< MemorySpace > > VectorType
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
std::string type_to_string(const T &t)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
static const unsigned int invalid_unsigned_int
boost::integer_range< IncrementableType > iota_view
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
const unsigned int n_inner_dofs
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
const unsigned int n_dofs
const unsigned int n_quadrature_points