Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
mapping_manifold.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2016 - 2018 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 #ifndef dealii_mapping_manifold_h
17 #define dealii_mapping_manifold_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #include <deal.II/base/array_view.h>
23 #include <deal.II/base/derivative_form.h>
24 #include <deal.II/base/quadrature.h>
25 
26 #include <deal.II/fe/mapping.h>
27 
28 #include <cmath>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 template <int, int>
33 class MappingQ;
34 
35 
38 
39 
59 template <int dim, int spacedim = dim>
60 class MappingManifold : public Mapping<dim, spacedim>
61 {
62 public:
66  MappingManifold() = default;
67 
72 
73  // for documentation, see the Mapping base class
74  virtual std::unique_ptr<Mapping<dim, spacedim>>
75  clone() const override;
76 
81  virtual bool
82  preserves_vertex_locations() const override;
83 
89  // for documentation, see the Mapping base class
90  virtual Point<spacedim>
93  const Point<dim> &p) const override;
94 
95  // for documentation, see the Mapping base class
96  virtual Point<dim>
99  const Point<spacedim> &p) const override;
100 
110  // for documentation, see the Mapping base class
111  virtual void
112  transform(const ArrayView<const Tensor<1, dim>> & input,
113  const MappingType type,
114  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
115  const ArrayView<Tensor<1, spacedim>> &output) const override;
116 
117  // for documentation, see the Mapping base class
118  virtual void
120  const MappingType type,
121  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
122  const ArrayView<Tensor<2, spacedim>> &output) const override;
123 
124  // for documentation, see the Mapping base class
125  virtual void
126  transform(const ArrayView<const Tensor<2, dim>> & input,
127  const MappingType type,
128  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
129  const ArrayView<Tensor<2, spacedim>> &output) const override;
130 
131  // for documentation, see the Mapping base class
132  virtual void
134  const MappingType type,
135  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
136  const ArrayView<Tensor<3, spacedim>> &output) const override;
137 
138  // for documentation, see the Mapping base class
139  virtual void
140  transform(const ArrayView<const Tensor<3, dim>> & input,
141  const MappingType type,
142  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
143  const ArrayView<Tensor<3, spacedim>> &output) const override;
144 
154 public:
166  class InternalData : public Mapping<dim, spacedim>::InternalDataBase
167  {
168  public:
172  InternalData() = default;
173 
182  void
183  initialize(const UpdateFlags update_flags,
184  const Quadrature<dim> &quadrature,
185  const unsigned int n_original_q_points);
186 
192  void
193  initialize_face(const UpdateFlags update_flags,
194  const Quadrature<dim> &quadrature,
195  const unsigned int n_original_q_points);
196 
197 
203  void
205 
209  void
211  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;
212 
216  virtual std::size_t
217  memory_consumption() const override;
218 
224  mutable std::vector<Point<spacedim>> vertices;
225 
232 
239 
240 
259  std::vector<std::vector<double>> cell_manifold_quadrature_weights;
260 
261 
276  mutable std::vector<double> vertex_weights;
277 
291  std::array<std::vector<Tensor<1, dim>>,
294 
304  mutable std::vector<DerivativeForm<1, dim, spacedim>> covariant;
305 
313  mutable std::vector<DerivativeForm<1, dim, spacedim>> contravariant;
314 
318  mutable std::vector<std::vector<Tensor<1, spacedim>>> aux;
319 
324  mutable std::vector<double> volume_elements;
325 
332  };
333 
334 
335  // documentation can be found in Mapping::requires_update_flags()
336  virtual UpdateFlags
337  requires_update_flags(const UpdateFlags update_flags) const override;
338 
339  // documentation can be found in Mapping::get_data()
340  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
341  get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
342 
343  // documentation can be found in Mapping::get_face_data()
344  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
345  get_face_data(const UpdateFlags flags,
346  const Quadrature<dim - 1> &quadrature) const override;
347 
348  // documentation can be found in Mapping::get_subface_data()
349  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
350  get_subface_data(const UpdateFlags flags,
351  const Quadrature<dim - 1> &quadrature) const override;
352 
353  // documentation can be found in Mapping::fill_fe_values()
356  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
357  const CellSimilarity::Similarity cell_similarity,
358  const Quadrature<dim> & quadrature,
359  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
361  &output_data) const override;
362 
363  // documentation can be found in Mapping::fill_fe_face_values()
364  virtual void
366  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
367  const unsigned int face_no,
368  const Quadrature<dim - 1> & quadrature,
369  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
371  &output_data) const override;
372 
373  // documentation can be found in Mapping::fill_fe_subface_values()
374  virtual void
376  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
377  const unsigned int face_no,
378  const unsigned int subface_no,
379  const Quadrature<dim - 1> & quadrature,
380  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
382  &output_data) const override;
383 
387 };
388 
389 
390 
393 /*----------------------------------------------------------------------*/
394 
395 #ifndef DOXYGEN
396 
397 template <int dim, int spacedim>
398 inline void
400  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
401 {
403  for (unsigned int i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
404  vertices[i] = cell->vertex(i);
405  this->cell = cell;
406 }
407 
408 
409 template <int dim, int spacedim>
410 inline void
413 {
414  cell_manifold_quadrature_weights.resize(
415  quad.size(), std::vector<double>(GeometryInfo<dim>::vertices_per_cell));
416  for (unsigned int q = 0; q < quad.size(); ++q)
417  {
418  for (unsigned int i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
419  {
420  cell_manifold_quadrature_weights[q][i] =
422  }
423  }
424 }
425 
426 
427 
428 template <int dim, int spacedim>
429 inline bool
431 {
432  return true;
433 }
434 
435 #endif // DOXYGEN
436 
437 /* -------------- declaration of explicit specializations ------------- */
438 
439 
440 DEAL_II_NAMESPACE_CLOSE
441 
442 #endif
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
virtual bool preserves_vertex_locations() const override
virtual void transform(const ArrayView< const Tensor< 1, dim >> &input, const MappingType type, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim >> &output) const override
virtual CellSimilarity::Similarity fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
SmartPointer< const Manifold< dim, spacedim > > manifold
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
MappingType
Definition: mapping.h:61
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim - 1)> unit_tangentials
const Point< dim > & point(const unsigned int i) const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
void store_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
virtual std::size_t memory_consumption() const override
Triangulation< dim, spacedim >::cell_iterator cell
std::vector< Point< spacedim > > vertices
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
std::vector< std::vector< Tensor< 1, spacedim > > > aux
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:57
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
void initialize(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
std::vector< double > volume_elements
unsigned int size() const
std::vector< double > vertex_weights
void compute_manifold_quadrature_weights(const Quadrature< dim > &quadrature)
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const Quadrature< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
void initialize_face(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
std::vector< std::vector< double > > cell_manifold_quadrature_weights
MappingManifold()=default