Reference documentation for deal.II version 9.1.1
|
#include <deal.II/fe/fe_values.h>
Public Member Functions | |
FESubfaceValues (const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &face_quadrature, const UpdateFlags update_flags) | |
FESubfaceValues (const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &face_quadrature, const UpdateFlags update_flags) | |
template<template< int, int > class DoFHandlerType, bool level_dof_access> | |
void | reinit (const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access >> &cell, const unsigned int face_no, const unsigned int subface_no) |
void | reinit (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no) |
const FESubfaceValues< dim, spacedim > & | get_present_fe_values () const |
Public Member Functions inherited from FEFaceValuesBase< dim, spacedim > | |
FEFaceValuesBase (const unsigned int n_q_points, const unsigned int dofs_per_cell, const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature) | |
const Tensor< 1, spacedim > & | boundary_form (const unsigned int i) const |
const std::vector< Tensor< 1, spacedim > > & | get_boundary_forms () const |
unsigned int | get_face_index () const |
const Quadrature< dim - 1 > & | get_quadrature () const |
std::size_t | memory_consumption () const |
Public Member Functions inherited from FEValuesBase< dim, spacedim > | |
FEValuesBase (const unsigned int n_q_points, const unsigned int dofs_per_cell, const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe) | |
FEValuesBase & | operator= (const FEValuesBase &)=delete |
FEValuesBase (const FEValuesBase &)=delete | |
virtual | ~FEValuesBase () override |
const double & | shape_value (const unsigned int function_no, const unsigned int point_no) const |
double | shape_value_component (const unsigned int function_no, const unsigned int point_no, const unsigned int component) const |
const Tensor< 1, spacedim > & | shape_grad (const unsigned int function_no, const unsigned int quadrature_point) const |
Tensor< 1, spacedim > | shape_grad_component (const unsigned int function_no, const unsigned int point_no, const unsigned int component) const |
const Tensor< 2, spacedim > & | shape_hessian (const unsigned int function_no, const unsigned int point_no) const |
Tensor< 2, spacedim > | shape_hessian_component (const unsigned int function_no, const unsigned int point_no, const unsigned int component) const |
const Tensor< 3, spacedim > & | shape_3rd_derivative (const unsigned int function_no, const unsigned int point_no) const |
Tensor< 3, spacedim > | shape_3rd_derivative_component (const unsigned int function_no, const unsigned int point_no, const unsigned int component) const |
template<class InputVector > | |
void | get_function_values (const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const |
template<class InputVector > | |
void | get_function_values (const InputVector &fe_function, std::vector< Vector< typename InputVector::value_type >> &values) const |
template<class InputVector > | |
void | get_function_values (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< typename InputVector::value_type > &values) const |
template<class InputVector > | |
void | get_function_values (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< Vector< typename InputVector::value_type >> &values) const |
template<class InputVector > | |
void | get_function_values (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, ArrayView< std::vector< typename InputVector::value_type >> values, const bool quadrature_points_fastest) const |
template<class InputVector > | |
void | get_function_gradients (const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type >> &gradients) const |
template<class InputVector > | |
void | get_function_gradients (const InputVector &fe_function, std::vector< std::vector< Tensor< 1, spacedim, typename InputVector::value_type >>> &gradients) const |
template<class InputVector > | |
void | get_function_gradients (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< Tensor< 1, spacedim, typename InputVector::value_type >> &gradients) const |
template<class InputVector > | |
void | get_function_gradients (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, ArrayView< std::vector< Tensor< 1, spacedim, typename InputVector::value_type >>> gradients, bool quadrature_points_fastest=false) const |
template<class InputVector > | |
void | get_function_hessians (const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type >> &hessians) const |
template<class InputVector > | |
void | get_function_hessians (const InputVector &fe_function, std::vector< std::vector< Tensor< 2, spacedim, typename InputVector::value_type >>> &hessians, bool quadrature_points_fastest=false) const |
template<class InputVector > | |
void | get_function_hessians (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< Tensor< 2, spacedim, typename InputVector::value_type >> &hessians) const |
template<class InputVector > | |
void | get_function_hessians (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, ArrayView< std::vector< Tensor< 2, spacedim, typename InputVector::value_type >>> hessians, bool quadrature_points_fastest=false) const |
template<class InputVector > | |
void | get_function_laplacians (const InputVector &fe_function, std::vector< typename InputVector::value_type > &laplacians) const |
template<class InputVector > | |
void | get_function_laplacians (const InputVector &fe_function, std::vector< Vector< typename InputVector::value_type >> &laplacians) const |
template<class InputVector > | |
void | get_function_laplacians (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< typename InputVector::value_type > &laplacians) const |
template<class InputVector > | |
void | get_function_laplacians (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< Vector< typename InputVector::value_type >> &laplacians) const |
template<class InputVector > | |
void | get_function_laplacians (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< std::vector< typename InputVector::value_type >> &laplacians, bool quadrature_points_fastest=false) const |
template<class InputVector > | |
void | get_function_third_derivatives (const InputVector &fe_function, std::vector< Tensor< 3, spacedim, typename InputVector::value_type >> &third_derivatives) const |
template<class InputVector > | |
void | get_function_third_derivatives (const InputVector &fe_function, std::vector< std::vector< Tensor< 3, spacedim, typename InputVector::value_type >>> &third_derivatives, bool quadrature_points_fastest=false) const |
template<class InputVector > | |
void | get_function_third_derivatives (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< Tensor< 3, spacedim, typename InputVector::value_type >> &third_derivatives) const |
template<class InputVector > | |
void | get_function_third_derivatives (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, ArrayView< std::vector< Tensor< 3, spacedim, typename InputVector::value_type >>> third_derivatives, bool quadrature_points_fastest=false) const |
const Point< spacedim > & | quadrature_point (const unsigned int q) const |
const std::vector< Point< spacedim > > & | get_quadrature_points () const |
double | JxW (const unsigned int quadrature_point) const |
const std::vector< double > & | get_JxW_values () const |
const DerivativeForm< 1, dim, spacedim > & | jacobian (const unsigned int quadrature_point) const |
const std::vector< DerivativeForm< 1, dim, spacedim > > & | get_jacobians () const |
const DerivativeForm< 2, dim, spacedim > & | jacobian_grad (const unsigned int quadrature_point) const |
const std::vector< DerivativeForm< 2, dim, spacedim > > & | get_jacobian_grads () const |
const Tensor< 3, spacedim > & | jacobian_pushed_forward_grad (const unsigned int quadrature_point) const |
const std::vector< Tensor< 3, spacedim > > & | get_jacobian_pushed_forward_grads () const |
const DerivativeForm< 3, dim, spacedim > & | jacobian_2nd_derivative (const unsigned int quadrature_point) const |
const std::vector< DerivativeForm< 3, dim, spacedim > > & | get_jacobian_2nd_derivatives () const |
const Tensor< 4, spacedim > & | jacobian_pushed_forward_2nd_derivative (const unsigned int quadrature_point) const |
const std::vector< Tensor< 4, spacedim > > & | get_jacobian_pushed_forward_2nd_derivatives () const |
const DerivativeForm< 4, dim, spacedim > & | jacobian_3rd_derivative (const unsigned int quadrature_point) const |
const std::vector< DerivativeForm< 4, dim, spacedim > > & | get_jacobian_3rd_derivatives () const |
const Tensor< 5, spacedim > & | jacobian_pushed_forward_3rd_derivative (const unsigned int quadrature_point) const |
const std::vector< Tensor< 5, spacedim > > & | get_jacobian_pushed_forward_3rd_derivatives () const |
const DerivativeForm< 1, spacedim, dim > & | inverse_jacobian (const unsigned int quadrature_point) const |
const std::vector< DerivativeForm< 1, spacedim, dim > > & | get_inverse_jacobians () const |
const Tensor< 1, spacedim > & | normal_vector (const unsigned int i) const |
const std::vector< Tensor< 1, spacedim > > & | get_all_normal_vectors () const |
const std::vector< Tensor< 1, spacedim > > & | get_normal_vectors () const |
const FEValuesViews::Scalar< dim, spacedim > & | operator[] (const FEValuesExtractors::Scalar &scalar) const |
const FEValuesViews::Vector< dim, spacedim > & | operator[] (const FEValuesExtractors::Vector &vector) const |
const FEValuesViews::SymmetricTensor< 2, dim, spacedim > & | operator[] (const FEValuesExtractors::SymmetricTensor< 2 > &tensor) const |
const FEValuesViews::Tensor< 2, dim, spacedim > & | operator[] (const FEValuesExtractors::Tensor< 2 > &tensor) const |
const Mapping< dim, spacedim > & | get_mapping () const |
const FiniteElement< dim, spacedim > & | get_fe () const |
UpdateFlags | get_update_flags () const |
const Triangulation< dim, spacedim >::cell_iterator | get_cell () const |
CellSimilarity::Similarity | get_cell_similarity () const |
std::size_t | memory_consumption () const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Public Member Functions | |
static ::ExceptionBase & | ExcReinitCalledWithBoundaryFace () |
static ::ExceptionBase & | ExcFaceHasNoSubfaces () |
Static Public Member Functions inherited from FEValuesBase< dim, spacedim > | |
static ::ExceptionBase & | ExcAccessToUninitializedField (std::string arg1) |
static ::ExceptionBase & | ExcFEDontMatch () |
static ::ExceptionBase & | ExcShapeFunctionNotPrimitive (int arg1) |
static ::ExceptionBase & | ExcFENotPrimitive () |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Static Public Attributes | |
static const unsigned int | dimension = dim |
static const unsigned int | space_dimension = spacedim |
static const unsigned int | integral_dimension = dim - 1 |
Static Public Attributes inherited from FEFaceValuesBase< dim, spacedim > | |
static const unsigned int | integral_dimension = dim - 1 |
Static Public Attributes inherited from FEValuesBase< dim, spacedim > | |
static const unsigned int | dimension = dim |
static const unsigned int | space_dimension = spacedim |
Private Member Functions | |
void | initialize (const UpdateFlags update_flags) |
void | do_reinit (const unsigned int face_no, const unsigned int subface_no) |
Additional Inherited Members | |
Public Attributes inherited from FEValuesBase< dim, spacedim > | |
const unsigned int | n_quadrature_points |
const unsigned int | dofs_per_cell |
Protected Member Functions inherited from FEValuesBase< dim, spacedim > | |
void | invalidate_present_cell () |
void | maybe_invalidate_previous_present_cell (const typename Triangulation< dim, spacedim >::cell_iterator &cell) |
UpdateFlags | compute_update_flags (const UpdateFlags update_flags) const |
void | check_cell_similarity (const typename Triangulation< dim, spacedim >::cell_iterator &cell) |
Protected Attributes inherited from FEFaceValuesBase< dim, spacedim > | |
unsigned int | present_face_index |
const Quadrature< dim - 1 > | quadrature |
Protected Attributes inherited from FEValuesBase< dim, spacedim > | |
std::unique_ptr< const CellIteratorBase > | present_cell |
boost::signals2::connection | tria_listener_refinement |
boost::signals2::connection | tria_listener_mesh_transform |
const SmartPointer< const Mapping< dim, spacedim >, FEValuesBase< dim, spacedim > > | mapping |
std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > | mapping_data |
::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > | mapping_output |
const SmartPointer< const FiniteElement< dim, spacedim >, FEValuesBase< dim, spacedim > > | fe |
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > | fe_data |
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > | finite_element_output |
UpdateFlags | update_flags |
CellSimilarity::Similarity | cell_similarity |
Finite element evaluated in quadrature points on a face.
This class adds the functionality of FEFaceValuesBase to FEValues; see there for more documentation.
This class is used for faces lying on a refinement edge. In this case, the neighboring cell is refined. To be able to compute differences between interior and exterior function values, the refinement of the neighboring cell must be simulated on this cell. This is achieved by applying a quadrature rule that simulates the refinement. The resulting data fields are split up to reflect the refinement structure of the neighbor: a subface number corresponds to the number of the child of the neighboring face.
FESubfaceValues< dim, spacedim >::FESubfaceValues | ( | const Mapping< dim, spacedim > & | mapping, |
const FiniteElement< dim, spacedim > & | fe, | ||
const Quadrature< dim - 1 > & | face_quadrature, | ||
const UpdateFlags | update_flags | ||
) |
Constructor. Gets cell independent data from mapping and finite element objects, matching the quadrature rule and update flags.
Definition at line 4855 of file fe_values.cc.
FESubfaceValues< dim, spacedim >::FESubfaceValues | ( | const FiniteElement< dim, spacedim > & | fe, |
const Quadrature< dim - 1 > & | face_quadrature, | ||
const UpdateFlags | update_flags | ||
) |
Constructor. This constructor is equivalent to the other one except that it makes the object use a \(Q_1\) mapping (i.e., an object of type MappingQGeneric(1)) implicitly.
Definition at line 4873 of file fe_values.cc.
void FESubfaceValues< dim, spacedim >::reinit | ( | const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, lda >> & | cell, |
const unsigned int | face_no, | ||
const unsigned int | subface_no | ||
) |
Reinitialize the gradients, Jacobi determinants, etc for the given cell of type "iterator into a DoFHandler object", and the finite element associated with this object. It is assumed that the finite element used by the given cell is also the one used by this FESubfaceValues object.
Definition at line 4940 of file fe_values.cc.
void FESubfaceValues< dim, spacedim >::reinit | ( | const typename Triangulation< dim, spacedim >::cell_iterator & | cell, |
const unsigned int | face_no, | ||
const unsigned int | subface_no | ||
) |
Reinitialize the gradients, Jacobi determinants, etc for the given subface on given cell of type "iterator into a Triangulation object", and the given finite element. Since iterators into triangulation alone only convey information about the geometry of a cell, but not about degrees of freedom possibly associated with this cell, you will not be able to call some functions of this class if they need information about degrees of freedom. These functions are, above all, the get_function_value/gradients/hessians/third_derivatives
functions. If you want to call these functions, you have to call the reinit
variants that take iterators into DoFHandler or other DoF handler type objects.
Definition at line 4988 of file fe_values.cc.
const FESubfaceValues<dim, spacedim>& FESubfaceValues< dim, spacedim >::get_present_fe_values | ( | ) | const |
Return a reference to this very object.
Though it seems that it is not very useful, this function is there to provide capability to the hp::FEValues class, in which case it provides the FEValues object for the present cell (remember that for hp finite elements, the actual FE object used may change from cell to cell, so we also need different FEValues objects for different cells; once you reinitialize the hp::FEValues object for a specific cell, it retrieves the FEValues object for the FE on that cell and returns it through a function of the same name as this one; this function here therefore only provides the same interface so that one can templatize on FEValues and hp::FEValues).
|
private |
Do work common to the two constructors.
Definition at line 4891 of file fe_values.cc.
|
private |
The reinit() functions do only that part of the work that requires knowledge of the type of iterator. After setting present_cell(), they pass on to this function, which does the real work, and which is independent of the actual type of the cell iterator.
Definition at line 5013 of file fe_values.cc.
|
static |
Dimension in which this object operates.
Definition at line 3835 of file fe_values.h.
|
static |
Dimension of the space in which this object operates.
Definition at line 3840 of file fe_values.h.
|
static |
Dimension of the object over which we integrate. For the present class, this is equal to dim-1
.
Definition at line 3846 of file fe_values.h.