16 #ifndef dealii__fe_base_h 17 #define dealii__fe_base_h 19 #include <deal.II/base/config.h> 20 #include <deal.II/base/exceptions.h> 21 #include <deal.II/base/subscriptor.h> 22 #include <deal.II/base/point.h> 23 #include <deal.II/base/tensor.h> 24 #include <deal.II/base/table.h> 25 #include <deal.II/base/vector_slice.h> 26 #include <deal.II/base/geometry_info.h> 27 #include <deal.II/lac/full_matrix.h> 28 #include <deal.II/lac/block_indices.h> 29 #include <deal.II/fe/fe_update_flags.h> 34 DEAL_II_NAMESPACE_OPEN
36 template<
int dim,
int spacedim>
class FESystem;
100 this_element_dominates,
101 other_element_dominates,
102 neither_element_dominates,
103 either_element_can_dominate,
221 static const unsigned int dimension = dim;
352 const unsigned int n_components,
353 const unsigned int degree,
360 unsigned int n_dofs_per_vertex ()
const;
365 unsigned int n_dofs_per_line ()
const;
370 unsigned int n_dofs_per_quad ()
const;
375 unsigned int n_dofs_per_hex ()
const;
381 unsigned int n_dofs_per_face ()
const;
387 unsigned int n_dofs_per_cell ()
const;
397 template <
int structdim>
398 unsigned int n_dofs_per_object ()
const;
405 unsigned int n_components ()
const;
412 unsigned int n_blocks ()
const;
427 bool is_primitive ()
const;
435 unsigned int tensor_degree ()
const;
458 void set_primitivity(
const bool value);
490 case this_element_dominates:
491 if ((d2 == this_element_dominates) ||
492 (d2 == either_element_can_dominate) ||
493 (d2 == no_requirements))
494 return this_element_dominates;
496 return neither_element_dominates;
498 case other_element_dominates:
499 if ((d2 == other_element_dominates) ||
500 (d2 == either_element_can_dominate) ||
501 (d2 == no_requirements))
502 return other_element_dominates;
504 return neither_element_dominates;
506 case neither_element_dominates:
507 return neither_element_dominates;
509 case either_element_can_dominate:
510 if (d2 == no_requirements)
511 return either_element_can_dominate;
515 case no_requirements:
520 Assert (
false, ExcInternalError());
523 return neither_element_dominates;
533 return dofs_per_vertex;
543 return dofs_per_line;
553 return dofs_per_quad;
573 return dofs_per_face;
583 return dofs_per_cell;
589 template <
int structdim>
597 return dofs_per_vertex;
599 return dofs_per_line;
601 return dofs_per_quad;
605 Assert (
false, ExcInternalError());
627 return cached_primitivity;
636 cached_primitivity = value;
645 return block_indices_data;
655 return block_indices_data.
size();
674 return ((space & conforming_space) == space);
681 DEAL_II_NAMESPACE_CLOSE
const unsigned int first_hex_index
bool is_primitive() const
static const unsigned int invalid_unsigned_int
const unsigned int components
unsigned int n_dofs_per_object() const
Auxiliary class aiding in the handling of block structures like in BlockVector or FESystem...
const unsigned int dofs_per_quad
const unsigned int degree
unsigned int n_dofs_per_vertex() const
Domination operator&(const Domination d1, const Domination d2)
unsigned int n_dofs_per_hex() const
const unsigned int dofs_per_line
unsigned int n_dofs_per_face() const
unsigned int n_blocks() const
const unsigned int first_face_line_index
const unsigned int first_quad_index
unsigned int n_dofs_per_line() const
const unsigned int dofs_per_hex
#define Assert(cond, exc)
const unsigned int dofs_per_cell
const BlockIndices block_indices_data
unsigned int n_dofs_per_quad() const
unsigned int n_components() const
unsigned int n_dofs_per_cell() const
const unsigned int first_face_quad_index
bool conforms(const Conformity) const
const Conformity conforming_space
const unsigned int dofs_per_face
const unsigned int first_line_index
const unsigned int dofs_per_vertex
void set_primitivity(const bool value)
unsigned int size() const
const BlockIndices & block_indices() const
unsigned int tensor_degree() const