Reference documentation for deal.II version Git 53084d8c6a 2020-10-21 00:34:33 -0400
\(\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 - 2020 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 
25 
26 template <int dim, int spacedim>
27 std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
29  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
30 {
31  std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell> vertices;
32  for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
33  {
34  vertices[i] = cell->vertex(i);
35  }
36  return vertices;
37 }
38 
39 
40 
41 template <int dim, int spacedim>
45  const bool map_center_of_reference_cell) const
46 {
47  if (map_center_of_reference_cell)
48  {
49  Point<dim> reference_center;
50  for (unsigned int d = 0; d < dim; ++d)
51  reference_center[d] = .5;
52  return transform_unit_to_real_cell(cell, reference_center);
53  }
54  else
55  {
56  const auto vertices = get_vertices(cell);
58  for (const auto &v : vertices)
59  center += v;
61  }
62 }
63 
64 
65 
66 template <int dim, int spacedim>
69  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
70 {
71  if (preserves_vertex_locations())
72  return cell->bounding_box();
73  else
74  return BoundingBox<spacedim>(get_vertices(cell));
75 }
76 
77 
78 
79 template <int dim, int spacedim>
80 void
83  const ArrayView<const Point<spacedim>> & real_points,
84  const ArrayView<Point<dim>> & unit_points) const
85 {
86  AssertDimension(real_points.size(), unit_points.size());
87  for (unsigned int i = 0; i < real_points.size(); ++i)
88  {
89  try
90  {
91  unit_points[i] = transform_real_to_unit_cell(cell, real_points[i]);
92  }
93  catch (typename Mapping<dim>::ExcTransformationFailed &)
94  {
95  unit_points[i] = Point<dim>();
96  unit_points[i][0] = std::numeric_limits<double>::infinity();
97  }
98  }
99 }
100 
101 
102 
103 template <int dim, int spacedim>
104 Point<dim - 1>
106  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
107  const unsigned int face_no,
108  const Point<spacedim> & p) const
109 {
110  // The function doesn't make physical sense for dim=1
111  Assert(dim > 1, ExcNotImplemented());
112  // Not implemented for higher dimensions
113  Assert(dim <= 3, ExcNotImplemented());
114 
115  Point<dim> unit_cell_pt = transform_real_to_unit_cell(cell, p);
116 
117  const unsigned int unit_normal_direction =
119 
120  if (dim == 2)
121  {
122  if (unit_normal_direction == 0)
123  return Point<dim - 1>{unit_cell_pt(1)};
124  else if (unit_normal_direction == 1)
125  return Point<dim - 1>{unit_cell_pt(0)};
126  }
127  else if (dim == 3)
128  {
129  if (unit_normal_direction == 0)
130  return Point<dim - 1>{unit_cell_pt(1), unit_cell_pt(2)};
131  else if (unit_normal_direction == 1)
132  return Point<dim - 1>{unit_cell_pt(0), unit_cell_pt(2)};
133  else if (unit_normal_direction == 2)
134  return Point<dim - 1>{unit_cell_pt(0), unit_cell_pt(1)};
135  }
136 
137  // We should never get here
138  Assert(false, ExcInternalError());
139  return {};
140 }
141 
142 /* ---------------------------- InternalDataBase --------------------------- */
143 
144 
145 template <int dim, int spacedim>
147  : update_each(update_default)
148 {}
149 
150 
151 
152 template <int dim, int spacedim>
153 std::size_t
155 {
156  return sizeof(*this);
157 }
158 
159 
160 /*------------------------------ InternalData ------------------------------*/
161 
162 
163 
164 // explicit instantiations
165 #include "mapping.inst"
166 
167 
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:105
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1580
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:43
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:81
virtual std::size_t memory_consumption() const
Definition: mapping.cc:154
virtual std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
Definition: mapping.cc:28
No update.
#define Assert(cond, exc)
Definition: exceptions.h:1423
Abstract base class for mapping classes.
Definition: mapping.h:301
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:369
Point< 3 > vertices[4]
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
virtual BoundingBox< spacedim > get_bounding_box(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
Definition: mapping.cc:68
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:368
Point< 3 > center
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcInternalError()