|
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
16 #ifndef dealii_mapping_q_generic_h
17 #define dealii_mapping_q_generic_h
134 template <
int dim,
int spacedim = dim>
151 virtual std::unique_ptr<Mapping<dim, spacedim>>
152 clone()
const override;
270 const unsigned int n_original_q_points);
280 const unsigned int n_original_q_points);
306 shape(
const unsigned int qpoint,
const unsigned int shape_nr)
const;
312 shape(
const unsigned int qpoint,
const unsigned int shape_nr);
318 derivative(
const unsigned int qpoint,
const unsigned int shape_nr)
const;
324 derivative(
const unsigned int qpoint,
const unsigned int shape_nr);
331 const unsigned int shape_nr)
const;
344 const unsigned int shape_nr)
const;
357 const unsigned int shape_nr)
const;
422 std::array<std::vector<Tensor<1, dim>>,
505 mutable std::vector<DerivativeForm<1, dim, spacedim>>
covariant;
519 mutable std::vector<std::vector<Tensor<1, spacedim>>>
aux;
546 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
550 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
555 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
567 &output_data)
const override;
573 const unsigned int face_no,
577 &output_data)
const override;
583 const unsigned int face_no,
584 const unsigned int subface_no,
588 &output_data)
const override;
676 virtual std::vector<Point<spacedim>>
746 template <
int dim,
int spacedim>
747 inline const double &
749 const unsigned int qpoint,
750 const unsigned int shape_nr)
const
758 template <
int dim,
int spacedim>
761 const unsigned int shape_nr)
763 AssertIndexRange(qpoint * n_shape_functions + shape_nr, shape_values.size());
764 return shape_values[qpoint * n_shape_functions + shape_nr];
768 template <
int dim,
int spacedim>
771 const unsigned int qpoint,
772 const unsigned int shape_nr)
const
775 shape_derivatives.size());
776 return shape_derivatives[qpoint * n_shape_functions + shape_nr];
781 template <
int dim,
int spacedim>
784 const unsigned int qpoint,
785 const unsigned int shape_nr)
788 shape_derivatives.size());
789 return shape_derivatives[qpoint * n_shape_functions + shape_nr];
793 template <
int dim,
int spacedim>
796 const unsigned int qpoint,
797 const unsigned int shape_nr)
const
800 shape_second_derivatives.size());
801 return shape_second_derivatives[qpoint * n_shape_functions + shape_nr];
805 template <
int dim,
int spacedim>
808 const unsigned int qpoint,
809 const unsigned int shape_nr)
812 shape_second_derivatives.size());
813 return shape_second_derivatives[qpoint * n_shape_functions + shape_nr];
816 template <
int dim,
int spacedim>
819 const unsigned int qpoint,
820 const unsigned int shape_nr)
const
823 shape_third_derivatives.size());
824 return shape_third_derivatives[qpoint * n_shape_functions + shape_nr];
828 template <
int dim,
int spacedim>
831 const unsigned int qpoint,
832 const unsigned int shape_nr)
835 shape_third_derivatives.size());
836 return shape_third_derivatives[qpoint * n_shape_functions + shape_nr];
840 template <
int dim,
int spacedim>
843 const unsigned int qpoint,
844 const unsigned int shape_nr)
const
847 shape_fourth_derivatives.size());
848 return shape_fourth_derivatives[qpoint * n_shape_functions + shape_nr];
852 template <
int dim,
int spacedim>
855 const unsigned int qpoint,
856 const unsigned int shape_nr)
859 shape_fourth_derivatives.size());
860 return shape_fourth_derivatives[qpoint * n_shape_functions + shape_nr];
865 template <
int dim,
int spacedim>
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
virtual bool preserves_vertex_locations() const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
std::vector< double > shape_values
std::vector< Tensor< 3, dim > > shape_third_derivatives
std::vector< Tensor< 1, dim > > shape_derivatives
AlignedVector< VectorizedArray< double > > hessians_quad
const Tensor< 2, dim > & second_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
const double & shape(const unsigned int qpoint, const unsigned int shape_nr) const
virtual void add_line_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim >> &a) const
std::vector< Tensor< 2, dim > > shape_second_derivatives
#define AssertIndexRange(index, range)
QGaussLobatto< 1 > line_support_points
Table< 2, double > support_point_weights_cell
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const Quadrature< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
AlignedVector< VectorizedArray< double > > values_dofs
const unsigned int polynomial_degree
InternalData(const unsigned int polynomial_degree)
void initialize(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
virtual void add_quad_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim >> &a) const
virtual void transform(const ArrayView< const Tensor< 1, dim >> &input, const MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim >> &output) const override
std::vector< std::vector< Tensor< 1, spacedim > > > aux
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
std::vector< double > volume_elements
Abstract base class for mapping classes.
bool tensor_product_quadrature
internal::MatrixFreeFunctions::ShapeInfo< VectorizedArray< double > > shape_info
#define DEAL_II_NAMESPACE_OPEN
std::vector< Table< 2, double > > support_point_weights_perimeter_to_interior
std::vector< Point< spacedim > > mapping_support_points
AlignedVector< VectorizedArray< double > > scratch
Triangulation< dim, spacedim >::cell_iterator cell_of_current_support_points
Point< dim > transform_real_to_unit_cell_internal(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_p_unit) const
const Tensor< 1, dim > & derivative(const unsigned int qpoint, const unsigned int shape_nr) const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
std::vector< Tensor< 4, dim > > shape_fourth_derivatives
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
unsigned int get_degree() const
AlignedVector< VectorizedArray< double > > gradients_quad
QGaussLobatto< 1 > line_support_points
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim - 1)> unit_tangentials
const Tensor< 3, dim > & third_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
const unsigned int n_shape_functions
AlignedVector< VectorizedArray< double > > values_quad
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
void initialize_face(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
const unsigned int polynomial_degree
#define DEAL_II_NAMESPACE_CLOSE
virtual CellSimilarity::Similarity fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
const Tensor< 4, dim > & fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
virtual std::size_t memory_consumption() const override
void compute_shape_function_values(const std::vector< Point< dim >> &unit_points)
MappingQGeneric(const unsigned int polynomial_degree)