16 #ifndef dealii_fe_interface_values_h
17 #define dealii_fe_interface_values_h
52 template <
int dim,
int spacedim = dim>
113 template <
class CellIteratorType>
115 reinit(
const CellIteratorType & cell,
116 const unsigned int face_no,
117 const unsigned int sub_face_no,
119 const unsigned int face_no_neighbor,
120 const unsigned int sub_face_no_neighbor);
133 template <
class CellIteratorType>
135 reinit(
const CellIteratorType &cell,
const unsigned int face_no);
185 JxW(
const unsigned int quadrature_point)
const;
192 const std::vector<double> &
204 const std::vector<Tensor<1, spacedim>> &
212 const std::vector<Point<spacedim>> &
239 std::vector<types::global_dof_index>
253 std::array<unsigned int, 2>
265 normal(
const unsigned int q_point_index)
const;
299 const unsigned int interface_dof_index,
300 const unsigned int q_point,
301 const unsigned int component = 0)
const;
318 jump(
const unsigned int interface_dof_index,
319 const unsigned int q_point,
320 const unsigned int component = 0)
const;
332 average(
const unsigned int interface_dof_index,
333 const unsigned int q_point,
334 const unsigned int component = 0)
const;
347 const unsigned int q_point,
348 const unsigned int component = 0)
const;
362 const unsigned int q_point,
363 const unsigned int component = 0)
const;
376 const unsigned int q_point,
377 const unsigned int component = 0)
const;
391 const unsigned int q_point,
392 const unsigned int component = 0)
const;
405 const unsigned int q_point,
406 const unsigned int component = 0)
const;
423 std::vector<std::array<unsigned int, 2>>
dofmap;
465 template <
int dim,
int spacedim>
471 : n_quadrature_points(quadrature.size())
472 , internal_fe_face_values(mapping, fe, quadrature, update_flags)
473 , internal_fe_subface_values(mapping, fe, quadrature, update_flags)
474 , internal_fe_face_values_neighbor(mapping, fe, quadrature, update_flags)
475 , internal_fe_subface_values_neighbor(mapping, fe, quadrature, update_flags)
476 , fe_face_values(nullptr)
477 , fe_face_values_neighbor(nullptr)
482 template <
int dim,
int spacedim>
487 : n_quadrature_points(quadrature.size())
496 , internal_fe_face_values_neighbor(
StaticMappingQ1<dim, spacedim>::mapping,
500 , internal_fe_subface_values_neighbor(
StaticMappingQ1<dim, spacedim>::mapping,
504 , fe_face_values(nullptr)
505 , fe_face_values_neighbor(nullptr)
510 template <
int dim,
int spacedim>
511 template <
class CellIteratorType>
514 const CellIteratorType & cell,
515 const unsigned int face_no,
516 const unsigned int sub_face_no,
518 const unsigned int face_no_neighbor,
519 const unsigned int sub_face_no_neighbor)
523 internal_fe_face_values.reinit(cell, face_no);
524 fe_face_values = &internal_fe_face_values;
528 internal_fe_subface_values.reinit(cell, face_no, sub_face_no);
529 fe_face_values = &internal_fe_subface_values;
533 internal_fe_face_values_neighbor.reinit(cell_neighbor, face_no_neighbor);
534 fe_face_values_neighbor = &internal_fe_face_values_neighbor;
538 internal_fe_subface_values_neighbor.reinit(cell_neighbor,
540 sub_face_no_neighbor);
541 fe_face_values_neighbor = &internal_fe_subface_values_neighbor;
547 std::vector<types::global_dof_index> v(
548 fe_face_values->get_fe().n_dofs_per_cell());
549 cell->get_active_or_mg_dof_indices(v);
550 std::vector<types::global_dof_index> v2(
551 fe_face_values_neighbor->get_fe().n_dofs_per_cell());
552 cell_neighbor->get_active_or_mg_dof_indices(v2);
558 std::map<types::global_dof_index, std::pair<unsigned int, unsigned int>>
560 std::pair<unsigned int, unsigned int> invalid_entry(
563 for (
unsigned int i = 0; i < v.size(); ++i)
566 auto result = tempmap.insert(std::make_pair(v[i], invalid_entry));
567 result.first->second.first = i;
570 for (
unsigned int i = 0; i < v2.size(); ++i)
573 auto result = tempmap.insert(std::make_pair(v2[i], invalid_entry));
574 result.first->second.second = i;
578 dofmap.resize(tempmap.size());
579 interface_dof_indices.resize(tempmap.size());
580 unsigned int idx = 0;
581 for (
auto &x : tempmap)
583 interface_dof_indices[idx] = x.first;
584 dofmap[idx] = {{x.second.first, x.second.second}};
592 template <
int dim,
int spacedim>
593 template <
class CellIteratorType>
596 const unsigned int face_no)
598 internal_fe_face_values.reinit(cell, face_no);
599 fe_face_values = &internal_fe_face_values;
600 fe_face_values_neighbor =
nullptr;
602 interface_dof_indices.resize(fe_face_values->get_fe().n_dofs_per_cell());
603 cell->get_active_or_mg_dof_indices(interface_dof_indices);
606 dofmap.resize(interface_dof_indices.size());
608 for (
unsigned int i = 0; i < interface_dof_indices.size(); ++i)
616 template <
int dim,
int spacedim>
620 Assert(fe_face_values !=
nullptr,
621 ExcMessage(
"This call requires a call to reinit() first."));
622 return fe_face_values->JxW(q);
627 template <
int dim,
int spacedim>
628 const std::vector<double> &
631 Assert(fe_face_values !=
nullptr,
632 ExcMessage(
"This call requires a call to reinit() first."));
633 return fe_face_values->get_JxW_values();
638 template <
int dim,
int spacedim>
639 const std::vector<Tensor<1, spacedim>> &
642 Assert(fe_face_values !=
nullptr,
643 ExcMessage(
"This call requires a call to reinit() first."));
644 return fe_face_values->get_normal_vectors();
649 template <
int dim,
int spacedim>
653 return internal_fe_face_values.get_quadrature();
658 template <
int dim,
int spacedim>
659 const std::vector<Point<spacedim>> &
662 Assert(fe_face_values !=
nullptr,
663 ExcMessage(
"This call requires a call to reinit() first."));
664 return fe_face_values->get_quadrature_points();
669 template <
int dim,
int spacedim>
673 return internal_fe_face_values.get_update_flags();
678 template <
int dim,
int spacedim>
683 interface_dof_indices.size() > 0,
685 "n_current_interface_dofs() is only available after a call to reinit()."));
686 return interface_dof_indices.size();
691 template <
int dim,
int spacedim>
695 return fe_face_values_neighbor ==
nullptr;
700 template <
int dim,
int spacedim>
701 std::vector<types::global_dof_index>
704 return interface_dof_indices;
709 template <
int dim,
int spacedim>
710 std::array<unsigned int, 2>
712 const unsigned int interface_dof_index)
const
715 return dofmap[interface_dof_index];
720 template <
int dim,
int spacedim>
723 const unsigned int cell_index)
const
727 cell_index == 0 || !at_boundary(),
729 "You are on a boundary, so you can only ask for the first FEFaceValues object."));
731 return (cell_index == 0) ? *fe_face_values : *fe_face_values_neighbor;
736 template <
int dim,
int spacedim>
740 return fe_face_values->normal_vector(q_point_index);
745 template <
int dim,
int spacedim>
748 const bool here_or_there,
749 const unsigned int interface_dof_index,
750 const unsigned int q_point,
751 const unsigned int component)
const
753 const auto dof_pair = dofmap[interface_dof_index];
756 return get_fe_face_values(0).shape_value_component(dof_pair[0],
760 return get_fe_face_values(1).shape_value_component(dof_pair[1],
769 template <
int dim,
int spacedim>
772 const unsigned int q_point,
773 const unsigned int component)
const
775 const auto dof_pair = dofmap[interface_dof_index];
780 value += get_fe_face_values(0).shape_value_component(dof_pair[0],
784 value -= get_fe_face_values(1).shape_value_component(dof_pair[1],
792 template <
int dim,
int spacedim>
795 const unsigned int interface_dof_index,
796 const unsigned int q_point,
797 const unsigned int component)
const
799 const auto dof_pair = dofmap[interface_dof_index];
802 return 1.0 * get_fe_face_values(0).shape_value_component(dof_pair[0],
809 value += 0.5 * get_fe_face_values(0).shape_value_component(dof_pair[0],
813 value += 0.5 * get_fe_face_values(1).shape_value_component(dof_pair[1],
822 template <
int dim,
int spacedim>
825 const unsigned int interface_dof_index,
826 const unsigned int q_point,
827 const unsigned int component)
const
829 const auto dof_pair = dofmap[interface_dof_index];
832 return get_fe_face_values(0).shape_grad_component(dof_pair[0],
839 value += 0.5 * get_fe_face_values(0).shape_grad_component(dof_pair[0],
843 value += 0.5 * get_fe_face_values(1).shape_grad_component(dof_pair[1],
852 template <
int dim,
int spacedim>
855 const unsigned int interface_dof_index,
856 const unsigned int q_point,
857 const unsigned int component)
const
859 const auto dof_pair = dofmap[interface_dof_index];
862 return get_fe_face_values(0).shape_hessian_component(dof_pair[0],
869 value += 0.5 * get_fe_face_values(0).shape_hessian_component(dof_pair[0],
873 value += 0.5 * get_fe_face_values(1).shape_hessian_component(dof_pair[1],
882 template <
int dim,
int spacedim>
885 const unsigned int interface_dof_index,
886 const unsigned int q_point,
887 const unsigned int component)
const
889 const auto dof_pair = dofmap[interface_dof_index];
892 return get_fe_face_values(0).shape_grad_component(dof_pair[0],
899 value += 1.0 * get_fe_face_values(0).shape_grad_component(dof_pair[0],
903 value += -1.0 * get_fe_face_values(1).shape_grad_component(dof_pair[1],
912 template <
int dim,
int spacedim>
915 const unsigned int interface_dof_index,
916 const unsigned int q_point,
917 const unsigned int component)
const
919 const auto dof_pair = dofmap[interface_dof_index];
922 return get_fe_face_values(0).shape_hessian_component(dof_pair[0],
929 value += 1.0 * get_fe_face_values(0).shape_hessian_component(dof_pair[0],
933 value += -1.0 * get_fe_face_values(1).shape_hessian_component(dof_pair[1],
941 template <
int dim,
int spacedim>
944 const unsigned int interface_dof_index,
945 const unsigned int q_point,
946 const unsigned int component)
const
948 const auto dof_pair = dofmap[interface_dof_index];
951 return get_fe_face_values(0).shape_3rd_derivative_component(dof_pair[0],
959 1.0 * get_fe_face_values(0).shape_3rd_derivative_component(dof_pair[0],
964 -1.0 * get_fe_face_values(1).shape_3rd_derivative_component(dof_pair[1],