Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
mapping_data_on_the_fly.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2014 - 2020 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.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16
17#ifndef dealii_matrix_free_mapping_data_on_the_fly_h
18#define dealii_matrix_free_mapping_data_on_the_fly_h
19
20
21#include <deal.II/base/config.h>
22
27
31
34
35
37
38
39namespace internal
40{
41 namespace MatrixFreeFunctions
42 {
59 template <int dim, typename Number, typename VectorizedArrayType>
61 {
62 static_assert(
63 std::is_same<Number, typename VectorizedArrayType::value_type>::value,
64 "Type of Number and of VectorizedArrayType do not match.");
65
66 public:
74 const Quadrature<1> &quadrature,
75 const UpdateFlags update_flags);
76
83 const UpdateFlags update_flags);
84
88 void
89 reinit(typename ::Triangulation<dim>::cell_iterator cell);
90
95 bool
97
101 typename ::Triangulation<dim>::cell_iterator
102 get_cell() const;
103
110 const ::FEValues<dim> &
112
121
125 const Quadrature<1> &
127
128 private:
133 typename ::Triangulation<dim>::cell_iterator present_cell;
134
140
145
150
157 };
158
159
160 /*-------------------------- Inline functions ---------------------------*/
161
162 template <int dim, typename Number, typename VectorizedArrayType>
164 MappingDataOnTheFly(const Mapping<dim> & mapping,
165 const Quadrature<1> &quadrature,
166 const UpdateFlags update_flags)
167 : fe_values(
168 mapping,
169 fe_dummy,
170 Quadrature<dim>(quadrature),
171 MappingInfo<dim, Number, VectorizedArrayType>::compute_update_flags(
172 update_flags))
173 , quadrature_1d(quadrature)
174 {
176 mapping_info_storage.descriptor[0].initialize(quadrature);
178 mapping_info_storage.JxW_values.resize(fe_values.n_quadrature_points);
179 mapping_info_storage.jacobians[0].resize(fe_values.n_quadrature_points);
180 if (update_flags & update_quadrature_points)
181 {
184 fe_values.n_quadrature_points);
185 }
186 if (fe_values.get_update_flags() & update_normal_vectors)
187 {
189 fe_values.n_quadrature_points);
191 fe_values.n_quadrature_points);
192 }
193 Assert(!(fe_values.get_update_flags() & update_jacobian_grads),
195 }
196
197
198
199 template <int dim, typename Number, typename VectorizedArrayType>
201 MappingDataOnTheFly(const Quadrature<1> &quadrature,
202 const UpdateFlags update_flags)
203 : MappingDataOnTheFly(::StaticMappingQ1<dim, dim>::mapping,
204 quadrature,
205 update_flags)
206 {}
207
208
209
210 template <int dim, typename Number, typename VectorizedArrayType>
211 inline void
213 typename ::Triangulation<dim>::cell_iterator cell)
214 {
215 if (present_cell == cell)
216 return;
217 present_cell = cell;
218 fe_values.reinit(present_cell);
219 for (unsigned int q = 0; q < fe_values.get_quadrature().size(); ++q)
220 {
221 if (fe_values.get_update_flags() & update_JxW_values)
222 mapping_info_storage.JxW_values[q] = fe_values.JxW(q);
223 if (fe_values.get_update_flags() & update_jacobians)
224 {
225 Tensor<2, dim> jac = fe_values.jacobian(q);
226 jac = invert(transpose(jac));
227 for (unsigned int d = 0; d < dim; ++d)
228 for (unsigned int e = 0; e < dim; ++e)
229 mapping_info_storage.jacobians[0][q][d][e] = jac[d][e];
230 }
231 if (fe_values.get_update_flags() & update_quadrature_points)
232 for (unsigned int d = 0; d < dim; ++d)
233 mapping_info_storage.quadrature_points[q][d] =
234 fe_values.quadrature_point(q)[d];
235 if (fe_values.get_update_flags() & update_normal_vectors)
236 {
237 for (unsigned int d = 0; d < dim; ++d)
238 mapping_info_storage.normal_vectors[q][d] =
239 fe_values.normal_vector(q)[d];
240 mapping_info_storage.normals_times_jacobians[0][q] =
241 mapping_info_storage.normal_vectors[q] *
242 mapping_info_storage.jacobians[0][q];
243 }
244 }
245 }
246
247
248
249 template <int dim, typename Number, typename VectorizedArrayType>
250 inline bool
252 const
253 {
254 return present_cell !=
255 typename ::Triangulation<dim>::cell_iterator();
256 }
257
258
259
260 template <int dim, typename Number, typename VectorizedArrayType>
261 inline typename ::Triangulation<dim>::cell_iterator
263 {
264 return fe_values.get_cell();
265 }
266
267
268
269 template <int dim, typename Number, typename VectorizedArrayType>
270 inline const ::FEValues<dim> &
272 {
273 return fe_values;
274 }
275
276
277
278 template <int dim, typename Number, typename VectorizedArrayType>
281 const
282 {
283 return mapping_info_storage;
284 }
285
286
287
288 template <int dim, typename Number, typename VectorizedArrayType>
289 inline const Quadrature<1> &
291 const
292 {
293 return quadrature_1d;
294 }
295
296
297 } // end of namespace MatrixFreeFunctions
298} // end of namespace internal
299
300
302
303#endif
void resize(const size_type new_size)
Abstract base class for mapping classes.
Definition: mapping.h:304
const MappingInfoStorage< dim, dim, Number, VectorizedArrayType > & get_data_storage() const
typename::Triangulation< dim >::cell_iterator present_cell
MappingInfoStorage< dim, dim, Number, VectorizedArrayType > mapping_info_storage
void reinit(typename ::Triangulation< dim >::cell_iterator cell)
typename::Triangulation< dim >::cell_iterator get_cell() const
MappingDataOnTheFly(const Quadrature< 1 > &quadrature, const UpdateFlags update_flags)
MappingDataOnTheFly(const Mapping< dim > &mapping, const Quadrature< 1 > &quadrature, const UpdateFlags update_flags)
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
UpdateFlags
@ update_jacobian_grads
Gradient of volume element.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_jacobians
Volume element.
@ update_quadrature_points
Transformed quadrature points.
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
Definition: exceptions.h:1465
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
AlignedVector< Point< spacedim, VectorizedArrayType > > quadrature_points
Definition: mapping_info.h:292
AlignedVector< unsigned int > quadrature_point_offsets
Definition: mapping_info.h:283
AlignedVector< Tensor< 1, spacedim, VectorizedArrayType > > normal_vectors
Definition: mapping_info.h:229
AlignedVector< VectorizedArrayType > JxW_values
Definition: mapping_info.h:222
AlignedVector< unsigned int > data_index_offsets
Definition: mapping_info.h:213
std::array< AlignedVector< Tensor< 1, spacedim, VectorizedArrayType > >, 2 > normals_times_jacobians
Definition: mapping_info.h:275
std::array< AlignedVector< Tensor< 2, spacedim, VectorizedArrayType > >, 2 > jacobians
Definition: mapping_info.h:243
std::vector< QuadratureDescriptor > descriptor
Definition: mapping_info.h:197
constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)