deal.II version GIT relicensing-2167-g9622207b8f 2024-11-21 12:40:00+00:00
\(\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
// 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
// non-hp data
2111
// 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
// 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 unsigned int active_fe_index = 0;
2343 unsigned int active_fe_index_neighbor = 0;
2344
2345 if (internal_fe_face_values)
2346 {
2347 if (sub_face_no == numbers::invalid_unsigned_int)
2348 {
2349 internal_fe_face_values->reinit(cell, face_no);
2350 fe_face_values = internal_fe_face_values.get();
2351 }
2352 else
2353 {
2354 internal_fe_subface_values->reinit(cell, face_no, sub_face_no);
2355 fe_face_values = internal_fe_subface_values.get();
2356 }
2357 if (sub_face_no_neighbor == numbers::invalid_unsigned_int)
2358 {
2359 internal_fe_face_values_neighbor->reinit(cell_neighbor,
2360 face_no_neighbor);
2361 fe_face_values_neighbor = internal_fe_face_values_neighbor.get();
2362 }
2363 else
2364 {
2365 internal_fe_subface_values_neighbor->reinit(cell_neighbor,
2366 face_no_neighbor,
2367 sub_face_no_neighbor);
2368 fe_face_values_neighbor = internal_fe_subface_values_neighbor.get();
2369 }
2370
2371 AssertDimension(fe_face_values->n_quadrature_points,
2372 fe_face_values_neighbor->n_quadrature_points);
2373
2374 const_cast<unsigned int &>(this->n_quadrature_points) =
2375 fe_face_values->n_quadrature_points;
2376 }
2377 else if (internal_hp_fe_face_values)
2378 {
2379 active_fe_index = fe_index_in;
2380 active_fe_index_neighbor =
2381 (fe_index_neighbor_in != numbers::invalid_unsigned_int) ?
2382 fe_index_neighbor_in :
2383 fe_index_in;
2384
2385 if (active_fe_index == numbers::invalid_unsigned_int)
2386 {
2387 if constexpr (is_dof_cell_accessor)
2388 active_fe_index = cell->active_fe_index();
2389 else
2390 active_fe_index = 0;
2391 }
2392
2393 if (active_fe_index_neighbor == numbers::invalid_unsigned_int)
2394 {
2395 if constexpr (is_dof_cell_accessor_neighbor)
2396 active_fe_index_neighbor = cell_neighbor->active_fe_index();
2397 else
2398 active_fe_index_neighbor = 0;
2399 }
2400
2401 unsigned int used_q_index = q_index;
2402 unsigned int used_mapping_index = mapping_index;
2403
2404 // First check. If there is only one element in a collection, and if none
2405 // had been specified explicitly, then that's clearly the one to take:
2406 if (used_q_index == numbers::invalid_unsigned_int)
2407 if (internal_hp_fe_face_values->get_quadrature_collection().size() == 1)
2408 used_q_index = 0;
2409
2410 if (used_mapping_index == numbers::invalid_unsigned_int)
2411 if (internal_hp_fe_face_values->get_mapping_collection().size() == 1)
2412 used_mapping_index = 0;
2413
2414 // Second check: See if the two quadrature objects are the same, because
2415 // in that case it does not matter which one we use. Unfortunately, we
2416 // currently have no way of testing that two mapping objects are the
2417 // same :-(
2418 if (used_q_index == numbers::invalid_unsigned_int)
2419 if (internal_hp_fe_face_values
2420 ->get_quadrature_collection()[active_fe_index] ==
2421 internal_hp_fe_face_values
2422 ->get_quadrature_collection()[active_fe_index_neighbor])
2423 used_q_index = active_fe_index;
2424
2425 // Third check, if the above did not already suffice. We see if we
2426 // can get somewhere via the dominated's finite element index.
2427 if ((used_q_index == numbers::invalid_unsigned_int) ||
2428 (used_mapping_index == numbers::invalid_unsigned_int))
2429 {
2430 const unsigned int dominated_fe_index =
2431 ((used_q_index == numbers::invalid_unsigned_int) ||
2432 (used_mapping_index == numbers::invalid_unsigned_int) ?
2433 internal_hp_fe_face_values->get_fe_collection()
2434 .find_dominated_fe(
2435 {active_fe_index, active_fe_index_neighbor}) :
2436 numbers::invalid_unsigned_int);
2437
2438 if (used_q_index == numbers::invalid_unsigned_int)
2439 {
2440 Assert(dominated_fe_index != numbers::invalid_fe_index,
2441 ExcMessage(
2442 "You called this function with 'q_index' left at its "
2443 "default value, but this can only work if one of "
2444 "the two finite elements adjacent to this face "
2445 "dominates the other. See the documentation "
2446 "of this function for more information of how "
2447 "to deal with this situation."));
2448 used_q_index = dominated_fe_index;
2449 }
2450
2451 if (used_mapping_index == numbers::invalid_unsigned_int)
2452 {
2453 Assert(dominated_fe_index != numbers::invalid_fe_index,
2454 ExcMessage(
2455 "You called this function with 'mapping_index' left "
2456 "at its default value, but this can only work if one "
2457 "of the two finite elements adjacent to this face "
2458 "dominates the other. See the documentation "
2459 "of this function for more information of how "
2460 "to deal with this situation."));
2461 used_mapping_index = dominated_fe_index;
2462 }
2463 }
2464
2465 // Same as if above, but when hp is enabled.
2466 if (sub_face_no == numbers::invalid_unsigned_int)
2467 {
2468 internal_hp_fe_face_values->reinit(
2469 cell, face_no, used_q_index, used_mapping_index, active_fe_index);
2470 fe_face_values = &const_cast<FEFaceValues<dim, spacedim> &>(
2471 internal_hp_fe_face_values->get_present_fe_values());
2472 }
2473 else
2474 {
2475 internal_hp_fe_subface_values->reinit(cell,
2476 face_no,
2477 sub_face_no,
2478 used_q_index,
2479 used_mapping_index,
2480 active_fe_index);
2481
2482 fe_face_values = &const_cast<FESubfaceValues<dim, spacedim> &>(
2483 internal_hp_fe_subface_values->get_present_fe_values());
2484 }
2485 if (sub_face_no_neighbor == numbers::invalid_unsigned_int)
2486 {
2487 internal_hp_fe_face_values_neighbor->reinit(cell_neighbor,
2488 face_no_neighbor,
2489 used_q_index,
2490 used_mapping_index,
2491 active_fe_index_neighbor);
2492
2493 fe_face_values_neighbor = &const_cast<FEFaceValues<dim, spacedim> &>(
2494 internal_hp_fe_face_values_neighbor->get_present_fe_values());
2495 }
2496 else
2497 {
2498 internal_hp_fe_subface_values_neighbor->reinit(
2499 cell_neighbor,
2500 face_no_neighbor,
2501 sub_face_no_neighbor,
2502 used_q_index,
2503 used_mapping_index,
2504 active_fe_index_neighbor);
2505
2506 fe_face_values_neighbor =
2507 &const_cast<FESubfaceValues<dim, spacedim> &>(
2508 internal_hp_fe_subface_values_neighbor->get_present_fe_values());
2509 }
2510
2511 AssertDimension(fe_face_values->n_quadrature_points,
2512 fe_face_values_neighbor->n_quadrature_points);
2513
2514 const_cast<unsigned int &>(this->n_quadrature_points) =
2515 fe_face_values->n_quadrature_points;
2516 }
2517
2518 // Set up dof mapping and remove duplicates (for continuous elements).
2519 if constexpr (is_dof_cell_accessor_neighbor && is_dof_cell_accessor)
2520 {
2521 // Get dof indices first:
2522 std::vector<types::global_dof_index> v(
2523 fe_face_values->get_fe().n_dofs_per_cell());
2524 cell->get_active_or_mg_dof_indices(v);
2525 std::vector<types::global_dof_index> v2(
2526 fe_face_values_neighbor->get_fe().n_dofs_per_cell());
2527 cell_neighbor->get_active_or_mg_dof_indices(v2);
2528
2529 // Fill a map from the global dof index to the left and right
2530 // local index.
2531 std::map<types::global_dof_index, std::pair<unsigned int, unsigned int>>
2532 tempmap;
2533 std::pair<unsigned int, unsigned int> invalid_entry(
2535
2536 for (unsigned int i = 0; i < v.size(); ++i)
2537 {
2538 // If not already existing, add an invalid entry:
2539 auto result = tempmap.insert(std::make_pair(v[i], invalid_entry));
2540 result.first->second.first = i;
2541 }
2542
2543 for (unsigned int i = 0; i < v2.size(); ++i)
2544 {
2545 // If not already existing, add an invalid entry:
2546 auto result = tempmap.insert(std::make_pair(v2[i], invalid_entry));
2547 result.first->second.second = i;
2548 }
2549
2550 // Transfer from the map to the sorted std::vectors.
2551 dofmap.resize(tempmap.size());
2552 interface_dof_indices.resize(tempmap.size());
2553 unsigned int idx = 0;
2554 for (auto &x : tempmap)
2555 {
2556 interface_dof_indices[idx] = x.first;
2557 dofmap[idx] = {{x.second.first, x.second.second}};
2558 ++idx;
2559 }
2560 }
2561 else
2562 {
2563 const unsigned int n_dofs_per_cell_1 = fe_face_values->dofs_per_cell;
2564 const unsigned int n_dofs_per_cell_2 =
2565 fe_face_values_neighbor->dofs_per_cell;
2566
2567 interface_dof_indices.resize(n_dofs_per_cell_1 + n_dofs_per_cell_2);
2568 dofmap.resize(n_dofs_per_cell_1 + n_dofs_per_cell_2);
2569
2570 for (unsigned int i = 0; i < n_dofs_per_cell_1; ++i)
2571 {
2572 interface_dof_indices[i] = numbers::invalid_dof_index;
2573 dofmap[i] = {{i, numbers::invalid_unsigned_int}};
2574 }
2575
2576 for (unsigned int i = 0; i < n_dofs_per_cell_2; ++i)
2577 {
2578 interface_dof_indices[i + n_dofs_per_cell_1] =
2580 dofmap[i + n_dofs_per_cell_1] = {{numbers::invalid_unsigned_int, i}};
2581 }
2582 }
2583}
2584
2585
2586
2587template <int dim, int spacedim>
2588template <typename CellIteratorType>
2589void
2590FEInterfaceValues<dim, spacedim>::reinit(const CellIteratorType &cell,
2591 const unsigned int face_no,
2592 const unsigned int q_index,
2593 const unsigned int mapping_index,
2594 const unsigned int fe_index)
2595{
2596 Assert(internal_fe_face_values || internal_hp_fe_face_values,
2598
2599 if (internal_fe_face_values)
2600 {
2601 Assert((q_index == 0 || q_index == numbers::invalid_unsigned_int),
2603 Assert((mapping_index == 0 ||
2604 mapping_index == numbers::invalid_unsigned_int),
2606 Assert((fe_index == 0 || fe_index == numbers::invalid_unsigned_int),
2608
2609 internal_fe_face_values->reinit(cell, face_no);
2610 fe_face_values = internal_fe_face_values.get();
2611 fe_face_values_neighbor = nullptr;
2612 }
2613 else if (internal_hp_fe_face_values)
2614 {
2615 internal_hp_fe_face_values->reinit(
2616 cell, face_no, q_index, mapping_index, fe_index);
2617 fe_face_values = &const_cast<FEFaceValues<dim> &>(
2618 internal_hp_fe_face_values->get_present_fe_values());
2619 fe_face_values_neighbor = nullptr;
2620 }
2621
2622 interface_dof_indices.resize(fe_face_values->get_fe().n_dofs_per_cell());
2623
2624 if constexpr (std::is_same_v<typename CellIteratorType::AccessorType,
2626 std::is_same_v<typename CellIteratorType::AccessorType,
2628 {
2629 cell->get_active_or_mg_dof_indices(interface_dof_indices);
2630 }
2631 else
2632 {
2633 for (auto &i : interface_dof_indices)
2635 }
2636
2637 dofmap.resize(interface_dof_indices.size());
2638 for (unsigned int i = 0; i < interface_dof_indices.size(); ++i)
2639 dofmap[i] = {{i, numbers::invalid_unsigned_int}};
2640}
2641
2642
2643
2644template <int dim, int spacedim>
2645inline double
2646FEInterfaceValues<dim, spacedim>::JxW(const unsigned int q) const
2647{
2648 Assert(fe_face_values != nullptr,
2649 ExcMessage("This call requires a call to reinit() first."));
2650 return fe_face_values->JxW(q);
2651}
2652
2653
2654
2655template <int dim, int spacedim>
2656const std::vector<double> &
2658{
2659 Assert(fe_face_values != nullptr,
2660 ExcMessage("This call requires a call to reinit() first."));
2661 return fe_face_values->get_JxW_values();
2662}
2663
2664
2665
2666template <int dim, int spacedim>
2667const std::vector<Tensor<1, spacedim>> &
2669{
2670 Assert(fe_face_values != nullptr,
2671 ExcMessage("This call requires a call to reinit() first."));
2672 return fe_face_values->get_normal_vectors();
2673}
2674
2675
2676
2677template <int dim, int spacedim>
2680{
2681 Assert(!has_hp_capabilities(), ExcOnlyAvailableWithoutHP());
2682 return internal_fe_face_values->get_mapping();
2683}
2684
2685
2686
2687template <int dim, int spacedim>
2690{
2691 Assert(!has_hp_capabilities(), ExcOnlyAvailableWithoutHP());
2692 return internal_fe_face_values->get_fe();
2693}
2694
2695
2696
2697template <int dim, int spacedim>
2698const Quadrature<dim - 1> &
2700{
2701 Assert(!has_hp_capabilities(), ExcOnlyAvailableWithoutHP());
2702 return internal_fe_face_values->get_quadrature();
2703}
2704
2705
2706
2707template <int dim, int spacedim>
2710{
2711 Assert(has_hp_capabilities(), ExcOnlyAvailableWithHP());
2712 return internal_hp_fe_face_values->get_mapping_collection();
2713}
2714
2715
2716
2717template <int dim, int spacedim>
2720{
2721 Assert(has_hp_capabilities(), ExcOnlyAvailableWithHP());
2722 return internal_hp_fe_face_values->get_fe_collection();
2723}
2724
2725
2726
2727template <int dim, int spacedim>
2728const hp::QCollection<dim - 1> &
2730{
2731 Assert(has_hp_capabilities(), ExcOnlyAvailableWithHP());
2732 return internal_hp_fe_face_values->get_quadrature_collection();
2733}
2734
2735
2736
2737template <int dim, int spacedim>
2738bool
2740{
2741 if (internal_hp_fe_face_values || internal_hp_fe_subface_values ||
2742 internal_hp_fe_face_values_neighbor ||
2743 internal_hp_fe_subface_values_neighbor)
2744 {
2745 Assert(!internal_fe_face_values, ExcInternalError());
2746 Assert(!internal_fe_subface_values, ExcInternalError());
2747 Assert(!internal_fe_face_values_neighbor, ExcInternalError());
2748 Assert(!internal_fe_subface_values_neighbor, ExcInternalError());
2749
2750 return true;
2751 }
2752
2753 Assert(internal_fe_face_values || internal_fe_subface_values ||
2754 internal_fe_face_values_neighbor ||
2755 internal_fe_subface_values_neighbor,
2757 Assert(!internal_hp_fe_face_values, ExcInternalError());
2758 Assert(!internal_hp_fe_subface_values, ExcInternalError());
2759 Assert(!internal_hp_fe_face_values_neighbor, ExcInternalError());
2760 Assert(!internal_hp_fe_subface_values_neighbor, ExcInternalError());
2761
2762 return false;
2763}
2764
2765
2766
2767template <int dim, int spacedim>
2770{
2771 return {0U, n_quadrature_points};
2772}
2773
2774
2775
2776template <int dim, int spacedim>
2777const std::vector<Point<spacedim>> &
2779{
2780 Assert(fe_face_values != nullptr,
2781 ExcMessage("This call requires a call to reinit() first."));
2782 return fe_face_values->get_quadrature_points();
2783}
2784
2785
2786
2787template <int dim, int spacedim>
2790{
2791 if (has_hp_capabilities())
2792 return internal_hp_fe_face_values->get_update_flags();
2793 else
2794 return internal_fe_face_values->get_update_flags();
2795}
2796
2797
2798
2799template <int dim, int spacedim>
2802{
2803 return get_fe_face_values(cell_index).get_cell();
2804}
2805
2806
2807
2808template <int dim, int spacedim>
2809inline unsigned int
2811 const unsigned int cell_index) const
2812{
2813 return get_fe_face_values(cell_index).get_face_number();
2814}
2815
2816
2817
2818template <int dim, int spacedim>
2819unsigned
2821{
2822 Assert(
2823 interface_dof_indices.size() > 0,
2824 ExcMessage(
2825 "n_current_interface_dofs() is only available after a call to reinit()."));
2826 return interface_dof_indices.size();
2827}
2828
2829
2830
2831template <int dim, int spacedim>
2834{
2835 return {0U, n_current_interface_dofs()};
2836}
2837
2838
2839
2840template <int dim, int spacedim>
2841bool
2843{
2844 return fe_face_values_neighbor == nullptr;
2845}
2846
2847
2848
2849template <int dim, int spacedim>
2850std::vector<types::global_dof_index>
2852{
2853 return interface_dof_indices;
2854}
2855
2856
2857
2858template <int dim, int spacedim>
2859std::array<unsigned int, 2>
2861 const unsigned int interface_dof_index) const
2862{
2863 AssertIndexRange(interface_dof_index, n_current_interface_dofs());
2864 return dofmap[interface_dof_index];
2865}
2866
2867
2868
2869template <int dim, int spacedim>
2872 const unsigned int cell_index) const
2873{
2875 Assert(
2876 cell_index == 0 || !at_boundary(),
2877 ExcMessage(
2878 "You are on a boundary, so you can only ask for the first FEFaceValues object."));
2879
2880 return (cell_index == 0) ? *fe_face_values : *fe_face_values_neighbor;
2881}
2882
2883
2884
2885template <int dim, int spacedim>
2887FEInterfaceValues<dim, spacedim>::normal(const unsigned int q_point_index) const
2888{
2889 return fe_face_values->normal_vector(q_point_index);
2890}
2891
2892
2893
2894template <int dim, int spacedim>
2895double
2897 const bool here_or_there,
2898 const unsigned int interface_dof_index,
2899 const unsigned int q_point,
2900 const unsigned int component) const
2901{
2902 const auto dof_pair = dofmap[interface_dof_index];
2903
2904 if (here_or_there && dof_pair[0] != numbers::invalid_unsigned_int)
2905 return get_fe_face_values(0).shape_value_component(dof_pair[0],
2906 q_point,
2907 component);
2908 if (!here_or_there && dof_pair[1] != numbers::invalid_unsigned_int)
2909 return get_fe_face_values(1).shape_value_component(dof_pair[1],
2910 q_point,
2911 component);
2912
2913 return 0.0;
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 double value = 0.0;
2928
2929 if (dof_pair[0] != numbers::invalid_unsigned_int)
2930 value += get_fe_face_values(0).shape_value_component(dof_pair[0],
2931 q_point,
2932 component);
2933 if (dof_pair[1] != numbers::invalid_unsigned_int)
2934 value -= get_fe_face_values(1).shape_value_component(dof_pair[1],
2935 q_point,
2936 component);
2937 return value;
2938}
2939
2940
2941
2942template <int dim, int spacedim>
2943double
2945 const unsigned int interface_dof_index,
2946 const unsigned int q_point,
2947 const unsigned int component) const
2948{
2949 const auto dof_pair = dofmap[interface_dof_index];
2950
2951 if (at_boundary())
2952 return get_fe_face_values(0).shape_value_component(dof_pair[0],
2953 q_point,
2954 component);
2955
2956 double value = 0.0;
2957
2958 if (dof_pair[0] != numbers::invalid_unsigned_int)
2959 value += 0.5 * get_fe_face_values(0).shape_value_component(dof_pair[0],
2960 q_point,
2961 component);
2962 if (dof_pair[1] != numbers::invalid_unsigned_int)
2963 value += 0.5 * get_fe_face_values(1).shape_value_component(dof_pair[1],
2964 q_point,
2965 component);
2966
2967 return value;
2968}
2969
2970
2971
2972template <int dim, int spacedim>
2975 const unsigned int interface_dof_index,
2976 const unsigned int q_point,
2977 const unsigned int component) const
2978{
2979 const auto dof_pair = dofmap[interface_dof_index];
2980
2981 if (at_boundary())
2982 return get_fe_face_values(0).shape_grad_component(dof_pair[0],
2983 q_point,
2984 component);
2985
2987
2988 if (dof_pair[0] != numbers::invalid_unsigned_int)
2989 value += 0.5 * get_fe_face_values(0).shape_grad_component(dof_pair[0],
2990 q_point,
2991 component);
2992 if (dof_pair[1] != numbers::invalid_unsigned_int)
2993 value += 0.5 * get_fe_face_values(1).shape_grad_component(dof_pair[1],
2994 q_point,
2995 component);
2996
2997 return value;
2998}
2999
3000
3001
3002template <int dim, int spacedim>
3005 const unsigned int interface_dof_index,
3006 const unsigned int q_point,
3007 const unsigned int component) const
3008{
3009 const auto dof_pair = dofmap[interface_dof_index];
3010
3011 if (at_boundary())
3012 return get_fe_face_values(0).shape_hessian_component(dof_pair[0],
3013 q_point,
3014 component);
3015
3017
3018 if (dof_pair[0] != numbers::invalid_unsigned_int)
3019 value += 0.5 * get_fe_face_values(0).shape_hessian_component(dof_pair[0],
3020 q_point,
3021 component);
3022 if (dof_pair[1] != numbers::invalid_unsigned_int)
3023 value += 0.5 * get_fe_face_values(1).shape_hessian_component(dof_pair[1],
3024 q_point,
3025 component);
3026
3027 return value;
3028}
3029
3030
3031
3032template <int dim, int spacedim>
3035 const unsigned int interface_dof_index,
3036 const unsigned int q_point,
3037 const unsigned int component) const
3038{
3039 const auto dof_pair = dofmap[interface_dof_index];
3040
3041 if (at_boundary())
3042 return get_fe_face_values(0).shape_grad_component(dof_pair[0],
3043 q_point,
3044 component);
3045
3047
3048 if (dof_pair[0] != numbers::invalid_unsigned_int)
3049 value += get_fe_face_values(0).shape_grad_component(dof_pair[0],
3050 q_point,
3051 component);
3052 if (dof_pair[1] != numbers::invalid_unsigned_int)
3053 value -= get_fe_face_values(1).shape_grad_component(dof_pair[1],
3054 q_point,
3055 component);
3056
3057 return value;
3058}
3059
3060
3061
3062template <int dim, int spacedim>
3065 const unsigned int interface_dof_index,
3066 const unsigned int q_point,
3067 const unsigned int component) const
3068{
3069 const auto dof_pair = dofmap[interface_dof_index];
3070
3071 if (at_boundary())
3072 return get_fe_face_values(0).shape_hessian_component(dof_pair[0],
3073 q_point,
3074 component);
3075
3077
3078 if (dof_pair[0] != numbers::invalid_unsigned_int)
3079 value += get_fe_face_values(0).shape_hessian_component(dof_pair[0],
3080 q_point,
3081 component);
3082 if (dof_pair[1] != numbers::invalid_unsigned_int)
3083 value -= get_fe_face_values(1).shape_hessian_component(dof_pair[1],
3084 q_point,
3085 component);
3086
3087 return value;
3088}
3089
3090
3091
3092template <int dim, int spacedim>
3095 const unsigned int interface_dof_index,
3096 const unsigned int q_point,
3097 const unsigned int component) const
3098{
3099 const auto dof_pair = dofmap[interface_dof_index];
3100
3101 if (at_boundary())
3102 return get_fe_face_values(0).shape_3rd_derivative_component(dof_pair[0],
3103 q_point,
3104 component);
3105
3107
3108 if (dof_pair[0] != numbers::invalid_unsigned_int)
3109 value += get_fe_face_values(0).shape_3rd_derivative_component(dof_pair[0],
3110 q_point,
3111 component);
3112 if (dof_pair[1] != numbers::invalid_unsigned_int)
3113 value -= get_fe_face_values(1).shape_3rd_derivative_component(dof_pair[1],
3114 q_point,
3115 component);
3116
3117 return value;
3118}
3119
3120
3121
3122template <int dim, int spacedim>
3123template <class InputVector>
3124void
3126 const InputVector &fe_function,
3127 std::vector<typename InputVector::value_type> &values) const
3128{
3129 AssertDimension(values.size(), n_quadrature_points);
3130
3132 this->operator[](scalar).get_jump_in_function_values(fe_function, values);
3133}
3134
3135
3136
3137template <int dim, int spacedim>
3138template <class InputVector>
3139void
3141 const InputVector &fe_function,
3143 const
3144{
3145 AssertDimension(gradients.size(), n_quadrature_points);
3146
3148 this->operator[](scalar).get_jump_in_function_gradients(fe_function,
3149 gradients);
3150}
3151
3152
3153
3154template <int dim, int spacedim>
3155template <class InputVector>
3156void
3158 const InputVector &fe_function,
3160 const
3161{
3162 AssertDimension(hessians.size(), n_quadrature_points);
3163
3165 this->operator[](scalar).get_jump_in_function_hessians(fe_function, hessians);
3166}
3167
3168
3169
3170template <int dim, int spacedim>
3171template <class InputVector>
3172void
3174 const InputVector &fe_function,
3176 &third_derivatives) const
3177{
3178 AssertDimension(third_derivatives.size(), n_quadrature_points);
3179
3181 this->operator[](scalar).get_jump_in_function_third_derivatives(
3182 fe_function, third_derivatives);
3183}
3184
3185
3186
3187template <int dim, int spacedim>
3188template <class InputVector>
3189void
3191 const InputVector &fe_function,
3192 std::vector<typename InputVector::value_type> &values) const
3193{
3194 AssertDimension(values.size(), n_quadrature_points);
3195
3197 this->operator[](scalar).get_average_of_function_values(fe_function, values);
3198}
3199
3200
3201
3202template <int dim, int spacedim>
3203template <class InputVector>
3204void
3206 const InputVector &fe_function,
3208 const
3209{
3210 AssertDimension(gradients.size(), n_quadrature_points);
3211
3213 this->operator[](scalar).get_average_of_function_gradients(fe_function,
3214 gradients);
3215}
3216
3217
3218
3219template <int dim, int spacedim>
3220template <class InputVector>
3221void
3223 const InputVector &fe_function,
3225 const
3226{
3227 AssertDimension(hessians.size(), n_quadrature_points);
3228
3230 this->operator[](scalar).get_average_of_function_hessians(fe_function,
3231 hessians);
3232}
3233
3234
3235
3236/*------------ Inline functions: FEInterfaceValues------------*/
3237template <int dim, int spacedim>
3240 const FEValuesExtractors::Scalar &scalar) const
3241{
3242 const unsigned int n_components =
3243 (this->has_hp_capabilities() ? this->get_fe_collection().n_components() :
3244 this->get_fe().n_components());
3245 (void)n_components;
3246 AssertIndexRange(scalar.component, n_components);
3247 return FEInterfaceViews::Scalar<dim, spacedim>(*this, scalar.component);
3248}
3249
3250
3251
3252template <int dim, int spacedim>
3255 const FEValuesExtractors::Vector &vector) const
3256{
3257 const unsigned int n_components =
3258 (this->has_hp_capabilities() ? this->get_fe_collection().n_components() :
3259 this->get_fe().n_components());
3260 const unsigned int n_vectors =
3263 0);
3264 (void)n_components;
3265 (void)n_vectors;
3266 AssertIndexRange(vector.first_vector_component, n_vectors);
3268 vector.first_vector_component);
3269}
3270
3271
3272
3273namespace FEInterfaceViews
3274{
3275 template <int dim, int spacedim>
3277 const FEInterfaceValues<dim, spacedim> &fe_interface)
3278 : fe_interface(&fe_interface)
3279 {}
3280
3281
3282
3283 template <int dim, int spacedim>
3284 Scalar<dim, spacedim>::Scalar(
3285 const FEInterfaceValues<dim, spacedim> &fe_interface,
3286 const unsigned int component)
3287 : Base<dim, spacedim>(fe_interface)
3288 , extractor(component)
3289 {}
3290
3291
3292
3293 template <int dim, int spacedim>
3294 template <class InputVector, class OutputVector>
3295 void
3296 Base<dim, spacedim>::get_local_dof_values(
3297 const InputVector &dof_values,
3298 OutputVector &local_dof_values) const
3299 {
3300 const auto &interface_dof_indices =
3301 this->fe_interface->get_interface_dof_indices();
3302
3303 AssertDimension(interface_dof_indices.size(), local_dof_values.size());
3304
3305 for (const unsigned int i : this->fe_interface->dof_indices())
3306 local_dof_values[i] = dof_values(interface_dof_indices[i]);
3307 }
3308
3309
3310
3311 template <int dim, int spacedim>
3312 typename Scalar<dim, spacedim>::value_type
3313 Scalar<dim, spacedim>::value(const bool here_or_there,
3314 const unsigned int interface_dof_index,
3315 const unsigned int q_point) const
3316 {
3317 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3318
3319 if (here_or_there && dof_pair[0] != numbers::invalid_unsigned_int)
3320 return (*(this->fe_interface->fe_face_values))[extractor].value(
3321 dof_pair[0], q_point);
3322
3323 if (!here_or_there && dof_pair[1] != numbers::invalid_unsigned_int)
3324 return (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3325 dof_pair[1], q_point);
3326
3327 return 0.0;
3328 }
3329
3330
3331
3332 template <int dim, int spacedim>
3333 typename Scalar<dim, spacedim>::value_type
3334 Scalar<dim, spacedim>::jump_in_values(const unsigned int interface_dof_index,
3335 const unsigned int q_point) const
3336 {
3337 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3338
3339 value_type value = 0.0;
3340
3341 if (dof_pair[0] != numbers::invalid_unsigned_int)
3342 value +=
3343 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3344 q_point);
3345
3346 if (dof_pair[1] != numbers::invalid_unsigned_int)
3347 value -=
3348 (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3349 dof_pair[1], q_point);
3350
3351 return value;
3352 }
3353
3354
3355
3356 template <int dim, int spacedim>
3357 typename Scalar<dim, spacedim>::value_type
3358 Scalar<dim, spacedim>::average_of_values(
3359 const unsigned int interface_dof_index,
3360 const unsigned int q_point) const
3361 {
3362 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3363
3364 if (this->fe_interface->at_boundary())
3365 return (*(this->fe_interface->fe_face_values))[extractor].value(
3366 dof_pair[0], q_point);
3367
3368 value_type value = 0.0;
3369
3370 if (dof_pair[0] != numbers::invalid_unsigned_int)
3371 value +=
3372 0.5 *
3373 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3374 q_point);
3375
3376 if (dof_pair[1] != numbers::invalid_unsigned_int)
3377 value +=
3378 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3379 dof_pair[1], q_point);
3380
3381 return value;
3382 }
3383
3384
3385
3386 template <int dim, int spacedim>
3387 typename Scalar<dim, spacedim>::gradient_type
3388 Scalar<dim, spacedim>::average_of_gradients(
3389 const unsigned int interface_dof_index,
3390 const unsigned int q_point) const
3391 {
3392 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3393
3394 if (this->fe_interface->at_boundary())
3395 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3396 dof_pair[0], q_point);
3397
3398 gradient_type value;
3399
3400 if (dof_pair[0] != numbers::invalid_unsigned_int)
3401 value +=
3402 0.5 *
3403 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3404 q_point);
3405
3406 if (dof_pair[1] != numbers::invalid_unsigned_int)
3407 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3408 .gradient(dof_pair[1], q_point);
3409
3410 return value;
3411 }
3412
3413
3414
3415 template <int dim, int spacedim>
3416 typename Scalar<dim, spacedim>::gradient_type
3417 Scalar<dim, spacedim>::jump_in_gradients(
3418 const unsigned int interface_dof_index,
3419 const unsigned int q_point) const
3420 {
3421 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3422
3423 if (this->fe_interface->at_boundary())
3424 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3425 dof_pair[0], q_point);
3426
3427 gradient_type value;
3428
3429 if (dof_pair[0] != numbers::invalid_unsigned_int)
3430 value +=
3431 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3432 q_point);
3433
3434 if (dof_pair[1] != numbers::invalid_unsigned_int)
3435 value -=
3436 (*(this->fe_interface->fe_face_values_neighbor))[extractor].gradient(
3437 dof_pair[1], q_point);
3438
3439 return value;
3440 }
3441
3442
3443
3444 template <int dim, int spacedim>
3445 typename Scalar<dim, spacedim>::hessian_type
3446 Scalar<dim, spacedim>::average_of_hessians(
3447 const unsigned int interface_dof_index,
3448 const unsigned int q_point) const
3449 {
3450 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3451
3452 if (this->fe_interface->at_boundary())
3453 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
3454 dof_pair[0], q_point);
3455
3456 hessian_type value;
3457
3458 if (dof_pair[0] != numbers::invalid_unsigned_int)
3459 value +=
3460 0.5 *
3461 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
3462 q_point);
3463
3464 if (dof_pair[1] != numbers::invalid_unsigned_int)
3465 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3466 .hessian(dof_pair[1], q_point);
3467
3468 return value;
3469 }
3470
3471
3472
3473 template <int dim, int spacedim>
3474 typename Scalar<dim, spacedim>::third_derivative_type
3475 Scalar<dim, spacedim>::jump_in_third_derivatives(
3476 const unsigned int interface_dof_index,
3477 const unsigned int q_point) const
3478 {
3479 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3480
3481 if (this->fe_interface->at_boundary())
3482 return (*(this->fe_interface->fe_face_values))[extractor]
3483 .third_derivative(dof_pair[0], q_point);
3484
3485 third_derivative_type value;
3486
3487 if (dof_pair[0] != numbers::invalid_unsigned_int)
3488 value +=
3489 (*(this->fe_interface->fe_face_values))[extractor].third_derivative(
3490 dof_pair[0], q_point);
3491
3492 if (dof_pair[1] != numbers::invalid_unsigned_int)
3493 value -= (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3494 .third_derivative(dof_pair[1], q_point);
3495
3496 return value;
3497 }
3498
3499
3500
3501 template <int dim, int spacedim>
3502 typename Scalar<dim, spacedim>::hessian_type
3503 Scalar<dim, spacedim>::jump_in_hessians(
3504 const unsigned int interface_dof_index,
3505 const unsigned int q_point) const
3506 {
3507 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3508
3509 if (this->fe_interface->at_boundary())
3510 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
3511 dof_pair[0], q_point);
3512
3513 hessian_type value;
3514
3515 if (dof_pair[0] != numbers::invalid_unsigned_int)
3516 value +=
3517 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
3518 q_point);
3519
3520 if (dof_pair[1] != numbers::invalid_unsigned_int)
3521 value -=
3522 (*(this->fe_interface->fe_face_values_neighbor))[extractor].hessian(
3523 dof_pair[1], q_point);
3524
3525 return value;
3526 }
3527
3528
3529
3530 template <int dim, int spacedim>
3531 template <class InputVector>
3532 void
3533 Scalar<dim, spacedim>::get_function_values_from_local_dof_values(
3534 const bool here_or_there,
3535 const InputVector &local_dof_values,
3536 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3537 const
3538 {
3539 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
3540
3541 for (const auto dof_index : this->fe_interface->dof_indices())
3542 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3543 {
3544 if (dof_index == 0)
3545 values[q_index] = 0.;
3546
3547 values[q_index] += local_dof_values[dof_index] *
3548 value(here_or_there, dof_index, q_index);
3549 }
3550 }
3551
3552
3553
3554 template <int dim, int spacedim>
3555 template <class InputVector>
3556 void
3557 Scalar<dim, spacedim>::get_function_values(
3558 const bool here_or_there,
3559 const InputVector &fe_function,
3560 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3561 const
3562 {
3563 std::vector<typename InputVector::value_type> local_dof_values(
3564 this->fe_interface->n_current_interface_dofs());
3565 this->get_local_dof_values(fe_function, local_dof_values);
3566
3567 get_function_values_from_local_dof_values(here_or_there,
3568 local_dof_values,
3569 values);
3570 }
3571
3572
3573
3574 template <int dim, int spacedim>
3575 template <class InputVector>
3576 void
3577 Scalar<dim, spacedim>::get_jump_in_function_values_from_local_dof_values(
3578 const InputVector &local_dof_values,
3579 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3580 const
3581 {
3582 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
3583
3584 for (const auto dof_index : this->fe_interface->dof_indices())
3585 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3586 {
3587 if (dof_index == 0)
3588 values[q_index] = 0.;
3589
3590 values[q_index] +=
3591 local_dof_values[dof_index] * jump_in_values(dof_index, q_index);
3592 }
3593 }
3594
3595
3596
3597 template <int dim, int spacedim>
3598 template <class InputVector>
3599 void
3600 Scalar<dim, spacedim>::get_jump_in_function_values(
3601 const InputVector &fe_function,
3602 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3603 const
3604 {
3605 std::vector<typename InputVector::value_type> local_dof_values(
3606 this->fe_interface->n_current_interface_dofs());
3607 this->get_local_dof_values(fe_function, local_dof_values);
3608
3609 get_jump_in_function_values_from_local_dof_values(local_dof_values, values);
3610 }
3611
3612
3613
3614 template <int dim, int spacedim>
3615 template <class InputVector>
3616 void
3617 Scalar<dim, spacedim>::get_jump_in_function_gradients_from_local_dof_values(
3618 const InputVector &local_dof_values,
3619 std::vector<solution_gradient_type<typename InputVector::value_type>>
3620 &gradients) const
3621 {
3622 AssertDimension(gradients.size(), this->fe_interface->n_quadrature_points);
3623
3624 for (const auto dof_index : this->fe_interface->dof_indices())
3625 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3626 {
3627 if (dof_index == 0)
3628 gradients[q_index] = 0.;
3629
3630 gradients[q_index] +=
3631 local_dof_values[dof_index] * jump_in_gradients(dof_index, q_index);
3632 }
3633 }
3634
3635
3636
3637 template <int dim, int spacedim>
3638 template <class InputVector>
3639 void
3640 Scalar<dim, spacedim>::get_jump_in_function_gradients(
3641 const InputVector &fe_function,
3642 std::vector<solution_gradient_type<typename InputVector::value_type>>
3643 &gradients) const
3644 {
3645 std::vector<typename InputVector::value_type> local_dof_values(
3646 this->fe_interface->n_current_interface_dofs());
3647 this->get_local_dof_values(fe_function, local_dof_values);
3648
3649 get_jump_in_function_gradients_from_local_dof_values(local_dof_values,
3650 gradients);
3651 }
3652
3653
3654
3655 template <int dim, int spacedim>
3656 template <class InputVector>
3657 void
3658 Scalar<dim, spacedim>::get_average_of_function_values_from_local_dof_values(
3659 const InputVector &local_dof_values,
3660 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3661 const
3662 {
3663 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
3664
3665 for (const auto dof_index : this->fe_interface->dof_indices())
3666 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3667 {
3668 if (dof_index == 0)
3669 values[q_index] = 0.;
3670
3671 values[q_index] +=
3672 local_dof_values[dof_index] * average_of_values(dof_index, q_index);
3673 }
3674 }
3675
3676
3677
3678 template <int dim, int spacedim>
3679 template <class InputVector>
3680 void
3681 Scalar<dim, spacedim>::get_average_of_function_values(
3682 const InputVector &fe_function,
3683 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3684 const
3685 {
3686 std::vector<typename InputVector::value_type> local_dof_values(
3687 this->fe_interface->n_current_interface_dofs());
3688 this->get_local_dof_values(fe_function, local_dof_values);
3689
3690 get_average_of_function_values_from_local_dof_values(local_dof_values,
3691 values);
3692 }
3693
3694
3695
3696 template <int dim, int spacedim>
3697 template <class InputVector>
3698 void
3699 Scalar<dim, spacedim>::
3700 get_average_of_function_gradients_from_local_dof_values(
3701 const InputVector &local_dof_values,
3702 std::vector<solution_gradient_type<typename InputVector::value_type>>
3703 &gradients) const
3704 {
3705 AssertDimension(gradients.size(), this->fe_interface->n_quadrature_points);
3706
3707 for (const auto dof_index : this->fe_interface->dof_indices())
3708 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3709 {
3710 if (dof_index == 0)
3711 gradients[q_index] = 0.;
3712
3713 gradients[q_index] += local_dof_values[dof_index] *
3714 average_of_gradients(dof_index, q_index);
3715 }
3716 }
3717
3718
3719
3720 template <int dim, int spacedim>
3721 template <class InputVector>
3722 void
3723 Scalar<dim, spacedim>::get_average_of_function_gradients(
3724 const InputVector &fe_function,
3725 std::vector<solution_gradient_type<typename InputVector::value_type>>
3726 &gradients) const
3727 {
3728 std::vector<typename InputVector::value_type> local_dof_values(
3729 this->fe_interface->n_current_interface_dofs());
3730 this->get_local_dof_values(fe_function, local_dof_values);
3731
3732 get_average_of_function_gradients_from_local_dof_values(local_dof_values,
3733 gradients);
3734 }
3735
3736
3737
3738 template <int dim, int spacedim>
3739 template <class InputVector>
3740 void
3741 Scalar<dim, spacedim>::get_jump_in_function_hessians_from_local_dof_values(
3742 const InputVector &local_dof_values,
3743 std::vector<solution_hessian_type<typename InputVector::value_type>>
3744 &hessians) const
3745 {
3746 AssertDimension(hessians.size(), this->fe_interface->n_quadrature_points);
3747
3748 for (const auto dof_index : this->fe_interface->dof_indices())
3749 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3750 {
3751 if (dof_index == 0)
3752 hessians[q_index] = 0.;
3753
3754 hessians[q_index] +=
3755 local_dof_values[dof_index] * jump_in_hessians(dof_index, q_index);
3756 }
3757 }
3758
3759
3760
3761 template <int dim, int spacedim>
3762 template <class InputVector>
3763 void
3764 Scalar<dim, spacedim>::get_jump_in_function_hessians(
3765 const InputVector &fe_function,
3766 std::vector<solution_hessian_type<typename InputVector::value_type>>
3767 &hessians) const
3768 {
3769 std::vector<typename InputVector::value_type> local_dof_values(
3770 this->fe_interface->n_current_interface_dofs());
3771 this->get_local_dof_values(fe_function, local_dof_values);
3772
3773 get_jump_in_function_hessians_from_local_dof_values(local_dof_values,
3774 hessians);
3775 }
3776
3777
3778
3779 template <int dim, int spacedim>
3780 template <class InputVector>
3781 void
3782 Scalar<dim, spacedim>::get_average_of_function_hessians_from_local_dof_values(
3783 const InputVector &local_dof_values,
3784 std::vector<solution_hessian_type<typename InputVector::value_type>>
3785 &hessians) const
3786 {
3787 AssertDimension(hessians.size(), this->fe_interface->n_quadrature_points);
3788
3789 for (const unsigned int dof_index : this->fe_interface->dof_indices())
3790 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3791 {
3792 if (dof_index == 0)
3793 hessians[q_index] = 0.;
3794
3795 hessians[q_index] += local_dof_values[dof_index] *
3796 average_of_hessians(dof_index, q_index);
3797 }
3798 }
3799
3800
3801
3802 template <int dim, int spacedim>
3803 template <class InputVector>
3804 void
3805 Scalar<dim, spacedim>::get_average_of_function_hessians(
3806 const InputVector &fe_function,
3807 std::vector<solution_hessian_type<typename InputVector::value_type>>
3808 &hessians) const
3809 {
3810 std::vector<typename InputVector::value_type> local_dof_values(
3811 this->fe_interface->n_current_interface_dofs());
3812 this->get_local_dof_values(fe_function, local_dof_values);
3813
3814 get_average_of_function_hessians_from_local_dof_values(local_dof_values,
3815 hessians);
3816 }
3817
3818
3819
3820 template <int dim, int spacedim>
3821 template <class InputVector>
3822 void
3823 Scalar<dim, spacedim>::
3824 get_jump_in_function_third_derivatives_from_local_dof_values(
3825 const InputVector &local_dof_values,
3826 std::vector<
3827 solution_third_derivative_type<typename InputVector::value_type>>
3828 &third_derivatives) const
3829 {
3830 AssertDimension(third_derivatives.size(),
3831 this->fe_interface->n_quadrature_points);
3832
3833 for (const unsigned int dof_index : this->fe_interface->dof_indices())
3834 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3835 {
3836 if (dof_index == 0)
3837 third_derivatives[q_index] = 0.;
3838
3839 third_derivatives[q_index] +=
3840 local_dof_values[dof_index] *
3841 jump_in_third_derivatives(dof_index, q_index);
3842 }
3843 }
3844
3845
3846
3847 template <int dim, int spacedim>
3848 template <class InputVector>
3849 void
3850 Scalar<dim, spacedim>::get_jump_in_function_third_derivatives(
3851 const InputVector &fe_function,
3852 std::vector<
3853 solution_third_derivative_type<typename InputVector::value_type>>
3854 &third_derivatives) const
3855 {
3856 std::vector<typename InputVector::value_type> local_dof_values(
3857 this->fe_interface->n_current_interface_dofs());
3858 this->get_local_dof_values(fe_function, local_dof_values);
3859
3860 get_jump_in_function_third_derivatives_from_local_dof_values(
3861 local_dof_values, third_derivatives);
3862 }
3863
3864
3865
3866 template <int dim, int spacedim>
3868 const FEInterfaceValues<dim, spacedim> &fe_interface,
3869 const unsigned int first_vector_component)
3870 : Base<dim, spacedim>(fe_interface)
3871 , extractor(first_vector_component)
3872 {}
3873
3874
3875
3876 template <int dim, int spacedim>
3878 Vector<dim, spacedim>::value(const bool here_or_there,
3879 const unsigned int interface_dof_index,
3880 const unsigned int q_point) const
3881 {
3882 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3883
3884 if (here_or_there && dof_pair[0] != numbers::invalid_unsigned_int)
3885 return (*(this->fe_interface->fe_face_values))[extractor].value(
3886 dof_pair[0], q_point);
3887
3888 if (!here_or_there && dof_pair[1] != numbers::invalid_unsigned_int)
3889 return (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3890 dof_pair[1], q_point);
3891
3892 return value_type();
3893 }
3894
3895
3896
3897 template <int dim, int spacedim>
3899 Vector<dim, spacedim>::jump_in_values(const unsigned int interface_dof_index,
3900 const unsigned int q_point) const
3901 {
3902 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3903
3904 value_type value;
3905
3906 if (dof_pair[0] != numbers::invalid_unsigned_int)
3907 value +=
3908 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3909 q_point);
3910
3911 if (dof_pair[1] != numbers::invalid_unsigned_int)
3912 value -=
3913 (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3914 dof_pair[1], q_point);
3915
3916 return value;
3917 }
3918
3919
3920
3921 template <int dim, int spacedim>
3924 const unsigned int interface_dof_index,
3925 const unsigned int q_point) const
3926 {
3927 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3928
3929 if (this->fe_interface->at_boundary())
3930 return (*(this->fe_interface->fe_face_values))[extractor].value(
3931 dof_pair[0], q_point);
3932
3933 value_type value;
3934
3935 if (dof_pair[0] != numbers::invalid_unsigned_int)
3936 value +=
3937 0.5 *
3938 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3939 q_point);
3940
3941 if (dof_pair[1] != numbers::invalid_unsigned_int)
3942 value +=
3943 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3944 dof_pair[1], q_point);
3945
3946 return value;
3947 }
3948
3949
3950
3951 template <int dim, int spacedim>
3954 const unsigned int interface_dof_index,
3955 const unsigned int q_point) const
3956 {
3957 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3958
3959 if (this->fe_interface->at_boundary())
3960 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3961 dof_pair[0], q_point);
3962
3963 gradient_type value;
3964
3965 if (dof_pair[0] != numbers::invalid_unsigned_int)
3966 value +=
3967 0.5 *
3968 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3969 q_point);
3970
3971 if (dof_pair[1] != numbers::invalid_unsigned_int)
3972 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3973 .gradient(dof_pair[1], q_point);
3974
3975 return value;
3976 }
3977
3978
3979
3980 template <int dim, int spacedim>
3983 const unsigned int interface_dof_index,
3984 const unsigned int q_point) const
3985 {
3986 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3987
3988 if (this->fe_interface->at_boundary())
3989 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3990 dof_pair[0], q_point);
3991
3992 gradient_type value;
3993
3994 if (dof_pair[0] != numbers::invalid_unsigned_int)
3995 value +=
3996 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3997 q_point);
3998
3999 if (dof_pair[1] != numbers::invalid_unsigned_int)
4000 value -=
4001 (*(this->fe_interface->fe_face_values_neighbor))[extractor].gradient(
4002 dof_pair[1], q_point);
4003
4004 return value;
4005 }
4006
4007
4008
4009 template <int dim, int spacedim>
4011 Vector<dim, spacedim>::jump_gradient(const unsigned int interface_dof_index,
4012 const unsigned int q_point) const
4013 {
4014 return jump_in_gradients(interface_dof_index, q_point);
4015 }
4016
4017
4018
4019 template <int dim, int spacedim>
4022 const unsigned int interface_dof_index,
4023 const unsigned int q_point) const
4024 {
4025 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
4026
4027 if (this->fe_interface->at_boundary())
4028 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
4029 dof_pair[0], q_point);
4030
4031 hessian_type value;
4032
4033 if (dof_pair[0] != numbers::invalid_unsigned_int)
4034 value +=
4035 0.5 *
4036 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
4037 q_point);
4038
4039 if (dof_pair[1] != numbers::invalid_unsigned_int)
4040 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
4041 .hessian(dof_pair[1], q_point);
4042
4043 return value;
4044 }
4045
4046
4047
4048 template <int dim, int spacedim>
4050 Vector<dim, spacedim>::average_hessian(const unsigned int interface_dof_index,
4051 const unsigned int q_point) const
4052 {
4053 return average_of_hessians(interface_dof_index, q_point);
4054 }
4055
4056
4057
4058 template <int dim, int spacedim>
4061 const unsigned int interface_dof_index,
4062 const unsigned int q_point) const
4063 {
4064 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
4065
4066 if (this->fe_interface->at_boundary())
4067 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
4068 dof_pair[0], q_point);
4069
4070 hessian_type value;
4071
4072 if (dof_pair[0] != numbers::invalid_unsigned_int)
4073 value +=
4074 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
4075 q_point);
4076
4077 if (dof_pair[1] != numbers::invalid_unsigned_int)
4078 value -=
4079 (*(this->fe_interface->fe_face_values_neighbor))[extractor].hessian(
4080 dof_pair[1], q_point);
4081
4082 return value;
4083 }
4084
4085
4086
4087 template <int dim, int spacedim>
4089 Vector<dim, spacedim>::jump_hessian(const unsigned int interface_dof_index,
4090 const unsigned int q_point) const
4091 {
4092 return jump_in_hessians(interface_dof_index, q_point);
4093 }
4094
4095
4096
4097 template <int dim, int spacedim>
4100 const unsigned int interface_dof_index,
4101 const unsigned int q_point) const
4102 {
4103 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
4104
4105 if (this->fe_interface->at_boundary())
4106 return (*(this->fe_interface->fe_face_values))[extractor]
4107 .third_derivative(dof_pair[0], q_point);
4108
4109 third_derivative_type value;
4110
4111 if (dof_pair[0] != numbers::invalid_unsigned_int)
4112 value +=
4113 (*(this->fe_interface->fe_face_values))[extractor].third_derivative(
4114 dof_pair[0], q_point);
4115
4116 if (dof_pair[1] != numbers::invalid_unsigned_int)
4117 value -= (*(this->fe_interface->fe_face_values_neighbor))[extractor]
4118 .third_derivative(dof_pair[1], q_point);
4119
4120 return value;
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 &local_dof_values,
4131 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4132 const
4133 {
4134 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
4135
4136 for (const auto dof_index : this->fe_interface->dof_indices())
4137 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4138 {
4139 if (dof_index == 0)
4140 values[q_index] = 0.;
4141
4142 values[q_index] += local_dof_values[dof_index] *
4143 value(here_or_there, dof_index, q_index);
4144 }
4145 }
4146
4147
4148
4149 template <int dim, int spacedim>
4150 template <class InputVector>
4151 void
4153 const bool here_or_there,
4154 const InputVector &fe_function,
4155 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4156 const
4157 {
4158 std::vector<typename InputVector::value_type> local_dof_values(
4159 this->fe_interface->n_current_interface_dofs());
4160 this->get_local_dof_values(fe_function, local_dof_values);
4161
4162 get_function_values_from_local_dof_values(here_or_there,
4163 local_dof_values,
4164 values);
4165 }
4166
4167
4168
4169 template <int dim, int spacedim>
4170 template <class InputVector>
4171 void
4173 const InputVector &local_dof_values,
4174 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4175 const
4176 {
4177 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
4178
4179 for (const auto dof_index : this->fe_interface->dof_indices())
4180 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4181 {
4182 if (dof_index == 0)
4183 values[q_index] = 0.;
4184
4185 values[q_index] +=
4186 local_dof_values[dof_index] * jump_in_values(dof_index, q_index);
4187 }
4188 }
4189
4190
4191
4192 template <int dim, int spacedim>
4193 template <class InputVector>
4194 void
4196 const InputVector &fe_function,
4197 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4198 const
4199 {
4200 std::vector<typename InputVector::value_type> local_dof_values(
4201 this->fe_interface->n_current_interface_dofs());
4202 this->get_local_dof_values(fe_function, local_dof_values);
4203
4204 get_jump_in_function_values_from_local_dof_values(local_dof_values, values);
4205 }
4206
4207
4208
4209 template <int dim, int spacedim>
4210 template <class InputVector>
4211 void
4213 const InputVector &local_dof_values,
4214 std::vector<solution_gradient_type<typename InputVector::value_type>>
4215 &gradients) const
4216 {
4217 AssertDimension(gradients.size(), this->fe_interface->n_quadrature_points);
4218
4219 for (const auto dof_index : this->fe_interface->dof_indices())
4220 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4221 {
4222 if (dof_index == 0)
4223 gradients[q_index] = 0.;
4224
4225 gradients[q_index] +=
4226 local_dof_values[dof_index] * jump_in_gradients(dof_index, q_index);
4227 }
4228 }
4229
4230
4231
4232 template <int dim, int spacedim>
4233 template <class InputVector>
4234 void
4236 const InputVector &fe_function,
4237 std::vector<solution_gradient_type<typename InputVector::value_type>>
4238 &gradients) const
4239 {
4240 std::vector<typename InputVector::value_type> local_dof_values(
4241 this->fe_interface->n_current_interface_dofs());
4242 this->get_local_dof_values(fe_function, local_dof_values);
4243
4244 get_jump_in_function_gradients_from_local_dof_values(local_dof_values,
4245 gradients);
4246 }
4247
4248
4249
4250 template <int dim, int spacedim>
4251 template <class InputVector>
4252 void
4254 const InputVector &local_dof_values,
4255 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4256 const
4257 {
4258 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
4259
4260 for (const auto dof_index : this->fe_interface->dof_indices())
4261 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4262 {
4263 if (dof_index == 0)
4264 values[q_index] = 0.;
4265
4266 values[q_index] +=
4267 local_dof_values[dof_index] * average_of_values(dof_index, q_index);
4268 }
4269 }
4270
4271
4272
4273 template <int dim, int spacedim>
4274 template <class InputVector>
4275 void
4277 const InputVector &fe_function,
4278 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4279 const
4280 {
4281 std::vector<typename InputVector::value_type> local_dof_values(
4282 this->fe_interface->n_current_interface_dofs());
4283 this->get_local_dof_values(fe_function, local_dof_values);
4284
4285 get_average_of_function_values_from_local_dof_values(local_dof_values,
4286 values);
4287 }
4288
4289
4290
4291 template <int dim, int spacedim>
4292 template <class InputVector>
4293 void
4296 const InputVector &local_dof_values,
4297 std::vector<solution_gradient_type<typename InputVector::value_type>>
4298 &gradients) const
4299 {
4300 AssertDimension(gradients.size(), this->fe_interface->n_quadrature_points);
4301
4302 for (const auto dof_index : this->fe_interface->dof_indices())
4303 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4304 {
4305 if (dof_index == 0)
4306 gradients[q_index] = 0.;
4307
4308 gradients[q_index] += local_dof_values[dof_index] *
4309 average_of_gradients(dof_index, q_index);
4310 }
4311 }
4312
4313
4314
4315 template <int dim, int spacedim>
4316 template <class InputVector>
4317 void
4319 const InputVector &fe_function,
4320 std::vector<solution_gradient_type<typename InputVector::value_type>>
4321 &gradients) const
4322 {
4323 std::vector<typename InputVector::value_type> local_dof_values(
4324 this->fe_interface->n_current_interface_dofs());
4325 this->get_local_dof_values(fe_function, local_dof_values);
4326
4327 get_average_of_function_gradients_from_local_dof_values(local_dof_values,
4328 gradients);
4329 }
4330
4331
4332
4333 template <int dim, int spacedim>
4334 template <class InputVector>
4335 void
4337 const InputVector &local_dof_values,
4338 std::vector<solution_hessian_type<typename InputVector::value_type>>
4339 &hessians) const
4340 {
4341 AssertDimension(hessians.size(), this->fe_interface->n_quadrature_points);
4342
4343 for (const auto dof_index : this->fe_interface->dof_indices())
4344 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4345 {
4346 if (dof_index == 0)
4347 hessians[q_index] = 0.;
4348
4349 hessians[q_index] +=
4350 local_dof_values[dof_index] * jump_in_hessians(dof_index, q_index);
4351 }
4352 }
4353
4354
4355
4356 template <int dim, int spacedim>
4357 template <class InputVector>
4358 void
4360 const InputVector &fe_function,
4361 std::vector<solution_hessian_type<typename InputVector::value_type>>
4362 &hessians) const
4363 {
4364 std::vector<typename InputVector::value_type> local_dof_values(
4365 this->fe_interface->n_current_interface_dofs());
4366 this->get_local_dof_values(fe_function, local_dof_values);
4367
4368 get_jump_in_function_hessians_from_local_dof_values(local_dof_values,
4369 hessians);
4370 }
4371
4372
4373
4374 template <int dim, int spacedim>
4375 template <class InputVector>
4376 void
4378 const InputVector &local_dof_values,
4379 std::vector<solution_hessian_type<typename InputVector::value_type>>
4380 &hessians) const
4381 {
4382 AssertDimension(hessians.size(), this->fe_interface->n_quadrature_points);
4383
4384 for (const unsigned int dof_index : this->fe_interface->dof_indices())
4385 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4386 {
4387 if (dof_index == 0)
4388 hessians[q_index] = 0.;
4389
4390 hessians[q_index] += local_dof_values[dof_index] *
4391 average_of_hessians(dof_index, q_index);
4392 }
4393 }
4394
4395
4396
4397 template <int dim, int spacedim>
4398 template <class InputVector>
4399 void
4401 const InputVector &fe_function,
4402 std::vector<solution_hessian_type<typename InputVector::value_type>>
4403 &hessians) const
4404 {
4405 std::vector<typename InputVector::value_type> local_dof_values(
4406 this->fe_interface->n_current_interface_dofs());
4407 this->get_local_dof_values(fe_function, local_dof_values);
4408
4409 get_average_of_function_hessians_from_local_dof_values(local_dof_values,
4410 hessians);
4411 }
4412
4413
4414
4415 template <int dim, int spacedim>
4416 template <class InputVector>
4417 void
4420 const InputVector &local_dof_values,
4421 std::vector<
4422 solution_third_derivative_type<typename InputVector::value_type>>
4423 &third_derivatives) const
4424 {
4425 AssertDimension(third_derivatives.size(),
4426 this->fe_interface->n_quadrature_points);
4427
4428 for (const unsigned int dof_index : this->fe_interface->dof_indices())
4429 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4430 {
4431 if (dof_index == 0)
4432 third_derivatives[q_index] = 0.;
4433
4434 third_derivatives[q_index] +=
4435 local_dof_values[dof_index] *
4436 jump_in_third_derivatives(dof_index, q_index);
4437 }
4438 }
4439
4440
4441
4442 template <int dim, int spacedim>
4443 template <class InputVector>
4444 void
4446 const InputVector &fe_function,
4447 std::vector<
4448 solution_third_derivative_type<typename InputVector::value_type>>
4449 &third_derivatives) const
4450 {
4451 std::vector<typename InputVector::value_type> local_dof_values(
4452 this->fe_interface->n_current_interface_dofs());
4453 this->get_local_dof_values(fe_function, local_dof_values);
4454
4455 get_jump_in_function_third_derivatives_from_local_dof_values(
4456 local_dof_values, third_derivatives);
4457 }
4458} // namespace FEInterfaceViews
4459
4460#endif // DOXYGEN
4461
4463
4464#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
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 FEValuesViews::Scalar< dim, spacedim >::gradient_type 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
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
typename ProductType< Number, third_derivative_type >::type solution_third_derivative_type
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)
typename FEValuesViews::Scalar< dim, spacedim >::hessian_type hessian_type
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
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
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
typename FEValuesViews::Scalar< dim, spacedim >::third_derivative_type third_derivative_type
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 ProductType< Number, hessian_type >::type solution_hessian_type
typename ProductType< Number, gradient_type >::type solution_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
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
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
typename ProductType< Number, gradient_type >::type solution_gradient_type
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
void get_average_of_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
value_type value(const bool here_or_there, const unsigned int interface_dof_index, const unsigned int q_point) const
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
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 ProductType< Number, hessian_type >::type solution_hessian_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
typename ProductType< Number, third_derivative_type >::type solution_third_derivative_type
typename FEValuesViews::Vector< dim, spacedim >::gradient_type gradient_type
gradient_type jump_in_gradients(const unsigned int interface_dof_index, const unsigned int q_point) const
typename FEValuesViews::Vector< dim, spacedim >::third_derivative_type third_derivative_type
hessian_type average_of_hessians(const unsigned int interface_dof_index, const unsigned int q_point) const
typename FEValuesViews::Vector< dim, spacedim >::hessian_type hessian_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:1124
Number value_type
Definition vector.h:137
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
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:489
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
std::size_t size
Definition mpi.cc:734
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
const types::global_dof_index invalid_dof_index
Definition types.h:252
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