Reference documentation for deal.II version 9.6.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\}}\)
Loading...
Searching...
No Matches
fe_update_flags.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1998 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_fe_update_flags_h
16#define dealii_fe_update_flags_h
17
18
19#include <deal.II/base/config.h>
20
21#include <deal.II/base/table.h>
22#include <deal.II/base/tensor.h>
23
24#include <vector>
25
26
28
29// Forward declaration
30#ifndef DOXYGEN
31template <int, int>
32class FiniteElement;
33#endif
34
64{
68
75 update_values = 0x0001,
77
83
89
95
102
129
136
143
149
154
160
167
174
180
186
216 update_rescale = 0x2000000,
218
226 // Direct data
234 // Transformation dependence
237 // Volume data
239 // Hermite needs several DOFs to be rescaled each time
242
243
249template <typename StreamType>
250inline StreamType &
251operator<<(StreamType &s, const UpdateFlags u)
252{
253 s << " UpdateFlags|";
254 if (u & update_values)
255 s << "values|";
256 if (u & update_gradients)
257 s << "gradients|";
258 if (u & update_hessians)
259 s << "hessians|";
261 s << "3rd_derivatives|";
263 s << "quadrature_points|";
264 if (u & update_JxW_values)
265 s << "JxW_values|";
266 if (u & update_normal_vectors)
267 s << "normal_vectors|";
268 if (u & update_jacobians)
269 s << "jacobians|";
271 s << "inverse_jacobians|";
272 if (u & update_jacobian_grads)
273 s << "jacobian_grads|";
275 s << "covariant_transformation|";
277 s << "contravariant_transformation|";
279 s << "transformation_values|";
281 s << "transformation_gradients|";
283 s << "jacobian_pushed_forward_grads|";
285 s << "jacobian_2nd_derivatives|";
287 s << "jacobian_pushed_forward_2nd_derivatives|";
289 s << "jacobian_3rd_derivatives|";
291 s << "jacobian_pushed_forward_3rd_derivatives|";
292
293 // TODO: check that 'u' really only has the flags set that are handled above
294 return s;
295}
296
297
307inline UpdateFlags
309{
310 return static_cast<UpdateFlags>(static_cast<unsigned int>(f1) |
311 static_cast<unsigned int>(f2));
312}
313
314
315
322inline UpdateFlags &
324{
325 f1 = f1 | f2;
326 return f1;
327}
328
329
339inline UpdateFlags
341{
342 return static_cast<UpdateFlags>(static_cast<unsigned int>(f1) &
343 static_cast<unsigned int>(f2));
344}
345
346
353inline UpdateFlags &
355{
356 f1 = f1 & f2;
357 return f1;
358}
359
360
361
395
396
397namespace internal
398{
399 namespace FEValuesImplementation
400 {
409 template <int dim, int spacedim = dim>
523 } // namespace FEValuesImplementation
524} // namespace internal
525
526
532
533#endif
void initialize(const unsigned int n_quadrature_points, const FiniteElement< dim, spacedim > &fe, const UpdateFlags flags)
Definition fe_values.cc:83
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
UpdateFlags & operator&=(UpdateFlags &f1, const UpdateFlags f2)
StreamType & operator<<(StreamType &s, const UpdateFlags u)
UpdateFlags & operator|=(UpdateFlags &f1, const UpdateFlags f2)
UpdateFlags operator|(const UpdateFlags f1, const UpdateFlags f2)
UpdateFlags operator&(const UpdateFlags f1, const UpdateFlags f2)
UpdateFlags
@ update_jacobian_pushed_forward_2nd_derivatives
@ update_volume_elements
Determinant of the Jacobian.
@ update_contravariant_transformation
Contravariant transformation.
@ update_jacobian_pushed_forward_grads
@ update_rescale
@ update_hessians
Second derivatives of shape functions.
@ update_jacobian_3rd_derivatives
@ update_values
Shape function values.
@ update_jacobian_grads
Gradient of volume element.
@ update_normal_vectors
Normal vectors.
@ update_3rd_derivatives
Third derivatives of shape functions.
@ update_JxW_values
Transformed quadrature weights.
@ update_covariant_transformation
Covariant transformation.
@ update_transformation_gradients
Shape function gradients of transformation.
@ update_jacobians
Volume element.
@ update_mapping
@ update_inverse_jacobians
Volume element.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
@ update_default
No update.
@ update_jacobian_pushed_forward_3rd_derivatives
@ update_transformation_values
Shape function values of transformation.
@ update_piola
Values needed for Piola transform.
@ update_boundary_forms
Outer normal vector, not normalized.
@ update_jacobian_2nd_derivatives