16#ifndef dealii_fe_coupling_values_h
17#define dealii_fe_coupling_values_h
34template <
int dim,
int spacedim>
168 template <
typename ViewType>
186 template <
typename Number>
194 template <
typename Number>
230 value(
const unsigned int shape_function,
const unsigned int q_point)
const;
250 const unsigned int q_point)
const;
266 template <
typename Number>
291 template <
class InputVector>
294 const InputVector &dof_values,
313 template <
typename Number>
339 template <
class InputVector>
342 const InputVector &dof_values,
370 template <
typename Number>
373 const unsigned int size,
374 const Number &exemplar_number)
const;
380 template <
typename InputVector>
388 template <
typename ValueType>
389 std::vector<ValueType> &
397 template <
typename ValueType>
400 std::vector<ValueType> &outer_values)
const;
730template <
int dim1,
int dim2 = dim1,
int spacedim = dim1>
829 template <
typename Extractor>
837 template <
typename Extractor>
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,
952 std::pair<unsigned int, unsigned int>
960 std::pair<unsigned int, unsigned int>
974 template <
typename Extractor>
985 template <
typename Extractor>
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)
1093 for (
const auto i : dof_renumbering)
1097 for (
const auto q : quadrature_renumbering)
1104 template <
typename ViewType>
1106 const RenumberingData &data)
1113 template <
typename ViewType>
1114 template <
typename Number>
1116 RenumberedView<ViewType>::get_unique_container_name(
1117 const std::string &prefix,
1118 const unsigned int size,
1119 const Number &exemplar_number)
const
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
1135 const auto inner_shape_function = data.dof_renumbering.empty() ?
1137 data.dof_renumbering[shape_function];
1138 const auto inner_q_point = data.quadrature_renumbering.empty() ?
1140 data.quadrature_renumbering[q_point];
1142 return value_type(0);
1147 return view.value(inner_shape_function, inner_q_point);
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
1161 const auto inner_shape_function = data.dof_renumbering.empty() ?
1163 data.dof_renumbering[shape_function];
1164 const auto inner_q_point = data.quadrature_renumbering.empty() ?
1166 data.quadrature_renumbering[q_point];
1168 return gradient_type();
1170 return view.gradient(inner_shape_function, inner_q_point);
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
1182 if (data.quadrature_renumbering.empty())
1184 return outer_values;
1189 get_unique_container_name(
"RenumberedView::outer_to_inner_values",
1190 data.n_inner_quadrature_points,
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;
1202 template <
typename ViewType>
1203 template <
typename VectorType>
1205 RenumberedView<ViewType>::outer_to_inner_dofs(
1206 const VectorType &outer_dofs)
const
1209 if (data.dof_renumbering.empty())
1216 get_unique_container_name(
"RenumberedView::outer_to_inner_dofs",
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)
1225 const auto inner_i = data.dof_renumbering[i];
1229 inner_dofs[inner_i] = outer_dofs[i];
1238 template <
typename ViewType>
1239 template <
typename ValueType>
1241 RenumberedView<ViewType>::inner_to_outer_values(
1242 const std::vector<ValueType> &inner_values,
1243 std::vector<ValueType> &outer_values)
const
1247 if (data.quadrature_renumbering.empty())
1252 for (
unsigned int i = 0; i < data.quadrature_renumbering.size(); ++i)
1254 outer_values[i] = inner_values[data.quadrature_renumbering[i]];
1260 template <
typename ViewType>
1261 template <
typename Number>
1263 RenumberedView<ViewType>::get_function_values(
1265 std::vector<solution_value_type<Number>> &values)
const
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);
1274 template <
typename ViewType>
1275 template <
typename InputVector>
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)
1282 const auto &inner_dof_values = outer_to_inner_dofs(dof_values);
1283 auto &inner_values = outer_to_inner_values(values);
1285 view.get_function_values_from_local_dof_values(inner_dof_values,
1287 inner_to_outer_values(inner_values, values);
1291 template <
typename ViewType>
1292 template <
typename Number>
1294 RenumberedView<ViewType>::get_function_gradients(
1296 std::vector<solution_gradient_type<Number>> &gradients)
const
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);
1305 template <
typename ViewType>
1306 template <
typename InputVector>
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>>
1313 const auto &inner_dof_values = outer_to_inner_dofs(dof_values);
1314 auto &inner_gradients = outer_to_inner_values(gradients);
1316 view.get_function_gradients_from_local_dof_values(inner_dof_values,
1318 inner_to_outer_values(inner_gradients, gradients);
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,
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);
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;
1351template <
int dim1,
int dim2,
int spacedim>
1353 : first_fe_values(nullptr)
1354 , second_fe_values(nullptr)
1356 , n_quadrature_points_(0)
1357 , n_coupling_dofs_(
numbers::invalid_unsigned_int)
1362template <
int dim1,
int dim2,
int spacedim>
1375template <
int dim1,
int dim2,
int spacedim>
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;
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;
1392 unsigned int first_n_inner_quadrature_points =
1394 unsigned int second_n_inner_quadrature_points =
1398 unsigned int first_n_dofs = 0;
1399 unsigned int second_n_dofs = 0;
1401 unsigned int first_n_quadrature_points = 0;
1402 unsigned int second_n_quadrature_points = 0;
1404 std::vector<unsigned int> first_dofs_map;
1405 std::vector<unsigned int> second_dofs_map;
1407 std::vector<unsigned int> first_quad_map;
1408 std::vector<unsigned int> second_quad_map;
1412 fe_values_2.
get_cell()->diameter()) *
1416 switch (dof_coupling_type)
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;
1432 first_n_dofs = n_coupling_dofs_;
1433 second_n_dofs = n_coupling_dofs_;
1435 first_dofs_map.resize(n_coupling_dofs_);
1436 second_dofs_map.resize(n_coupling_dofs_);
1438 unsigned int idx = 0;
1439 for (
const unsigned int &i : fe_values_1.dof_indices())
1441 first_dofs_map[idx] = i;
1444 for (
const unsigned int &i : fe_values_2.dof_indices())
1447 second_dofs_map[idx++] = i;
1457 switch (quadrature_coupling_type)
1464 n_quadrature_points_ =
1465 quadrature_points_1.size() * quadrature_points_2.size();
1467 first_n_quadrature_points = n_quadrature_points_;
1468 second_n_quadrature_points = n_quadrature_points_;
1470 first_quad_map.resize(n_quadrature_points_);
1471 second_quad_map.resize(n_quadrature_points_);
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())
1477 first_quad_map[idx] = i;
1478 second_quad_map[idx] = j;
1487 ExcMessage(
"The two FEValuesBase objects must have the same "
1488 "number of quadrature points"));
1492 first_n_quadrature_points = n_quadrature_points_;
1493 second_n_quadrature_points = n_quadrature_points_;
1495 first_quad_map.clear();
1496 second_quad_map.clear();
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"));
1509 for (
const unsigned int &i : fe_values_1.quadrature_point_indices())
1511 Assert(quadrature_points_1[i].distance(quadrature_points_2[i]) <
1514 "The two FEValuesBase objects must have the same "
1515 "quadrature points"));
1520 first_n_quadrature_points = n_quadrature_points_;
1521 second_n_quadrature_points = n_quadrature_points_;
1523 first_quad_map.clear();
1524 second_quad_map.clear();
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"));
1539 first_n_quadrature_points = n_quadrature_points_;
1540 second_n_quadrature_points = n_quadrature_points_;
1543 first_quad_map.clear();
1547 for (
const unsigned int &i : fe_values_1.quadrature_point_indices())
1550 for (
const unsigned int &j :
1551 fe_values_2.quadrature_point_indices())
1552 if (quadrature_points_1[i].distance(quadrature_points_2[j]) <
1556 second_quad_map[i] = j;
1561 "The two FEValuesBase objects must have the same "
1562 "quadrature points, even if not in the same order."));
1572 for (
const unsigned int &i : fe_values_1.quadrature_point_indices())
1574 for (
const unsigned int &j :
1575 fe_values_2.quadrature_point_indices())
1576 if (quadrature_points_1[i].distance(quadrature_points_2[j]) <
1579 first_quad_map.emplace_back(i);
1580 second_quad_map.emplace_back(j);
1584 n_quadrature_points_ = first_quad_map.size();
1586 first_n_quadrature_points = n_quadrature_points_;
1587 second_n_quadrature_points = n_quadrature_points_;
1594 first_renumbering_data = std::make_unique<FEValuesViews::RenumberingData>(
1596 first_n_inner_quadrature_points,
1600 second_renumbering_data = std::make_unique<FEValuesViews::RenumberingData>(
1601 second_n_inner_dofs,
1602 second_n_inner_quadrature_points,
1608template <
int dim1,
int dim2,
int spacedim>
1614 const auto first_q = first_renumbering_data->quadrature_renumbering.empty() ?
1616 first_renumbering_data->quadrature_renumbering[q];
1621 const auto second_q =
1622 second_renumbering_data->quadrature_renumbering.empty() ?
1624 second_renumbering_data->quadrature_renumbering[q];
1625 return first_fe_values->JxW(first_q) * second_fe_values->JxW(second_q);
1628 return first_fe_values->JxW(first_q);
1633template <
int dim1,
int dim2,
int spacedim>
1637 return {0U, n_quadrature_points_};
1642template <
int dim1,
int dim2,
int spacedim>
1648 "Dofs are independent. You cannot ask for coupling dofs."));
1649 return {0U, n_coupling_dofs_};
1654template <
int dim1,
int dim2,
int spacedim>
1658 return {0U, n_first_dofs()};
1663template <
int dim1,
int dim2,
int spacedim>
1667 return {0U, n_second_dofs()};
1672template <
int dim1,
int dim2,
int spacedim>
1678 "Dofs are independent. You cannot ask for coupling dofs."));
1679 return n_coupling_dofs_;
1684template <
int dim1,
int dim2,
int spacedim>
1688 return first_renumbering_data->n_dofs;
1693template <
int dim1,
int dim2,
int spacedim>
1697 return second_renumbering_data->n_dofs;
1702template <
int dim1,
int dim2,
int spacedim>
1706 return n_quadrature_points_;
1711template <
int dim1,
int dim2,
int spacedim>
1714 const unsigned int quadrature_point)
const
1717 const auto first_q = first_fe_values->quadrature_point(
1718 first_renumbering_data->quadrature_renumbering.empty() ?
1720 first_renumbering_data->quadrature_renumbering[quadrature_point]);
1722 const auto second_q = second_fe_values->quadrature_point(
1723 second_renumbering_data->quadrature_renumbering.empty() ?
1725 second_renumbering_data->quadrature_renumbering[quadrature_point]);
1727 return {first_q, second_q};
1732template <
int dim1,
int dim2,
int spacedim>
1733std::pair<unsigned int, unsigned int>
1736 const unsigned int quadrature_point)
const
1739 const auto first_id =
1740 first_renumbering_data->quadrature_renumbering.empty() ?
1742 first_renumbering_data->quadrature_renumbering[quadrature_point];
1744 const auto second_id =
1745 second_renumbering_data->quadrature_renumbering.empty() ?
1747 second_renumbering_data->quadrature_renumbering[quadrature_point];
1748 return std::make_pair(first_id, second_id);
1753template <
int dim1,
int dim2,
int spacedim>
1754std::pair<unsigned int, unsigned int>
1756 const unsigned int coupling_dof_index)
const
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];
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);
1773template <
int dim1,
int dim2,
int spacedim>
1774template <
typename Extractor>
1777 const Extractor &extractor)
const
1784template <
int dim1,
int dim2,
int spacedim>
1785template <
typename Extractor>
1788 const Extractor &extractor)
const
1795template <
int dim1,
int dim2,
int spacedim>
1796template <
typename Extractor>
1804 (*first_fe_values)[extractor.
extractor], *first_renumbering_data);
1809template <
int dim1,
int dim2,
int spacedim>
1810template <
typename Extractor>
1818 (*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_
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
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)
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)
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