Reference documentation for deal.II version 9.0.0
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 at
12 // the top level of the deal.II distribution.
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 #include <deal.II/base/table.h>
22 #include <deal.II/base/derivative_form.h>
23 #include <deal.II/base/point.h>
24 #include <deal.II/base/tensor.h>
25 
26 #include <vector>
27 
28 
29 DEAL_II_NAMESPACE_OPEN
30 
31 template <int,int> class FiniteElement;
32 
33 
36 
63 {
67 
74  update_values = 0x0001,
76 
80  update_gradients = 0x0002,
82 
86  update_hessians = 0x0004,
88 
94 
101 
117 
124 
139 
145 
150 
156 
163 
170 
176 
182 
218 
226  // Direct data
239  // Transformation dependence
244  // Volume data
246 };
247 
248 
254 template <class StreamType>
255 inline
256 StreamType &operator << (StreamType &s,
257  const UpdateFlags u)
258 {
259  s << " UpdateFlags|";
260  if (u & update_values) s << "values|";
261  if (u & update_gradients) s << "gradients|";
262  if (u & update_hessians) s << "hessians|";
263  if (u & update_3rd_derivatives) s << "3rd_derivatives|";
264  if (u & update_quadrature_points) s << "quadrature_points|";
265  if (u & update_JxW_values) s << "JxW_values|";
266  if (u & update_normal_vectors) s << "normal_vectors|";
267  if (u & update_jacobians) s << "jacobians|";
268  if (u & update_inverse_jacobians) s << "inverse_jacobians|";
269  if (u & update_jacobian_grads) s << "jacobian_grads|";
270  if (u & update_covariant_transformation) s << "covariant_transformation|";
271  if (u & update_contravariant_transformation) s << "contravariant_transformation|";
272  if (u & update_transformation_values) s << "transformation_values|";
273  if (u & update_transformation_gradients) s << "transformation_gradients|";
274  if (u & update_jacobian_pushed_forward_grads) s << "jacobian_pushed_forward_grads|";
275  if (u & update_jacobian_2nd_derivatives) s << "jacobian_2nd_derivatives|";
276  if (u & update_jacobian_pushed_forward_2nd_derivatives) s << "jacobian_pushed_forward_2nd_derivatives|";
277  if (u &update_jacobian_3rd_derivatives) s << "jacobian_3rd_derivatives|";
278  if (u & update_jacobian_pushed_forward_3rd_derivatives) 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
297  const UpdateFlags f2)
298 {
299  return static_cast<UpdateFlags> (
300  static_cast<unsigned int> (f1) |
301  static_cast<unsigned int> (f2));
302 }
303 
304 
305 
306 
313 inline
314 UpdateFlags &
316  const UpdateFlags f2)
317 {
318  f1 = f1 | f2;
319  return f1;
320 }
321 
322 
332 inline
335  const UpdateFlags f2)
336 {
337  return static_cast<UpdateFlags> (
338  static_cast<unsigned int> (f1) &
339  static_cast<unsigned int> (f2));
340 }
341 
342 
349 inline
350 UpdateFlags &
352  const UpdateFlags f2)
353 {
354  f1 = f1 & f2;
355  return f1;
356 }
357 
358 
359 
370 namespace CellSimilarity
371 {
373  {
391  };
392 }
393 
394 
395 namespace internal
396 {
397  namespace FEValuesImplementation
398  {
411  template <int dim, int spacedim=dim>
413  {
414  public:
418  void initialize (const unsigned int n_quadrature_points,
419  const UpdateFlags flags);
420 
425  std::size_t memory_consumption () const;
426 
440  std::vector<double> JxW_values;
441 
445  std::vector< DerivativeForm<1,dim,spacedim> > jacobians;
446 
451  std::vector<DerivativeForm<2,dim,spacedim> > jacobian_grads;
452 
456  std::vector<DerivativeForm<1,spacedim,dim> > inverse_jacobians;
457 
462  std::vector<Tensor<3,spacedim> > jacobian_pushed_forward_grads;
463 
468  std::vector<DerivativeForm<3,dim,spacedim> > jacobian_2nd_derivatives;
469 
474  std::vector<Tensor<4,spacedim> > jacobian_pushed_forward_2nd_derivatives;
475 
480  std::vector<DerivativeForm<4,dim,spacedim> > jacobian_3rd_derivatives;
481 
486  std::vector<Tensor<5,spacedim> > jacobian_pushed_forward_3rd_derivatives;
487 
493  std::vector<Point<spacedim> > quadrature_points;
494 
498  std::vector<Tensor<1,spacedim> > normal_vectors;
499 
503  std::vector<Tensor<1,spacedim> > boundary_forms;
504  };
505 
506 
515  template <int dim, int spacedim=dim>
517  {
518  public:
522  void initialize (const unsigned int n_quadrature_points,
523  const FiniteElement<dim,spacedim> &fe,
524  const UpdateFlags flags);
525 
530  std::size_t memory_consumption () const;
531 
550  typedef ::Table<2,double> ShapeVector;
551 
556  typedef ::Table<2,Tensor<1,spacedim> > GradientVector;
557 
561  typedef ::Table<2,Tensor<2,spacedim> > HessianVector;
562 
566  typedef ::Table<2,Tensor<3,spacedim> > ThirdDerivativeVector;
567 
574 
581 
588 
595 
625  std::vector<unsigned int> shape_function_to_row_table;
626  };
627  }
628 }
629 
630 
635 DEAL_II_NAMESPACE_CLOSE
636 
637 #endif
Transformed quadrature weights.
Shape function values.
std::vector< Tensor< 1, spacedim > > normal_vectors
Contravariant transformation.
std::vector< DerivativeForm< 1, dim, spacedim > > jacobians
std::vector< Tensor< 4, spacedim > > jacobian_pushed_forward_2nd_derivatives
std::vector< DerivativeForm< 4, dim, spacedim > > jacobian_3rd_derivatives
Volume element.
Outer normal vector, not normalized.
Determinant of the Jacobian.
Transformed quadrature points.
::Table< 2, Tensor< 3, spacedim > > ThirdDerivativeVector
std::vector< Tensor< 5, spacedim > > jacobian_pushed_forward_3rd_derivatives
UpdateFlags operator|(const UpdateFlags f1, const UpdateFlags f2)
Shape function gradients of transformation.
std::vector< DerivativeForm< 2, dim, spacedim > > jacobian_grads
No update.
Third derivatives of shape functions.
UpdateFlags
void initialize(const unsigned int n_quadrature_points, const UpdateFlags flags)
Definition: fe_values.cc:2587
Shape function values of transformation.
void initialize(const unsigned int n_quadrature_points, const FiniteElement< dim, spacedim > &fe, const UpdateFlags flags)
Definition: fe_values.cc:2662
Second derivatives of shape functions.
Gradient of volume element.
std::vector< Tensor< 3, spacedim > > jacobian_pushed_forward_grads
std::vector< DerivativeForm< 3, dim, spacedim > > jacobian_2nd_derivatives
std::vector< Point< spacedim > > quadrature_points
OutputOperator< VectorType > & operator<<(OutputOperator< VectorType > &out, unsigned int step)
Definition: operator.h:163
UpdateFlags & operator|=(UpdateFlags &f1, const UpdateFlags f2)
std::vector< DerivativeForm< 1, spacedim, dim > > inverse_jacobians
Shape function gradients.
Normal vectors.
UpdateFlags operator&(const UpdateFlags f1, const UpdateFlags f2)
std::vector< Tensor< 1, spacedim > > boundary_forms
UpdateFlags & operator&=(UpdateFlags &f1, const UpdateFlags f2)
Values needed for Piola transform.
Covariant transformation.