Reference documentation for deal.II version 8.1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
fe_update_flags.h
1 // ---------------------------------------------------------------------
2 // @f$Id: fe_update_flags.h 30036 2013-07-18 16:55:32Z maier @f$
3 //
4 // Copyright (C) 1998 - 2013 by the deal.II authors
5 //
6 // This file is part of the deal.II library.
7 //
8 // The deal.II library is free software; you can use it, redistribute
9 // it, and/or modify it under the terms of the GNU Lesser General
10 // Public License as published by the Free Software Foundation; either
11 // version 2.1 of the License, or (at your option) any later version.
12 // The full text of the license can be found in the file LICENSE at
13 // the top level of the deal.II distribution.
14 //
15 // ---------------------------------------------------------------------
16 
17 #ifndef __deal2__fe_update_flags_h
18 #define __deal2__fe_update_flags_h
19 
20 
21 #include <deal.II/base/config.h>
22 
24 
27 
84 {
88 
95  update_values = 0x0001,
97 
103 
107  update_hessians = 0x0004,
109 
116 
121 
128 
143 
149 
154 
160 
167 
174 
180 
186 
212 
216 };
217 
218 
224 template <class STREAM>
225 inline
226 STREAM &operator << (STREAM &s, UpdateFlags u)
227 {
228  s << " UpdateFlags|";
229  if (u & update_values) s << "values|";
230  if (u & update_gradients) s << "gradients|";
231  if (u & update_hessians) s << "hessians|";
232  if (u & update_quadrature_points) s << "quadrature_points|";
233  if (u & update_JxW_values) s << "JxW_values|";
234  if (u & update_normal_vectors) s << "normal_vectors|";
235  if (u & update_jacobians) s << "jacobians|";
236  if (u & update_inverse_jacobians) s << "inverse_jacobians|";
237  if (u & update_jacobian_grads) s << "jacobian_grads|";
238  if (u & update_covariant_transformation) s << "covariant_transformation|";
239  if (u & update_contravariant_transformation) s << "contravariant_transformation|";
240  if (u & update_transformation_values) s << "transformation_values|";
241  if (u & update_transformation_gradients) s << "transformation_gradients|";
242  if (u & update_support_points) s << "support_points|";
243  if (u & update_support_jacobians) s << "support_jacobians|";
244  if (u & update_support_inverse_jacobians) s << "support_inverse_jacobians|";
245 
246 //TODO: check that 'u' really only has the flags set that are handled above
247  return s;
248 }
249 
250 
261 inline
264 {
265  return static_cast<UpdateFlags> (
266  static_cast<unsigned int> (f1) |
267  static_cast<unsigned int> (f2));
268 }
269 
270 
271 
272 
279 inline
280 UpdateFlags &
282 {
283  f1 = f1 | f2;
284  return f1;
285 }
286 
287 
298 inline
301 {
302  return static_cast<UpdateFlags> (
303  static_cast<unsigned int> (f1) &
304  static_cast<unsigned int> (f2));
305 }
306 
307 
314 inline
315 UpdateFlags &
317 {
318  f1 = f1 & f2;
319  return f1;
320 }
321 
322 
323 
334 namespace CellSimilarity
335 {
336  enum Similarity
337  {
338  none,
339  translation,
340  inverted_translation,
341  invalid_next_cell
342  };
343 }
344 
345 
350 DEAL_II_NAMESPACE_CLOSE
351 
352 #endif
Transformed quadrature weights.
Shape function values.
UpdateFlags operator&(UpdateFlags f1, UpdateFlags f2)
Contravariant transformation.
Volume element.
Outer normal vector, not normalized.
Determinant of the Jacobian.
Transformed quadrature points.
UpdateFlags & operator&=(UpdateFlags &f1, UpdateFlags f2)
Shape function gradients of transformation.
No update.
UpdateFlags
Shape function values of transformation.
UpdateFlags operator|(UpdateFlags f1, UpdateFlags f2)
Second derivatives of shape functions.
Gradient of volume element.
OS & operator<<(OS &o, const Event &e)
Definition: event.h:304
Shape function gradients.
Normal vectors.
UpdateFlags & operator|=(UpdateFlags &f1, UpdateFlags f2)
Values needed for Piola transform.
Covariant transformation.