Reference documentation for deal.II version 9.3.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.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 2021 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 
18 
19 #include <deal.II/fe/mapping.h>
20 
21 #include <deal.II/grid/tria.h>
22 
24 #include <boost/geometry.hpp>
26 
28 
29 
30 template <int dim, int spacedim>
31 boost::container::small_vector<Point<spacedim>,
34  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
35 {
36  boost::container::small_vector<Point<spacedim>,
38  vertices;
39  for (const unsigned int i : cell->vertex_indices())
40  vertices.push_back(cell->vertex(i));
41 
42  return vertices;
43 }
44 
45 
46 
47 template <int dim, int spacedim>
51  const bool map_center_of_reference_cell) const
52 {
53  if (map_center_of_reference_cell)
54  {
55  Point<dim> reference_center;
56  for (unsigned int d = 0; d < dim; ++d)
57  reference_center[d] = .5;
58  return transform_unit_to_real_cell(cell, reference_center);
59  }
60  else
61  {
62  const auto vertices = get_vertices(cell);
64  for (const auto &v : vertices)
65  center += v;
67  }
68 }
69 
70 
71 
72 template <int dim, int spacedim>
75  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
76 {
77  if (preserves_vertex_locations())
78  return cell->bounding_box();
79  else
80  return BoundingBox<spacedim>(get_vertices(cell));
81 }
82 
83 
84 
85 template <int dim, int spacedim>
86 void
89  const ArrayView<const Point<spacedim>> & real_points,
90  const ArrayView<Point<dim>> & unit_points) const
91 {
92  AssertDimension(real_points.size(), unit_points.size());
93  for (unsigned int i = 0; i < real_points.size(); ++i)
94  {
95  try
96  {
97  unit_points[i] = transform_real_to_unit_cell(cell, real_points[i]);
98  }
99  catch (typename Mapping<dim>::ExcTransformationFailed &)
100  {
101  unit_points[i] = Point<dim>();
102  unit_points[i][0] = std::numeric_limits<double>::infinity();
103  }
104  }
105 }
106 
107 
108 
109 template <int dim, int spacedim>
110 Point<dim - 1>
112  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
113  const unsigned int face_no,
114  const Point<spacedim> & p) const
115 {
116  // The function doesn't make physical sense for dim=1
117  Assert(dim > 1, ExcNotImplemented());
118  // Not implemented for higher dimensions
119  Assert(dim <= 3, ExcNotImplemented());
120 
121  Point<dim> unit_cell_pt = transform_real_to_unit_cell(cell, p);
122 
123  const unsigned int unit_normal_direction =
125 
126  if (dim == 2)
127  {
128  if (unit_normal_direction == 0)
129  return Point<dim - 1>{unit_cell_pt(1)};
130  else if (unit_normal_direction == 1)
131  return Point<dim - 1>{unit_cell_pt(0)};
132  }
133  else if (dim == 3)
134  {
135  if (unit_normal_direction == 0)
136  return Point<dim - 1>{unit_cell_pt(1), unit_cell_pt(2)};
137  else if (unit_normal_direction == 1)
138  return Point<dim - 1>{unit_cell_pt(0), unit_cell_pt(2)};
139  else if (unit_normal_direction == 2)
140  return Point<dim - 1>{unit_cell_pt(0), unit_cell_pt(1)};
141  }
142 
143  // We should never get here
144  Assert(false, ExcInternalError());
145  return {};
146 }
147 
148 
149 
150 #ifndef DOXYGEN
151 template <int dim, int spacedim>
152 void
154  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
155  const unsigned int face_no,
156  const hp::QCollection<dim - 1> & quadrature,
157  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
159  &output_data) const
160 {
161  // base class version, implement overridden function in derived classes
162  AssertDimension(quadrature.size(), 1);
163  fill_fe_face_values(cell, face_no, quadrature[0], internal_data, output_data);
164 }
165 
166 
167 
168 template <int dim, int spacedim>
169 void
171  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
172  const unsigned int face_no,
173  const Quadrature<dim - 1> & quadrature,
174  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
176  &output_data) const
177 {
178  Assert(false,
179  ExcMessage("Use of a deprecated interface, please implement "
180  "fill_fe_face_values taking a hp::QCollection argument"));
181  (void)cell;
182  (void)face_no;
183  (void)quadrature;
184  (void)internal_data;
185  (void)output_data;
186 }
187 
188 
189 
190 template <int dim, int spacedim>
191 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
193  const UpdateFlags update_flags,
194  const hp::QCollection<dim - 1> &quadrature) const
195 {
196  // base class version, implement overridden function in derived classes
197  return get_face_data(update_flags, quadrature[0]);
198 }
199 
200 
201 
202 template <int dim, int spacedim>
203 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
205  const UpdateFlags update_flags,
206  const Quadrature<dim - 1> &quadrature) const
207 {
208  Assert(false,
209  ExcMessage("Use of a deprecated interface, please implement "
210  "fill_fe_face_values taking a hp::QCollection argument"));
211  (void)update_flags;
212  (void)quadrature;
213 
214  return std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>();
215 }
216 #endif
217 
218 /* ---------------------------- InternalDataBase --------------------------- */
219 
220 
221 template <int dim, int spacedim>
223  : update_each(update_default)
224 {}
225 
226 
227 
228 template <int dim, int spacedim>
229 std::size_t
231 {
232  return sizeof(*this);
233 }
234 
235 
236 /* ------------------------------ Global functions ------------------------- */
237 
238 template <int dim, int spacedim>
241 {
242  const auto &reference_cells = triangulation.get_reference_cells();
243  Assert(reference_cells.size() == 1,
244  ExcMessage(
245  "This function can only work for triangulations that "
246  "use only a single cell type -- for example, only triangles "
247  "or only quadrilaterals. For mixed meshes, there is no "
248  "single linear mapping object that can be used for all "
249  "cells of the triangulation. The triangulation you are "
250  "passing to this function uses multiple cell types."));
251 
252  return reference_cells.front()
253  .template get_default_linear_mapping<dim, spacedim>();
254 }
255 
256 
257 
258 // explicit instantiations
259 #include "mapping.inst"
260 
261 
Point< dim - 1 > project_real_point_to_unit_point_on_face(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Point< spacedim > &p) const
Definition: mapping.cc:111
Point< 1 > transform_real_to_unit_cell(const std::array< Point< spacedim >, GeometryInfo< 1 >::vertices_per_cell > &vertices, const Point< spacedim > &p)
const std::vector< ReferenceCell > & get_reference_cells() const
Definition: tria.cc:13503
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
unsigned int size() const
Definition: collection.h:109
virtual Point< spacedim > get_center(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const bool map_center_of_reference_cell=true) const
Definition: mapping.cc:49
virtual void transform_points_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< spacedim >> &real_points, const ArrayView< Point< dim >> &unit_points) const
Definition: mapping.cc:87
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition: mapping.cc:240
virtual std::size_t memory_consumption() const
Definition: mapping.cc:230
virtual std::unique_ptr< InternalDataBase > get_face_data(const UpdateFlags update_flags, const hp::QCollection< dim - 1 > &quadrature) const
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition: config.h:408
static ::ExceptionBase & ExcMessage(std::string arg1)
No update.
#define Assert(cond, exc)
Definition: exceptions.h:1465
UpdateFlags
Abstract base class for mapping classes.
Definition: mapping.h:303
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:395
Point< 3 > vertices[4]
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
virtual BoundingBox< spacedim > get_bounding_box(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
Definition: mapping.cc:74
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition: config.h:445
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:394
Point< 3 > center
virtual boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
Definition: mapping.cc:33
static ::ExceptionBase & ExcNotImplemented()
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
static ::ExceptionBase & ExcInternalError()