Flow123d  release_2.2.0-914-gf1a3a4f
update_flags.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file update_flags.hh
15  * @brief Enum type UpdateFlags indicates which quantities are to be
16  * recomputed on each finite element cell.
17  * @author Jan Stebel
18  */
19 
20 #ifndef UPDATE_FLAGS_HH_
21 #define UPDATE_FLAGS_HH_
22 
23 /**
24  * \enum UpdateFlags
25  * @brief Enum type UpdateFlags indicates which quantities are to be
26  * recomputed on each finite element cell.
27  *
28  * Selecting these flags in a restrictive way is crucial for the
29  * efficiency of FEValues::reinit() and FESideValues::reinit().
30  * Therefore, only the flags actually
31  * needed should be selected. It is the responsibility of the involved
32  * Mapping and FiniteElement to add additional flags according to
33  * their own requirements.
34 
35  * By default, all flags are off, i.e. no reinitialization will be
36  * done.
37  *
38  * You can select more than one flag by concatenation
39  * using the bitwise or operator|(UpdateFlags,UpdateFlags).
40  *
41  * <h3>Generating the actual flags</h3>
42  *
43  * When given a set of UpdateFlags @p flags, the FEValues object must
44  * determine, which values will have to be computed once only for the
45  * reference cell and which values will have to be updated for each
46  * cell. Here, it is important to note that in many cases, the
47  * FiniteElement will require additional updates from the Mapping. To
48  * this end, the following auxiliary functions have been implemented:
49  *
50  * FiniteElement::update_each() determine the values required by
51  * the FiniteElement on each cell. The same function exists in Mapping.
52  *
53  * FEValuesBase::update_each() is used to compute the union
54  * of all values to be computed ever. It does this by first adding to
55  * the flags set by the user all flags added by the FiniteElement.
56  * This new set of flags is then given to the Mapping and all flags
57  * required there are added.
58  *
59  * This union of all flags is given to Mapping::fill_fe_values() and
60  * FiniteElement::fill_fe_values(), where the quantities indicated by
61  * the flags are computed.
62  *
63  * The flags finally stored in FEValues then are the union of all the
64  * flags required by the user, by FiniteElement and by Mapping, for
65  * computation once or on each cell.
66  */
68 {
69  //! No update
71  //! Shape function values
72  /**
73  * Compute the values of the
74  * shape functions at the
75  * quadrature points on the
76  * real space cell. For the
77  * usual Lagrange elements,
78  * these values are equal to
79  * the values of the shape
80  * functions at the quadrature
81  * points on the unit cell, but
82  * they are different for more
83  * complicated elements, such
84  * as Raviart-Thomas
85  * elements.
86  */
87  update_values = 0x0001,
88  //! Shape function gradients
89  /**
90  * Compute the gradients of the
91  * shape functions in
92  * coordinates of the real
93  * cell.
94  */
95  update_gradients = 0x0002,
96  //! Transformed quadrature points
97  /**
98  * Compute the quadrature
99  * points transformed into real
100  * cell coordinates.
101  */
103  //! Transformed quadrature weights
104  /**
105  * Compute the quadrature
106  * weights on the real cell,
107  * i.e. the weights of the
108  * quadrature rule multiplied
109  * with the determinant of the
110  * Jacobian of the
111  * transformation from
112  * reference to real cell.
113  */
115  //! Normal vectors
116  /**
117  * Compute the normal vectors,
118  * either for a face or for a
119  * cell of codimension
120  * one. Setting this flag for
121  * any other object will raise
122  * an error.
123  */
125  //! Volume element
126  /**
127  * Compute the Jacobian of the
128  * transformation from the
129  * reference cell to the real
130  * cell.
131  */
133  //! Volume element
134  /**
135  * Compute the inverse
136  * Jacobian of the
137  * transformation from the
138  * reference cell to the real
139  * cell.
140  */
142  //! Determinant of the Jacobian
143  /**
144  * Compute the volume element
145  * in each quadrature point.
146  */
148  //! Transformed quadrature weight for cell sides
149  /**
150  * Same as update_JxW_values but for quadratures living
151  * on a side of the cell.
152  */
154 };
155 
156 
157 /**
158  * Output operator which outputs update flags as a set of or'd text values.
159  *
160  * @ref UpdateFlags
161  */
162 template <class STREAM>
163 inline
164 STREAM& operator << (STREAM& s, UpdateFlags u)
165 {
166  s << " UpdateFlags|";
167  if (u & update_values) s << "values|";
168  if (u & update_gradients) s << "gradients|";
169  if (u & update_quadrature_points) s << "quadrature_points|";
170  if (u & update_JxW_values) s << "JxW_values|";
171  if (u & update_normal_vectors) s << "normal_vectors|";
172  if (u & update_jacobians) s << "jacobians|";
173  if (u & update_inverse_jacobians) s << "inverse_jacobians|";
174  return s;
175 }
176 
177 
178 /**
179  * Global operator which returns an object in which all bits are set
180  * which are either set in the first or the second argument. This
181  * operator exists since if it did not then the result of the bit-or
182  * <tt>operator |</tt> would be an integer which would in turn trigger
183  * a compiler warning when we tried to assign it to an object of type
184  * UpdateFlags.
185  *
186  * @ref UpdateFlags
187  */
188 inline
191 {
192  return static_cast<UpdateFlags> (
193  static_cast<unsigned int> (f1) |
194  static_cast<unsigned int> (f2));
195 }
196 
197 
198 
199 
200 /**
201  * Global operator which sets the bits from the second argument also
202  * in the first one.
203  *
204  * @ref UpdateFlags
205  */
206 inline
207 UpdateFlags &
209 {
210  f1 = f1 | f2;
211  return f1;
212 }
213 
214 
215 /**
216  * Global operator which returns an object in which all bits are set
217  * which are set in the first as well as the second argument. This
218  * operator exists since if it did not then the result of the bit-and
219  * <tt>operator &</tt> would be an integer which would in turn trigger
220  * a compiler warning when we tried to assign it to an object of type
221  * UpdateFlags.
222  *
223  * @ref UpdateFlags
224  */
225 inline
228 {
229  return static_cast<UpdateFlags> (
230  static_cast<unsigned int> (f1) &
231  static_cast<unsigned int> (f2));
232 }
233 
234 
235 /**
236  * Global operator which clears all the bits in the first argument if
237  * they are not also set in the second argument.
238  *
239  * @ref UpdateFlags
240  */
241 inline
242 UpdateFlags &
244 {
245  f1 = f1 & f2;
246  return f1;
247 }
248 
249 
250 
251 
252 
253 
254 
255 
256 #endif /* UPDATE_FLAGS_HH_ */
Shape function values.
Definition: update_flags.hh:87
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Definition: update_flags.hh:67
No update.
Definition: update_flags.hh:70
Transformed quadrature weight for cell sides.
Determinant of the Jacobian.
UpdateFlags & operator&=(UpdateFlags &f1, UpdateFlags f2)
UpdateFlags operator|(UpdateFlags f1, UpdateFlags f2)
Volume element.
UpdateFlags & operator|=(UpdateFlags &f1, UpdateFlags f2)
Transformed quadrature points.
Shape function gradients.
Definition: update_flags.hh:95
STREAM & operator<<(STREAM &s, UpdateFlags u)
Normal vectors.
UpdateFlags operator&(UpdateFlags f1, UpdateFlags f2)
Transformed quadrature weights.