16#ifndef dealii_fe_interface_values_h
17#define dealii_fe_interface_values_h
31template <
int dim,
int spacedim>
44 template <
int dim,
int spacedim = dim>
63 template <
class InputVector,
class OutputVector>
66 OutputVector & local_dof_values)
const;
74 template <
int dim,
int spacedim = dim>
110 template <
typename Number>
119 template <
typename Number>
129 template <
typename Number>
139 template <
typename Number>
147 const unsigned int component);
175 const unsigned int interface_dof_index,
176 const unsigned int q_point)
const;
197 const unsigned int q_point)
const;
206 jump(
const unsigned int interface_dof_index,
207 const unsigned int q_point)
const;
221 const unsigned int q_point)
const;
231 const unsigned int q_point)
const;
245 const unsigned int q_point)
const;
255 const unsigned int q_point)
const;
269 const unsigned int q_point)
const;
279 const unsigned int q_point)
const;
300 const unsigned int q_point)
const;
310 const unsigned int q_point)
const;
319 average(
const unsigned int interface_dof_index,
320 const unsigned int q_point)
const;
334 const unsigned int q_point)
const;
344 const unsigned int q_point)
const;
359 const unsigned int q_point)
const;
369 const unsigned int q_point)
const;
397 template <
class InputVector>
400 const bool here_or_there,
401 const InputVector &fe_function,
427 template <
class InputVector>
430 const bool here_or_there,
431 const InputVector &local_dof_values,
455 template <
class InputVector>
458 const InputVector &fe_function,
469 template <
class InputVector>
472 const InputVector &local_dof_values,
489 template <
class InputVector>
492 const InputVector &fe_function,
503 template <
class InputVector>
506 const InputVector &local_dof_values,
523 template <
class InputVector>
526 const InputVector &fe_function,
537 template <
class InputVector>
540 const InputVector &local_dof_values,
558 template <
class InputVector>
561 const InputVector &fe_function,
564 &third_derivatives)
const;
573 template <
class InputVector>
576 const InputVector &local_dof_values,
579 &third_derivatives)
const;
601 template <
class InputVector>
604 const InputVector &fe_function,
615 template <
class InputVector>
618 const InputVector &local_dof_values,
635 template <
class InputVector>
638 const InputVector &fe_function,
649 template <
class InputVector>
652 const InputVector &local_dof_values,
669 template <
class InputVector>
672 const InputVector &fe_function,
683 template <
class InputVector>
686 const InputVector &local_dof_values,
704 template <
int dim,
int spacedim = dim>
743 template <
typename Number>
752 template <
typename Number>
762 template <
typename Number>
772 template <
typename Number>
780 const unsigned int first_vector_component);
808 const unsigned int interface_dof_index,
809 const unsigned int q_point)
const;
829 const unsigned int q_point)
const;
838 jump(
const unsigned int interface_dof_index,
839 const unsigned int q_point)
const;
852 const unsigned int q_point)
const;
861 const unsigned int q_point)
const;
875 const unsigned int q_point)
const;
884 const unsigned int q_point)
const;
898 const unsigned int q_point)
const;
908 const unsigned int q_point)
const;
928 const unsigned int q_point)
const;
937 average(
const unsigned int interface_dof_index,
938 const unsigned int q_point)
const;
951 const unsigned int q_point)
const;
961 const unsigned int q_point)
const;
976 const unsigned int q_point)
const;
985 const unsigned int q_point)
const;
1013 template <
class InputVector>
1016 const bool here_or_there,
1017 const InputVector &fe_function,
1043 template <
class InputVector>
1046 const bool here_or_there,
1047 const InputVector &local_dof_values,
1071 template <
class InputVector>
1074 const InputVector &fe_function,
1085 template <
class InputVector>
1088 const InputVector &local_dof_values,
1105 template <
class InputVector>
1108 const InputVector &fe_function,
1119 template <
class InputVector>
1122 const InputVector &local_dof_values,
1139 template <
class InputVector>
1142 const InputVector &fe_function,
1153 template <
class InputVector>
1156 const InputVector &local_dof_values,
1174 template <
class InputVector>
1177 const InputVector &fe_function,
1180 &third_derivatives)
const;
1189 template <
class InputVector>
1192 const InputVector &local_dof_values,
1195 &third_derivatives)
const;
1217 template <
class InputVector>
1220 const InputVector &fe_function,
1231 template <
class InputVector>
1234 const InputVector &local_dof_values,
1251 template <
class InputVector>
1254 const InputVector &fe_function,
1265 template <
class InputVector>
1268 const InputVector &local_dof_values,
1285 template <
class InputVector>
1288 const InputVector &fe_function,
1299 template <
class InputVector>
1302 const InputVector &local_dof_values,
1325 template <
int dim,
int spacedim,
typename Extractor>
1336 template <
int dim,
int spacedim>
1339 using type = typename ::FEInterfaceViews::Scalar<dim, spacedim>;
1349 template <
int dim,
int spacedim>
1352 using type = typename ::FEInterfaceViews::Vector<dim, spacedim>;
1363 template <
int dim,
int spacedim,
typename Extractor>
1364 using View = typename ::internal::FEInterfaceViews::
1365 ViewType<dim, spacedim, Extractor>::type;
1394template <
int dim,
int spacedim = dim>
1464 template <
class CellIteratorType,
class CellNeighborIteratorType>
1467 const unsigned int face_no,
1468 const unsigned int sub_face_no,
1469 const CellNeighborIteratorType &cell_neighbor,
1470 const unsigned int face_no_neighbor,
1471 const unsigned int sub_face_no_neighbor);
1484 template <
class CellIteratorType>
1486 reinit(
const CellIteratorType &cell,
const unsigned int face_no);
1567 JxW(
const unsigned int quadrature_point)
const;
1574 const std::vector<double> &
1586 const std::vector<Tensor<1, spacedim>> &
1604 const std::vector<Point<spacedim>> &
1660 std::vector<types::global_dof_index>
1675 std::array<unsigned int, 2>
1687 normal(
const unsigned int q_point_index)
const;
1721 const unsigned int interface_dof_index,
1722 const unsigned int q_point,
1723 const unsigned int component = 0)
const;
1754 const unsigned int q_point,
1755 const unsigned int component = 0)
const;
1764 jump(
const unsigned int interface_dof_index,
1765 const unsigned int q_point,
1766 const unsigned int component = 0)
const;
1783 const unsigned int q_point,
1784 const unsigned int component = 0)
const;
1794 const unsigned int q_point,
1795 const unsigned int component = 0)
const;
1813 const unsigned int q_point,
1814 const unsigned int component = 0)
const;
1824 const unsigned int q_point,
1825 const unsigned int component = 0)
const;
1842 const unsigned int q_point,
1843 const unsigned int component = 0)
const;
1853 const unsigned int q_point,
1854 const unsigned int component = 0)
const;
1880 const unsigned int q_point,
1881 const unsigned int component = 0)
const;
1891 const unsigned int q_point,
1892 const unsigned int component = 0)
const;
1909 const unsigned int q_point,
1910 const unsigned int component = 0)
const;
1920 const unsigned int q_point,
1921 const unsigned int component = 0)
const;
1939 const unsigned int q_point,
1940 const unsigned int component = 0)
const;
1950 const unsigned int q_point,
1951 const unsigned int component = 0)
const;
1972 template <
class InputVector>
1975 const InputVector & fe_function,
1976 std::vector<typename InputVector::value_type> &values)
const;
1986 template <
class InputVector>
1989 const InputVector &fe_function,
2000 template <
class InputVector>
2003 const InputVector &fe_function,
2015 template <
class InputVector>
2018 const InputVector &fe_function,
2020 &third_derivatives)
const;
2037 template <
class InputVector>
2040 const InputVector & fe_function,
2041 std::vector<typename InputVector::value_type> &values)
const;
2050 template <
class InputVector>
2053 const InputVector &fe_function,
2064 template <
class InputVector>
2067 const InputVector &fe_function,
2115 std::vector<std::array<unsigned int, 2>>
dofmap;
2165template <
int dim,
int spacedim>
2171 : n_quadrature_points(quadrature.size())
2172 , internal_fe_face_values(mapping, fe, quadrature, update_flags)
2173 , internal_fe_subface_values(mapping, fe, quadrature, update_flags)
2174 , internal_fe_face_values_neighbor(mapping, fe, quadrature, update_flags)
2175 , internal_fe_subface_values_neighbor(mapping, fe, quadrature, update_flags)
2176 , fe_face_values(nullptr)
2177 , fe_face_values_neighbor(nullptr)
2182template <
int dim,
int spacedim>
2188 : n_quadrature_points(quadrature.max_n_quadrature_points())
2189 , internal_fe_face_values(mapping, fe, quadrature, update_flags)
2190 , internal_fe_subface_values(mapping, fe, quadrature, update_flags)
2191 , internal_fe_face_values_neighbor(mapping, fe, quadrature[0], update_flags)
2192 , internal_fe_subface_values_neighbor(mapping,
2196 , fe_face_values(nullptr)
2197 , fe_face_values_neighbor(nullptr)
2202template <
int dim,
int spacedim>
2207 : n_quadrature_points(quadrature.size())
2208 , internal_fe_face_values(
2213 , internal_fe_subface_values(
2218 , internal_fe_face_values_neighbor(
2223 , internal_fe_subface_values_neighbor(
2228 , fe_face_values(nullptr)
2229 , fe_face_values_neighbor(nullptr)
2234template <
int dim,
int spacedim>
2235template <
class CellIteratorType,
class CellNeighborIteratorType>
2238 const CellIteratorType & cell,
2239 const unsigned int face_no,
2240 const unsigned int sub_face_no,
2241 const CellNeighborIteratorType &cell_neighbor,
2242 const unsigned int face_no_neighbor,
2243 const unsigned int sub_face_no_neighbor)
2247 internal_fe_face_values.reinit(cell, face_no);
2248 fe_face_values = &internal_fe_face_values;
2252 internal_fe_subface_values.reinit(cell, face_no, sub_face_no);
2253 fe_face_values = &internal_fe_subface_values;
2257 internal_fe_face_values_neighbor.reinit(cell_neighbor, face_no_neighbor);
2258 fe_face_values_neighbor = &internal_fe_face_values_neighbor;
2262 internal_fe_subface_values_neighbor.reinit(cell_neighbor,
2264 sub_face_no_neighbor);
2265 fe_face_values_neighbor = &internal_fe_subface_values_neighbor;
2269 fe_face_values_neighbor->n_quadrature_points);
2271 const_cast<unsigned int &
>(this->n_quadrature_points) =
2272 fe_face_values->n_quadrature_points;
2277 std::vector<types::global_dof_index> v(
2278 fe_face_values->get_fe().n_dofs_per_cell());
2279 cell->get_active_or_mg_dof_indices(v);
2280 std::vector<types::global_dof_index> v2(
2281 fe_face_values_neighbor->get_fe().n_dofs_per_cell());
2282 cell_neighbor->get_active_or_mg_dof_indices(v2);
2286 std::map<types::global_dof_index, std::pair<unsigned int, unsigned int>>
2288 std::pair<unsigned int, unsigned int> invalid_entry(
2291 for (
unsigned int i = 0; i < v.size(); ++i)
2294 auto result = tempmap.insert(std::make_pair(v[i], invalid_entry));
2295 result.first->second.first = i;
2298 for (
unsigned int i = 0; i < v2.size(); ++i)
2301 auto result = tempmap.insert(std::make_pair(v2[i], invalid_entry));
2302 result.first->second.second = i;
2306 dofmap.resize(tempmap.size());
2307 interface_dof_indices.resize(tempmap.size());
2308 unsigned int idx = 0;
2309 for (
auto &x : tempmap)
2311 interface_dof_indices[idx] = x.first;
2312 dofmap[idx] = {{x.second.first, x.second.second}};
2320template <
int dim,
int spacedim>
2321template <
class CellIteratorType>
2324 const unsigned int face_no)
2326 internal_fe_face_values.reinit(cell, face_no);
2327 fe_face_values = &internal_fe_face_values;
2328 fe_face_values_neighbor =
nullptr;
2330 interface_dof_indices.resize(fe_face_values->get_fe().n_dofs_per_cell());
2331 cell->get_active_or_mg_dof_indices(interface_dof_indices);
2333 dofmap.resize(interface_dof_indices.size());
2335 for (
unsigned int i = 0; i < interface_dof_indices.size(); ++i)
2343template <
int dim,
int spacedim>
2347 Assert(fe_face_values !=
nullptr,
2348 ExcMessage(
"This call requires a call to reinit() first."));
2349 return fe_face_values->JxW(q);
2354template <
int dim,
int spacedim>
2355const std::vector<double> &
2358 Assert(fe_face_values !=
nullptr,
2359 ExcMessage(
"This call requires a call to reinit() first."));
2360 return fe_face_values->get_JxW_values();
2365template <
int dim,
int spacedim>
2366const std::vector<Tensor<1, spacedim>> &
2369 Assert(fe_face_values !=
nullptr,
2370 ExcMessage(
"This call requires a call to reinit() first."));
2371 return fe_face_values->get_normal_vectors();
2376template <
int dim,
int spacedim>
2380 return internal_fe_face_values.get_mapping();
2385template <
int dim,
int spacedim>
2389 return internal_fe_face_values.get_fe();
2394template <
int dim,
int spacedim>
2398 return internal_fe_face_values.get_quadrature();
2403template <
int dim,
int spacedim>
2407 return {0U, n_quadrature_points};
2412template <
int dim,
int spacedim>
2413const std::vector<Point<spacedim>> &
2416 Assert(fe_face_values !=
nullptr,
2417 ExcMessage(
"This call requires a call to reinit() first."));
2418 return fe_face_values->get_quadrature_points();
2423template <
int dim,
int spacedim>
2427 return internal_fe_face_values.get_update_flags();
2432template <
int dim,
int spacedim>
2436 return get_fe_face_values(
cell_index).get_cell();
2441template <
int dim,
int spacedim>
2446 return get_fe_face_values(
cell_index).get_face_number();
2451template <
int dim,
int spacedim>
2456 interface_dof_indices.size() > 0,
2458 "n_current_interface_dofs() is only available after a call to reinit()."));
2459 return interface_dof_indices.size();
2464template <
int dim,
int spacedim>
2468 return {0U, n_current_interface_dofs()};
2473template <
int dim,
int spacedim>
2477 return fe_face_values_neighbor ==
nullptr;
2482template <
int dim,
int spacedim>
2483std::vector<types::global_dof_index>
2486 return interface_dof_indices;
2491template <
int dim,
int spacedim>
2492std::array<unsigned int, 2>
2494 const unsigned int interface_dof_index)
const
2497 return dofmap[interface_dof_index];
2502template <
int dim,
int spacedim>
2511 "You are on a boundary, so you can only ask for the first FEFaceValues object."));
2513 return (
cell_index == 0) ? *fe_face_values : *fe_face_values_neighbor;
2518template <
int dim,
int spacedim>
2522 return fe_face_values->normal_vector(q_point_index);
2527template <
int dim,
int spacedim>
2530 const bool here_or_there,
2531 const unsigned int interface_dof_index,
2532 const unsigned int q_point,
2533 const unsigned int component)
const
2535 const auto dof_pair = dofmap[interface_dof_index];
2538 return get_fe_face_values(0).shape_value_component(dof_pair[0],
2542 return get_fe_face_values(1).shape_value_component(dof_pair[1],
2551template <
int dim,
int spacedim>
2554 const unsigned int interface_dof_index,
2555 const unsigned int q_point,
2556 const unsigned int component)
const
2558 const auto dof_pair = dofmap[interface_dof_index];
2563 value += get_fe_face_values(0).shape_value_component(dof_pair[0],
2567 value -= get_fe_face_values(1).shape_value_component(dof_pair[1],
2575template <
int dim,
int spacedim>
2578 const unsigned int q_point,
2579 const unsigned int component)
const
2581 return jump_in_shape_values(interface_dof_index, q_point, component);
2586template <
int dim,
int spacedim>
2589 const unsigned int interface_dof_index,
2590 const unsigned int q_point,
2591 const unsigned int component)
const
2593 const auto dof_pair = dofmap[interface_dof_index];
2596 return get_fe_face_values(0).shape_value_component(dof_pair[0],
2603 value += 0.5 * get_fe_face_values(0).shape_value_component(dof_pair[0],
2607 value += 0.5 * get_fe_face_values(1).shape_value_component(dof_pair[1],
2616template <
int dim,
int spacedim>
2619 const unsigned int interface_dof_index,
2620 const unsigned int q_point,
2621 const unsigned int component)
const
2623 return average_of_shape_values(interface_dof_index, q_point, component);
2628template <
int dim,
int spacedim>
2631 const unsigned int interface_dof_index,
2632 const unsigned int q_point,
2633 const unsigned int component)
const
2635 const auto dof_pair = dofmap[interface_dof_index];
2638 return get_fe_face_values(0).shape_grad_component(dof_pair[0],
2645 value += 0.5 * get_fe_face_values(0).shape_grad_component(dof_pair[0],
2649 value += 0.5 * get_fe_face_values(1).shape_grad_component(dof_pair[1],
2658template <
int dim,
int spacedim>
2661 const unsigned int interface_dof_index,
2662 const unsigned int q_point,
2663 const unsigned int component)
const
2665 return average_of_shape_gradients(interface_dof_index, q_point, component);
2670template <
int dim,
int spacedim>
2673 const unsigned int interface_dof_index,
2674 const unsigned int q_point,
2675 const unsigned int component)
const
2677 const auto dof_pair = dofmap[interface_dof_index];
2680 return get_fe_face_values(0).shape_hessian_component(dof_pair[0],
2687 value += 0.5 * get_fe_face_values(0).shape_hessian_component(dof_pair[0],
2691 value += 0.5 * get_fe_face_values(1).shape_hessian_component(dof_pair[1],
2700template <
int dim,
int spacedim>
2703 const unsigned int interface_dof_index,
2704 const unsigned int q_point,
2705 const unsigned int component)
const
2707 return average_of_shape_hessians(interface_dof_index, q_point, component);
2712template <
int dim,
int spacedim>
2715 const unsigned int interface_dof_index,
2716 const unsigned int q_point,
2717 const unsigned int component)
const
2719 const auto dof_pair = dofmap[interface_dof_index];
2722 return get_fe_face_values(0).shape_grad_component(dof_pair[0],
2729 value += get_fe_face_values(0).shape_grad_component(dof_pair[0],
2733 value -= get_fe_face_values(1).shape_grad_component(dof_pair[1],
2742template <
int dim,
int spacedim>
2745 const unsigned int interface_dof_index,
2746 const unsigned int q_point,
2747 const unsigned int component)
const
2749 return jump_in_shape_gradients(interface_dof_index, q_point, component);
2754template <
int dim,
int spacedim>
2757 const unsigned int interface_dof_index,
2758 const unsigned int q_point,
2759 const unsigned int component)
const
2761 const auto dof_pair = dofmap[interface_dof_index];
2764 return get_fe_face_values(0).shape_hessian_component(dof_pair[0],
2771 value += get_fe_face_values(0).shape_hessian_component(dof_pair[0],
2775 value -= get_fe_face_values(1).shape_hessian_component(dof_pair[1],
2784template <
int dim,
int spacedim>
2787 const unsigned int interface_dof_index,
2788 const unsigned int q_point,
2789 const unsigned int component)
const
2791 return jump_in_shape_hessians(interface_dof_index, q_point, component);
2796template <
int dim,
int spacedim>
2799 const unsigned int interface_dof_index,
2800 const unsigned int q_point,
2801 const unsigned int component)
const
2803 const auto dof_pair = dofmap[interface_dof_index];
2806 return get_fe_face_values(0).shape_3rd_derivative_component(dof_pair[0],
2813 value += get_fe_face_values(0).shape_3rd_derivative_component(dof_pair[0],
2817 value -= get_fe_face_values(1).shape_3rd_derivative_component(dof_pair[1],
2826template <
int dim,
int spacedim>
2829 const unsigned int interface_dof_index,
2830 const unsigned int q_point,
2831 const unsigned int component)
const
2833 return jump_in_shape_3rd_derivatives(interface_dof_index, q_point, component);
2838template <
int dim,
int spacedim>
2839template <
class InputVector>
2842 const InputVector & fe_function,
2843 std::vector<typename InputVector::value_type> &values)
const
2848 this->operator[](scalar).get_jump_in_function_values(fe_function, values);
2853template <
int dim,
int spacedim>
2854template <
class InputVector>
2857 const InputVector &fe_function,
2864 this->operator[](scalar).get_jump_in_function_gradients(fe_function,
2870template <
int dim,
int spacedim>
2871template <
class InputVector>
2874 const InputVector &fe_function,
2881 this->operator[](scalar).get_jump_in_function_hessians(fe_function, hessians);
2886template <
int dim,
int spacedim>
2887template <
class InputVector>
2890 const InputVector &fe_function,
2892 &third_derivatives)
const
2897 this->operator[](scalar).get_jump_in_function_third_derivatives(
2898 fe_function, third_derivatives);
2903template <
int dim,
int spacedim>
2904template <
class InputVector>
2907 const InputVector & fe_function,
2908 std::vector<typename InputVector::value_type> &values)
const
2913 this->operator[](scalar).get_average_of_function_values(fe_function, values);
2918template <
int dim,
int spacedim>
2919template <
class InputVector>
2922 const InputVector &fe_function,
2929 this->operator[](scalar).get_average_of_function_gradients(fe_function,
2935template <
int dim,
int spacedim>
2936template <
class InputVector>
2939 const InputVector &fe_function,
2946 this->operator[](scalar).get_average_of_function_hessians(fe_function,
2953template <
int dim,
int spacedim>
2964template <
int dim,
int spacedim>
2970 const unsigned int n_vectors =
2984 template <
int dim,
int spacedim>
2987 : fe_interface(&fe_interface)
2992 template <
int dim,
int spacedim>
2993 Scalar<dim, spacedim>::Scalar(
2995 const unsigned int component)
2996 : Base<dim, spacedim>(fe_interface)
2997 , extractor(component)
3002 template <
int dim,
int spacedim>
3003 template <
class InputVector,
class OutputVector>
3005 Base<dim, spacedim>::get_local_dof_values(
3006 const InputVector &dof_values,
3007 OutputVector & local_dof_values)
const
3009 const auto &interface_dof_indices =
3010 this->fe_interface->get_interface_dof_indices();
3012 AssertDimension(interface_dof_indices.size(), local_dof_values.size());
3014 for (
unsigned int i : this->fe_interface->dof_indices())
3015 local_dof_values[i] = dof_values(interface_dof_indices[i]);
3020 template <
int dim,
int spacedim>
3021 typename Scalar<dim, spacedim>::value_type
3022 Scalar<dim, spacedim>::value(
const bool here_or_there,
3023 const unsigned int interface_dof_index,
3024 const unsigned int q_point)
const
3026 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3029 return (*(this->fe_interface->fe_face_values))[extractor].value(
3030 dof_pair[0], q_point);
3033 return (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3034 dof_pair[1], q_point);
3041 template <
int dim,
int spacedim>
3042 typename Scalar<dim, spacedim>::value_type
3043 Scalar<dim, spacedim>::jump_in_values(
const unsigned int interface_dof_index,
3044 const unsigned int q_point)
const
3046 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3048 value_type
value = 0.0;
3052 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3057 (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3058 dof_pair[1], q_point);
3065 template <
int dim,
int spacedim>
3066 typename Scalar<dim, spacedim>::value_type
3067 Scalar<dim, spacedim>::jump(
const unsigned int interface_dof_index,
3068 const unsigned int q_point)
const
3070 return jump_in_values(interface_dof_index, q_point);
3075 template <
int dim,
int spacedim>
3076 typename Scalar<dim, spacedim>::value_type
3077 Scalar<dim, spacedim>::average_of_values(
3078 const unsigned int interface_dof_index,
3079 const unsigned int q_point)
const
3081 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3083 if (this->fe_interface->at_boundary())
3084 return (*(this->fe_interface->fe_face_values))[extractor].value(
3085 dof_pair[0], q_point);
3087 value_type
value = 0.0;
3092 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3097 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3098 dof_pair[1], q_point);
3105 template <
int dim,
int spacedim>
3106 typename Scalar<dim, spacedim>::value_type
3107 Scalar<dim, spacedim>::average(
const unsigned int interface_dof_index,
3108 const unsigned int q_point)
const
3110 return average_of_values(interface_dof_index, q_point);
3115 template <
int dim,
int spacedim>
3116 typename Scalar<dim, spacedim>::gradient_type
3117 Scalar<dim, spacedim>::average_of_gradients(
3118 const unsigned int interface_dof_index,
3119 const unsigned int q_point)
const
3121 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3123 if (this->fe_interface->at_boundary())
3124 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3125 dof_pair[0], q_point);
3127 gradient_type
value;
3132 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3136 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3137 .gradient(dof_pair[1], q_point);
3143 template <
int dim,
int spacedim>
3144 typename Scalar<dim, spacedim>::gradient_type
3145 Scalar<dim, spacedim>::average_gradient(
3146 const unsigned int interface_dof_index,
3147 const unsigned int q_point)
const
3149 return average_of_gradients(interface_dof_index, q_point);
3154 template <
int dim,
int spacedim>
3155 typename Scalar<dim, spacedim>::gradient_type
3156 Scalar<dim, spacedim>::jump_in_gradients(
3157 const unsigned int interface_dof_index,
3158 const unsigned int q_point)
const
3160 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3162 if (this->fe_interface->at_boundary())
3163 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3164 dof_pair[0], q_point);
3166 gradient_type
value;
3170 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3175 (*(this->fe_interface->fe_face_values_neighbor))[extractor].gradient(
3176 dof_pair[1], q_point);
3183 template <
int dim,
int spacedim>
3184 typename Scalar<dim, spacedim>::gradient_type
3185 Scalar<dim, spacedim>::jump_gradient(
const unsigned int interface_dof_index,
3186 const unsigned int q_point)
const
3188 return jump_in_gradients(interface_dof_index, q_point);
3193 template <
int dim,
int spacedim>
3194 typename Scalar<dim, spacedim>::hessian_type
3195 Scalar<dim, spacedim>::average_of_hessians(
3196 const unsigned int interface_dof_index,
3197 const unsigned int q_point)
const
3199 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3201 if (this->fe_interface->at_boundary())
3202 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
3203 dof_pair[0], q_point);
3210 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
3214 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3215 .hessian(dof_pair[1], q_point);
3222 template <
int dim,
int spacedim>
3223 typename Scalar<dim, spacedim>::hessian_type
3224 Scalar<dim, spacedim>::average_hessian(
const unsigned int interface_dof_index,
3225 const unsigned int q_point)
const
3227 return average_of_hessians(interface_dof_index, q_point);
3232 template <
int dim,
int spacedim>
3233 typename Scalar<dim, spacedim>::third_derivative_type
3234 Scalar<dim, spacedim>::jump_in_third_derivatives(
3235 const unsigned int interface_dof_index,
3236 const unsigned int q_point)
const
3238 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3240 if (this->fe_interface->at_boundary())
3241 return (*(this->fe_interface->fe_face_values))[extractor]
3242 .third_derivative(dof_pair[0], q_point);
3244 third_derivative_type
value;
3248 (*(this->fe_interface->fe_face_values))[extractor].third_derivative(
3249 dof_pair[0], q_point);
3252 value -= (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3253 .third_derivative(dof_pair[1], q_point);
3260 template <
int dim,
int spacedim>
3261 typename Scalar<dim, spacedim>::third_derivative_type
3262 Scalar<dim, spacedim>::jump_3rd_derivative(
3263 const unsigned int interface_dof_index,
3264 const unsigned int q_point)
const
3266 return jump_in_third_derivatives(interface_dof_index, q_point);
3271 template <
int dim,
int spacedim>
3272 typename Scalar<dim, spacedim>::hessian_type
3273 Scalar<dim, spacedim>::jump_in_hessians(
3274 const unsigned int interface_dof_index,
3275 const unsigned int q_point)
const
3277 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3279 if (this->fe_interface->at_boundary())
3280 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
3281 dof_pair[0], q_point);
3287 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
3292 (*(this->fe_interface->fe_face_values_neighbor))[extractor].hessian(
3293 dof_pair[1], q_point);
3300 template <
int dim,
int spacedim>
3301 typename Scalar<dim, spacedim>::hessian_type
3302 Scalar<dim, spacedim>::jump_hessian(
const unsigned int interface_dof_index,
3303 const unsigned int q_point)
const
3305 return jump_in_hessians(interface_dof_index, q_point);
3310 template <
int dim,
int spacedim>
3311 template <
class InputVector>
3313 Scalar<dim, spacedim>::get_function_values_from_local_dof_values(
3314 const bool here_or_there,
3315 const InputVector &local_dof_values,
3316 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3321 for (
const auto dof_index : this->fe_interface->dof_indices())
3322 for (
const auto q_index : this->fe_interface->quadrature_point_indices())
3327 values[q_index] += local_dof_values[dof_index] *
3328 value(here_or_there, dof_index, q_index);
3334 template <
int dim,
int spacedim>
3335 template <
class InputVector>
3337 Scalar<dim, spacedim>::get_function_values(
3338 const bool here_or_there,
3339 const InputVector &fe_function,
3340 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3343 std::vector<typename InputVector::value_type> local_dof_values(
3344 this->fe_interface->n_current_interface_dofs());
3345 this->get_local_dof_values(fe_function, local_dof_values);
3347 get_function_values_from_local_dof_values(here_or_there,
3354 template <
int dim,
int spacedim>
3355 template <
class InputVector>
3357 Scalar<dim, spacedim>::get_jump_in_function_values_from_local_dof_values(
3358 const InputVector &local_dof_values,
3359 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3364 for (
const auto dof_index : this->fe_interface->dof_indices())
3365 for (
const auto q_index : this->fe_interface->quadrature_point_indices())
3371 local_dof_values[dof_index] * jump_in_values(dof_index, q_index);
3377 template <
int dim,
int spacedim>
3378 template <
class InputVector>
3380 Scalar<dim, spacedim>::get_jump_in_function_values(
3381 const InputVector &fe_function,
3382 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3385 std::vector<typename InputVector::value_type> local_dof_values(
3386 this->fe_interface->n_current_interface_dofs());
3387 this->get_local_dof_values(fe_function, local_dof_values);
3389 get_jump_in_function_values_from_local_dof_values(local_dof_values, values);
3394 template <
int dim,
int spacedim>
3395 template <
class InputVector>
3397 Scalar<dim, spacedim>::get_jump_in_function_gradients_from_local_dof_values(
3398 const InputVector &local_dof_values,
3399 std::vector<solution_gradient_type<typename InputVector::value_type>>
3404 for (
const auto dof_index : this->fe_interface->dof_indices())
3405 for (
const auto q_index : this->fe_interface->quadrature_point_indices())
3411 local_dof_values[dof_index] * jump_in_gradients(dof_index, q_index);
3417 template <
int dim,
int spacedim>
3418 template <
class InputVector>
3420 Scalar<dim, spacedim>::get_jump_in_function_gradients(
3421 const InputVector &fe_function,
3422 std::vector<solution_gradient_type<typename InputVector::value_type>>
3425 std::vector<typename InputVector::value_type> local_dof_values(
3426 this->fe_interface->n_current_interface_dofs());
3427 this->get_local_dof_values(fe_function, local_dof_values);
3429 get_jump_in_function_gradients_from_local_dof_values(local_dof_values,
3435 template <
int dim,
int spacedim>
3436 template <
class InputVector>
3438 Scalar<dim, spacedim>::get_average_of_function_values_from_local_dof_values(
3439 const InputVector &local_dof_values,
3440 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3445 for (
const auto dof_index : this->fe_interface->dof_indices())
3446 for (
const auto q_index : this->fe_interface->quadrature_point_indices())
3452 local_dof_values[dof_index] * average_of_values(dof_index, q_index);
3458 template <
int dim,
int spacedim>
3459 template <
class InputVector>
3461 Scalar<dim, spacedim>::get_average_of_function_values(
3462 const InputVector &fe_function,
3463 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3466 std::vector<typename InputVector::value_type> local_dof_values(
3467 this->fe_interface->n_current_interface_dofs());
3468 this->get_local_dof_values(fe_function, local_dof_values);
3470 get_average_of_function_values_from_local_dof_values(local_dof_values,
3476 template <
int dim,
int spacedim>
3477 template <
class InputVector>
3479 Scalar<dim, spacedim>::
3480 get_average_of_function_gradients_from_local_dof_values(
3481 const InputVector &local_dof_values,
3482 std::vector<solution_gradient_type<typename InputVector::value_type>>
3487 for (
const auto dof_index : this->fe_interface->dof_indices())
3488 for (
const auto q_index : this->fe_interface->quadrature_point_indices())
3493 gradients[q_index] += local_dof_values[dof_index] *
3494 average_of_gradients(dof_index, q_index);
3500 template <
int dim,
int spacedim>
3501 template <
class InputVector>
3503 Scalar<dim, spacedim>::get_average_of_function_gradients(
3504 const InputVector &fe_function,
3505 std::vector<solution_gradient_type<typename InputVector::value_type>>
3508 std::vector<typename InputVector::value_type> local_dof_values(
3509 this->fe_interface->n_current_interface_dofs());
3510 this->get_local_dof_values(fe_function, local_dof_values);
3512 get_average_of_function_gradients_from_local_dof_values(local_dof_values,
3518 template <
int dim,
int spacedim>
3519 template <
class InputVector>
3521 Scalar<dim, spacedim>::get_jump_in_function_hessians_from_local_dof_values(
3522 const InputVector &local_dof_values,
3523 std::vector<solution_hessian_type<typename InputVector::value_type>>
3528 for (
const auto dof_index : this->fe_interface->dof_indices())
3529 for (
const auto q_index : this->fe_interface->quadrature_point_indices())
3535 local_dof_values[dof_index] * jump_in_hessians(dof_index, q_index);
3541 template <
int dim,
int spacedim>
3542 template <
class InputVector>
3544 Scalar<dim, spacedim>::get_jump_in_function_hessians(
3545 const InputVector &fe_function,
3546 std::vector<solution_hessian_type<typename InputVector::value_type>>
3549 std::vector<typename InputVector::value_type> local_dof_values(
3550 this->fe_interface->n_current_interface_dofs());
3551 this->get_local_dof_values(fe_function, local_dof_values);
3553 get_jump_in_function_hessians_from_local_dof_values(local_dof_values,
3559 template <
int dim,
int spacedim>
3560 template <
class InputVector>
3562 Scalar<dim, spacedim>::get_average_of_function_hessians_from_local_dof_values(
3563 const InputVector &local_dof_values,
3564 std::vector<solution_hessian_type<typename InputVector::value_type>>
3569 for (
unsigned int dof_index : this->fe_interface->dof_indices())
3570 for (
const auto q_index : this->fe_interface->quadrature_point_indices())
3575 hessians[q_index] += local_dof_values[dof_index] *
3576 average_of_hessians(dof_index, q_index);
3582 template <
int dim,
int spacedim>
3583 template <
class InputVector>
3585 Scalar<dim, spacedim>::get_average_of_function_hessians(
3586 const InputVector &fe_function,
3587 std::vector<solution_hessian_type<typename InputVector::value_type>>
3590 std::vector<typename InputVector::value_type> local_dof_values(
3591 this->fe_interface->n_current_interface_dofs());
3592 this->get_local_dof_values(fe_function, local_dof_values);
3594 get_average_of_function_hessians_from_local_dof_values(local_dof_values,
3600 template <
int dim,
int spacedim>
3601 template <
class InputVector>
3603 Scalar<dim, spacedim>::
3604 get_jump_in_function_third_derivatives_from_local_dof_values(
3605 const InputVector &local_dof_values,
3607 solution_third_derivative_type<typename InputVector::value_type>>
3608 &third_derivatives)
const
3611 this->fe_interface->n_quadrature_points);
3613 for (
unsigned int dof_index : this->fe_interface->dof_indices())
3614 for (
const auto q_index : this->fe_interface->quadrature_point_indices())
3617 third_derivatives[q_index] = 0.;
3619 third_derivatives[q_index] +=
3620 local_dof_values[dof_index] *
3621 jump_in_third_derivatives(dof_index, q_index);
3627 template <
int dim,
int spacedim>
3628 template <
class InputVector>
3630 Scalar<dim, spacedim>::get_jump_in_function_third_derivatives(
3631 const InputVector &fe_function,
3633 solution_third_derivative_type<typename InputVector::value_type>>
3634 &third_derivatives)
const
3636 std::vector<typename InputVector::value_type> local_dof_values(
3637 this->fe_interface->n_current_interface_dofs());
3638 this->get_local_dof_values(fe_function, local_dof_values);
3640 get_jump_in_function_third_derivatives_from_local_dof_values(
3641 local_dof_values, third_derivatives);
3646 template <
int dim,
int spacedim>
3649 const unsigned int first_vector_component)
3650 : Base<dim, spacedim>(fe_interface)
3651 , extractor(first_vector_component)
3656 template <
int dim,
int spacedim>
3659 const unsigned int interface_dof_index,
3660 const unsigned int q_point)
const
3662 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3665 return (*(this->fe_interface->fe_face_values))[extractor].value(
3666 dof_pair[0], q_point);
3669 return (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3670 dof_pair[1], q_point);
3672 return value_type();
3677 template <
int dim,
int spacedim>
3680 const unsigned int q_point)
const
3682 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3688 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3693 (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3694 dof_pair[1], q_point);
3701 template <
int dim,
int spacedim>
3704 const unsigned int q_point)
const
3706 return jump_in_values(interface_dof_index, q_point);
3711 template <
int dim,
int spacedim>
3714 const unsigned int interface_dof_index,
3715 const unsigned int q_point)
const
3717 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3719 if (this->fe_interface->at_boundary())
3720 return (*(this->fe_interface->fe_face_values))[extractor].value(
3721 dof_pair[0], q_point);
3728 (*(this->fe_interface->fe_face_values))[extractor].value(dof_pair[0],
3733 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor].value(
3734 dof_pair[1], q_point);
3741 template <
int dim,
int spacedim>
3744 const unsigned int q_point)
const
3746 return average_of_values(interface_dof_index, q_point);
3751 template <
int dim,
int spacedim>
3754 const unsigned int interface_dof_index,
3755 const unsigned int q_point)
const
3757 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3759 if (this->fe_interface->at_boundary())
3760 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3761 dof_pair[0], q_point);
3763 gradient_type
value;
3768 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3772 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3773 .gradient(dof_pair[1], q_point);
3780 template <
int dim,
int spacedim>
3783 const unsigned int interface_dof_index,
3784 const unsigned int q_point)
const
3786 return average_of_gradients(interface_dof_index, q_point);
3791 template <
int dim,
int spacedim>
3794 const unsigned int interface_dof_index,
3795 const unsigned int q_point)
const
3797 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3799 if (this->fe_interface->at_boundary())
3800 return (*(this->fe_interface->fe_face_values))[extractor].gradient(
3801 dof_pair[0], q_point);
3803 gradient_type
value;
3807 (*(this->fe_interface->fe_face_values))[extractor].gradient(dof_pair[0],
3812 (*(this->fe_interface->fe_face_values_neighbor))[extractor].gradient(
3813 dof_pair[1], q_point);
3820 template <
int dim,
int spacedim>
3823 const unsigned int q_point)
const
3825 return jump_in_gradients(interface_dof_index, q_point);
3830 template <
int dim,
int spacedim>
3833 const unsigned int interface_dof_index,
3834 const unsigned int q_point)
const
3836 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3838 if (this->fe_interface->at_boundary())
3839 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
3840 dof_pair[0], q_point);
3847 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
3851 value += 0.5 * (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3852 .hessian(dof_pair[1], q_point);
3859 template <
int dim,
int spacedim>
3862 const unsigned int q_point)
const
3864 return average_of_hessians(interface_dof_index, q_point);
3869 template <
int dim,
int spacedim>
3872 const unsigned int interface_dof_index,
3873 const unsigned int q_point)
const
3875 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3877 if (this->fe_interface->at_boundary())
3878 return (*(this->fe_interface->fe_face_values))[extractor].hessian(
3879 dof_pair[0], q_point);
3885 (*(this->fe_interface->fe_face_values))[extractor].hessian(dof_pair[0],
3890 (*(this->fe_interface->fe_face_values_neighbor))[extractor].hessian(
3891 dof_pair[1], q_point);
3898 template <
int dim,
int spacedim>
3901 const unsigned int q_point)
const
3903 return jump_in_hessians(interface_dof_index, q_point);
3908 template <
int dim,
int spacedim>
3911 const unsigned int interface_dof_index,
3912 const unsigned int q_point)
const
3914 const auto dof_pair = this->fe_interface->dofmap[interface_dof_index];
3916 if (this->fe_interface->at_boundary())
3917 return (*(this->fe_interface->fe_face_values))[extractor]
3918 .third_derivative(dof_pair[0], q_point);
3920 third_derivative_type
value;
3924 (*(this->fe_interface->fe_face_values))[extractor].third_derivative(
3925 dof_pair[0], q_point);
3928 value -= (*(this->fe_interface->fe_face_values_neighbor))[extractor]
3929 .third_derivative(dof_pair[1], q_point);
3936 template <
int dim,
int spacedim>
3939 const unsigned int interface_dof_index,
3940 const unsigned int q_point)
const
3942 return jump_in_third_derivatives(interface_dof_index, q_point);
3947 template <
int dim,
int spacedim>
3948 template <
class InputVector>
3951 const bool here_or_there,
3952 const InputVector &local_dof_values,
3953 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3958 for (
const auto dof_index : this->fe_interface->dof_indices())
3959 for (
const auto q_index : this->fe_interface->quadrature_point_indices())
3964 values[q_index] += local_dof_values[dof_index] *
3965 value(here_or_there, dof_index, q_index);
3971 template <
int dim,
int spacedim>
3972 template <
class InputVector>
3975 const bool here_or_there,
3976 const InputVector &fe_function,
3977 std::vector<solution_value_type<typename InputVector::value_type>> &values)
3980 std::vector<typename InputVector::value_type> local_dof_values(
3981 this->fe_interface->n_current_interface_dofs());
3982 this->get_local_dof_values(fe_function, local_dof_values);
3984 get_function_values_from_local_dof_values(here_or_there,
3991 template <
int dim,
int spacedim>
3992 template <
class InputVector>
3995 const InputVector &local_dof_values,
3996 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4001 for (
const auto dof_index : this->fe_interface->dof_indices())
4002 for (
const auto q_index : this->fe_interface->quadrature_point_indices())
4008 local_dof_values[dof_index] * jump_in_values(dof_index, q_index);
4014 template <
int dim,
int spacedim>
4015 template <
class InputVector>
4018 const InputVector &fe_function,
4019 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4022 std::vector<typename InputVector::value_type> local_dof_values(
4023 this->fe_interface->n_current_interface_dofs());
4024 this->get_local_dof_values(fe_function, local_dof_values);
4026 get_jump_in_function_values_from_local_dof_values(local_dof_values, values);
4031 template <
int dim,
int spacedim>
4032 template <
class InputVector>
4035 const InputVector &local_dof_values,
4036 std::vector<solution_gradient_type<typename InputVector::value_type>>
4041 for (
const auto dof_index : this->fe_interface->dof_indices())
4042 for (
const auto q_index : this->fe_interface->quadrature_point_indices())
4048 local_dof_values[dof_index] * jump_in_gradients(dof_index, q_index);
4054 template <
int dim,
int spacedim>
4055 template <
class InputVector>
4058 const InputVector &fe_function,
4059 std::vector<solution_gradient_type<typename InputVector::value_type>>
4062 std::vector<typename InputVector::value_type> local_dof_values(
4063 this->fe_interface->n_current_interface_dofs());
4064 this->get_local_dof_values(fe_function, local_dof_values);
4066 get_jump_in_function_gradients_from_local_dof_values(local_dof_values,
4072 template <
int dim,
int spacedim>
4073 template <
class InputVector>
4076 const InputVector &local_dof_values,
4077 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4082 for (
const auto dof_index : this->fe_interface->dof_indices())
4083 for (
const auto q_index : this->fe_interface->quadrature_point_indices())
4089 local_dof_values[dof_index] * average_of_values(dof_index, q_index);
4095 template <
int dim,
int spacedim>
4096 template <
class InputVector>
4099 const InputVector &fe_function,
4100 std::vector<solution_value_type<typename InputVector::value_type>> &values)
4103 std::vector<typename InputVector::value_type> local_dof_values(
4104 this->fe_interface->n_current_interface_dofs());
4105 this->get_local_dof_values(fe_function, local_dof_values);
4107 get_average_of_function_values_from_local_dof_values(local_dof_values,
4113 template <
int dim,
int spacedim>
4114 template <
class InputVector>
4118 const InputVector &local_dof_values,
4119 std::vector<solution_gradient_type<typename InputVector::value_type>>
4124 for (
const auto dof_index : this->fe_interface->dof_indices())
4125 for (
const auto q_index : this->fe_interface->quadrature_point_indices())
4130 gradients[q_index] += local_dof_values[dof_index] *
4131 average_of_gradients(dof_index, q_index);
4137 template <
int dim,
int spacedim>
4138 template <
class InputVector>
4141 const InputVector &fe_function,
4142 std::vector<solution_gradient_type<typename InputVector::value_type>>
4145 std::vector<typename InputVector::value_type> local_dof_values(
4146 this->fe_interface->n_current_interface_dofs());
4147 this->get_local_dof_values(fe_function, local_dof_values);
4149 get_average_of_function_gradients_from_local_dof_values(local_dof_values,
4155 template <
int dim,
int spacedim>
4156 template <
class InputVector>
4159 const InputVector &local_dof_values,
4160 std::vector<solution_hessian_type<typename InputVector::value_type>>
4165 for (
const auto dof_index : this->fe_interface->dof_indices())
4166 for (
const auto q_index : this->fe_interface->quadrature_point_indices())
4172 local_dof_values[dof_index] * jump_in_hessians(dof_index, q_index);
4178 template <
int dim,
int spacedim>
4179 template <
class InputVector>
4182 const InputVector &fe_function,
4183 std::vector<solution_hessian_type<typename InputVector::value_type>>
4186 std::vector<typename InputVector::value_type> local_dof_values(
4187 this->fe_interface->n_current_interface_dofs());
4188 this->get_local_dof_values(fe_function, local_dof_values);
4190 get_jump_in_function_hessians_from_local_dof_values(local_dof_values,
4196 template <
int dim,
int spacedim>
4197 template <
class InputVector>
4200 const InputVector &local_dof_values,
4201 std::vector<solution_hessian_type<typename InputVector::value_type>>
4206 for (
unsigned int dof_index : this->fe_interface->dof_indices())
4207 for (
const auto q_index : this->fe_interface->quadrature_point_indices())
4212 hessians[q_index] += local_dof_values[dof_index] *
4213 average_of_hessians(dof_index, q_index);
4219 template <
int dim,
int spacedim>
4220 template <
class InputVector>
4223 const InputVector &fe_function,
4224 std::vector<solution_hessian_type<typename InputVector::value_type>>
4227 std::vector<typename InputVector::value_type> local_dof_values(
4228 this->fe_interface->n_current_interface_dofs());
4229 this->get_local_dof_values(fe_function, local_dof_values);
4231 get_average_of_function_hessians_from_local_dof_values(local_dof_values,
4237 template <
int dim,
int spacedim>
4238 template <
class InputVector>
4242 const InputVector &local_dof_values,
4244 solution_third_derivative_type<typename InputVector::value_type>>
4245 &third_derivatives)
const
4248 this->fe_interface->n_quadrature_points);
4250 for (
unsigned int dof_index : this->fe_interface->dof_indices())
4251 for (
const auto q_index : this->fe_interface->quadrature_point_indices())
4254 third_derivatives[q_index] = 0.;
4256 third_derivatives[q_index] +=
4257 local_dof_values[dof_index] *
4258 jump_in_third_derivatives(dof_index, q_index);
4264 template <
int dim,
int spacedim>
4265 template <
class InputVector>
4268 const InputVector &fe_function,
4270 solution_third_derivative_type<typename InputVector::value_type>>
4271 &third_derivatives)
const
4273 std::vector<typename InputVector::value_type> local_dof_values(
4274 this->fe_interface->n_current_interface_dofs());
4275 this->get_local_dof_values(fe_function, local_dof_values);
4277 get_jump_in_function_third_derivatives_from_local_dof_values(
4278 local_dof_values, third_derivatives);
UpdateFlags get_update_flags() 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
Tensor< 2, spacedim > average_hessian(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) 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
const Triangulation< dim, spacedim >::cell_iterator get_cell(const unsigned int cell_index) const
FESubfaceValues< dim, spacedim > internal_fe_subface_values_neighbor
FEInterfaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature, const UpdateFlags update_flags)
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
const unsigned int n_quadrature_points
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
double average(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
double shape_value(const bool here_or_there, const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
const FEFaceValuesBase< dim, spacedim > & get_fe_face_values(const unsigned int cell_index) const
Tensor< 2, spacedim > average_of_shape_hessians(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
double jump(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
Tensor< 3, spacedim > jump_3rd_derivative(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
FEInterfaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const hp::QCollection< dim - 1 > &quadrature, const UpdateFlags update_flags)
FEFaceValuesBase< dim, spacedim > * fe_face_values
std::vector< types::global_dof_index > interface_dof_indices
Tensor< 1, spacedim > jump_gradient(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
std::vector< std::array< unsigned int, 2 > > dofmap
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
FESubfaceValues< dim, spacedim > internal_fe_subface_values
Tensor< 2, spacedim > jump_hessian(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) 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< 1, spacedim > average_gradient(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
void reinit(const CellIteratorType &cell, const unsigned int face_no)
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)
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_cxx20::ranges::iota_view< unsigned int, unsigned int > quadrature_point_indices() const
double jump_in_shape_values(const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
const FEInterfaceViews::Scalar< dim, spacedim > operator[](const FEValuesExtractors::Scalar &scalar) 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
FEFaceValues< dim, spacedim > internal_fe_face_values
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
FEFaceValues< dim, spacedim > internal_fe_face_values_neighbor
const Quadrature< dim - 1 > & get_quadrature() const
const FEInterfaceViews::Vector< dim, spacedim > operator[](const FEValuesExtractors::Vector &vector) 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 average_gradient(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 jump_hessian(const unsigned int interface_dof_index, const unsigned int q_point) const
gradient_type jump_gradient(const unsigned int interface_dof_index, const unsigned int q_point) const
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
value_type jump(const unsigned int interface_dof_index, const unsigned int q_point) 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
value_type average(const unsigned int interface_dof_index, const unsigned int q_point) 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
value_type average_value(const unsigned int interface_dof_index, const unsigned int q_point) const
void get_average_of_function_values(const InputVector &fe_function, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
hessian_type average_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_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
third_derivative_type jump_3rd_derivative(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
value_type average(const unsigned int interface_dof_index, const unsigned int q_point) const
Vector(const FEInterfaceValues< dim, spacedim > &fe_interface, const unsigned int first_vector_component)
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
third_derivative_type jump_3rd_derivative(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
value_type jump_in_values(const unsigned int interface_dof_index, const unsigned int q_point) const
value_type jump(const unsigned int interface_dof_index, const unsigned int q_point) const
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
gradient_type average_gradient(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 FEValuesViews::Vector< dim, spacedim >::hessian_type hessian_type
unsigned int n_components() const
Abstract base class for mapping classes.
#define DEAL_II_DEPRECATED
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
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)
static const unsigned int invalid_unsigned_int
boost::integer_range< IncrementableType > iota_view
typename internal::ProductTypeImpl< typename std::decay< T >::type, typename std::decay< U >::type >::type type