Reference documentation for deal.II version 9.4.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 Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const hp::QCollection< dim - 1 > &face_quadrature, const UpdateFlags update_flags) | |
FESubfaceValues (const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &face_quadrature, const UpdateFlags update_flags) | |
FESubfaceValues (const FiniteElement< dim, spacedim > &fe, const hp::QCollection< dim - 1 > &face_quadrature, const UpdateFlags update_flags) | |
template<bool level_dof_access> | |
void | reinit (const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell, const unsigned int face_no, const unsigned int subface_no) |
template<bool level_dof_access> | |
void | reinit (const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell, const typename Triangulation< dim, spacedim >::face_iterator &face, const typename Triangulation< dim, spacedim >::face_iterator &subface) |
void | reinit (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no) |
void | reinit (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::face_iterator &face, const typename Triangulation< dim, spacedim >::face_iterator &subface) |
const FESubfaceValues< dim, spacedim > & | get_present_fe_values () const |
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_number () const |
unsigned int | get_face_index () const |
const Quadrature< dim - 1 > & | get_quadrature () const |
std::size_t | memory_consumption () const |
Access to shape function values | |
These fields are filled by the finite element. | |
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 |
Access to values of global finite element fields | |
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 |
Access to derivatives of global finite element fields | |
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, const bool quadrature_points_fastest=false) const |
Access to second derivatives | |
Hessian matrices and Laplacians of global finite element fields | |
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, const 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, const 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, const bool quadrature_points_fastest=false) const |
Access to third derivatives of global finite element fields | |
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, const 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, const bool quadrature_points_fastest=false) const |
Cell degrees of freedom | |
std_cxx20::ranges::iota_view< unsigned int, unsigned int > | dof_indices () const |
std_cxx20::ranges::iota_view< unsigned int, unsigned int > | dof_indices_starting_at (const unsigned int start_dof_index) const |
std_cxx20::ranges::iota_view< unsigned int, unsigned int > | dof_indices_ending_at (const unsigned int end_dof_index) const |
Geometry of the cell | |
std_cxx20::ranges::iota_view< unsigned int, unsigned int > | quadrature_point_indices () 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_normal_vectors () const |
Extractors Methods to extract individual components | |
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 |
Static Public Member Functions | |
static ::ExceptionBase & | ExcReinitCalledWithBoundaryFace () |
static ::ExceptionBase & | ExcFaceHasNoSubfaces () |
Public Attributes | |
const unsigned int | n_quadrature_points |
const unsigned int | max_n_quadrature_points |
const unsigned int | dofs_per_cell |
Static Public Attributes | |
static constexpr unsigned int | dimension = dim |
static constexpr unsigned int | space_dimension = spacedim |
static constexpr unsigned int | integral_dimension = dim - 1 |
Protected Attributes | |
unsigned int | present_face_no |
unsigned int | present_face_index |
const hp::QCollection< dim - 1 > | quadrature |
Private Member Functions | |
void | initialize (const UpdateFlags update_flags) |
void | do_reinit (const unsigned int face_no, const unsigned int subface_no) |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
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 ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
void | check_no_subscribers () const noexcept |
std::atomic< unsigned int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
static std::mutex | mutex |
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.
Definition at line 4469 of file fe_values.h.
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.
FESubfaceValues< dim, spacedim >::FESubfaceValues | ( | const Mapping< dim, spacedim > & | mapping, |
const FiniteElement< dim, spacedim > & | fe, | ||
const hp::QCollection< dim - 1 > & | face_quadrature, | ||
const UpdateFlags | update_flags | ||
) |
Like the function above, but taking a collection of quadrature rules.
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 MappingQ(1)) implicitly.
FESubfaceValues< dim, spacedim >::FESubfaceValues | ( | const FiniteElement< dim, spacedim > & | fe, |
const hp::QCollection< dim - 1 > & | face_quadrature, | ||
const UpdateFlags | update_flags | ||
) |
Like the function above, but taking a collection of quadrature rules.
void FESubfaceValues< dim, spacedim >::reinit | ( | const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > & | 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.
void FESubfaceValues< dim, spacedim >::reinit | ( | const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > & | cell, |
const typename Triangulation< dim, spacedim >::face_iterator & | face, | ||
const typename Triangulation< dim, spacedim >::face_iterator & | subface | ||
) |
Alternative reinitialization function that takes, as arguments, iterators to the face and subface instead of their numbers.
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 a given cell of type "iterator into a Triangulation object", and the given finite element. Since iterators into a 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.
void FESubfaceValues< dim, spacedim >::reinit | ( | const typename Triangulation< dim, spacedim >::cell_iterator & | cell, |
const typename Triangulation< dim, spacedim >::face_iterator & | face, | ||
const typename Triangulation< dim, spacedim >::face_iterator & | subface | ||
) |
Reinitialize the gradients, Jacobi determinants, etc for the given subface on a given cell of type "iterator into a Triangulation object", and the given finite element. Since iterators into a 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.
This does the same thing as the previous function but takes iterators instead of numbers as arguments.
face
and subface
must correspond to a face (and a subface of that face) of cell
. 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.
|
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.
|
inherited |
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.
|
inherited |
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.
|
inherited |
Return the number of the face selected the last time the reinit() function was called.
|
inherited |
Return the index of the face selected the last time the reinit() function was called.
|
inherited |
Return a reference to the copy of the quadrature formula stored by this object.
|
inherited |
Determine an estimate for the memory consumption (in bytes) of this object.
|
inherited |
Value of a shape function at a quadrature point on the cell, face or subface selected the last time the reinit
function of the derived class was called.
If the shape function is vector-valued, then this returns the only non- zero component. If the shape function has more than one non-zero component (i.e. it is not primitive), then throw an exception of type ExcShapeFunctionNotPrimitive. In that case, use the shape_value_component() function.
function_no | Number of the shape function to be evaluated. Note that this number runs from zero to dofs_per_cell, even in the case of an FEFaceValues or FESubfaceValues object. |
point_no | Number of the quadrature point at which function is to be evaluated |
update_values
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.
|
inherited |
Compute one vector component of the value of a shape function at a quadrature point. If the finite element is scalar, then only component zero is allowed and the return value equals that of the shape_value() function. If the finite element is vector valued but all shape functions are primitive (i.e. they are non-zero in only one component), then the value returned by shape_value() equals that of this function for exactly one component. This function is therefore only of greater interest if the shape function is not primitive, but then it is necessary since the other function cannot be used.
function_no | Number of the shape function to be evaluated. |
point_no | Number of the quadrature point at which function is to be evaluated. |
component | vector component to be evaluated. |
update_values
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.
|
inherited |
Compute the gradient of the function_no
th shape function at the quadrature_point
th quadrature point with respect to real cell coordinates. If you want to get the derivative in one of the coordinate directions, use the appropriate function of the Tensor class to extract one component of the Tensor returned by this function. Since only a reference to the gradient's value is returned, there should be no major performance drawback.
If the shape function is vector-valued, then this returns the only non- zero component. If the shape function has more than one non-zero component (i.e. it is not primitive), then it will throw an exception of type ExcShapeFunctionNotPrimitive. In that case, use the shape_grad_component() function.
The same holds for the arguments of this function as for the shape_value() function.
function_no | Number of the shape function to be evaluated. |
quadrature_point | Number of the quadrature point at which function is to be evaluated. |
update_gradients
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.
|
inherited |
Return one vector component of the gradient of a shape function at a quadrature point. If the finite element is scalar, then only component zero is allowed and the return value equals that of the shape_grad() function. If the finite element is vector valued but all shape functions are primitive (i.e. they are non-zero in only one component), then the value returned by shape_grad() equals that of this function for exactly one component. This function is therefore only of greater interest if the shape function is not primitive, but then it is necessary since the other function cannot be used.
The same holds for the arguments of this function as for the shape_value_component() function.
update_gradients
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.
|
inherited |
Second derivatives of the function_no
th shape function at the point_no
th quadrature point with respect to real cell coordinates. If you want to get the derivatives in one of the coordinate directions, use the appropriate function of the Tensor class to extract one component. Since only a reference to the hessian values is returned, there should be no major performance drawback.
If the shape function is vector-valued, then this returns the only non- zero component. If the shape function has more than one non-zero component (i.e. it is not primitive), then throw an exception of type ExcShapeFunctionNotPrimitive. In that case, use the shape_hessian_component() function.
The same holds for the arguments of this function as for the shape_value() function.
update_hessians
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.
|
inherited |
Return one vector component of the hessian of a shape function at a quadrature point. If the finite element is scalar, then only component zero is allowed and the return value equals that of the shape_hessian() function. If the finite element is vector valued but all shape functions are primitive (i.e. they are non-zero in only one component), then the value returned by shape_hessian() equals that of this function for exactly one component. This function is therefore only of greater interest if the shape function is not primitive, but then it is necessary since the other function cannot be used.
The same holds for the arguments of this function as for the shape_value_component() function.
update_hessians
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.
|
inherited |
Third derivatives of the function_no
th shape function at the point_no
th quadrature point with respect to real cell coordinates. If you want to get the 3rd derivatives in one of the coordinate directions, use the appropriate function of the Tensor class to extract one component. Since only a reference to the 3rd derivative values is returned, there should be no major performance drawback.
If the shape function is vector-valued, then this returns the only non- zero component. If the shape function has more than one non-zero component (i.e. it is not primitive), then throw an exception of type ExcShapeFunctionNotPrimitive. In that case, use the shape_3rdderivative_component() function.
The same holds for the arguments of this function as for the shape_value() function.
update_3rd_derivatives
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.
|
inherited |
Return one vector component of the third derivative of a shape function at a quadrature point. If the finite element is scalar, then only component zero is allowed and the return value equals that of the shape_3rdderivative() function. If the finite element is vector valued but all shape functions are primitive (i.e. they are non-zero in only one component), then the value returned by shape_3rdderivative() equals that of this function for exactly one component. This function is therefore only of greater interest if the shape function is not primitive, but then it is necessary since the other function cannot be used.
The same holds for the arguments of this function as for the shape_value_component() function.
update_3rd_derivatives
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.
|
inherited |
Return the values of a finite element function restricted to the current cell, face or subface selected the last time the reinit
function of the derived class was called, at the quadrature points.
If the present cell is not active then values are interpolated to the current cell and point values are computed from that.
This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. To get values of multi- component elements, there is another get_function_values() below, returning a vector of vectors of results.
[in] | fe_function | A vector of values that describes (globally) the finite element function that this function should evaluate at the quadrature points of the current cell. |
[out] | values | The values of the function specified by fe_function at the quadrature points of the current cell. The object is assume to already have the correct size. The data type stored by this output vector must be what you get when you multiply the values of shape function times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument). This happens to be equal to the type of the elements of the solution vector. |
values[q]
will contain the value of the field described by fe_function at the \(q\)th quadrature point.update_values
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 3356 of file fe_values.cc.
|
inherited |
This function does the same as the other get_function_values(), but applied to multi-component (vector-valued) elements. The meaning of the arguments is as explained there.
values[q]
is a vector of values of the field described by fe_function at the \(q\)th quadrature point. The size of the vector accessed by values[q]
equals the number of components of the finite element, i.e. values[q](c)
returns the value of the \(c\)th vector component at the \(q\)th quadrature point.update_values
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 3404 of file fe_values.cc.
|
inherited |
Generate function values from an arbitrary vector. This function does in essence the same as the first function of this name above, except that it does not make the assumption that the input vector corresponds to a DoFHandler that describes the unknowns of a finite element field (and for which we would then assume that fe_function.size() == dof_handler.n_dofs()
). Rather, the nodal values corresponding to the current cell are elements of an otherwise arbitrary vector, and these elements are indexed by the second argument to this function. What the rest of the fe_function
input argument corresponds to is of no consequence to this function.
Given this, the function above corresponds to passing fe_function
as first argument to the current function, and using the local_dof_indices
array that results from the following call as second argument to the current function:
(See DoFCellAccessor::get_dof_indices() for more information.)
Likewise, the function above is equivalent to calling
and then calling the current function with local_dof_values
as first argument, and an array with indices {0,...,fe.dofs_per_cell-1}
as second argument.
The point of the current function is that one sometimes wants to evaluate finite element functions at quadrature points with nodal values that are not stored in a global vector – for example, one could modify these local values first, such as by applying a limiter or by ensuring that all nodal values are positive, before evaluating the finite element field that corresponds to these local values on the current cell. Another application is where one wants to postprocess the solution on a cell into a different finite element space on every cell, without actually creating a corresponding DoFHandler – in that case, all one would compute is a local representation of that postprocessed function, characterized by its nodal values; this function then allows the evaluation of that representation at quadrature points.
[in] | fe_function | A vector of nodal values. This vector can have an arbitrary size, as long as all elements index by indices can actually be accessed. |
[in] | indices | A vector of indices into fe_function . This vector must have length equal to the number of degrees of freedom on the current cell, and must identify elements in fe_function in the order in which degrees of freedom are indexed on the reference cell. |
[out] | values | A vector of values of the given finite element field, at the quadrature points on the current object. |
update_values
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 3380 of file fe_values.cc.
|
inherited |
Generate vector function values from an arbitrary vector.
This function corresponds to the previous one, just for the vector-valued case.
update_values
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 3431 of file fe_values.cc.
|
inherited |
Generate vector function values from an arbitrary vector. This function is similar to the previous one, but the indices
vector may also be a multiple of the number of dofs per cell. Then, the vectors in value
should allow for the same multiple of the components of the finite element.
Depending on the value of the last argument, the outer vector of values
has either the length of the quadrature rule (quadrature_points_fastest == false
) or the length of components to be filled quadrature_points_fastest == true
. If p
is the current quadrature point number and i
is the vector component of the solution desired, the access to values
is values[p][i]
if quadrature_points_fastest == false
, and values[i][p]
otherwise.
Since this function allows for fairly general combinations of argument sizes, be aware that the checks on the arguments may not detect errors.
update_values
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 3462 of file fe_values.cc.
|
inherited |
Compute the gradients of a finite element at the quadrature points of a cell. This function is the equivalent of the corresponding get_function_values() function (see there for more information) but evaluates the finite element field's gradient instead of its value.
This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. There is a corresponding function of the same name for vector-valued finite elements.
[in] | fe_function | A vector of values that describes (globally) the finite element function that this function should evaluate at the quadrature points of the current cell. |
[out] | gradients | The gradients of the function specified by fe_function at the quadrature points of the current cell. The gradients are computed in real space (as opposed to on the unit cell). The object is assume to already have the correct size. The data type stored by this output vector must be what you get when you multiply the gradients of shape function times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument). |
gradients[q]
will contain the gradient of the field described by fe_function at the \(q\)th quadrature point. gradients[q][d]
represents the derivative in coordinate direction \(d\) at quadrature point \(q\).update_gradients
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 3495 of file fe_values.cc.
|
inherited |
This function does the same as the other get_function_gradients(), but applied to multi-component (vector-valued) elements. The meaning of the arguments is as explained there.
gradients[q]
is a vector of gradients of the field described by fe_function at the \(q\)th quadrature point. The size of the vector accessed by gradients[q]
equals the number of components of the finite element, i.e. gradients[q][c]
returns the gradient of the \(c\)th vector component at the \(q\)th quadrature point. Consequently, gradients[q][c][d]
is the derivative in coordinate direction \(d\) of the \(c\)th vector component of the vector field at quadrature point \(q\) of the current cell.update_gradients
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 3545 of file fe_values.cc.
|
inherited |
This function relates to the first of the get_function_gradients() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.
update_gradients
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 3520 of file fe_values.cc.
|
inherited |
This function relates to the first of the get_function_gradients() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.
update_gradients
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 3573 of file fe_values.cc.
|
inherited |
Compute the tensor of second derivatives of a finite element at the quadrature points of a cell. This function is the equivalent of the corresponding get_function_values() function (see there for more information) but evaluates the finite element field's second derivatives instead of its value.
This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. There is a corresponding function of the same name for vector-valued finite elements.
[in] | fe_function | A vector of values that describes (globally) the finite element function that this function should evaluate at the quadrature points of the current cell. |
[out] | hessians | The Hessians of the function specified by fe_function at the quadrature points of the current cell. The Hessians are computed in real space (as opposed to on the unit cell). The object is assume to already have the correct size. The data type stored by this output vector must be what you get when you multiply the Hessians of shape function times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument). |
hessians[q]
will contain the Hessian of the field described by fe_function at the \(q\)th quadrature point. hessians[q][i][j]
represents the \((i,j)\)th component of the matrix of second derivatives at quadrature point \(q\).update_hessians
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 3606 of file fe_values.cc.
|
inherited |
This function does the same as the other get_function_hessians(), but applied to multi-component (vector-valued) elements. The meaning of the arguments is as explained there.
hessians[q]
is a vector of Hessians of the field described by fe_function at the \(q\)th quadrature point. The size of the vector accessed by hessians[q]
equals the number of components of the finite element, i.e. hessians[q][c]
returns the Hessian of the \(c\)th vector component at the \(q\)th quadrature point. Consequently, hessians[q][c][i][j]
is the \((i,j)\)th component of the matrix of second derivatives of the \(c\)th vector component of the vector field at quadrature point \(q\) of the current cell.update_hessians
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 3656 of file fe_values.cc.
|
inherited |
This function relates to the first of the get_function_hessians() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.
update_hessians
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 3631 of file fe_values.cc.
|
inherited |
This function relates to the first of the get_function_hessians() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.
update_hessians
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 3686 of file fe_values.cc.
|
inherited |
Compute the (scalar) Laplacian (i.e. the trace of the tensor of second derivatives) of a finite element at the quadrature points of a cell. This function is the equivalent of the corresponding get_function_values() function (see there for more information) but evaluates the finite element field's second derivatives instead of its value.
This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. There is a corresponding function of the same name for vector-valued finite elements.
[in] | fe_function | A vector of values that describes (globally) the finite element function that this function should evaluate at the quadrature points of the current cell. |
[out] | laplacians | The Laplacians of the function specified by fe_function at the quadrature points of the current cell. The Laplacians are computed in real space (as opposed to on the unit cell). The object is assume to already have the correct size. The data type stored by this output vector must be what you get when you multiply the Laplacians of shape function times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument). This happens to be equal to the type of the elements of the input vector. |
laplacians[q]
will contain the Laplacian of the field described by fe_function at the \(q\)th quadrature point.laplacians[q]=trace(hessians[q])
, where hessians
would be the output of the get_function_hessians() function.update_hessians
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 3717 of file fe_values.cc.
|
inherited |
This function does the same as the other get_function_laplacians(), but applied to multi-component (vector-valued) elements. The meaning of the arguments is as explained there.
laplacians[q]
is a vector of Laplacians of the field described by fe_function at the \(q\)th quadrature point. The size of the vector accessed by laplacians[q]
equals the number of components of the finite element, i.e. laplacians[q][c]
returns the Laplacian of the \(c\)th vector component at the \(q\)th quadrature point.laplacians[q][c]=trace(hessians[q][c])
, where hessians
would be the output of the get_function_hessians() function.update_hessians
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 3765 of file fe_values.cc.
|
inherited |
This function relates to the first of the get_function_laplacians() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.
update_hessians
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 3741 of file fe_values.cc.
|
inherited |
This function relates to the first of the get_function_laplacians() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.
update_hessians
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 3791 of file fe_values.cc.
|
inherited |
This function relates to the first of the get_function_laplacians() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.
update_hessians
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 3822 of file fe_values.cc.
|
inherited |
Compute the tensor of third derivatives of a finite element at the quadrature points of a cell. This function is the equivalent of the corresponding get_function_values() function (see there for more information) but evaluates the finite element field's third derivatives instead of its value.
This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. There is a corresponding function of the same name for vector-valued finite elements.
[in] | fe_function | A vector of values that describes (globally) the finite element function that this function should evaluate at the quadrature points of the current cell. |
[out] | third_derivatives | The third derivatives of the function specified by fe_function at the quadrature points of the current cell. The third derivatives are computed in real space (as opposed to on the unit cell). The object is assumed to already have the correct size. The data type stored by this output vector must be what you get when you multiply the third derivatives of shape function times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument). |
third_derivatives[q]
will contain the third derivatives of the field described by fe_function at the \(q\)th quadrature point. third_derivatives[q][i][j][k]
represents the \((i,j,k)\)th component of the 3rd order tensor of third derivatives at quadrature point \(q\).update_3rd_derivatives
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 3852 of file fe_values.cc.
|
inherited |
This function does the same as the other get_function_third_derivatives(), but applied to multi-component (vector- valued) elements. The meaning of the arguments is as explained there.
third_derivatives[q]
is a vector of third derivatives of the field described by fe_function at the \(q\)th quadrature point. The size of the vector accessed by third_derivatives[q]
equals the number of components of the finite element, i.e. third_derivatives[q][c]
returns the third derivative of the \(c\)th vector component at the \(q\)th quadrature point. Consequently, third_derivatives[q][c][i][j][k]
is the \((i,j,k)\)th component of the tensor of third derivatives of the \(c\)th vector component of the vector field at quadrature point \(q\) of the current cell.update_3rd_derivatives
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 3904 of file fe_values.cc.
|
inherited |
This function relates to the first of the get_function_third_derivatives() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.
update_3rd_derivatives
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 3878 of file fe_values.cc.
|
inherited |
This function relates to the first of the get_function_third_derivatives() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.
update_3rd_derivatives
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 3934 of file fe_values.cc.
|
inherited |
Return an object that can be thought of as an array containing all indices from zero (inclusive) to dofs_per_cell
(exclusive). This allows one to write code using range-based for
loops of the following kind:
Here, we are looping over all degrees of freedom on all cells, with i
and j
taking on all valid indices for cell degrees of freedom, as defined by the finite element passed to fe_values
.
|
inherited |
Return an object that can be thought of as an array containing all indices from start_dof_index
(inclusive) to dofs_per_cell
(exclusive). This allows one to write code using range-based for
loops of the following kind:
Here, we are looping over all local degrees of freedom on all cells, with i
taking on all valid indices for cell degrees of freedom, as defined by the finite element passed to fe_values
, and j
taking on a specified subset of i
's range, starting at i
itself and ending at the number of cell degrees of freedom. In this way, we can construct the upper half and the diagonal of a stiffness matrix contribution (assuming it is symmetric, and that only one half of it needs to be computed), for example.
start_dof_index
is equal to the number of DoFs in the cell, then the returned index range is empty.
|
inherited |
Return an object that can be thought of as an array containing all indices from zero (inclusive) to end_dof_index
(inclusive). This allows one to write code using range-based for
loops of the following kind:
Here, we are looping over all local degrees of freedom on all cells, with i
taking on all valid indices for cell degrees of freedom, as defined by the finite element passed to fe_values
, and j
taking on a specified subset of i
's range, starting at zero and ending at i
itself. In this way, we can construct the lower half and the diagonal of a stiffness matrix contribution (assuming it is symmetric, and that only one half of it needs to be computed), for example.
end_dof_index
is equal to zero, then the returned index range is empty.
|
inherited |
Return an object that can be thought of as an array containing all indices from zero to n_quadrature_points
. This allows to write code using range-based for
loops of the following kind:
Here, we are looping over all quadrature points on all cells, with q_point
taking on all valid indices for quadrature points, as defined by the quadrature rule passed to fe_values
.
|
inherited |
Position of the q
th quadrature point in real space.
update_quadrature_points
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.
|
inherited |
Return a reference to the vector of quadrature points in real space.
update_quadrature_points
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.
|
inherited |
Mapped quadrature weight. If this object refers to a volume evaluation (i.e. the derived class is of type FEValues), then this is the Jacobi determinant times the weight of the *i
th unit quadrature point.
For surface evaluations (i.e. classes FEFaceValues or FESubfaceValues), it is the mapped surface element times the weight of the quadrature point.
You can think of the quantity returned by this function as the volume or surface element \(dx, ds\) in the integral that we implement here by quadrature.
update_JxW_values
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.
|
inherited |
Return a reference to the array holding the values returned by JxW().
|
inherited |
Return the Jacobian of the transformation at the specified quadrature point, i.e. \(J_{ij}=dx_i/d\hat x_j\)
update_jacobians
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.
|
inherited |
Return a reference to the array holding the values returned by jacobian().
update_jacobians
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.
|
inherited |
Return the second derivative of the transformation from unit to real cell, i.e. the first derivative of the Jacobian, at the specified quadrature point, i.e. \(G_{ijk}=dJ_{jk}/d\hat x_i\).
update_jacobian_grads
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.
|
inherited |
Return a reference to the array holding the values returned by jacobian_grads().
update_jacobian_grads
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.
|
inherited |
Return the second derivative of the transformation from unit to real cell, i.e. the first derivative of the Jacobian, at the specified quadrature point, pushed forward to the real cell coordinates, i.e. \(G_{ijk}=dJ_{iJ}/d\hat x_K (J_{jJ})^{-1} (J_{kK})^{-1}\).
update_jacobian_pushed_forward_grads
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.
|
inherited |
Return a reference to the array holding the values returned by jacobian_pushed_forward_grads().
update_jacobian_pushed_forward_grads
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.
|
inherited |
Return the third derivative of the transformation from unit to real cell, i.e. the second derivative of the Jacobian, at the specified quadrature point, i.e. \(G_{ijkl}=\frac{d^2J_{ij}}{d\hat x_k d\hat x_l}\).
update_jacobian_2nd_derivatives
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.
|
inherited |
Return a reference to the array holding the values returned by jacobian_2nd_derivatives().
update_jacobian_2nd_derivatives
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.
|
inherited |
Return the third derivative of the transformation from unit to real cell, i.e. the second derivative of the Jacobian, at the specified quadrature point, pushed forward to the real cell coordinates, i.e. \(G_{ijkl}=\frac{d^2J_{iJ}}{d\hat x_K d\hat x_L} (J_{jJ})^{-1} (J_{kK})^{-1}(J_{lL})^{-1}\).
update_jacobian_pushed_forward_2nd_derivatives
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.
|
inherited |
Return a reference to the array holding the values returned by jacobian_pushed_forward_2nd_derivatives().
update_jacobian_pushed_forward_2nd_derivatives
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.
|
inherited |
Return the fourth derivative of the transformation from unit to real cell, i.e. the third derivative of the Jacobian, at the specified quadrature point, i.e. \(G_{ijklm}=\frac{d^2J_{ij}}{d\hat x_k d\hat x_l d\hat x_m}\).
update_jacobian_3rd_derivatives
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.
|
inherited |
Return a reference to the array holding the values returned by jacobian_3rd_derivatives().
update_jacobian_3rd_derivatives
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.
|
inherited |
Return the fourth derivative of the transformation from unit to real cell, i.e. the third derivative of the Jacobian, at the specified quadrature point, pushed forward to the real cell coordinates, i.e. \(G_{ijklm}=\frac{d^3J_{iJ}}{d\hat x_K d\hat x_L d\hat x_M} (J_{jJ})^{-1} (J_{kK})^{-1} (J_{lL})^{-1} (J_{mM})^{-1}\).
update_jacobian_pushed_forward_3rd_derivatives
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.
|
inherited |
Return a reference to the array holding the values returned by jacobian_pushed_forward_3rd_derivatives().
update_jacobian_pushed_forward_2nd_derivatives
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.
|
inherited |
Return the inverse Jacobian of the transformation at the specified quadrature point, i.e. \(J_{ij}=d\hat x_i/dx_j\)
update_inverse_jacobians
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.
|
inherited |
Return a reference to the array holding the values returned by inverse_jacobian().
update_inverse_jacobians
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.
|
inherited |
Return the normal vector at a quadrature point. If you call this function for a face (i.e., when using a FEFaceValues or FESubfaceValues object), then this function returns the outward normal vector to the cell at the i
th quadrature point of the face.
In contrast, if you call this function for a cell of codimension one (i.e., when using a FEValues<dim,spacedim>
object with spacedim>dim
), then this function returns the normal vector to the cell – in other words, an approximation to the normal vector to the manifold in which the triangulation is embedded. There are of course two normal directions to a manifold in that case, and this function returns the "up" direction as induced by the numbering of the vertices.
The length of the vector is normalized to one.
update_normal_vectors
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.
|
inherited |
Return the normal vectors at all quadrature points represented by this object. See the normal_vector() function for what the normal vectors represent.
update_normal_vectors
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 3973 of file fe_values.cc.
|
inherited |
Create a view of the current FEValues object that represents a particular scalar component of the possibly vector-valued finite element. The concept of views is explained in the documentation of the namespace FEValuesViews and in particular in the Handling vector valued problems module.
|
inherited |
Create a view of the current FEValues object that represents a set of dim
scalar components (i.e. a vector) of the vector-valued finite element. The concept of views is explained in the documentation of the namespace FEValuesViews and in particular in the Handling vector valued problems module.
|
inherited |
Create a view of the current FEValues object that represents a set of (dim*dim + dim)/2
scalar components (i.e. a symmetric 2nd order tensor) of the vector-valued finite element. The concept of views is explained in the documentation of the namespace FEValuesViews and in particular in the Handling vector valued problems module.
|
inherited |
Create a view of the current FEValues object that represents a set of (dim*dim)
scalar components (i.e. a 2nd order tensor) of the vector-valued finite element. The concept of views is explained in the documentation of the namespace FEValuesViews and in particular in the Handling vector valued problems module.
|
inherited |
Constant reference to the selected mapping object.
|
inherited |
Constant reference to the selected finite element object.
|
inherited |
Return the update flags set for this object.
|
inherited |
Return a triangulation iterator to the current cell.
Definition at line 3964 of file fe_values.cc.
|
inherited |
Return the relation of the current cell to the previous cell. This allows re-use of some cell data (like local matrices for equations with constant coefficients) if the result is CellSimilarity::translation
.
Definition at line 4142 of file fe_values.cc.
|
protectedinherited |
A function that is connected to the triangulation in order to reset the stored 'present_cell' iterator to an invalid one whenever the triangulation is changed and the iterator consequently becomes invalid.
Definition at line 4026 of file fe_values.cc.
|
protectedinherited |
This function is called by the various reinit() functions in derived classes. Given the cell indicated by the argument, test whether we have to throw away the previously stored present_cell argument because it would require us to compare cells from different triangulations. In checking all this, also make sure that we have tria_listener connected to the triangulation to which we will set present_cell right after calling this function.
Definition at line 4044 of file fe_values.cc.
|
protectedinherited |
Initialize some update flags. Called from the initialize
functions of derived classes, which are in turn called from their constructors.
Basically, this function finds out using the finite element and mapping object already stored which flags need to be set to compute everything the user wants, as expressed through the flags passed as argument.
Definition at line 4007 of file fe_values.cc.
|
inlineprotectedinherited |
A function that checks whether the new cell is similar to the one previously used. Then, a significant amount of the data can be reused, e.g. the derivatives of the basis functions in real space, shape_grad.
Definition at line 4087 of file fe_values.cc.
|
staticconstexpr |
Dimension in which this object operates.
Definition at line 4475 of file fe_values.h.
|
staticconstexpr |
Dimension of the space in which this object operates.
Definition at line 4480 of file fe_values.h.
|
staticconstexpr |
Dimension of the object over which we integrate. For the present class, this is equal to dim-1
.
Definition at line 4486 of file fe_values.h.
|
protectedinherited |
Number of the face selected the last time the reinit() function was called.
Definition at line 4273 of file fe_values.h.
|
protectedinherited |
Index of the face selected the last time the reinit() function was called.
Definition at line 4279 of file fe_values.h.
|
protectedinherited |
Store a copy of the quadrature formula here.
Definition at line 4284 of file fe_values.h.
|
inherited |
Number of quadrature points of the current object. Its value is initialized by the value of max_n_quadrature_points and is updated, e.g., if FEFaceValues::reinit() is called for a new cell/face.
Definition at line 2432 of file fe_values.h.
|
inherited |
Maximum number of quadrature points. This value might be different from n_quadrature_points, e.g., if a QCollection with different face quadrature rules has been passed to initialize FEFaceValues.
This is mostly useful to initialize arrays to allocate the maximum amount of memory that may be used when re-sizing later on to a the current number of quadrature points given by n_quadrature_points.
Definition at line 2443 of file fe_values.h.
|
inherited |
Number of shape functions per cell. If we use this base class to evaluate a finite element on faces of cells, this is still the number of degrees of freedom per cell, not per face.
Definition at line 2450 of file fe_values.h.
|
protectedinherited |
Store the cell selected last time the reinit() function was called. This is necessary for the get_function_*
functions as well as the functions of same name in the extractor classes.
Definition at line 3894 of file fe_values.h.
|
protectedinherited |
A signal connection we use to ensure we get informed whenever the triangulation changes by refinement. We need to know about that because it invalidates all cell iterators and, as part of that, the 'present_cell' iterator we keep around between subsequent calls to reinit() in order to compute the cell similarity.
Definition at line 3903 of file fe_values.h.
|
protectedinherited |
A signal connection we use to ensure we get informed whenever the triangulation changes by mesh transformations. We need to know about that because it invalidates all cell iterators and, as part of that, the 'present_cell' iterator we keep around between subsequent calls to reinit() in order to compute the cell similarity.
Definition at line 3912 of file fe_values.h.
|
protectedinherited |
A pointer to the mapping object associated with this FEValues object.
Definition at line 3939 of file fe_values.h.
|
protectedinherited |
A pointer to the internal data object of mapping, obtained from Mapping::get_data(), Mapping::get_face_data(), or Mapping::get_subface_data().
Definition at line 3947 of file fe_values.h.
|
protectedinherited |
An object into which the Mapping::fill_fe_values() and similar functions place their output.
Definition at line 3954 of file fe_values.h.
|
protectedinherited |
A pointer to the finite element object associated with this FEValues object.
Definition at line 3963 of file fe_values.h.
|
protectedinherited |
A pointer to the internal data object of finite element, obtained from FiniteElement::get_data(), Mapping::get_face_data(), or FiniteElement::get_subface_data().
Definition at line 3971 of file fe_values.h.
|
protectedinherited |
An object into which the FiniteElement::fill_fe_values() and similar functions place their output.
Definition at line 3979 of file fe_values.h.
|
protectedinherited |
Original update flags handed to the constructor of FEValues.
Definition at line 3985 of file fe_values.h.
|
protectedinherited |
An enum variable that can store different states of the current cell in comparison to the previously visited cell. If wanted, additional states can be checked here and used in one of the methods used during reinit.
Definition at line 4003 of file fe_values.h.
|
privateinherited |
A cache for all possible FEValuesViews objects.
Definition at line 4018 of file fe_values.h.