Reference documentation for deal.II version 9.1.0-pre
fe_update_flags.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2018 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/derivative_form.h>
23 #include <deal.II/base/point.h>
24 #include <deal.II/base/table.h>
25 #include <deal.II/base/tensor.h>
26 
27 #include <vector>
28 
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 template <int, int>
33 class FiniteElement;
34 
35 
38 
65 {
69 
76  update_values = 0x0001,
78 
82  update_gradients = 0x0002,
84 
88  update_hessians = 0x0004,
90 
96 
103 
119 
126 
141 
147 
152 
158 
165 
172 
178 
184 
220 
228  // Direct data
236  // Transformation dependence
239  // Volume data
241 };
242 
243 
249 template <class StreamType>
250 inline StreamType &
251 operator<<(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|";
260  if (u & update_3rd_derivatives)
261  s << "3rd_derivatives|";
262  if (u & update_quadrature_points)
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|";
270  if (u & update_inverse_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 
307 inline UpdateFlags
308 operator|(const UpdateFlags f1, const UpdateFlags f2)
309 {
310  return static_cast<UpdateFlags>(static_cast<unsigned int>(f1) |
311  static_cast<unsigned int>(f2));
312 }
313 
314 
315 
322 inline UpdateFlags &
324 {
325  f1 = f1 | f2;
326  return f1;
327 }
328 
329 
339 inline UpdateFlags operator&(const UpdateFlags f1, const UpdateFlags f2)
340 {
341  return static_cast<UpdateFlags>(static_cast<unsigned int>(f1) &
342  static_cast<unsigned int>(f2));
343 }
344 
345 
352 inline UpdateFlags &
354 {
355  f1 = f1 & f2;
356  return f1;
357 }
358 
359 
360 
371 namespace CellSimilarity
372 {
374  {
392  };
393 }
394 
395 
396 namespace internal
397 {
398  namespace FEValuesImplementation
399  {
412  template <int dim, int spacedim = dim>
414  {
415  public:
419  void
420  initialize(const unsigned int n_quadrature_points,
421  const UpdateFlags flags);
422 
427  std::size_t
428  memory_consumption() const;
429 
443  std::vector<double> JxW_values;
444 
448  std::vector<DerivativeForm<1, dim, spacedim>> jacobians;
449 
454  std::vector<DerivativeForm<2, dim, spacedim>> jacobian_grads;
455 
459  std::vector<DerivativeForm<1, spacedim, dim>> inverse_jacobians;
460 
465  std::vector<Tensor<3, spacedim>> jacobian_pushed_forward_grads;
466 
471  std::vector<DerivativeForm<3, dim, spacedim>> jacobian_2nd_derivatives;
472 
477  std::vector<Tensor<4, spacedim>> jacobian_pushed_forward_2nd_derivatives;
478 
483  std::vector<DerivativeForm<4, dim, spacedim>> jacobian_3rd_derivatives;
484 
489  std::vector<Tensor<5, spacedim>> jacobian_pushed_forward_3rd_derivatives;
490 
496  std::vector<Point<spacedim>> quadrature_points;
497 
501  std::vector<Tensor<1, spacedim>> normal_vectors;
502 
506  std::vector<Tensor<1, spacedim>> boundary_forms;
507  };
508 
509 
518  template <int dim, int spacedim = dim>
520  {
521  public:
525  void
526  initialize(const unsigned int n_quadrature_points,
528  const UpdateFlags flags);
529 
534  std::size_t
535  memory_consumption() const;
536 
556 
562 
567 
572 
579 
586 
593 
600 
630  std::vector<unsigned int> shape_function_to_row_table;
631  };
632  } // namespace FEValuesImplementation
633 } // namespace internal
634 
635 
640 DEAL_II_NAMESPACE_CLOSE
641 
642 #endif
Transformed quadrature weights.
Shape function values.
Contravariant transformation.
std::vector< Tensor< 1, spacedim > > boundary_forms
Volume element.
std::vector< DerivativeForm< 3, dim, spacedim > > jacobian_2nd_derivatives
Outer normal vector, not normalized.
Determinant of the Jacobian.
std::vector< DerivativeForm< 4, dim, spacedim > > jacobian_3rd_derivatives
Transformed quadrature points.
std::vector< Tensor< 4, spacedim > > jacobian_pushed_forward_2nd_derivatives
UpdateFlags & operator&=(UpdateFlags &f1, const UpdateFlags f2)
UpdateFlags operator&(const UpdateFlags f1, const UpdateFlags f2)
UpdateFlags operator|(const UpdateFlags f1, const UpdateFlags f2)
Shape function gradients of transformation.
No update.
Third derivatives of shape functions.
UpdateFlags
::Table< 2, Tensor< 3, spacedim >> ThirdDerivativeVector
Shape function values of transformation.
Second derivatives of shape functions.
Gradient of volume element.
std::vector< Tensor< 3, spacedim > > jacobian_pushed_forward_grads
std::vector< Point< spacedim > > quadrature_points
OutputOperator< VectorType > & operator<<(OutputOperator< VectorType > &out, unsigned int step)
Definition: operator.h:170
UpdateFlags & operator|=(UpdateFlags &f1, const UpdateFlags f2)
Shape function gradients.
Normal vectors.
std::vector< DerivativeForm< 1, dim, spacedim > > jacobians
std::vector< Tensor< 5, spacedim > > jacobian_pushed_forward_3rd_derivatives
std::vector< DerivativeForm< 2, dim, spacedim > > jacobian_grads
std::vector< DerivativeForm< 1, spacedim, dim > > inverse_jacobians
Values needed for Piola transform.
Covariant transformation.
std::vector< Tensor< 1, spacedim > > normal_vectors