Reference documentation for deal.II version 9.6.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 - 2024 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
175 value(const bool here_or_there,
176 const unsigned int interface_dof_index,
177 const unsigned int q_point) const;
178
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
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
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
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
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
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
728 value(const bool here_or_there,
729 const unsigned int interface_dof_index,
730 const unsigned int q_point) const;
731
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
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
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
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
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
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
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
1583 bool
1585
1597 double
1598 JxW(const unsigned int quadrature_point) const;
1599
1605 const std::vector<double> &
1607
1617 const std::vector<Tensor<1, spacedim>> &
1619
1629
1635 const std::vector<Point<spacedim>> &
1637
1648 unsigned
1650
1680
1691 std::vector<types::global_dof_index>
1693
1706 std::array<unsigned int, 2>
1707 interface_dof_to_dof_indices(const unsigned int interface_dof_index) const;
1708
1718 normal(const unsigned int q_point_index) const;
1719
1750 double
1751 shape_value(const bool here_or_there,
1752 const unsigned int interface_dof_index,
1753 const unsigned int q_point,
1754 const unsigned int component = 0) const;
1755
1783 double
1784 jump_in_shape_values(const unsigned int interface_dof_index,
1785 const unsigned int q_point,
1786 const unsigned int component = 0) const;
1787
1802 jump_in_shape_gradients(const unsigned int interface_dof_index,
1803 const unsigned int q_point,
1804 const unsigned int component = 0) const;
1805
1821 jump_in_shape_hessians(const unsigned int interface_dof_index,
1822 const unsigned int q_point,
1823 const unsigned int component = 0) const;
1824
1839 jump_in_shape_3rd_derivatives(const unsigned int interface_dof_index,
1840 const unsigned int q_point,
1841 const unsigned int component = 0) const;
1842
1865 double
1866 average_of_shape_values(const unsigned int interface_dof_index,
1867 const unsigned int q_point,
1868 const unsigned int component = 0) const;
1869
1884 average_of_shape_gradients(const unsigned int interface_dof_index,
1885 const unsigned int q_point,
1886 const unsigned int component = 0) const;
1887
1903 average_of_shape_hessians(const unsigned int interface_dof_index,
1904 const unsigned int q_point,
1905 const unsigned int component = 0) const;
1906
1926 template <class InputVector>
1927 void
1929 const InputVector &fe_function,
1930 std::vector<typename InputVector::value_type> &values) const;
1931
1940 template <class InputVector>
1941 void
1943 const InputVector &fe_function,
1945 &gradients) const;
1946
1954 template <class InputVector>
1955 void
1957 const InputVector &fe_function,
1959 &hessians) const;
1960
1969 template <class InputVector>
1970 void
1972 const InputVector &fe_function,
1974 &third_derivatives) const;
1975
1991 template <class InputVector>
1992 void
1994 const InputVector &fe_function,
1995 std::vector<typename InputVector::value_type> &values) const;
1996
2004 template <class InputVector>
2005 void
2007 const InputVector &fe_function,
2009 &gradients) const;
2010
2018 template <class InputVector>
2019 void
2021 const InputVector &fe_function,
2023 &hessians) const;
2024
2044
2053
2058private:
2062 std::vector<types::global_dof_index> interface_dof_indices;
2063
2069 std::vector<std::array<unsigned int, 2>> dofmap;
2070
2076
2083
2087 // non-hp data
2088
2092 std::unique_ptr<FEFaceValues<dim, spacedim>> internal_fe_face_values;
2093
2097 std::unique_ptr<FESubfaceValues<dim, spacedim>> internal_fe_subface_values;
2098
2102 std::unique_ptr<FEFaceValues<dim, spacedim>> internal_fe_face_values_neighbor;
2103
2107 std::unique_ptr<FESubfaceValues<dim, spacedim>>
2109
2110 // non-hp data
2111
2115 // hp data
2116
2121 std::unique_ptr<hp::FEFaceValues<dim, spacedim>> internal_hp_fe_face_values;
2122
2127 std::unique_ptr<hp::FESubfaceValues<dim, spacedim>>
2129
2134 std::unique_ptr<hp::FEFaceValues<dim, spacedim>>
2136
2141 std::unique_ptr<hp::FESubfaceValues<dim, spacedim>>
2143
2151 "The current function doesn't make sense when used with a "
2152 "FEInterfaceValues object with hp-capabilities.");
2153
2161 "The current function doesn't make sense when used with a "
2162 "FEInterfaceValues object without hp-capabilities.");
2163
2164 // hp data
2165
2166 /*
2167 * Make the view classes friends of this class, since they
2168 * access internal data.
2169 */
2170 template <int, int>
2172 template <int, int>
2174};
2175
2176
2177
2178#ifndef DOXYGEN
2179
2180/*---------------------- Inline functions ---------------------*/
2181
2182template <int dim, int spacedim>
2184 const Mapping<dim, spacedim> &mapping,
2186 const Quadrature<dim - 1> &quadrature,
2187 const UpdateFlags update_flags)
2188 : n_quadrature_points(quadrature.size())
2189 , fe_face_values(nullptr)
2190 , fe_face_values_neighbor(nullptr)
2191 , internal_fe_face_values(
2192 std::make_unique<FEFaceValues<dim, spacedim>>(mapping,
2193 fe,
2194 quadrature,
2195 update_flags))
2196 , internal_fe_subface_values(
2197 std::make_unique<FESubfaceValues<dim, spacedim>>(mapping,
2198 fe,
2199 quadrature,
2200 update_flags))
2201 , internal_fe_face_values_neighbor(
2202 std::make_unique<FEFaceValues<dim, spacedim>>(mapping,
2203 fe,
2204 quadrature,
2205 update_flags))
2206 , internal_fe_subface_values_neighbor(
2207 std::make_unique<FESubfaceValues<dim, spacedim>>(mapping,
2208 fe,
2209 quadrature,
2210 update_flags))
2211{}
2212
2213
2214
2215template <int dim, int spacedim>
2218 const Quadrature<dim - 1> &quadrature,
2219 const UpdateFlags update_flags)
2221 fe.reference_cell().template get_default_linear_mapping<dim, spacedim>(),
2222 fe,
2223 quadrature,
2224 update_flags)
2225{}
2226
2227
2228
2229template <int dim, int spacedim>
2231 const Mapping<dim, spacedim> &mapping,
2233 const hp::QCollection<dim - 1> &quadrature,
2234 const UpdateFlags update_flags)
2235 : n_quadrature_points(quadrature.max_n_quadrature_points())
2236 , fe_face_values(nullptr)
2237 , fe_face_values_neighbor(nullptr)
2238 , internal_fe_face_values(
2239 std::make_unique<FEFaceValues<dim, spacedim>>(mapping,
2240 fe,
2241 quadrature,
2242 update_flags))
2243 , internal_fe_subface_values(
2244 std::make_unique<FESubfaceValues<dim, spacedim>>(mapping,
2245 fe,
2246 quadrature,
2247 update_flags))
2248 , internal_fe_face_values_neighbor(
2249 std::make_unique<FEFaceValues<dim, spacedim>>(mapping,
2250 fe,
2251 quadrature[0],
2252 update_flags))
2253 , internal_fe_subface_values_neighbor(
2254 std::make_unique<FESubfaceValues<dim, spacedim>>(mapping,
2255 fe,
2256 quadrature[0],
2257 update_flags))
2258{}
2259
2260
2261
2262template <int dim, int spacedim>
2264 const hp::MappingCollection<dim, spacedim> &mapping_collection,
2265 const hp::FECollection<dim, spacedim> &fe_collection,
2266 const hp::QCollection<dim - 1> &quadrature_collection,
2267 const UpdateFlags update_flags)
2268 : n_quadrature_points(quadrature_collection.max_n_quadrature_points())
2269 , fe_face_values(nullptr)
2270 , fe_face_values_neighbor(nullptr)
2271 , internal_hp_fe_face_values(
2272 std::make_unique<hp::FEFaceValues<dim, spacedim>>(mapping_collection,
2273 fe_collection,
2274 quadrature_collection,
2275 update_flags))
2276 , internal_hp_fe_subface_values(
2277 std::make_unique<hp::FESubfaceValues<dim, spacedim>>(
2278 mapping_collection,
2279 fe_collection,
2280 quadrature_collection,
2281 update_flags))
2282 , internal_hp_fe_face_values_neighbor(
2283 std::make_unique<hp::FEFaceValues<dim, spacedim>>(mapping_collection,
2284 fe_collection,
2285 quadrature_collection,
2286 update_flags))
2287 , internal_hp_fe_subface_values_neighbor(
2288 std::make_unique<hp::FESubfaceValues<dim, spacedim>>(
2289 mapping_collection,
2290 fe_collection,
2291 quadrature_collection,
2292 update_flags))
2293{
2294 AssertDimension(dim, spacedim);
2295}
2296
2297
2298
2299template <int dim, int spacedim>
2301 const hp::FECollection<dim, spacedim> &fe_collection,
2302 const hp::QCollection<dim - 1> &quadrature_collection,
2303 const UpdateFlags update_flags)
2304 : FEInterfaceValues(fe_collection.get_reference_cell_default_linear_mapping(),
2305 fe_collection,
2306 quadrature_collection,
2307 update_flags)
2308{}
2309
2310
2311
2312template <int dim, int spacedim>
2313template <typename CellIteratorType, typename CellNeighborIteratorType>
2314void
2316 const CellIteratorType &cell,
2317 const unsigned int face_no,
2318 const unsigned int sub_face_no,
2319 const CellNeighborIteratorType &cell_neighbor,
2320 const unsigned int face_no_neighbor,
2321 const unsigned int sub_face_no_neighbor,
2322 const unsigned int q_index,
2323 const unsigned int mapping_index,
2324 const unsigned int fe_index_in,
2325 const unsigned int fe_index_neighbor_in)
2326{
2327 Assert(internal_fe_face_values || internal_hp_fe_face_values,
2329
2330 constexpr bool is_dof_cell_accessor =
2331 std::is_same_v<DoFCellAccessor<dim, spacedim, true>,
2332 typename CellIteratorType::AccessorType> ||
2333 std::is_same_v<DoFCellAccessor<dim, spacedim, false>,
2334 typename CellIteratorType::AccessorType>;
2335
2336 constexpr bool is_dof_cell_accessor_neighbor =
2337 std::is_same_v<DoFCellAccessor<dim, spacedim, true>,
2338 typename CellNeighborIteratorType::AccessorType> ||
2339 std::is_same_v<DoFCellAccessor<dim, spacedim, false>,
2340 typename CellNeighborIteratorType::AccessorType>;
2341
2342 if (internal_fe_face_values)
2343 {
2344 if (sub_face_no == numbers::invalid_unsigned_int)
2345 {
2346 internal_fe_face_values->reinit(cell, face_no);
2347 fe_face_values = internal_fe_face_values.get();
2348 }
2349 else
2350 {
2351 internal_fe_subface_values->reinit(cell, face_no, sub_face_no);
2352 fe_face_values = internal_fe_subface_values.get();
2353 }
2354 if (sub_face_no_neighbor == numbers::invalid_unsigned_int)
2355 {
2356 internal_fe_face_values_neighbor->reinit(cell_neighbor,
2357 face_no_neighbor);
2358 fe_face_values_neighbor = internal_fe_face_values_neighbor.get();
2359 }
2360 else
2361 {
2362 internal_fe_subface_values_neighbor->reinit(cell_neighbor,
2363 face_no_neighbor,
2364 sub_face_no_neighbor);
2365 fe_face_values_neighbor = internal_fe_subface_values_neighbor.get();
2366 }
2367
2368 AssertDimension(fe_face_values->n_quadrature_points,
2369 fe_face_values_neighbor->n_quadrature_points);
2370
2371 const_cast<unsigned int &>(this->n_quadrature_points) =
2372 fe_face_values->n_quadrature_points;
2373 }
2374 else if (internal_hp_fe_face_values)
2375 {
2376 unsigned int active_fe_index = fe_index_in;
2377 unsigned int active_fe_index_neighbor =
2378 (fe_index_neighbor_in != numbers::invalid_unsigned_int) ?
2379 fe_index_neighbor_in :
2380 fe_index_in;
2381
2382 if (active_fe_index == numbers::invalid_unsigned_int)
2383 {
2384 if constexpr (is_dof_cell_accessor)
2385 active_fe_index = cell->active_fe_index();
2386 else
2387 active_fe_index = 0;
2388 }
2389
2390 if (active_fe_index_neighbor == numbers::invalid_unsigned_int)
2391 {
2392 if constexpr (is_dof_cell_accessor_neighbor)
2393 active_fe_index_neighbor = cell_neighbor->active_fe_index();
2394 else
2395 active_fe_index_neighbor = 0;
2396 }
2397
2398 unsigned int used_q_index = q_index;
2399 unsigned int used_mapping_index = mapping_index;
2400
2401 // First check. If there is only one element in a collection, and if none
2402 // had been specified explicitly, then that's clearly the one to take:
2403 if (used_q_index == numbers::invalid_unsigned_int)
2404 if (internal_hp_fe_face_values->get_quadrature_collection().size() == 1)
2405 used_q_index = 0;
2406
2407 if (used_mapping_index == numbers::invalid_unsigned_int)
2408 if (internal_hp_fe_face_values->get_mapping_collection().size() == 1)
2409 used_mapping_index = 0;
2410
2411 // Second check: See if the two quadrature objects are the same, because
2412 // in that case it does not matter which one we use. Unfortunately, we
2413 // currently have no way of testing that two mapping objects are the
2414 // same :-(
2415 if (used_q_index == numbers::invalid_unsigned_int)
2416 if (internal_hp_fe_face_values
2417 ->get_quadrature_collection()[active_fe_index] ==
2418 internal_hp_fe_face_values
2419 ->get_quadrature_collection()[active_fe_index_neighbor])
2420 used_q_index = active_fe_index;
2421
2422 // Third check, if the above did not already suffice. We see if we
2423 // can get somewhere via the dominated's finite element index.
2424 const unsigned int dominated_fe_index =
2425 ((used_q_index == numbers::invalid_unsigned_int) ||
2426 (used_mapping_index == numbers::invalid_unsigned_int) ?
2427 internal_hp_fe_face_values->get_fe_collection().find_dominated_fe(
2428 {active_fe_index, active_fe_index_neighbor}) :
2429 numbers::invalid_unsigned_int);
2430
2431 if (used_q_index == numbers::invalid_unsigned_int)
2432 {
2433 Assert(dominated_fe_index != numbers::invalid_fe_index,
2434 ExcMessage(
2435 "You called this function with 'q_index' left at its "
2436 "default value, but this can only work if one of "
2437 "the two finite elements adjacent to this face "
2438 "dominates the other. See the documentation "
2439 "of this function for more information of how "
2440 "to deal with this situation."));
2441 used_q_index = dominated_fe_index;
2442 }
2443
2444 if (used_mapping_index == numbers::invalid_unsigned_int)
2445 {
2446 Assert(dominated_fe_index != numbers::invalid_fe_index,
2447 ExcMessage(
2448 "You called this function with 'mapping_index' left "
2449 "at its default value, but this can only work if one "
2450 "of the two finite elements adjacent to this face "
2451 "dominates the other. See the documentation "
2452 "of this function for more information of how "
2453 "to deal with this situation."));
2454 used_mapping_index = dominated_fe_index;
2455 }
2456
2457 // Same as if above, but when hp is enabled.
2458 if (sub_face_no == numbers::invalid_unsigned_int)
2459 {
2460 internal_hp_fe_face_values->reinit(
2461 cell, face_no, used_q_index, used_mapping_index, active_fe_index);
2462 fe_face_values = &const_cast<FEFaceValues<dim, spacedim> &>(
2463 internal_hp_fe_face_values->get_present_fe_values());
2464 }
2465 else
2466 {
2467 internal_hp_fe_subface_values->reinit(cell,
2468 face_no,
2469 sub_face_no,
2470 used_q_index,
2471 used_mapping_index,
2472 active_fe_index);
2473
2474 fe_face_values = &const_cast<FESubfaceValues<dim, spacedim> &>(
2475 internal_hp_fe_subface_values->get_present_fe_values());
2476 }
2477 if (sub_face_no_neighbor == numbers::invalid_unsigned_int)
2478 {
2479 internal_hp_fe_face_values_neighbor->reinit(cell_neighbor,
2480 face_no_neighbor,
2481 used_q_index,
2482 used_mapping_index,
2483 active_fe_index_neighbor);
2484
2485 fe_face_values_neighbor = &const_cast<FEFaceValues<dim, spacedim> &>(
2486 internal_hp_fe_face_values_neighbor->get_present_fe_values());
2487 }
2488 else
2489 {
2490 internal_hp_fe_subface_values_neighbor->reinit(
2491 cell_neighbor,
2492 face_no_neighbor,
2493 sub_face_no_neighbor,
2494 used_q_index,
2495 used_mapping_index,
2496 active_fe_index_neighbor);
2497
2498 fe_face_values_neighbor =
2499 &const_cast<FESubfaceValues<dim, spacedim> &>(
2500 internal_hp_fe_subface_values_neighbor->get_present_fe_values());
2501 }
2502
2503 AssertDimension(fe_face_values->n_quadrature_points,
2504 fe_face_values_neighbor->n_quadrature_points);
2505
2506 const_cast<unsigned int &>(this->n_quadrature_points) =
2507 fe_face_values->n_quadrature_points;
2508 }
2509
2510 // Set up dof mapping and remove duplicates (for continuous elements).
2511 if constexpr (is_dof_cell_accessor_neighbor && is_dof_cell_accessor)
2512 {
2513 // Get dof indices first:
2514 std::vector<types::global_dof_index> v(
2515 fe_face_values->get_fe().n_dofs_per_cell());
2516 cell->get_active_or_mg_dof_indices(v);
2517 std::vector<types::global_dof_index> v2(
2518 fe_face_values_neighbor->get_fe().n_dofs_per_cell());
2519 cell_neighbor->get_active_or_mg_dof_indices(v2);
2520
2521 // Fill a map from the global dof index to the left and right
2522 // local index.
2523 std::map<types::global_dof_index, std::pair<unsigned int, unsigned int>>
2524 tempmap;
2525 std::pair<unsigned int, unsigned int> invalid_entry(
2527
2528 for (unsigned int i = 0; i < v.size(); ++i)
2529 {
2530 // If not already existing, add an invalid entry:
2531 auto result = tempmap.insert(std::make_pair(v[i], invalid_entry));
2532 result.first->second.first = i;
2533 }
2534
2535 for (unsigned int i = 0; i < v2.size(); ++i)
2536 {
2537 // If not already existing, add an invalid entry:
2538 auto result = tempmap.insert(std::make_pair(v2[i], invalid_entry));
2539 result.first->second.second = i;
2540 }
2541
2542 // Transfer from the map to the sorted std::vectors.
2543 dofmap.resize(tempmap.size());
2544 interface_dof_indices.resize(tempmap.size());
2545 unsigned int idx = 0;
2546 for (auto &x : tempmap)
2547 {
2548 interface_dof_indices[idx] = x.first;
2549 dofmap[idx] = {{x.second.first, x.second.second}};
2550 ++idx;
2551 }
2552 }
2553}
2554
2555
2556
2557template <int dim, int spacedim>
2558template <typename CellIteratorType>
2559void
2560FEInterfaceValues<dim, spacedim>::reinit(const CellIteratorType &cell,
2561 const unsigned int face_no,
2562 const unsigned int q_index,
2563 const unsigned int mapping_index,
2564 const unsigned int fe_index)
2565{
2566 Assert(internal_fe_face_values || internal_hp_fe_face_values,
2568
2569 if (internal_fe_face_values)
2570 {
2571 Assert((q_index == 0 || q_index == numbers::invalid_unsigned_int),
2573 Assert((mapping_index == 0 ||
2574 mapping_index == numbers::invalid_unsigned_int),
2576 Assert((fe_index == 0 || fe_index == numbers::invalid_unsigned_int),
2578
2579 internal_fe_face_values->reinit(cell, face_no);
2580 fe_face_values = internal_fe_face_values.get();
2581 fe_face_values_neighbor = nullptr;
2582 }
2583 else if (internal_hp_fe_face_values)
2584 {
2585 internal_hp_fe_face_values->reinit(
2586 cell, face_no, q_index, mapping_index, fe_index);
2587 fe_face_values = &const_cast<FEFaceValues<dim> &>(
2588 internal_hp_fe_face_values->get_present_fe_values());
2589 fe_face_values_neighbor = nullptr;
2590 }
2591
2592 if constexpr (std::is_same_v<typename CellIteratorType::AccessorType,
2594 std::is_same_v<typename CellIteratorType::AccessorType,
2596 {
2597 if (internal_fe_face_values)
2598 {
2599 interface_dof_indices.resize(
2600 fe_face_values->get_fe().n_dofs_per_cell());
2601 cell->get_active_or_mg_dof_indices(interface_dof_indices);
2602 }
2603 else if (internal_hp_fe_face_values)
2604 {
2605 interface_dof_indices.resize(
2606 fe_face_values->get_fe().n_dofs_per_cell());
2607 cell->get_active_or_mg_dof_indices(interface_dof_indices);
2608 }
2609
2610 dofmap.resize(interface_dof_indices.size());
2611 for (unsigned int i = 0; i < interface_dof_indices.size(); ++i)
2612 {
2613 dofmap[i] = {{i, numbers::invalid_unsigned_int}};
2614 }
2615 }
2616}
2617
2618
2619
2620template <int dim, int spacedim>
2621inline double
2622FEInterfaceValues<dim, spacedim>::JxW(const unsigned int q) const
2623{
2624 Assert(fe_face_values != nullptr,
2625 ExcMessage("This call requires a call to reinit() first."));
2626 return fe_face_values->JxW(q);
2627}
2628
2629
2630
2631template <int dim, int spacedim>
2632const std::vector<double> &
2634{
2635 Assert(fe_face_values != nullptr,
2636 ExcMessage("This call requires a call to reinit() first."));
2637 return fe_face_values->get_JxW_values();
2638}
2639
2640
2641
2642template <int dim, int spacedim>
2643const std::vector<Tensor<1, spacedim>> &
2645{
2646 Assert(fe_face_values != nullptr,
2647 ExcMessage("This call requires a call to reinit() first."));
2648 return fe_face_values->get_normal_vectors();
2649}
2650
2651
2652
2653template <int dim, int spacedim>
2656{
2657 Assert(!has_hp_capabilities(), ExcOnlyAvailableWithoutHP());
2658 return internal_fe_face_values->get_mapping();
2659}
2660
2661
2662
2663template <int dim, int spacedim>
2666{
2667 Assert(!has_hp_capabilities(), ExcOnlyAvailableWithoutHP());
2668 return internal_fe_face_values->get_fe();
2669}
2670
2671
2672
2673template <int dim, int spacedim>
2674const Quadrature<dim - 1> &
2676{
2677 Assert(!has_hp_capabilities(), ExcOnlyAvailableWithoutHP());
2678 return internal_fe_face_values->get_quadrature();
2679}
2680
2681
2682
2683template <int dim, int spacedim>
2686{
2687 Assert(has_hp_capabilities(), ExcOnlyAvailableWithHP());
2688 return internal_hp_fe_face_values->get_mapping_collection();
2689}
2690
2691
2692
2693template <int dim, int spacedim>
2696{
2697 Assert(has_hp_capabilities(), ExcOnlyAvailableWithHP());
2698 return internal_hp_fe_face_values->get_fe_collection();
2699}
2700
2701
2702
2703template <int dim, int spacedim>
2704const hp::QCollection<dim - 1> &
2706{
2707 Assert(has_hp_capabilities(), ExcOnlyAvailableWithHP());
2708 return internal_hp_fe_face_values->get_quadrature_collection();
2709}
2710
2711
2712
2713template <int dim, int spacedim>
2714bool
2716{
2717 if (internal_hp_fe_face_values || internal_hp_fe_subface_values ||
2718 internal_hp_fe_face_values_neighbor ||
2719 internal_hp_fe_subface_values_neighbor)
2720 {
2721 Assert(!internal_fe_face_values, ExcInternalError());
2722 Assert(!internal_fe_subface_values, ExcInternalError());
2723 Assert(!internal_fe_face_values_neighbor, ExcInternalError());
2724 Assert(!internal_fe_subface_values_neighbor, ExcInternalError());
2725
2726 return true;
2727 }
2728
2729 Assert(internal_fe_face_values || internal_fe_subface_values ||
2730 internal_fe_face_values_neighbor ||
2731 internal_fe_subface_values_neighbor,
2733 Assert(!internal_hp_fe_face_values, ExcInternalError());
2734 Assert(!internal_hp_fe_subface_values, ExcInternalError());
2735 Assert(!internal_hp_fe_face_values_neighbor, ExcInternalError());
2736 Assert(!internal_hp_fe_subface_values_neighbor, ExcInternalError());
2737
2738 return false;
2739}
2740
2741
2742
2743template <int dim, int spacedim>
2746{
2747 return {0U, n_quadrature_points};
2748}
2749
2750
2751
2752template <int dim, int spacedim>
2753const std::vector<Point<spacedim>> &
2755{
2756 Assert(fe_face_values != nullptr,
2757 ExcMessage("This call requires a call to reinit() first."));
2758 return fe_face_values->get_quadrature_points();
2759}
2760
2761
2762
2763template <int dim, int spacedim>
2766{
2767 if (has_hp_capabilities())
2768 return internal_hp_fe_face_values->get_update_flags();
2769 else
2770 return internal_fe_face_values->get_update_flags();
2771}
2772
2773
2774
2775template <int dim, int spacedim>
2778{
2779 return get_fe_face_values(cell_index).get_cell();
2780}
2781
2782
2783
2784template <int dim, int spacedim>
2785inline unsigned int
2787 const unsigned int cell_index) const
2788{
2789 return get_fe_face_values(cell_index).get_face_number();
2790}
2791
2792
2793
2794template <int dim, int spacedim>
2795unsigned
2797{
2798 Assert(
2799 interface_dof_indices.size() > 0,
2800 ExcMessage(
2801 "n_current_interface_dofs() is only available after a call to reinit()."));
2802 return interface_dof_indices.size();
2803}
2804
2805
2806
2807template <int dim, int spacedim>
2810{
2811 return {0U, n_current_interface_dofs()};
2812}
2813
2814
2815
2816template <int dim, int spacedim>
2817bool
2819{
2820 return fe_face_values_neighbor == nullptr;
2821}
2822
2823
2824
2825template <int dim, int spacedim>
2826std::vector<types::global_dof_index>
2828{
2829 return interface_dof_indices;
2830}
2831
2832
2833
2834template <int dim, int spacedim>
2835std::array<unsigned int, 2>
2837 const unsigned int interface_dof_index) const
2838{
2839 AssertIndexRange(interface_dof_index, n_current_interface_dofs());
2840 return dofmap[interface_dof_index];
2841}
2842
2843
2844
2845template <int dim, int spacedim>
2848 const unsigned int cell_index) const
2849{
2851 Assert(
2852 cell_index == 0 || !at_boundary(),
2853 ExcMessage(
2854 "You are on a boundary, so you can only ask for the first FEFaceValues object."));
2855
2856 return (cell_index == 0) ? *fe_face_values : *fe_face_values_neighbor;
2857}
2858
2859
2860
2861template <int dim, int spacedim>
2863FEInterfaceValues<dim, spacedim>::normal(const unsigned int q_point_index) const
2864{
2865 return fe_face_values->normal_vector(q_point_index);
2866}
2867
2868
2869
2870template <int dim, int spacedim>
2871double
2873 const bool here_or_there,
2874 const unsigned int interface_dof_index,
2875 const unsigned int q_point,
2876 const unsigned int component) const
2877{
2878 const auto dof_pair = dofmap[interface_dof_index];
2879
2880 if (here_or_there && dof_pair[0] != numbers::invalid_unsigned_int)
2881 return get_fe_face_values(0).shape_value_component(dof_pair[0],
2882 q_point,
2883 component);
2884 if (!here_or_there && dof_pair[1] != numbers::invalid_unsigned_int)
2885 return get_fe_face_values(1).shape_value_component(dof_pair[1],
2886 q_point,
2887 component);
2888
2889 return 0.0;
2890}
2891
2892
2893
2894template <int dim, int spacedim>
2895double
2897 const unsigned int interface_dof_index,
2898 const unsigned int q_point,
2899 const unsigned int component) const
2900{
2901 const auto dof_pair = dofmap[interface_dof_index];
2902
2903 double value = 0.0;
2904
2905 if (dof_pair[0] != numbers::invalid_unsigned_int)
2906 value += get_fe_face_values(0).shape_value_component(dof_pair[0],
2907 q_point,
2908 component);
2909 if (dof_pair[1] != numbers::invalid_unsigned_int)
2910 value -= get_fe_face_values(1).shape_value_component(dof_pair[1],
2911 q_point,
2912 component);
2913 return value;
2914}
2915
2916
2917
2918template <int dim, int spacedim>
2919double
2921 const unsigned int interface_dof_index,
2922 const unsigned int q_point,
2923 const unsigned int component) const
2924{
2925 const auto dof_pair = dofmap[interface_dof_index];
2926
2927 if (at_boundary())
2928 return get_fe_face_values(0).shape_value_component(dof_pair[0],
2929 q_point,
2930 component);
2931
2932 double value = 0.0;
2933
2934 if (dof_pair[0] != numbers::invalid_unsigned_int)
2935 value += 0.5 * get_fe_face_values(0).shape_value_component(dof_pair[0],
2936 q_point,
2937 component);
2938 if (dof_pair[1] != numbers::invalid_unsigned_int)
2939 value += 0.5 * get_fe_face_values(1).shape_value_component(dof_pair[1],
2940 q_point,
2941 component);
2942
2943 return value;
2944}
2945
2946
2947
2948template <int dim, int spacedim>
2951 const unsigned int interface_dof_index,
2952 const unsigned int q_point,
2953 const unsigned int component) const
2954{
2955 const auto dof_pair = dofmap[interface_dof_index];
2956
2957 if (at_boundary())
2958 return get_fe_face_values(0).shape_grad_component(dof_pair[0],
2959 q_point,
2960 component);
2961
2963
2964 if (dof_pair[0] != numbers::invalid_unsigned_int)
2965 value += 0.5 * get_fe_face_values(0).shape_grad_component(dof_pair[0],
2966 q_point,
2967 component);
2968 if (dof_pair[1] != numbers::invalid_unsigned_int)
2969 value += 0.5 * get_fe_face_values(1).shape_grad_component(dof_pair[1],
2970 q_point,
2971 component);
2972
2973 return value;
2974}
2975
2976
2977
2978template <int dim, int spacedim>
2981 const unsigned int interface_dof_index,
2982 const unsigned int q_point,
2983 const unsigned int component) const
2984{
2985 const auto dof_pair = dofmap[interface_dof_index];
2986
2987 if (at_boundary())
2988 return get_fe_face_values(0).shape_hessian_component(dof_pair[0],
2989 q_point,
2990 component);
2991
2993
2994 if (dof_pair[0] != numbers::invalid_unsigned_int)
2995 value += 0.5 * get_fe_face_values(0).shape_hessian_component(dof_pair[0],
2996 q_point,
2997 component);
2998 if (dof_pair[1] != numbers::invalid_unsigned_int)
2999 value += 0.5 * get_fe_face_values(1).shape_hessian_component(dof_pair[1],
3000 q_point,
3001 component);
3002
3003 return value;
3004}
3005
3006
3007
3008template <int dim, int spacedim>
3011 const unsigned int interface_dof_index,
3012 const unsigned int q_point,
3013 const unsigned int component) const
3014{
3015 const auto dof_pair = dofmap[interface_dof_index];
3016
3017 if (at_boundary())
3018 return get_fe_face_values(0).shape_grad_component(dof_pair[0],
3019 q_point,
3020 component);
3021
3023
3024 if (dof_pair[0] != numbers::invalid_unsigned_int)
3025 value += get_fe_face_values(0).shape_grad_component(dof_pair[0],
3026 q_point,
3027 component);
3028 if (dof_pair[1] != numbers::invalid_unsigned_int)
3029 value -= get_fe_face_values(1).shape_grad_component(dof_pair[1],
3030 q_point,
3031 component);
3032
3033 return value;
3034}
3035
3036
3037
3038template <int dim, int spacedim>
3041 const unsigned int interface_dof_index,
3042 const unsigned int q_point,
3043 const unsigned int component) const
3044{
3045 const auto dof_pair = dofmap[interface_dof_index];
3046
3047 if (at_boundary())
3048 return get_fe_face_values(0).shape_hessian_component(dof_pair[0],
3049 q_point,
3050 component);
3051
3053
3054 if (dof_pair[0] != numbers::invalid_unsigned_int)
3055 value += get_fe_face_values(0).shape_hessian_component(dof_pair[0],
3056 q_point,
3057 component);
3058 if (dof_pair[1] != numbers::invalid_unsigned_int)
3059 value -= get_fe_face_values(1).shape_hessian_component(dof_pair[1],
3060 q_point,
3061 component);
3062
3063 return value;
3064}
3065
3066
3067
3068template <int dim, int spacedim>
3071 const unsigned int interface_dof_index,
3072 const unsigned int q_point,
3073 const unsigned int component) const
3074{
3075 const auto dof_pair = dofmap[interface_dof_index];
3076
3077 if (at_boundary())
3078 return get_fe_face_values(0).shape_3rd_derivative_component(dof_pair[0],
3079 q_point,
3080 component);
3081
3083
3084 if (dof_pair[0] != numbers::invalid_unsigned_int)
3085 value += get_fe_face_values(0).shape_3rd_derivative_component(dof_pair[0],
3086 q_point,
3087 component);
3088 if (dof_pair[1] != numbers::invalid_unsigned_int)
3089 value -= get_fe_face_values(1).shape_3rd_derivative_component(dof_pair[1],
3090 q_point,
3091 component);
3092
3093 return value;
3094}
3095
3096
3097
3098template <int dim, int spacedim>
3099template <class InputVector>
3100void
3102 const InputVector &fe_function,
3103 std::vector<typename InputVector::value_type> &values) const
3104{
3105 AssertDimension(values.size(), n_quadrature_points);
3106
3108 this->operator[](scalar).get_jump_in_function_values(fe_function, values);
3109}
3110
3111
3112
3113template <int dim, int spacedim>
3114template <class InputVector>
3115void
3117 const InputVector &fe_function,
3119 const
3120{
3121 AssertDimension(gradients.size(), n_quadrature_points);
3122
3124 this->operator[](scalar).get_jump_in_function_gradients(fe_function,
3125 gradients);
3126}
3127
3128
3129
3130template <int dim, int spacedim>
3131template <class InputVector>
3132void
3134 const InputVector &fe_function,
3136 const
3137{
3138 AssertDimension(hessians.size(), n_quadrature_points);
3139
3141 this->operator[](scalar).get_jump_in_function_hessians(fe_function, hessians);
3142}
3143
3144
3145
3146template <int dim, int spacedim>
3147template <class InputVector>
3148void
3150 const InputVector &fe_function,
3152 &third_derivatives) const
3153{
3154 AssertDimension(third_derivatives.size(), n_quadrature_points);
3155
3157 this->operator[](scalar).get_jump_in_function_third_derivatives(
3158 fe_function, third_derivatives);
3159}
3160
3161
3162
3163template <int dim, int spacedim>
3164template <class InputVector>
3165void
3167 const InputVector &fe_function,
3168 std::vector<typename InputVector::value_type> &values) const
3169{
3170 AssertDimension(values.size(), n_quadrature_points);
3171
3173 this->operator[](scalar).get_average_of_function_values(fe_function, values);
3174}
3175
3176
3177
3178template <int dim, int spacedim>
3179template <class InputVector>
3180void
3182 const InputVector &fe_function,
3184 const
3185{
3186 AssertDimension(gradients.size(), n_quadrature_points);
3187
3189 this->operator[](scalar).get_average_of_function_gradients(fe_function,
3190 gradients);
3191}
3192
3193
3194
3195template <int dim, int spacedim>
3196template <class InputVector>
3197void
3199 const InputVector &fe_function,
3201 const
3202{
3203 AssertDimension(hessians.size(), n_quadrature_points);
3204
3206 this->operator[](scalar).get_average_of_function_hessians(fe_function,
3207 hessians);
3208}
3209
3210
3211
3212/*------------ Inline functions: FEInterfaceValues------------*/
3213template <int dim, int spacedim>
3216 const FEValuesExtractors::Scalar &scalar) const
3217{
3218 const unsigned int n_components =
3219 (this->has_hp_capabilities() ? this->get_fe_collection().n_components() :
3220 this->get_fe().n_components());
3221 (void)n_components;
3222 AssertIndexRange(scalar.component, n_components);
3223 return FEInterfaceViews::Scalar<dim, spacedim>(*this, scalar.component);
3224}
3225
3226
3227
3228template <int dim, int spacedim>
3231 const FEValuesExtractors::Vector &vector) const
3232{
3233 const unsigned int n_components =
3234 (this->has_hp_capabilities() ? this->get_fe_collection().n_components() :
3235 this->get_fe().n_components());
3236 const unsigned int n_vectors =
3239 0);
3240 (void)n_components;
3241 (void)n_vectors;
3242 AssertIndexRange(vector.first_vector_component, n_vectors);
3244 vector.first_vector_component);
3245}
3246
3247
3248
3249namespace FEInterfaceViews
3250{
3251 template <int dim, int spacedim>
3253 const FEInterfaceValues<dim, spacedim> &fe_interface)
3254 : fe_interface(&fe_interface)
3255 {}
3256
3257
3258
3259 template <int dim, int spacedim>
3260 Scalar<dim, spacedim>::Scalar(
3261 const FEInterfaceValues<dim, spacedim> &fe_interface,
3262 const unsigned int component)
3263 : Base<dim, spacedim>(fe_interface)
3264 , extractor(component)
3265 {}
3266
3267
3268
3269 template <int dim, int spacedim>
3270 template <class InputVector, class OutputVector>
3271 void
3272 Base<dim, spacedim>::get_local_dof_values(
3273 const InputVector &dof_values,
3274 OutputVector &local_dof_values) const
3275 {
3276 const auto &interface_dof_indices =
3277 this->fe_interface->get_interface_dof_indices();
3278
3279 AssertDimension(interface_dof_indices.size(), local_dof_values.size());
3280
3281 for (const unsigned int i : this->fe_interface->dof_indices())
3282 local_dof_values[i] = dof_values(interface_dof_indices[i]);
3283 }
3284
3285
3286
3287 template <int dim, int spacedim>
3288 typename Scalar<dim, spacedim>::value_type
3289 Scalar<dim, spacedim>::value(const bool here_or_there,
3290 const unsigned int interface_dof_index,
3291 const unsigned int q_point) const
3292 {
3293 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3294
3295 if (here_or_there && dof_pair[0] != numbers::invalid_unsigned_int)
3296 return (*(this->fe_interface->fe_face_values))[extractor].value(
3297 dof_pair[0], q_point);
3298
3299 if (!here_or_there && dof_pair[1] != numbers::invalid_unsigned_int)
3300 return (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3301 dof_pair[1], q_point);
3302
3303 return 0.0;
3304 }
3305
3306
3307
3308 template <int dim, int spacedim>
3309 typename Scalar<dim, spacedim>::value_type
3310 Scalar<dim, spacedim>::jump_in_values(const unsigned int interface_dof_index,
3311 const unsigned int q_point) const
3312 {
3313 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3314
3315 value_type value = 0.0;
3316
3317 if (dof_pair[0] != numbers::invalid_unsigned_int)
3318 value +=
3319 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3320 q_point);
3321
3322 if (dof_pair[1] != numbers::invalid_unsigned_int)
3323 value -=
3324 (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3325 dof_pair[1], q_point);
3326
3327 return value;
3328 }
3329
3330
3331
3332 template <int dim, int spacedim>
3333 typename Scalar<dim, spacedim>::value_type
3334 Scalar<dim, spacedim>::average_of_values(
3335 const unsigned int interface_dof_index,
3336 const unsigned int q_point) const
3337 {
3338 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3339
3340 if (this->fe_interface->at_boundary())
3341 return (*(this->fe_interface->fe_face_values))[extractor].value(
3342 dof_pair[0], q_point);
3343
3344 value_type value = 0.0;
3345
3346 if (dof_pair[0] != numbers::invalid_unsigned_int)
3347 value +=
3348 0.5 *
3349 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3350 q_point);
3351
3352 if (dof_pair[1] != numbers::invalid_unsigned_int)
3353 value +=
3354 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3355 dof_pair[1], q_point);
3356
3357 return value;
3358 }
3359
3360
3361
3362 template <int dim, int spacedim>
3363 typename Scalar<dim, spacedim>::gradient_type
3364 Scalar<dim, spacedim>::average_of_gradients(
3365 const unsigned int interface_dof_index,
3366 const unsigned int q_point) const
3367 {
3368 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3369
3370 if (this->fe_interface->at_boundary())
3371 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3372 dof_pair[0], q_point);
3373
3374 gradient_type value;
3375
3376 if (dof_pair[0] != numbers::invalid_unsigned_int)
3377 value +=
3378 0.5 *
3379 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3380 q_point);
3381
3382 if (dof_pair[1] != numbers::invalid_unsigned_int)
3383 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3384 .gradient(dof_pair[1], q_point);
3385
3386 return value;
3387 }
3388
3389
3390
3391 template <int dim, int spacedim>
3392 typename Scalar<dim, spacedim>::gradient_type
3393 Scalar<dim, spacedim>::jump_in_gradients(
3394 const unsigned int interface_dof_index,
3395 const unsigned int q_point) const
3396 {
3397 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3398
3399 if (this->fe_interface->at_boundary())
3400 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3401 dof_pair[0], q_point);
3402
3403 gradient_type value;
3404
3405 if (dof_pair[0] != numbers::invalid_unsigned_int)
3406 value +=
3407 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3408 q_point);
3409
3410 if (dof_pair[1] != numbers::invalid_unsigned_int)
3411 value -=
3412 (*(this->fe_interface->fe_face_values_neighbor))[extractor].gradient(
3413 dof_pair[1], q_point);
3414
3415 return value;
3416 }
3417
3418
3419
3420 template <int dim, int spacedim>
3421 typename Scalar<dim, spacedim>::hessian_type
3422 Scalar<dim, spacedim>::average_of_hessians(
3423 const unsigned int interface_dof_index,
3424 const unsigned int q_point) const
3425 {
3426 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3427
3428 if (this->fe_interface->at_boundary())
3429 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
3430 dof_pair[0], q_point);
3431
3432 hessian_type value;
3433
3434 if (dof_pair[0] != numbers::invalid_unsigned_int)
3435 value +=
3436 0.5 *
3437 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
3438 q_point);
3439
3440 if (dof_pair[1] != numbers::invalid_unsigned_int)
3441 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3442 .hessian(dof_pair[1], q_point);
3443
3444 return value;
3445 }
3446
3447
3448
3449 template <int dim, int spacedim>
3450 typename Scalar<dim, spacedim>::third_derivative_type
3451 Scalar<dim, spacedim>::jump_in_third_derivatives(
3452 const unsigned int interface_dof_index,
3453 const unsigned int q_point) const
3454 {
3455 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3456
3457 if (this->fe_interface->at_boundary())
3458 return (*(this->fe_interface->fe_face_values))[extractor]
3459 .third_derivative(dof_pair[0], q_point);
3460
3461 third_derivative_type value;
3462
3463 if (dof_pair[0] != numbers::invalid_unsigned_int)
3464 value +=
3465 (*(this->fe_interface->fe_face_values))[extractor].third_derivative(
3466 dof_pair[0], q_point);
3467
3468 if (dof_pair[1] != numbers::invalid_unsigned_int)
3469 value -= (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3470 .third_derivative(dof_pair[1], q_point);
3471
3472 return value;
3473 }
3474
3475
3476
3477 template <int dim, int spacedim>
3478 typename Scalar<dim, spacedim>::hessian_type
3479 Scalar<dim, spacedim>::jump_in_hessians(
3480 const unsigned int interface_dof_index,
3481 const unsigned int q_point) const
3482 {
3483 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3484
3485 if (this->fe_interface->at_boundary())
3486 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
3487 dof_pair[0], q_point);
3488
3489 hessian_type value;
3490
3491 if (dof_pair[0] != numbers::invalid_unsigned_int)
3492 value +=
3493 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
3494 q_point);
3495
3496 if (dof_pair[1] != numbers::invalid_unsigned_int)
3497 value -=
3498 (*(this->fe_interface->fe_face_values_neighbor))[extractor].hessian(
3499 dof_pair[1], q_point);
3500
3501 return value;
3502 }
3503
3504
3505
3506 template <int dim, int spacedim>
3507 template <class InputVector>
3508 void
3509 Scalar<dim, spacedim>::get_function_values_from_local_dof_values(
3510 const bool here_or_there,
3511 const InputVector &local_dof_values,
3512 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3513 const
3514 {
3515 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
3516
3517 for (const auto dof_index : this->fe_interface->dof_indices())
3518 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3519 {
3520 if (dof_index == 0)
3521 values[q_index] = 0.;
3522
3523 values[q_index] += local_dof_values[dof_index] *
3524 value(here_or_there, dof_index, q_index);
3525 }
3526 }
3527
3528
3529
3530 template <int dim, int spacedim>
3531 template <class InputVector>
3532 void
3533 Scalar<dim, spacedim>::get_function_values(
3534 const bool here_or_there,
3535 const InputVector &fe_function,
3536 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3537 const
3538 {
3539 std::vector<typename InputVector::value_type> local_dof_values(
3540 this->fe_interface->n_current_interface_dofs());
3541 this->get_local_dof_values(fe_function, local_dof_values);
3542
3543 get_function_values_from_local_dof_values(here_or_there,
3544 local_dof_values,
3545 values);
3546 }
3547
3548
3549
3550 template <int dim, int spacedim>
3551 template <class InputVector>
3552 void
3553 Scalar<dim, spacedim>::get_jump_in_function_values_from_local_dof_values(
3554 const InputVector &local_dof_values,
3555 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3556 const
3557 {
3558 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
3559
3560 for (const auto dof_index : this->fe_interface->dof_indices())
3561 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3562 {
3563 if (dof_index == 0)
3564 values[q_index] = 0.;
3565
3566 values[q_index] +=
3567 local_dof_values[dof_index] * jump_in_values(dof_index, q_index);
3568 }
3569 }
3570
3571
3572
3573 template <int dim, int spacedim>
3574 template <class InputVector>
3575 void
3576 Scalar<dim, spacedim>::get_jump_in_function_values(
3577 const InputVector &fe_function,
3578 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3579 const
3580 {
3581 std::vector<typename InputVector::value_type> local_dof_values(
3582 this->fe_interface->n_current_interface_dofs());
3583 this->get_local_dof_values(fe_function, local_dof_values);
3584
3585 get_jump_in_function_values_from_local_dof_values(local_dof_values, values);
3586 }
3587
3588
3589
3590 template <int dim, int spacedim>
3591 template <class InputVector>
3592 void
3593 Scalar<dim, spacedim>::get_jump_in_function_gradients_from_local_dof_values(
3594 const InputVector &local_dof_values,
3595 std::vector<solution_gradient_type<typename InputVector::value_type>>
3596 &gradients) const
3597 {
3598 AssertDimension(gradients.size(), this->fe_interface->n_quadrature_points);
3599
3600 for (const auto dof_index : this->fe_interface->dof_indices())
3601 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3602 {
3603 if (dof_index == 0)
3604 gradients[q_index] = 0.;
3605
3606 gradients[q_index] +=
3607 local_dof_values[dof_index] * jump_in_gradients(dof_index, q_index);
3608 }
3609 }
3610
3611
3612
3613 template <int dim, int spacedim>
3614 template <class InputVector>
3615 void
3616 Scalar<dim, spacedim>::get_jump_in_function_gradients(
3617 const InputVector &fe_function,
3618 std::vector<solution_gradient_type<typename InputVector::value_type>>
3619 &gradients) const
3620 {
3621 std::vector<typename InputVector::value_type> local_dof_values(
3622 this->fe_interface->n_current_interface_dofs());
3623 this->get_local_dof_values(fe_function, local_dof_values);
3624
3625 get_jump_in_function_gradients_from_local_dof_values(local_dof_values,
3626 gradients);
3627 }
3628
3629
3630
3631 template <int dim, int spacedim>
3632 template <class InputVector>
3633 void
3634 Scalar<dim, spacedim>::get_average_of_function_values_from_local_dof_values(
3635 const InputVector &local_dof_values,
3636 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3637 const
3638 {
3639 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
3640
3641 for (const auto dof_index : this->fe_interface->dof_indices())
3642 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3643 {
3644 if (dof_index == 0)
3645 values[q_index] = 0.;
3646
3647 values[q_index] +=
3648 local_dof_values[dof_index] * average_of_values(dof_index, q_index);
3649 }
3650 }
3651
3652
3653
3654 template <int dim, int spacedim>
3655 template <class InputVector>
3656 void
3657 Scalar<dim, spacedim>::get_average_of_function_values(
3658 const InputVector &fe_function,
3659 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3660 const
3661 {
3662 std::vector<typename InputVector::value_type> local_dof_values(
3663 this->fe_interface->n_current_interface_dofs());
3664 this->get_local_dof_values(fe_function, local_dof_values);
3665
3666 get_average_of_function_values_from_local_dof_values(local_dof_values,
3667 values);
3668 }
3669
3670
3671
3672 template <int dim, int spacedim>
3673 template <class InputVector>
3674 void
3675 Scalar<dim, spacedim>::
3676 get_average_of_function_gradients_from_local_dof_values(
3677 const InputVector &local_dof_values,
3678 std::vector<solution_gradient_type<typename InputVector::value_type>>
3679 &gradients) const
3680 {
3681 AssertDimension(gradients.size(), this->fe_interface->n_quadrature_points);
3682
3683 for (const auto dof_index : this->fe_interface->dof_indices())
3684 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3685 {
3686 if (dof_index == 0)
3687 gradients[q_index] = 0.;
3688
3689 gradients[q_index] += local_dof_values[dof_index] *
3690 average_of_gradients(dof_index, q_index);
3691 }
3692 }
3693
3694
3695
3696 template <int dim, int spacedim>
3697 template <class InputVector>
3698 void
3699 Scalar<dim, spacedim>::get_average_of_function_gradients(
3700 const InputVector &fe_function,
3701 std::vector<solution_gradient_type<typename InputVector::value_type>>
3702 &gradients) const
3703 {
3704 std::vector<typename InputVector::value_type> local_dof_values(
3705 this->fe_interface->n_current_interface_dofs());
3706 this->get_local_dof_values(fe_function, local_dof_values);
3707
3708 get_average_of_function_gradients_from_local_dof_values(local_dof_values,
3709 gradients);
3710 }
3711
3712
3713
3714 template <int dim, int spacedim>
3715 template <class InputVector>
3716 void
3717 Scalar<dim, spacedim>::get_jump_in_function_hessians_from_local_dof_values(
3718 const InputVector &local_dof_values,
3719 std::vector<solution_hessian_type<typename InputVector::value_type>>
3720 &hessians) const
3721 {
3722 AssertDimension(hessians.size(), this->fe_interface->n_quadrature_points);
3723
3724 for (const auto dof_index : this->fe_interface->dof_indices())
3725 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3726 {
3727 if (dof_index == 0)
3728 hessians[q_index] = 0.;
3729
3730 hessians[q_index] +=
3731 local_dof_values[dof_index] * jump_in_hessians(dof_index, q_index);
3732 }
3733 }
3734
3735
3736
3737 template <int dim, int spacedim>
3738 template <class InputVector>
3739 void
3740 Scalar<dim, spacedim>::get_jump_in_function_hessians(
3741 const InputVector &fe_function,
3742 std::vector<solution_hessian_type<typename InputVector::value_type>>
3743 &hessians) const
3744 {
3745 std::vector<typename InputVector::value_type> local_dof_values(
3746 this->fe_interface->n_current_interface_dofs());
3747 this->get_local_dof_values(fe_function, local_dof_values);
3748
3749 get_jump_in_function_hessians_from_local_dof_values(local_dof_values,
3750 hessians);
3751 }
3752
3753
3754
3755 template <int dim, int spacedim>
3756 template <class InputVector>
3757 void
3758 Scalar<dim, spacedim>::get_average_of_function_hessians_from_local_dof_values(
3759 const InputVector &local_dof_values,
3760 std::vector<solution_hessian_type<typename InputVector::value_type>>
3761 &hessians) const
3762 {
3763 AssertDimension(hessians.size(), this->fe_interface->n_quadrature_points);
3764
3765 for (const unsigned int dof_index : this->fe_interface->dof_indices())
3766 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3767 {
3768 if (dof_index == 0)
3769 hessians[q_index] = 0.;
3770
3771 hessians[q_index] += local_dof_values[dof_index] *
3772 average_of_hessians(dof_index, q_index);
3773 }
3774 }
3775
3776
3777
3778 template <int dim, int spacedim>
3779 template <class InputVector>
3780 void
3781 Scalar<dim, spacedim>::get_average_of_function_hessians(
3782 const InputVector &fe_function,
3783 std::vector<solution_hessian_type<typename InputVector::value_type>>
3784 &hessians) const
3785 {
3786 std::vector<typename InputVector::value_type> local_dof_values(
3787 this->fe_interface->n_current_interface_dofs());
3788 this->get_local_dof_values(fe_function, local_dof_values);
3789
3790 get_average_of_function_hessians_from_local_dof_values(local_dof_values,
3791 hessians);
3792 }
3793
3794
3795
3796 template <int dim, int spacedim>
3797 template <class InputVector>
3798 void
3799 Scalar<dim, spacedim>::
3800 get_jump_in_function_third_derivatives_from_local_dof_values(
3801 const InputVector &local_dof_values,
3802 std::vector<
3803 solution_third_derivative_type<typename InputVector::value_type>>
3804 &third_derivatives) const
3805 {
3806 AssertDimension(third_derivatives.size(),
3807 this->fe_interface->n_quadrature_points);
3808
3809 for (const unsigned int dof_index : this->fe_interface->dof_indices())
3810 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3811 {
3812 if (dof_index == 0)
3813 third_derivatives[q_index] = 0.;
3814
3815 third_derivatives[q_index] +=
3816 local_dof_values[dof_index] *
3817 jump_in_third_derivatives(dof_index, q_index);
3818 }
3819 }
3820
3821
3822
3823 template <int dim, int spacedim>
3824 template <class InputVector>
3825 void
3826 Scalar<dim, spacedim>::get_jump_in_function_third_derivatives(
3827 const InputVector &fe_function,
3828 std::vector<
3829 solution_third_derivative_type<typename InputVector::value_type>>
3830 &third_derivatives) const
3831 {
3832 std::vector<typename InputVector::value_type> local_dof_values(
3833 this->fe_interface->n_current_interface_dofs());
3834 this->get_local_dof_values(fe_function, local_dof_values);
3835
3836 get_jump_in_function_third_derivatives_from_local_dof_values(
3837 local_dof_values, third_derivatives);
3838 }
3839
3840
3841
3842 template <int dim, int spacedim>
3844 const FEInterfaceValues<dim, spacedim> &fe_interface,
3845 const unsigned int first_vector_component)
3846 : Base<dim, spacedim>(fe_interface)
3847 , extractor(first_vector_component)
3848 {}
3849
3850
3851
3852 template <int dim, int spacedim>
3854 Vector<dim, spacedim>::value(const bool here_or_there,
3855 const unsigned int interface_dof_index,
3856 const unsigned int q_point) const
3857 {
3858 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3859
3860 if (here_or_there && dof_pair[0] != numbers::invalid_unsigned_int)
3861 return (*(this->fe_interface->fe_face_values))[extractor].value(
3862 dof_pair[0], q_point);
3863
3864 if (!here_or_there && dof_pair[1] != numbers::invalid_unsigned_int)
3865 return (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3866 dof_pair[1], q_point);
3867
3868 return value_type();
3869 }
3870
3871
3872
3873 template <int dim, int spacedim>
3875 Vector<dim, spacedim>::jump_in_values(const unsigned int interface_dof_index,
3876 const unsigned int q_point) const
3877 {
3878 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3879
3880 value_type value;
3881
3882 if (dof_pair[0] != numbers::invalid_unsigned_int)
3883 value +=
3884 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3885 q_point);
3886
3887 if (dof_pair[1] != numbers::invalid_unsigned_int)
3888 value -=
3889 (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3890 dof_pair[1], q_point);
3891
3892 return value;
3893 }
3894
3895
3896
3897 template <int dim, int spacedim>
3900 const unsigned int interface_dof_index,
3901 const unsigned int q_point) const
3902 {
3903 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3904
3905 if (this->fe_interface->at_boundary())
3906 return (*(this->fe_interface->fe_face_values))[extractor].value(
3907 dof_pair[0], q_point);
3908
3909 value_type value;
3910
3911 if (dof_pair[0] != numbers::invalid_unsigned_int)
3912 value +=
3913 0.5 *
3914 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3915 q_point);
3916
3917 if (dof_pair[1] != numbers::invalid_unsigned_int)
3918 value +=
3919 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3920 dof_pair[1], q_point);
3921
3922 return value;
3923 }
3924
3925
3926
3927 template <int dim, int spacedim>
3930 const unsigned int interface_dof_index,
3931 const unsigned int q_point) const
3932 {
3933 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3934
3935 if (this->fe_interface->at_boundary())
3936 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3937 dof_pair[0], q_point);
3938
3939 gradient_type value;
3940
3941 if (dof_pair[0] != numbers::invalid_unsigned_int)
3942 value +=
3943 0.5 *
3944 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3945 q_point);
3946
3947 if (dof_pair[1] != numbers::invalid_unsigned_int)
3948 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3949 .gradient(dof_pair[1], q_point);
3950
3951 return value;
3952 }
3953
3954
3955
3956 template <int dim, int spacedim>
3959 const unsigned int interface_dof_index,
3960 const unsigned int q_point) const
3961 {
3962 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3963
3964 if (this->fe_interface->at_boundary())
3965 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3966 dof_pair[0], q_point);
3967
3968 gradient_type value;
3969
3970 if (dof_pair[0] != numbers::invalid_unsigned_int)
3971 value +=
3972 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3973 q_point);
3974
3975 if (dof_pair[1] != numbers::invalid_unsigned_int)
3976 value -=
3977 (*(this->fe_interface->fe_face_values_neighbor))[extractor].gradient(
3978 dof_pair[1], q_point);
3979
3980 return value;
3981 }
3982
3983
3984
3985 template <int dim, int spacedim>
3987 Vector<dim, spacedim>::jump_gradient(const unsigned int interface_dof_index,
3988 const unsigned int q_point) const
3989 {
3990 return jump_in_gradients(interface_dof_index, q_point);
3991 }
3992
3993
3994
3995 template <int dim, int spacedim>
3998 const unsigned int interface_dof_index,
3999 const unsigned int q_point) const
4000 {
4001 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
4002
4003 if (this->fe_interface->at_boundary())
4004 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
4005 dof_pair[0], q_point);
4006
4007 hessian_type value;
4008
4009 if (dof_pair[0] != numbers::invalid_unsigned_int)
4010 value +=
4011 0.5 *
4012 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
4013 q_point);
4014
4015 if (dof_pair[1] != numbers::invalid_unsigned_int)
4016 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
4017 .hessian(dof_pair[1], q_point);
4018
4019 return value;
4020 }
4021
4022
4023
4024 template <int dim, int spacedim>
4026 Vector<dim, spacedim>::average_hessian(const unsigned int interface_dof_index,
4027 const unsigned int q_point) const
4028 {
4029 return average_of_hessians(interface_dof_index, q_point);
4030 }
4031
4032
4033
4034 template <int dim, int spacedim>
4037 const unsigned int interface_dof_index,
4038 const unsigned int q_point) const
4039 {
4040 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
4041
4042 if (this->fe_interface->at_boundary())
4043 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
4044 dof_pair[0], q_point);
4045
4046 hessian_type value;
4047
4048 if (dof_pair[0] != numbers::invalid_unsigned_int)
4049 value +=
4050 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
4051 q_point);
4052
4053 if (dof_pair[1] != numbers::invalid_unsigned_int)
4054 value -=
4055 (*(this->fe_interface->fe_face_values_neighbor))[extractor].hessian(
4056 dof_pair[1], q_point);
4057
4058 return value;
4059 }
4060
4061
4062
4063 template <int dim, int spacedim>
4065 Vector<dim, spacedim>::jump_hessian(const unsigned int interface_dof_index,
4066 const unsigned int q_point) const
4067 {
4068 return jump_in_hessians(interface_dof_index, q_point);
4069 }
4070
4071
4072
4073 template <int dim, int spacedim>
4076 const unsigned int interface_dof_index,
4077 const unsigned int q_point) const
4078 {
4079 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
4080
4081 if (this->fe_interface->at_boundary())
4082 return (*(this->fe_interface->fe_face_values))[extractor]
4083 .third_derivative(dof_pair[0], q_point);
4084
4085 third_derivative_type value;
4086
4087 if (dof_pair[0] != numbers::invalid_unsigned_int)
4088 value +=
4089 (*(this->fe_interface->fe_face_values))[extractor].third_derivative(
4090 dof_pair[0], q_point);
4091
4092 if (dof_pair[1] != numbers::invalid_unsigned_int)
4093 value -= (*(this->fe_interface->fe_face_values_neighbor))[extractor]
4094 .third_derivative(dof_pair[1], q_point);
4095
4096 return value;
4097 }
4098
4099
4100
4101 template <int dim, int spacedim>
4102 template <class InputVector>
4103 void
4105 const bool here_or_there,
4106 const InputVector &local_dof_values,
4107 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4108 const
4109 {
4110 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
4111
4112 for (const auto dof_index : this->fe_interface->dof_indices())
4113 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4114 {
4115 if (dof_index == 0)
4116 values[q_index] = 0.;
4117
4118 values[q_index] += local_dof_values[dof_index] *
4119 value(here_or_there, dof_index, q_index);
4120 }
4121 }
4122
4123
4124
4125 template <int dim, int spacedim>
4126 template <class InputVector>
4127 void
4129 const bool here_or_there,
4130 const InputVector &fe_function,
4131 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4132 const
4133 {
4134 std::vector<typename InputVector::value_type> local_dof_values(
4135 this->fe_interface->n_current_interface_dofs());
4136 this->get_local_dof_values(fe_function, local_dof_values);
4137
4138 get_function_values_from_local_dof_values(here_or_there,
4139 local_dof_values,
4140 values);
4141 }
4142
4143
4144
4145 template <int dim, int spacedim>
4146 template <class InputVector>
4147 void
4149 const InputVector &local_dof_values,
4150 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4151 const
4152 {
4153 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
4154
4155 for (const auto dof_index : this->fe_interface->dof_indices())
4156 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4157 {
4158 if (dof_index == 0)
4159 values[q_index] = 0.;
4160
4161 values[q_index] +=
4162 local_dof_values[dof_index] * jump_in_values(dof_index, q_index);
4163 }
4164 }
4165
4166
4167
4168 template <int dim, int spacedim>
4169 template <class InputVector>
4170 void
4172 const InputVector &fe_function,
4173 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4174 const
4175 {
4176 std::vector<typename InputVector::value_type> local_dof_values(
4177 this->fe_interface->n_current_interface_dofs());
4178 this->get_local_dof_values(fe_function, local_dof_values);
4179
4180 get_jump_in_function_values_from_local_dof_values(local_dof_values, values);
4181 }
4182
4183
4184
4185 template <int dim, int spacedim>
4186 template <class InputVector>
4187 void
4189 const InputVector &local_dof_values,
4190 std::vector<solution_gradient_type<typename InputVector::value_type>>
4191 &gradients) const
4192 {
4193 AssertDimension(gradients.size(), this->fe_interface->n_quadrature_points);
4194
4195 for (const auto dof_index : this->fe_interface->dof_indices())
4196 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4197 {
4198 if (dof_index == 0)
4199 gradients[q_index] = 0.;
4200
4201 gradients[q_index] +=
4202 local_dof_values[dof_index] * jump_in_gradients(dof_index, q_index);
4203 }
4204 }
4205
4206
4207
4208 template <int dim, int spacedim>
4209 template <class InputVector>
4210 void
4212 const InputVector &fe_function,
4213 std::vector<solution_gradient_type<typename InputVector::value_type>>
4214 &gradients) const
4215 {
4216 std::vector<typename InputVector::value_type> local_dof_values(
4217 this->fe_interface->n_current_interface_dofs());
4218 this->get_local_dof_values(fe_function, local_dof_values);
4219
4220 get_jump_in_function_gradients_from_local_dof_values(local_dof_values,
4221 gradients);
4222 }
4223
4224
4225
4226 template <int dim, int spacedim>
4227 template <class InputVector>
4228 void
4230 const InputVector &local_dof_values,
4231 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4232 const
4233 {
4234 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
4235
4236 for (const auto dof_index : this->fe_interface->dof_indices())
4237 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4238 {
4239 if (dof_index == 0)
4240 values[q_index] = 0.;
4241
4242 values[q_index] +=
4243 local_dof_values[dof_index] * average_of_values(dof_index, q_index);
4244 }
4245 }
4246
4247
4248
4249 template <int dim, int spacedim>
4250 template <class InputVector>
4251 void
4253 const InputVector &fe_function,
4254 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4255 const
4256 {
4257 std::vector<typename InputVector::value_type> local_dof_values(
4258 this->fe_interface->n_current_interface_dofs());
4259 this->get_local_dof_values(fe_function, local_dof_values);
4260
4261 get_average_of_function_values_from_local_dof_values(local_dof_values,
4262 values);
4263 }
4264
4265
4266
4267 template <int dim, int spacedim>
4268 template <class InputVector>
4269 void
4272 const InputVector &local_dof_values,
4273 std::vector<solution_gradient_type<typename InputVector::value_type>>
4274 &gradients) const
4275 {
4276 AssertDimension(gradients.size(), this->fe_interface->n_quadrature_points);
4277
4278 for (const auto dof_index : this->fe_interface->dof_indices())
4279 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4280 {
4281 if (dof_index == 0)
4282 gradients[q_index] = 0.;
4283
4284 gradients[q_index] += local_dof_values[dof_index] *
4285 average_of_gradients(dof_index, q_index);
4286 }
4287 }
4288
4289
4290
4291 template <int dim, int spacedim>
4292 template <class InputVector>
4293 void
4295 const InputVector &fe_function,
4296 std::vector<solution_gradient_type<typename InputVector::value_type>>
4297 &gradients) const
4298 {
4299 std::vector<typename InputVector::value_type> local_dof_values(
4300 this->fe_interface->n_current_interface_dofs());
4301 this->get_local_dof_values(fe_function, local_dof_values);
4302
4303 get_average_of_function_gradients_from_local_dof_values(local_dof_values,
4304 gradients);
4305 }
4306
4307
4308
4309 template <int dim, int spacedim>
4310 template <class InputVector>
4311 void
4313 const InputVector &local_dof_values,
4314 std::vector<solution_hessian_type<typename InputVector::value_type>>
4315 &hessians) const
4316 {
4317 AssertDimension(hessians.size(), this->fe_interface->n_quadrature_points);
4318
4319 for (const auto dof_index : this->fe_interface->dof_indices())
4320 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4321 {
4322 if (dof_index == 0)
4323 hessians[q_index] = 0.;
4324
4325 hessians[q_index] +=
4326 local_dof_values[dof_index] * jump_in_hessians(dof_index, q_index);
4327 }
4328 }
4329
4330
4331
4332 template <int dim, int spacedim>
4333 template <class InputVector>
4334 void
4336 const InputVector &fe_function,
4337 std::vector<solution_hessian_type<typename InputVector::value_type>>
4338 &hessians) const
4339 {
4340 std::vector<typename InputVector::value_type> local_dof_values(
4341 this->fe_interface->n_current_interface_dofs());
4342 this->get_local_dof_values(fe_function, local_dof_values);
4343
4344 get_jump_in_function_hessians_from_local_dof_values(local_dof_values,
4345 hessians);
4346 }
4347
4348
4349
4350 template <int dim, int spacedim>
4351 template <class InputVector>
4352 void
4354 const InputVector &local_dof_values,
4355 std::vector<solution_hessian_type<typename InputVector::value_type>>
4356 &hessians) const
4357 {
4358 AssertDimension(hessians.size(), this->fe_interface->n_quadrature_points);
4359
4360 for (const unsigned int dof_index : this->fe_interface->dof_indices())
4361 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4362 {
4363 if (dof_index == 0)
4364 hessians[q_index] = 0.;
4365
4366 hessians[q_index] += local_dof_values[dof_index] *
4367 average_of_hessians(dof_index, q_index);
4368 }
4369 }
4370
4371
4372
4373 template <int dim, int spacedim>
4374 template <class InputVector>
4375 void
4377 const InputVector &fe_function,
4378 std::vector<solution_hessian_type<typename InputVector::value_type>>
4379 &hessians) const
4380 {
4381 std::vector<typename InputVector::value_type> local_dof_values(
4382 this->fe_interface->n_current_interface_dofs());
4383 this->get_local_dof_values(fe_function, local_dof_values);
4384
4385 get_average_of_function_hessians_from_local_dof_values(local_dof_values,
4386 hessians);
4387 }
4388
4389
4390
4391 template <int dim, int spacedim>
4392 template <class InputVector>
4393 void
4396 const InputVector &local_dof_values,
4397 std::vector<
4398 solution_third_derivative_type<typename InputVector::value_type>>
4399 &third_derivatives) const
4400 {
4401 AssertDimension(third_derivatives.size(),
4402 this->fe_interface->n_quadrature_points);
4403
4404 for (const unsigned int dof_index : this->fe_interface->dof_indices())
4405 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4406 {
4407 if (dof_index == 0)
4408 third_derivatives[q_index] = 0.;
4409
4410 third_derivatives[q_index] +=
4411 local_dof_values[dof_index] *
4412 jump_in_third_derivatives(dof_index, q_index);
4413 }
4414 }
4415
4416
4417
4418 template <int dim, int spacedim>
4419 template <class InputVector>
4420 void
4422 const InputVector &fe_function,
4423 std::vector<
4424 solution_third_derivative_type<typename InputVector::value_type>>
4425 &third_derivatives) const
4426 {
4427 std::vector<typename InputVector::value_type> local_dof_values(
4428 this->fe_interface->n_current_interface_dofs());
4429 this->get_local_dof_values(fe_function, local_dof_values);
4430
4431 get_jump_in_function_third_derivatives_from_local_dof_values(
4432 local_dof_values, third_derivatives);
4433 }
4434} // namespace FEInterfaceViews
4435
4436#endif // DOXYGEN
4437
4439
4440#endif
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell, const unsigned int face_no)
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
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
std::unique_ptr< FESubfaceValues< dim, spacedim > > internal_fe_subface_values_neighbor
std::vector< std::array< unsigned int, 2 > > dofmap
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
FEInterfaceViews::Vector< dim, spacedim > operator[](const FEValuesExtractors::Vector &vector) 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 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
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell, const unsigned int face_no, const unsigned int subface_no)
Abstract base class for mapping classes.
Definition mapping.h:318
friend class Vector
Definition vector.h:1108
Number value_type
Definition vector.h:138
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
unsigned int cell_index
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcOnlyAvailableWithHP()
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:494
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcOnlyAvailableWithoutHP()
static ::ExceptionBase & ExcMessage(std::string arg1)
UpdateFlags
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition mapping.cc:294
typename ::internal::FEInterfaceViews:: ViewType< dim, spacedim, Extractor >::type View
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
Definition hp.h:117
const types::fe_index invalid_fe_index
Definition types.h:243
static const unsigned int invalid_unsigned_int
Definition types.h:220
boost::integer_range< IncrementableType > iota_view
Definition iota_view.h:45
STL namespace.
typename internal::ProductTypeImpl< std::decay_t< T >, std::decay_t< U > >::type type