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