16 #ifndef dealii_fe_interface_values_h
17 #define dealii_fe_interface_values_h
29 template <
int dim,
int spacedim>
42 template <
int dim,
int spacedim = dim>
63 template <
int dim,
int spacedim = dim>
96 const unsigned int component);
119 const unsigned int interface_dof_index,
120 const unsigned int q_point)
const;
129 jump(
const unsigned int interface_dof_index,
130 const unsigned int q_point)
const;
139 average(
const unsigned int interface_dof_index,
140 const unsigned int q_point)
const;
150 const unsigned int q_point)
const;
160 const unsigned int q_point)
const;
171 const unsigned int q_point)
const;
181 const unsigned int q_point)
const;
191 const unsigned int q_point)
const;
205 template <
int dim,
int spacedim = dim>
242 const unsigned int first_vector_component);
265 const unsigned int interface_dof_index,
266 const unsigned int q_point)
const;
274 jump(
const unsigned int interface_dof_index,
275 const unsigned int q_point)
const;
283 average(
const unsigned int interface_dof_index,
284 const unsigned int q_point)
const;
293 const unsigned int q_point)
const;
302 const unsigned int q_point)
const;
313 const unsigned int q_point)
const;
323 const unsigned int q_point)
const;
333 const unsigned int q_point)
const;
369 template <
int dim,
int spacedim = dim>
439 template <
class CellIteratorType>
442 const unsigned int face_no,
443 const unsigned int sub_face_no,
445 const unsigned int face_no_neighbor,
446 const unsigned int sub_face_no_neighbor);
459 template <
class CellIteratorType>
461 reinit(
const CellIteratorType &cell,
const unsigned int face_no);
523 JxW(
const unsigned int quadrature_point)
const;
530 const std::vector<double> &
542 const std::vector<Tensor<1, spacedim>> &
550 const std::vector<Point<spacedim>> &
576 std::vector<types::global_dof_index>
590 std::array<unsigned int, 2>
602 normal(
const unsigned int q_point_index)
const;
636 const unsigned int interface_dof_index,
637 const unsigned int q_point,
638 const unsigned int component = 0)
const;
655 jump(
const unsigned int interface_dof_index,
656 const unsigned int q_point,
657 const unsigned int component = 0)
const;
669 average(
const unsigned int interface_dof_index,
670 const unsigned int q_point,
671 const unsigned int component = 0)
const;
684 const unsigned int q_point,
685 const unsigned int component = 0)
const;
699 const unsigned int q_point,
700 const unsigned int component = 0)
const;
713 const unsigned int q_point,
714 const unsigned int component = 0)
const;
728 const unsigned int q_point,
729 const unsigned int component = 0)
const;
742 const unsigned int q_point,
743 const unsigned int component = 0)
const;
778 std::vector<std::array<unsigned int, 2>>
dofmap;
828 template <
int dim,
int spacedim>
834 : n_quadrature_points(quadrature.size())
835 , internal_fe_face_values(mapping, fe, quadrature, update_flags)
836 , internal_fe_subface_values(mapping, fe, quadrature, update_flags)
837 , internal_fe_face_values_neighbor(mapping, fe, quadrature, update_flags)
838 , internal_fe_subface_values_neighbor(mapping, fe, quadrature, update_flags)
839 , fe_face_values(nullptr)
840 , fe_face_values_neighbor(nullptr)
845 template <
int dim,
int spacedim>
851 : n_quadrature_points(quadrature.max_n_quadrature_points())
852 , internal_fe_face_values(mapping, fe, quadrature, update_flags)
853 , internal_fe_subface_values(mapping, fe, quadrature, update_flags)
854 , internal_fe_face_values_neighbor(mapping, fe, quadrature[0], update_flags)
855 , internal_fe_subface_values_neighbor(mapping,
859 , fe_face_values(nullptr)
860 , fe_face_values_neighbor(nullptr)
865 template <
int dim,
int spacedim>
870 : n_quadrature_points(quadrature.size())
871 , internal_fe_face_values(
876 , internal_fe_subface_values(
881 , internal_fe_face_values_neighbor(
886 , internal_fe_subface_values_neighbor(
891 , fe_face_values(nullptr)
892 , fe_face_values_neighbor(nullptr)
897 template <
int dim,
int spacedim>
898 template <
class CellIteratorType>
901 const CellIteratorType & cell,
902 const unsigned int face_no,
903 const unsigned int sub_face_no,
905 const unsigned int face_no_neighbor,
906 const unsigned int sub_face_no_neighbor)
910 internal_fe_face_values.reinit(cell, face_no);
911 fe_face_values = &internal_fe_face_values;
915 internal_fe_subface_values.reinit(cell, face_no, sub_face_no);
916 fe_face_values = &internal_fe_subface_values;
920 internal_fe_face_values_neighbor.reinit(cell_neighbor, face_no_neighbor);
921 fe_face_values_neighbor = &internal_fe_face_values_neighbor;
925 internal_fe_subface_values_neighbor.reinit(cell_neighbor,
927 sub_face_no_neighbor);
928 fe_face_values_neighbor = &internal_fe_subface_values_neighbor;
932 fe_face_values_neighbor->n_quadrature_points);
934 const_cast<unsigned int &
>(this->n_quadrature_points) =
935 fe_face_values->n_quadrature_points;
940 std::vector<types::global_dof_index> v(
941 fe_face_values->get_fe().n_dofs_per_cell());
942 cell->get_active_or_mg_dof_indices(v);
943 std::vector<types::global_dof_index> v2(
944 fe_face_values_neighbor->get_fe().n_dofs_per_cell());
945 cell_neighbor->get_active_or_mg_dof_indices(v2);
949 std::map<types::global_dof_index, std::pair<unsigned int, unsigned int>>
951 std::pair<unsigned int, unsigned int> invalid_entry(
954 for (
unsigned int i = 0; i < v.size(); ++i)
957 auto result = tempmap.insert(std::make_pair(v[i], invalid_entry));
958 result.first->second.first = i;
961 for (
unsigned int i = 0; i < v2.size(); ++i)
964 auto result = tempmap.insert(std::make_pair(v2[i], invalid_entry));
965 result.first->second.second = i;
969 dofmap.resize(tempmap.size());
970 interface_dof_indices.resize(tempmap.size());
971 unsigned int idx = 0;
972 for (
auto &x : tempmap)
974 interface_dof_indices[idx] = x.first;
975 dofmap[idx] = {{x.second.first, x.second.second}};
983 template <
int dim,
int spacedim>
984 template <
class CellIteratorType>
987 const unsigned int face_no)
989 internal_fe_face_values.reinit(cell, face_no);
990 fe_face_values = &internal_fe_face_values;
991 fe_face_values_neighbor =
nullptr;
993 interface_dof_indices.resize(fe_face_values->get_fe().n_dofs_per_cell());
994 cell->get_active_or_mg_dof_indices(interface_dof_indices);
996 dofmap.resize(interface_dof_indices.size());
998 for (
unsigned int i = 0; i < interface_dof_indices.size(); ++i)
1006 template <
int dim,
int spacedim>
1010 Assert(fe_face_values !=
nullptr,
1011 ExcMessage(
"This call requires a call to reinit() first."));
1012 return fe_face_values->JxW(q);
1017 template <
int dim,
int spacedim>
1018 const std::vector<double> &
1021 Assert(fe_face_values !=
nullptr,
1022 ExcMessage(
"This call requires a call to reinit() first."));
1023 return fe_face_values->get_JxW_values();
1028 template <
int dim,
int spacedim>
1029 const std::vector<Tensor<1, spacedim>> &
1032 Assert(fe_face_values !=
nullptr,
1033 ExcMessage(
"This call requires a call to reinit() first."));
1034 return fe_face_values->get_normal_vectors();
1039 template <
int dim,
int spacedim>
1043 return internal_fe_face_values.get_mapping();
1048 template <
int dim,
int spacedim>
1052 return internal_fe_face_values.get_fe();
1057 template <
int dim,
int spacedim>
1061 return internal_fe_face_values.get_quadrature();
1066 template <
int dim,
int spacedim>
1067 const std::vector<Point<spacedim>> &
1070 Assert(fe_face_values !=
nullptr,
1071 ExcMessage(
"This call requires a call to reinit() first."));
1072 return fe_face_values->get_quadrature_points();
1077 template <
int dim,
int spacedim>
1081 return internal_fe_face_values.get_update_flags();
1086 template <
int dim,
int spacedim>
1091 interface_dof_indices.size() > 0,
1093 "n_current_interface_dofs() is only available after a call to reinit()."));
1094 return interface_dof_indices.size();
1099 template <
int dim,
int spacedim>
1103 return fe_face_values_neighbor ==
nullptr;
1108 template <
int dim,
int spacedim>
1109 std::vector<types::global_dof_index>
1112 return interface_dof_indices;
1117 template <
int dim,
int spacedim>
1118 std::array<unsigned int, 2>
1120 const unsigned int interface_dof_index)
const
1123 return dofmap[interface_dof_index];
1128 template <
int dim,
int spacedim>
1137 "You are on a boundary, so you can only ask for the first FEFaceValues object."));
1139 return (
cell_index == 0) ? *fe_face_values : *fe_face_values_neighbor;
1144 template <
int dim,
int spacedim>
1148 return fe_face_values->normal_vector(q_point_index);
1153 template <
int dim,
int spacedim>
1156 const bool here_or_there,
1157 const unsigned int interface_dof_index,
1158 const unsigned int q_point,
1159 const unsigned int component)
const
1161 const auto dof_pair = dofmap[interface_dof_index];
1164 return get_fe_face_values(0).shape_value_component(dof_pair[0],
1168 return get_fe_face_values(1).shape_value_component(dof_pair[1],
1177 template <
int dim,
int spacedim>
1180 const unsigned int q_point,
1181 const unsigned int component)
const
1183 const auto dof_pair = dofmap[interface_dof_index];
1188 value += get_fe_face_values(0).shape_value_component(dof_pair[0],
1192 value -= get_fe_face_values(1).shape_value_component(dof_pair[1],
1200 template <
int dim,
int spacedim>
1203 const unsigned int interface_dof_index,
1204 const unsigned int q_point,
1205 const unsigned int component)
const
1207 const auto dof_pair = dofmap[interface_dof_index];
1210 return get_fe_face_values(0).shape_value_component(dof_pair[0],
1217 value += 0.5 * get_fe_face_values(0).shape_value_component(dof_pair[0],
1221 value += 0.5 * get_fe_face_values(1).shape_value_component(dof_pair[1],
1230 template <
int dim,
int spacedim>
1233 const unsigned int interface_dof_index,
1234 const unsigned int q_point,
1235 const unsigned int component)
const
1237 const auto dof_pair = dofmap[interface_dof_index];
1240 return get_fe_face_values(0).shape_grad_component(dof_pair[0],
1247 value += 0.5 * get_fe_face_values(0).shape_grad_component(dof_pair[0],
1251 value += 0.5 * get_fe_face_values(1).shape_grad_component(dof_pair[1],
1260 template <
int dim,
int spacedim>
1263 const unsigned int interface_dof_index,
1264 const unsigned int q_point,
1265 const unsigned int component)
const
1267 const auto dof_pair = dofmap[interface_dof_index];
1270 return get_fe_face_values(0).shape_hessian_component(dof_pair[0],
1277 value += 0.5 * get_fe_face_values(0).shape_hessian_component(dof_pair[0],
1281 value += 0.5 * get_fe_face_values(1).shape_hessian_component(dof_pair[1],
1290 template <
int dim,
int spacedim>
1293 const unsigned int interface_dof_index,
1294 const unsigned int q_point,
1295 const unsigned int component)
const
1297 const auto dof_pair = dofmap[interface_dof_index];
1300 return get_fe_face_values(0).shape_grad_component(dof_pair[0],
1307 value += get_fe_face_values(0).shape_grad_component(dof_pair[0],
1311 value -= get_fe_face_values(1).shape_grad_component(dof_pair[1],
1320 template <
int dim,
int spacedim>
1323 const unsigned int interface_dof_index,
1324 const unsigned int q_point,
1325 const unsigned int component)
const
1327 const auto dof_pair = dofmap[interface_dof_index];
1330 return get_fe_face_values(0).shape_hessian_component(dof_pair[0],
1337 value += get_fe_face_values(0).shape_hessian_component(dof_pair[0],
1341 value -= get_fe_face_values(1).shape_hessian_component(dof_pair[1],
1350 template <
int dim,
int spacedim>
1353 const unsigned int interface_dof_index,
1354 const unsigned int q_point,
1355 const unsigned int component)
const
1357 const auto dof_pair = dofmap[interface_dof_index];
1360 return get_fe_face_values(0).shape_3rd_derivative_component(dof_pair[0],
1367 value += get_fe_face_values(0).shape_3rd_derivative_component(dof_pair[0],
1371 value -= get_fe_face_values(1).shape_3rd_derivative_component(dof_pair[1],
1381 template <
int dim,
int spacedim>
1392 template <
int dim,
int spacedim>
1398 const unsigned int n_vectors =
1412 template <
int dim,
int spacedim>
1415 : fe_interface(&fe_interface)
1420 template <
int dim,
int spacedim>
1423 const unsigned int component)
1424 : Base<dim, spacedim>(fe_interface)
1425 , extractor(component)
1430 template <
int dim,
int spacedim>
1433 const unsigned int interface_dof_index,
1434 const unsigned int q_point)
const
1436 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1439 return (*(this->fe_interface->fe_face_values))[extractor].value(
1440 dof_pair[0], q_point);
1443 return (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
1444 dof_pair[1], q_point);
1451 template <
int dim,
int spacedim>
1453 Scalar<dim, spacedim>::jump(
const unsigned int interface_dof_index,
1454 const unsigned int q_point)
const
1456 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1458 value_type
value = 0.0;
1462 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
1467 (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
1468 dof_pair[1], q_point);
1475 template <
int dim,
int spacedim>
1477 Scalar<dim, spacedim>::average(
const unsigned int interface_dof_index,
1478 const unsigned int q_point)
const
1480 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1482 if (this->fe_interface->at_boundary())
1483 return (*(this->fe_interface->fe_face_values))[extractor].value(
1484 dof_pair[0], q_point);
1486 value_type
value = 0.0;
1491 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
1496 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
1497 dof_pair[1], q_point);
1504 template <
int dim,
int spacedim>
1506 Scalar<dim, spacedim>::average_gradient(
1507 const unsigned int interface_dof_index,
1508 const unsigned int q_point)
const
1510 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1512 if (this->fe_interface->at_boundary())
1513 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
1514 dof_pair[0], q_point);
1516 gradient_type
value;
1521 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
1525 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
1526 .gradient(dof_pair[1], q_point);
1533 template <
int dim,
int spacedim>
1535 Scalar<dim, spacedim>::jump_gradient(
const unsigned int interface_dof_index,
1536 const unsigned int q_point)
const
1538 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1540 if (this->fe_interface->at_boundary())
1541 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
1542 dof_pair[0], q_point);
1544 gradient_type
value;
1548 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
1553 (*(this->fe_interface->fe_face_values_neighbor))[extractor].gradient(
1554 dof_pair[1], q_point);
1561 template <
int dim,
int spacedim>
1563 Scalar<dim, spacedim>::average_hessian(
const unsigned int interface_dof_index,
1564 const unsigned int q_point)
const
1566 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1568 if (this->fe_interface->at_boundary())
1569 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
1570 dof_pair[0], q_point);
1577 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
1581 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
1582 .hessian(dof_pair[1], q_point);
1589 template <
int dim,
int spacedim>
1591 Scalar<dim, spacedim>::jump_3rd_derivative(
1592 const unsigned int interface_dof_index,
1593 const unsigned int q_point)
const
1595 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1597 if (this->fe_interface->at_boundary())
1598 return (*(this->fe_interface->fe_face_values))[extractor]
1599 .third_derivative(dof_pair[0], q_point);
1601 third_derivative_type
value;
1605 (*(this->fe_interface->fe_face_values))[extractor].third_derivative(
1606 dof_pair[0], q_point);
1609 value -= (*(this->fe_interface->fe_face_values_neighbor))[extractor]
1610 .third_derivative(dof_pair[1], q_point);
1617 template <
int dim,
int spacedim>
1619 Scalar<dim, spacedim>::jump_hessian(
const unsigned int interface_dof_index,
1620 const unsigned int q_point)
const
1622 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1624 if (this->fe_interface->at_boundary())
1625 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
1626 dof_pair[0], q_point);
1632 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
1637 (*(this->fe_interface->fe_face_values_neighbor))[extractor].hessian(
1638 dof_pair[1], q_point);
1645 template <
int dim,
int spacedim>
1648 const unsigned int first_vector_component)
1649 : Base<dim, spacedim>(fe_interface)
1650 , extractor(first_vector_component)
1655 template <
int dim,
int spacedim>
1658 const unsigned int interface_dof_index,
1659 const unsigned int q_point)
const
1661 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1664 return (*(this->fe_interface->fe_face_values))[extractor].value(
1665 dof_pair[0], q_point);
1668 return (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
1669 dof_pair[1], q_point);
1671 return value_type();
1676 template <
int dim,
int spacedim>
1679 const unsigned int q_point)
const
1681 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1687 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
1692 (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
1693 dof_pair[1], q_point);
1700 template <
int dim,
int spacedim>
1703 const unsigned int q_point)
const
1705 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1707 if (this->fe_interface->at_boundary())
1708 return (*(this->fe_interface->fe_face_values))[extractor].value(
1709 dof_pair[0], q_point);
1716 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
1721 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
1722 dof_pair[1], q_point);
1729 template <
int dim,
int spacedim>
1732 const unsigned int interface_dof_index,
1733 const unsigned int q_point)
const
1735 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1737 if (this->fe_interface->at_boundary())
1738 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
1739 dof_pair[0], q_point);
1741 gradient_type
value;
1746 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
1750 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
1751 .gradient(dof_pair[1], q_point);
1758 template <
int dim,
int spacedim>
1761 const unsigned int q_point)
const
1763 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1765 if (this->fe_interface->at_boundary())
1766 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
1767 dof_pair[0], q_point);
1769 gradient_type
value;
1773 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
1778 (*(this->fe_interface->fe_face_values_neighbor))[extractor].gradient(
1779 dof_pair[1], q_point);
1786 template <
int dim,
int spacedim>
1789 const unsigned int q_point)
const
1791 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1793 if (this->fe_interface->at_boundary())
1794 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
1795 dof_pair[0], q_point);
1802 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
1806 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
1807 .hessian(dof_pair[1], q_point);
1814 template <
int dim,
int spacedim>
1817 const unsigned int q_point)
const
1819 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1821 if (this->fe_interface->at_boundary())
1822 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
1823 dof_pair[0], q_point);
1829 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
1834 (*(this->fe_interface->fe_face_values_neighbor))[extractor].hessian(
1835 dof_pair[1], q_point);
1842 template <
int dim,
int spacedim>
1845 const unsigned int interface_dof_index,
1846 const unsigned int q_point)
const
1848 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1850 if (this->fe_interface->at_boundary())
1851 return (*(this->fe_interface->fe_face_values))[extractor]
1852 .third_derivative(dof_pair[0], q_point);
1854 third_derivative_type
value;
1858 (*(this->fe_interface->fe_face_values))[extractor].third_derivative(
1859 dof_pair[0], q_point);
1862 value -= (*(this->fe_interface->fe_face_values_neighbor))[extractor]
1863 .third_derivative(dof_pair[1], q_point);
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors() const
std::vector< types::global_dof_index > get_interface_dof_indices() const
UpdateFlags get_update_flags() const
Tensor< 1, spacedim > normal(const unsigned int q_point_index) const
FEFaceValuesBase< dim, spacedim > * fe_face_values_neighbor
Tensor< 2, spacedim > average_hessian(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
const FEInterfaceViews::Scalar< dim, spacedim > operator[](const FEValuesExtractors::Scalar &scalar) const
FESubfaceValues< dim, spacedim > internal_fe_subface_values_neighbor
FEInterfaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature, const UpdateFlags update_flags)
const unsigned int n_quadrature_points
const Mapping< dim, spacedim > & get_mapping() const
Tensor< 1, spacedim > jump_gradient(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
std::array< unsigned int, 2 > interface_dof_to_dof_indices(const unsigned int interface_dof_index) const
unsigned n_current_interface_dofs() const
Tensor< 1, spacedim > average_gradient(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
double average(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
double shape_value(const bool here_or_there, const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
const FEFaceValuesBase< dim, spacedim > & get_fe_face_values(const unsigned int cell_index) const
const std::vector< Point< spacedim > > & get_quadrature_points() const
double jump(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
FEInterfaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const hp::QCollection< dim - 1 > &quadrature, const UpdateFlags update_flags)
FEFaceValuesBase< dim, spacedim > * fe_face_values
std::vector< types::global_dof_index > interface_dof_indices
const std::vector< double > & get_JxW_values() const
const FEInterfaceViews::Vector< dim, spacedim > operator[](const FEValuesExtractors::Vector &vector) const
std::vector< std::array< unsigned int, 2 > > dofmap
double JxW(const unsigned int quadrature_point) const
FESubfaceValues< dim, spacedim > internal_fe_subface_values
const Quadrature< dim - 1 > & get_quadrature() const
void reinit(const CellIteratorType &cell, const unsigned int face_no)
FEInterfaceValues(const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature, const UpdateFlags update_flags)
Tensor< 3, spacedim > jump_3rd_derivative(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
Tensor< 2, spacedim > jump_hessian(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
const FiniteElement< dim, spacedim > & get_fe() const
FEFaceValues< dim, spacedim > internal_fe_face_values
FEFaceValues< dim, spacedim > internal_fe_face_values_neighbor
void reinit(const CellIteratorType &cell, const unsigned int face_no, const unsigned int sub_face_no, const typename identity< CellIteratorType >::type &cell_neighbor, const unsigned int face_no_neighbor, const unsigned int sub_face_no_neighbor)
Base(const FEInterfaceValues< dim, spacedim > &fe_interface)
const FEInterfaceValues< dim, spacedim > * fe_interface
gradient_type average_gradient(const unsigned int interface_dof_index, const unsigned int q_point) const
hessian_type jump_hessian(const unsigned int interface_dof_index, const unsigned int q_point) const
gradient_type jump_gradient(const unsigned int interface_dof_index, const unsigned int q_point) const
const FEValuesExtractors::Scalar extractor
value_type jump(const unsigned int interface_dof_index, const unsigned int q_point) const
typename FEValuesViews::Scalar< dim, spacedim >::gradient_type gradient_type
value_type value(const bool here_or_there, const unsigned int interface_dof_index, const unsigned int q_point) const
value_type average(const unsigned int interface_dof_index, const unsigned int q_point) const
hessian_type average_hessian(const unsigned int interface_dof_index, const unsigned int q_point) const
Scalar(const FEInterfaceValues< dim, spacedim > &fe_interface, const unsigned int component)
typename FEValuesViews::Scalar< dim, spacedim >::hessian_type hessian_type
third_derivative_type jump_3rd_derivative(const unsigned int interface_dof_index, const unsigned int q_point) const
typename FEValuesViews::Scalar< dim, spacedim >::third_derivative_type third_derivative_type
const FEValuesExtractors::Vector extractor
hessian_type average_hessian(const unsigned int interface_dof_index, const unsigned int q_point) const
value_type average(const unsigned int interface_dof_index, const unsigned int q_point) const
Vector(const FEInterfaceValues< dim, spacedim > &fe_interface, const unsigned int first_vector_component)
hessian_type jump_hessian(const unsigned int interface_dof_index, const unsigned int q_point) const
typename FEValuesViews::Vector< dim, spacedim >::value_type value_type
value_type value(const bool here_or_there, const unsigned int interface_dof_index, const unsigned int q_point) const
third_derivative_type jump_3rd_derivative(const unsigned int interface_dof_index, const unsigned int q_point) const
value_type jump(const unsigned int interface_dof_index, const unsigned int q_point) const
gradient_type jump_gradient(const unsigned int interface_dof_index, const unsigned int q_point) const
typename FEValuesViews::Vector< dim, spacedim >::gradient_type gradient_type
typename FEValuesViews::Vector< dim, spacedim >::third_derivative_type third_derivative_type
gradient_type average_gradient(const unsigned int interface_dof_index, const unsigned int q_point) const
typename FEValuesViews::Vector< dim, spacedim >::hessian_type hessian_type
value_type value(const unsigned int shape_function, const unsigned int q_point) const
::Tensor< 1, spacedim > gradient_type
::Tensor< 2, spacedim > hessian_type
::Tensor< 3, spacedim > third_derivative_type
unsigned int n_components() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcMessage(std::string arg1)
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
static const unsigned int invalid_unsigned_int