Reference documentation for deal.II version 8.5.1
fe_update_flags.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2016 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 
106 
113 
128 
134 
139 
145 
152 
159 
165 
171 
207 
215  // Direct data
228  // Transformation dependence
233  // Volume data
235 };
236 
237 
243 template <class StreamType>
244 inline
245 StreamType &operator << (StreamType &s, UpdateFlags u)
246 {
247  s << " UpdateFlags|";
248  if (u & update_values) s << "values|";
249  if (u & update_gradients) s << "gradients|";
250  if (u & update_hessians) s << "hessians|";
251  if (u & update_3rd_derivatives) s << "3rd_derivatives|";
252  if (u & update_quadrature_points) s << "quadrature_points|";
253  if (u & update_JxW_values) s << "JxW_values|";
254  if (u & update_normal_vectors) s << "normal_vectors|";
255  if (u & update_jacobians) s << "jacobians|";
256  if (u & update_inverse_jacobians) s << "inverse_jacobians|";
257  if (u & update_jacobian_grads) s << "jacobian_grads|";
258  if (u & update_covariant_transformation) s << "covariant_transformation|";
259  if (u & update_contravariant_transformation) s << "contravariant_transformation|";
260  if (u & update_transformation_values) s << "transformation_values|";
261  if (u & update_transformation_gradients) s << "transformation_gradients|";
262  if (u & update_jacobian_pushed_forward_grads) s << "jacobian_pushed_forward_grads|";
263  if (u & update_jacobian_2nd_derivatives) s << "jacobian_2nd_derivatives|";
264  if (u & update_jacobian_pushed_forward_2nd_derivatives) s << "jacobian_pushed_forward_2nd_derivatives|";
265  if (u &update_jacobian_3rd_derivatives) s << "jacobian_3rd_derivatives|";
266  if (u & update_jacobian_pushed_forward_3rd_derivatives) s << "jacobian_pushed_forward_3rd_derivatives|";
267 
268 //TODO: check that 'u' really only has the flags set that are handled above
269  return s;
270 }
271 
272 
282 inline
285 {
286  return static_cast<UpdateFlags> (
287  static_cast<unsigned int> (f1) |
288  static_cast<unsigned int> (f2));
289 }
290 
291 
292 
293 
300 inline
301 UpdateFlags &
303 {
304  f1 = f1 | f2;
305  return f1;
306 }
307 
308 
318 inline
321 {
322  return static_cast<UpdateFlags> (
323  static_cast<unsigned int> (f1) &
324  static_cast<unsigned int> (f2));
325 }
326 
327 
334 inline
335 UpdateFlags &
337 {
338  f1 = f1 & f2;
339  return f1;
340 }
341 
342 
343 
354 namespace CellSimilarity
355 {
357  {
375  };
376 }
377 
378 
379 namespace internal
380 {
381  namespace FEValues
382  {
395  template <int dim, int spacedim=dim>
397  {
398  public:
402  void initialize (const unsigned int n_quadrature_points,
403  const UpdateFlags flags);
404 
409  std::size_t memory_consumption () const;
410 
424  std::vector<double> JxW_values;
425 
429  std::vector< DerivativeForm<1,dim,spacedim> > jacobians;
430 
435  std::vector<DerivativeForm<2,dim,spacedim> > jacobian_grads;
436 
440  std::vector<DerivativeForm<1,spacedim,dim> > inverse_jacobians;
441 
446  std::vector<Tensor<3,spacedim> > jacobian_pushed_forward_grads;
447 
452  std::vector<DerivativeForm<3,dim,spacedim> > jacobian_2nd_derivatives;
453 
458  std::vector<Tensor<4,spacedim> > jacobian_pushed_forward_2nd_derivatives;
459 
464  std::vector<DerivativeForm<4,dim,spacedim> > jacobian_3rd_derivatives;
465 
470  std::vector<Tensor<5,spacedim> > jacobian_pushed_forward_3rd_derivatives;
471 
477  std::vector<Point<spacedim> > quadrature_points;
478 
482  std::vector<Tensor<1,spacedim> > normal_vectors;
483 
487  std::vector<Tensor<1,spacedim> > boundary_forms;
488  };
489 
490 
499  template <int dim, int spacedim=dim>
501  {
502  public:
506  void initialize (const unsigned int n_quadrature_points,
507  const FiniteElement<dim,spacedim> &fe,
508  const UpdateFlags flags);
509 
514  std::size_t memory_consumption () const;
515 
534  typedef ::Table<2,double> ShapeVector;
535 
540  typedef ::Table<2,Tensor<1,spacedim> > GradientVector;
541 
545  typedef ::Table<2,Tensor<2,spacedim> > HessianVector;
546 
550  typedef ::Table<2,Tensor<3,spacedim> > ThirdDerivativeVector;
551 
558 
565 
572 
579 
609  std::vector<unsigned int> shape_function_to_row_table;
610  };
611  }
612 }
613 
614 
619 DEAL_II_NAMESPACE_CLOSE
620 
621 #endif
Transformed quadrature weights.
Shape function values.
UpdateFlags & operator&=(UpdateFlags &f1, UpdateFlags f2)
Contravariant transformation.
Volume element.
std::size_t memory_consumption() const
Definition: fe_values.cc:2160
Outer normal vector, not normalized.
Determinant of the Jacobian.
Transformed quadrature points.
std::vector< unsigned int > shape_function_to_row_table
void initialize(const unsigned int n_quadrature_points, const FiniteElement< dim, spacedim > &fe, const UpdateFlags flags)
Definition: fe_values.cc:2181
std::vector< Tensor< 1, spacedim > > boundary_forms
std::vector< Tensor< 4, spacedim > > jacobian_pushed_forward_2nd_derivatives
::Table< 2, Tensor< 2, spacedim > > HessianVector
std::vector< DerivativeForm< 2, dim, spacedim > > jacobian_grads
Shape function gradients of transformation.
std::vector< Tensor< 3, spacedim > > jacobian_pushed_forward_grads
No update.
::Table< 2, Tensor< 3, spacedim > > ThirdDerivativeVector
Third derivatives of shape functions.
UpdateFlags
::Table< 2, Tensor< 1, spacedim > > GradientVector
void initialize(const unsigned int n_quadrature_points, const UpdateFlags flags)
Definition: fe_values.cc:2105
std::vector< DerivativeForm< 1, dim, spacedim > > jacobians
std::vector< Tensor< 1, spacedim > > normal_vectors
Shape function values of transformation.
UpdateFlags operator|(UpdateFlags f1, UpdateFlags f2)
Second derivatives of shape functions.
Gradient of volume element.
UpdateFlags operator&(UpdateFlags f1, UpdateFlags f2)
std::vector< Tensor< 5, spacedim > > jacobian_pushed_forward_3rd_derivatives
OutputOperator< VectorType > & operator<<(OutputOperator< VectorType > &out, unsigned int step)
Definition: operator.h:154
std::vector< DerivativeForm< 4, dim, spacedim > > jacobian_3rd_derivatives
Shape function gradients.
std::vector< Point< spacedim > > quadrature_points
Normal vectors.
Definition: fe.h:30
UpdateFlags & operator|=(UpdateFlags &f1, UpdateFlags f2)
std::vector< DerivativeForm< 1, spacedim, dim > > inverse_jacobians
std::vector< DerivativeForm< 3, dim, spacedim > > jacobian_2nd_derivatives
Values needed for Piola transform.
Covariant transformation.