Reference documentation for deal.II version 9.2.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\}}\)
mapping_data_on_the_fly.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2014 - 2019 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 
28 #include <deal.II/fe/fe_nothing.h>
29 #include <deal.II/fe/fe_values.h>
30 #include <deal.II/fe/mapping_q1.h>
31 
34 
35 
37 
38 
39 namespace internal
40 {
41  namespace MatrixFreeFunctions
42  {
61  template <int dim, typename Number, typename VectorizedArrayType>
63  {
64  static_assert(
66  "Type of Number and of VectorizedArrayType do not match.");
67 
68  public:
75  MappingDataOnTheFly(const Mapping<dim> & mapping,
76  const Quadrature<1> &quadrature,
77  const UpdateFlags update_flags);
78 
84  MappingDataOnTheFly(const Quadrature<1> &quadrature,
85  const UpdateFlags update_flags);
86 
90  void
91  reinit(typename ::Triangulation<dim>::cell_iterator cell);
92 
97  bool
98  is_initialized() const;
99 
103  typename ::Triangulation<dim>::cell_iterator
104  get_cell() const;
105 
112  const ::FEValues<dim> &
113  get_fe_values() const;
114 
122  get_data_storage() const;
123 
127  const Quadrature<1> &
128  get_quadrature() const;
129 
130  private:
135  typename ::Triangulation<dim>::cell_iterator present_cell;
136 
142 
147 
152 
159  };
160 
161 
162  /*-------------------------- Inline functions ---------------------------*/
163 
164  template <int dim, typename Number, typename VectorizedArrayType>
167  const Quadrature<1> &quadrature,
168  const UpdateFlags update_flags)
169  : fe_values(
170  mapping,
171  fe_dummy,
172  Quadrature<dim>(quadrature),
173  MappingInfo<dim, Number, VectorizedArrayType>::compute_update_flags(
174  update_flags))
175  , quadrature_1d(quadrature)
176  {
178  mapping_info_storage.descriptor[0].initialize(quadrature);
182  if (update_flags & update_quadrature_points)
183  {
187  }
189  {
194  }
197  }
198 
199 
200 
201  template <int dim, typename Number, typename VectorizedArrayType>
204  const UpdateFlags update_flags)
205  : MappingDataOnTheFly(::StaticMappingQ1<dim, dim>::mapping,
206  quadrature,
207  update_flags)
208  {}
209 
210 
211 
212  template <int dim, typename Number, typename VectorizedArrayType>
213  inline void
215  typename ::Triangulation<dim>::cell_iterator cell)
216  {
217  if (present_cell == cell)
218  return;
219  present_cell = cell;
221  for (unsigned int q = 0; q < fe_values.get_quadrature().size(); ++q)
222  {
226  {
228  jac = invert(transpose(jac));
229  for (unsigned int d = 0; d < dim; ++d)
230  for (unsigned int e = 0; e < dim; ++e)
231  mapping_info_storage.jacobians[0][q][d][e] = jac[d][e];
232  }
234  for (unsigned int d = 0; d < dim; ++d)
238  {
239  for (unsigned int d = 0; d < dim; ++d)
245  }
246  }
247  }
248 
249 
250 
251  template <int dim, typename Number, typename VectorizedArrayType>
252  inline bool
254  const
255  {
256  return present_cell !=
257  typename ::Triangulation<dim>::cell_iterator();
258  }
259 
260 
261 
262  template <int dim, typename Number, typename VectorizedArrayType>
263  inline typename ::Triangulation<dim>::cell_iterator
265  {
266  return fe_values.get_cell();
267  }
268 
269 
270 
271  template <int dim, typename Number, typename VectorizedArrayType>
272  inline const ::FEValues<dim> &
274  {
275  return fe_values;
276  }
277 
278 
279 
280  template <int dim, typename Number, typename VectorizedArrayType>
283  const
284  {
285  return mapping_info_storage;
286  }
287 
288 
289 
290  template <int dim, typename Number, typename VectorizedArrayType>
291  inline const Quadrature<1> &
293  const
294  {
295  return quadrature_1d;
296  }
297 
298 
299  } // end of namespace MatrixFreeFunctions
300 } // end of namespace internal
301 
302 
304 
305 #endif
Transformed quadrature weights.
typename ::Triangulation< dim >::cell_iterator get_cell() const
MappingInfoStorage< dim, dim, Number, VectorizedArrayType > mapping_info_storage
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
AlignedVector< unsigned int > data_index_offsets
Definition: mapping_info.h:196
Volume element.
AlignedVector< Tensor< 2, spacedim, VectorizedArrayType > > jacobians[2]
Definition: mapping_info.h:225
UpdateFlags get_update_flags() const
Transformed quadrature points.
const DerivativeForm< 1, dim, spacedim > & jacobian(const unsigned int quadrature_point) const
void resize(const size_type size_in)
typename ::Triangulation< dim >::cell_iterator present_cell
const Point< spacedim > & quadrature_point(const unsigned int q) const
AlignedVector< Point< spacedim, VectorizedArrayType > > quadrature_points
Definition: mapping_info.h:272
AlignedVector< unsigned int > quadrature_point_offsets
Definition: mapping_info.h:263
std::vector< QuadratureDescriptor > descriptor
Definition: mapping_info.h:188
#define Assert(cond, exc)
Definition: exceptions.h:1419
UpdateFlags
const MappingInfoStorage< dim, dim, Number, VectorizedArrayType > & get_data_storage() const
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
Gradient of volume element.
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access >> &cell)
Definition: fe_values.cc:4533
const Quadrature< dim > & get_quadrature() const
const unsigned int n_quadrature_points
Definition: fe_values.h:2101
void reinit(typename ::Triangulation< dim >::cell_iterator cell)
MappingDataOnTheFly(const Mapping< dim > &mapping, const Quadrature< 1 > &quadrature, const UpdateFlags update_flags)
double JxW(const unsigned int quadrature_point) const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
AlignedVector< Tensor< 1, spacedim, VectorizedArrayType > > normals_times_jacobians[2]
Definition: mapping_info.h:255
AlignedVector< VectorizedArrayType > JxW_values
Definition: mapping_info.h:205
Normal vectors.
static ::ExceptionBase & ExcNotImplemented()
AlignedVector< Tensor< 1, spacedim, VectorizedArrayType > > normal_vectors
Definition: mapping_info.h:212
const Triangulation< dim, spacedim >::cell_iterator get_cell() const
Definition: fe_values.cc:4190
constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
static const bool value
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const