23 template <
int dim,
int spacedim>
33 , cell_quadrature(quadrature)
34 , face_quadrature(face_quadrature)
35 , hp_capability_enabled(false)
36 , cell_update_flags(update_flags)
37 , neighbor_cell_update_flags(update_flags)
38 , face_update_flags(face_update_flags)
39 , neighbor_face_update_flags(face_update_flags)
40 , local_dof_indices(fe.n_dofs_per_cell())
41 , neighbor_dof_indices(fe.n_dofs_per_cell())
46 template <
int dim,
int spacedim>
58 , cell_quadrature(quadrature)
59 , face_quadrature(face_quadrature)
60 , hp_capability_enabled(false)
61 , cell_update_flags(update_flags)
62 , neighbor_cell_update_flags(neighbor_update_flags)
63 , face_update_flags(face_update_flags)
64 , neighbor_face_update_flags(neighbor_face_update_flags)
65 , local_dof_indices(fe.n_dofs_per_cell())
66 , neighbor_dof_indices(fe.n_dofs_per_cell())
71 template <
int dim,
int spacedim>
89 template <
int dim,
int spacedim>
103 neighbor_update_flags,
106 neighbor_face_update_flags)
111 template <
int dim,
int spacedim>
119 : mapping_collection(&mapping_collection)
120 , fe_collection(&fe_collection)
121 , cell_quadrature_collection(cell_quadrature_collection)
122 , face_quadrature_collection(face_quadrature_collection)
123 , hp_capability_enabled(true)
124 , cell_update_flags(cell_update_flags)
125 , neighbor_cell_update_flags(cell_update_flags)
126 , face_update_flags(face_update_flags)
127 , neighbor_face_update_flags(face_update_flags)
135 template <
int dim,
int spacedim>
145 : mapping_collection(&mapping_collection)
146 , fe_collection(&fe_collection)
147 , cell_quadrature_collection(cell_quadrature_collection)
148 , face_quadrature_collection(face_quadrature_collection)
149 , hp_capability_enabled(true)
150 , cell_update_flags(cell_update_flags)
151 , neighbor_cell_update_flags(neighbor_cell_update_flags)
152 , face_update_flags(face_update_flags)
153 , neighbor_face_update_flags(neighbor_face_update_flags)
161 template <
int dim,
int spacedim>
168 :
ScratchData(fe_collection.get_reference_cell_default_linear_mapping(),
170 cell_quadrature_collection,
172 face_quadrature_collection,
178 template <
int dim,
int spacedim>
187 :
ScratchData(fe_collection.get_reference_cell_default_linear_mapping(),
189 cell_quadrature_collection,
191 neighbor_cell_update_flags,
192 face_quadrature_collection,
194 neighbor_face_update_flags)
199 template <
int dim,
int spacedim>
202 : mapping(scratch.mapping)
204 , cell_quadrature(scratch.cell_quadrature)
205 , face_quadrature(scratch.face_quadrature)
206 , mapping_collection(scratch.mapping_collection)
207 , fe_collection(scratch.fe_collection)
208 , cell_quadrature_collection(scratch.cell_quadrature_collection)
209 , face_quadrature_collection(scratch.face_quadrature_collection)
210 , hp_capability_enabled(scratch.hp_capability_enabled)
211 , cell_update_flags(scratch.cell_update_flags)
212 , neighbor_cell_update_flags(scratch.neighbor_cell_update_flags)
213 , face_update_flags(scratch.face_update_flags)
214 , neighbor_face_update_flags(scratch.neighbor_face_update_flags)
215 , local_dof_indices(scratch.local_dof_indices)
216 , neighbor_dof_indices(scratch.neighbor_dof_indices)
217 , user_data_storage(scratch.user_data_storage)
218 , internal_data_storage(scratch.internal_data_storage)
223 template <
int dim,
int spacedim>
228 if (hp_capability_enabled ==
false)
231 fe_values = std::make_unique<FEValues<dim, spacedim>>(
232 *mapping, *fe, cell_quadrature, cell_update_flags);
235 local_dof_indices.resize(fe_values->dofs_per_cell);
236 cell->get_dof_indices(local_dof_indices);
237 current_fe_values = fe_values.get();
243 hp_fe_values = std::make_unique<hp::FEValues<dim, spacedim>>(
246 cell_quadrature_collection,
249 hp_fe_values->
reinit(cell);
250 const auto &fe_values = hp_fe_values->get_present_fe_values();
253 (*fe_collection)[cell->active_fe_index()].n_dofs_per_cell(),
254 fe_values.dofs_per_cell);
255 local_dof_indices.resize(fe_values.dofs_per_cell);
256 cell->get_dof_indices(local_dof_indices);
258 current_fe_values = &fe_values;
265 template <
int dim,
int spacedim>
269 const unsigned int face_no)
271 if (hp_capability_enabled ==
false)
274 fe_face_values = std::make_unique<FEFaceValues<dim, spacedim>>(
275 *mapping, *fe, face_quadrature, face_update_flags);
277 fe_face_values->
reinit(cell, face_no);
278 local_dof_indices.resize(fe->n_dofs_per_cell());
279 cell->get_dof_indices(local_dof_indices);
280 current_fe_values = fe_face_values.get();
281 return *fe_face_values;
285 return reinit(cell, cell, face_no);
291 template <
int dim,
int spacedim>
297 const unsigned int face_no)
299 Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
301 if (!hp_fe_face_values)
302 hp_fe_face_values = std::make_unique<hp::FEFaceValues<dim, spacedim>>(
305 face_quadrature_collection,
308 if (neighbor_cell == cell)
310 hp_fe_face_values->reinit(cell, face_no);
319 unsigned int dominated_fe_index = fe_collection->find_dominated_fe(
320 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
329 hp_fe_face_values->reinit(cell,
335 const auto &fe_face_values = hp_fe_face_values->get_present_fe_values();
336 const auto &fe = (*fe_collection)[cell->active_fe_index()];
338 local_dof_indices.resize(fe.n_dofs_per_cell());
339 cell->get_dof_indices(local_dof_indices);
341 current_fe_values = &fe_face_values;
342 return fe_face_values;
347 template <
int dim,
int spacedim>
351 const unsigned int face_no,
352 const unsigned int subface_no)
356 if (hp_capability_enabled ==
false)
358 if (!fe_subface_values)
360 std::make_unique<FESubfaceValues<dim, spacedim>>(
361 *mapping, *fe, face_quadrature, face_update_flags);
362 fe_subface_values->reinit(cell, face_no, subface_no);
363 local_dof_indices.resize(fe->n_dofs_per_cell());
364 cell->get_dof_indices(local_dof_indices);
366 current_fe_values = fe_subface_values.get();
367 return *fe_subface_values;
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)
389 Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
393 if (!hp_fe_subface_values)
394 hp_fe_subface_values =
395 std::make_unique<hp::FESubfaceValues<dim, spacedim>>(
398 face_quadrature_collection,
401 if (neighbor_cell == cell)
403 hp_fe_subface_values->reinit(cell, face_no, subface_no);
412 unsigned int dominated_fe_index = fe_collection->find_dominated_fe(
413 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
422 hp_fe_subface_values->reinit(cell,
429 const auto &fe_subface_values =
430 hp_fe_subface_values->get_present_fe_values();
431 const auto &fe = (*fe_collection)[cell->active_fe_index()];
433 local_dof_indices.resize(fe.n_dofs_per_cell());
434 cell->get_dof_indices(local_dof_indices);
436 current_fe_values = &fe_subface_values;
437 return fe_subface_values;
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)
477 if (hp_capability_enabled ==
false)
479 if (!interface_fe_values)
480 interface_fe_values =
481 std::make_unique<FEInterfaceValues<dim, spacedim>>(
482 *mapping, *fe, face_quadrature, face_update_flags);
484 interface_fe_values->
reinit(cell,
489 sub_face_no_neighbor);
493 if (!interface_fe_values)
494 interface_fe_values =
495 std::make_unique<FEInterfaceValues<dim, spacedim>>(
498 face_quadrature_collection,
506 unsigned int dominated_fe_index = fe_collection->find_dominated_fe(
507 {cell->active_fe_index(), cell_neighbor->active_fe_index()});
516 interface_fe_values->reinit(cell,
521 sub_face_no_neighbor,
526 current_fe_values = &interface_fe_values->get_fe_face_values(0);
527 current_neighbor_fe_values = &interface_fe_values->get_fe_face_values(1);
529 local_dof_indices = interface_fe_values->get_interface_dof_indices();
530 return *interface_fe_values;
535 template <
int dim,
int spacedim>
540 if (hp_capability_enabled ==
false)
542 if (!neighbor_fe_values)
543 neighbor_fe_values = std::make_unique<FEValues<dim, spacedim>>(
544 *mapping, *fe, cell_quadrature, neighbor_cell_update_flags);
546 neighbor_fe_values->
reinit(cell);
547 cell->get_dof_indices(neighbor_dof_indices);
548 current_neighbor_fe_values = neighbor_fe_values.get();
549 return *neighbor_fe_values;
553 if (!neighbor_hp_fe_values)
554 neighbor_hp_fe_values = std::make_unique<hp::FEValues<dim, spacedim>>(
557 cell_quadrature_collection,
558 neighbor_cell_update_flags);
560 neighbor_hp_fe_values->
reinit(cell);
561 const auto &neighbor_fe_values =
562 neighbor_hp_fe_values->get_present_fe_values();
565 (*fe_collection)[cell->active_fe_index()].n_dofs_per_cell(),
566 neighbor_fe_values.dofs_per_cell);
567 neighbor_dof_indices.resize(neighbor_fe_values.dofs_per_cell);
568 cell->get_dof_indices(neighbor_dof_indices);
570 current_neighbor_fe_values = &neighbor_fe_values;
571 return neighbor_fe_values;
577 template <
int dim,
int spacedim>
581 const unsigned int face_no)
583 if (hp_capability_enabled ==
false)
585 if (!neighbor_fe_face_values)
586 neighbor_fe_face_values =
587 std::make_unique<FEFaceValues<dim, spacedim>>(
588 *mapping, *fe, face_quadrature, neighbor_face_update_flags);
589 neighbor_fe_face_values->
reinit(cell, face_no);
590 cell->get_dof_indices(neighbor_dof_indices);
591 current_neighbor_fe_values = neighbor_fe_face_values.get();
592 return *neighbor_fe_face_values;
596 return reinit_neighbor(cell, cell, face_no);
602 template <
int dim,
int spacedim>
608 const unsigned int face_no)
610 Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
612 if (!neighbor_hp_fe_face_values)
613 neighbor_hp_fe_face_values =
614 std::make_unique<hp::FEFaceValues<dim, spacedim>>(
617 face_quadrature_collection,
618 neighbor_face_update_flags);
620 if (neighbor_cell == cell)
622 neighbor_hp_fe_face_values->reinit(neighbor_cell, face_no);
631 unsigned int dominated_fe_index = fe_collection->find_dominated_fe(
632 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
641 neighbor_hp_fe_face_values->reinit(neighbor_cell,
647 const auto &neighbor_fe_face_values =
648 neighbor_hp_fe_face_values->get_present_fe_values();
649 const auto &neighbor_fe =
650 (*fe_collection)[neighbor_cell->active_fe_index()];
652 neighbor_dof_indices.resize(neighbor_fe.n_dofs_per_cell());
653 neighbor_cell->get_dof_indices(neighbor_dof_indices);
655 current_neighbor_fe_values = &neighbor_fe_face_values;
656 return neighbor_fe_face_values;
661 template <
int dim,
int spacedim>
665 const unsigned int face_no,
666 const unsigned int subface_no)
670 if (hp_capability_enabled ==
false)
672 if (!neighbor_fe_subface_values)
673 neighbor_fe_subface_values =
674 std::make_unique<FESubfaceValues<dim, spacedim>>(
675 *mapping, *fe, face_quadrature, neighbor_face_update_flags);
676 neighbor_fe_subface_values->reinit(cell, face_no, subface_no);
677 cell->get_dof_indices(neighbor_dof_indices);
678 current_neighbor_fe_values = neighbor_fe_subface_values.get();
679 return *neighbor_fe_subface_values;
683 return reinit_neighbor(cell, cell, face_no, subface_no);
687 return reinit_neighbor(cell, face_no);
692 template <
int dim,
int spacedim>
698 const unsigned int face_no,
699 const unsigned int subface_no)
701 Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
705 if (!neighbor_hp_fe_subface_values)
706 neighbor_hp_fe_subface_values =
707 std::make_unique<hp::FESubfaceValues<dim, spacedim>>(
710 face_quadrature_collection,
711 neighbor_face_update_flags);
713 if (neighbor_cell == cell)
715 neighbor_hp_fe_subface_values->reinit(neighbor_cell,
726 unsigned int dominated_fe_index = fe_collection->find_dominated_fe(
727 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
736 neighbor_hp_fe_subface_values->reinit(neighbor_cell,
743 const auto &neighbor_fe_subface_values =
744 neighbor_hp_fe_subface_values->get_present_fe_values();
745 const auto &neighbor_fe =
746 (*fe_collection)[neighbor_cell->active_fe_index()];
748 neighbor_dof_indices.resize(neighbor_fe.n_dofs_per_cell());
749 neighbor_cell->get_dof_indices(neighbor_dof_indices);
751 current_neighbor_fe_values = &neighbor_fe_subface_values;
752 return neighbor_fe_subface_values;
756 return reinit_neighbor(cell, neighbor_cell, face_no);
762 template <
int dim,
int spacedim>
766 Assert(current_fe_values !=
nullptr,
767 ExcMessage(
"You have to initialize the cache using one of the "
768 "reinit functions first!"));
769 return *current_fe_values;
774 template <
int dim,
int spacedim>
778 Assert(interface_fe_values !=
nullptr,
779 ExcMessage(
"You have to initialize the cache using one of the "
780 "reinit functions first!"));
781 return *interface_fe_values;
786 template <
int dim,
int spacedim>
790 Assert(current_neighbor_fe_values !=
nullptr,
791 ExcMessage(
"You have to initialize the cache using one of the "
792 "reinit functions first!"));
793 return *current_neighbor_fe_values;
798 template <
int dim,
int spacedim>
799 const std::vector<Point<spacedim>> &
802 return get_current_fe_values().get_quadrature_points();
807 template <
int dim,
int spacedim>
808 const std::vector<double> &
811 return get_current_fe_values().get_JxW_values();
816 template <
int dim,
int spacedim>
817 const std::vector<double> &
820 return get_current_neighbor_fe_values().get_JxW_values();
825 template <
int dim,
int spacedim>
826 const std::vector<Tensor<1, spacedim>> &
829 return get_current_fe_values().get_normal_vectors();
834 template <
int dim,
int spacedim>
835 const std::vector<Tensor<1, spacedim>> &
838 return get_current_neighbor_fe_values().get_normal_vectors();
843 template <
int dim,
int spacedim>
844 const std::vector<types::global_dof_index> &
847 return local_dof_indices;
852 template <
int dim,
int spacedim>
856 return local_dof_indices.size();
861 template <
int dim,
int spacedim>
862 const std::vector<types::global_dof_index> &
865 return neighbor_dof_indices;
870 template <
int dim,
int spacedim>
874 return neighbor_dof_indices.size();
879 template <
int dim,
int spacedim>
883 return user_data_storage;
888 template <
int dim,
int spacedim>
892 return user_data_storage;
897 template <
int dim,
int spacedim>
901 Assert(hp_capability_enabled ==
false, ExcOnlyAvailableWithoutHP());
908 template <
int dim,
int spacedim>
912 Assert(hp_capability_enabled ==
false, ExcOnlyAvailableWithoutHP());
919 template <
int dim,
int spacedim>
923 Assert(hp_capability_enabled ==
false, ExcOnlyAvailableWithoutHP());
924 return cell_quadrature;
929 template <
int dim,
int spacedim>
933 Assert(hp_capability_enabled ==
false, ExcOnlyAvailableWithoutHP());
934 return face_quadrature;
939 template <
int dim,
int spacedim>
943 Assert(hp_capability_enabled ==
true, ExcOnlyAvailableWithHP());
945 return *mapping_collection;
950 template <
int dim,
int spacedim>
954 Assert(hp_capability_enabled ==
true, ExcOnlyAvailableWithHP());
956 return *fe_collection;
961 template <
int dim,
int spacedim>
965 Assert(hp_capability_enabled ==
true, ExcOnlyAvailableWithHP());
966 return cell_quadrature_collection;
971 template <
int dim,
int spacedim>
975 Assert(hp_capability_enabled ==
true, ExcOnlyAvailableWithHP());
976 return face_quadrature_collection;
981 template <
int dim,
int spacedim>
985 return hp_capability_enabled;
990 template <
int dim,
int spacedim>
994 return cell_update_flags;
999 template <
int dim,
int spacedim>
1003 return neighbor_cell_update_flags;
1008 template <
int dim,
int spacedim>
1012 return face_update_flags;
1017 template <
int dim,
int spacedim>
1021 return neighbor_face_update_flags;
1031#include "scratch_data.inst"
void reinit(const Triangulation< dim, spacedim > &tria)
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &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, 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)
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
Abstract base class for mapping classes.
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
unsigned int n_neighbor_dofs_per_cell() const
const std::vector< double > & get_JxW_values() const
const Mapping< dim, spacedim > & get_mapping() const
const std::vector< double > & get_neighbor_JxW_values() const
ObserverPointer< const hp::FECollection< dim, spacedim > > fe_collection
const FiniteElement< dim, spacedim > & get_fe() const
const hp::QCollection< dim > & get_cell_quadrature_collection() const
std::vector< types::global_dof_index > neighbor_dof_indices
const std::vector< types::global_dof_index > & get_neighbor_dof_indices() const
UpdateFlags get_face_update_flags() const
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::vector< types::global_dof_index > local_dof_indices
const std::vector< Point< spacedim > > & get_quadrature_points() const
const FEValuesBase< dim, spacedim > & get_current_neighbor_fe_values() const
const std::vector< types::global_dof_index > & get_local_dof_indices() const
const FEInterfaceValues< dim, spacedim > & get_current_interface_fe_values() const
UpdateFlags get_cell_update_flags() const
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)
const FEValues< dim, spacedim > & reinit(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
UpdateFlags get_neighbor_face_update_flags() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcMessage(std::string arg1)
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