Reference documentation for deal.II version 9.6.0
\(\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\}}\)
Loading...
Searching...
No Matches
mapping_data_on_the_fly.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2014 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
16#ifndef dealii_matrix_free_mapping_data_on_the_fly_h
17#define dealii_matrix_free_mapping_data_on_the_fly_h
18
19
20#include <deal.II/base/config.h>
21
26
30
33
34
36
37
38namespace internal
39{
40 namespace MatrixFreeFunctions
41 {
58 template <int dim, typename Number>
60 {
61 public:
66
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
131
135 const Quadrature<1> &
137
138 private:
143 typename ::Triangulation<dim>::cell_iterator present_cell;
144
149 std::unique_ptr<FE_Nothing<dim>> fe_dummy;
150
154 std::unique_ptr<::FEValues<dim>> fe_values;
155
160
166 };
167
168
169 /*-------------------------- Inline functions ---------------------------*/
170
171 template <int dim, typename Number>
173 const Mapping<dim> &mapping,
174 const Quadrature<1> &quadrature,
175 const UpdateFlags update_flags)
176 : fe_dummy(std::make_unique<FE_Nothing<dim>>())
177 , fe_values(std::make_unique<::FEValues<dim>>(
178 mapping,
179 *fe_dummy,
180 Quadrature<dim>(quadrature),
181 MappingInfoStorage<dim, dim, Number>::compute_update_flags(
182 update_flags)))
183 , quadrature_1d(quadrature)
184 {
186 mapping_info_storage.descriptor[0].initialize(quadrature);
188 mapping_info_storage.JxW_values.resize(fe_values->n_quadrature_points);
189 mapping_info_storage.jacobians[0].resize(fe_values->n_quadrature_points);
190 if (update_flags & update_quadrature_points)
191 {
194 fe_values->n_quadrature_points);
195 }
196 if (fe_values->get_update_flags() & update_normal_vectors)
197 {
199 fe_values->n_quadrature_points);
201 fe_values->n_quadrature_points);
202 }
203 Assert(!(fe_values->get_update_flags() & update_jacobian_grads),
205 }
206
207
208
209 template <int dim, typename Number>
211 const Quadrature<1> &quadrature,
212 const UpdateFlags update_flags)
213 : MappingDataOnTheFly(::StaticMappingQ1<dim, dim>::mapping,
214 quadrature,
215 update_flags)
216 {}
217
218
219
220 template <int dim, typename Number>
221 inline void
223 typename ::Triangulation<dim>::cell_iterator cell)
224 {
225 if (present_cell == cell)
226 return;
227 present_cell = cell;
228 fe_values->reinit(present_cell);
229 for (unsigned int q = 0; q < fe_values->get_quadrature().size(); ++q)
230 {
231 if (fe_values->get_update_flags() & update_JxW_values)
232 mapping_info_storage.JxW_values[q] = fe_values->JxW(q);
233 if (fe_values->get_update_flags() & update_jacobians)
234 {
235 Tensor<2, dim> jac = fe_values->jacobian(q);
236 jac = invert(transpose(jac));
237 for (unsigned int d = 0; d < dim; ++d)
238 for (unsigned int e = 0; e < dim; ++e)
239 mapping_info_storage.jacobians[0][q][d][e] = jac[d][e];
240 }
241 if (fe_values->get_update_flags() & update_quadrature_points)
242 for (unsigned int d = 0; d < dim; ++d)
243 mapping_info_storage.quadrature_points[q][d] =
244 fe_values->quadrature_point(q)[d];
245 if (fe_values->get_update_flags() & update_normal_vectors)
246 {
247 for (unsigned int d = 0; d < dim; ++d)
248 mapping_info_storage.normal_vectors[q][d] =
249 fe_values->normal_vector(q)[d];
250 mapping_info_storage.normals_times_jacobians[0][q] =
251 mapping_info_storage.normal_vectors[q] *
252 mapping_info_storage.jacobians[0][q];
253 }
254 }
255 }
256
257
258
259 template <int dim, typename Number>
260 inline bool
262 {
263 return present_cell !=
264 typename ::Triangulation<dim>::cell_iterator();
265 }
266
267
268
269 template <int dim, typename Number>
270 inline typename ::Triangulation<dim>::cell_iterator
272 {
273 return fe_values->get_cell();
274 }
275
276
277
278 template <int dim, typename Number>
279 inline const ::FEValues<dim> &
281 {
282 return *fe_values;
283 }
284
285
286
287 template <int dim, typename Number>
290 {
291 return mapping_info_storage;
292 }
293
294
295
296 template <int dim, typename Number>
299 {
300 return mapping_info_storage;
301 }
302
303
304
305 template <int dim, typename Number>
306 inline const Quadrature<1> &
308 {
309 return quadrature_1d;
310 }
311
312
313 } // end of namespace MatrixFreeFunctions
314} // end of namespace internal
315
316
318
319#endif
void resize(const size_type new_size)
Abstract base class for mapping classes.
Definition mapping.h:318
const MappingInfoStorage< dim, dim, Number > & get_data_storage() const
void reinit(typename ::Triangulation< dim >::cell_iterator cell)
MappingDataOnTheFly(const Quadrature< 1 > &quadrature, const UpdateFlags update_flags)
MappingInfoStorage< dim, dim, Number > & get_data_storage()
typename::Triangulation< dim >::cell_iterator present_cell
MappingDataOnTheFly(const Mapping< dim > &mapping, const Quadrature< 1 > &quadrature, const UpdateFlags update_flags)
MappingInfoStorage< dim, dim, Number > mapping_info_storage
typename::Triangulation< dim >::cell_iterator get_cell() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
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.
STL namespace.
AlignedVector< Tensor< 1, spacedim, Number > > normal_vectors
std::array< AlignedVector< Tensor< 1, spacedim, Number > >, 2 > normals_times_jacobians
AlignedVector< Point< spacedim, Number > > quadrature_points
std::array< AlignedVector< Tensor< 2, spacedim, Number > >, 2 > jacobians
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)