Reference documentation for deal.II version 9.1.0-pre
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/derivative_form.h>
23 #include <deal.II/base/qprojector.h>
24 #include <deal.II/base/quadrature_lib.h>
25 #include <deal.II/base/table.h>
26 
27 #include <deal.II/dofs/dof_accessor.h>
28 
29 #include <deal.II/fe/mapping.h>
30 
31 #include <deal.II/grid/tria_iterator.h>
32 
33 #include <cmath>
34 
35 DEAL_II_NAMESPACE_OPEN
36 
37 template <int, int>
38 class MappingQ;
39 
40 
43 
44 
64 template <int dim, int spacedim = dim>
65 class MappingManifold : public Mapping<dim, spacedim>
66 {
67 public:
71  MappingManifold() = default;
72 
77 
78  // for documentation, see the Mapping base class
79  virtual std::unique_ptr<Mapping<dim, spacedim>>
80  clone() const override;
81 
86  virtual bool
87  preserves_vertex_locations() const override;
88 
94  // for documentation, see the Mapping base class
95  virtual Point<spacedim>
98  const Point<dim> &p) const override;
99 
100  // for documentation, see the Mapping base class
101  virtual Point<dim>
103  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
104  const Point<spacedim> &p) const override;
105 
115  // for documentation, see the Mapping base class
116  virtual void
117  transform(const ArrayView<const Tensor<1, dim>> & input,
118  const MappingType type,
119  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
120  const ArrayView<Tensor<1, spacedim>> &output) const override;
121 
122  // for documentation, see the Mapping base class
123  virtual void
125  const MappingType type,
126  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
127  const ArrayView<Tensor<2, spacedim>> &output) const override;
128 
129  // for documentation, see the Mapping base class
130  virtual void
131  transform(const ArrayView<const Tensor<2, dim>> & input,
132  const MappingType type,
133  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
134  const ArrayView<Tensor<2, spacedim>> &output) const override;
135 
136  // for documentation, see the Mapping base class
137  virtual void
139  const MappingType type,
140  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
141  const ArrayView<Tensor<3, spacedim>> &output) const override;
142 
143  // for documentation, see the Mapping base class
144  virtual void
145  transform(const ArrayView<const Tensor<3, dim>> & input,
146  const MappingType type,
147  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
148  const ArrayView<Tensor<3, spacedim>> &output) const override;
149 
159 public:
171  class InternalData : public Mapping<dim, spacedim>::InternalDataBase
172  {
173  public:
177  InternalData() = default;
178 
187  void
188  initialize(const UpdateFlags update_flags,
189  const Quadrature<dim> &quadrature,
190  const unsigned int n_original_q_points);
191 
197  void
198  initialize_face(const UpdateFlags update_flags,
199  const Quadrature<dim> &quadrature,
200  const unsigned int n_original_q_points);
201 
202 
208  void
210 
214  void
216  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;
217 
221  virtual std::size_t
222  memory_consumption() const override;
223 
229  mutable std::vector<Point<spacedim>> vertices;
230 
237 
244 
245 
264  std::vector<std::vector<double>> cell_manifold_quadrature_weights;
265 
266 
281  mutable std::vector<double> vertex_weights;
282 
296  std::array<std::vector<Tensor<1, dim>>,
299 
309  mutable std::vector<DerivativeForm<1, dim, spacedim>> covariant;
310 
318  mutable std::vector<DerivativeForm<1, dim, spacedim>> contravariant;
319 
323  mutable std::vector<std::vector<Tensor<1, spacedim>>> aux;
324 
329  mutable std::vector<double> volume_elements;
330 
337  };
338 
339 
340  // documentation can be found in Mapping::requires_update_flags()
341  virtual UpdateFlags
342  requires_update_flags(const UpdateFlags update_flags) const override;
343 
344  // documentation can be found in Mapping::get_data()
345  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
346  get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
347 
348  // documentation can be found in Mapping::get_face_data()
349  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
350  get_face_data(const UpdateFlags flags,
351  const Quadrature<dim - 1> &quadrature) const override;
352 
353  // documentation can be found in Mapping::get_subface_data()
354  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
355  get_subface_data(const UpdateFlags flags,
356  const Quadrature<dim - 1> &quadrature) const override;
357 
358  // documentation can be found in Mapping::fill_fe_values()
361  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
362  const CellSimilarity::Similarity cell_similarity,
363  const Quadrature<dim> & quadrature,
364  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
366  &output_data) const override;
367 
368  // documentation can be found in Mapping::fill_fe_face_values()
369  virtual void
371  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
372  const unsigned int face_no,
373  const Quadrature<dim - 1> & quadrature,
374  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
376  &output_data) const override;
377 
378  // documentation can be found in Mapping::fill_fe_subface_values()
379  virtual void
381  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
382  const unsigned int face_no,
383  const unsigned int subface_no,
384  const Quadrature<dim - 1> & quadrature,
385  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
387  &output_data) const override;
388 
392 };
393 
394 
395 
398 /*----------------------------------------------------------------------*/
399 
400 #ifndef DOXYGEN
401 
402 template <int dim, int spacedim>
403 inline void
406 {
408  for (unsigned int i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
409  vertices[i] = cell->vertex(i);
410  this->cell = cell;
411 }
412 
413 
414 template <int dim, int spacedim>
415 inline void
418 {
420  quad.size(), std::vector<double>(GeometryInfo<dim>::vertices_per_cell));
421  for (unsigned int q = 0; q < quad.size(); ++q)
422  {
423  for (unsigned int i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
424  {
427  }
428  }
429 }
430 
431 
432 
433 template <int dim, int spacedim>
434 inline bool
436 {
437  return true;
438 }
439 
440 #endif // DOXYGEN
441 
442 /* -------------- declaration of explicit specializations ------------- */
443 
444 
445 DEAL_II_NAMESPACE_CLOSE
446 
447 #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 std::unique_ptr< Mapping< dim, spacedim > > clone() 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
SmartPointer< const Manifold< dim, spacedim > > manifold
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
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
virtual std::size_t memory_consumption() const override
Triangulation< dim, spacedim >::cell_iterator cell
const Point< dim > & point(const unsigned int i) const
std::vector< Point< spacedim > > vertices
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
void store_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
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
std::vector< std::vector< Tensor< 1, spacedim > > > aux
unsigned int size() const
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:57
void initialize(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
std::vector< double > volume_elements
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 Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) 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)
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
std::vector< std::vector< double > > cell_manifold_quadrature_weights
MappingManifold()=default