23 template <
int dim,
int spacedim>
28 const UpdateFlags &update_flags,
46 template <
int dim,
int spacedim>
51 const UpdateFlags &update_flags,
52 const UpdateFlags &neighbor_update_flags,
71 template <
int dim,
int spacedim>
75 const UpdateFlags &update_flags,
89 template <
int dim,
int spacedim>
93 const UpdateFlags &update_flags,
94 const UpdateFlags &neighbor_update_flags,
103 neighbor_update_flags,
111 template <
int dim,
int spacedim>
135 template <
int dim,
int spacedim>
161 template <
int dim,
int spacedim>
178 template <
int dim,
int spacedim>
199 template <
int dim,
int spacedim>
223 template <
int dim,
int spacedim>
231 fe_values = std::make_unique<FEValues<dim, spacedim>>(
243 hp_fe_values = std::make_unique<hp::FEValues<dim, spacedim>>(
265 template <
int dim,
int spacedim>
269 const unsigned int face_no)
285 return reinit(cell, cell, face_no);
291 template <
int dim,
int spacedim>
297 const unsigned int face_no)
308 if (neighbor_cell == cell)
319 unsigned int dominated_fe_index =
fe_collection->find_dominated_fe(
320 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
336 const auto &
fe = (*fe_collection)[cell->active_fe_index()];
347 template <
int dim,
int spacedim>
351 const unsigned int face_no,
352 const unsigned int subface_no)
360 std::make_unique<FESubfaceValues<dim, spacedim>>(
371 return reinit(cell, cell, face_no, subface_no);
375 return reinit(cell, face_no);
380 template <
int dim,
int spacedim>
386 const unsigned int face_no,
387 const unsigned int subface_no)
395 std::make_unique<hp::FESubfaceValues<dim, spacedim>>(
401 if (neighbor_cell == cell)
412 unsigned int dominated_fe_index =
fe_collection->find_dominated_fe(
413 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
431 const auto &
fe = (*fe_collection)[cell->active_fe_index()];
441 return reinit(cell, neighbor_cell, face_no);
447 template <
int dim,
int spacedim>
451 const unsigned int face_no,
454 const unsigned int face_no_neighbor)
466 template <
int dim,
int spacedim>
470 const unsigned int face_no,
471 const unsigned int sub_face_no,
474 const unsigned int face_no_neighbor,
475 const unsigned int sub_face_no_neighbor)
481 std::make_unique<FEInterfaceValues<dim, spacedim>>(
489 sub_face_no_neighbor);
495 std::make_unique<FEInterfaceValues<dim, spacedim>>(
506 unsigned int dominated_fe_index =
fe_collection->find_dominated_fe(
507 {cell->active_fe_index(), cell_neighbor->active_fe_index()});
521 sub_face_no_neighbor,
535 template <
int dim,
int spacedim>
577 template <
int dim,
int spacedim>
581 const unsigned int face_no)
587 std::make_unique<FEFaceValues<dim, spacedim>>(
602 template <
int dim,
int spacedim>
608 const unsigned int face_no)
614 std::make_unique<hp::FEFaceValues<dim, spacedim>>(
620 if (neighbor_cell == cell)
631 unsigned int dominated_fe_index =
fe_collection->find_dominated_fe(
632 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
649 const auto &neighbor_fe =
650 (*fe_collection)[neighbor_cell->active_fe_index()];
661 template <
int dim,
int spacedim>
665 const unsigned int face_no,
666 const unsigned int subface_no)
674 std::make_unique<FESubfaceValues<dim, spacedim>>(
692 template <
int dim,
int spacedim>
698 const unsigned int face_no,
699 const unsigned int subface_no)
707 std::make_unique<hp::FESubfaceValues<dim, spacedim>>(
713 if (neighbor_cell == cell)
726 unsigned int dominated_fe_index =
fe_collection->find_dominated_fe(
727 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
745 const auto &neighbor_fe =
746 (*fe_collection)[neighbor_cell->active_fe_index()];
762 template <
int dim,
int spacedim>
767 ExcMessage(
"You have to initialize the cache using one of the "
768 "reinit functions first!"));
774 template <
int dim,
int spacedim>
779 ExcMessage(
"You have to initialize the cache using one of the "
780 "reinit functions first!"));
786 template <
int dim,
int spacedim>
791 ExcMessage(
"You have to initialize the cache using one of the "
792 "reinit functions first!"));
798 template <
int dim,
int spacedim>
799 const std::vector<Point<spacedim>> &
807 template <
int dim,
int spacedim>
808 const std::vector<double> &
816 template <
int dim,
int spacedim>
817 const std::vector<double> &
825 template <
int dim,
int spacedim>
826 const std::vector<Tensor<1, spacedim>> &
834 template <
int dim,
int spacedim>
835 const std::vector<Tensor<1, spacedim>> &
843 template <
int dim,
int spacedim>
844 const std::vector<types::global_dof_index> &
852 template <
int dim,
int spacedim>
861 template <
int dim,
int spacedim>
862 const std::vector<types::global_dof_index> &
870 template <
int dim,
int spacedim>
879 template <
int dim,
int spacedim>
888 template <
int dim,
int spacedim>
897 template <
int dim,
int spacedim>
908 template <
int dim,
int spacedim>
919 template <
int dim,
int spacedim>
929 template <
int dim,
int spacedim>
939 template <
int dim,
int spacedim>
950 template <
int dim,
int spacedim>
961 template <
int dim,
int spacedim>
971 template <
int dim,
int spacedim>
981 template <
int dim,
int spacedim>
990 template <
int dim,
int spacedim>
999 template <
int dim,
int spacedim>
1008 template <
int dim,
int spacedim>
1017 template <
int dim,
int spacedim>
1031#include "scratch_data.inst"
Abstract base class for mapping classes.
GeneralDataStorage user_data_storage
hp::QCollection< dim - 1 > face_quadrature_collection
bool has_hp_capabilities() const
const hp::MappingCollection< dim, spacedim > & get_mapping_collection() const
const FEValuesBase< dim, spacedim > & get_current_fe_values() const
GeneralDataStorage & get_general_data_storage()
const Quadrature< dim - 1 > & get_face_quadrature() const
const std::vector< Tensor< 1, spacedim > > & get_neighbor_normal_vectors()
unsigned int n_dofs_per_cell() const
std::unique_ptr< hp::FESubfaceValues< dim, spacedim > > hp_fe_subface_values
unsigned int n_neighbor_dofs_per_cell() const
std::unique_ptr< FEValues< dim, spacedim > > neighbor_fe_values
std::unique_ptr< FEFaceValues< dim, spacedim > > fe_face_values
std::unique_ptr< FEValues< dim, spacedim > > fe_values
UpdateFlags cell_update_flags
const std::vector< double > & get_JxW_values() const
const Mapping< dim, spacedim > & get_mapping() const
std::unique_ptr< FEInterfaceValues< dim, spacedim > > interface_fe_values
const std::vector< double > & get_neighbor_JxW_values() const
Quadrature< dim > cell_quadrature
const FiniteElement< dim, spacedim > & get_fe() const
SmartPointer< const hp::FECollection< dim, spacedim > > fe_collection
std::unique_ptr< FEFaceValues< dim, spacedim > > neighbor_fe_face_values
const hp::QCollection< dim > & get_cell_quadrature_collection() const
std::vector< types::global_dof_index > neighbor_dof_indices
SmartPointer< const FiniteElement< dim, spacedim > > fe
UpdateFlags face_update_flags
const std::vector< types::global_dof_index > & get_neighbor_dof_indices() const
UpdateFlags get_face_update_flags() const
std::unique_ptr< hp::FEValues< dim, spacedim > > hp_fe_values
SmartPointer< const hp::MappingCollection< dim, spacedim > > mapping_collection
UpdateFlags neighbor_cell_update_flags
bool hp_capability_enabled
UpdateFlags get_neighbor_cell_update_flags() const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors() const
const FEValues< dim, spacedim > & reinit_neighbor(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
std::unique_ptr< hp::FEFaceValues< dim, spacedim > > neighbor_hp_fe_face_values
std::vector< types::global_dof_index > local_dof_indices
SmartPointer< const Mapping< dim, spacedim > > mapping
SmartPointer< const FEValuesBase< dim, spacedim > > current_fe_values
const std::vector< Point< spacedim > > & get_quadrature_points() const
std::unique_ptr< hp::FEValues< dim, spacedim > > neighbor_hp_fe_values
const FEValuesBase< dim, spacedim > & get_current_neighbor_fe_values() const
const std::vector< types::global_dof_index > & get_local_dof_indices() const
SmartPointer< const FEValuesBase< dim, spacedim > > current_neighbor_fe_values
Quadrature< dim - 1 > face_quadrature
const FEInterfaceValues< dim, spacedim > & get_current_interface_fe_values() const
hp::QCollection< dim > cell_quadrature_collection
UpdateFlags get_cell_update_flags() const
UpdateFlags neighbor_face_update_flags
const Quadrature< dim > & get_cell_quadrature() const
const hp::QCollection< dim - 1 > & get_face_quadrature_collection() const
ScratchData(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags &update_flags, const Quadrature< dim - 1 > &face_quadrature=Quadrature< dim - 1 >(), const UpdateFlags &face_update_flags=update_default)
std::unique_ptr< FESubfaceValues< dim, spacedim > > neighbor_fe_subface_values
const FEValues< dim, spacedim > & reinit(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
std::unique_ptr< hp::FESubfaceValues< dim, spacedim > > neighbor_hp_fe_subface_values
std::unique_ptr< FESubfaceValues< dim, spacedim > > fe_subface_values
GeneralDataStorage internal_data_storage
UpdateFlags get_neighbor_face_update_flags() const
std::unique_ptr< hp::FEFaceValues< dim, spacedim > > hp_fe_face_values
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcOnlyAvailableWithoutHP()
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcOnlyAvailableWithHP()
typename ActiveSelector::active_cell_iterator active_cell_iterator
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
const types::fe_index invalid_fe_index
static const unsigned int invalid_unsigned_int