Reference documentation for deal.II version GIT a2efd28e10 2023-02-01 14:45: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 
122 
129 
136 
142 
147 
153 
160 
167 
173 
179 
207 
215  // Direct data
223  // Transformation dependence
226  // Volume data
228 };
229 
230 
236 template <class StreamType>
237 inline StreamType &
238 operator<<(StreamType &s, const UpdateFlags u)
239 {
240  s << " UpdateFlags|";
241  if (u & update_values)
242  s << "values|";
243  if (u & update_gradients)
244  s << "gradients|";
245  if (u & update_hessians)
246  s << "hessians|";
247  if (u & update_3rd_derivatives)
248  s << "3rd_derivatives|";
249  if (u & update_quadrature_points)
250  s << "quadrature_points|";
251  if (u & update_JxW_values)
252  s << "JxW_values|";
253  if (u & update_normal_vectors)
254  s << "normal_vectors|";
255  if (u & update_jacobians)
256  s << "jacobians|";
257  if (u & update_inverse_jacobians)
258  s << "inverse_jacobians|";
259  if (u & update_jacobian_grads)
260  s << "jacobian_grads|";
262  s << "covariant_transformation|";
264  s << "contravariant_transformation|";
266  s << "transformation_values|";
268  s << "transformation_gradients|";
270  s << "jacobian_pushed_forward_grads|";
272  s << "jacobian_2nd_derivatives|";
274  s << "jacobian_pushed_forward_2nd_derivatives|";
276  s << "jacobian_3rd_derivatives|";
278  s << "jacobian_pushed_forward_3rd_derivatives|";
279 
280  // TODO: check that 'u' really only has the flags set that are handled above
281  return s;
282 }
283 
284 
294 inline UpdateFlags
295 operator|(const UpdateFlags f1, const UpdateFlags f2)
296 {
297  return static_cast<UpdateFlags>(static_cast<unsigned int>(f1) |
298  static_cast<unsigned int>(f2));
299 }
300 
301 
302 
309 inline UpdateFlags &
311 {
312  f1 = f1 | f2;
313  return f1;
314 }
315 
316 
326 inline UpdateFlags
327 operator&(const UpdateFlags f1, const UpdateFlags f2)
328 {
329  return static_cast<UpdateFlags>(static_cast<unsigned int>(f1) &
330  static_cast<unsigned int>(f2));
331 }
332 
333 
340 inline UpdateFlags &
342 {
343  f1 = f1 & f2;
344  return f1;
345 }
346 
347 
348 
359 namespace CellSimilarity
360 {
362  {
380  };
381 }
382 
383 
384 namespace internal
385 {
386  namespace FEValuesImplementation
387  {
396  template <int dim, int spacedim = dim>
398  {
399  public:
403  void
404  initialize(const unsigned int n_quadrature_points,
406  const UpdateFlags flags);
407 
412  std::size_t
413  memory_consumption() const;
414 
434 
440 
445 
450 
457 
464 
471 
478 
508  std::vector<unsigned int> shape_function_to_row_table;
509  };
510  } // namespace FEValuesImplementation
511 } // namespace internal
512 
513 
519 
520 #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:2704
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:461
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:462
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