Reference documentation for deal.II version 9.3.2
\(\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\}}\)
fe_interface_values.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2018 - 2021 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_interface_values_h
17 #define dealii_fe_interface_values_h
18 
19 #include <deal.II/base/config.h>
20 
21 #include <deal.II/fe/fe_values.h>
22 #include <deal.II/fe/mapping_q1.h>
23 
25 
27 
28 #ifndef DOXYGEN
29 template <int dim, int spacedim>
30 class FEInterfaceValues;
31 #endif
32 
38 {
42  template <int dim, int spacedim = dim>
43  class Base
44  {
45  public:
50 
51  protected:
56  };
57 
58 
59 
63  template <int dim, int spacedim = dim>
64  class Scalar : public Base<dim, spacedim>
65  {
66  public:
70  using value_type = double;
71 
76  using gradient_type =
78 
82  using hessian_type =
84 
91 
96  const unsigned int component);
97 
117  value_type
118  value(const bool here_or_there,
119  const unsigned int interface_dof_index,
120  const unsigned int q_point) const;
121 
128  value_type
129  jump(const unsigned int interface_dof_index,
130  const unsigned int q_point) const;
131 
138  value_type
139  average(const unsigned int interface_dof_index,
140  const unsigned int q_point) const;
141 
149  average_gradient(const unsigned int interface_dof_index,
150  const unsigned int q_point) const;
151 
159  jump_gradient(const unsigned int interface_dof_index,
160  const unsigned int q_point) const;
161 
170  average_hessian(const unsigned int interface_dof_index,
171  const unsigned int q_point) const;
172 
180  jump_hessian(const unsigned int interface_dof_index,
181  const unsigned int q_point) const;
182 
190  jump_3rd_derivative(const unsigned int interface_dof_index,
191  const unsigned int q_point) const;
192 
193  private:
198  };
199 
200 
201 
205  template <int dim, int spacedim = dim>
206  class Vector : public Base<dim, spacedim>
207  {
208  public:
212  using value_type =
214 
221 
227  using hessian_type =
229 
237 
242  const unsigned int first_vector_component);
243 
263  value_type
264  value(const bool here_or_there,
265  const unsigned int interface_dof_index,
266  const unsigned int q_point) const;
267 
273  value_type
274  jump(const unsigned int interface_dof_index,
275  const unsigned int q_point) const;
276 
282  value_type
283  average(const unsigned int interface_dof_index,
284  const unsigned int q_point) const;
285 
292  average_gradient(const unsigned int interface_dof_index,
293  const unsigned int q_point) const;
294 
301  jump_gradient(const unsigned int interface_dof_index,
302  const unsigned int q_point) const;
303 
312  average_hessian(const unsigned int interface_dof_index,
313  const unsigned int q_point) const;
314 
322  jump_hessian(const unsigned int interface_dof_index,
323  const unsigned int q_point) const;
324 
332  jump_3rd_derivative(const unsigned int interface_dof_index,
333  const unsigned int q_point) const;
334 
335  private:
340  };
341 } // namespace FEInterfaceViews
342 
343 
344 
369 template <int dim, int spacedim = dim>
371 {
372 public:
376  const unsigned int n_quadrature_points;
377 
385  const Quadrature<dim - 1> & quadrature,
386  const UpdateFlags update_flags);
387 
395  const hp::QCollection<dim - 1> & quadrature,
396  const UpdateFlags update_flags);
397 
405  const Quadrature<dim - 1> & quadrature,
406  const UpdateFlags update_flags);
407 
439  template <class CellIteratorType>
440  void
441  reinit(const CellIteratorType & cell,
442  const unsigned int face_no,
443  const unsigned int sub_face_no,
444  const typename identity<CellIteratorType>::type &cell_neighbor,
445  const unsigned int face_no_neighbor,
446  const unsigned int sub_face_no_neighbor);
447 
459  template <class CellIteratorType>
460  void
461  reinit(const CellIteratorType &cell, const unsigned int face_no);
462 
471  get_fe_face_values(const unsigned int cell_index) const;
472 
476  const Mapping<dim, spacedim> &
477  get_mapping() const;
478 
483  get_fe() const;
484 
488  const Quadrature<dim - 1> &
489  get_quadrature() const;
490 
496 
508  bool
509  at_boundary() const;
510 
522  double
523  JxW(const unsigned int quadrature_point) const;
524 
530  const std::vector<double> &
531  get_JxW_values() const;
532 
542  const std::vector<Tensor<1, spacedim>> &
544 
550  const std::vector<Point<spacedim>> &
552 
563  unsigned
565 
576  std::vector<types::global_dof_index>
578 
590  std::array<unsigned int, 2>
591  interface_dof_to_dof_indices(const unsigned int interface_dof_index) const;
592 
602  normal(const unsigned int q_point_index) const;
603 
634  double
635  shape_value(const bool here_or_there,
636  const unsigned int interface_dof_index,
637  const unsigned int q_point,
638  const unsigned int component = 0) const;
639 
654  double
655  jump(const unsigned int interface_dof_index,
656  const unsigned int q_point,
657  const unsigned int component = 0) const;
658 
668  double
669  average(const unsigned int interface_dof_index,
670  const unsigned int q_point,
671  const unsigned int component = 0) const;
672 
683  average_gradient(const unsigned int interface_dof_index,
684  const unsigned int q_point,
685  const unsigned int component = 0) const;
686 
698  average_hessian(const unsigned int interface_dof_index,
699  const unsigned int q_point,
700  const unsigned int component = 0) const;
701 
712  jump_gradient(const unsigned int interface_dof_index,
713  const unsigned int q_point,
714  const unsigned int component = 0) const;
715 
727  jump_hessian(const unsigned int interface_dof_index,
728  const unsigned int q_point,
729  const unsigned int component = 0) const;
730 
741  jump_3rd_derivative(const unsigned int interface_dof_index,
742  const unsigned int q_point,
743  const unsigned int component = 0) const;
744 
753 
762 
767 private:
771  std::vector<types::global_dof_index> interface_dof_indices;
772 
778  std::vector<std::array<unsigned int, 2>> dofmap;
779 
784 
789 
794 
799 
805 
812 
813  /* Make the view classes friends of this class, since they access internal
814  * data.
815  */
816  template <int, int>
818  template <int, int>
820 };
821 
822 
823 
824 #ifndef DOXYGEN
825 
826 /*---------------------- Inline functions ---------------------*/
827 
828 template <int dim, int spacedim>
830  const Mapping<dim, spacedim> & mapping,
832  const Quadrature<dim - 1> & quadrature,
833  const UpdateFlags update_flags)
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)
841 {}
842 
843 
844 
845 template <int dim, int spacedim>
847  const Mapping<dim, spacedim> & mapping,
849  const hp::QCollection<dim - 1> & quadrature,
850  const UpdateFlags update_flags)
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,
856  fe,
857  quadrature[0],
858  update_flags)
859  , fe_face_values(nullptr)
860  , fe_face_values_neighbor(nullptr)
861 {}
862 
863 
864 
865 template <int dim, int spacedim>
868  const Quadrature<dim - 1> & quadrature,
869  const UpdateFlags update_flags)
870  : n_quadrature_points(quadrature.size())
871  , internal_fe_face_values(
872  fe.reference_cell().template get_default_linear_mapping<dim, spacedim>(),
873  fe,
874  quadrature,
875  update_flags)
876  , internal_fe_subface_values(
877  fe.reference_cell().template get_default_linear_mapping<dim, spacedim>(),
878  fe,
879  quadrature,
880  update_flags)
881  , internal_fe_face_values_neighbor(
882  fe.reference_cell().template get_default_linear_mapping<dim, spacedim>(),
883  fe,
884  quadrature,
885  update_flags)
886  , internal_fe_subface_values_neighbor(
887  fe.reference_cell().template get_default_linear_mapping<dim, spacedim>(),
888  fe,
889  quadrature,
890  update_flags)
891  , fe_face_values(nullptr)
892  , fe_face_values_neighbor(nullptr)
893 {}
894 
895 
896 
897 template <int dim, int spacedim>
898 template <class CellIteratorType>
899 void
901  const CellIteratorType & cell,
902  const unsigned int face_no,
903  const unsigned int sub_face_no,
904  const typename identity<CellIteratorType>::type &cell_neighbor,
905  const unsigned int face_no_neighbor,
906  const unsigned int sub_face_no_neighbor)
907 {
908  if (sub_face_no == numbers::invalid_unsigned_int)
909  {
910  internal_fe_face_values.reinit(cell, face_no);
911  fe_face_values = &internal_fe_face_values;
912  }
913  else
914  {
915  internal_fe_subface_values.reinit(cell, face_no, sub_face_no);
916  fe_face_values = &internal_fe_subface_values;
917  }
918  if (sub_face_no_neighbor == numbers::invalid_unsigned_int)
919  {
920  internal_fe_face_values_neighbor.reinit(cell_neighbor, face_no_neighbor);
921  fe_face_values_neighbor = &internal_fe_face_values_neighbor;
922  }
923  else
924  {
925  internal_fe_subface_values_neighbor.reinit(cell_neighbor,
926  face_no_neighbor,
927  sub_face_no_neighbor);
928  fe_face_values_neighbor = &internal_fe_subface_values_neighbor;
929  }
930 
931  AssertDimension(fe_face_values->n_quadrature_points,
932  fe_face_values_neighbor->n_quadrature_points);
933 
934  const_cast<unsigned int &>(this->n_quadrature_points) =
935  fe_face_values->n_quadrature_points;
936 
937  // Set up dof mapping and remove duplicates (for continuous elements).
938  {
939  // Get dof indices first:
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);
946 
947  // Fill a map from the global dof index to the left and right
948  // local index.
949  std::map<types::global_dof_index, std::pair<unsigned int, unsigned int>>
950  tempmap;
951  std::pair<unsigned int, unsigned int> invalid_entry(
953 
954  for (unsigned int i = 0; i < v.size(); ++i)
955  {
956  // If not already existing, add an invalid entry:
957  auto result = tempmap.insert(std::make_pair(v[i], invalid_entry));
958  result.first->second.first = i;
959  }
960 
961  for (unsigned int i = 0; i < v2.size(); ++i)
962  {
963  // If not already existing, add an invalid entry:
964  auto result = tempmap.insert(std::make_pair(v2[i], invalid_entry));
965  result.first->second.second = i;
966  }
967 
968  // Transfer from the map to the sorted std::vectors.
969  dofmap.resize(tempmap.size());
970  interface_dof_indices.resize(tempmap.size());
971  unsigned int idx = 0;
972  for (auto &x : tempmap)
973  {
974  interface_dof_indices[idx] = x.first;
975  dofmap[idx] = {{x.second.first, x.second.second}};
976  ++idx;
977  }
978  }
979 }
980 
981 
982 
983 template <int dim, int spacedim>
984 template <class CellIteratorType>
985 void
986 FEInterfaceValues<dim, spacedim>::reinit(const CellIteratorType &cell,
987  const unsigned int face_no)
988 {
989  internal_fe_face_values.reinit(cell, face_no);
990  fe_face_values = &internal_fe_face_values;
991  fe_face_values_neighbor = nullptr;
992 
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);
995 
996  dofmap.resize(interface_dof_indices.size());
997 
998  for (unsigned int i = 0; i < interface_dof_indices.size(); ++i)
999  {
1000  dofmap[i] = {{i, numbers::invalid_unsigned_int}};
1001  }
1002 }
1003 
1004 
1005 
1006 template <int dim, int spacedim>
1007 inline double
1008 FEInterfaceValues<dim, spacedim>::JxW(const unsigned int q) const
1009 {
1010  Assert(fe_face_values != nullptr,
1011  ExcMessage("This call requires a call to reinit() first."));
1012  return fe_face_values->JxW(q);
1013 }
1014 
1015 
1016 
1017 template <int dim, int spacedim>
1018 const std::vector<double> &
1020 {
1021  Assert(fe_face_values != nullptr,
1022  ExcMessage("This call requires a call to reinit() first."));
1023  return fe_face_values->get_JxW_values();
1024 }
1025 
1026 
1027 
1028 template <int dim, int spacedim>
1029 const std::vector<Tensor<1, spacedim>> &
1031 {
1032  Assert(fe_face_values != nullptr,
1033  ExcMessage("This call requires a call to reinit() first."));
1034  return fe_face_values->get_normal_vectors();
1035 }
1036 
1037 
1038 
1039 template <int dim, int spacedim>
1040 const Mapping<dim, spacedim> &
1042 {
1043  return internal_fe_face_values.get_mapping();
1044 }
1045 
1046 
1047 
1048 template <int dim, int spacedim>
1051 {
1052  return internal_fe_face_values.get_fe();
1053 }
1054 
1055 
1056 
1057 template <int dim, int spacedim>
1058 const Quadrature<dim - 1> &
1060 {
1061  return internal_fe_face_values.get_quadrature();
1062 }
1063 
1064 
1065 
1066 template <int dim, int spacedim>
1067 const std::vector<Point<spacedim>> &
1069 {
1070  Assert(fe_face_values != nullptr,
1071  ExcMessage("This call requires a call to reinit() first."));
1072  return fe_face_values->get_quadrature_points();
1073 }
1074 
1075 
1076 
1077 template <int dim, int spacedim>
1080 {
1081  return internal_fe_face_values.get_update_flags();
1082 }
1083 
1084 
1085 
1086 template <int dim, int spacedim>
1087 unsigned
1089 {
1090  Assert(
1091  interface_dof_indices.size() > 0,
1092  ExcMessage(
1093  "n_current_interface_dofs() is only available after a call to reinit()."));
1094  return interface_dof_indices.size();
1095 }
1096 
1097 
1098 
1099 template <int dim, int spacedim>
1100 bool
1102 {
1103  return fe_face_values_neighbor == nullptr;
1104 }
1105 
1106 
1107 
1108 template <int dim, int spacedim>
1109 std::vector<types::global_dof_index>
1111 {
1112  return interface_dof_indices;
1113 }
1114 
1115 
1116 
1117 template <int dim, int spacedim>
1118 std::array<unsigned int, 2>
1120  const unsigned int interface_dof_index) const
1121 {
1122  AssertIndexRange(interface_dof_index, n_current_interface_dofs());
1123  return dofmap[interface_dof_index];
1124 }
1125 
1126 
1127 
1128 template <int dim, int spacedim>
1131  const unsigned int cell_index) const
1132 {
1134  Assert(
1135  cell_index == 0 || !at_boundary(),
1136  ExcMessage(
1137  "You are on a boundary, so you can only ask for the first FEFaceValues object."));
1138 
1139  return (cell_index == 0) ? *fe_face_values : *fe_face_values_neighbor;
1140 }
1141 
1142 
1143 
1144 template <int dim, int spacedim>
1146 FEInterfaceValues<dim, spacedim>::normal(const unsigned int q_point_index) const
1147 {
1148  return fe_face_values->normal_vector(q_point_index);
1149 }
1150 
1151 
1152 
1153 template <int dim, int spacedim>
1154 double
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
1160 {
1161  const auto dof_pair = dofmap[interface_dof_index];
1162 
1163  if (here_or_there && dof_pair[0] != numbers::invalid_unsigned_int)
1164  return get_fe_face_values(0).shape_value_component(dof_pair[0],
1165  q_point,
1166  component);
1167  if (!here_or_there && dof_pair[1] != numbers::invalid_unsigned_int)
1168  return get_fe_face_values(1).shape_value_component(dof_pair[1],
1169  q_point,
1170  component);
1171 
1172  return 0.0;
1173 }
1174 
1175 
1176 
1177 template <int dim, int spacedim>
1178 double
1179 FEInterfaceValues<dim, spacedim>::jump(const unsigned int interface_dof_index,
1180  const unsigned int q_point,
1181  const unsigned int component) const
1182 {
1183  const auto dof_pair = dofmap[interface_dof_index];
1184 
1185  double value = 0.0;
1186 
1187  if (dof_pair[0] != numbers::invalid_unsigned_int)
1188  value += get_fe_face_values(0).shape_value_component(dof_pair[0],
1189  q_point,
1190  component);
1191  if (dof_pair[1] != numbers::invalid_unsigned_int)
1192  value -= get_fe_face_values(1).shape_value_component(dof_pair[1],
1193  q_point,
1194  component);
1195  return value;
1196 }
1197 
1198 
1199 
1200 template <int dim, int spacedim>
1201 double
1203  const unsigned int interface_dof_index,
1204  const unsigned int q_point,
1205  const unsigned int component) const
1206 {
1207  const auto dof_pair = dofmap[interface_dof_index];
1208 
1209  if (at_boundary())
1210  return get_fe_face_values(0).shape_value_component(dof_pair[0],
1211  q_point,
1212  component);
1213 
1214  double value = 0.0;
1215 
1216  if (dof_pair[0] != numbers::invalid_unsigned_int)
1217  value += 0.5 * get_fe_face_values(0).shape_value_component(dof_pair[0],
1218  q_point,
1219  component);
1220  if (dof_pair[1] != numbers::invalid_unsigned_int)
1221  value += 0.5 * get_fe_face_values(1).shape_value_component(dof_pair[1],
1222  q_point,
1223  component);
1224 
1225  return value;
1226 }
1227 
1228 
1229 
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
1236 {
1237  const auto dof_pair = dofmap[interface_dof_index];
1238 
1239  if (at_boundary())
1240  return get_fe_face_values(0).shape_grad_component(dof_pair[0],
1241  q_point,
1242  component);
1243 
1245 
1246  if (dof_pair[0] != numbers::invalid_unsigned_int)
1247  value += 0.5 * get_fe_face_values(0).shape_grad_component(dof_pair[0],
1248  q_point,
1249  component);
1250  if (dof_pair[1] != numbers::invalid_unsigned_int)
1251  value += 0.5 * get_fe_face_values(1).shape_grad_component(dof_pair[1],
1252  q_point,
1253  component);
1254 
1255  return value;
1256 }
1257 
1258 
1259 
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
1266 {
1267  const auto dof_pair = dofmap[interface_dof_index];
1268 
1269  if (at_boundary())
1270  return get_fe_face_values(0).shape_hessian_component(dof_pair[0],
1271  q_point,
1272  component);
1273 
1275 
1276  if (dof_pair[0] != numbers::invalid_unsigned_int)
1277  value += 0.5 * get_fe_face_values(0).shape_hessian_component(dof_pair[0],
1278  q_point,
1279  component);
1280  if (dof_pair[1] != numbers::invalid_unsigned_int)
1281  value += 0.5 * get_fe_face_values(1).shape_hessian_component(dof_pair[1],
1282  q_point,
1283  component);
1284 
1285  return value;
1286 }
1287 
1288 
1289 
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
1296 {
1297  const auto dof_pair = dofmap[interface_dof_index];
1298 
1299  if (at_boundary())
1300  return get_fe_face_values(0).shape_grad_component(dof_pair[0],
1301  q_point,
1302  component);
1303 
1305 
1306  if (dof_pair[0] != numbers::invalid_unsigned_int)
1307  value += get_fe_face_values(0).shape_grad_component(dof_pair[0],
1308  q_point,
1309  component);
1310  if (dof_pair[1] != numbers::invalid_unsigned_int)
1311  value -= get_fe_face_values(1).shape_grad_component(dof_pair[1],
1312  q_point,
1313  component);
1314 
1315  return value;
1316 }
1317 
1318 
1319 
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
1326 {
1327  const auto dof_pair = dofmap[interface_dof_index];
1328 
1329  if (at_boundary())
1330  return get_fe_face_values(0).shape_hessian_component(dof_pair[0],
1331  q_point,
1332  component);
1333 
1335 
1336  if (dof_pair[0] != numbers::invalid_unsigned_int)
1337  value += get_fe_face_values(0).shape_hessian_component(dof_pair[0],
1338  q_point,
1339  component);
1340  if (dof_pair[1] != numbers::invalid_unsigned_int)
1341  value -= get_fe_face_values(1).shape_hessian_component(dof_pair[1],
1342  q_point,
1343  component);
1344 
1345  return value;
1346 }
1347 
1348 
1349 
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
1356 {
1357  const auto dof_pair = dofmap[interface_dof_index];
1358 
1359  if (at_boundary())
1360  return get_fe_face_values(0).shape_3rd_derivative_component(dof_pair[0],
1361  q_point,
1362  component);
1363 
1365 
1366  if (dof_pair[0] != numbers::invalid_unsigned_int)
1367  value += get_fe_face_values(0).shape_3rd_derivative_component(dof_pair[0],
1368  q_point,
1369  component);
1370  if (dof_pair[1] != numbers::invalid_unsigned_int)
1371  value -= get_fe_face_values(1).shape_3rd_derivative_component(dof_pair[1],
1372  q_point,
1373  component);
1374 
1375  return value;
1376 }
1377 
1378 
1379 
1380 /*------------ Inline functions: FEInterfaceValues------------*/
1381 template <int dim, int spacedim>
1384  operator[](const FEValuesExtractors::Scalar &scalar) const
1385 {
1386  AssertIndexRange(scalar.component, this->get_fe().n_components());
1387  return FEInterfaceViews::Scalar<dim, spacedim>(*this, scalar.component);
1388 }
1389 
1390 
1391 
1392 template <int dim, int spacedim>
1395  operator[](const FEValuesExtractors::Vector &vector) const
1396 {
1397  const FiniteElement<dim, spacedim> &fe = this->get_fe();
1398  const unsigned int n_vectors =
1401  0);
1402  (void)n_vectors;
1403  AssertIndexRange(vector.first_vector_component, n_vectors);
1405  vector.first_vector_component);
1406 }
1407 
1408 
1409 
1410 namespace FEInterfaceViews
1411 {
1412  template <int dim, int spacedim>
1414  const FEInterfaceValues<dim, spacedim> &fe_interface)
1415  : fe_interface(&fe_interface)
1416  {}
1417 
1418 
1419 
1420  template <int dim, int spacedim>
1422  const FEInterfaceValues<dim, spacedim> &fe_interface,
1423  const unsigned int component)
1424  : Base<dim, spacedim>(fe_interface)
1425  , extractor(component)
1426  {}
1427 
1428 
1429 
1430  template <int dim, int spacedim>
1432  Scalar<dim, spacedim>::value(const bool here_or_there,
1433  const unsigned int interface_dof_index,
1434  const unsigned int q_point) const
1435  {
1436  const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1437 
1438  if (here_or_there && dof_pair[0] != numbers::invalid_unsigned_int)
1439  return (*(this->fe_interface->fe_face_values))[extractor].value(
1440  dof_pair[0], q_point);
1441 
1442  if (!here_or_there && dof_pair[1] != numbers::invalid_unsigned_int)
1443  return (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
1444  dof_pair[1], q_point);
1445 
1446  return 0.0;
1447  }
1448 
1449 
1450 
1451  template <int dim, int spacedim>
1453  Scalar<dim, spacedim>::jump(const unsigned int interface_dof_index,
1454  const unsigned int q_point) const
1455  {
1456  const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1457 
1458  value_type value = 0.0;
1459 
1460  if (dof_pair[0] != numbers::invalid_unsigned_int)
1461  value +=
1462  (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
1463  q_point);
1464 
1465  if (dof_pair[1] != numbers::invalid_unsigned_int)
1466  value -=
1467  (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
1468  dof_pair[1], q_point);
1469 
1470  return value;
1471  }
1472 
1473 
1474 
1475  template <int dim, int spacedim>
1477  Scalar<dim, spacedim>::average(const unsigned int interface_dof_index,
1478  const unsigned int q_point) const
1479  {
1480  const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1481 
1482  if (this->fe_interface->at_boundary())
1483  return (*(this->fe_interface->fe_face_values))[extractor].value(
1484  dof_pair[0], q_point);
1485 
1486  value_type value = 0.0;
1487 
1488  if (dof_pair[0] != numbers::invalid_unsigned_int)
1489  value +=
1490  0.5 *
1491  (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
1492  q_point);
1493 
1494  if (dof_pair[1] != numbers::invalid_unsigned_int)
1495  value +=
1496  0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
1497  dof_pair[1], q_point);
1498 
1499  return value;
1500  }
1501 
1502 
1503 
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
1509  {
1510  const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1511 
1512  if (this->fe_interface->at_boundary())
1513  return (*(this->fe_interface->fe_face_values))[extractor].gradient(
1514  dof_pair[0], q_point);
1515 
1516  gradient_type value;
1517 
1518  if (dof_pair[0] != numbers::invalid_unsigned_int)
1519  value +=
1520  0.5 *
1521  (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
1522  q_point);
1523 
1524  if (dof_pair[1] != numbers::invalid_unsigned_int)
1525  value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
1526  .gradient(dof_pair[1], q_point);
1527 
1528  return value;
1529  }
1530 
1531 
1532 
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
1537  {
1538  const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1539 
1540  if (this->fe_interface->at_boundary())
1541  return (*(this->fe_interface->fe_face_values))[extractor].gradient(
1542  dof_pair[0], q_point);
1543 
1544  gradient_type value;
1545 
1546  if (dof_pair[0] != numbers::invalid_unsigned_int)
1547  value +=
1548  (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
1549  q_point);
1550 
1551  if (dof_pair[1] != numbers::invalid_unsigned_int)
1552  value -=
1553  (*(this->fe_interface->fe_face_values_neighbor))[extractor].gradient(
1554  dof_pair[1], q_point);
1555 
1556  return value;
1557  }
1558 
1559 
1560 
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
1565  {
1566  const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1567 
1568  if (this->fe_interface->at_boundary())
1569  return (*(this->fe_interface->fe_face_values))[extractor].hessian(
1570  dof_pair[0], q_point);
1571 
1572  hessian_type value;
1573 
1574  if (dof_pair[0] != numbers::invalid_unsigned_int)
1575  value +=
1576  0.5 *
1577  (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
1578  q_point);
1579 
1580  if (dof_pair[1] != numbers::invalid_unsigned_int)
1581  value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
1582  .hessian(dof_pair[1], q_point);
1583 
1584  return value;
1585  }
1586 
1587 
1588 
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
1594  {
1595  const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1596 
1597  if (this->fe_interface->at_boundary())
1598  return (*(this->fe_interface->fe_face_values))[extractor]
1599  .third_derivative(dof_pair[0], q_point);
1600 
1601  third_derivative_type value;
1602 
1603  if (dof_pair[0] != numbers::invalid_unsigned_int)
1604  value +=
1605  (*(this->fe_interface->fe_face_values))[extractor].third_derivative(
1606  dof_pair[0], q_point);
1607 
1608  if (dof_pair[1] != numbers::invalid_unsigned_int)
1609  value -= (*(this->fe_interface->fe_face_values_neighbor))[extractor]
1610  .third_derivative(dof_pair[1], q_point);
1611 
1612  return value;
1613  }
1614 
1615 
1616 
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
1621  {
1622  const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1623 
1624  if (this->fe_interface->at_boundary())
1625  return (*(this->fe_interface->fe_face_values))[extractor].hessian(
1626  dof_pair[0], q_point);
1627 
1628  hessian_type value;
1629 
1630  if (dof_pair[0] != numbers::invalid_unsigned_int)
1631  value +=
1632  (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
1633  q_point);
1634 
1635  if (dof_pair[1] != numbers::invalid_unsigned_int)
1636  value -=
1637  (*(this->fe_interface->fe_face_values_neighbor))[extractor].hessian(
1638  dof_pair[1], q_point);
1639 
1640  return value;
1641  }
1642 
1643 
1644 
1645  template <int dim, int spacedim>
1647  const FEInterfaceValues<dim, spacedim> &fe_interface,
1648  const unsigned int first_vector_component)
1649  : Base<dim, spacedim>(fe_interface)
1650  , extractor(first_vector_component)
1651  {}
1652 
1653 
1654 
1655  template <int dim, int spacedim>
1657  Vector<dim, spacedim>::value(const bool here_or_there,
1658  const unsigned int interface_dof_index,
1659  const unsigned int q_point) const
1660  {
1661  const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1662 
1663  if (here_or_there && dof_pair[0] != numbers::invalid_unsigned_int)
1664  return (*(this->fe_interface->fe_face_values))[extractor].value(
1665  dof_pair[0], q_point);
1666 
1667  if (!here_or_there && dof_pair[1] != numbers::invalid_unsigned_int)
1668  return (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
1669  dof_pair[1], q_point);
1670 
1671  return value_type();
1672  }
1673 
1674 
1675 
1676  template <int dim, int spacedim>
1678  Vector<dim, spacedim>::jump(const unsigned int interface_dof_index,
1679  const unsigned int q_point) const
1680  {
1681  const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1682 
1683  value_type value;
1684 
1685  if (dof_pair[0] != numbers::invalid_unsigned_int)
1686  value +=
1687  (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
1688  q_point);
1689 
1690  if (dof_pair[1] != numbers::invalid_unsigned_int)
1691  value -=
1692  (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
1693  dof_pair[1], q_point);
1694 
1695  return value;
1696  }
1697 
1698 
1699 
1700  template <int dim, int spacedim>
1702  Vector<dim, spacedim>::average(const unsigned int interface_dof_index,
1703  const unsigned int q_point) const
1704  {
1705  const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1706 
1707  if (this->fe_interface->at_boundary())
1708  return (*(this->fe_interface->fe_face_values))[extractor].value(
1709  dof_pair[0], q_point);
1710 
1711  value_type value;
1712 
1713  if (dof_pair[0] != numbers::invalid_unsigned_int)
1714  value +=
1715  0.5 *
1716  (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
1717  q_point);
1718 
1719  if (dof_pair[1] != numbers::invalid_unsigned_int)
1720  value +=
1721  0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
1722  dof_pair[1], q_point);
1723 
1724  return value;
1725  }
1726 
1727 
1728 
1729  template <int dim, int spacedim>
1732  const unsigned int interface_dof_index,
1733  const unsigned int q_point) const
1734  {
1735  const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1736 
1737  if (this->fe_interface->at_boundary())
1738  return (*(this->fe_interface->fe_face_values))[extractor].gradient(
1739  dof_pair[0], q_point);
1740 
1741  gradient_type value;
1742 
1743  if (dof_pair[0] != numbers::invalid_unsigned_int)
1744  value +=
1745  0.5 *
1746  (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
1747  q_point);
1748 
1749  if (dof_pair[1] != numbers::invalid_unsigned_int)
1750  value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
1751  .gradient(dof_pair[1], q_point);
1752 
1753  return value;
1754  }
1755 
1756 
1757 
1758  template <int dim, int spacedim>
1760  Vector<dim, spacedim>::jump_gradient(const unsigned int interface_dof_index,
1761  const unsigned int q_point) const
1762  {
1763  const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1764 
1765  if (this->fe_interface->at_boundary())
1766  return (*(this->fe_interface->fe_face_values))[extractor].gradient(
1767  dof_pair[0], q_point);
1768 
1769  gradient_type value;
1770 
1771  if (dof_pair[0] != numbers::invalid_unsigned_int)
1772  value +=
1773  (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
1774  q_point);
1775 
1776  if (dof_pair[1] != numbers::invalid_unsigned_int)
1777  value -=
1778  (*(this->fe_interface->fe_face_values_neighbor))[extractor].gradient(
1779  dof_pair[1], q_point);
1780 
1781  return value;
1782  }
1783 
1784 
1785 
1786  template <int dim, int spacedim>
1788  Vector<dim, spacedim>::average_hessian(const unsigned int interface_dof_index,
1789  const unsigned int q_point) const
1790  {
1791  const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1792 
1793  if (this->fe_interface->at_boundary())
1794  return (*(this->fe_interface->fe_face_values))[extractor].hessian(
1795  dof_pair[0], q_point);
1796 
1797  hessian_type value;
1798 
1799  if (dof_pair[0] != numbers::invalid_unsigned_int)
1800  value +=
1801  0.5 *
1802  (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
1803  q_point);
1804 
1805  if (dof_pair[1] != numbers::invalid_unsigned_int)
1806  value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
1807  .hessian(dof_pair[1], q_point);
1808 
1809  return value;
1810  }
1811 
1812 
1813 
1814  template <int dim, int spacedim>
1816  Vector<dim, spacedim>::jump_hessian(const unsigned int interface_dof_index,
1817  const unsigned int q_point) const
1818  {
1819  const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1820 
1821  if (this->fe_interface->at_boundary())
1822  return (*(this->fe_interface->fe_face_values))[extractor].hessian(
1823  dof_pair[0], q_point);
1824 
1825  hessian_type value;
1826 
1827  if (dof_pair[0] != numbers::invalid_unsigned_int)
1828  value +=
1829  (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
1830  q_point);
1831 
1832  if (dof_pair[1] != numbers::invalid_unsigned_int)
1833  value -=
1834  (*(this->fe_interface->fe_face_values_neighbor))[extractor].hessian(
1835  dof_pair[1], q_point);
1836 
1837  return value;
1838  }
1839 
1840 
1841 
1842  template <int dim, int spacedim>
1845  const unsigned int interface_dof_index,
1846  const unsigned int q_point) const
1847  {
1848  const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
1849 
1850  if (this->fe_interface->at_boundary())
1851  return (*(this->fe_interface->fe_face_values))[extractor]
1852  .third_derivative(dof_pair[0], q_point);
1853 
1854  third_derivative_type value;
1855 
1856  if (dof_pair[0] != numbers::invalid_unsigned_int)
1857  value +=
1858  (*(this->fe_interface->fe_face_values))[extractor].third_derivative(
1859  dof_pair[0], q_point);
1860 
1861  if (dof_pair[1] != numbers::invalid_unsigned_int)
1862  value -= (*(this->fe_interface->fe_face_values_neighbor))[extractor]
1863  .third_derivative(dof_pair[1], q_point);
1864 
1865  return value;
1866  }
1867 } // namespace FEInterfaceViews
1868 
1869 #endif // DOXYGEN
1870 
1872 
1873 #endif
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
bool at_boundary() 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
Definition: fe_values.h:161
::Tensor< 2, spacedim > hessian_type
Definition: fe_values.h:168
::Tensor< 3, spacedim > third_derivative_type
Definition: fe_values.h:175
unsigned int n_components() const
Definition: vector.h:110
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:396
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:397
UpdateFlags
unsigned int cell_index
Definition: grid_tools.cc:1092
#define Assert(cond, exc)
Definition: exceptions.h:1465
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
static ::ExceptionBase & ExcMessage(std::string arg1)
friend class Vector
Definition: vector.h:1033
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition: mapping.cc:240
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
static const unsigned int invalid_unsigned_int
Definition: types.h:196