|
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
16 #ifndef dealii_fe_update_flags_h
17 #define dealii_fe_update_flags_h
254 template <
class StreamType>
258 s <<
" UpdateFlags|";
266 s <<
"3rd_derivatives|";
268 s <<
"quadrature_points|";
272 s <<
"normal_vectors|";
276 s <<
"inverse_jacobians|";
278 s <<
"jacobian_grads|";
280 s <<
"covariant_transformation|";
282 s <<
"contravariant_transformation|";
284 s <<
"transformation_values|";
286 s <<
"transformation_gradients|";
288 s <<
"jacobian_pushed_forward_grads|";
290 s <<
"jacobian_2nd_derivatives|";
292 s <<
"jacobian_pushed_forward_2nd_derivatives|";
294 s <<
"jacobian_3rd_derivatives|";
296 s <<
"jacobian_pushed_forward_3rd_derivatives|";
315 return static_cast<UpdateFlags>(
static_cast<unsigned int>(f1) |
316 static_cast<unsigned int>(f2));
346 return static_cast<UpdateFlags>(
static_cast<unsigned int>(f1) &
347 static_cast<unsigned int>(f2));
403 namespace FEValuesImplementation
417 template <
int dim,
int spacedim = dim>
425 initialize(
const unsigned int n_quadrature_points,
453 std::vector<DerivativeForm<1, dim, spacedim>>
jacobians;
523 template <
int dim,
int spacedim = dim>
531 initialize(
const unsigned int n_quadrature_points,
@ update_quadrature_points
Transformed quadrature points.
@ update_jacobian_pushed_forward_3rd_derivatives
@ update_3rd_derivatives
Third derivatives of shape functions.
@ update_default
No update.
@ update_transformation_values
Shape function values of transformation.
UpdateFlags operator|(const UpdateFlags f1, const UpdateFlags f2)
UpdateFlags & operator|=(UpdateFlags &f1, const UpdateFlags f2)
@ update_values
Shape function values.
@ update_jacobian_pushed_forward_grads
@ update_piola
Values needed for Piola transform.
@ update_cell_normal_vectors
@ update_face_normal_vectors
@ update_covariant_transformation
Covariant transformation.
#define DEAL_II_NAMESPACE_OPEN
@ update_jacobians
Volume element.
@ update_gradients
Shape function gradients.
#define DEAL_II_DEPRECATED
@ update_jacobian_pushed_forward_2nd_derivatives
@ update_second_derivatives
@ update_jacobian_3rd_derivatives
UpdateFlags & operator&=(UpdateFlags &f1, const UpdateFlags f2)
UpdateFlags operator&(const UpdateFlags f1, const UpdateFlags f2)
@ update_jacobian_grads
Gradient of volume element.
@ update_JxW_values
Transformed quadrature weights.
@ update_normal_vectors
Normal vectors.
@ update_hessians
Second derivatives of shape functions.
@ update_transformation_gradients
Shape function gradients of transformation.
OutputOperator< VectorType > & operator<<(OutputOperator< VectorType > &out, unsigned int step)
@ update_inverse_jacobians
Volume element.
#define DEAL_II_NAMESPACE_CLOSE
@ update_volume_elements
Determinant of the Jacobian.
@ update_contravariant_transformation
Contravariant transformation.
@ update_boundary_forms
Outer normal vector, not normalized.
@ update_jacobian_2nd_derivatives