16 #include <deal.II/base/memory_consumption.h> 17 #include <deal.II/fe/mapping.h> 18 #include <deal.II/fe/fe.h> 19 #include <deal.II/fe/fe_values.h> 20 #include <deal.II/base/quadrature.h> 21 #include <deal.II/base/qprojector.h> 22 #include <deal.II/grid/tria.h> 23 #include <deal.II/grid/tria_iterator.h> 24 #include <deal.II/dofs/dof_accessor.h> 25 #include <deal.II/grid/tria_boundary.h> 32 DEAL_II_NAMESPACE_OPEN
38 template <
int dim,
int spacedim>
45 template <
int dim,
int spacedim>
54 template <
int dim,
int spacedim>
57 const std::vector<bool> &r_i_a_f,
58 const std::vector<ComponentMask> &nonzero_c)
69 std::make_pair(
std::make_pair(0U, 0U), 0U)),
86 std::bind (
std::not_equal_to<unsigned int>(),
134 ref<RefinementCase<dim>::isotropic_refinement+1; ++ref)
149 template <
int dim,
int spacedim>
150 std::pair<std::unique_ptr<FiniteElement<dim, spacedim> >,
unsigned int>
153 return {this->clone(), multiplicity};
158 template <
int dim,
int spacedim>
163 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
169 template <
int dim,
int spacedim>
173 const unsigned int)
const 175 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
181 template <
int dim,
int spacedim>
186 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
192 template <
int dim,
int spacedim>
196 const unsigned int)
const 198 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
204 template <
int dim,
int spacedim>
209 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
215 template <
int dim,
int spacedim>
219 const unsigned int)
const 221 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
227 template <
int dim,
int spacedim>
232 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
238 template <
int dim,
int spacedim>
242 const unsigned int)
const 244 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
250 template <
int dim,
int spacedim>
255 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
261 template <
int dim,
int spacedim>
265 const unsigned int)
const 267 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
272 template <
int dim,
int spacedim>
275 const bool isotropic_restriction_only,
276 const bool isotropic_prolongation_only)
279 ref_case <= RefinementCase<dim>::isotropic_refinement; ++ref_case)
283 for (
unsigned int i=0; i<nc; ++i)
285 if (this->restriction[ref_case-1][i].m() != this->dofs_per_cell
288 this->restriction[ref_case-1][i].reinit (this->dofs_per_cell,
289 this->dofs_per_cell);
290 if (this->prolongation[ref_case-1][i].m() != this->dofs_per_cell
293 this->prolongation[ref_case-1][i].reinit (this->dofs_per_cell,
294 this->dofs_per_cell);
300 template <
int dim,
int spacedim>
308 ExcMessage(
"Restriction matrices are only available for refined cells!"));
314 Assert (restriction[refinement_case-1][child].n() == this->dofs_per_cell, ExcProjectionVoid());
315 return restriction[refinement_case-1][child];
320 template <
int dim,
int spacedim>
328 ExcMessage(
"Prolongation matrices are only available for refined cells!"));
336 Assert (prolongation[refinement_case-1][child].n() == this->dofs_per_cell, ExcEmbeddingVoid());
337 return prolongation[refinement_case-1][child];
342 template <
int dim,
int spacedim>
349 return first_block_of_base(component_to_base_table[index].first.first)
350 + component_to_base_table[index].second;
354 template <
int dim,
int spacedim>
372 template <
int dim,
int spacedim>
391 template <
int dim,
int spacedim>
398 this->n_components());
414 template <
int dim,
int spacedim>
426 std::vector<bool> component_mask (this->
n_components(),
false);
428 if (block_mask[component_to_block_index(c)] ==
true)
429 component_mask[c] =
true;
431 return component_mask;
436 template <
int dim,
int spacedim>
443 return block_mask(component_mask(scalar));
447 template <
int dim,
int spacedim>
454 return block_mask(component_mask(vector));
458 template <
int dim,
int spacedim>
465 return block_mask(component_mask(sym_tensor));
470 template <
int dim,
int spacedim>
491 std::vector<bool> block_mask (this->n_blocks(),
false);
494 const unsigned int block = component_to_block_index(c);
495 if (component_mask[c] ==
true)
496 block_mask[block] =
true;
504 (component_to_block_index(c) == block))
506 Assert (component_mask[c] == block_mask[block],
507 ExcMessage (
"The component mask argument given to this function " 508 "is not a mask where the individual components belonging " 509 "to one block of the finite element are either all " 510 "selected or not selected. You can't call this function " 511 "with a component mask that splits blocks."));
522 template <
int dim,
int spacedim>
526 const unsigned int face,
527 const bool face_orientation,
528 const bool face_flip,
529 const bool face_rotation)
const 531 Assert (face_index < this->dofs_per_face,
549 if ((face_orientation !=
true) || (face_flip !=
false) || (face_rotation !=
false))
550 Assert ((this->dofs_per_line <= 1) && (this->dofs_per_quad <= 1),
551 ExcMessage (
"The function in this base class can not handle this case. " 552 "Rather, the derived class you are using must provide " 553 "an overloaded version but apparently hasn't done so. See " 554 "the documentation of this function for more information."));
558 if (face_index < this->first_face_line_index)
563 const unsigned int face_vertex = face_index / this->dofs_per_vertex;
564 const unsigned int dof_index_on_vertex = face_index % this->dofs_per_vertex;
572 * this->dofs_per_vertex
574 dof_index_on_vertex);
576 else if (face_index < this->first_face_quad_index)
581 const unsigned int index = face_index - this->first_face_line_index;
583 const unsigned int face_line = index / this->dofs_per_line;
584 const unsigned int dof_index_on_line = index % this->dofs_per_line;
586 return (this->first_line_index
591 * this->dofs_per_line
601 const unsigned int index = face_index - this->first_face_quad_index;
603 return (this->first_quad_index
604 + face * this->dofs_per_quad
612 template <
int dim,
int spacedim>
615 const bool face_orientation,
616 const bool face_flip,
617 const bool face_rotation)
const 637 Assert (adjust_quad_dof_index_for_face_orientation_table.n_elements()==8*this->dofs_per_quad,
639 return index+adjust_quad_dof_index_for_face_orientation_table(index,4*face_orientation+2*face_flip+face_rotation);
644 template <
int dim,
int spacedim>
647 const bool line_orientation)
const 657 Assert (adjust_line_dof_index_for_line_orientation_table.size()==this->dofs_per_line,
659 if (line_orientation)
662 return index+adjust_line_dof_index_for_line_orientation_table[index];
667 template <
int dim,
int spacedim>
672 ref_case<RefinementCase<dim>::isotropic_refinement+1; ++ref_case)
673 for (
unsigned int c=0;
678 Assert ((prolongation[ref_case-1][c].m() == this->dofs_per_cell) ||
679 (prolongation[ref_case-1][c].m() == 0),
681 Assert ((prolongation[ref_case-1][c].n() == this->dofs_per_cell) ||
682 (prolongation[ref_case-1][c].n() == 0),
684 if ((prolongation[ref_case-1][c].m() == 0) ||
685 (prolongation[ref_case-1][c].n() == 0))
693 template <
int dim,
int spacedim>
698 ref_case<RefinementCase<dim>::isotropic_refinement+1; ++ref_case)
699 for (
unsigned int c=0;
704 Assert ((restriction[ref_case-1][c].m() == this->dofs_per_cell) ||
705 (restriction[ref_case-1][c].m() == 0),
707 Assert ((restriction[ref_case-1][c].n() == this->dofs_per_cell) ||
708 (restriction[ref_case-1][c].n() == 0),
710 if ((restriction[ref_case-1][c].m() == 0) ||
711 (restriction[ref_case-1][c].n() == 0))
719 template <
int dim,
int spacedim>
725 for (
unsigned int c=0;
730 Assert ((prolongation[ref_case-1][c].m() == this->dofs_per_cell) ||
731 (prolongation[ref_case-1][c].m() == 0),
733 Assert ((prolongation[ref_case-1][c].n() == this->dofs_per_cell) ||
734 (prolongation[ref_case-1][c].n() == 0),
736 if ((prolongation[ref_case-1][c].m() == 0) ||
737 (prolongation[ref_case-1][c].n() == 0))
745 template <
int dim,
int spacedim>
751 for (
unsigned int c=0;
756 Assert ((restriction[ref_case-1][c].m() == this->dofs_per_cell) ||
757 (restriction[ref_case-1][c].m() == 0),
759 Assert ((restriction[ref_case-1][c].n() == this->dofs_per_cell) ||
760 (restriction[ref_case-1][c].n() == 0),
762 if ((restriction[ref_case-1][c].m() == 0) ||
763 (restriction[ref_case-1][c].n() == 0))
771 template <
int dim,
int spacedim>
776 return (this->dofs_per_face == 0) || (interface_constraints.m() != 0);
783 template <
int dim,
int spacedim>
792 template <
int dim,
int spacedim>
798 ExcMessage(
"Constraints for this element are only implemented " 799 "for the case that faces are refined isotropically " 800 "(which is always the case in 2d, and in 3d requires " 801 "that the neighboring cell of a coarse cell presents " 802 "exactly four children on the common face)."));
803 Assert ((this->dofs_per_face == 0) || (interface_constraints.m() != 0),
804 ExcMessage (
"The finite element for which you try to obtain " 805 "hanging node constraints does not appear to " 809 Assert ((interface_constraints.m()==0) && (interface_constraints.n()==0),
810 ExcWrongInterfaceMatrixSize(interface_constraints.m(),
811 interface_constraints.n()));
813 return interface_constraints;
818 template <
int dim,
int spacedim>
828 2*this->dofs_per_line,
829 this->dofs_per_face);
832 12*this->dofs_per_line +
833 4*this->dofs_per_quad,
834 this->dofs_per_face);
844 template <
int dim,
int spacedim>
855 ExcInterpolationNotImplemented()));
860 template <
int dim,
int spacedim>
871 ExcInterpolationNotImplemented()));
876 template <
int dim,
int spacedim>
888 ExcInterpolationNotImplemented()));
893 template <
int dim,
int spacedim>
894 std::vector<std::pair<unsigned int, unsigned int> >
899 return std::vector<std::pair<unsigned int, unsigned int> > ();
904 template <
int dim,
int spacedim>
905 std::vector<std::pair<unsigned int, unsigned int> >
910 return std::vector<std::pair<unsigned int, unsigned int> > ();
915 template <
int dim,
int spacedim>
916 std::vector<std::pair<unsigned int, unsigned int> >
921 return std::vector<std::pair<unsigned int, unsigned int> > ();
926 template <
int dim,
int spacedim>
937 template <
int dim,
int spacedim>
943 return ((
typeid(*
this) ==
typeid(f))
955 template <
int dim,
int spacedim>
959 return !(*
this == f);
964 template <
int dim,
int spacedim>
965 const std::vector<Point<dim> > &
972 Assert ((unit_support_points.size() == 0) ||
973 (unit_support_points.size() == this->dofs_per_cell),
975 return unit_support_points;
980 template <
int dim,
int spacedim>
984 return (unit_support_points.size() != 0);
989 template <
int dim,
int spacedim>
990 const std::vector<Point<dim> > &
995 return ((generalized_support_points.size() == 0)
996 ? unit_support_points
997 : generalized_support_points);
1002 template <
int dim,
int spacedim>
1006 return (get_generalized_support_points().size() != 0);
1011 template <
int dim,
int spacedim>
1015 Assert (index < this->dofs_per_cell,
1017 Assert (unit_support_points.size() == this->dofs_per_cell,
1018 ExcFEHasNoSupportPoints ());
1019 return unit_support_points[index];
1024 template <
int dim,
int spacedim>
1025 const std::vector<
Point<dim-1> > &
1032 Assert ((unit_face_support_points.size() == 0) ||
1033 (unit_face_support_points.size() == this->dofs_per_face),
1035 return unit_face_support_points;
1040 template <
int dim,
int spacedim>
1044 return (unit_face_support_points.size() != 0);
1049 template <
int dim,
int spacedim>
1050 const std::vector<
Point<dim-1> > &
1057 return ((generalized_face_support_points.size() == 0)
1058 ? unit_face_support_points
1059 : generalized_face_support_points);
1064 template <
int dim,
int spacedim>
1068 return (generalized_face_support_points.size() != 0);
1073 template <
int dim,
int spacedim>
1077 Assert (index < this->dofs_per_face,
1079 Assert (unit_face_support_points.size() == this->dofs_per_face,
1080 ExcFEHasNoSupportPoints ());
1081 return unit_face_support_points[index];
1086 template <
int dim,
int spacedim>
1090 const unsigned int)
const 1097 template <
int dim,
int spacedim>
1104 const unsigned int n_total_components = this->
n_components();
1106 ExcMessage(
"The given ComponentMask has the wrong size."));
1110 ExcMessage(
"You need at least one selected component."));
1116 for (
unsigned int c=0; c<n_total_components; ++c)
1117 Assert((c<first_selected && (!mask[c]))
1119 (c>=first_selected && c<first_selected+n_selected && mask[c])
1121 (c>=first_selected+n_selected && !mask[c]),
1122 ExcMessage(
"Error: the given ComponentMask is not contiguous!"));
1125 return get_sub_fe(first_selected, n_selected);
1130 template <
int dim,
int spacedim>
1134 const unsigned int n_selected_components)
const 1139 ExcMessage(
"You can only select a whole FiniteElement, not a part of one."));
1141 (void)first_component;
1142 (void)n_selected_components;
1148 template <
int dim,
int spacedim>
1149 std::pair<Table<2,bool>, std::vector<unsigned int> >
1153 return std::pair<Table<2,bool>, std::vector<unsigned int> >
1160 template <
int dim,
int spacedim>
1164 std::vector<double> &)
const 1166 Assert (has_generalized_support_points(),
1167 ExcMessage (
"The element for which you are calling the current " 1168 "function does not have generalized support points (see " 1169 "the glossary for a definition of generalized support " 1170 "points). Consequently, the current function can not " 1171 "be defined and is not implemented by the element."));
1177 template <
int dim,
int spacedim>
1197 template <
int dim,
int spacedim>
1198 std::vector<unsigned int>
1200 const std::vector<ComponentMask> &nonzero_components)
1202 std::vector<unsigned int> retval (nonzero_components.size());
1203 for (
unsigned int i=0; i<nonzero_components.size(); ++i)
1204 retval[i] = nonzero_components[i].n_selected_components();
1212 template <
int dim,
int spacedim>
1213 std::unique_ptr<typename FiniteElement<dim,spacedim>::InternalDataBase>
1219 return get_data (flags, mapping,
1226 template <
int dim,
int spacedim>
1227 std::unique_ptr<typename FiniteElement<dim,spacedim>::InternalDataBase>
1233 return get_data (flags, mapping,
1240 template <
int dim,
int spacedim>
1258 DEAL_II_NAMESPACE_CLOSE
bool represents_the_all_selected_mask() const
virtual std::unique_ptr< InternalDataBase > get_face_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
static const unsigned int invalid_unsigned_int
BlockMask block_mask(const FEValuesExtractors::Scalar &scalar) const
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
std::vector< std::vector< FullMatrix< double > > > restriction
const unsigned int components
#define AssertDimension(dim1, dim2)
virtual std::size_t memory_consumption() const
const std::vector< Point< dim-1 > > & get_unit_face_support_points() const
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual bool operator==(const FiniteElement< dim, spacedim > &fe) const
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const
FullMatrix< double > interface_constraints
bool operator!=(const FiniteElement< dim, spacedim > &) const
const std::vector< Point< dim-1 > > & get_generalized_face_support_points() const
FiniteElement(const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
const FiniteElement< dim, spacedim > & get_sub_fe(const ComponentMask &mask) const
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual Point< dim > unit_support_point(const unsigned int index) const
bool restriction_is_implemented() const
const std::vector< Point< dim > > & get_generalized_support_points() const
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const
virtual Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
#define AssertIndexRange(index, range)
const unsigned int dofs_per_quad
Table< 2, int > adjust_quad_dof_index_for_face_orientation_table
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
bool isotropic_prolongation_is_implemented() const
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > component_to_base_table
bool has_generalized_support_points() const
bool represents_the_all_selected_mask() const
ComponentMask component_mask(const FEValuesExtractors::Scalar &scalar) const
#define AssertThrow(cond, exc)
static std::vector< unsigned int > compute_n_nonzero_components(const std::vector< ComponentMask > &nonzero_components)
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > face_system_to_base_table
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
bool is_primitive() const
const std::vector< unsigned int > n_nonzero_components_table
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
const unsigned int dofs_per_line
bool has_face_support_points() const
bool has_support_points() const
unsigned int component_to_block_index(const unsigned int component) const
static ::ExceptionBase & ExcMessage(std::string arg1)
const FullMatrix< double > & constraints(const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual std::size_t memory_consumption() const
std::vector< std::vector< FullMatrix< double > > > prolongation
#define Assert(cond, exc)
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
bool constraints_are_implemented(const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
unsigned int size() const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
virtual std::string get_name() const =0
bool has_generalized_face_support_points() const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) 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
virtual Point< dim-1 > unit_face_support_point(const unsigned int index) const
const unsigned int dofs_per_cell
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
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
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
unsigned int n_components() const
std::vector< std::pair< unsigned int, unsigned int > > system_to_component_table
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
virtual bool hp_constraints_are_implemented() const
unsigned int first_selected_component(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
const unsigned int dofs_per_face
bool prolongation_is_implemented() const
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
const bool cached_primitivity
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const
const std::vector< ComponentMask > nonzero_components
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
static ::ExceptionBase & ExcNotImplemented()
std::vector< std::pair< unsigned int, unsigned int > > face_system_to_component_table
bool isotropic_restriction_is_implemented() const
BlockIndices base_to_block_indices
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const
unsigned int size() const
unsigned int adjust_quad_dof_index_for_face_orientation(const unsigned int index, const bool face_orientation, const bool face_flip, const bool face_rotation) const
std::vector< int > adjust_line_dof_index_for_line_orientation_table
TableIndices< 2 > interface_constraints_size() const
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
void fill(InputIterator entries, const bool C_style_indexing=true)
const std::vector< bool > restriction_is_additive_flags
std::pair< std::unique_ptr< FiniteElement< dim, spacedim > >, unsigned int > operator^(const unsigned int multiplicity) 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
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > system_to_base_table
static ::ExceptionBase & ExcInternalError()