Reference documentation for deal.II version 9.0.0
|
#include <deal.II/fe/fe_values.h>
Public Member Functions | |
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 () | |
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 VectorSlice< const std::vector< 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 VectorSlice< const std::vector< 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 VectorSlice< const std::vector< types::global_dof_index > > &indices, VectorSlice< std::vector< 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 VectorSlice< const std::vector< 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 VectorSlice< const std::vector< types::global_dof_index > > &indices, VectorSlice< std::vector< 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 VectorSlice< const std::vector< 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 VectorSlice< const std::vector< types::global_dof_index > > &indices, VectorSlice< std::vector< 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 VectorSlice< const std::vector< 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 VectorSlice< const std::vector< 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 VectorSlice< const std::vector< 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 VectorSlice< const std::vector< 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 VectorSlice< const std::vector< types::global_dof_index > > &indices, VectorSlice< std::vector< 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 (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Public Attributes | |
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 |
Protected Attributes | |
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 |
Additional Inherited Members | |
Static Public Member Functions inherited from FEValuesBase< dim, spacedim > | |
static ::ExceptionBase & | ExcAccessToUninitializedField (char *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) |
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) |
Extend the interface of FEValuesBase to values that only make sense when evaluating something on the surface of a cell. All the data that is available in the interior of cells is also available here.
See FEValuesBase
Definition at line 3252 of file fe_values.h.
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 | ||
) |
Constructor. Call the constructor of the base class and set up the arrays of this class with the right sizes. Actually filling these arrays is a duty of the derived class's constructors.
n_faces_or_subfaces
is the number of faces or subfaces that this object is to store. The actual number depends on the derived class, for FEFaceValues it is 2*dim
, while for the FESubfaceValues class it is 2*dim*(1<<(dim-1))
, i.e. the number of faces times the number of subfaces per face.
Definition at line 4176 of file fe_values.cc.
const Tensor<1,spacedim>& FEFaceValuesBase< dim, spacedim >::boundary_form | ( | const unsigned int | i | ) | const |
Boundary form of the transformation of the cell at the i
th quadrature point. See GlossBoundaryForm.
update_boundary_forms
flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information. const std::vector< Tensor< 1, spacedim > > & FEFaceValuesBase< dim, spacedim >::get_boundary_forms | ( | ) | const |
Return the list of outward normal vectors times the Jacobian of the surface mapping.
update_boundary_forms
flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information. Definition at line 4196 of file fe_values.cc.
unsigned int FEFaceValuesBase< dim, spacedim >::get_face_index | ( | ) | const |
Return the index of the face selected the last time the reinit() function was called.
const Quadrature<dim-1>& FEFaceValuesBase< dim, spacedim >::get_quadrature | ( | ) | const |
Return a reference to the copy of the quadrature formula stored by this object.
std::size_t FEFaceValuesBase< dim, spacedim >::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object.
Definition at line 4207 of file fe_values.cc.
|
static |
Dimension of the object over which we integrate. For the present class, this is equal to dim-1
.
Definition at line 3259 of file fe_values.h.
|
protected |
Index of the face selected the last time the reinit() function was called.
Definition at line 3320 of file fe_values.h.
|
protected |
Store a copy of the quadrature formula here.
Definition at line 3325 of file fe_values.h.