Reference documentation for deal.II version 9.0.0
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 at
12 // the top level of the deal.II distribution.
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 #include <deal.II/base/derivative_form.h>
22 #include <deal.II/base/table.h>
23 #include <deal.II/base/quadrature_lib.h>
24 #include <deal.II/base/qprojector.h>
25 #include <deal.II/grid/tria_iterator.h>
26 #include <deal.II/dofs/dof_accessor.h>
27 #include <deal.II/fe/mapping.h>
28 
29 #include <cmath>
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 template <int,int> 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
75  std::unique_ptr<Mapping<dim,spacedim>> clone () const override;
76 
81  virtual
82  bool preserves_vertex_locations () const override;
83 
89  // for documentation, see the Mapping base class
90  virtual
93  const Point<dim> &p) const override;
94 
95  // for documentation, see the Mapping base class
96  virtual
99  const Point<spacedim> &p) const override;
100 
110  // for documentation, see the Mapping base class
111  virtual
112  void
113  transform (const ArrayView<const Tensor<1,dim> > &input,
114  const MappingType type,
115  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
116  const ArrayView<Tensor<1,spacedim> > &output) const override;
117 
118  // for documentation, see the Mapping base class
119  virtual
120  void
122  const MappingType type,
123  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
124  const ArrayView<Tensor<2,spacedim> > &output) const override;
125 
126  // for documentation, see the Mapping base class
127  virtual
128  void
129  transform (const ArrayView<const Tensor<2, dim> > &input,
130  const MappingType type,
131  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
132  const ArrayView<Tensor<2,spacedim> > &output) const override;
133 
134  // for documentation, see the Mapping base class
135  virtual
136  void
138  const MappingType type,
139  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
140  const ArrayView<Tensor<3,spacedim> > &output) const override;
141 
142  // for documentation, see the Mapping base class
143  virtual
144  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 
220  virtual
221  std::size_t memory_consumption () const override;
222 
228  mutable std::vector<Point<spacedim> > vertices;
229 
236 
243 
244 
263  std::vector<std::vector<double> > cell_manifold_quadrature_weights;
264 
265 
280  mutable std::vector<double> vertex_weights;
281 
295  std::array<std::vector<Tensor<1,dim> >, GeometryInfo<dim>::faces_per_cell *(dim-1)> unit_tangentials;
296 
306  mutable std::vector<DerivativeForm<1,dim, spacedim > > covariant;
307 
315  mutable std::vector<DerivativeForm<1,dim,spacedim> > contravariant;
316 
320  mutable std::vector<std::vector<Tensor<1,spacedim> > > aux;
321 
326  mutable std::vector<double> volume_elements;
327 
334  };
335 
336 
337  // documentation can be found in Mapping::requires_update_flags()
338  virtual
340  requires_update_flags (const UpdateFlags update_flags) const override;
341 
342  // documentation can be found in Mapping::get_data()
343  virtual
344  std::unique_ptr<typename Mapping<dim,spacedim>::InternalDataBase>
345  get_data (const UpdateFlags,
346  const Quadrature<dim> &quadrature) const override;
347 
348  // documentation can be found in Mapping::get_face_data()
349  virtual
350  std::unique_ptr<typename Mapping<dim,spacedim>::InternalDataBase>
351  get_face_data (const UpdateFlags flags,
352  const Quadrature<dim-1>& quadrature) const override;
353 
354  // documentation can be found in Mapping::get_subface_data()
355  virtual
356  std::unique_ptr<typename Mapping<dim,spacedim>::InternalDataBase>
357  get_subface_data (const UpdateFlags flags,
358  const Quadrature<dim-1>& quadrature) const override;
359 
360  // documentation can be found in Mapping::fill_fe_values()
361  virtual
364  const CellSimilarity::Similarity cell_similarity,
365  const Quadrature<dim> &quadrature,
366  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
368 
369  // documentation can be found in Mapping::fill_fe_face_values()
370  virtual
371  void
373  const unsigned int face_no,
374  const Quadrature<dim-1> &quadrature,
375  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
377 
378  // documentation can be found in Mapping::fill_fe_subface_values()
379  virtual
380  void
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 
391 };
392 
393 
394 
397 /*----------------------------------------------------------------------*/
398 
399 #ifndef DOXYGEN
400 
401 template <int dim, int spacedim>
402 inline
403 void
405 {
407  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
408  vertices[i] = cell->vertex(i);
409  this->cell = cell;
410 }
411 
412 
413 template <int dim, int spacedim>
414 inline
415 void
417 {
418  cell_manifold_quadrature_weights.resize(quad.size(), std::vector<double>(GeometryInfo<dim>::vertices_per_cell));
419  for (unsigned int q=0; q<quad.size(); ++q)
420  {
421  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
422  {
423  cell_manifold_quadrature_weights[q][i] = GeometryInfo<dim>::d_linear_shape_function(quad.point(q), i);
424  }
425  }
426 }
427 
428 
429 
430 template <int dim, int spacedim>
431 inline
432 bool
434 {
435  return true;
436 }
437 
438 #endif // DOXYGEN
439 
440 /* -------------- declaration of explicit specializations ------------- */
441 
442 
443 DEAL_II_NAMESPACE_CLOSE
444 
445 #endif
virtual bool preserves_vertex_locations() 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
MappingType
Definition: mapping.h:51
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim-1)> unit_tangentials
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
const Point< dim > & point(const unsigned int i) const
void store_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
Triangulation< dim, spacedim >::cell_iterator cell
virtual std::size_t memory_consumption() const override
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
std::vector< std::vector< Tensor< 1, spacedim > > > aux
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:46
void initialize(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
std::vector< double > volume_elements
unsigned int size() const
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
std::vector< double > vertex_weights
void compute_manifold_quadrature_weights(const Quadrature< dim > &quadrature)
SmartPointer< const Manifold< dim, spacedim > > manifold
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
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 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
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
MappingManifold()=default