24 template <
int dim,
int spacedim>
34 , cell_quadrature(quadrature)
35 , face_quadrature(face_quadrature)
36 , hp_capability_enabled(false)
37 , cell_update_flags(update_flags)
38 , neighbor_cell_update_flags(update_flags)
39 , face_update_flags(face_update_flags)
40 , neighbor_face_update_flags(face_update_flags)
41 , local_dof_indices(fe.n_dofs_per_cell())
42 , neighbor_dof_indices(fe.n_dofs_per_cell())
47 template <
int dim,
int spacedim>
59 , cell_quadrature(quadrature)
60 , face_quadrature(face_quadrature)
61 , hp_capability_enabled(false)
62 , cell_update_flags(update_flags)
63 , neighbor_cell_update_flags(neighbor_update_flags)
64 , face_update_flags(face_update_flags)
65 , neighbor_face_update_flags(neighbor_face_update_flags)
66 , local_dof_indices(fe.n_dofs_per_cell())
67 , neighbor_dof_indices(fe.n_dofs_per_cell())
72 template <
int dim,
int spacedim>
90 template <
int dim,
int spacedim>
104 neighbor_update_flags,
107 neighbor_face_update_flags)
112 template <
int dim,
int spacedim>
120 : mapping_collection(&mapping_collection)
121 , fe_collection(&fe_collection)
122 , cell_quadrature_collection(cell_quadrature_collection)
123 , face_quadrature_collection(face_quadrature_collection)
124 , hp_capability_enabled(true)
125 , cell_update_flags(cell_update_flags)
126 , neighbor_cell_update_flags(cell_update_flags)
127 , face_update_flags(face_update_flags)
128 , neighbor_face_update_flags(face_update_flags)
136 template <
int dim,
int spacedim>
146 : mapping_collection(&mapping_collection)
147 , fe_collection(&fe_collection)
148 , cell_quadrature_collection(cell_quadrature_collection)
149 , face_quadrature_collection(face_quadrature_collection)
150 , hp_capability_enabled(true)
151 , cell_update_flags(cell_update_flags)
152 , neighbor_cell_update_flags(neighbor_cell_update_flags)
153 , face_update_flags(face_update_flags)
154 , neighbor_face_update_flags(neighbor_face_update_flags)
162 template <
int dim,
int spacedim>
169 :
ScratchData(fe_collection.get_reference_cell_default_linear_mapping(),
171 cell_quadrature_collection,
173 face_quadrature_collection,
179 template <
int dim,
int spacedim>
188 :
ScratchData(fe_collection.get_reference_cell_default_linear_mapping(),
190 cell_quadrature_collection,
192 neighbor_cell_update_flags,
193 face_quadrature_collection,
195 neighbor_face_update_flags)
200 template <
int dim,
int spacedim>
203 : mapping(scratch.mapping)
205 , cell_quadrature(scratch.cell_quadrature)
206 , face_quadrature(scratch.face_quadrature)
207 , mapping_collection(scratch.mapping_collection)
208 , fe_collection(scratch.fe_collection)
209 , cell_quadrature_collection(scratch.cell_quadrature_collection)
210 , face_quadrature_collection(scratch.face_quadrature_collection)
211 , hp_capability_enabled(scratch.hp_capability_enabled)
212 , cell_update_flags(scratch.cell_update_flags)
213 , neighbor_cell_update_flags(scratch.neighbor_cell_update_flags)
214 , face_update_flags(scratch.face_update_flags)
215 , neighbor_face_update_flags(scratch.neighbor_face_update_flags)
216 , local_dof_indices(scratch.local_dof_indices)
217 , neighbor_dof_indices(scratch.neighbor_dof_indices)
218 , user_data_storage(scratch.user_data_storage)
219 , internal_data_storage(scratch.internal_data_storage)
224 template <
int dim,
int spacedim>
229 if (hp_capability_enabled ==
false)
232 fe_values = std::make_unique<FEValues<dim, spacedim>>(
233 *mapping, *fe, cell_quadrature, cell_update_flags);
235 fe_values->reinit(cell);
236 local_dof_indices.resize(fe_values->dofs_per_cell);
237 cell->get_dof_indices(local_dof_indices);
238 current_fe_values = fe_values.get();
244 hp_fe_values = std::make_unique<hp::FEValues<dim, spacedim>>(
247 cell_quadrature_collection,
250 hp_fe_values->
reinit(cell);
251 const auto &fe_values = hp_fe_values->get_present_fe_values();
254 (*fe_collection)[cell->active_fe_index()].n_dofs_per_cell(),
255 fe_values.dofs_per_cell);
256 local_dof_indices.resize(fe_values.dofs_per_cell);
257 cell->get_dof_indices(local_dof_indices);
259 current_fe_values = &fe_values;
266 template <
int dim,
int spacedim>
270 const unsigned int face_no)
272 if (hp_capability_enabled ==
false)
275 fe_face_values = std::make_unique<FEFaceValues<dim, spacedim>>(
276 *mapping, *fe, face_quadrature, face_update_flags);
278 fe_face_values->reinit(cell, face_no);
279 local_dof_indices.resize(fe->n_dofs_per_cell());
280 cell->get_dof_indices(local_dof_indices);
281 current_fe_values = fe_face_values.get();
282 return *fe_face_values;
286 return reinit(cell, cell, face_no);
292 template <
int dim,
int spacedim>
298 const unsigned int face_no)
300 Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
302 if (!hp_fe_face_values)
303 hp_fe_face_values = std::make_unique<hp::FEFaceValues<dim, spacedim>>(
306 face_quadrature_collection,
309 if (neighbor_cell == cell)
311 hp_fe_face_values->reinit(cell, face_no);
320 const unsigned int dominated_fe_index =
321 fe_collection->find_dominated_fe(
322 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
324 hp_fe_face_values->reinit(cell,
330 const auto &fe_face_values = hp_fe_face_values->get_present_fe_values();
331 const auto &fe = (*fe_collection)[cell->active_fe_index()];
333 local_dof_indices.resize(fe.n_dofs_per_cell());
334 cell->get_dof_indices(local_dof_indices);
336 current_fe_values = &fe_face_values;
337 return fe_face_values;
342 template <
int dim,
int spacedim>
346 const unsigned int face_no,
347 const unsigned int subface_no)
351 if (hp_capability_enabled ==
false)
353 if (!fe_subface_values)
355 std::make_unique<FESubfaceValues<dim, spacedim>>(
356 *mapping, *fe, face_quadrature, face_update_flags);
357 fe_subface_values->reinit(cell, face_no, subface_no);
358 local_dof_indices.resize(fe->n_dofs_per_cell());
359 cell->get_dof_indices(local_dof_indices);
361 current_fe_values = fe_subface_values.get();
362 return *fe_subface_values;
366 return reinit(cell, cell, face_no, subface_no);
370 return reinit(cell, face_no);
375 template <
int dim,
int spacedim>
381 const unsigned int face_no,
382 const unsigned int subface_no)
384 Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
388 if (!hp_fe_subface_values)
389 hp_fe_subface_values =
390 std::make_unique<hp::FESubfaceValues<dim, spacedim>>(
393 face_quadrature_collection,
396 if (neighbor_cell == cell)
398 hp_fe_subface_values->reinit(cell, face_no, subface_no);
407 const unsigned int dominated_fe_index =
408 fe_collection->find_dominated_fe(
409 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
411 hp_fe_subface_values->reinit(cell,
418 const auto &fe_subface_values =
419 hp_fe_subface_values->get_present_fe_values();
420 const auto &fe = (*fe_collection)[cell->active_fe_index()];
422 local_dof_indices.resize(fe.n_dofs_per_cell());
423 cell->get_dof_indices(local_dof_indices);
425 current_fe_values = &fe_subface_values;
426 return fe_subface_values;
430 return reinit(cell, neighbor_cell, face_no);
436 template <
int dim,
int spacedim>
440 const unsigned int face_no,
441 const unsigned int sub_face_no,
444 const unsigned int face_no_neighbor,
445 const unsigned int sub_face_no_neighbor)
447 if (hp_capability_enabled ==
false)
449 if (!interface_fe_values)
450 interface_fe_values =
451 std::make_unique<FEInterfaceValues<dim, spacedim>>(
452 *mapping, *fe, face_quadrature, face_update_flags);
453 interface_fe_values->reinit(cell,
458 sub_face_no_neighbor);
460 current_fe_values = &interface_fe_values->get_fe_face_values(0);
461 current_neighbor_fe_values =
462 &interface_fe_values->get_fe_face_values(1);
464 cell_neighbor->get_dof_indices(neighbor_dof_indices);
465 local_dof_indices = interface_fe_values->get_interface_dof_indices();
466 return *interface_fe_values;
471 return *interface_fe_values;
477 template <
int dim,
int spacedim>
482 if (hp_capability_enabled ==
false)
484 if (!neighbor_fe_values)
485 neighbor_fe_values = std::make_unique<FEValues<dim, spacedim>>(
486 *mapping, *fe, cell_quadrature, neighbor_cell_update_flags);
488 neighbor_fe_values->reinit(cell);
489 cell->get_dof_indices(neighbor_dof_indices);
490 current_neighbor_fe_values = neighbor_fe_values.get();
491 return *neighbor_fe_values;
495 if (!neighbor_hp_fe_values)
496 neighbor_hp_fe_values = std::make_unique<hp::FEValues<dim, spacedim>>(
499 cell_quadrature_collection,
500 neighbor_cell_update_flags);
502 neighbor_hp_fe_values->
reinit(cell);
503 const auto &neighbor_fe_values =
504 neighbor_hp_fe_values->get_present_fe_values();
507 (*fe_collection)[cell->active_fe_index()].n_dofs_per_cell(),
508 neighbor_fe_values.dofs_per_cell);
509 neighbor_dof_indices.resize(neighbor_fe_values.dofs_per_cell);
510 cell->get_dof_indices(neighbor_dof_indices);
512 current_neighbor_fe_values = &neighbor_fe_values;
513 return neighbor_fe_values;
519 template <
int dim,
int spacedim>
523 const unsigned int face_no)
525 if (hp_capability_enabled ==
false)
527 if (!neighbor_fe_face_values)
528 neighbor_fe_face_values =
529 std::make_unique<FEFaceValues<dim, spacedim>>(
530 *mapping, *fe, face_quadrature, neighbor_face_update_flags);
531 neighbor_fe_face_values->reinit(cell, face_no);
532 cell->get_dof_indices(neighbor_dof_indices);
533 current_neighbor_fe_values = neighbor_fe_face_values.get();
534 return *neighbor_fe_face_values;
538 return reinit_neighbor(cell, cell, face_no);
544 template <
int dim,
int spacedim>
550 const unsigned int face_no)
552 Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
554 if (!neighbor_hp_fe_face_values)
555 neighbor_hp_fe_face_values =
556 std::make_unique<hp::FEFaceValues<dim, spacedim>>(
559 face_quadrature_collection,
560 neighbor_face_update_flags);
562 if (neighbor_cell == cell)
564 neighbor_hp_fe_face_values->reinit(neighbor_cell, face_no);
573 const unsigned int dominated_fe_index =
574 fe_collection->find_dominated_fe(
575 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
577 neighbor_hp_fe_face_values->reinit(neighbor_cell,
583 const auto &neighbor_fe_face_values =
584 neighbor_hp_fe_face_values->get_present_fe_values();
585 const auto &neighbor_fe =
586 (*fe_collection)[neighbor_cell->active_fe_index()];
588 neighbor_dof_indices.resize(neighbor_fe.n_dofs_per_cell());
589 neighbor_cell->get_dof_indices(neighbor_dof_indices);
591 current_neighbor_fe_values = &neighbor_fe_face_values;
592 return neighbor_fe_face_values;
597 template <
int dim,
int spacedim>
601 const unsigned int face_no,
602 const unsigned int subface_no)
606 if (hp_capability_enabled ==
false)
608 if (!neighbor_fe_subface_values)
609 neighbor_fe_subface_values =
610 std::make_unique<FESubfaceValues<dim, spacedim>>(
611 *mapping, *fe, face_quadrature, neighbor_face_update_flags);
612 neighbor_fe_subface_values->reinit(cell, face_no, subface_no);
613 cell->get_dof_indices(neighbor_dof_indices);
614 current_neighbor_fe_values = neighbor_fe_subface_values.get();
615 return *neighbor_fe_subface_values;
619 return reinit_neighbor(cell, cell, face_no, subface_no);
623 return reinit_neighbor(cell, face_no);
628 template <
int dim,
int spacedim>
634 const unsigned int face_no,
635 const unsigned int subface_no)
637 Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
641 if (!neighbor_hp_fe_subface_values)
642 neighbor_hp_fe_subface_values =
643 std::make_unique<hp::FESubfaceValues<dim, spacedim>>(
646 face_quadrature_collection,
647 neighbor_face_update_flags);
649 if (neighbor_cell == cell)
651 neighbor_hp_fe_subface_values->reinit(neighbor_cell,
662 const unsigned int dominated_fe_index =
663 fe_collection->find_dominated_fe(
664 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
666 neighbor_hp_fe_subface_values->reinit(neighbor_cell,
673 const auto &neighbor_fe_subface_values =
674 neighbor_hp_fe_subface_values->get_present_fe_values();
675 const auto &neighbor_fe =
676 (*fe_collection)[neighbor_cell->active_fe_index()];
678 neighbor_dof_indices.resize(neighbor_fe.n_dofs_per_cell());
679 neighbor_cell->get_dof_indices(neighbor_dof_indices);
681 current_neighbor_fe_values = &neighbor_fe_subface_values;
682 return neighbor_fe_subface_values;
686 return reinit_neighbor(cell, neighbor_cell, face_no);
692 template <
int dim,
int spacedim>
696 Assert(current_fe_values !=
nullptr,
697 ExcMessage(
"You have to initialize the cache using one of the "
698 "reinit functions first!"));
699 return *current_fe_values;
704 template <
int dim,
int spacedim>
708 Assert(interface_fe_values !=
nullptr,
709 ExcMessage(
"You have to initialize the cache using one of the "
710 "reinit functions first!"));
711 return *interface_fe_values;
716 template <
int dim,
int spacedim>
720 Assert(current_neighbor_fe_values !=
nullptr,
721 ExcMessage(
"You have to initialize the cache using one of the "
722 "reinit functions first!"));
723 return *current_neighbor_fe_values;
728 template <
int dim,
int spacedim>
729 const std::vector<Point<spacedim>> &
732 return get_current_fe_values().get_quadrature_points();
737 template <
int dim,
int spacedim>
738 const std::vector<double> &
741 return get_current_fe_values().get_JxW_values();
746 template <
int dim,
int spacedim>
747 const std::vector<double> &
750 return get_current_neighbor_fe_values().get_JxW_values();
755 template <
int dim,
int spacedim>
756 const std::vector<Tensor<1, spacedim>> &
759 return get_current_fe_values().get_normal_vectors();
764 template <
int dim,
int spacedim>
765 const std::vector<Tensor<1, spacedim>> &
768 return get_current_neighbor_fe_values().get_normal_vectors();
773 template <
int dim,
int spacedim>
774 const std::vector<types::global_dof_index> &
777 return local_dof_indices;
782 template <
int dim,
int spacedim>
786 return local_dof_indices.size();
791 template <
int dim,
int spacedim>
792 const std::vector<types::global_dof_index> &
795 return neighbor_dof_indices;
800 template <
int dim,
int spacedim>
804 return neighbor_dof_indices.size();
809 template <
int dim,
int spacedim>
813 return user_data_storage;
818 template <
int dim,
int spacedim>
822 return user_data_storage;
827 template <
int dim,
int spacedim>
831 Assert(hp_capability_enabled ==
false, ExcOnlyAvailableWithoutHP());
838 template <
int dim,
int spacedim>
842 Assert(hp_capability_enabled ==
false, ExcOnlyAvailableWithoutHP());
849 template <
int dim,
int spacedim>
853 Assert(hp_capability_enabled ==
false, ExcOnlyAvailableWithoutHP());
854 return cell_quadrature;
859 template <
int dim,
int spacedim>
863 Assert(hp_capability_enabled ==
false, ExcOnlyAvailableWithoutHP());
864 return face_quadrature;
869 template <
int dim,
int spacedim>
873 Assert(hp_capability_enabled ==
true, ExcOnlyAvailableWithHP());
875 return *mapping_collection;
880 template <
int dim,
int spacedim>
884 Assert(hp_capability_enabled ==
true, ExcOnlyAvailableWithHP());
886 return *fe_collection;
891 template <
int dim,
int spacedim>
895 Assert(hp_capability_enabled ==
true, ExcOnlyAvailableWithHP());
896 return cell_quadrature_collection;
901 template <
int dim,
int spacedim>
905 Assert(hp_capability_enabled ==
true, ExcOnlyAvailableWithHP());
906 return face_quadrature_collection;
911 template <
int dim,
int spacedim>
915 return hp_capability_enabled;
920 template <
int dim,
int spacedim>
924 return cell_update_flags;
929 template <
int dim,
int spacedim>
933 return neighbor_cell_update_flags;
938 template <
int dim,
int spacedim>
942 return face_update_flags;
947 template <
int dim,
int spacedim>
951 return neighbor_face_update_flags;
961#include "scratch_data.inst"
void reinit(const Triangulation< dim, spacedim > &tria)
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
const FiniteElement< dim, spacedim > & get_fe() const
SmartPointer< const hp::FECollection< dim, spacedim > > fe_collection
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)
#define AssertThrow(cond, exc)
typename ActiveSelector::active_cell_iterator active_cell_iterator
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
static const unsigned int invalid_unsigned_int