Reference documentation for deal.II version 9.2.0
|
#include <deal.II/base/config.h>
#include <deal.II/base/derivative_form.h>
#include <deal.II/base/point.h>
#include <deal.II/base/table.h>
#include <deal.II/base/tensor.h>
#include <vector>
Go to the source code of this file.
Classes | |
class | internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > |
class | internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > |
Namespaces | |
CellSimilarity | |
internal | |
internal::FEValuesImplementation | |
Functions | |
template<class StreamType > | |
StreamType & | operator<< (StreamType &s, const UpdateFlags u) |
UpdateFlags | operator| (const UpdateFlags f1, const UpdateFlags f2) |
UpdateFlags & | operator|= (UpdateFlags &f1, const UpdateFlags f2) |
UpdateFlags | operator& (const UpdateFlags f1, const UpdateFlags f2) |
UpdateFlags & | operator&= (UpdateFlags &f1, const UpdateFlags f2) |
enum UpdateFlags |
The enum type given to the constructors of FEValues, FEFaceValues and FESubfaceValues, telling those objects which data will be needed on each mesh cell.
Selecting these flags in a restrictive way is crucial for the efficiency of FEValues::reinit(), FEFaceValues::reinit() and FESubfaceValues::reinit(). Therefore, only the flags actually needed should be selected. It is the responsibility of the involved Mapping and FiniteElement to add additional flags according to their own requirements. For instance, most finite elements will add update_covariant_transformation if update_gradients is selected. By default, all flags are off, i.e. no reinitialization will be done.
You can select more than one flag by concatenation using the bitwise or operator|(UpdateFlags,UpdateFlags).
More information on the use of this type both in user code as well as internally can be found in the documentation modules on The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues and How Mapping, FiniteElement, and FEValues work together.
Enumerator | |
---|---|
update_default | No update. |
update_values | Shape function values. Compute the values of the shape functions at the quadrature points on the real space cell. For the usual Lagrange elements, these values are equal to the values of the shape functions at the quadrature points on the unit cell, but they are different for more complicated elements, such as FE_RaviartThomas elements. |
update_gradients | Shape function gradients. Compute the gradients of the shape functions in coordinates of the real cell. |
update_hessians | Second derivatives of shape functions. Compute the second derivatives of the shape functions in coordinates of the real cell. |
update_3rd_derivatives | Third derivatives of shape functions. Compute the third derivatives of the shape functions in coordinates of the real cell |
update_boundary_forms | Outer normal vector, not normalized. Vector product of tangential vectors, yielding a normal vector with a length corresponding to the surface element; may be more efficient than computing both. |
update_quadrature_points | Transformed quadrature points. Compute the quadrature points location in real cell coordinates. FEValues objects take the quadrature point locations on the reference cell as an argument of the constructor (via the Quadrature object). For most finite elements, knowing the location of quadrature points on the reference cell is all that is necessary to evaluate shape functions, evaluate the mapping, and other things. On the other hand, if you want to evaluate a right hand side function \(f(\mathbf x_q)\) at quadrature point locations \(\mathbf x_q\) on the real cell, you need to pass this flag to the FEValues constructor to make sure you can later access them. In the context of DataPostprocessor, DataPostprocessorInputs::CommonInputs::evaluation_points will be updated. |
update_JxW_values | Transformed quadrature weights. Compute the quadrature weights on the real cell, i.e. the weights of the quadrature rule multiplied with the determinant of the Jacobian of the transformation from reference to real cell. |
update_normal_vectors | Normal vectors. Compute the normal vectors, either for a face or for a cell of codimension one. Setting this flag for any other object will raise an error. |
update_face_normal_vectors | |
update_cell_normal_vectors | |
update_jacobians | Volume element. Compute the Jacobian of the transformation from the reference cell to the real cell. |
update_jacobian_grads | Gradient of volume element. Compute the derivatives of the Jacobian of the transformation. |
update_inverse_jacobians | Volume element. Compute the inverse Jacobian of the transformation from the reference cell to the real cell. |
update_covariant_transformation | Covariant transformation. Compute all values the Mapping needs to perform a contravariant transformation of vectors. For special mappings like MappingCartesian this may be simpler than update_inverse_jacobians. |
update_contravariant_transformation | Contravariant transformation. Compute all values the Mapping needs to perform a contravariant transformation of vectors. For special mappings like MappingCartesian this may be simpler than update_jacobians. |
update_transformation_values | Shape function values of transformation. Compute the shape function values of the transformation defined by the Mapping. |
update_transformation_gradients | Shape function gradients of transformation. Compute the shape function gradients of the transformation defined by the Mapping. |
update_volume_elements | Determinant of the Jacobian. Compute the volume element in each quadrature point. |
update_jacobian_pushed_forward_grads | Compute the derivatives of the Jacobian of the transformation pushed forward to the real cell coordinates. |
update_jacobian_2nd_derivatives | Compute the second derivatives of the Jacobian of the transformation. |
update_jacobian_pushed_forward_2nd_derivatives | Compute the second derivatives of the Jacobian of the transformation pushed forward to the real cell coordinates. |
update_jacobian_3rd_derivatives | Compute the third derivatives of the Jacobian of the transformation. |
update_jacobian_pushed_forward_3rd_derivatives | Compute the third derivatives of the Jacobian of the transformation pushed forward to the real cell coordinates. |
update_q_points | |
update_second_derivatives | |
update_piola | Values needed for Piola transform. Combination of the flags needed for Piola transform of Hdiv elements. |
update_mapping | Combination of the flags that require a mapping calculation |
Definition at line 66 of file fe_update_flags.h.
|
inline |
Output operator which outputs update flags as a set of or'd text values.
The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues
Definition at line 256 of file fe_update_flags.h.
|
inline |
Global operator which returns an object in which all bits are set which are either set in the first or the second argument. This operator exists since if it did not then the result of the bit-or operator |
would be an integer which would in turn trigger a compiler warning when we tried to assign it to an object of type UpdateFlags.
The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues
Definition at line 313 of file fe_update_flags.h.
|
inline |
Global operator which sets the bits from the second argument also in the first one.
The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues
Definition at line 328 of file fe_update_flags.h.
|
inline |
Global operator which returns an object in which all bits are set which are set in the first as well as the second argument. This operator exists since if it did not then the result of the bit-and operator &
would be an integer which would in turn trigger a compiler warning when we tried to assign it to an object of type UpdateFlags.
The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues
Definition at line 344 of file fe_update_flags.h.
|
inline |
Global operator which clears all the bits in the first argument if they are not also set in the second argument.
The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues
Definition at line 358 of file fe_update_flags.h.