15#ifndef dealii_fe_interface_values_h
16#define dealii_fe_interface_values_h
32template <
int dim,
int spacedim>
45 template <
int dim,
int spacedim = dim>
64 template <
class InputVector,
class OutputVector>
67 OutputVector &local_dof_values)
const;
75 template <
int dim,
int spacedim = dim>
111 template <
typename Number>
120 template <
typename Number>
130 template <
typename Number>
140 template <
typename Number>
148 const unsigned int component);
176 const unsigned int interface_dof_index,
177 const unsigned int q_point)
const;
198 const unsigned int q_point)
const;
212 const unsigned int q_point)
const;
226 const unsigned int q_point)
const;
240 const unsigned int q_point)
const;
261 const unsigned int q_point)
const;
275 const unsigned int q_point)
const;
290 const unsigned int q_point)
const;
318 template <
class InputVector>
321 const bool here_or_there,
322 const InputVector &fe_function,
348 template <
class InputVector>
351 const bool here_or_there,
352 const InputVector &local_dof_values,
376 template <
class InputVector>
379 const InputVector &fe_function,
390 template <
class InputVector>
393 const InputVector &local_dof_values,
410 template <
class InputVector>
413 const InputVector &fe_function,
424 template <
class InputVector>
427 const InputVector &local_dof_values,
444 template <
class InputVector>
447 const InputVector &fe_function,
458 template <
class InputVector>
461 const InputVector &local_dof_values,
479 template <
class InputVector>
482 const InputVector &fe_function,
485 &third_derivatives)
const;
494 template <
class InputVector>
497 const InputVector &local_dof_values,
500 &third_derivatives)
const;
522 template <
class InputVector>
525 const InputVector &fe_function,
536 template <
class InputVector>
539 const InputVector &local_dof_values,
556 template <
class InputVector>
559 const InputVector &fe_function,
570 template <
class InputVector>
573 const InputVector &local_dof_values,
590 template <
class InputVector>
593 const InputVector &fe_function,
604 template <
class InputVector>
607 const InputVector &local_dof_values,
625 template <
int dim,
int spacedim = dim>
664 template <
typename Number>
673 template <
typename Number>
683 template <
typename Number>
693 template <
typename Number>
701 const unsigned int first_vector_component);
729 const unsigned int interface_dof_index,
730 const unsigned int q_point)
const;
750 const unsigned int q_point)
const;
763 const unsigned int q_point)
const;
772 const unsigned int q_point)
const;
786 const unsigned int q_point)
const;
795 const unsigned int q_point)
const;
809 const unsigned int q_point)
const;
829 const unsigned int q_point)
const;
842 const unsigned int q_point)
const;
857 const unsigned int q_point)
const;
866 const unsigned int q_point)
const;
894 template <
class InputVector>
897 const bool here_or_there,
898 const InputVector &fe_function,
924 template <
class InputVector>
927 const bool here_or_there,
928 const InputVector &local_dof_values,
952 template <
class InputVector>
955 const InputVector &fe_function,
966 template <
class InputVector>
969 const InputVector &local_dof_values,
986 template <
class InputVector>
989 const InputVector &fe_function,
1000 template <
class InputVector>
1003 const InputVector &local_dof_values,
1020 template <
class InputVector>
1023 const InputVector &fe_function,
1034 template <
class InputVector>
1037 const InputVector &local_dof_values,
1055 template <
class InputVector>
1058 const InputVector &fe_function,
1061 &third_derivatives)
const;
1070 template <
class InputVector>
1073 const InputVector &local_dof_values,
1076 &third_derivatives)
const;
1098 template <
class InputVector>
1101 const InputVector &fe_function,
1112 template <
class InputVector>
1115 const InputVector &local_dof_values,
1132 template <
class InputVector>
1135 const InputVector &fe_function,
1146 template <
class InputVector>
1149 const InputVector &local_dof_values,
1166 template <
class InputVector>
1169 const InputVector &fe_function,
1180 template <
class InputVector>
1183 const InputVector &local_dof_values,
1206 template <
int dim,
int spacedim,
typename Extractor>
1217 template <
int dim,
int spacedim>
1220 using type = typename ::FEInterfaceViews::Scalar<dim, spacedim>;
1230 template <
int dim,
int spacedim>
1233 using type = typename ::FEInterfaceViews::Vector<dim, spacedim>;
1244 template <
int dim,
int spacedim,
typename Extractor>
1245 using View = typename ::internal::FEInterfaceViews::
1246 ViewType<dim, spacedim, Extractor>::type;
1275template <
int dim,
int spacedim = dim>
1447 template <
typename CellIteratorType,
typename CellNeighborIteratorType>
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,
1486 template <
typename CellIteratorType>
1489 const unsigned int face_no,
1598 JxW(
const unsigned int quadrature_point)
const;
1605 const std::vector<double> &
1617 const std::vector<Tensor<1, spacedim>> &
1635 const std::vector<Point<spacedim>> &
1691 std::vector<types::global_dof_index>
1706 std::array<unsigned int, 2>
1718 normal(
const unsigned int q_point_index)
const;
1752 const unsigned int interface_dof_index,
1753 const unsigned int q_point,
1754 const unsigned int component = 0)
const;
1785 const unsigned int q_point,
1786 const unsigned int component = 0)
const;
1803 const unsigned int q_point,
1804 const unsigned int component = 0)
const;
1822 const unsigned int q_point,
1823 const unsigned int component = 0)
const;
1840 const unsigned int q_point,
1841 const unsigned int component = 0)
const;
1867 const unsigned int q_point,
1868 const unsigned int component = 0)
const;
1885 const unsigned int q_point,
1886 const unsigned int component = 0)
const;
1904 const unsigned int q_point,
1905 const unsigned int component = 0)
const;
1926 template <
class InputVector>
1929 const InputVector &fe_function,
1930 std::vector<typename InputVector::value_type> &values)
const;
1940 template <
class InputVector>
1943 const InputVector &fe_function,
1954 template <
class InputVector>
1957 const InputVector &fe_function,
1969 template <
class InputVector>
1972 const InputVector &fe_function,
1974 &third_derivatives)
const;
1991 template <
class InputVector>
1994 const InputVector &fe_function,
1995 std::vector<typename InputVector::value_type> &values)
const;
2004 template <
class InputVector>
2007 const InputVector &fe_function,
2018 template <
class InputVector>
2021 const InputVector &fe_function,
2069 std::vector<std::array<unsigned int, 2>>
dofmap;
2107 std::unique_ptr<FESubfaceValues<dim, spacedim>>
2127 std::unique_ptr<hp::FESubfaceValues<dim, spacedim>>
2134 std::unique_ptr<hp::FEFaceValues<dim, spacedim>>
2141 std::unique_ptr<hp::FESubfaceValues<dim, spacedim>>
2151 "The current function doesn't make sense when used with a "
2152 "FEInterfaceValues object with hp-capabilities.");
2161 "The current function doesn't make sense when used with a "
2162 "FEInterfaceValues object without hp-capabilities.");
2182template <
int dim,
int spacedim>
2188 : n_quadrature_points(quadrature.size())
2189 , fe_face_values(nullptr)
2190 , fe_face_values_neighbor(nullptr)
2191 , internal_fe_face_values(
2196 , internal_fe_subface_values(
2201 , internal_fe_face_values_neighbor(
2206 , internal_fe_subface_values_neighbor(
2215template <
int dim,
int spacedim>
2229template <
int dim,
int spacedim>
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(
2243 , internal_fe_subface_values(
2248 , internal_fe_face_values_neighbor(
2253 , internal_fe_subface_values_neighbor(
2262template <
int dim,
int spacedim>
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(
2274 quadrature_collection,
2276 , internal_hp_fe_subface_values(
2280 quadrature_collection,
2282 , internal_hp_fe_face_values_neighbor(
2285 quadrature_collection,
2287 , internal_hp_fe_subface_values_neighbor(
2291 quadrature_collection,
2299template <
int dim,
int spacedim>
2306 quadrature_collection,
2312template <
int dim,
int spacedim>
2313template <
typename CellIteratorType,
typename CellNeighborIteratorType>
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)
2327 Assert(internal_fe_face_values || internal_hp_fe_face_values,
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>;
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>;
2342 if (internal_fe_face_values)
2346 internal_fe_face_values->reinit(cell, face_no);
2347 fe_face_values = internal_fe_face_values.get();
2351 internal_fe_subface_values->reinit(cell, face_no, sub_face_no);
2352 fe_face_values = internal_fe_subface_values.get();
2356 internal_fe_face_values_neighbor->reinit(cell_neighbor,
2358 fe_face_values_neighbor = internal_fe_face_values_neighbor.get();
2362 internal_fe_subface_values_neighbor->reinit(cell_neighbor,
2364 sub_face_no_neighbor);
2365 fe_face_values_neighbor = internal_fe_subface_values_neighbor.get();
2369 fe_face_values_neighbor->n_quadrature_points);
2371 const_cast<unsigned int &
>(this->n_quadrature_points) =
2372 fe_face_values->n_quadrature_points;
2374 else if (internal_hp_fe_face_values)
2376 unsigned int active_fe_index = fe_index_in;
2377 unsigned int active_fe_index_neighbor =
2379 fe_index_neighbor_in :
2384 if constexpr (is_dof_cell_accessor)
2385 active_fe_index = cell->active_fe_index();
2387 active_fe_index = 0;
2392 if constexpr (is_dof_cell_accessor_neighbor)
2393 active_fe_index_neighbor = cell_neighbor->active_fe_index();
2395 active_fe_index_neighbor = 0;
2398 unsigned int used_q_index = q_index;
2399 unsigned int used_mapping_index = mapping_index;
2404 if (internal_hp_fe_face_values->get_quadrature_collection().size() == 1)
2408 if (internal_hp_fe_face_values->get_mapping_collection().size() == 1)
2409 used_mapping_index = 0;
2416 if (internal_hp_fe_face_values
2417 ->get_quadrature_collection()[active_fe_index] ==
2418 internal_hp_fe_face_values
2419 ->get_quadrature_collection()[active_fe_index_neighbor])
2420 used_q_index = active_fe_index;
2424 const unsigned int dominated_fe_index =
2427 internal_hp_fe_face_values->get_fe_collection().find_dominated_fe(
2428 {active_fe_index, active_fe_index_neighbor}) :
2429 numbers::invalid_unsigned_int);
2435 "You called this function with 'q_index' left at its "
2436 "default value, but this can only work if one of "
2437 "the two finite elements adjacent to this face "
2438 "dominates the other. See the documentation "
2439 "of this function for more information of how "
2440 "to deal with this situation."));
2441 used_q_index = dominated_fe_index;
2448 "You called this function with 'mapping_index' left "
2449 "at its default value, but this can only work if one "
2450 "of the two finite elements adjacent to this face "
2451 "dominates the other. See the documentation "
2452 "of this function for more information of how "
2453 "to deal with this situation."));
2454 used_mapping_index = dominated_fe_index;
2460 internal_hp_fe_face_values->reinit(
2461 cell, face_no, used_q_index, used_mapping_index, active_fe_index);
2463 internal_hp_fe_face_values->get_present_fe_values());
2467 internal_hp_fe_subface_values->
reinit(cell,
2475 internal_hp_fe_subface_values->get_present_fe_values());
2479 internal_hp_fe_face_values_neighbor->
reinit(cell_neighbor,
2483 active_fe_index_neighbor);
2486 internal_hp_fe_face_values_neighbor->get_present_fe_values());
2490 internal_hp_fe_subface_values_neighbor->
reinit(
2493 sub_face_no_neighbor,
2496 active_fe_index_neighbor);
2498 fe_face_values_neighbor =
2500 internal_hp_fe_subface_values_neighbor->get_present_fe_values());
2504 fe_face_values_neighbor->n_quadrature_points);
2506 const_cast<unsigned int &
>(this->n_quadrature_points) =
2507 fe_face_values->n_quadrature_points;
2511 if constexpr (is_dof_cell_accessor_neighbor && is_dof_cell_accessor)
2514 std::vector<types::global_dof_index> v(
2515 fe_face_values->get_fe().n_dofs_per_cell());
2516 cell->get_active_or_mg_dof_indices(v);
2517 std::vector<types::global_dof_index> v2(
2518 fe_face_values_neighbor->get_fe().n_dofs_per_cell());
2519 cell_neighbor->get_active_or_mg_dof_indices(v2);
2523 std::map<types::global_dof_index, std::pair<unsigned int, unsigned int>>
2525 std::pair<unsigned int, unsigned int> invalid_entry(
2528 for (
unsigned int i = 0; i < v.size(); ++i)
2531 auto result = tempmap.insert(std::make_pair(v[i], invalid_entry));
2532 result.first->second.first = i;
2535 for (
unsigned int i = 0; i < v2.size(); ++i)
2538 auto result = tempmap.insert(std::make_pair(v2[i], invalid_entry));
2539 result.first->second.second = i;
2543 dofmap.resize(tempmap.size());
2544 interface_dof_indices.resize(tempmap.size());
2545 unsigned int idx = 0;
2546 for (
auto &x : tempmap)
2548 interface_dof_indices[idx] = x.first;
2549 dofmap[idx] = {{x.second.first, x.second.second}};
2557template <
int dim,
int spacedim>
2558template <
typename CellIteratorType>
2561 const unsigned int face_no,
2562 const unsigned int q_index,
2563 const unsigned int mapping_index,
2564 const unsigned int fe_index)
2566 Assert(internal_fe_face_values || internal_hp_fe_face_values,
2569 if (internal_fe_face_values)
2573 Assert((mapping_index == 0 ||
2579 internal_fe_face_values->reinit(cell, face_no);
2580 fe_face_values = internal_fe_face_values.get();
2581 fe_face_values_neighbor =
nullptr;
2583 else if (internal_hp_fe_face_values)
2585 internal_hp_fe_face_values->reinit(
2586 cell, face_no, q_index, mapping_index, fe_index);
2588 internal_hp_fe_face_values->get_present_fe_values());
2589 fe_face_values_neighbor =
nullptr;
2592 if constexpr (std::is_same_v<
typename CellIteratorType::AccessorType,
2594 std::is_same_v<
typename CellIteratorType::AccessorType,
2597 if (internal_fe_face_values)
2599 interface_dof_indices.resize(
2600 fe_face_values->get_fe().n_dofs_per_cell());
2601 cell->get_active_or_mg_dof_indices(interface_dof_indices);
2603 else if (internal_hp_fe_face_values)
2605 interface_dof_indices.resize(
2606 fe_face_values->get_fe().n_dofs_per_cell());
2607 cell->get_active_or_mg_dof_indices(interface_dof_indices);
2610 dofmap.resize(interface_dof_indices.size());
2611 for (
unsigned int i = 0; i < interface_dof_indices.size(); ++i)
2620template <
int dim,
int spacedim>
2624 Assert(fe_face_values !=
nullptr,
2625 ExcMessage(
"This call requires a call to reinit() first."));
2626 return fe_face_values->JxW(q);
2631template <
int dim,
int spacedim>
2632const std::vector<double> &
2635 Assert(fe_face_values !=
nullptr,
2636 ExcMessage(
"This call requires a call to reinit() first."));
2637 return fe_face_values->get_JxW_values();
2642template <
int dim,
int spacedim>
2643const std::vector<Tensor<1, spacedim>> &
2646 Assert(fe_face_values !=
nullptr,
2647 ExcMessage(
"This call requires a call to reinit() first."));
2648 return fe_face_values->get_normal_vectors();
2653template <
int dim,
int spacedim>
2657 Assert(!has_hp_capabilities(), ExcOnlyAvailableWithoutHP());
2658 return internal_fe_face_values->get_mapping();
2663template <
int dim,
int spacedim>
2667 Assert(!has_hp_capabilities(), ExcOnlyAvailableWithoutHP());
2668 return internal_fe_face_values->get_fe();
2673template <
int dim,
int spacedim>
2677 Assert(!has_hp_capabilities(), ExcOnlyAvailableWithoutHP());
2678 return internal_fe_face_values->get_quadrature();
2683template <
int dim,
int spacedim>
2687 Assert(has_hp_capabilities(), ExcOnlyAvailableWithHP());
2688 return internal_hp_fe_face_values->get_mapping_collection();
2693template <
int dim,
int spacedim>
2697 Assert(has_hp_capabilities(), ExcOnlyAvailableWithHP());
2698 return internal_hp_fe_face_values->get_fe_collection();
2703template <
int dim,
int spacedim>
2707 Assert(has_hp_capabilities(), ExcOnlyAvailableWithHP());
2708 return internal_hp_fe_face_values->get_quadrature_collection();
2713template <
int dim,
int spacedim>
2717 if (internal_hp_fe_face_values || internal_hp_fe_subface_values ||
2718 internal_hp_fe_face_values_neighbor ||
2719 internal_hp_fe_subface_values_neighbor)
2729 Assert(internal_fe_face_values || internal_fe_subface_values ||
2730 internal_fe_face_values_neighbor ||
2731 internal_fe_subface_values_neighbor,
2743template <
int dim,
int spacedim>
2747 return {0U, n_quadrature_points};
2752template <
int dim,
int spacedim>
2753const std::vector<Point<spacedim>> &
2756 Assert(fe_face_values !=
nullptr,
2757 ExcMessage(
"This call requires a call to reinit() first."));
2758 return fe_face_values->get_quadrature_points();
2763template <
int dim,
int spacedim>
2767 if (has_hp_capabilities())
2768 return internal_hp_fe_face_values->get_update_flags();
2770 return internal_fe_face_values->get_update_flags();
2775template <
int dim,
int spacedim>
2779 return get_fe_face_values(
cell_index).get_cell();
2784template <
int dim,
int spacedim>
2789 return get_fe_face_values(
cell_index).get_face_number();
2794template <
int dim,
int spacedim>
2799 interface_dof_indices.size() > 0,
2801 "n_current_interface_dofs() is only available after a call to reinit()."));
2802 return interface_dof_indices.size();
2807template <
int dim,
int spacedim>
2811 return {0U, n_current_interface_dofs()};
2816template <
int dim,
int spacedim>
2820 return fe_face_values_neighbor ==
nullptr;
2825template <
int dim,
int spacedim>
2826std::vector<types::global_dof_index>
2829 return interface_dof_indices;
2834template <
int dim,
int spacedim>
2835std::array<unsigned int, 2>
2837 const unsigned int interface_dof_index)
const
2840 return dofmap[interface_dof_index];
2845template <
int dim,
int spacedim>
2854 "You are on a boundary, so you can only ask for the first FEFaceValues object."));
2856 return (
cell_index == 0) ? *fe_face_values : *fe_face_values_neighbor;
2861template <
int dim,
int spacedim>
2865 return fe_face_values->normal_vector(q_point_index);
2870template <
int dim,
int spacedim>
2873 const bool here_or_there,
2874 const unsigned int interface_dof_index,
2875 const unsigned int q_point,
2876 const unsigned int component)
const
2878 const auto dof_pair = dofmap[interface_dof_index];
2881 return get_fe_face_values(0).shape_value_component(dof_pair[0],
2885 return get_fe_face_values(1).shape_value_component(dof_pair[1],
2894template <
int dim,
int spacedim>
2897 const unsigned int interface_dof_index,
2898 const unsigned int q_point,
2899 const unsigned int component)
const
2901 const auto dof_pair = dofmap[interface_dof_index];
2906 value += get_fe_face_values(0).shape_value_component(dof_pair[0],
2910 value -= get_fe_face_values(1).shape_value_component(dof_pair[1],
2918template <
int dim,
int spacedim>
2921 const unsigned int interface_dof_index,
2922 const unsigned int q_point,
2923 const unsigned int component)
const
2925 const auto dof_pair = dofmap[interface_dof_index];
2928 return get_fe_face_values(0).shape_value_component(dof_pair[0],
2935 value += 0.5 * get_fe_face_values(0).shape_value_component(dof_pair[0],
2939 value += 0.5 * get_fe_face_values(1).shape_value_component(dof_pair[1],
2948template <
int dim,
int spacedim>
2951 const unsigned int interface_dof_index,
2952 const unsigned int q_point,
2953 const unsigned int component)
const
2955 const auto dof_pair = dofmap[interface_dof_index];
2958 return get_fe_face_values(0).shape_grad_component(dof_pair[0],
2965 value += 0.5 * get_fe_face_values(0).shape_grad_component(dof_pair[0],
2969 value += 0.5 * get_fe_face_values(1).shape_grad_component(dof_pair[1],
2978template <
int dim,
int spacedim>
2981 const unsigned int interface_dof_index,
2982 const unsigned int q_point,
2983 const unsigned int component)
const
2985 const auto dof_pair = dofmap[interface_dof_index];
2988 return get_fe_face_values(0).shape_hessian_component(dof_pair[0],
2995 value += 0.5 * get_fe_face_values(0).shape_hessian_component(dof_pair[0],
2999 value += 0.5 * get_fe_face_values(1).shape_hessian_component(dof_pair[1],
3008template <
int dim,
int spacedim>
3011 const unsigned int interface_dof_index,
3012 const unsigned int q_point,
3013 const unsigned int component)
const
3015 const auto dof_pair = dofmap[interface_dof_index];
3018 return get_fe_face_values(0).shape_grad_component(dof_pair[0],
3025 value += get_fe_face_values(0).shape_grad_component(dof_pair[0],
3029 value -= get_fe_face_values(1).shape_grad_component(dof_pair[1],
3038template <
int dim,
int spacedim>
3041 const unsigned int interface_dof_index,
3042 const unsigned int q_point,
3043 const unsigned int component)
const
3045 const auto dof_pair = dofmap[interface_dof_index];
3048 return get_fe_face_values(0).shape_hessian_component(dof_pair[0],
3055 value += get_fe_face_values(0).shape_hessian_component(dof_pair[0],
3059 value -= get_fe_face_values(1).shape_hessian_component(dof_pair[1],
3068template <
int dim,
int spacedim>
3071 const unsigned int interface_dof_index,
3072 const unsigned int q_point,
3073 const unsigned int component)
const
3075 const auto dof_pair = dofmap[interface_dof_index];
3078 return get_fe_face_values(0).shape_3rd_derivative_component(dof_pair[0],
3085 value += get_fe_face_values(0).shape_3rd_derivative_component(dof_pair[0],
3089 value -= get_fe_face_values(1).shape_3rd_derivative_component(dof_pair[1],
3098template <
int dim,
int spacedim>
3099template <
class InputVector>
3102 const InputVector &fe_function,
3103 std::vector<typename InputVector::value_type> &values)
const
3108 this->operator[](scalar).get_jump_in_function_values(fe_function, values);
3113template <
int dim,
int spacedim>
3114template <
class InputVector>
3117 const InputVector &fe_function,
3124 this->operator[](scalar).get_jump_in_function_gradients(fe_function,
3130template <
int dim,
int spacedim>
3131template <
class InputVector>
3134 const InputVector &fe_function,
3141 this->operator[](scalar).get_jump_in_function_hessians(fe_function, hessians);
3146template <
int dim,
int spacedim>
3147template <
class InputVector>
3150 const InputVector &fe_function,
3152 &third_derivatives)
const
3157 this->operator[](scalar).get_jump_in_function_third_derivatives(
3158 fe_function, third_derivatives);
3163template <
int dim,
int spacedim>
3164template <
class InputVector>
3167 const InputVector &fe_function,
3168 std::vector<typename InputVector::value_type> &values)
const
3173 this->operator[](scalar).get_average_of_function_values(fe_function, values);
3178template <
int dim,
int spacedim>
3179template <
class InputVector>
3182 const InputVector &fe_function,
3189 this->operator[](scalar).get_average_of_function_gradients(fe_function,
3195template <
int dim,
int spacedim>
3196template <
class InputVector>
3199 const InputVector &fe_function,
3206 this->operator[](scalar).get_average_of_function_hessians(fe_function,
3213template <
int dim,
int spacedim>
3218 const unsigned int n_components =
3219 (this->has_hp_capabilities() ? this->get_fe_collection().n_components() :
3220 this->get_fe().n_components());
3228template <
int dim,
int spacedim>
3233 const unsigned int n_components =
3234 (this->has_hp_capabilities() ? this->get_fe_collection().n_components() :
3235 this->get_fe().n_components());
3236 const unsigned int n_vectors =
3251 template <
int dim,
int spacedim>
3254 : fe_interface(&fe_interface)
3259 template <
int dim,
int spacedim>
3260 Scalar<dim, spacedim>::Scalar(
3262 const unsigned int component)
3263 : Base<dim, spacedim>(fe_interface)
3264 , extractor(component)
3269 template <
int dim,
int spacedim>
3270 template <
class InputVector,
class OutputVector>
3272 Base<dim, spacedim>::get_local_dof_values(
3273 const InputVector &dof_values,
3274 OutputVector &local_dof_values)
const
3276 const auto &interface_dof_indices =
3277 this->fe_interface->get_interface_dof_indices();
3279 AssertDimension(interface_dof_indices.size(), local_dof_values.size());
3281 for (
const unsigned int i : this->fe_interface->dof_indices())
3282 local_dof_values[i] = dof_values(interface_dof_indices[i]);
3287 template <
int dim,
int spacedim>
3288 typename Scalar<dim, spacedim>::value_type
3289 Scalar<dim, spacedim>::value(
const bool here_or_there,
3290 const unsigned int interface_dof_index,
3291 const unsigned int q_point)
const
3293 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3296 return (*(this->fe_interface->fe_face_values))[extractor].value(
3297 dof_pair[0], q_point);
3300 return (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3301 dof_pair[1], q_point);
3308 template <
int dim,
int spacedim>
3309 typename Scalar<dim, spacedim>::value_type
3310 Scalar<dim, spacedim>::jump_in_values(
const unsigned int interface_dof_index,
3311 const unsigned int q_point)
const
3313 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3315 value_type
value = 0.0;
3319 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3324 (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3325 dof_pair[1], q_point);
3332 template <
int dim,
int spacedim>
3333 typename Scalar<dim, spacedim>::value_type
3334 Scalar<dim, spacedim>::average_of_values(
3335 const unsigned int interface_dof_index,
3336 const unsigned int q_point)
const
3338 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3340 if (this->fe_interface->at_boundary())
3341 return (*(this->fe_interface->fe_face_values))[extractor].value(
3342 dof_pair[0], q_point);
3344 value_type
value = 0.0;
3349 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3354 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3355 dof_pair[1], q_point);
3362 template <
int dim,
int spacedim>
3363 typename Scalar<dim, spacedim>::gradient_type
3364 Scalar<dim, spacedim>::average_of_gradients(
3365 const unsigned int interface_dof_index,
3366 const unsigned int q_point)
const
3368 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3370 if (this->fe_interface->at_boundary())
3371 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3372 dof_pair[0], q_point);
3374 gradient_type
value;
3379 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3383 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3384 .gradient(dof_pair[1], q_point);
3391 template <
int dim,
int spacedim>
3392 typename Scalar<dim, spacedim>::gradient_type
3393 Scalar<dim, spacedim>::jump_in_gradients(
3394 const unsigned int interface_dof_index,
3395 const unsigned int q_point)
const
3397 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3399 if (this->fe_interface->at_boundary())
3400 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3401 dof_pair[0], q_point);
3403 gradient_type
value;
3407 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3412 (*(this->fe_interface->fe_face_values_neighbor))[extractor].gradient(
3413 dof_pair[1], q_point);
3420 template <
int dim,
int spacedim>
3421 typename Scalar<dim, spacedim>::hessian_type
3422 Scalar<dim, spacedim>::average_of_hessians(
3423 const unsigned int interface_dof_index,
3424 const unsigned int q_point)
const
3426 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3428 if (this->fe_interface->at_boundary())
3429 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
3430 dof_pair[0], q_point);
3437 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
3441 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3442 .hessian(dof_pair[1], q_point);
3449 template <
int dim,
int spacedim>
3450 typename Scalar<dim, spacedim>::third_derivative_type
3451 Scalar<dim, spacedim>::jump_in_third_derivatives(
3452 const unsigned int interface_dof_index,
3453 const unsigned int q_point)
const
3455 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3457 if (this->fe_interface->at_boundary())
3458 return (*(this->fe_interface->fe_face_values))[extractor]
3459 .third_derivative(dof_pair[0], q_point);
3461 third_derivative_type
value;
3465 (*(this->fe_interface->fe_face_values))[extractor].third_derivative(
3466 dof_pair[0], q_point);
3469 value -= (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3470 .third_derivative(dof_pair[1], q_point);
3477 template <
int dim,
int spacedim>
3478 typename Scalar<dim, spacedim>::hessian_type
3479 Scalar<dim, spacedim>::jump_in_hessians(
3480 const unsigned int interface_dof_index,
3481 const unsigned int q_point)
const
3483 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3485 if (this->fe_interface->at_boundary())
3486 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
3487 dof_pair[0], q_point);
3493 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
3498 (*(this->fe_interface->fe_face_values_neighbor))[extractor].hessian(
3499 dof_pair[1], q_point);
3506 template <
int dim,
int spacedim>
3507 template <
class InputVector>
3509 Scalar<dim, spacedim>::get_function_values_from_local_dof_values(
3510 const bool here_or_there,
3511 const InputVector &local_dof_values,
3512 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3515 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
3517 for (
const auto dof_index : this->fe_interface->dof_indices())
3518 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3521 values[q_index] = 0.;
3523 values[q_index] += local_dof_values[dof_index] *
3524 value(here_or_there, dof_index, q_index);
3530 template <
int dim,
int spacedim>
3531 template <
class InputVector>
3533 Scalar<dim, spacedim>::get_function_values(
3534 const bool here_or_there,
3535 const InputVector &fe_function,
3536 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3539 std::vector<typename InputVector::value_type> local_dof_values(
3540 this->fe_interface->n_current_interface_dofs());
3541 this->get_local_dof_values(fe_function, local_dof_values);
3543 get_function_values_from_local_dof_values(here_or_there,
3550 template <
int dim,
int spacedim>
3551 template <
class InputVector>
3553 Scalar<dim, spacedim>::get_jump_in_function_values_from_local_dof_values(
3554 const InputVector &local_dof_values,
3555 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3558 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
3560 for (
const auto dof_index : this->fe_interface->dof_indices())
3561 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3564 values[q_index] = 0.;
3567 local_dof_values[dof_index] * jump_in_values(dof_index, q_index);
3573 template <
int dim,
int spacedim>
3574 template <
class InputVector>
3576 Scalar<dim, spacedim>::get_jump_in_function_values(
3577 const InputVector &fe_function,
3578 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3581 std::vector<typename InputVector::value_type> local_dof_values(
3582 this->fe_interface->n_current_interface_dofs());
3583 this->get_local_dof_values(fe_function, local_dof_values);
3585 get_jump_in_function_values_from_local_dof_values(local_dof_values, values);
3590 template <
int dim,
int spacedim>
3591 template <
class InputVector>
3593 Scalar<dim, spacedim>::get_jump_in_function_gradients_from_local_dof_values(
3594 const InputVector &local_dof_values,
3595 std::vector<solution_gradient_type<typename InputVector::value_type>>
3600 for (
const auto dof_index : this->fe_interface->dof_indices())
3601 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3607 local_dof_values[dof_index] * jump_in_gradients(dof_index, q_index);
3613 template <
int dim,
int spacedim>
3614 template <
class InputVector>
3616 Scalar<dim, spacedim>::get_jump_in_function_gradients(
3617 const InputVector &fe_function,
3618 std::vector<solution_gradient_type<typename InputVector::value_type>>
3621 std::vector<typename InputVector::value_type> local_dof_values(
3622 this->fe_interface->n_current_interface_dofs());
3623 this->get_local_dof_values(fe_function, local_dof_values);
3625 get_jump_in_function_gradients_from_local_dof_values(local_dof_values,
3631 template <
int dim,
int spacedim>
3632 template <
class InputVector>
3634 Scalar<dim, spacedim>::get_average_of_function_values_from_local_dof_values(
3635 const InputVector &local_dof_values,
3636 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3639 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
3641 for (
const auto dof_index : this->fe_interface->dof_indices())
3642 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3645 values[q_index] = 0.;
3648 local_dof_values[dof_index] * average_of_values(dof_index, q_index);
3654 template <
int dim,
int spacedim>
3655 template <
class InputVector>
3657 Scalar<dim, spacedim>::get_average_of_function_values(
3658 const InputVector &fe_function,
3659 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3662 std::vector<typename InputVector::value_type> local_dof_values(
3663 this->fe_interface->n_current_interface_dofs());
3664 this->get_local_dof_values(fe_function, local_dof_values);
3666 get_average_of_function_values_from_local_dof_values(local_dof_values,
3672 template <
int dim,
int spacedim>
3673 template <
class InputVector>
3675 Scalar<dim, spacedim>::
3676 get_average_of_function_gradients_from_local_dof_values(
3677 const InputVector &local_dof_values,
3678 std::vector<solution_gradient_type<typename InputVector::value_type>>
3683 for (
const auto dof_index : this->fe_interface->dof_indices())
3684 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3689 gradients[q_index] += local_dof_values[dof_index] *
3690 average_of_gradients(dof_index, q_index);
3696 template <
int dim,
int spacedim>
3697 template <
class InputVector>
3699 Scalar<dim, spacedim>::get_average_of_function_gradients(
3700 const InputVector &fe_function,
3701 std::vector<solution_gradient_type<typename InputVector::value_type>>
3704 std::vector<typename InputVector::value_type> local_dof_values(
3705 this->fe_interface->n_current_interface_dofs());
3706 this->get_local_dof_values(fe_function, local_dof_values);
3708 get_average_of_function_gradients_from_local_dof_values(local_dof_values,
3714 template <
int dim,
int spacedim>
3715 template <
class InputVector>
3717 Scalar<dim, spacedim>::get_jump_in_function_hessians_from_local_dof_values(
3718 const InputVector &local_dof_values,
3719 std::vector<solution_hessian_type<typename InputVector::value_type>>
3724 for (
const auto dof_index : this->fe_interface->dof_indices())
3725 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3731 local_dof_values[dof_index] * jump_in_hessians(dof_index, q_index);
3737 template <
int dim,
int spacedim>
3738 template <
class InputVector>
3740 Scalar<dim, spacedim>::get_jump_in_function_hessians(
3741 const InputVector &fe_function,
3742 std::vector<solution_hessian_type<typename InputVector::value_type>>
3745 std::vector<typename InputVector::value_type> local_dof_values(
3746 this->fe_interface->n_current_interface_dofs());
3747 this->get_local_dof_values(fe_function, local_dof_values);
3749 get_jump_in_function_hessians_from_local_dof_values(local_dof_values,
3755 template <
int dim,
int spacedim>
3756 template <
class InputVector>
3758 Scalar<dim, spacedim>::get_average_of_function_hessians_from_local_dof_values(
3759 const InputVector &local_dof_values,
3760 std::vector<solution_hessian_type<typename InputVector::value_type>>
3765 for (
const unsigned int dof_index : this->fe_interface->dof_indices())
3766 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3771 hessians[q_index] += local_dof_values[dof_index] *
3772 average_of_hessians(dof_index, q_index);
3778 template <
int dim,
int spacedim>
3779 template <
class InputVector>
3781 Scalar<dim, spacedim>::get_average_of_function_hessians(
3782 const InputVector &fe_function,
3783 std::vector<solution_hessian_type<typename InputVector::value_type>>
3786 std::vector<typename InputVector::value_type> local_dof_values(
3787 this->fe_interface->n_current_interface_dofs());
3788 this->get_local_dof_values(fe_function, local_dof_values);
3790 get_average_of_function_hessians_from_local_dof_values(local_dof_values,
3796 template <
int dim,
int spacedim>
3797 template <
class InputVector>
3799 Scalar<dim, spacedim>::
3800 get_jump_in_function_third_derivatives_from_local_dof_values(
3801 const InputVector &local_dof_values,
3803 solution_third_derivative_type<typename InputVector::value_type>>
3804 &third_derivatives)
const
3807 this->fe_interface->n_quadrature_points);
3809 for (
const unsigned int dof_index : this->fe_interface->dof_indices())
3810 for (const auto q_index : this->fe_interface->quadrature_point_indices())
3813 third_derivatives[q_index] = 0.;
3815 third_derivatives[q_index] +=
3816 local_dof_values[dof_index] *
3817 jump_in_third_derivatives(dof_index, q_index);
3823 template <
int dim,
int spacedim>
3824 template <
class InputVector>
3826 Scalar<dim, spacedim>::get_jump_in_function_third_derivatives(
3827 const InputVector &fe_function,
3829 solution_third_derivative_type<typename InputVector::value_type>>
3830 &third_derivatives)
const
3832 std::vector<typename InputVector::value_type> local_dof_values(
3833 this->fe_interface->n_current_interface_dofs());
3834 this->get_local_dof_values(fe_function, local_dof_values);
3836 get_jump_in_function_third_derivatives_from_local_dof_values(
3837 local_dof_values, third_derivatives);
3842 template <
int dim,
int spacedim>
3845 const unsigned int first_vector_component)
3846 : Base<dim, spacedim>(fe_interface)
3847 , extractor(first_vector_component)
3852 template <
int dim,
int spacedim>
3855 const unsigned int interface_dof_index,
3856 const unsigned int q_point)
const
3858 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3861 return (*(this->fe_interface->fe_face_values))[extractor].value(
3862 dof_pair[0], q_point);
3865 return (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3866 dof_pair[1], q_point);
3868 return value_type();
3873 template <
int dim,
int spacedim>
3876 const unsigned int q_point)
const
3878 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3884 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3889 (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3890 dof_pair[1], q_point);
3897 template <
int dim,
int spacedim>
3900 const unsigned int interface_dof_index,
3901 const unsigned int q_point)
const
3903 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3905 if (this->fe_interface->at_boundary())
3906 return (*(this->fe_interface->fe_face_values))[extractor].value(
3907 dof_pair[0], q_point);
3914 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3919 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3920 dof_pair[1], q_point);
3927 template <
int dim,
int spacedim>
3930 const unsigned int interface_dof_index,
3931 const unsigned int q_point)
const
3933 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3935 if (this->fe_interface->at_boundary())
3936 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3937 dof_pair[0], q_point);
3939 gradient_type
value;
3944 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3948 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3949 .gradient(dof_pair[1], q_point);
3956 template <
int dim,
int spacedim>
3959 const unsigned int interface_dof_index,
3960 const unsigned int q_point)
const
3962 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3964 if (this->fe_interface->at_boundary())
3965 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3966 dof_pair[0], q_point);
3968 gradient_type
value;
3972 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3977 (*(this->fe_interface->fe_face_values_neighbor))[extractor].gradient(
3978 dof_pair[1], q_point);
3985 template <
int dim,
int spacedim>
3988 const unsigned int q_point)
const
3990 return jump_in_gradients(interface_dof_index, q_point);
3995 template <
int dim,
int spacedim>
3998 const unsigned int interface_dof_index,
3999 const unsigned int q_point)
const
4001 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
4003 if (this->fe_interface->at_boundary())
4004 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
4005 dof_pair[0], q_point);
4012 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
4016 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
4017 .hessian(dof_pair[1], q_point);
4024 template <
int dim,
int spacedim>
4027 const unsigned int q_point)
const
4029 return average_of_hessians(interface_dof_index, q_point);
4034 template <
int dim,
int spacedim>
4037 const unsigned int interface_dof_index,
4038 const unsigned int q_point)
const
4040 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
4042 if (this->fe_interface->at_boundary())
4043 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
4044 dof_pair[0], q_point);
4050 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
4055 (*(this->fe_interface->fe_face_values_neighbor))[extractor].hessian(
4056 dof_pair[1], q_point);
4063 template <
int dim,
int spacedim>
4066 const unsigned int q_point)
const
4068 return jump_in_hessians(interface_dof_index, q_point);
4073 template <
int dim,
int spacedim>
4076 const unsigned int interface_dof_index,
4077 const unsigned int q_point)
const
4079 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
4081 if (this->fe_interface->at_boundary())
4082 return (*(this->fe_interface->fe_face_values))[extractor]
4083 .third_derivative(dof_pair[0], q_point);
4085 third_derivative_type
value;
4089 (*(this->fe_interface->fe_face_values))[extractor].third_derivative(
4090 dof_pair[0], q_point);
4093 value -= (*(this->fe_interface->fe_face_values_neighbor))[extractor]
4094 .third_derivative(dof_pair[1], q_point);
4101 template <
int dim,
int spacedim>
4102 template <
class InputVector>
4105 const bool here_or_there,
4106 const InputVector &local_dof_values,
4107 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4110 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
4112 for (
const auto dof_index : this->fe_interface->dof_indices())
4113 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4116 values[q_index] = 0.;
4118 values[q_index] += local_dof_values[dof_index] *
4119 value(here_or_there, dof_index, q_index);
4125 template <
int dim,
int spacedim>
4126 template <
class InputVector>
4129 const bool here_or_there,
4130 const InputVector &fe_function,
4131 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4134 std::vector<typename InputVector::value_type> local_dof_values(
4135 this->fe_interface->n_current_interface_dofs());
4136 this->get_local_dof_values(fe_function, local_dof_values);
4138 get_function_values_from_local_dof_values(here_or_there,
4145 template <
int dim,
int spacedim>
4146 template <
class InputVector>
4149 const InputVector &local_dof_values,
4150 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4153 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
4155 for (
const auto dof_index : this->fe_interface->dof_indices())
4156 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4159 values[q_index] = 0.;
4162 local_dof_values[dof_index] * jump_in_values(dof_index, q_index);
4168 template <
int dim,
int spacedim>
4169 template <
class InputVector>
4172 const InputVector &fe_function,
4173 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4176 std::vector<typename InputVector::value_type> local_dof_values(
4177 this->fe_interface->n_current_interface_dofs());
4178 this->get_local_dof_values(fe_function, local_dof_values);
4180 get_jump_in_function_values_from_local_dof_values(local_dof_values, values);
4185 template <
int dim,
int spacedim>
4186 template <
class InputVector>
4189 const InputVector &local_dof_values,
4190 std::vector<solution_gradient_type<typename InputVector::value_type>>
4195 for (
const auto dof_index : this->fe_interface->dof_indices())
4196 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4202 local_dof_values[dof_index] * jump_in_gradients(dof_index, q_index);
4208 template <
int dim,
int spacedim>
4209 template <
class InputVector>
4212 const InputVector &fe_function,
4213 std::vector<solution_gradient_type<typename InputVector::value_type>>
4216 std::vector<typename InputVector::value_type> local_dof_values(
4217 this->fe_interface->n_current_interface_dofs());
4218 this->get_local_dof_values(fe_function, local_dof_values);
4220 get_jump_in_function_gradients_from_local_dof_values(local_dof_values,
4226 template <
int dim,
int spacedim>
4227 template <
class InputVector>
4230 const InputVector &local_dof_values,
4231 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4234 AssertDimension(values.size(), this->fe_interface->n_quadrature_points);
4236 for (
const auto dof_index : this->fe_interface->dof_indices())
4237 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4240 values[q_index] = 0.;
4243 local_dof_values[dof_index] * average_of_values(dof_index, q_index);
4249 template <
int dim,
int spacedim>
4250 template <
class InputVector>
4253 const InputVector &fe_function,
4254 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4257 std::vector<typename InputVector::value_type> local_dof_values(
4258 this->fe_interface->n_current_interface_dofs());
4259 this->get_local_dof_values(fe_function, local_dof_values);
4261 get_average_of_function_values_from_local_dof_values(local_dof_values,
4267 template <
int dim,
int spacedim>
4268 template <
class InputVector>
4272 const InputVector &local_dof_values,
4273 std::vector<solution_gradient_type<typename InputVector::value_type>>
4278 for (
const auto dof_index : this->fe_interface->dof_indices())
4279 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4284 gradients[q_index] += local_dof_values[dof_index] *
4285 average_of_gradients(dof_index, q_index);
4291 template <
int dim,
int spacedim>
4292 template <
class InputVector>
4295 const InputVector &fe_function,
4296 std::vector<solution_gradient_type<typename InputVector::value_type>>
4299 std::vector<typename InputVector::value_type> local_dof_values(
4300 this->fe_interface->n_current_interface_dofs());
4301 this->get_local_dof_values(fe_function, local_dof_values);
4303 get_average_of_function_gradients_from_local_dof_values(local_dof_values,
4309 template <
int dim,
int spacedim>
4310 template <
class InputVector>
4313 const InputVector &local_dof_values,
4314 std::vector<solution_hessian_type<typename InputVector::value_type>>
4319 for (
const auto dof_index : this->fe_interface->dof_indices())
4320 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4326 local_dof_values[dof_index] * jump_in_hessians(dof_index, q_index);
4332 template <
int dim,
int spacedim>
4333 template <
class InputVector>
4336 const InputVector &fe_function,
4337 std::vector<solution_hessian_type<typename InputVector::value_type>>
4340 std::vector<typename InputVector::value_type> local_dof_values(
4341 this->fe_interface->n_current_interface_dofs());
4342 this->get_local_dof_values(fe_function, local_dof_values);
4344 get_jump_in_function_hessians_from_local_dof_values(local_dof_values,
4350 template <
int dim,
int spacedim>
4351 template <
class InputVector>
4354 const InputVector &local_dof_values,
4355 std::vector<solution_hessian_type<typename InputVector::value_type>>
4360 for (
const unsigned int dof_index : this->fe_interface->dof_indices())
4361 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4366 hessians[q_index] += local_dof_values[dof_index] *
4367 average_of_hessians(dof_index, q_index);
4373 template <
int dim,
int spacedim>
4374 template <
class InputVector>
4377 const InputVector &fe_function,
4378 std::vector<solution_hessian_type<typename InputVector::value_type>>
4381 std::vector<typename InputVector::value_type> local_dof_values(
4382 this->fe_interface->n_current_interface_dofs());
4383 this->get_local_dof_values(fe_function, local_dof_values);
4385 get_average_of_function_hessians_from_local_dof_values(local_dof_values,
4391 template <
int dim,
int spacedim>
4392 template <
class InputVector>
4396 const InputVector &local_dof_values,
4398 solution_third_derivative_type<typename InputVector::value_type>>
4399 &third_derivatives)
const
4402 this->fe_interface->n_quadrature_points);
4404 for (
const unsigned int dof_index : this->fe_interface->dof_indices())
4405 for (const auto q_index : this->fe_interface->quadrature_point_indices())
4408 third_derivatives[q_index] = 0.;
4410 third_derivatives[q_index] +=
4411 local_dof_values[dof_index] *
4412 jump_in_third_derivatives(dof_index, q_index);
4418 template <
int dim,
int spacedim>
4419 template <
class InputVector>
4422 const InputVector &fe_function,
4424 solution_third_derivative_type<typename InputVector::value_type>>
4425 &third_derivatives)
const
4427 std::vector<typename InputVector::value_type> local_dof_values(
4428 this->fe_interface->n_current_interface_dofs());
4429 this->get_local_dof_values(fe_function, local_dof_values);
4431 get_jump_in_function_third_derivatives_from_local_dof_values(
4432 local_dof_values, third_derivatives);
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
const std::vector< double > & get_JxW_values() const
void get_average_of_function_hessians(const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type > > &hessians) const
FEInterfaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature, const UpdateFlags update_flags)
void reinit(const CellIteratorType &cell, const unsigned int face_no, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
const hp::FECollection< dim, spacedim > & get_fe_collection() const
Tensor< 3, spacedim > jump_in_shape_3rd_derivatives(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
FEInterfaceValues(const hp::MappingCollection< dim, spacedim > &mapping_collection, const hp::FECollection< dim, spacedim > &fe_collection, const hp::QCollection< dim - 1 > &quadrature_collection, const UpdateFlags update_flags)
std::unique_ptr< hp::FESubfaceValues< dim, spacedim > > internal_hp_fe_subface_values
const unsigned int n_quadrature_points
bool has_hp_capabilities() const
void reinit(const CellIteratorType &cell, const unsigned int face_no, const unsigned int sub_face_no, const CellNeighborIteratorType &cell_neighbor, const unsigned int face_no_neighbor, const unsigned int sub_face_no_neighbor, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int, const unsigned int fe_index_neighbor=numbers::invalid_unsigned_int)
FEInterfaceViews::Scalar< dim, spacedim > operator[](const FEValuesExtractors::Scalar &scalar) const
unsigned n_current_interface_dofs() const
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors() const
const FiniteElement< dim, spacedim > & get_fe() const
double average_of_shape_values(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
const std::vector< Point< spacedim > > & get_quadrature_points() const
std::unique_ptr< FESubfaceValues< dim, spacedim > > internal_fe_subface_values
double shape_value(const bool here_or_there, const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
const FEFaceValuesBase< dim, spacedim > & get_fe_face_values(const unsigned int cell_index) const
Tensor< 2, spacedim > average_of_shape_hessians(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
std::unique_ptr< hp::FESubfaceValues< dim, spacedim > > internal_hp_fe_subface_values_neighbor
FEInterfaceValues(const hp::FECollection< dim, spacedim > &fe_collection, const hp::QCollection< dim - 1 > &quadrature_collection, const UpdateFlags update_flags)
const hp::QCollection< dim - 1 > & get_quadrature_collection() const
FEInterfaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const hp::QCollection< dim - 1 > &quadrature, const UpdateFlags update_flags)
FEFaceValuesBase< dim, spacedim > * fe_face_values
std::vector< types::global_dof_index > interface_dof_indices
std::unique_ptr< FESubfaceValues< dim, spacedim > > internal_fe_subface_values_neighbor
std::vector< std::array< unsigned int, 2 > > dofmap
Triangulation< dim, spacedim >::cell_iterator get_cell(const unsigned int cell_index) const
void get_jump_in_function_hessians(const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type > > &hessians) const
double JxW(const unsigned int quadrature_point) const
Tensor< 1, spacedim > normal(const unsigned int q_point_index) const
void get_average_of_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
Tensor< 1, spacedim > jump_in_shape_gradients(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
Tensor< 2, spacedim > jump_in_shape_hessians(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
std::unique_ptr< hp::FEFaceValues< dim, spacedim > > internal_hp_fe_face_values_neighbor
std_cxx20::ranges::iota_view< unsigned int, unsigned int > quadrature_point_indices() const
FEInterfaceViews::Vector< dim, spacedim > operator[](const FEValuesExtractors::Vector &vector) const
double jump_in_shape_values(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
unsigned int get_face_number(const unsigned int cell_index) const
std::vector< types::global_dof_index > get_interface_dof_indices() const
FEInterfaceValues(const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature, const UpdateFlags update_flags)
std::array< unsigned int, 2 > interface_dof_to_dof_indices(const unsigned int interface_dof_index) const
std::unique_ptr< hp::FEFaceValues< dim, spacedim > > internal_hp_fe_face_values
std::unique_ptr< FEFaceValues< dim, spacedim > > internal_fe_face_values_neighbor
void get_jump_in_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
void get_jump_in_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &gradients) const
const Mapping< dim, spacedim > & get_mapping() const
const Quadrature< dim - 1 > & get_quadrature() const
void get_local_dof_values(const InputVector &dof_values, OutputVector &local_dof_values) const
Base(const FEInterfaceValues< dim, spacedim > &fe_interface)
const FEInterfaceValues< dim, spacedim > * fe_interface
void get_jump_in_function_third_derivatives(const InputVector &fe_function, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &third_derivatives) const
void get_jump_in_function_gradients_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
typename FEValuesViews::Scalar< dim, spacedim >::third_derivative_type third_derivative_type
gradient_type jump_in_gradients(const unsigned int interface_dof_index, const unsigned int q_point) const
void get_jump_in_function_hessians(const InputVector &fe_function, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
hessian_type jump_in_hessians(const unsigned int interface_dof_index, const unsigned int q_point) const
const FEValuesExtractors::Scalar extractor
void get_jump_in_function_third_derivatives_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &third_derivatives) const
typename ProductType< Number, gradient_type >::type solution_gradient_type
hessian_type average_of_hessians(const unsigned int interface_dof_index, const unsigned int q_point) const
void get_jump_in_function_gradients(const InputVector &fe_function, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
value_type value(const bool here_or_there, const unsigned int interface_dof_index, const unsigned int q_point) const
void get_average_of_function_hessians_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
void get_average_of_function_values_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
typename ProductType< Number, third_derivative_type >::type solution_third_derivative_type
void get_average_of_function_values(const InputVector &fe_function, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
typename ProductType< Number, value_type >::type solution_value_type
void get_average_of_function_gradients(const InputVector &fe_function, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
void get_jump_in_function_values_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
void get_jump_in_function_values(const InputVector &fe_function, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
Scalar(const FEInterfaceValues< dim, spacedim > &fe_interface, const unsigned int component)
third_derivative_type jump_in_third_derivatives(const unsigned int interface_dof_index, const unsigned int q_point) const
gradient_type average_of_gradients(const unsigned int interface_dof_index, const unsigned int q_point) const
void get_average_of_function_hessians(const InputVector &fe_function, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
typename ProductType< Number, hessian_type >::type solution_hessian_type
void get_average_of_function_gradients_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
value_type average_of_values(const unsigned int interface_dof_index, const unsigned int q_point) const
void get_jump_in_function_hessians_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
value_type jump_in_values(const unsigned int interface_dof_index, const unsigned int q_point) const
typename FEValuesViews::Scalar< dim, spacedim >::hessian_type hessian_type
void get_function_values_from_local_dof_values(const bool here_or_there, const InputVector &local_dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
void get_function_values(const bool here_or_there, const InputVector &fe_function, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
typename FEValuesViews::Scalar< dim, spacedim >::gradient_type gradient_type
const FEValuesExtractors::Vector extractor
void get_average_of_function_gradients_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
void get_average_of_function_hessians(const InputVector &fe_function, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
void get_jump_in_function_values(const InputVector &fe_function, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
typename FEValuesViews::Vector< dim, spacedim >::value_type value_type
void get_average_of_function_values_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
void get_function_values_from_local_dof_values(const bool here_or_there, const InputVector &local_dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
void get_average_of_function_hessians_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
hessian_type average_hessian(const unsigned int interface_dof_index, const unsigned int q_point) const
typename FEValuesViews::Vector< dim, spacedim >::third_derivative_type third_derivative_type
Vector(const FEInterfaceValues< dim, spacedim > &fe_interface, const unsigned int first_vector_component)
void get_average_of_function_gradients(const InputVector &fe_function, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
void get_function_values(const bool here_or_there, const InputVector &fe_function, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
hessian_type jump_hessian(const unsigned int interface_dof_index, const unsigned int q_point) const
typename ProductType< Number, value_type >::type solution_value_type
typename FEValuesViews::Vector< dim, spacedim >::hessian_type hessian_type
void get_average_of_function_values(const InputVector &fe_function, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
value_type value(const bool here_or_there, const unsigned int interface_dof_index, const unsigned int q_point) const
typename ProductType< Number, hessian_type >::type solution_hessian_type
void get_jump_in_function_third_derivatives_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &third_derivatives) const
void get_jump_in_function_gradients(const InputVector &fe_function, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
third_derivative_type jump_in_third_derivatives(const unsigned int interface_dof_index, const unsigned int q_point) const
void get_jump_in_function_gradients_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
void get_jump_in_function_hessians_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
typename ProductType< Number, third_derivative_type >::type solution_third_derivative_type
void get_jump_in_function_hessians(const InputVector &fe_function, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
value_type jump_in_values(const unsigned int interface_dof_index, const unsigned int q_point) const
typename FEValuesViews::Vector< dim, spacedim >::gradient_type gradient_type
void get_jump_in_function_third_derivatives(const InputVector &fe_function, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &third_derivatives) const
gradient_type jump_gradient(const unsigned int interface_dof_index, const unsigned int q_point) const
value_type average_of_values(const unsigned int interface_dof_index, const unsigned int q_point) const
hessian_type jump_in_hessians(const unsigned int interface_dof_index, const unsigned int q_point) const
void get_jump_in_function_values_from_local_dof_values(const InputVector &local_dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
gradient_type average_of_gradients(const unsigned int interface_dof_index, const unsigned int q_point) const
gradient_type jump_in_gradients(const unsigned int interface_dof_index, const unsigned int q_point) const
hessian_type average_of_hessians(const unsigned int interface_dof_index, const unsigned int q_point) const
typename ProductType< Number, gradient_type >::type solution_gradient_type
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell, const unsigned int face_no, const unsigned int subface_no)
Abstract base class for mapping classes.
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcOnlyAvailableWithHP()
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcOnlyAvailableWithoutHP()
static ::ExceptionBase & ExcMessage(std::string arg1)
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
typename ::internal::FEInterfaceViews:: ViewType< dim, spacedim, Extractor >::type View
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
const types::fe_index invalid_fe_index
static const unsigned int invalid_unsigned_int
boost::integer_range< IncrementableType > iota_view
typename internal::ProductTypeImpl< std::decay_t< T >, std::decay_t< U > >::type type