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> 31 DEAL_II_NAMESPACE_OPEN
37 template <
int dim,
int spacedim>
44 template <
int dim,
int spacedim>
50 template <
int dim,
int spacedim>
59 template <
int dim,
int spacedim>
62 const std::vector<bool> &r_i_a_f,
63 const std::vector<ComponentMask> &nonzero_c)
74 std::make_pair(
std::make_pair(0U, 0U), 0U)),
139 ref<RefinementCase<dim>::isotropic_refinement+1; ++ref)
154 template <
int dim,
int spacedim>
161 template <
int dim,
int spacedim>
166 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
172 template <
int dim,
int spacedim>
176 const unsigned int)
const 178 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
184 template <
int dim,
int spacedim>
189 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
195 template <
int dim,
int spacedim>
199 const unsigned int)
const 201 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
207 template <
int dim,
int spacedim>
212 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
218 template <
int dim,
int spacedim>
222 const unsigned int)
const 224 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
230 template <
int dim,
int spacedim>
235 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
241 template <
int dim,
int spacedim>
245 const unsigned int)
const 247 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
253 template <
int dim,
int spacedim>
258 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
264 template <
int dim,
int spacedim>
268 const unsigned int)
const 270 AssertThrow(
false, ExcUnitShapeValuesDoNotExist());
275 template <
int dim,
int spacedim>
278 const bool isotropic_restriction_only,
279 const bool isotropic_prolongation_only)
282 ref_case <= RefinementCase<dim>::isotropic_refinement; ++ref_case)
286 for (
unsigned int i=0; i<nc; ++i)
288 if (this->restriction[ref_case-1][i].m() != this->dofs_per_cell
291 this->restriction[ref_case-1][i].reinit (this->dofs_per_cell,
292 this->dofs_per_cell);
293 if (this->prolongation[ref_case-1][i].m() != this->dofs_per_cell
296 this->prolongation[ref_case-1][i].reinit (this->dofs_per_cell,
297 this->dofs_per_cell);
303 template <
int dim,
int spacedim>
311 ExcMessage(
"Restriction matrices are only available for refined cells!"));
317 Assert (restriction[refinement_case-1][child].n() == this->dofs_per_cell, ExcProjectionVoid());
318 return restriction[refinement_case-1][child];
323 template <
int dim,
int spacedim>
331 ExcMessage(
"Prolongation matrices are only available for refined cells!"));
339 Assert (prolongation[refinement_case-1][child].n() == this->dofs_per_cell, ExcEmbeddingVoid());
340 return prolongation[refinement_case-1][child];
345 template <
int dim,
int spacedim>
352 return first_block_of_base(component_to_base_table[index].first.first)
353 + component_to_base_table[index].second;
357 template <
int dim,
int spacedim>
375 template <
int dim,
int spacedim>
394 template <
int dim,
int spacedim>
401 this->n_components());
417 template <
int dim,
int spacedim>
429 std::vector<bool> component_mask (this->
n_components(),
false);
431 if (block_mask[component_to_block_index(c)] ==
true)
432 component_mask[c] =
true;
434 return component_mask;
439 template <
int dim,
int spacedim>
446 return block_mask(component_mask(scalar));
450 template <
int dim,
int spacedim>
457 return block_mask(component_mask(vector));
461 template <
int dim,
int spacedim>
468 return block_mask(component_mask(sym_tensor));
473 template <
int dim,
int spacedim>
494 std::vector<bool> block_mask (this->n_blocks(),
false);
497 const unsigned int block = component_to_block_index(c);
498 if (component_mask[c] ==
true)
499 block_mask[block] =
true;
507 (component_to_block_index(c) == block))
509 Assert (component_mask[c] == block_mask[block],
510 ExcMessage (
"The component mask argument given to this function " 511 "is not a mask where the individual components belonging " 512 "to one block of the finite element are either all " 513 "selected or not selected. You can't call this function " 514 "with a component mask that splits blocks."));
525 template <
int dim,
int spacedim>
529 const unsigned int face,
530 const bool face_orientation,
531 const bool face_flip,
532 const bool face_rotation)
const 534 Assert (face_index < this->dofs_per_face,
552 if ((face_orientation !=
true) || (face_flip !=
false) || (face_rotation !=
false))
553 Assert ((this->dofs_per_line <= 1) && (this->dofs_per_quad <= 1),
554 ExcMessage (
"The function in this base class can not handle this case. " 555 "Rather, the derived class you are using must provide " 556 "an overloaded version but apparently hasn't done so. See " 557 "the documentation of this function for more information."));
561 if (face_index < this->first_face_line_index)
566 const unsigned int face_vertex = face_index / this->dofs_per_vertex;
567 const unsigned int dof_index_on_vertex = face_index % this->dofs_per_vertex;
575 * this->dofs_per_vertex
577 dof_index_on_vertex);
579 else if (face_index < this->first_face_quad_index)
584 const unsigned int index = face_index - this->first_face_line_index;
586 const unsigned int face_line = index / this->dofs_per_line;
587 const unsigned int dof_index_on_line = index % this->dofs_per_line;
589 return (this->first_line_index
594 * this->dofs_per_line
604 const unsigned int index = face_index - this->first_face_quad_index;
606 return (this->first_quad_index
607 + face * this->dofs_per_quad
615 template <
int dim,
int spacedim>
618 const bool face_orientation,
619 const bool face_flip,
620 const bool face_rotation)
const 640 Assert (adjust_quad_dof_index_for_face_orientation_table.n_elements()==8*this->dofs_per_quad,
642 return index+adjust_quad_dof_index_for_face_orientation_table(index,4*face_orientation+2*face_flip+face_rotation);
647 template <
int dim,
int spacedim>
650 const bool line_orientation)
const 660 Assert (adjust_line_dof_index_for_line_orientation_table.size()==this->dofs_per_line,
662 if (line_orientation)
665 return index+adjust_line_dof_index_for_line_orientation_table[index];
670 template <
int dim,
int spacedim>
675 ref_case<RefinementCase<dim>::isotropic_refinement+1; ++ref_case)
676 for (
unsigned int c=0;
681 Assert ((prolongation[ref_case-1][c].m() == this->dofs_per_cell) ||
682 (prolongation[ref_case-1][c].m() == 0),
684 Assert ((prolongation[ref_case-1][c].n() == this->dofs_per_cell) ||
685 (prolongation[ref_case-1][c].n() == 0),
687 if ((prolongation[ref_case-1][c].m() == 0) ||
688 (prolongation[ref_case-1][c].n() == 0))
696 template <
int dim,
int spacedim>
701 ref_case<RefinementCase<dim>::isotropic_refinement+1; ++ref_case)
702 for (
unsigned int c=0;
707 Assert ((restriction[ref_case-1][c].m() == this->dofs_per_cell) ||
708 (restriction[ref_case-1][c].m() == 0),
710 Assert ((restriction[ref_case-1][c].n() == this->dofs_per_cell) ||
711 (restriction[ref_case-1][c].n() == 0),
713 if ((restriction[ref_case-1][c].m() == 0) ||
714 (restriction[ref_case-1][c].n() == 0))
722 template <
int dim,
int spacedim>
728 for (
unsigned int c=0;
733 Assert ((prolongation[ref_case-1][c].m() == this->dofs_per_cell) ||
734 (prolongation[ref_case-1][c].m() == 0),
736 Assert ((prolongation[ref_case-1][c].n() == this->dofs_per_cell) ||
737 (prolongation[ref_case-1][c].n() == 0),
739 if ((prolongation[ref_case-1][c].m() == 0) ||
740 (prolongation[ref_case-1][c].n() == 0))
748 template <
int dim,
int spacedim>
754 for (
unsigned int c=0;
759 Assert ((restriction[ref_case-1][c].m() == this->dofs_per_cell) ||
760 (restriction[ref_case-1][c].m() == 0),
762 Assert ((restriction[ref_case-1][c].n() == this->dofs_per_cell) ||
763 (restriction[ref_case-1][c].n() == 0),
765 if ((restriction[ref_case-1][c].m() == 0) ||
766 (restriction[ref_case-1][c].n() == 0))
774 template <
int dim,
int spacedim>
779 return (this->dofs_per_face == 0) || (interface_constraints.m() != 0);
786 template <
int dim,
int spacedim>
795 template <
int dim,
int spacedim>
801 ExcMessage(
"Constraints for this element are only implemented " 802 "for the case that faces are refined isotropically " 803 "(which is always the case in 2d, and in 3d requires " 804 "that the neighboring cell of a coarse cell presents " 805 "exactly four children on the common face)."));
806 Assert ((this->dofs_per_face == 0) || (interface_constraints.m() != 0),
807 ExcMessage (
"The finite element for which you try to obtain " 808 "hanging node constraints does not appear to " 812 Assert ((interface_constraints.m()==0) && (interface_constraints.n()==0),
813 ExcWrongInterfaceMatrixSize(interface_constraints.m(),
814 interface_constraints.n()));
816 return interface_constraints;
821 template <
int dim,
int spacedim>
831 2*this->dofs_per_line,
832 this->dofs_per_face);
835 12*this->dofs_per_line +
836 4*this->dofs_per_quad,
837 this->dofs_per_face);
847 template <
int dim,
int spacedim>
859 ExcInterpolationNotImplemented());
864 template <
int dim,
int spacedim>
876 ExcInterpolationNotImplemented());
881 template <
int dim,
int spacedim>
893 typename FEE::ExcInterpolationNotImplemented());
898 template <
int dim,
int spacedim>
899 std::vector<std::pair<unsigned int, unsigned int> >
904 return std::vector<std::pair<unsigned int, unsigned int> > ();
909 template <
int dim,
int spacedim>
910 std::vector<std::pair<unsigned int, unsigned int> >
915 return std::vector<std::pair<unsigned int, unsigned int> > ();
920 template <
int dim,
int spacedim>
921 std::vector<std::pair<unsigned int, unsigned int> >
926 return std::vector<std::pair<unsigned int, unsigned int> > ();
931 template <
int dim,
int spacedim>
942 template <
int dim,
int spacedim>
953 template <
int dim,
int spacedim>
954 const std::vector<Point<dim> > &
961 Assert ((unit_support_points.size() == 0) ||
962 (unit_support_points.size() == this->dofs_per_cell),
964 return unit_support_points;
969 template <
int dim,
int spacedim>
973 return (unit_support_points.size() != 0);
978 template <
int dim,
int spacedim>
979 const std::vector<Point<dim> > &
986 return ((generalized_support_points.size() == 0)
987 ? unit_support_points
988 : generalized_support_points);
993 template <
int dim,
int spacedim>
997 return (get_generalized_support_points().size() != 0);
1002 template <
int dim,
int spacedim>
1006 Assert (index < this->dofs_per_cell,
1008 Assert (unit_support_points.size() == this->dofs_per_cell,
1009 ExcFEHasNoSupportPoints ());
1010 return unit_support_points[index];
1015 template <
int dim,
int spacedim>
1016 const std::vector<
Point<dim-1> > &
1023 Assert ((unit_face_support_points.size() == 0) ||
1024 (unit_face_support_points.size() == this->dofs_per_face),
1026 return unit_face_support_points;
1031 template <
int dim,
int spacedim>
1035 return (unit_face_support_points.size() != 0);
1040 template <
int dim,
int spacedim>
1041 const std::vector<
Point<dim-1> > &
1048 return ((generalized_face_support_points.size() == 0)
1049 ? unit_face_support_points
1050 : generalized_face_support_points);
1055 template <
int dim,
int spacedim>
1059 return (generalized_face_support_points.size() != 0);
1064 template <
int dim,
int spacedim>
1068 Assert (index < this->dofs_per_face,
1070 Assert (unit_face_support_points.size() == this->dofs_per_face,
1071 ExcFEHasNoSupportPoints ());
1072 return unit_face_support_points[index];
1077 template <
int dim,
int spacedim>
1081 const unsigned int)
const 1088 template <
int dim,
int spacedim>
1089 std::pair<Table<2,bool>, std::vector<unsigned int> >
1093 return std::pair<Table<2,bool>, std::vector<unsigned int> >
1100 template <
int dim,
int spacedim>
1104 std::vector<double> &)
const 1106 Assert (has_generalized_support_points(),
1107 ExcMessage (
"The element for which you are calling the current " 1108 "function does not have generalized support points (see " 1109 "the glossary for a definition of generalized support " 1110 "points). Consequently, the current function can not " 1111 "be defined and is not implemented by the element."));
1116 template <
int dim,
int spacedim>
1119 std::vector<double> &local_dofs,
1120 const std::vector<double> &values)
const 1122 Assert (has_support_points(), ExcFEHasNoSupportPoints());
1123 Assert (values.size() == unit_support_points.size(),
1125 Assert (local_dofs.size() == this->dofs_per_cell,
1130 std::copy(values.begin(), values.end(), local_dofs.begin());
1136 template <
int dim,
int spacedim>
1139 std::vector<double> &local_dofs,
1141 const unsigned int offset)
const 1143 Assert (has_support_points(), ExcFEHasNoSupportPoints());
1144 Assert (values.size() == unit_support_points.size(),
1146 Assert (local_dofs.size() == this->dofs_per_cell,
1151 for (
unsigned int i=0; i<this->dofs_per_cell; ++i)
1153 const std::pair<unsigned int, unsigned int> index
1154 = this->system_to_component_index(i);
1155 local_dofs[i] = values[i](offset+index.first);
1162 template <
int dim,
int spacedim>
1165 std::vector<double> &local_dofs,
1166 const VectorSlice<
const std::vector<std::vector<double> > > &values)
const 1168 Assert (has_support_points(), ExcFEHasNoSupportPoints());
1169 Assert (values[0].size() == unit_support_points.size(),
1171 Assert (local_dofs.size() == this->dofs_per_cell,
1176 for (
unsigned int i=0; i<this->dofs_per_cell; ++i)
1178 const std::pair<unsigned int, unsigned int> index
1179 = this->system_to_component_index(i);
1180 local_dofs[i] = values[index.first][i];
1186 template <
int dim,
int spacedim>
1206 template <
int dim,
int spacedim>
1207 std::vector<unsigned int>
1209 const std::vector<ComponentMask> &nonzero_components)
1211 std::vector<unsigned int> retval (nonzero_components.size());
1212 for (
unsigned int i=0; i<nonzero_components.size(); ++i)
1213 retval[i] = nonzero_components[i].n_selected_components();
1221 template <
int dim,
int spacedim>
1228 return get_data (flags, mapping,
1235 template <
int dim,
int spacedim>
1242 return get_data (flags, mapping,
1249 template <
int dim,
int spacedim>
1267 DEAL_II_NAMESPACE_CLOSE
bool represents_the_all_selected_mask() 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 Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const
FullMatrix< double > interface_constraints
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)
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const 1
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
virtual void convert_generalized_support_point_values_to_nodal_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
bool has_face_support_points() const
bool operator==(const FiniteElement< dim, spacedim > &) 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
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
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
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 InternalDataBase * get_face_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const
virtual Point< dim-1 > unit_face_support_point(const unsigned int index) const
virtual ~InternalDataBase()
const unsigned int dofs_per_cell
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
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
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 bool hp_constraints_are_implemented() 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()
virtual InternalDataBase * get_subface_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const
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
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()