Reference documentation for deal.II version 8.5.1
mapping_data_on_the_fly.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2014 - 2016 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 at
12 // the top level of the deal.II distribution.
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 #include <deal.II/base/exceptions.h>
23 #include <deal.II/base/subscriptor.h>
24 #include <deal.II/base/vectorization.h>
25 #include <deal.II/base/aligned_vector.h>
26 #include <deal.II/matrix_free/shape_info.h>
27 #include <deal.II/matrix_free/mapping_info.h>
28 #include <deal.II/fe/fe_values.h>
29 #include <deal.II/fe/fe_nothing.h>
30 
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 
35 namespace internal
36 {
37  namespace MatrixFreeFunctions
38  {
57  template <int dim, typename Number=double>
59  {
60  public:
67  MappingDataOnTheFly (const Mapping<dim> &mapping,
68  const Quadrature<1> &quadrature,
69  const UpdateFlags update_flags);
70 
76  MappingDataOnTheFly (const Quadrature<1> &quadrature,
77  const UpdateFlags update_flags);
78 
82  void reinit(typename ::Triangulation<dim>::cell_iterator cell);
83 
88  bool is_initialized() const;
89 
93  typename ::Triangulation<dim>::cell_iterator get_cell () const;
94 
101  const ::FEValues<dim> &get_fe_values () const;
102 
109  get_inverse_jacobians() const;
110 
117  get_JxW_values() const;
118 
125  get_quadrature_points() const;
126 
133  get_normal_vectors() const;
134 
138  const Quadrature<1> &
139  get_quadrature () const;
140 
141  private:
146  typename ::Triangulation<dim>::cell_iterator present_cell;
147 
153 
158 
163 
168 
173 
178 
183  };
184 
185 
186  /*----------------------- Inline functions ----------------------------------*/
187 
188  template <int dim, typename Number>
189  inline
191  const Quadrature<1> &quadrature,
192  const UpdateFlags update_flags)
193  :
194  fe_values(mapping, fe_dummy, Quadrature<dim>(quadrature),
195  internal::MatrixFreeFunctions::MappingInfo<dim,Number>::compute_update_flags(update_flags)),
196  quadrature_1d(quadrature),
197  inverse_jacobians(fe_values.get_quadrature().size()),
198  jxw_values(fe_values.get_quadrature().size()),
199  quadrature_points(fe_values.get_quadrature().size()),
200  normal_vectors(fe_values.get_quadrature().size())
201  {
204  }
205 
206 
207 
208  template <int dim, typename Number>
209  inline
211  const UpdateFlags update_flags)
212  :
213  fe_values(fe_dummy, Quadrature<dim>(quadrature),
214  internal::MatrixFreeFunctions::MappingInfo<dim,Number>::compute_update_flags(update_flags)),
215  quadrature_1d(quadrature),
216  inverse_jacobians(fe_values.get_quadrature().size()),
217  jxw_values(fe_values.get_quadrature().size()),
218  quadrature_points(fe_values.get_quadrature().size()),
219  normal_vectors(fe_values.get_quadrature().size())
220  {
223  }
224 
225 
226 
227  template <int dim, typename Number>
228  inline
229  void
230  MappingDataOnTheFly<dim,Number>::reinit(typename ::Triangulation<dim>::cell_iterator cell)
231  {
232  if (present_cell == cell)
233  return;
234  present_cell = cell;
235  fe_values.reinit(present_cell);
236  for (unsigned int q=0; q<fe_values.get_quadrature().size(); ++q)
237  {
238  if (fe_values.get_update_flags() & update_inverse_jacobians)
239  for (unsigned int d=0; d<dim; ++d)
240  for (unsigned int e=0; e<dim; ++e)
241  inverse_jacobians[q][d][e] = fe_values.inverse_jacobian(q)[e][d];
242  if (fe_values.get_update_flags() & update_quadrature_points)
243  for (unsigned int d=0; d<dim; ++d)
244  quadrature_points[q][d] = fe_values.quadrature_point(q)[d];
245  if (fe_values.get_update_flags() & update_normal_vectors)
246  for (unsigned int d=0; d<dim; ++d)
247  normal_vectors[q][d] = fe_values.normal_vector(q)[d];
248  if (fe_values.get_update_flags() & update_JxW_values)
249  jxw_values[q] = fe_values.JxW(q);
250  }
251  }
252 
253 
254 
255  template <int dim, typename Number>
256  inline
257  bool
259  {
260  return present_cell != typename ::Triangulation<dim>::cell_iterator();
261  }
262 
263 
264 
265  template <int dim, typename Number>
266  inline
267  typename ::Triangulation<dim>::cell_iterator
269  {
270  return fe_values.get_cell();
271  }
272 
273 
274 
275  template <int dim, typename Number>
276  inline
277  const ::FEValues<dim> &
279  {
280  return fe_values;
281  }
282 
283 
284 
285  template <int dim, typename Number>
286  inline
289  {
290  return inverse_jacobians;
291  }
292 
293 
294 
295  template <int dim, typename Number>
296  inline
299  {
300  return normal_vectors;
301  }
302 
303 
304 
305  template <int dim, typename Number>
306  inline
309  {
310  return quadrature_points;
311  }
312 
313 
314 
315  template <int dim, typename Number>
316  inline
319  {
320  return jxw_values;
321  }
322 
323 
324 
325  template <int dim, typename Number>
326  inline
327  const Quadrature<1> &
329  {
330  return quadrature_1d;
331  }
332 
333 
334  } // end of namespace MatrixFreeFunctions
335 } // end of namespace internal
336 
337 
338 DEAL_II_NAMESPACE_CLOSE
339 
340 #endif
Transformed quadrature weights.
typename ::Triangulation< dim >::cell_iterator present_cell
const AlignedVector< VectorizedArray< Number > > & get_JxW_values() const
const AlignedVector< Tensor< 2, dim, VectorizedArray< Number > > > & get_inverse_jacobians() const
UpdateFlags get_update_flags() const
Transformed quadrature points.
typename ::Triangulation< dim >::cell_iterator get_cell() const
AlignedVector< Point< dim, VectorizedArray< Number > > > quadrature_points
AlignedVector< Tensor< 2, dim, VectorizedArray< Number > > > inverse_jacobians
#define Assert(cond, exc)
Definition: exceptions.h:313
UpdateFlags
const AlignedVector< Tensor< 1, dim, VectorizedArray< Number > > > & get_normal_vectors() const
AlignedVector< Tensor< 1, dim, VectorizedArray< Number > > > normal_vectors
AlignedVector< VectorizedArray< Number > > jxw_values
Gradient of volume element.
const AlignedVector< Point< dim, VectorizedArray< Number > > > & get_quadrature_points() const
Normal vectors.
static ::ExceptionBase & ExcNotImplemented()
void reinit(typename ::Triangulation< dim >::cell_iterator cell)
MappingDataOnTheFly(const Mapping< dim > &mapping, const Quadrature< 1 > &quadrature, const UpdateFlags update_flags)