Reference documentation for deal.II version GIT b065060f03 2023-05-30 23:50:02+00:00
\(\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\}}\)
fe_update_flags.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2021 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_fe_update_flags_h
17 #define dealii_fe_update_flags_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #include <deal.II/base/table.h>
23 #include <deal.II/base/tensor.h>
24 
25 #include <vector>
26 
27 
29 
30 // Forward declaration
31 #ifndef DOXYGEN
32 template <int, int>
33 class FiniteElement;
34 #endif
35 
65 {
69 
76  update_values = 0x0001,
78 
82  update_gradients = 0x0002,
84 
88  update_hessians = 0x0004,
90 
96 
103 
130 
137 
144 
150 
155 
161 
168 
175 
181 
187 
215 
223  // Direct data
231  // Transformation dependence
234  // Volume data
236 };
237 
238 
244 template <class StreamType>
245 inline StreamType &
246 operator<<(StreamType &s, const UpdateFlags u)
247 {
248  s << " UpdateFlags|";
249  if (u & update_values)
250  s << "values|";
251  if (u & update_gradients)
252  s << "gradients|";
253  if (u & update_hessians)
254  s << "hessians|";
255  if (u & update_3rd_derivatives)
256  s << "3rd_derivatives|";
257  if (u & update_quadrature_points)
258  s << "quadrature_points|";
259  if (u & update_JxW_values)
260  s << "JxW_values|";
261  if (u & update_normal_vectors)
262  s << "normal_vectors|";
263  if (u & update_jacobians)
264  s << "jacobians|";
265  if (u & update_inverse_jacobians)
266  s << "inverse_jacobians|";
267  if (u & update_jacobian_grads)
268  s << "jacobian_grads|";
270  s << "covariant_transformation|";
272  s << "contravariant_transformation|";
274  s << "transformation_values|";
276  s << "transformation_gradients|";
278  s << "jacobian_pushed_forward_grads|";
280  s << "jacobian_2nd_derivatives|";
282  s << "jacobian_pushed_forward_2nd_derivatives|";
284  s << "jacobian_3rd_derivatives|";
286  s << "jacobian_pushed_forward_3rd_derivatives|";
287 
288  // TODO: check that 'u' really only has the flags set that are handled above
289  return s;
290 }
291 
292 
302 inline UpdateFlags
303 operator|(const UpdateFlags f1, const UpdateFlags f2)
304 {
305  return static_cast<UpdateFlags>(static_cast<unsigned int>(f1) |
306  static_cast<unsigned int>(f2));
307 }
308 
309 
310 
317 inline UpdateFlags &
319 {
320  f1 = f1 | f2;
321  return f1;
322 }
323 
324 
334 inline UpdateFlags
335 operator&(const UpdateFlags f1, const UpdateFlags f2)
336 {
337  return static_cast<UpdateFlags>(static_cast<unsigned int>(f1) &
338  static_cast<unsigned int>(f2));
339 }
340 
341 
348 inline UpdateFlags &
350 {
351  f1 = f1 & f2;
352  return f1;
353 }
354 
355 
356 
367 namespace CellSimilarity
368 {
370  {
388  };
389 }
390 
391 
392 namespace internal
393 {
394  namespace FEValuesImplementation
395  {
404  template <int dim, int spacedim = dim>
406  {
407  public:
411  void
412  initialize(const unsigned int n_quadrature_points,
414  const UpdateFlags flags);
415 
420  std::size_t
421  memory_consumption() const;
422 
442 
448 
453 
458 
465 
472 
479 
486 
516  std::vector<unsigned int> shape_function_to_row_table;
517  };
518  } // namespace FEValuesImplementation
519 } // namespace internal
520 
521 
527 
528 #endif
OutputOperator< VectorType > & operator<<(OutputOperator< VectorType > &out, unsigned int step)
Definition: operator.h:165
void initialize(const unsigned int n_quadrature_points, const FiniteElement< dim, spacedim > &fe, const UpdateFlags flags)
Definition: fe_values.cc:2702
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:474
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:475
IteratorRange< FilteredIterator< BaseIterator > > operator|(IteratorRange< BaseIterator > i, const Predicate &p)
UpdateFlags & operator|=(UpdateFlags &f1, const UpdateFlags f2)
UpdateFlags operator&(const UpdateFlags f1, const UpdateFlags f2)
UpdateFlags
UpdateFlags & operator&=(UpdateFlags &f1, const UpdateFlags f2)
@ update_jacobian_pushed_forward_2nd_derivatives
@ update_volume_elements
Determinant of the Jacobian.
@ update_contravariant_transformation
Contravariant transformation.
@ update_jacobian_pushed_forward_grads
@ 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