40template <
int dim,
int spacedim>
47template <
int dim,
int spacedim>
56template <
int dim,
int spacedim>
59 const std::vector<bool> & r_i_a_f,
60 const std::vector<ComponentMask> &nonzero_c)
66 std::make_pair(
std::make_pair(0U, 0U), 0U))
76 nonzero_c.size() == 1 ?
83 return n_components != 1U;
84 }) == n_nonzero_components_table.end())
86 Assert(restriction_is_additive_flags.size() == this->n_dofs_per_cell(),
88 this->n_dofs_per_cell()));
90 for (
unsigned int i = 0; i < nonzero_components.size(); ++i)
92 Assert(nonzero_components[i].size() == this->n_components(),
94 Assert(nonzero_components[i].n_selected_components() >= 1,
97 Assert(n_nonzero_components_table[i] <= this->n_components(),
104 if (this->is_primitive())
106 system_to_component_table.resize(this->n_dofs_per_cell());
107 for (
unsigned int j = 0; j < this->n_dofs_per_cell(); ++j)
108 system_to_component_table[j] = std::pair<unsigned, unsigned>(0, j);
110 face_system_to_component_table.resize(this->n_unique_faces());
111 for (
unsigned int f = 0; f < this->n_unique_faces(); ++f)
113 face_system_to_component_table[f].resize(this->n_dofs_per_face(f));
114 for (
unsigned int j = 0; j < this->n_dofs_per_face(f); ++j)
115 face_system_to_component_table[f][j] =
116 std::pair<unsigned, unsigned>(0, j);
120 for (
unsigned int j = 0; j < this->n_dofs_per_cell(); ++j)
121 system_to_base_table[j] = std::make_pair(std::make_pair(0U, 0U), j);
123 face_system_to_base_table.resize(this->n_unique_faces());
124 for (
unsigned int f = 0; f < this->n_unique_faces(); ++f)
126 face_system_to_base_table[f].resize(this->n_dofs_per_face(f));
127 for (
unsigned int j = 0; j < this->n_dofs_per_face(f); ++j)
128 face_system_to_base_table[f][j] =
129 std::make_pair(std::make_pair(0U, 0U), j);
133 base_to_block_indices.reinit(1, 1);
140 ref < RefinementCase<dim>::isotropic_refinement + 1;
154 adjust_quad_dof_index_for_face_orientation_table.resize(
155 this->n_unique_quads());
157 for (
unsigned int f = 0; f < this->n_unique_quads(); ++f)
159 adjust_quad_dof_index_for_face_orientation_table[f] =
165 adjust_quad_dof_index_for_face_orientation_table[f].fill(0);
169 unit_face_support_points.resize(this->n_unique_faces());
170 generalized_face_support_points.resize(this->n_unique_faces());
175template <
int dim,
int spacedim>
176std::pair<std::unique_ptr<FiniteElement<dim, spacedim>>,
unsigned int>
179 return {this->clone(), multiplicity};
184template <
int dim,
int spacedim>
189 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
195template <
int dim,
int spacedim>
199 const unsigned int)
const
201 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
207template <
int dim,
int spacedim>
212 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
218template <
int dim,
int spacedim>
222 const unsigned int)
const
224 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
230template <
int dim,
int spacedim>
235 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
241template <
int dim,
int spacedim>
246 const unsigned int)
const
248 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
254template <
int dim,
int spacedim>
259 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
265template <
int dim,
int spacedim>
270 const unsigned int)
const
272 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
278template <
int dim,
int spacedim>
283 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
289template <
int dim,
int spacedim>
294 const unsigned int)
const
296 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
301template <
int dim,
int spacedim>
304 const bool isotropic_restriction_only,
305 const bool isotropic_prolongation_only)
308 ref_case <= RefinementCase<dim>::isotropic_refinement;
311 const unsigned int nc =
314 for (
unsigned int i = 0; i < nc; ++i)
316 if (this->restriction[ref_case - 1][i].m() !=
317 this->n_dofs_per_cell() &&
318 (!isotropic_restriction_only ||
320 this->restriction[ref_case - 1][i].reinit(this->n_dofs_per_cell(),
321 this->n_dofs_per_cell());
322 if (this->prolongation[ref_case - 1][i].m() !=
323 this->n_dofs_per_cell() &&
324 (!isotropic_prolongation_only ||
326 this->prolongation[ref_case - 1][i].reinit(this->n_dofs_per_cell(),
327 this->n_dofs_per_cell());
333template <
int dim,
int spacedim>
336 const unsigned int child,
343 "Restriction matrices are only available for refined cells!"));
349 Assert(restriction[refinement_case - 1][child].n() == this->n_dofs_per_cell(),
350 ExcProjectionVoid());
351 return restriction[refinement_case - 1][child];
356template <
int dim,
int spacedim>
359 const unsigned int child,
366 "Prolongation matrices are only available for refined cells!"));
374 Assert(prolongation[refinement_case - 1][child].n() ==
375 this->n_dofs_per_cell(),
377 return prolongation[refinement_case - 1][child];
382template <
int dim,
int spacedim>
385 const unsigned int index)
const
389 return first_block_of_base(component_to_base_table[index].
first.first) +
390 component_to_base_table[
index].second;
394template <
int dim,
int spacedim>
406 std::vector<bool>
mask(this->n_components(),
false);
412template <
int dim,
int spacedim>
418 this->n_components());
425 std::vector<bool>
mask(this->n_components(),
false);
434template <
int dim,
int spacedim>
441 this->n_components());
448 std::vector<bool>
mask(this->n_components(),
false);
459template <
int dim,
int spacedim>
470 std::vector<bool> component_mask(this->n_components(),
false);
471 for (
unsigned int c = 0; c < this->n_components(); ++c)
472 if (block_mask[component_to_block_index(c)] ==
true)
473 component_mask[c] =
true;
475 return component_mask;
480template <
int dim,
int spacedim>
487 return block_mask(component_mask(scalar));
491template <
int dim,
int spacedim>
498 return block_mask(component_mask(vector));
502template <
int dim,
int spacedim>
509 return block_mask(component_mask(sym_tensor));
514template <
int dim,
int spacedim>
535 std::vector<bool> block_mask(this->
n_blocks(),
false);
536 for (
unsigned int c = 0; c < this->n_components();)
538 const unsigned int block = component_to_block_index(c);
539 if (component_mask[c] ==
true)
540 block_mask[block] =
true;
546 while ((c < this->n_components()) &&
547 (component_to_block_index(c) == block))
549 Assert(component_mask[c] == block_mask[block],
551 "The component mask argument given to this function "
552 "is not a mask where the individual components belonging "
553 "to one block of the finite element are either all "
554 "selected or not selected. You can't call this function "
555 "with a component mask that splits blocks."));
566template <
int dim,
int spacedim>
569 const unsigned int face,
570 const bool face_orientation,
571 const bool face_flip,
572 const bool face_rotation)
const
590 if ((face_orientation !=
true) || (face_flip !=
false) ||
591 (face_rotation !=
false))
592 Assert((this->n_dofs_per_line() <= 1) && (this->n_dofs_per_quad(face) <= 1),
594 "The function in this base class can not handle this case. "
595 "Rather, the derived class you are using must provide "
596 "an overloaded version but apparently hasn't done so. See "
597 "the documentation of this function for more information."));
601 if (face_index < this->get_first_face_line_index(face))
606 const unsigned int face_vertex = face_index / this->n_dofs_per_vertex();
607 const unsigned int dof_index_on_vertex =
608 face_index % this->n_dofs_per_vertex();
615 (face_orientation ? 1 : 0) + (face_rotation ? 2 : 0) +
616 (face_flip ? 4 : 0)) *
617 this->n_dofs_per_vertex() +
618 dof_index_on_vertex);
620 else if (face_index < this->get_first_face_quad_index(face))
625 const unsigned int index =
626 face_index - this->get_first_face_line_index(face);
628 const unsigned int face_line =
index / this->n_dofs_per_line();
629 const unsigned int dof_index_on_line =
index % this->n_dofs_per_line();
632 this->get_first_line_index() +
635 (face_orientation ? 1 : 0) +
636 (face_rotation ? 2 : 0) +
637 (face_flip ? 4 : 0)) *
638 this->n_dofs_per_line() +
647 const unsigned int index =
648 face_index - this->get_first_face_quad_index(face);
650 return (this->get_first_quad_index(face) + index);
656template <
int dim,
int spacedim>
659 const unsigned int index,
660 const unsigned int face,
661 const bool face_orientation,
662 const bool face_flip,
663 const bool face_rotation)
const
683 Assert(adjust_quad_dof_index_for_face_orientation_table
684 [this->n_unique_quads() == 1 ? 0 : face]
689 this->n_dofs_per_quad(face),
692 adjust_quad_dof_index_for_face_orientation_table
693 [this->n_unique_quads() == 1 ? 0 : face](
index,
694 (face_orientation ? 4 : 0) +
695 (face_flip ? 2 : 0) +
696 (face_rotation ? 1 : 0));
701template <
int dim,
int spacedim>
704 const unsigned int index,
705 const bool line_orientation)
const
715 Assert(adjust_line_dof_index_for_line_orientation_table.size() ==
716 this->n_dofs_per_line(),
718 if (line_orientation)
721 return index + adjust_line_dof_index_for_line_orientation_table[
index];
726template <
int dim,
int spacedim>
731 ref_case < RefinementCase<dim>::isotropic_refinement + 1;
733 for (
unsigned int c = 0;
739 Assert((prolongation[ref_case - 1][c].m() == this->n_dofs_per_cell()) ||
740 (prolongation[ref_case - 1][c].m() == 0),
742 Assert((prolongation[ref_case - 1][c].n() == this->n_dofs_per_cell()) ||
743 (prolongation[ref_case - 1][c].n() == 0),
745 if ((prolongation[ref_case - 1][c].m() == 0) ||
746 (prolongation[ref_case - 1][c].n() == 0))
754template <
int dim,
int spacedim>
759 ref_case < RefinementCase<dim>::isotropic_refinement + 1;
761 for (
unsigned int c = 0;
767 Assert((restriction[ref_case - 1][c].m() == this->n_dofs_per_cell()) ||
768 (restriction[ref_case - 1][c].m() == 0),
770 Assert((restriction[ref_case - 1][c].n() == this->n_dofs_per_cell()) ||
771 (restriction[ref_case - 1][c].n() == 0),
773 if ((restriction[ref_case - 1][c].m() == 0) ||
774 (restriction[ref_case - 1][c].n() == 0))
782template <
int dim,
int spacedim>
789 for (
unsigned int c = 0;
795 Assert((prolongation[ref_case - 1][c].m() == this->n_dofs_per_cell()) ||
796 (prolongation[ref_case - 1][c].m() == 0),
798 Assert((prolongation[ref_case - 1][c].n() == this->n_dofs_per_cell()) ||
799 (prolongation[ref_case - 1][c].n() == 0),
801 if ((prolongation[ref_case - 1][c].m() == 0) ||
802 (prolongation[ref_case - 1][c].n() == 0))
810template <
int dim,
int spacedim>
817 for (
unsigned int c = 0;
823 Assert((restriction[ref_case - 1][c].m() == this->n_dofs_per_cell()) ||
824 (restriction[ref_case - 1][c].m() == 0),
826 Assert((restriction[ref_case - 1][c].n() == this->n_dofs_per_cell()) ||
827 (restriction[ref_case - 1][c].n() == 0),
829 if ((restriction[ref_case - 1][c].m() == 0) ||
830 (restriction[ref_case - 1][c].n() == 0))
838template <
int dim,
int spacedim>
846 const unsigned int face_no = 0;
849 return (this->n_dofs_per_face(face_no) == 0) ||
850 (interface_constraints.m() != 0);
857template <
int dim,
int spacedim>
866template <
int dim,
int spacedim>
874 const unsigned int face_no = 0;
879 ExcMessage(
"Constraints for this element are only implemented "
880 "for the case that faces are refined isotropically "
881 "(which is always the case in 2d, and in 3d requires "
882 "that the neighboring cell of a coarse cell presents "
883 "exactly four children on the common face)."));
884 Assert((this->n_dofs_per_face(face_no) == 0) ||
885 (interface_constraints.m() != 0),
886 ExcMessage(
"The finite element for which you try to obtain "
887 "hanging node constraints does not appear to "
891 Assert((interface_constraints.m() == 0) && (interface_constraints.n() == 0),
892 ExcWrongInterfaceMatrixSize(interface_constraints.m(),
893 interface_constraints.n()));
895 return interface_constraints;
900template <
int dim,
int spacedim>
907 const unsigned int face_no = 0;
914 return {this->n_dofs_per_vertex() + 2 * this->n_dofs_per_line(),
915 this->n_dofs_per_face(face_no)};
917 return {5 * this->n_dofs_per_vertex() + 12 * this->n_dofs_per_line() +
918 4 * this->n_dofs_per_quad(face_no),
919 this->n_dofs_per_face(face_no)};
928template <
int dim,
int spacedim>
944template <
int dim,
int spacedim>
949 const unsigned int)
const
961template <
int dim,
int spacedim>
967 const unsigned int)
const
979template <
int dim,
int spacedim>
980std::vector<std::pair<unsigned int, unsigned int>>
985 return std::vector<std::pair<unsigned int, unsigned int>>();
990template <
int dim,
int spacedim>
991std::vector<std::pair<unsigned int, unsigned int>>
996 return std::vector<std::pair<unsigned int, unsigned int>>();
1001template <
int dim,
int spacedim>
1002std::vector<std::pair<unsigned int, unsigned int>>
1005 const unsigned int)
const
1008 return std::vector<std::pair<unsigned int, unsigned int>>();
1013template <
int dim,
int spacedim>
1017 const unsigned int)
const
1025template <
int dim,
int spacedim>
1032 return ((
typeid(*
this) ==
typeid(f)) && (this->get_name() == f.
get_name()) &&
1040template <
int dim,
int spacedim>
1045 return !(*
this == f);
1050template <
int dim,
int spacedim>
1051const std::vector<Point<dim>> &
1066template <
int dim,
int spacedim>
1075template <
int dim,
int spacedim>
1076const std::vector<Point<dim>> &
1081 return ((generalized_support_points.size() == 0) ?
1082 unit_support_points :
1083 generalized_support_points);
1088template <
int dim,
int spacedim>
1092 return (get_generalized_support_points().size() != 0);
1097template <
int dim,
int spacedim>
1103 ExcFEHasNoSupportPoints());
1109template <
int dim,
int spacedim>
1110const std::vector<
Point<dim - 1>> &
1112 const unsigned int face_no)
const
1118 Assert((unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
1120 (unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
1121 .size() == this->n_dofs_per_face(face_no)),
1123 return unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no];
1128template <
int dim,
int spacedim>
1131 const unsigned int face_no)
const
1133 return (unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
1139template <
int dim,
int spacedim>
1142 const unsigned int index,
1143 const unsigned int face_no)
const
1146 Assert(unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
1147 .size() == this->n_dofs_per_face(face_no),
1148 ExcFEHasNoSupportPoints());
1149 return unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
1155template <
int dim,
int spacedim>
1158 const unsigned int)
const
1165template <
int dim,
int spacedim>
1171 const unsigned int n_total_components = this->n_components();
1172 Assert((n_total_components ==
mask.size()) || (
mask.size() == 0),
1173 ExcMessage(
"The given ComponentMask has the wrong size."));
1175 const unsigned int n_selected =
1176 mask.n_selected_components(n_total_components);
1178 ExcMessage(
"You need at least one selected component."));
1180 const unsigned int first_selected =
1181 mask.first_selected_component(n_total_components);
1185 for (
unsigned int c = 0; c < n_total_components; ++c)
1186 Assert((c < first_selected && (!mask[c])) ||
1187 (c >= first_selected && c < first_selected + n_selected &&
1189 (c >= first_selected + n_selected && !mask[c]),
1190 ExcMessage(
"Error: the given ComponentMask is not contiguous!"));
1193 return get_sub_fe(first_selected, n_selected);
1198template <
int dim,
int spacedim>
1201 const unsigned int first_component,
1202 const unsigned int n_selected_components)
const
1206 Assert(first_component == 0 && n_selected_components == this->n_components(),
1208 "You can only select a whole FiniteElement, not a part of one."));
1210 (void)first_component;
1211 (void)n_selected_components;
1217template <
int dim,
int spacedim>
1218std::pair<Table<2, bool>, std::vector<unsigned int>>
1222 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
1224 std::vector<unsigned int>(this->n_components()));
1229template <
int dim,
int spacedim>
1234 std::vector<double> &)
const
1236 Assert(has_generalized_support_points(),
1237 ExcMessage(
"The element for which you are calling the current "
1238 "function does not have generalized support points (see "
1239 "the glossary for a definition of generalized support "
1240 "points). Consequently, the current function can not "
1241 "be defined and is not implemented by the element."));
1247template <
int dim,
int spacedim>
1268template <
int dim,
int spacedim>
1269std::vector<unsigned int>
1271 const std::vector<ComponentMask> &nonzero_components)
1273 std::vector<unsigned int> retval(nonzero_components.size());
1274 for (
unsigned int i = 0; i < nonzero_components.size(); ++i)
1275 retval[i] = nonzero_components[i].n_selected_components();
1284template <
int dim,
int spacedim>
1285std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1294 return get_data(flags,
1303template <
int dim,
int spacedim>
1304std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1313 return get_data(flags,
1322template <
int dim,
int spacedim>
1326 const unsigned int face_no,
1330 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
1340 fill_fe_face_values(cell,
1352template <
int dim,
int spacedim>
1356 const unsigned int face_no,
1360 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
1369 ExcMessage(
"Use of a deprecated interface, please implement "
1370 "fill_fe_face_values taking a hp::QCollection argument"));
1375 (void)mapping_internal;
1384template <
int dim,
int spacedim>
1385std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1394 return get_data(flags,
1403template <
int dim,
int spacedim>
bool represents_the_all_selected_mask() const
unsigned int size() const
bool represents_the_all_selected_mask() const
unsigned int size() const
const unsigned int components
unsigned int n_dofs_per_cell() const
unsigned int n_dofs_per_line() const
unsigned int n_components() const
virtual std::size_t memory_consumption() const
bool constraints_are_implemented(const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
bool isotropic_prolongation_is_implemented() const
virtual std::string get_name() const =0
virtual Point< dim > unit_support_point(const unsigned int index) const
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const
const std::vector< unsigned int > n_nonzero_components_table
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
bool prolongation_is_implemented() const
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const
virtual std::unique_ptr< InternalDataBase > get_face_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const hp::QCollection< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
virtual Point< dim - 1 > unit_face_support_point(const unsigned int index, const unsigned int face_no=0) const
static std::vector< unsigned int > compute_n_nonzero_components(const std::vector< ComponentMask > &nonzero_components)
virtual std::unique_ptr< InternalDataBase > get_subface_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
bool has_face_support_points(const unsigned int face_no=0) const
const bool cached_primitivity
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other, const unsigned int face_no=0) const
ComponentMask component_mask(const FEValuesExtractors::Scalar &scalar) const
bool operator!=(const FiniteElement< dim, spacedim > &) const
bool has_support_points() const
bool has_generalized_support_points() const
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const
virtual bool operator==(const FiniteElement< dim, spacedim > &fe) const
unsigned int adjust_line_dof_index_for_line_orientation(const unsigned int index, const bool line_orientation) const
const std::vector< Point< dim > > & get_unit_support_points() const
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
const std::vector< Point< dim - 1 > > & get_unit_face_support_points(const unsigned int face_no=0) const
bool isotropic_restriction_is_implemented() const
std::vector< int > adjust_line_dof_index_for_line_orientation_table
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false) const
const FiniteElement< dim, spacedim > & get_sub_fe(const ComponentMask &mask) const
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const
unsigned int component_to_block_index(const unsigned int component) const
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
BlockMask block_mask(const FEValuesExtractors::Scalar &scalar) const
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > system_to_base_table
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
FiniteElement(const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
virtual Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
const FullMatrix< double > & constraints(const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
std::pair< std::unique_ptr< FiniteElement< dim, spacedim > >, unsigned int > operator^(const unsigned int multiplicity) const
const std::vector< bool > restriction_is_additive_flags
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > component_to_base_table
TableIndices< 2 > interface_constraints_size() const
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
bool restriction_is_implemented() const
FullMatrix< double > interface_constraints
const std::vector< Point< dim > > & get_generalized_support_points() const
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const
virtual std::size_t memory_consumption() const
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const
unsigned int adjust_quad_dof_index_for_face_orientation(const unsigned int index, const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation) const
const std::vector< ComponentMask > nonzero_components
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual bool hp_constraints_are_implemented() const
Abstract base class for mapping classes.
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
@ update_default
No update.
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
@ neither_element_dominates
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
std::enable_if< IsBlockVector< VectorType >::value, unsignedint >::type n_blocks(const VectorType &vector)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
constexpr const ReferenceCell Quadrilateral
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
std::vector< Point< dim > > unit_support_points(const std::vector< Point< 1 > > &line_support_points, const std::vector< unsigned int > &renumbering)
static const unsigned int invalid_unsigned_int
static unsigned int n_children(const RefinementCase< dim > &refinement_case)