deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
fe_interface_values.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2019 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_fe_interface_values_h
16#define dealii_fe_interface_values_h
17
18#include <deal.II/base/config.h>
19
21
23#include <deal.II/fe/mapping.h>
24
28
30
31#ifndef DOXYGEN
32template <int dim, int spacedim>
34#endif
35
41{
45 template <int dim, int spacedim = dim>
46 class Base
47 {
48 public:
53
54 protected:
59
64 template <class InputVector, class OutputVector>
65 void
66 get_local_dof_values(const InputVector &dof_values,
67 OutputVector &local_dof_values) const;
68 };
69
70
71
75 template <int dim, int spacedim = dim>
76 class Scalar : public Base<dim, spacedim>
77 {
78 public:
82 using value_type = double;
83
90
97
104
111 template <typename Number>
113
120 template <typename Number>
123
130 template <typename Number>
133
140 template <typename Number>
143
148 const unsigned int component);
149
154
175 value(const bool here_or_there,
176 const unsigned int interface_dof_index,
177 const unsigned int q_point) const;
178
180
185
197 jump_in_values(const unsigned int interface_dof_index,
198 const unsigned int q_point) const;
199
211 jump_in_gradients(const unsigned int interface_dof_index,
212 const unsigned int q_point) const;
213
225 jump_in_hessians(const unsigned int interface_dof_index,
226 const unsigned int q_point) const;
227
239 jump_in_third_derivatives(const unsigned int interface_dof_index,
240 const unsigned int q_point) const;
241
243
248
260 average_of_values(const unsigned int interface_dof_index,
261 const unsigned int q_point) const;
262
274 average_of_gradients(const unsigned int interface_dof_index,
275 const unsigned int q_point) const;
276
289 average_of_hessians(const unsigned int interface_dof_index,
290 const unsigned int q_point) const;
291
293
298
318 template <class InputVector>
319 void
321 const bool here_or_there,
322 const InputVector &fe_function,
324 &values) const;
325
348 template <class InputVector>
349 void
351 const bool here_or_there,
352 const InputVector &local_dof_values,
354 &values) const;
355
357
362
376 template <class InputVector>
377 void
379 const InputVector &fe_function,
381 &values) const;
382
390 template <class InputVector>
391 void
393 const InputVector &local_dof_values,
395 &values) const;
396
410 template <class InputVector>
411 void
413 const InputVector &fe_function,
415 &gradients) const;
416
424 template <class InputVector>
425 void
427 const InputVector &local_dof_values,
429 &gradients) const;
430
444 template <class InputVector>
445 void
447 const InputVector &fe_function,
449 &hessians) const;
450
458 template <class InputVector>
459 void
461 const InputVector &local_dof_values,
463 &hessians) const;
464
479 template <class InputVector>
480 void
482 const InputVector &fe_function,
483 std::vector<
485 &third_derivatives) const;
486
494 template <class InputVector>
495 void
497 const InputVector &local_dof_values,
498 std::vector<
500 &third_derivatives) const;
501
503
508
522 template <class InputVector>
523 void
525 const InputVector &fe_function,
527 &values) const;
528
536 template <class InputVector>
537 void
539 const InputVector &local_dof_values,
541 &values) const;
542
556 template <class InputVector>
557 void
559 const InputVector &fe_function,
561 &gradients) const;
562
570 template <class InputVector>
571 void
573 const InputVector &local_dof_values,
575 &gradients) const;
576
590 template <class InputVector>
591 void
593 const InputVector &fe_function,
595 &hessians) const;
596
604 template <class InputVector>
605 void
607 const InputVector &local_dof_values,
609 &hessians) const;
610
612
613 private:
618 };
619
620
621
625 template <int dim, int spacedim = dim>
626 class Vector : public Base<dim, spacedim>
627 {
628 public:
634
641
649
657
664 template <typename Number>
666
673 template <typename Number>
676
683 template <typename Number>
686
693 template <typename Number>
696
701 const unsigned int first_vector_component);
702
707
728 value(const bool here_or_there,
729 const unsigned int interface_dof_index,
730 const unsigned int q_point) const;
731
733
738
749 jump_in_values(const unsigned int interface_dof_index,
750 const unsigned int q_point) const;
751
762 jump_in_gradients(const unsigned int interface_dof_index,
763 const unsigned int q_point) const;
764
771 jump_gradient(const unsigned int interface_dof_index,
772 const unsigned int q_point) const;
773
785 jump_in_hessians(const unsigned int interface_dof_index,
786 const unsigned int q_point) const;
787
794 jump_hessian(const unsigned int interface_dof_index,
795 const unsigned int q_point) const;
796
808 jump_in_third_derivatives(const unsigned int interface_dof_index,
809 const unsigned int q_point) const;
810
812
817
828 average_of_values(const unsigned int interface_dof_index,
829 const unsigned int q_point) const;
830
841 average_of_gradients(const unsigned int interface_dof_index,
842 const unsigned int q_point) const;
843
856 average_of_hessians(const unsigned int interface_dof_index,
857 const unsigned int q_point) const;
858
865 average_hessian(const unsigned int interface_dof_index,
866 const unsigned int q_point) const;
867
869
874
894 template <class InputVector>
895 void
897 const bool here_or_there,
898 const InputVector &fe_function,
900 &values) const;
901
924 template <class InputVector>
925 void
927 const bool here_or_there,
928 const InputVector &local_dof_values,
930 &values) const;
931
933
938
952 template <class InputVector>
953 void
955 const InputVector &fe_function,
957 &values) const;
958
966 template <class InputVector>
967 void
969 const InputVector &local_dof_values,
971 &values) const;
972
986 template <class InputVector>
987 void
989 const InputVector &fe_function,
991 &gradients) const;
992
1000 template <class InputVector>
1001 void
1003 const InputVector &local_dof_values,
1005 &gradients) const;
1006
1020 template <class InputVector>
1021 void
1023 const InputVector &fe_function,
1025 &hessians) const;
1026
1034 template <class InputVector>
1035 void
1037 const InputVector &local_dof_values,
1039 &hessians) const;
1040
1055 template <class InputVector>
1056 void
1058 const InputVector &fe_function,
1059 std::vector<
1061 &third_derivatives) const;
1062
1070 template <class InputVector>
1071 void
1073 const InputVector &local_dof_values,
1074 std::vector<
1076 &third_derivatives) const;
1077
1079
1084
1098 template <class InputVector>
1099 void
1101 const InputVector &fe_function,
1103 &values) const;
1104
1112 template <class InputVector>
1113 void
1115 const InputVector &local_dof_values,
1117 &values) const;
1118
1132 template <class InputVector>
1133 void
1135 const InputVector &fe_function,
1137 &gradients) const;
1138
1146 template <class InputVector>
1147 void
1149 const InputVector &local_dof_values,
1151 &gradients) const;
1152
1166 template <class InputVector>
1167 void
1169 const InputVector &fe_function,
1171 &hessians) const;
1172
1180 template <class InputVector>
1181 void
1183 const InputVector &local_dof_values,
1185 &hessians) const;
1186
1188
1189 private:
1194 };
1195} // namespace FEInterfaceViews
1196
1197
1198namespace internal
1199{
1201 {
1206 template <int dim, int spacedim, typename Extractor>
1208 {};
1209
1217 template <int dim, int spacedim>
1218 struct ViewType<dim, spacedim, FEValuesExtractors::Scalar>
1219 {
1220 using type = typename ::FEInterfaceViews::Scalar<dim, spacedim>;
1221 };
1222
1230 template <int dim, int spacedim>
1231 struct ViewType<dim, spacedim, FEValuesExtractors::Vector>
1232 {
1233 using type = typename ::FEInterfaceViews::Vector<dim, spacedim>;
1234 };
1235 } // namespace FEInterfaceViews
1236} // namespace internal
1237
1238namespace FEInterfaceViews
1239{
1244 template <int dim, int spacedim, typename Extractor>
1245 using View = typename ::internal::FEInterfaceViews::
1246 ViewType<dim, spacedim, Extractor>::type;
1247} // namespace FEInterfaceViews
1248
1249
1250
1275template <int dim, int spacedim = dim>
1277{
1278public:
1282 const unsigned int n_quadrature_points;
1283
1291 const Quadrature<dim - 1> &quadrature,
1292 const UpdateFlags update_flags);
1293
1301 const hp::QCollection<dim - 1> &quadrature,
1302 const UpdateFlags update_flags);
1303
1311 const Quadrature<dim - 1> &quadrature,
1312 const UpdateFlags update_flags);
1313
1319 const hp::MappingCollection<dim, spacedim> &mapping_collection,
1320 const hp::FECollection<dim, spacedim> &fe_collection,
1321 const hp::QCollection<dim - 1> &quadrature_collection,
1322 const UpdateFlags update_flags);
1323
1328 const hp::QCollection<dim - 1> &quadrature_collection,
1329 const UpdateFlags update_flags);
1330
1447 template <typename CellIteratorType, typename CellNeighborIteratorType>
1448 void
1449 reinit(const CellIteratorType &cell,
1450 const unsigned int face_no,
1451 const unsigned int sub_face_no,
1452 const CellNeighborIteratorType &cell_neighbor,
1453 const unsigned int face_no_neighbor,
1454 const unsigned int sub_face_no_neighbor,
1455 const unsigned int q_index = numbers::invalid_unsigned_int,
1456 const unsigned int mapping_index = numbers::invalid_unsigned_int,
1457 const unsigned int fe_index = numbers::invalid_unsigned_int,
1458 const unsigned int fe_index_neighbor = numbers::invalid_unsigned_int);
1459
1486 template <typename CellIteratorType>
1487 void
1488 reinit(const CellIteratorType &cell,
1489 const unsigned int face_no,
1490 const unsigned int q_index = numbers::invalid_unsigned_int,
1491 const unsigned int mapping_index = numbers::invalid_unsigned_int,
1492 const unsigned int fe_index = numbers::invalid_unsigned_int);
1493
1502 get_fe_face_values(const unsigned int cell_index) const;
1503
1509
1514 get_fe() const;
1515
1519 const Quadrature<dim - 1> &
1521
1527
1533
1537 const hp::QCollection<dim - 1> &
1539
1544 bool
1546
1550 UpdateFlags
1552
1560 get_cell(const unsigned int cell_index) const;
1561
1569 unsigned int
1570 get_face_number(const unsigned int cell_index) const;
1571
1576
1583 bool
1585
1597 double
1598 JxW(const unsigned int quadrature_point) const;
1599
1605 const std::vector<double> &
1607
1618 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT("Use the function normal_vector().")
1619 Tensor<1, spacedim>
1620 normal(const unsigned int q_point_index) const;
1621
1630 Tensor<1, spacedim>
1631 normal_vector(const unsigned int q_point_index) const;
1632
1642 const std::vector<Tensor<1, spacedim>> &
1644
1652 std_cxx20::ranges::iota_view<unsigned int, unsigned int>
1654
1661 const Point<spacedim> &
1662 quadrature_point(const unsigned int q_point) const;
1663
1669 const std::vector<Point<spacedim>> &
1671
1682 unsigned
1684
1712 std_cxx20::ranges::iota_view<unsigned int, unsigned int>
1714
1725 std::vector<types::global_dof_index>
1727
1740 std::array<unsigned int, 2>
1741 interface_dof_to_dof_indices(const unsigned int interface_dof_index) const;
1742
1746
1751
1773 double
1774 shape_value(const bool here_or_there,
1775 const unsigned int interface_dof_index,
1776 const unsigned int q_point,
1777 const unsigned int component = 0) const;
1778
1788 Tensor<1, spacedim>
1789 shape_grad(const bool here_or_there,
1790 const unsigned int interface_dof_index,
1791 const unsigned int q_point,
1792 const unsigned int component = 0) const;
1793
1797
1802
1821 double
1822 jump_in_shape_values(const unsigned int interface_dof_index,
1823 const unsigned int q_point,
1824 const unsigned int component = 0) const;
1825
1839 Tensor<1, spacedim>
1840 jump_in_shape_gradients(const unsigned int interface_dof_index,
1841 const unsigned int q_point,
1842 const unsigned int component = 0) const;
1843
1858 Tensor<2, spacedim>
1859 jump_in_shape_hessians(const unsigned int interface_dof_index,
1860 const unsigned int q_point,
1861 const unsigned int component = 0) const;
1862
1876 Tensor<3, spacedim>
1877 jump_in_shape_3rd_derivatives(const unsigned int interface_dof_index,
1878 const unsigned int q_point,
1879 const unsigned int component = 0) const;
1880
1884
1889
1903 double
1904 average_of_shape_values(const unsigned int interface_dof_index,
1905 const unsigned int q_point,
1906 const unsigned int component = 0) const;
1907
1921 Tensor<1, spacedim>
1922 average_of_shape_gradients(const unsigned int interface_dof_index,
1923 const unsigned int q_point,
1924 const unsigned int component = 0) const;
1925
1940 Tensor<2, spacedim>
1941 average_of_shape_hessians(const unsigned int interface_dof_index,
1942 const unsigned int q_point,
1943 const unsigned int component = 0) const;
1944
1948
1949
1950
1955
1964 template <class InputVector>
1965 void
1967 const InputVector &fe_function,
1968 std::vector<typename InputVector::value_type> &values) const;
1969
1978 template <class InputVector>
1979 void
1981 const InputVector &fe_function,
1982 std::vector<Tensor<1, spacedim, typename InputVector::value_type>>
1983 &gradients) const;
1984
1992 template <class InputVector>
1993 void
1995 const InputVector &fe_function,
1996 std::vector<Tensor<2, spacedim, typename InputVector::value_type>>
1997 &hessians) const;
1998
2007 template <class InputVector>
2008 void
2010 const InputVector &fe_function,
2011 std::vector<Tensor<3, spacedim, typename InputVector::value_type>>
2012 &third_derivatives) const;
2013
2015
2020
2029 template <class InputVector>
2030 void
2032 const InputVector &fe_function,
2033 std::vector<typename InputVector::value_type> &values) const;
2034
2042 template <class InputVector>
2043 void
2045 const InputVector &fe_function,
2046 std::vector<Tensor<1, spacedim, typename InputVector::value_type>>
2047 &gradients) const;
2048
2056 template <class InputVector>
2057 void
2059 const InputVector &fe_function,
2060 std::vector<Tensor<2, spacedim, typename InputVector::value_type>>
2061 &hessians) const;
2062
2066
2067
2068
2073
2080 FEInterfaceViews::Scalar<dim, spacedim>
2081 operator[](const FEValuesExtractors::Scalar &scalar) const;
2082
2089 FEInterfaceViews::Vector<dim, spacedim>
2090 operator[](const FEValuesExtractors::Vector &vector) const;
2091
2095
2096private:
2100 std::vector<types::global_dof_index> interface_dof_indices;
2101
2107 std::vector<std::array<unsigned int, 2>> dofmap;
2108
2114
2121 // non-hp data
2126
2130 std::unique_ptr<FEFaceValues<dim, spacedim>> internal_fe_face_values;
2131
2136
2141
2145 std::unique_ptr<FESubfaceValues<dim, spacedim>>
2147 // non-hp data
2149 // hp data
2154
2159 std::unique_ptr<hp::FEFaceValues<dim, spacedim>> internal_hp_fe_face_values;
2160
2165 std::unique_ptr<hp::FESubfaceValues<dim, spacedim>>
2167
2172 std::unique_ptr<hp::FEFaceValues<dim, spacedim>>
2174
2179 std::unique_ptr<hp::FESubfaceValues<dim, spacedim>>
2181
2189 "The current function doesn't make sense when used with a "
2190 "FEInterfaceValues object with hp-capabilities.");
2191
2199 "The current function doesn't make sense when used with a "
2200 "FEInterfaceValues object without hp-capabilities.");
2201 // hp data
2203
2204 /*
2205 * Make the view classes friends of this class, since they
2206 * access internal data.
2207 */
2208 template <int, int>
2210 template <int, int>
2212};
2213
2214
2215
2216#ifndef DOXYGEN
2217
2218/*---------------------- Inline functions ---------------------*/
2219
2220template <int dim, int spacedim>
2224 const Quadrature<dim - 1> &quadrature,
2225 const UpdateFlags update_flags)
2226 : n_quadrature_points(quadrature.size())
2227 , fe_face_values(nullptr)
2228 , fe_face_values_neighbor(nullptr)
2230 std::make_unique<FEFaceValues<dim, spacedim>>(mapping,
2231 fe,
2232 quadrature,
2233 update_flags))
2235 std::make_unique<FESubfaceValues<dim, spacedim>>(mapping,
2236 fe,
2237 quadrature,
2238 update_flags))
2240 std::make_unique<FEFaceValues<dim, spacedim>>(mapping,
2241 fe,
2242 quadrature,
2243 update_flags))
2245 std::make_unique<FESubfaceValues<dim, spacedim>>(mapping,
2246 fe,
2247 quadrature,
2248 update_flags))
2249{}
2250
2251
2252
2253template <int dim, int spacedim>
2256 const Quadrature<dim - 1> &quadrature,
2257 const UpdateFlags update_flags)
2259 fe.reference_cell().template get_default_linear_mapping<dim, spacedim>(),
2260 fe,
2261 quadrature,
2262 update_flags)
2263{}
2264
2265
2266
2267template <int dim, int spacedim>
2271 const hp::QCollection<dim - 1> &quadrature,
2272 const UpdateFlags update_flags)
2273 : n_quadrature_points(quadrature.max_n_quadrature_points())
2274 , fe_face_values(nullptr)
2275 , fe_face_values_neighbor(nullptr)
2276 , internal_fe_face_values(
2277 std::make_unique<FEFaceValues<dim, spacedim>>(mapping,
2278 fe,
2279 quadrature,
2280 update_flags))
2281 , internal_fe_subface_values(
2282 std::make_unique<FESubfaceValues<dim, spacedim>>(mapping,
2283 fe,
2284 quadrature,
2285 update_flags))
2286 , internal_fe_face_values_neighbor(
2287 std::make_unique<FEFaceValues<dim, spacedim>>(mapping,
2288 fe,
2289 quadrature[0],
2290 update_flags))
2291 , internal_fe_subface_values_neighbor(
2292 std::make_unique<FESubfaceValues<dim, spacedim>>(mapping,
2293 fe,
2294 quadrature[0],
2295 update_flags))
2296{}
2297
2298
2299
2300template <int dim, int spacedim>
2302 const hp::MappingCollection<dim, spacedim> &mapping_collection,
2303 const hp::FECollection<dim, spacedim> &fe_collection,
2304 const hp::QCollection<dim - 1> &quadrature_collection,
2305 const UpdateFlags update_flags)
2306 : n_quadrature_points(quadrature_collection.max_n_quadrature_points())
2307 , fe_face_values(nullptr)
2308 , fe_face_values_neighbor(nullptr)
2309 , internal_hp_fe_face_values(
2310 std::make_unique<hp::FEFaceValues<dim, spacedim>>(mapping_collection,
2311 fe_collection,
2312 quadrature_collection,
2313 update_flags))
2314 , internal_hp_fe_subface_values(
2315 std::make_unique<hp::FESubfaceValues<dim, spacedim>>(
2316 mapping_collection,
2317 fe_collection,
2318 quadrature_collection,
2319 update_flags))
2320 , internal_hp_fe_face_values_neighbor(
2321 std::make_unique<hp::FEFaceValues<dim, spacedim>>(mapping_collection,
2322 fe_collection,
2323 quadrature_collection,
2324 update_flags))
2325 , internal_hp_fe_subface_values_neighbor(
2326 std::make_unique<hp::FESubfaceValues<dim, spacedim>>(
2327 mapping_collection,
2328 fe_collection,
2329 quadrature_collection,
2330 update_flags))
2331{
2332 AssertDimension(dim, spacedim);
2333}
2334
2335
2336
2337template <int dim, int spacedim>
2339 const hp::FECollection<dim, spacedim> &fe_collection,
2340 const hp::QCollection<dim - 1> &quadrature_collection,
2341 const UpdateFlags update_flags)
2342 : FEInterfaceValues(fe_collection.get_reference_cell_default_linear_mapping(),
2343 fe_collection,
2344 quadrature_collection,
2345 update_flags)
2346{}
2347
2348
2349
2350template <int dim, int spacedim>
2351template <typename CellIteratorType, typename CellNeighborIteratorType>
2352void
2354 const CellIteratorType &cell,
2355 const unsigned int face_no,
2356 const unsigned int sub_face_no,
2357 const CellNeighborIteratorType &cell_neighbor,
2358 const unsigned int face_no_neighbor,
2359 const unsigned int sub_face_no_neighbor,
2360 const unsigned int q_index,
2361 const unsigned int mapping_index,
2362 const unsigned int fe_index_in,
2363 const unsigned int fe_index_neighbor_in)
2364{
2365 Assert(internal_fe_face_values || internal_hp_fe_face_values,
2367
2368 constexpr bool is_dof_cell_accessor =
2369 std::is_same_v<DoFCellAccessor<dim, spacedim, true>,
2370 typename CellIteratorType::AccessorType> ||
2371 std::is_same_v<DoFCellAccessor<dim, spacedim, false>,
2372 typename CellIteratorType::AccessorType>;
2373
2374 constexpr bool is_dof_cell_accessor_neighbor =
2375 std::is_same_v<DoFCellAccessor<dim, spacedim, true>,
2376 typename CellNeighborIteratorType::AccessorType> ||
2377 std::is_same_v<DoFCellAccessor<dim, spacedim, false>,
2378 typename CellNeighborIteratorType::AccessorType>;
2379
2380 unsigned int active_fe_index = 0;
2381 unsigned int active_fe_index_neighbor = 0;
2382
2383 if (internal_fe_face_values)
2384 {
2385 if (sub_face_no == numbers::invalid_unsigned_int)
2386 {
2387 internal_fe_face_values->reinit(cell, face_no);
2388 fe_face_values = internal_fe_face_values.get();
2389 }
2390 else
2391 {
2392 internal_fe_subface_values->reinit(cell, face_no, sub_face_no);
2393 fe_face_values = internal_fe_subface_values.get();
2394 }
2395 if (sub_face_no_neighbor == numbers::invalid_unsigned_int)
2396 {
2397 internal_fe_face_values_neighbor->reinit(cell_neighbor,
2398 face_no_neighbor);
2399 fe_face_values_neighbor = internal_fe_face_values_neighbor.get();
2400 }
2401 else
2402 {
2403 internal_fe_subface_values_neighbor->reinit(cell_neighbor,
2404 face_no_neighbor,
2405 sub_face_no_neighbor);
2406 fe_face_values_neighbor = internal_fe_subface_values_neighbor.get();
2407 }
2408
2409 AssertDimension(fe_face_values->n_quadrature_points,
2410 fe_face_values_neighbor->n_quadrature_points);
2411
2412 const_cast<unsigned int &>(this->n_quadrature_points) =
2413 fe_face_values->n_quadrature_points;
2414 }
2415 else if (internal_hp_fe_face_values)
2416 {
2417 active_fe_index = fe_index_in;
2418 active_fe_index_neighbor =
2419 (fe_index_neighbor_in != numbers::invalid_unsigned_int) ?
2420 fe_index_neighbor_in :
2421 fe_index_in;
2422
2423 if (active_fe_index == numbers::invalid_unsigned_int)
2424 {
2425 if constexpr (is_dof_cell_accessor)
2426 active_fe_index = cell->active_fe_index();
2427 else
2428 active_fe_index = 0;
2429 }
2430
2431 if (active_fe_index_neighbor == numbers::invalid_unsigned_int)
2432 {
2433 if constexpr (is_dof_cell_accessor_neighbor)
2434 active_fe_index_neighbor = cell_neighbor->active_fe_index();
2435 else
2436 active_fe_index_neighbor = 0;
2437 }
2438
2439 unsigned int used_q_index = q_index;
2440 unsigned int used_mapping_index = mapping_index;
2441
2442 // First check. If there is only one element in a collection, and if none
2443 // had been specified explicitly, then that's clearly the one to take:
2444 if (used_q_index == numbers::invalid_unsigned_int)
2445 if (internal_hp_fe_face_values->get_quadrature_collection().size() == 1)
2446 used_q_index = 0;
2447
2448 if (used_mapping_index == numbers::invalid_unsigned_int)
2449 if (internal_hp_fe_face_values->get_mapping_collection().size() == 1)
2450 used_mapping_index = 0;
2451
2452 // Second check: See if the two quadrature objects are the same, because
2453 // in that case it does not matter which one we use. Unfortunately, we
2454 // currently have no way of testing that two mapping objects are the
2455 // same :-(
2456 if (used_q_index == numbers::invalid_unsigned_int)
2457 if (internal_hp_fe_face_values
2458 ->get_quadrature_collection()[active_fe_index] ==
2459 internal_hp_fe_face_values
2460 ->get_quadrature_collection()[active_fe_index_neighbor])
2461 used_q_index = active_fe_index;
2462
2463 // Third check, if the above did not already suffice. We see if we
2464 // can get somewhere via the dominated's finite element index.
2465 if ((used_q_index == numbers::invalid_unsigned_int) ||
2466 (used_mapping_index == numbers::invalid_unsigned_int))
2467 {
2468 const unsigned int dominated_fe_index =
2469 ((used_q_index == numbers::invalid_unsigned_int) ||
2470 (used_mapping_index == numbers::invalid_unsigned_int) ?
2471 internal_hp_fe_face_values->get_fe_collection()
2472 .find_dominated_fe(
2473 {active_fe_index, active_fe_index_neighbor}) :
2475
2476 if (used_q_index == numbers::invalid_unsigned_int)
2477 {
2478 Assert(dominated_fe_index != numbers::invalid_fe_index,
2479 ExcMessage(
2480 "You called this function with 'q_index' left at its "
2481 "default value, but this can only work if one of "
2482 "the two finite elements adjacent to this face "
2483 "dominates the other. See the documentation "
2484 "of this function for more information of how "
2485 "to deal with this situation."));
2486 used_q_index = dominated_fe_index;
2487 }
2488
2489 if (used_mapping_index == numbers::invalid_unsigned_int)
2490 {
2491 Assert(dominated_fe_index != numbers::invalid_fe_index,
2492 ExcMessage(
2493 "You called this function with 'mapping_index' left "
2494 "at its default value, but this can only work if one "
2495 "of the two finite elements adjacent to this face "
2496 "dominates the other. See the documentation "
2497 "of this function for more information of how "
2498 "to deal with this situation."));
2499 used_mapping_index = dominated_fe_index;
2500 }
2501 }
2502
2503 // Same as if above, but when hp is enabled.
2504 if (sub_face_no == numbers::invalid_unsigned_int)
2505 {
2506 internal_hp_fe_face_values->reinit(
2507 cell, face_no, used_q_index, used_mapping_index, active_fe_index);
2508 fe_face_values = &const_cast<FEFaceValues<dim, spacedim> &>(
2509 internal_hp_fe_face_values->get_present_fe_values());
2510 }
2511 else
2512 {
2513 internal_hp_fe_subface_values->reinit(cell,
2514 face_no,
2515 sub_face_no,
2516 used_q_index,
2517 used_mapping_index,
2518 active_fe_index);
2519
2520 fe_face_values = &const_cast<FESubfaceValues<dim, spacedim> &>(
2521 internal_hp_fe_subface_values->get_present_fe_values());
2522 }
2523 if (sub_face_no_neighbor == numbers::invalid_unsigned_int)
2524 {
2525 internal_hp_fe_face_values_neighbor->reinit(cell_neighbor,
2526 face_no_neighbor,
2527 used_q_index,
2528 used_mapping_index,
2529 active_fe_index_neighbor);
2530
2531 fe_face_values_neighbor = &const_cast<FEFaceValues<dim, spacedim> &>(
2532 internal_hp_fe_face_values_neighbor->get_present_fe_values());
2533 }
2534 else
2535 {
2536 internal_hp_fe_subface_values_neighbor->reinit(
2537 cell_neighbor,
2538 face_no_neighbor,
2539 sub_face_no_neighbor,
2540 used_q_index,
2541 used_mapping_index,
2542 active_fe_index_neighbor);
2543
2544 fe_face_values_neighbor =
2545 &const_cast<FESubfaceValues<dim, spacedim> &>(
2546 internal_hp_fe_subface_values_neighbor->get_present_fe_values());
2547 }
2548
2549 AssertDimension(fe_face_values->n_quadrature_points,
2550 fe_face_values_neighbor->n_quadrature_points);
2551
2552 const_cast<unsigned int &>(this->n_quadrature_points) =
2553 fe_face_values->n_quadrature_points;
2554 }
2555
2556 // Set up dof mapping and remove duplicates (for continuous elements).
2557 if constexpr (is_dof_cell_accessor_neighbor && is_dof_cell_accessor)
2558 {
2559 // Get dof indices first:
2560 std::vector<types::global_dof_index> v(
2561 fe_face_values->get_fe().n_dofs_per_cell());
2562 cell->get_active_or_mg_dof_indices(v);
2563 std::vector<types::global_dof_index> v2(
2564 fe_face_values_neighbor->get_fe().n_dofs_per_cell());
2565 cell_neighbor->get_active_or_mg_dof_indices(v2);
2566
2567 // Fill a map from the global dof index to the left and right
2568 // local index.
2569 std::map<types::global_dof_index, std::pair<unsigned int, unsigned int>>
2570 tempmap;
2571 std::pair<unsigned int, unsigned int> invalid_entry(
2573
2574 for (unsigned int i = 0; i < v.size(); ++i)
2575 {
2576 // If not already existing, add an invalid entry:
2577 auto result = tempmap.insert(std::make_pair(v[i], invalid_entry));
2578 result.first->second.first = i;
2579 }
2580
2581 for (unsigned int i = 0; i < v2.size(); ++i)
2582 {
2583 // If not already existing, add an invalid entry:
2584 auto result = tempmap.insert(std::make_pair(v2[i], invalid_entry));
2585 result.first->second.second = i;
2586 }
2587
2588 // Transfer from the map to the sorted std::vectors.
2589 dofmap.resize(tempmap.size());
2590 interface_dof_indices.resize(tempmap.size());
2591 unsigned int idx = 0;
2592 for (auto &x : tempmap)
2593 {
2594 interface_dof_indices[idx] = x.first;
2595 dofmap[idx] = {{x.second.first, x.second.second}};
2596 ++idx;
2597 }
2598 }
2599 else
2600 {
2601 const unsigned int n_dofs_per_cell_1 = fe_face_values->dofs_per_cell;
2602 const unsigned int n_dofs_per_cell_2 =
2603 fe_face_values_neighbor->dofs_per_cell;
2604
2605 interface_dof_indices.resize(n_dofs_per_cell_1 + n_dofs_per_cell_2);
2606 dofmap.resize(n_dofs_per_cell_1 + n_dofs_per_cell_2);
2607
2608 for (unsigned int i = 0; i < n_dofs_per_cell_1; ++i)
2609 {
2610 interface_dof_indices[i] = numbers::invalid_dof_index;
2611 dofmap[i] = {{i, numbers::invalid_unsigned_int}};
2612 }
2613
2614 for (unsigned int i = 0; i < n_dofs_per_cell_2; ++i)
2615 {
2616 interface_dof_indices[i + n_dofs_per_cell_1] =
2618 dofmap[i + n_dofs_per_cell_1] = {{numbers::invalid_unsigned_int, i}};
2619 }
2620 }
2621}
2622
2623
2624
2625template <int dim, int spacedim>
2626template <typename CellIteratorType>
2627void
2628FEInterfaceValues<dim, spacedim>::reinit(const CellIteratorType &cell,
2629 const unsigned int face_no,
2630 const unsigned int q_index,
2631 const unsigned int mapping_index,
2632 const unsigned int fe_index)
2633{
2634 Assert(internal_fe_face_values || internal_hp_fe_face_values,
2636
2637 if (internal_fe_face_values)
2638 {
2639 Assert((q_index == 0 || q_index == numbers::invalid_unsigned_int),
2641 Assert((mapping_index == 0 ||
2642 mapping_index == numbers::invalid_unsigned_int),
2644 Assert((fe_index == 0 || fe_index == numbers::invalid_unsigned_int),
2646
2647 internal_fe_face_values->reinit(cell, face_no);
2648 fe_face_values = internal_fe_face_values.get();
2649 fe_face_values_neighbor = nullptr;
2650 }
2651 else if (internal_hp_fe_face_values)
2652 {
2653 internal_hp_fe_face_values->reinit(
2654 cell, face_no, q_index, mapping_index, fe_index);
2655 fe_face_values = &const_cast<FEFaceValues<dim> &>(
2656 internal_hp_fe_face_values->get_present_fe_values());
2657 fe_face_values_neighbor = nullptr;
2658 }
2659
2660 interface_dof_indices.resize(fe_face_values->get_fe().n_dofs_per_cell());
2661
2662 if constexpr (std::is_same_v<typename CellIteratorType::AccessorType,
2664 std::is_same_v<typename CellIteratorType::AccessorType,
2666 {
2667 cell->get_active_or_mg_dof_indices(interface_dof_indices);
2668 }
2669 else
2670 {
2671 for (auto &i : interface_dof_indices)
2673 }
2674
2675 dofmap.resize(interface_dof_indices.size());
2676 for (unsigned int i = 0; i < interface_dof_indices.size(); ++i)
2677 dofmap[i] = {{i, numbers::invalid_unsigned_int}};
2678}
2679
2680
2681
2682template <int dim, int spacedim>
2683inline double
2684FEInterfaceValues<dim, spacedim>::JxW(const unsigned int q) const
2685{
2686 Assert(fe_face_values != nullptr,
2687 ExcMessage("This call requires a call to reinit() first."));
2688 return fe_face_values->JxW(q);
2689}
2690
2691
2692
2693template <int dim, int spacedim>
2694const std::vector<double> &
2696{
2697 Assert(fe_face_values != nullptr,
2698 ExcMessage("This call requires a call to reinit() first."));
2699 return fe_face_values->get_JxW_values();
2700}
2701
2702
2703
2704template <int dim, int spacedim>
2705const std::vector<Tensor<1, spacedim>> &
2707{
2708 Assert(fe_face_values != nullptr,
2709 ExcMessage("This call requires a call to reinit() first."));
2710 return fe_face_values->get_normal_vectors();
2711}
2712
2713
2714
2715template <int dim, int spacedim>
2718{
2719 Assert(!has_hp_capabilities(), ExcOnlyAvailableWithoutHP());
2720 return internal_fe_face_values->get_mapping();
2721}
2722
2723
2724
2725template <int dim, int spacedim>
2728{
2729 Assert(!has_hp_capabilities(), ExcOnlyAvailableWithoutHP());
2730 return internal_fe_face_values->get_fe();
2731}
2732
2733
2734
2735template <int dim, int spacedim>
2736const Quadrature<dim - 1> &
2738{
2739 Assert(!has_hp_capabilities(), ExcOnlyAvailableWithoutHP());
2740 return internal_fe_face_values->get_quadrature();
2741}
2742
2743
2744
2745template <int dim, int spacedim>
2748{
2749 Assert(has_hp_capabilities(), ExcOnlyAvailableWithHP());
2750 return internal_hp_fe_face_values->get_mapping_collection();
2751}
2752
2753
2754
2755template <int dim, int spacedim>
2758{
2759 Assert(has_hp_capabilities(), ExcOnlyAvailableWithHP());
2760 return internal_hp_fe_face_values->get_fe_collection();
2761}
2762
2763
2764
2765template <int dim, int spacedim>
2766const hp::QCollection<dim - 1> &
2768{
2769 Assert(has_hp_capabilities(), ExcOnlyAvailableWithHP());
2770 return internal_hp_fe_face_values->get_quadrature_collection();
2771}
2772
2773
2774
2775template <int dim, int spacedim>
2776bool
2778{
2779 if (internal_hp_fe_face_values || internal_hp_fe_subface_values ||
2780 internal_hp_fe_face_values_neighbor ||
2781 internal_hp_fe_subface_values_neighbor)
2782 {
2783 Assert(!internal_fe_face_values, ExcInternalError());
2784 Assert(!internal_fe_subface_values, ExcInternalError());
2785 Assert(!internal_fe_face_values_neighbor, ExcInternalError());
2786 Assert(!internal_fe_subface_values_neighbor, ExcInternalError());
2787
2788 return true;
2789 }
2790
2791 Assert(internal_fe_face_values || internal_fe_subface_values ||
2792 internal_fe_face_values_neighbor ||
2793 internal_fe_subface_values_neighbor,
2795 Assert(!internal_hp_fe_face_values, ExcInternalError());
2796 Assert(!internal_hp_fe_subface_values, ExcInternalError());
2797 Assert(!internal_hp_fe_face_values_neighbor, ExcInternalError());
2798 Assert(!internal_hp_fe_subface_values_neighbor, ExcInternalError());
2799
2800 return false;
2801}
2802
2803
2804
2805template <int dim, int spacedim>
2808{
2810 0U, n_quadrature_points);
2811}
2812
2813
2814
2815template <int dim, int spacedim>
2816const Point<spacedim> &
2818 const unsigned int q_point) const
2819{
2820 Assert(fe_face_values != nullptr,
2821 ExcMessage("This call requires a call to reinit() first."));
2822 return fe_face_values->quadrature_point(q_point);
2823}
2824
2825
2826
2827template <int dim, int spacedim>
2828const std::vector<Point<spacedim>> &
2830{
2831 Assert(fe_face_values != nullptr,
2832 ExcMessage("This call requires a call to reinit() first."));
2833 return fe_face_values->get_quadrature_points();
2834}
2835
2836
2837
2838template <int dim, int spacedim>
2841{
2842 if (has_hp_capabilities())
2843 return internal_hp_fe_face_values->get_update_flags();
2844 else
2845 return internal_fe_face_values->get_update_flags();
2846}
2847
2848
2849
2850template <int dim, int spacedim>
2852FEInterfaceValues<dim, spacedim>::get_cell(const unsigned int cell_index) const
2853{
2854 return get_fe_face_values(cell_index).get_cell();
2855}
2856
2857
2858
2859template <int dim, int spacedim>
2860inline unsigned int
2862 const unsigned int cell_index) const
2863{
2864 return get_fe_face_values(cell_index).get_face_number();
2865}
2866
2867
2868
2869template <int dim, int spacedim>
2870unsigned
2872{
2873 Assert(
2874 interface_dof_indices.size() > 0,
2875 ExcMessage(
2876 "n_current_interface_dofs() is only available after a call to reinit()."));
2877 return interface_dof_indices.size();
2878}
2879
2880
2881
2882template <int dim, int spacedim>
2885{
2886 return {0U, n_current_interface_dofs()};
2887}
2888
2889
2890
2891template <int dim, int spacedim>
2892bool
2894{
2895 return fe_face_values_neighbor == nullptr;
2896}
2897
2898
2899
2900template <int dim, int spacedim>
2901std::vector<types::global_dof_index>
2903{
2904 return interface_dof_indices;
2905}
2906
2907
2908
2909template <int dim, int spacedim>
2910std::array<unsigned int, 2>
2912 const unsigned int interface_dof_index) const
2913{
2914 AssertIndexRange(interface_dof_index, n_current_interface_dofs());
2915 return dofmap[interface_dof_index];
2916}
2917
2918
2919
2920template <int dim, int spacedim>
2923 const unsigned int cell_index) const
2924{
2925 AssertIndexRange(cell_index, 2);
2926 Assert(
2927 cell_index == 0 || !at_boundary(),
2928 ExcMessage(
2929 "You are on a boundary, so you can only ask for the first FEFaceValues object."));
2930
2931 return (cell_index == 0) ? *fe_face_values : *fe_face_values_neighbor;
2932}
2933
2934
2935
2936template <int dim, int spacedim>
2938FEInterfaceValues<dim, spacedim>::normal(const unsigned int q_point_index) const
2939{
2940 return normal_vector(q_point_index);
2941}
2942
2943
2944
2945template <int dim, int spacedim>
2948 const unsigned int q_point_index) const
2949{
2950 return fe_face_values->normal_vector(q_point_index);
2951}
2952
2953
2954
2955template <int dim, int spacedim>
2956double
2958 const bool here_or_there,
2959 const unsigned int interface_dof_index,
2960 const unsigned int q_point,
2961 const unsigned int component) const
2962{
2963 const auto dof_pair = dofmap[interface_dof_index];
2964
2965 if (here_or_there && dof_pair[0] != numbers::invalid_unsigned_int)
2966 return get_fe_face_values(0).shape_value_component(dof_pair[0],
2967 q_point,
2968 component);
2969 if (!here_or_there && dof_pair[1] != numbers::invalid_unsigned_int)
2970 return get_fe_face_values(1).shape_value_component(dof_pair[1],
2971 q_point,
2972 component);
2973
2974 return 0.0;
2975}
2976
2977
2978template <int dim, int spacedim>
2981 const bool here_or_there,
2982 const unsigned int interface_dof_index,
2983 const unsigned int q_point,
2984 const unsigned int component) const
2985{
2986 const auto dof_pair = dofmap[interface_dof_index];
2987
2989
2990 if (here_or_there && dof_pair[0] != numbers::invalid_unsigned_int)
2991 value = get_fe_face_values(0).shape_grad_component(dof_pair[0],
2992 q_point,
2993 component);
2994 if (!here_or_there && dof_pair[1] != numbers::invalid_unsigned_int)
2995 value = get_fe_face_values(1).shape_grad_component(dof_pair[1],
2996 q_point,
2997 component);
2998
2999 return value;
3000}
3001
3002template <int dim, int spacedim>
3003double
3005 const unsigned int interface_dof_index,
3006 const unsigned int q_point,
3007 const unsigned int component) const
3008{
3009 const auto dof_pair = dofmap[interface_dof_index];
3010
3011 double value = 0.0;
3012
3013 if (dof_pair[0] != numbers::invalid_unsigned_int)
3014 value += get_fe_face_values(0).shape_value_component(dof_pair[0],
3015 q_point,
3016 component);
3017 if (dof_pair[1] != numbers::invalid_unsigned_int)
3018 value -= get_fe_face_values(1).shape_value_component(dof_pair[1],
3019 q_point,
3020 component);
3021 return value;
3022}
3023
3024
3025
3026template <int dim, int spacedim>
3027double
3029 const unsigned int interface_dof_index,
3030 const unsigned int q_point,
3031 const unsigned int component) const
3032{
3033 const auto dof_pair = dofmap[interface_dof_index];
3034
3035 if (at_boundary())
3036 return get_fe_face_values(0).shape_value_component(dof_pair[0],
3037 q_point,
3038 component);
3039
3040 double value = 0.0;
3041
3042 if (dof_pair[0] != numbers::invalid_unsigned_int)
3043 value += 0.5 * get_fe_face_values(0).shape_value_component(dof_pair[0],
3044 q_point,
3045 component);
3046 if (dof_pair[1] != numbers::invalid_unsigned_int)
3047 value += 0.5 * get_fe_face_values(1).shape_value_component(dof_pair[1],
3048 q_point,
3049 component);
3050
3051 return value;
3052}
3053
3054
3055
3056template <int dim, int spacedim>
3059 const unsigned int interface_dof_index,
3060 const unsigned int q_point,
3061 const unsigned int component) const
3062{
3063 const auto dof_pair = dofmap[interface_dof_index];
3064
3065 if (at_boundary())
3066 return get_fe_face_values(0).shape_grad_component(dof_pair[0],
3067 q_point,
3068 component);
3069
3071
3072 if (dof_pair[0] != numbers::invalid_unsigned_int)
3073 value += 0.5 * get_fe_face_values(0).shape_grad_component(dof_pair[0],
3074 q_point,
3075 component);
3076 if (dof_pair[1] != numbers::invalid_unsigned_int)
3077 value += 0.5 * get_fe_face_values(1).shape_grad_component(dof_pair[1],
3078 q_point,
3079 component);
3080
3081 return value;
3082}
3083
3084
3085
3086template <int dim, int spacedim>
3089 const unsigned int interface_dof_index,
3090 const unsigned int q_point,
3091 const unsigned int component) const
3092{
3093 const auto dof_pair = dofmap[interface_dof_index];
3094
3095 if (at_boundary())
3096 return get_fe_face_values(0).shape_hessian_component(dof_pair[0],
3097 q_point,
3098 component);
3099
3101
3102 if (dof_pair[0] != numbers::invalid_unsigned_int)
3103 value += 0.5 * get_fe_face_values(0).shape_hessian_component(dof_pair[0],
3104 q_point,
3105 component);
3106 if (dof_pair[1] != numbers::invalid_unsigned_int)
3107 value += 0.5 * get_fe_face_values(1).shape_hessian_component(dof_pair[1],
3108 q_point,
3109 component);
3110
3111 return value;
3112}
3113
3114
3115
3116template <int dim, int spacedim>
3119 const unsigned int interface_dof_index,
3120 const unsigned int q_point,
3121 const unsigned int component) const
3122{
3123 const auto dof_pair = dofmap[interface_dof_index];
3124
3125 if (at_boundary())
3126 return get_fe_face_values(0).shape_grad_component(dof_pair[0],
3127 q_point,
3128 component);
3129
3131
3132 if (dof_pair[0] != numbers::invalid_unsigned_int)
3133 value += get_fe_face_values(0).shape_grad_component(dof_pair[0],
3134 q_point,
3135 component);
3136 if (dof_pair[1] != numbers::invalid_unsigned_int)
3137 value -= get_fe_face_values(1).shape_grad_component(dof_pair[1],
3138 q_point,
3139 component);
3140
3141 return value;
3142}
3143
3144
3145
3146template <int dim, int spacedim>
3149 const unsigned int interface_dof_index,
3150 const unsigned int q_point,
3151 const unsigned int component) const
3152{
3153 const auto dof_pair = dofmap[interface_dof_index];
3154
3155 if (at_boundary())
3156 return get_fe_face_values(0).shape_hessian_component(dof_pair[0],
3157 q_point,
3158 component);
3159
3161
3162 if (dof_pair[0] != numbers::invalid_unsigned_int)
3163 value += get_fe_face_values(0).shape_hessian_component(dof_pair[0],
3164 q_point,
3165 component);
3166 if (dof_pair[1] != numbers::invalid_unsigned_int)
3167 value -= get_fe_face_values(1).shape_hessian_component(dof_pair[1],
3168 q_point,
3169 component);
3170
3171 return value;
3172}
3173
3174
3175
3176template <int dim, int spacedim>
3179 const unsigned int interface_dof_index,
3180 const unsigned int q_point,
3181 const unsigned int component) const
3182{
3183 const auto dof_pair = dofmap[interface_dof_index];
3184
3185 if (at_boundary())
3186 return get_fe_face_values(0).shape_3rd_derivative_component(dof_pair[0],
3187 q_point,
3188 component);
3189
3191
3192 if (dof_pair[0] != numbers::invalid_unsigned_int)
3193 value += get_fe_face_values(0).shape_3rd_derivative_component(dof_pair[0],
3194 q_point,
3195 component);
3196 if (dof_pair[1] != numbers::invalid_unsigned_int)
3197 value -= get_fe_face_values(1).shape_3rd_derivative_component(dof_pair[1],
3198 q_point,
3199 component);
3200
3201 return value;
3202}
3203
3204
3205
3206template <int dim, int spacedim>
3207template <class InputVector>
3208void
3210 const InputVector &fe_function,
3211 std::vector<typename InputVector::value_type> &values) const
3212{
3213 AssertDimension(values.size(), n_quadrature_points);
3214
3216 this->operator[](scalar).get_jump_in_function_values(fe_function, values);
3217}
3218
3219
3220
3221template <int dim, int spacedim>
3222template <class InputVector>
3223void
3225 const InputVector &fe_function,
3227 const
3228{
3229 AssertDimension(gradients.size(), n_quadrature_points);
3230
3232 this->operator[](scalar).get_jump_in_function_gradients(fe_function,
3233 gradients);
3234}
3235
3236
3237
3238template <int dim, int spacedim>
3239template <class InputVector>
3240void
3242 const InputVector &fe_function,
3244 const
3245{
3246 AssertDimension(hessians.size(), n_quadrature_points);
3247
3249 this->operator[](scalar).get_jump_in_function_hessians(fe_function, hessians);
3250}
3251
3252
3253
3254template <int dim, int spacedim>
3255template <class InputVector>
3256void
3258 const InputVector &fe_function,
3260 &third_derivatives) const
3261{
3262 AssertDimension(third_derivatives.size(), n_quadrature_points);
3263
3265 this->operator[](scalar).get_jump_in_function_third_derivatives(
3266 fe_function, third_derivatives);
3267}
3268
3269
3270
3271template <int dim, int spacedim>
3272template <class InputVector>
3273void
3275 const InputVector &fe_function,
3276 std::vector<typename InputVector::value_type> &values) const
3277{
3278 AssertDimension(values.size(), n_quadrature_points);
3279
3281 this->operator[](scalar).get_average_of_function_values(fe_function, values);
3282}
3283
3284
3285
3286template <int dim, int spacedim>
3287template <class InputVector>
3288void
3290 const InputVector &fe_function,
3292 const
3293{
3294 AssertDimension(gradients.size(), n_quadrature_points);
3295
3297 this->operator[](scalar).get_average_of_function_gradients(fe_function,
3298 gradients);
3299}
3300
3301
3302
3303template <int dim, int spacedim>
3304template <class InputVector>
3305void
3307 const InputVector &fe_function,
3309 const
3310{
3311 AssertDimension(hessians.size(), n_quadrature_points);
3312
3314 this->operator[](scalar).get_average_of_function_hessians(fe_function,
3315 hessians);
3316}
3317
3318
3319
3320/*------------ Inline functions: FEInterfaceValues------------*/
3321template <int dim, int spacedim>
3324 const FEValuesExtractors::Scalar &scalar) const
3325{
3326 const unsigned int n_components =
3327 (this->has_hp_capabilities() ? this->get_fe_collection().n_components() :
3328 this->get_fe().n_components());
3329 (void)n_components;
3330 AssertIndexRange(scalar.component, n_components);
3331 return FEInterfaceViews::Scalar<dim, spacedim>(*this, scalar.component);
3332}
3333
3334
3335
3336template <int dim, int spacedim>
3339 const FEValuesExtractors::Vector &vector) const
3340{
3341 const unsigned int n_components =
3342 (this->has_hp_capabilities() ? this->get_fe_collection().n_components() :
3343 this->get_fe().n_components());
3344 const unsigned int n_vectors =
3347 0);
3348 (void)n_components;
3349 (void)n_vectors;
3350 AssertIndexRange(vector.first_vector_component, n_vectors);
3352 vector.first_vector_component);
3353}
3354
3355
3356
3357namespace FEInterfaceViews
3358{
3359 template <int dim, int spacedim>
3361 const FEInterfaceValues<dim, spacedim> &fe_interface)
3362 : fe_interface(&fe_interface)
3363 {}
3364
3365
3366
3367 template <int dim, int spacedim>
3368 Scalar<dim, spacedim>::Scalar(
3369 const FEInterfaceValues<dim, spacedim> &fe_interface,
3370 const unsigned int component)
3371 : Base<dim, spacedim>(fe_interface)
3372 , extractor(component)
3373 {}
3374
3375
3376
3377 template <int dim, int spacedim>
3378 template <class InputVector, class OutputVector>
3379 void
3380 Base<dim, spacedim>::get_local_dof_values(
3381 const InputVector &dof_values,
3382 OutputVector &local_dof_values) const
3383 {
3384 const auto &interface_dof_indices =
3385 this->fe_interface->get_interface_dof_indices();
3386
3387 AssertDimension(interface_dof_indices.size(), local_dof_values.size());
3388
3389 for (const unsigned int i : this->fe_interface->dof_indices())
3390 local_dof_values[i] = dof_values(interface_dof_indices[i]);
3391 }
3392
3393
3394
3395 template <int dim, int spacedim>
3396 typename Scalar<dim, spacedim>::value_type
3397 Scalar<dim, spacedim>::value(const bool here_or_there,
3398 const unsigned int interface_dof_index,
3399 const unsigned int q_point) const
3400 {
3401 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3402
3403 if (here_or_there && dof_pair[0] != numbers::invalid_unsigned_int)
3404 return (*(this->fe_interface->fe_face_values))[extractor].value(
3405 dof_pair[0], q_point);
3406
3407 if (!here_or_there && dof_pair[1] != numbers::invalid_unsigned_int)
3408 return (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3409 dof_pair[1], q_point);
3410
3411 return 0.0;
3412 }
3413
3414
3415
3416 template <int dim, int spacedim>
3417 typename Scalar<dim, spacedim>::value_type
3418 Scalar<dim, spacedim>::jump_in_values(const unsigned int interface_dof_index,
3419 const unsigned int q_point) const
3420 {
3421 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3422
3423 value_type value = 0.0;
3424
3425 if (dof_pair[0] != numbers::invalid_unsigned_int)
3426 value +=
3427 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3428 q_point);
3429
3430 if (dof_pair[1] != numbers::invalid_unsigned_int)
3431 value -=
3432 (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3433 dof_pair[1], q_point);
3434
3435 return value;
3436 }
3437
3438
3439
3440 template <int dim, int spacedim>
3441 typename Scalar<dim, spacedim>::value_type
3442 Scalar<dim, spacedim>::average_of_values(
3443 const unsigned int interface_dof_index,
3444 const unsigned int q_point) const
3445 {
3446 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3447
3448 if (this->fe_interface->at_boundary())
3449 return (*(this->fe_interface->fe_face_values))[extractor].value(
3450 dof_pair[0], q_point);
3451
3452 value_type value = 0.0;
3453
3454 if (dof_pair[0] != numbers::invalid_unsigned_int)
3455 value +=
3456 0.5 *
3457 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3458 q_point);
3459
3460 if (dof_pair[1] != numbers::invalid_unsigned_int)
3461 value +=
3462 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3463 dof_pair[1], q_point);
3464
3465 return value;
3466 }
3467
3468
3469
3470 template <int dim, int spacedim>
3471 typename Scalar<dim, spacedim>::gradient_type
3472 Scalar<dim, spacedim>::average_of_gradients(
3473 const unsigned int interface_dof_index,
3474 const unsigned int q_point) const
3475 {
3476 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3477
3478 if (this->fe_interface->at_boundary())
3479 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3480 dof_pair[0], q_point);
3481
3482 gradient_type value;
3483
3484 if (dof_pair[0] != numbers::invalid_unsigned_int)
3485 value +=
3486 0.5 *
3487 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3488 q_point);
3489
3490 if (dof_pair[1] != numbers::invalid_unsigned_int)
3491 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3492 .gradient(dof_pair[1], q_point);
3493
3494 return value;
3495 }
3496
3497
3498
3499 template <int dim, int spacedim>
3500 typename Scalar<dim, spacedim>::gradient_type
3501 Scalar<dim, spacedim>::jump_in_gradients(
3502 const unsigned int interface_dof_index,
3503 const unsigned int q_point) const
3504 {
3505 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3506
3507 if (this->fe_interface->at_boundary())
3508 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3509 dof_pair[0], q_point);
3510
3511 gradient_type value;
3512
3513 if (dof_pair[0] != numbers::invalid_unsigned_int)
3514 value +=
3515 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3516 q_point);
3517
3518 if (dof_pair[1] != numbers::invalid_unsigned_int)
3519 value -=
3520 (*(this->fe_interface->fe_face_values_neighbor))[extractor].gradient(
3521 dof_pair[1], q_point);
3522
3523 return value;
3524 }
3525
3526
3527
3528 template <int dim, int spacedim>
3529 typename Scalar<dim, spacedim>::hessian_type
3530 Scalar<dim, spacedim>::average_of_hessians(
3531 const unsigned int interface_dof_index,
3532 const unsigned int q_point) const
3533 {
3534 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3535
3536 if (this->fe_interface->at_boundary())
3537 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
3538 dof_pair[0], q_point);
3539
3540 hessian_type value;
3541
3542 if (dof_pair[0] != numbers::invalid_unsigned_int)
3543 value +=
3544 0.5 *
3545 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
3546 q_point);
3547
3548 if (dof_pair[1] != numbers::invalid_unsigned_int)
3549 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3550 .hessian(dof_pair[1], q_point);
3551
3552 return value;
3553 }
3554
3555
3556
3557 template <int dim, int spacedim>
3558 typename Scalar<dim, spacedim>::third_derivative_type
3559 Scalar<dim, spacedim>::jump_in_third_derivatives(
3560 const unsigned int interface_dof_index,
3561 const unsigned int q_point) const
3562 {
3563 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3564
3565 if (this->fe_interface->at_boundary())
3566 return (*(this->fe_interface->fe_face_values))[extractor]
3567 .third_derivative(dof_pair[0], q_point);
3568
3569 third_derivative_type value;
3570
3571 if (dof_pair[0] != numbers::invalid_unsigned_int)
3572 value +=
3573 (*(this->fe_interface->fe_face_values))[extractor].third_derivative(
3574 dof_pair[0], q_point);
3575
3576 if (dof_pair[1] != numbers::invalid_unsigned_int)
3577 value -= (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3578 .third_derivative(dof_pair[1], q_point);
3579
3580 return value;
3581 }
3582
3583
3584
3585 template <int dim, int spacedim>
3586 typename Scalar<dim, spacedim>::hessian_type
3587 Scalar<dim, spacedim>::jump_in_hessians(
3588 const unsigned int interface_dof_index,
3589 const unsigned int q_point) const
3590 {
3591 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3592
3593 if (this->fe_interface->at_boundary())
3594 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
3595 dof_pair[0], q_point);
3596
3597 hessian_type value;
3598
3599 if (dof_pair[0] != numbers::invalid_unsigned_int)
3600 value +=
3601 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
3602 q_point);
3603
3604 if (dof_pair[1] != numbers::invalid_unsigned_int)
3605 value -=
3606 (*(this->fe_interface->fe_face_values_neighbor))[extractor].hessian(
3607 dof_pair[1], q_point);
3608
3609 return value;
3610 }
3611
3612
3613
3614 template <int dim, int spacedim>
3615 template <class InputVector>
3616 void
3617 Scalar<dim, spacedim>::get_function_values_from_local_dof_values(
3618 const bool here_or_there,
3619 const InputVector &local_dof_values,
3620 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3621 const
3622 {
3623 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
3624
3625 for (const auto dof_index : this->fe_interface->dof_indices())
3626 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3627 {
3628 if (dof_index == 0)
3629 values[q_index] = 0.;
3630
3631 values[q_index] += local_dof_values[dof_index] *
3632 value(here_or_there, dof_index, q_index);
3633 }
3634 }
3635
3636
3637
3638 template <int dim, int spacedim>
3639 template <class InputVector>
3640 void
3641 Scalar<dim, spacedim>::get_function_values(
3642 const bool here_or_there,
3643 const InputVector &fe_function,
3644 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3645 const
3646 {
3647 std::vector<typename InputVector::value_type> local_dof_values(
3648 this->fe_interface->n_current_interface_dofs());
3649 this->get_local_dof_values(fe_function, local_dof_values);
3650
3651 get_function_values_from_local_dof_values(here_or_there,
3652 local_dof_values,
3653 values);
3654 }
3655
3656
3657
3658 template <int dim, int spacedim>
3659 template <class InputVector>
3660 void
3661 Scalar<dim, spacedim>::get_jump_in_function_values_from_local_dof_values(
3662 const InputVector &local_dof_values,
3663 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3664 const
3665 {
3666 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
3667
3668 for (const auto dof_index : this->fe_interface->dof_indices())
3669 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3670 {
3671 if (dof_index == 0)
3672 values[q_index] = 0.;
3673
3674 values[q_index] +=
3675 local_dof_values[dof_index] * jump_in_values(dof_index, q_index);
3676 }
3677 }
3678
3679
3680
3681 template <int dim, int spacedim>
3682 template <class InputVector>
3683 void
3684 Scalar<dim, spacedim>::get_jump_in_function_values(
3685 const InputVector &fe_function,
3686 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3687 const
3688 {
3689 std::vector<typename InputVector::value_type> local_dof_values(
3690 this->fe_interface->n_current_interface_dofs());
3691 this->get_local_dof_values(fe_function, local_dof_values);
3692
3693 get_jump_in_function_values_from_local_dof_values(local_dof_values, values);
3694 }
3695
3696
3697
3698 template <int dim, int spacedim>
3699 template <class InputVector>
3700 void
3701 Scalar<dim, spacedim>::get_jump_in_function_gradients_from_local_dof_values(
3702 const InputVector &local_dof_values,
3703 std::vector<solution_gradient_type<typename InputVector::value_type>>
3704 &gradients) const
3705 {
3706 AssertDimension(gradients.size(), this->fe_interface->n_quadrature_points);
3707
3708 for (const auto dof_index : this->fe_interface->dof_indices())
3709 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3710 {
3711 if (dof_index == 0)
3712 gradients[q_index] = 0.;
3713
3714 gradients[q_index] +=
3715 local_dof_values[dof_index] * jump_in_gradients(dof_index, q_index);
3716 }
3717 }
3718
3719
3720
3721 template <int dim, int spacedim>
3722 template <class InputVector>
3723 void
3724 Scalar<dim, spacedim>::get_jump_in_function_gradients(
3725 const InputVector &fe_function,
3726 std::vector<solution_gradient_type<typename InputVector::value_type>>
3727 &gradients) const
3728 {
3729 std::vector<typename InputVector::value_type> local_dof_values(
3730 this->fe_interface->n_current_interface_dofs());
3731 this->get_local_dof_values(fe_function, local_dof_values);
3732
3733 get_jump_in_function_gradients_from_local_dof_values(local_dof_values,
3734 gradients);
3735 }
3736
3737
3738
3739 template <int dim, int spacedim>
3740 template <class InputVector>
3741 void
3742 Scalar<dim, spacedim>::get_average_of_function_values_from_local_dof_values(
3743 const InputVector &local_dof_values,
3744 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3745 const
3746 {
3747 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
3748
3749 for (const auto dof_index : this->fe_interface->dof_indices())
3750 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3751 {
3752 if (dof_index == 0)
3753 values[q_index] = 0.;
3754
3755 values[q_index] +=
3756 local_dof_values[dof_index] * average_of_values(dof_index, q_index);
3757 }
3758 }
3759
3760
3761
3762 template <int dim, int spacedim>
3763 template <class InputVector>
3764 void
3765 Scalar<dim, spacedim>::get_average_of_function_values(
3766 const InputVector &fe_function,
3767 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3768 const
3769 {
3770 std::vector<typename InputVector::value_type> local_dof_values(
3771 this->fe_interface->n_current_interface_dofs());
3772 this->get_local_dof_values(fe_function, local_dof_values);
3773
3774 get_average_of_function_values_from_local_dof_values(local_dof_values,
3775 values);
3776 }
3777
3778
3779
3780 template <int dim, int spacedim>
3781 template <class InputVector>
3782 void
3783 Scalar<dim, spacedim>::
3784 get_average_of_function_gradients_from_local_dof_values(
3785 const InputVector &local_dof_values,
3786 std::vector<solution_gradient_type<typename InputVector::value_type>>
3787 &gradients) const
3788 {
3789 AssertDimension(gradients.size(), this->fe_interface->n_quadrature_points);
3790
3791 for (const auto dof_index : this->fe_interface->dof_indices())
3792 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3793 {
3794 if (dof_index == 0)
3795 gradients[q_index] = 0.;
3796
3797 gradients[q_index] += local_dof_values[dof_index] *
3798 average_of_gradients(dof_index, q_index);
3799 }
3800 }
3801
3802
3803
3804 template <int dim, int spacedim>
3805 template <class InputVector>
3806 void
3807 Scalar<dim, spacedim>::get_average_of_function_gradients(
3808 const InputVector &fe_function,
3809 std::vector<solution_gradient_type<typename InputVector::value_type>>
3810 &gradients) const
3811 {
3812 std::vector<typename InputVector::value_type> local_dof_values(
3813 this->fe_interface->n_current_interface_dofs());
3814 this->get_local_dof_values(fe_function, local_dof_values);
3815
3816 get_average_of_function_gradients_from_local_dof_values(local_dof_values,
3817 gradients);
3818 }
3819
3820
3821
3822 template <int dim, int spacedim>
3823 template <class InputVector>
3824 void
3825 Scalar<dim, spacedim>::get_jump_in_function_hessians_from_local_dof_values(
3826 const InputVector &local_dof_values,
3827 std::vector<solution_hessian_type<typename InputVector::value_type>>
3828 &hessians) const
3829 {
3830 AssertDimension(hessians.size(), this->fe_interface->n_quadrature_points);
3831
3832 for (const auto dof_index : this->fe_interface->dof_indices())
3833 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3834 {
3835 if (dof_index == 0)
3836 hessians[q_index] = 0.;
3837
3838 hessians[q_index] +=
3839 local_dof_values[dof_index] * jump_in_hessians(dof_index, q_index);
3840 }
3841 }
3842
3843
3844
3845 template <int dim, int spacedim>
3846 template <class InputVector>
3847 void
3848 Scalar<dim, spacedim>::get_jump_in_function_hessians(
3849 const InputVector &fe_function,
3850 std::vector<solution_hessian_type<typename InputVector::value_type>>
3851 &hessians) const
3852 {
3853 std::vector<typename InputVector::value_type> local_dof_values(
3854 this->fe_interface->n_current_interface_dofs());
3855 this->get_local_dof_values(fe_function, local_dof_values);
3856
3857 get_jump_in_function_hessians_from_local_dof_values(local_dof_values,
3858 hessians);
3859 }
3860
3861
3862
3863 template <int dim, int spacedim>
3864 template <class InputVector>
3865 void
3866 Scalar<dim, spacedim>::get_average_of_function_hessians_from_local_dof_values(
3867 const InputVector &local_dof_values,
3868 std::vector<solution_hessian_type<typename InputVector::value_type>>
3869 &hessians) const
3870 {
3871 AssertDimension(hessians.size(), this->fe_interface->n_quadrature_points);
3872
3873 for (const unsigned int dof_index : this->fe_interface->dof_indices())
3874 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3875 {
3876 if (dof_index == 0)
3877 hessians[q_index] = 0.;
3878
3879 hessians[q_index] += local_dof_values[dof_index] *
3880 average_of_hessians(dof_index, q_index);
3881 }
3882 }
3883
3884
3885
3886 template <int dim, int spacedim>
3887 template <class InputVector>
3888 void
3889 Scalar<dim, spacedim>::get_average_of_function_hessians(
3890 const InputVector &fe_function,
3891 std::vector<solution_hessian_type<typename InputVector::value_type>>
3892 &hessians) const
3893 {
3894 std::vector<typename InputVector::value_type> local_dof_values(
3895 this->fe_interface->n_current_interface_dofs());
3896 this->get_local_dof_values(fe_function, local_dof_values);
3897
3898 get_average_of_function_hessians_from_local_dof_values(local_dof_values,
3899 hessians);
3900 }
3901
3902
3903
3904 template <int dim, int spacedim>
3905 template <class InputVector>
3906 void
3907 Scalar<dim, spacedim>::
3908 get_jump_in_function_third_derivatives_from_local_dof_values(
3909 const InputVector &local_dof_values,
3910 std::vector<
3911 solution_third_derivative_type<typename InputVector::value_type>>
3912 &third_derivatives) const
3913 {
3914 AssertDimension(third_derivatives.size(),
3915 this->fe_interface->n_quadrature_points);
3916
3917 for (const unsigned int dof_index : this->fe_interface->dof_indices())
3918 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3919 {
3920 if (dof_index == 0)
3921 third_derivatives[q_index] = 0.;
3922
3923 third_derivatives[q_index] +=
3924 local_dof_values[dof_index] *
3925 jump_in_third_derivatives(dof_index, q_index);
3926 }
3927 }
3928
3929
3930
3931 template <int dim, int spacedim>
3932 template <class InputVector>
3933 void
3934 Scalar<dim, spacedim>::get_jump_in_function_third_derivatives(
3935 const InputVector &fe_function,
3936 std::vector<
3937 solution_third_derivative_type<typename InputVector::value_type>>
3938 &third_derivatives) const
3939 {
3940 std::vector<typename InputVector::value_type> local_dof_values(
3941 this->fe_interface->n_current_interface_dofs());
3942 this->get_local_dof_values(fe_function, local_dof_values);
3943
3944 get_jump_in_function_third_derivatives_from_local_dof_values(
3945 local_dof_values, third_derivatives);
3946 }
3947
3948
3949
3950 template <int dim, int spacedim>
3952 const FEInterfaceValues<dim, spacedim> &fe_interface,
3953 const unsigned int first_vector_component)
3954 : Base<dim, spacedim>(fe_interface)
3955 , extractor(first_vector_component)
3956 {}
3957
3958
3959
3960 template <int dim, int spacedim>
3962 Vector<dim, spacedim>::value(const bool here_or_there,
3963 const unsigned int interface_dof_index,
3964 const unsigned int q_point) const
3965 {
3966 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3967
3968 if (here_or_there && dof_pair[0] != numbers::invalid_unsigned_int)
3969 return (*(this->fe_interface->fe_face_values))[extractor].value(
3970 dof_pair[0], q_point);
3971
3972 if (!here_or_there && dof_pair[1] != numbers::invalid_unsigned_int)
3973 return (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3974 dof_pair[1], q_point);
3975
3976 return value_type();
3977 }
3978
3979
3980
3981 template <int dim, int spacedim>
3983 Vector<dim, spacedim>::jump_in_values(const unsigned int interface_dof_index,
3984 const unsigned int q_point) const
3985 {
3986 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3987
3988 value_type value;
3989
3990 if (dof_pair[0] != numbers::invalid_unsigned_int)
3991 value +=
3992 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3993 q_point);
3994
3995 if (dof_pair[1] != numbers::invalid_unsigned_int)
3996 value -=
3997 (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3998 dof_pair[1], q_point);
3999
4000 return value;
4001 }
4002
4003
4004
4005 template <int dim, int spacedim>
4008 const unsigned int interface_dof_index,
4009 const unsigned int q_point) const
4010 {
4011 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
4012
4013 if (this->fe_interface->at_boundary())
4014 return (*(this->fe_interface->fe_face_values))[extractor].value(
4015 dof_pair[0], q_point);
4016
4017 value_type value;
4018
4019 if (dof_pair[0] != numbers::invalid_unsigned_int)
4020 value +=
4021 0.5 *
4022 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
4023 q_point);
4024
4025 if (dof_pair[1] != numbers::invalid_unsigned_int)
4026 value +=
4027 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
4028 dof_pair[1], q_point);
4029
4030 return value;
4031 }
4032
4033
4034
4035 template <int dim, int spacedim>
4038 const unsigned int interface_dof_index,
4039 const unsigned int q_point) const
4040 {
4041 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
4042
4043 if (this->fe_interface->at_boundary())
4044 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
4045 dof_pair[0], q_point);
4046
4047 gradient_type value;
4048
4049 if (dof_pair[0] != numbers::invalid_unsigned_int)
4050 value +=
4051 0.5 *
4052 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
4053 q_point);
4054
4055 if (dof_pair[1] != numbers::invalid_unsigned_int)
4056 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
4057 .gradient(dof_pair[1], q_point);
4058
4059 return value;
4060 }
4061
4062
4063
4064 template <int dim, int spacedim>
4067 const unsigned int interface_dof_index,
4068 const unsigned int q_point) const
4069 {
4070 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
4071
4072 if (this->fe_interface->at_boundary())
4073 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
4074 dof_pair[0], q_point);
4075
4076 gradient_type value;
4077
4078 if (dof_pair[0] != numbers::invalid_unsigned_int)
4079 value +=
4080 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
4081 q_point);
4082
4083 if (dof_pair[1] != numbers::invalid_unsigned_int)
4084 value -=
4085 (*(this->fe_interface->fe_face_values_neighbor))[extractor].gradient(
4086 dof_pair[1], q_point);
4087
4088 return value;
4089 }
4090
4091
4092
4093 template <int dim, int spacedim>
4095 Vector<dim, spacedim>::jump_gradient(const unsigned int interface_dof_index,
4096 const unsigned int q_point) const
4097 {
4098 return jump_in_gradients(interface_dof_index, q_point);
4099 }
4100
4101
4102
4103 template <int dim, int spacedim>
4106 const unsigned int interface_dof_index,
4107 const unsigned int q_point) const
4108 {
4109 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
4110
4111 if (this->fe_interface->at_boundary())
4112 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
4113 dof_pair[0], q_point);
4114
4115 hessian_type value;
4116
4117 if (dof_pair[0] != numbers::invalid_unsigned_int)
4118 value +=
4119 0.5 *
4120 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
4121 q_point);
4122
4123 if (dof_pair[1] != numbers::invalid_unsigned_int)
4124 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
4125 .hessian(dof_pair[1], q_point);
4126
4127 return value;
4128 }
4129
4130
4131
4132 template <int dim, int spacedim>
4134 Vector<dim, spacedim>::average_hessian(const unsigned int interface_dof_index,
4135 const unsigned int q_point) const
4136 {
4137 return average_of_hessians(interface_dof_index, q_point);
4138 }
4139
4140
4141
4142 template <int dim, int spacedim>
4145 const unsigned int interface_dof_index,
4146 const unsigned int q_point) const
4147 {
4148 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
4149
4150 if (this->fe_interface->at_boundary())
4151 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
4152 dof_pair[0], q_point);
4153
4154 hessian_type value;
4155
4156 if (dof_pair[0] != numbers::invalid_unsigned_int)
4157 value +=
4158 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
4159 q_point);
4160
4161 if (dof_pair[1] != numbers::invalid_unsigned_int)
4162 value -=
4163 (*(this->fe_interface->fe_face_values_neighbor))[extractor].hessian(
4164 dof_pair[1], q_point);
4165
4166 return value;
4167 }
4168
4169
4170
4171 template <int dim, int spacedim>
4173 Vector<dim, spacedim>::jump_hessian(const unsigned int interface_dof_index,
4174 const unsigned int q_point) const
4175 {
4176 return jump_in_hessians(interface_dof_index, q_point);
4177 }
4178
4179
4180
4181 template <int dim, int spacedim>
4184 const unsigned int interface_dof_index,
4185 const unsigned int q_point) const
4186 {
4187 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
4188
4189 if (this->fe_interface->at_boundary())
4190 return (*(this->fe_interface->fe_face_values))[extractor]
4191 .third_derivative(dof_pair[0], q_point);
4192
4193 third_derivative_type value;
4194
4195 if (dof_pair[0] != numbers::invalid_unsigned_int)
4196 value +=
4197 (*(this->fe_interface->fe_face_values))[extractor].third_derivative(
4198 dof_pair[0], q_point);
4199
4200 if (dof_pair[1] != numbers::invalid_unsigned_int)
4201 value -= (*(this->fe_interface->fe_face_values_neighbor))[extractor]
4202 .third_derivative(dof_pair[1], q_point);
4203
4204 return value;
4205 }
4206
4207
4208
4209 template <int dim, int spacedim>
4210 template <class InputVector>
4211 void
4213 const bool here_or_there,
4214 const InputVector &local_dof_values,
4215 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4216 const
4217 {
4218 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
4219
4220 for (const auto dof_index : this->fe_interface->dof_indices())
4221 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4222 {
4223 if (dof_index == 0)
4224 values[q_index] = 0.;
4225
4226 values[q_index] += local_dof_values[dof_index] *
4227 value(here_or_there, dof_index, q_index);
4228 }
4229 }
4230
4231
4232
4233 template <int dim, int spacedim>
4234 template <class InputVector>
4235 void
4237 const bool here_or_there,
4238 const InputVector &fe_function,
4239 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4240 const
4241 {
4242 std::vector<typename InputVector::value_type> local_dof_values(
4243 this->fe_interface->n_current_interface_dofs());
4244 this->get_local_dof_values(fe_function, local_dof_values);
4245
4246 get_function_values_from_local_dof_values(here_or_there,
4247 local_dof_values,
4248 values);
4249 }
4250
4251
4252
4253 template <int dim, int spacedim>
4254 template <class InputVector>
4255 void
4257 const InputVector &local_dof_values,
4258 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4259 const
4260 {
4261 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
4262
4263 for (const auto dof_index : this->fe_interface->dof_indices())
4264 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4265 {
4266 if (dof_index == 0)
4267 values[q_index] = 0.;
4268
4269 values[q_index] +=
4270 local_dof_values[dof_index] * jump_in_values(dof_index, q_index);
4271 }
4272 }
4273
4274
4275
4276 template <int dim, int spacedim>
4277 template <class InputVector>
4278 void
4280 const InputVector &fe_function,
4281 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4282 const
4283 {
4284 std::vector<typename InputVector::value_type> local_dof_values(
4285 this->fe_interface->n_current_interface_dofs());
4286 this->get_local_dof_values(fe_function, local_dof_values);
4287
4288 get_jump_in_function_values_from_local_dof_values(local_dof_values, values);
4289 }
4290
4291
4292
4293 template <int dim, int spacedim>
4294 template <class InputVector>
4295 void
4297 const InputVector &local_dof_values,
4298 std::vector<solution_gradient_type<typename InputVector::value_type>>
4299 &gradients) const
4300 {
4301 AssertDimension(gradients.size(), this->fe_interface->n_quadrature_points);
4302
4303 for (const auto dof_index : this->fe_interface->dof_indices())
4304 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4305 {
4306 if (dof_index == 0)
4307 gradients[q_index] = 0.;
4308
4309 gradients[q_index] +=
4310 local_dof_values[dof_index] * jump_in_gradients(dof_index, q_index);
4311 }
4312 }
4313
4314
4315
4316 template <int dim, int spacedim>
4317 template <class InputVector>
4318 void
4320 const InputVector &fe_function,
4321 std::vector<solution_gradient_type<typename InputVector::value_type>>
4322 &gradients) const
4323 {
4324 std::vector<typename InputVector::value_type> local_dof_values(
4325 this->fe_interface->n_current_interface_dofs());
4326 this->get_local_dof_values(fe_function, local_dof_values);
4327
4328 get_jump_in_function_gradients_from_local_dof_values(local_dof_values,
4329 gradients);
4330 }
4331
4332
4333
4334 template <int dim, int spacedim>
4335 template <class InputVector>
4336 void
4338 const InputVector &local_dof_values,
4339 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4340 const
4341 {
4342 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
4343
4344 for (const auto dof_index : this->fe_interface->dof_indices())
4345 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4346 {
4347 if (dof_index == 0)
4348 values[q_index] = 0.;
4349
4350 values[q_index] +=
4351 local_dof_values[dof_index] * average_of_values(dof_index, q_index);
4352 }
4353 }
4354
4355
4356
4357 template <int dim, int spacedim>
4358 template <class InputVector>
4359 void
4361 const InputVector &fe_function,
4362 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4363 const
4364 {
4365 std::vector<typename InputVector::value_type> local_dof_values(
4366 this->fe_interface->n_current_interface_dofs());
4367 this->get_local_dof_values(fe_function, local_dof_values);
4368
4369 get_average_of_function_values_from_local_dof_values(local_dof_values,
4370 values);
4371 }
4372
4373
4374
4375 template <int dim, int spacedim>
4376 template <class InputVector>
4377 void
4380 const InputVector &local_dof_values,
4381 std::vector<solution_gradient_type<typename InputVector::value_type>>
4382 &gradients) const
4383 {
4384 AssertDimension(gradients.size(), this->fe_interface->n_quadrature_points);
4385
4386 for (const auto dof_index : this->fe_interface->dof_indices())
4387 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4388 {
4389 if (dof_index == 0)
4390 gradients[q_index] = 0.;
4391
4392 gradients[q_index] += local_dof_values[dof_index] *
4393 average_of_gradients(dof_index, q_index);
4394 }
4395 }
4396
4397
4398
4399 template <int dim, int spacedim>
4400 template <class InputVector>
4401 void
4403 const InputVector &fe_function,
4404 std::vector<solution_gradient_type<typename InputVector::value_type>>
4405 &gradients) const
4406 {
4407 std::vector<typename InputVector::value_type> local_dof_values(
4408 this->fe_interface->n_current_interface_dofs());
4409 this->get_local_dof_values(fe_function, local_dof_values);
4410
4411 get_average_of_function_gradients_from_local_dof_values(local_dof_values,
4412 gradients);
4413 }
4414
4415
4416
4417 template <int dim, int spacedim>
4418 template <class InputVector>
4419 void
4421 const InputVector &local_dof_values,
4422 std::vector<solution_hessian_type<typename InputVector::value_type>>
4423 &hessians) const
4424 {
4425 AssertDimension(hessians.size(), this->fe_interface->n_quadrature_points);
4426
4427 for (const auto dof_index : this->fe_interface->dof_indices())
4428 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4429 {
4430 if (dof_index == 0)
4431 hessians[q_index] = 0.;
4432
4433 hessians[q_index] +=
4434 local_dof_values[dof_index] * jump_in_hessians(dof_index, q_index);
4435 }
4436 }
4437
4438
4439
4440 template <int dim, int spacedim>
4441 template <class InputVector>
4442 void
4444 const InputVector &fe_function,
4445 std::vector<solution_hessian_type<typename InputVector::value_type>>
4446 &hessians) const
4447 {
4448 std::vector<typename InputVector::value_type> local_dof_values(
4449 this->fe_interface->n_current_interface_dofs());
4450 this->get_local_dof_values(fe_function, local_dof_values);
4451
4452 get_jump_in_function_hessians_from_local_dof_values(local_dof_values,
4453 hessians);
4454 }
4455
4456
4457
4458 template <int dim, int spacedim>
4459 template <class InputVector>
4460 void
4462 const InputVector &local_dof_values,
4463 std::vector<solution_hessian_type<typename InputVector::value_type>>
4464 &hessians) const
4465 {
4466 AssertDimension(hessians.size(), this->fe_interface->n_quadrature_points);
4467
4468 for (const unsigned int dof_index : this->fe_interface->dof_indices())
4469 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4470 {
4471 if (dof_index == 0)
4472 hessians[q_index] = 0.;
4473
4474 hessians[q_index] += local_dof_values[dof_index] *
4475 average_of_hessians(dof_index, q_index);
4476 }
4477 }
4478
4479
4480
4481 template <int dim, int spacedim>
4482 template <class InputVector>
4483 void
4485 const InputVector &fe_function,
4486 std::vector<solution_hessian_type<typename InputVector::value_type>>
4487 &hessians) const
4488 {
4489 std::vector<typename InputVector::value_type> local_dof_values(
4490 this->fe_interface->n_current_interface_dofs());
4491 this->get_local_dof_values(fe_function, local_dof_values);
4492
4493 get_average_of_function_hessians_from_local_dof_values(local_dof_values,
4494 hessians);
4495 }
4496
4497
4498
4499 template <int dim, int spacedim>
4500 template <class InputVector>
4501 void
4504 const InputVector &local_dof_values,
4505 std::vector<
4506 solution_third_derivative_type<typename InputVector::value_type>>
4507 &third_derivatives) const
4508 {
4509 AssertDimension(third_derivatives.size(),
4510 this->fe_interface->n_quadrature_points);
4511
4512 for (const unsigned int dof_index : this->fe_interface->dof_indices())
4513 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4514 {
4515 if (dof_index == 0)
4516 third_derivatives[q_index] = 0.;
4517
4518 third_derivatives[q_index] +=
4519 local_dof_values[dof_index] *
4520 jump_in_third_derivatives(dof_index, q_index);
4521 }
4522 }
4523
4524
4525
4526 template <int dim, int spacedim>
4527 template <class InputVector>
4528 void
4530 const InputVector &fe_function,
4531 std::vector<
4532 solution_third_derivative_type<typename InputVector::value_type>>
4533 &third_derivatives) const
4534 {
4535 std::vector<typename InputVector::value_type> local_dof_values(
4536 this->fe_interface->n_current_interface_dofs());
4537 this->get_local_dof_values(fe_function, local_dof_values);
4538
4539 get_jump_in_function_third_derivatives_from_local_dof_values(
4540 local_dof_values, third_derivatives);
4541 }
4542} // namespace FEInterfaceViews
4543
4544#endif // DOXYGEN
4545
4547
4548#endif
friend class FEInterfaceViews::Scalar
std::unique_ptr< FEFaceValues< dim, spacedim > > internal_fe_face_values
UpdateFlags get_update_flags() const
const hp::MappingCollection< dim, spacedim > & get_mapping_collection() const
Tensor< 1, spacedim > average_of_shape_gradients(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
void get_average_of_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &gradients) const
void get_jump_in_function_third_derivatives(const InputVector &fe_function, std::vector< Tensor< 3, spacedim, typename InputVector::value_type > > &third_derivatives) const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > dof_indices() const
FEFaceValuesBase< dim, spacedim > * fe_face_values_neighbor
friend class FEInterfaceViews::Vector
bool at_boundary() const
const std::vector< double > & get_JxW_values() const
void get_average_of_function_hessians(const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type > > &hessians) const
FEInterfaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature, const UpdateFlags update_flags)
void reinit(const CellIteratorType &cell, const unsigned int face_no, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
const hp::FECollection< dim, spacedim > & get_fe_collection() const
Tensor< 3, spacedim > jump_in_shape_3rd_derivatives(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
FEInterfaceValues(const hp::MappingCollection< dim, spacedim > &mapping_collection, const hp::FECollection< dim, spacedim > &fe_collection, const hp::QCollection< dim - 1 > &quadrature_collection, const UpdateFlags update_flags)
std::unique_ptr< hp::FESubfaceValues< dim, spacedim > > internal_hp_fe_subface_values
const unsigned int n_quadrature_points
bool has_hp_capabilities() const
void reinit(const CellIteratorType &cell, const unsigned int face_no, const unsigned int sub_face_no, const CellNeighborIteratorType &cell_neighbor, const unsigned int face_no_neighbor, const unsigned int sub_face_no_neighbor, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int, const unsigned int fe_index_neighbor=numbers::invalid_unsigned_int)
FEInterfaceViews::Scalar< dim, spacedim > operator[](const FEValuesExtractors::Scalar &scalar) const
unsigned n_current_interface_dofs() const
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors() const
const FiniteElement< dim, spacedim > & get_fe() const
double average_of_shape_values(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
const std::vector< Point< spacedim > > & get_quadrature_points() const
std::unique_ptr< FESubfaceValues< dim, spacedim > > internal_fe_subface_values
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
Tensor< 2, spacedim > average_of_shape_hessians(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
std::unique_ptr< hp::FESubfaceValues< dim, spacedim > > internal_hp_fe_subface_values_neighbor
FEInterfaceValues(const hp::FECollection< dim, spacedim > &fe_collection, const hp::QCollection< dim - 1 > &quadrature_collection, const UpdateFlags update_flags)
const hp::QCollection< dim - 1 > & get_quadrature_collection() 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
Tensor< 1, spacedim > normal_vector(const unsigned int q_point_index) const
std::unique_ptr< FESubfaceValues< dim, spacedim > > internal_fe_subface_values_neighbor
std::vector< std::array< unsigned int, 2 > > dofmap
Tensor< 1, spacedim > shape_grad(const bool here_or_there, const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
Triangulation< dim, spacedim >::cell_iterator get_cell(const unsigned int cell_index) const
void get_jump_in_function_hessians(const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type > > &hessians) const
double JxW(const unsigned int quadrature_point) const
Tensor< 1, spacedim > normal(const unsigned int q_point_index) const
void get_average_of_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
Tensor< 1, spacedim > jump_in_shape_gradients(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
Tensor< 2, spacedim > jump_in_shape_hessians(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
std::unique_ptr< hp::FEFaceValues< dim, spacedim > > internal_hp_fe_face_values_neighbor
std_cxx20::ranges::iota_view< unsigned int, unsigned int > quadrature_point_indices() const
double jump_in_shape_values(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
unsigned int get_face_number(const unsigned int cell_index) const
std::vector< types::global_dof_index > get_interface_dof_indices() const
FEInterfaceValues(const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature, const UpdateFlags update_flags)
std::array< unsigned int, 2 > interface_dof_to_dof_indices(const unsigned int interface_dof_index) const
std::unique_ptr< hp::FEFaceValues< dim, spacedim > > internal_hp_fe_face_values
std::unique_ptr< FEFaceValues< dim, spacedim > > internal_fe_face_values_neighbor
void get_jump_in_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
void get_jump_in_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &gradients) const
const Mapping< dim, spacedim > & get_mapping() const
const Point< spacedim > & quadrature_point(const unsigned int q_point) const
const Quadrature< dim - 1 > & get_quadrature() const
void get_local_dof_values(const InputVector &dof_values, OutputVector &local_dof_values) const
Base(const FEInterfaceValues< dim, spacedim > &fe_interface)
const FEInterfaceValues< dim, spacedim > * fe_interface
void get_jump_in_function_third_derivatives(const InputVector &fe_function, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &third_derivatives) const
void get_jump_in_function_gradients_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
typename FEValuesViews::Scalar< dim, spacedim >::third_derivative_type third_derivative_type
gradient_type jump_in_gradients(const unsigned int interface_dof_index, const unsigned int q_point) const
void get_jump_in_function_hessians(const InputVector &fe_function, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
hessian_type jump_in_hessians(const unsigned int interface_dof_index, const unsigned int q_point) const
const FEValuesExtractors::Scalar extractor
void get_jump_in_function_third_derivatives_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &third_derivatives) const
typename ProductType< Number, gradient_type >::type solution_gradient_type
hessian_type average_of_hessians(const unsigned int interface_dof_index, const unsigned int q_point) const
void get_jump_in_function_gradients(const InputVector &fe_function, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
value_type value(const bool here_or_there, const unsigned int interface_dof_index, const unsigned int q_point) const
void get_average_of_function_hessians_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
void get_average_of_function_values_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
typename ProductType< Number, third_derivative_type >::type solution_third_derivative_type
void get_average_of_function_values(const InputVector &fe_function, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
typename ProductType< Number, value_type >::type solution_value_type
void get_average_of_function_gradients(const InputVector &fe_function, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
void get_jump_in_function_values_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
void get_jump_in_function_values(const InputVector &fe_function, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
Scalar(const FEInterfaceValues< dim, spacedim > &fe_interface, const unsigned int component)
third_derivative_type jump_in_third_derivatives(const unsigned int interface_dof_index, const unsigned int q_point) const
gradient_type average_of_gradients(const unsigned int interface_dof_index, const unsigned int q_point) const
void get_average_of_function_hessians(const InputVector &fe_function, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
typename ProductType< Number, hessian_type >::type solution_hessian_type
void get_average_of_function_gradients_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
value_type average_of_values(const unsigned int interface_dof_index, const unsigned int q_point) const
void get_jump_in_function_hessians_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
value_type jump_in_values(const unsigned int interface_dof_index, const unsigned int q_point) const
typename FEValuesViews::Scalar< dim, spacedim >::hessian_type hessian_type
void get_function_values_from_local_dof_values(const bool here_or_there, const InputVector &local_dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
void get_function_values(const bool here_or_there, const InputVector &fe_function, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
typename FEValuesViews::Scalar< dim, spacedim >::gradient_type gradient_type
const FEValuesExtractors::Vector extractor
void get_average_of_function_gradients_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
void get_average_of_function_hessians(const InputVector &fe_function, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
void get_jump_in_function_values(const InputVector &fe_function, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
typename FEValuesViews::Vector< dim, spacedim >::value_type value_type
void get_average_of_function_values_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
void get_function_values_from_local_dof_values(const bool here_or_there, const InputVector &local_dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
void get_average_of_function_hessians_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
hessian_type average_hessian(const unsigned int interface_dof_index, const unsigned int q_point) const
typename FEValuesViews::Vector< dim, spacedim >::third_derivative_type third_derivative_type
Vector(const FEInterfaceValues< dim, spacedim > &fe_interface, const unsigned int first_vector_component)
void get_average_of_function_gradients(const InputVector &fe_function, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
void get_function_values(const bool here_or_there, const InputVector &fe_function, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
hessian_type jump_hessian(const unsigned int interface_dof_index, const unsigned int q_point) const
typename ProductType< Number, value_type >::type solution_value_type
typename FEValuesViews::Vector< dim, spacedim >::hessian_type hessian_type
void get_average_of_function_values(const InputVector &fe_function, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
value_type value(const bool here_or_there, const unsigned int interface_dof_index, const unsigned int q_point) const
typename ProductType< Number, hessian_type >::type solution_hessian_type
void get_jump_in_function_third_derivatives_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &third_derivatives) const
void get_jump_in_function_gradients(const InputVector &fe_function, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
third_derivative_type jump_in_third_derivatives(const unsigned int interface_dof_index, const unsigned int q_point) const
void get_jump_in_function_gradients_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
void get_jump_in_function_hessians_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
typename ProductType< Number, third_derivative_type >::type solution_third_derivative_type
void get_jump_in_function_hessians(const InputVector &fe_function, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
value_type jump_in_values(const unsigned int interface_dof_index, const unsigned int q_point) const
typename FEValuesViews::Vector< dim, spacedim >::gradient_type gradient_type
void get_jump_in_function_third_derivatives(const InputVector &fe_function, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &third_derivatives) const
gradient_type jump_gradient(const unsigned int interface_dof_index, const unsigned int q_point) const
value_type average_of_values(const unsigned int interface_dof_index, const unsigned int q_point) const
hessian_type jump_in_hessians(const unsigned int interface_dof_index, const unsigned int q_point) const
void get_jump_in_function_values_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
gradient_type average_of_gradients(const unsigned int interface_dof_index, const unsigned int q_point) const
gradient_type jump_in_gradients(const unsigned int interface_dof_index, const unsigned int q_point) const
hessian_type average_of_hessians(const unsigned int interface_dof_index, const unsigned int q_point) const
typename ProductType< Number, gradient_type >::type solution_gradient_type
::Tensor< 1, spacedim > gradient_type
::Tensor< 2, spacedim > hessian_type
::Tensor< 3, spacedim > third_derivative_type
::Tensor< 4, spacedim > third_derivative_type
::Tensor< 2, spacedim > gradient_type
::Tensor< 3, spacedim > hessian_type
::Tensor< 1, spacedim > value_type
Abstract base class for mapping classes.
Definition mapping.h:320
Definition point.h:113
static constexpr unsigned int n_independent_components
Definition tensor.h:498
friend class Vector
Definition vector.h:1124
Number value_type
Definition vector.h:137
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcOnlyAvailableWithHP()
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcOnlyAvailableWithoutHP()
static ::ExceptionBase & ExcMessage(std::string arg1)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition tria.h:1557
const bool IsBlockVector< VectorType >::value
UpdateFlags
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition mapping.cc:316
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
Definition mapping_q1.h:104
typename ::internal::FEInterfaceViews:: ViewType< dim, spacedim, Extractor >::type View
constexpr char U
Definition hp.h:117
constexpr types::global_dof_index invalid_dof_index
Definition types.h:269
constexpr unsigned int invalid_unsigned_int
Definition types.h:238
constexpr types::fe_index invalid_fe_index
Definition types.h:260
boost::integer_range< IncrementableType > iota_view
Definition iota_view.h:45
STL namespace.
Definition types.h:32
typename internal::ProductTypeImpl< std::decay_t< T >, std::decay_t< U > >::type type