Reference documentation for deal.II version GIT relicensing-439-g5fda5c893d 2024-04-20 06:50:02+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
mapping.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2001 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
17
18#include <deal.II/fe/mapping.h>
19
20#include <deal.II/grid/tria.h>
21
22#ifdef DEAL_II_BOOST_HAS_BROKEN_HEADER_DEPRECATIONS
23# define BOOST_ALLOW_DEPRECATED_HEADERS
24#endif
25#include <boost/geometry.hpp>
26#ifdef DEAL_II_BOOST_HAS_BROKEN_HEADER_DEPRECATIONS
27# undef BOOST_ALLOW_DEPRECATED_HEADERS
28#endif
29
30#include <limits>
31
33#ifndef DOXYGEN
34
35template <int dim, int spacedim>
36boost::container::small_vector<Point<spacedim>,
39 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
40{
41 boost::container::small_vector<Point<spacedim>,
44 for (const unsigned int i : cell->vertex_indices())
45 vertices.push_back(cell->vertex(i));
46
47 return vertices;
48}
49
50
51
52template <int dim, int spacedim>
53boost::container::small_vector<Point<spacedim>,
57 const unsigned int face_no) const
58{
59 boost::container::small_vector<Point<spacedim>,
61 face_vertices;
62
63 const auto &cell_vertices = get_vertices(cell);
64 const auto &reference_cell = cell->reference_cell();
65 const auto face_orientation = cell->combined_face_orientation(face_no);
66
67 for (const unsigned int v :
68 reference_cell.face_reference_cell(face_no).vertex_indices())
69 {
70 face_vertices.push_back(
71 cell_vertices[reference_cell.face_to_cell_vertices(
72 face_no, v, face_orientation)]);
73 }
74
75 return face_vertices;
76}
77
78
79
80template <int dim, int spacedim>
84 const bool map_barycenter_of_reference_cell) const
85{
86 if (map_barycenter_of_reference_cell)
87 {
88 return transform_unit_to_real_cell(
89 cell, cell->reference_cell().template barycenter<dim>());
90 }
91 else
92 {
93 const auto vertices = get_vertices(cell);
95 for (const auto &v : vertices)
96 center += v;
97 return center / cell->n_vertices();
98 }
99}
100
101
102
103template <int dim, int spacedim>
106 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
107{
108 if (preserves_vertex_locations())
109 return cell->bounding_box();
110 else
111 return BoundingBox<spacedim>(get_vertices(cell));
112}
113
114
115
116template <int dim, int spacedim>
117void
123{
125}
126
127
128
129template <int dim, int spacedim>
130void
133 const ArrayView<const Point<spacedim>> &real_points,
134 const ArrayView<Point<dim>> &unit_points) const
135{
136 AssertDimension(real_points.size(), unit_points.size());
137 for (unsigned int i = 0; i < real_points.size(); ++i)
138 {
139 try
140 {
141 unit_points[i] = transform_real_to_unit_cell(cell, real_points[i]);
142 }
144 {
145 unit_points[i] = Point<dim>();
146 unit_points[i][0] = std::numeric_limits<double>::infinity();
147 }
148 }
149}
150
151
152
153template <int dim, int spacedim>
154Point<dim - 1>
157 const unsigned int face_no,
158 const Point<spacedim> &p) const
159{
160 // The function doesn't make physical sense for dim=1
161 Assert(dim > 1, ExcNotImplemented());
162 // Not implemented for higher dimensions
163 Assert(dim <= 3, ExcNotImplemented());
164
165 Point<dim> unit_cell_pt = transform_real_to_unit_cell(cell, p);
166
167 const unsigned int unit_normal_direction =
169
170 if (dim == 2)
171 {
172 if (unit_normal_direction == 0)
173 return Point<dim - 1>{unit_cell_pt[1]};
174 else if (unit_normal_direction == 1)
175 return Point<dim - 1>{unit_cell_pt[0]};
176 }
177 else if (dim == 3)
178 {
179 if (unit_normal_direction == 0)
180 return Point<dim - 1>{unit_cell_pt[1], unit_cell_pt[2]};
181 else if (unit_normal_direction == 1)
182 return Point<dim - 1>{unit_cell_pt[0], unit_cell_pt[2]};
183 else if (unit_normal_direction == 2)
184 return Point<dim - 1>{unit_cell_pt[0], unit_cell_pt[1]};
185 }
186
187 // We should never get here
189 return {};
190}
191
192
193
194# ifndef DOXYGEN
195template <int dim, int spacedim>
196void
199 const unsigned int face_no,
200 const hp::QCollection<dim - 1> &quadrature,
201 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
203 &output_data) const
204{
205 // base class version, implement overridden function in derived classes
206 AssertDimension(quadrature.size(), 1);
207 fill_fe_face_values(cell, face_no, quadrature[0], internal_data, output_data);
208}
209
210
211
212template <int dim, int spacedim>
213void
216 const unsigned int face_no,
217 const Quadrature<dim - 1> &quadrature,
218 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
220 &output_data) const
221{
222 Assert(false,
223 ExcMessage("Use of a deprecated interface, please implement "
224 "fill_fe_face_values taking a hp::QCollection argument"));
225 (void)cell;
226 (void)face_no;
227 (void)quadrature;
228 (void)internal_data;
229 (void)output_data;
230}
231
232
233
234template <int dim, int spacedim>
235std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
237 const UpdateFlags update_flags,
238 const hp::QCollection<dim - 1> &quadrature) const
239{
240 // base class version, implement overridden function in derived classes
241 return get_face_data(update_flags, quadrature[0]);
242}
243
244
245
246template <int dim, int spacedim>
247std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
249 const UpdateFlags update_flags,
250 const Quadrature<dim - 1> &quadrature) const
251{
252 Assert(false,
253 ExcMessage("Use of a deprecated interface, please implement "
254 "fill_fe_face_values taking a hp::QCollection argument"));
255 (void)update_flags;
256 (void)quadrature;
257
258 return std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>();
259}
260# endif
261
262/* ---------------------------- InternalDataBase --------------------------- */
263
264
265template <int dim, int spacedim>
267 : update_each(update_default)
268{}
269
270
271
272template <int dim, int spacedim>
273void
275 const Quadrature<dim> &)
276{
278}
279
280
281
282template <int dim, int spacedim>
283std::size_t
285{
286 return sizeof(*this);
287}
288#endif
289
290/* ------------------------------ Global functions ------------------------- */
291
292template <int dim, int spacedim>
295{
296 const auto &reference_cells = triangulation.get_reference_cells();
297 Assert(reference_cells.size() == 1,
299 "This function can only work for triangulations that "
300 "use only a single cell type -- for example, only triangles "
301 "or only quadrilaterals. For mixed meshes, there is no "
302 "single linear mapping object that can be used for all "
303 "cells of the triangulation. The triangulation you are "
304 "passing to this function uses multiple cell types."));
305
306 return reference_cells.front()
307 .template get_default_linear_mapping<dim, spacedim>();
308}
309
310
311
312// explicit instantiations
313#include "mapping.inst"
314
315
virtual std::size_t memory_consumption() const
virtual void reinit(const UpdateFlags update_flags, const Quadrature< dim > &quadrature)
Abstract base class for mapping classes.
Definition mapping.h:318
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
virtual BoundingBox< spacedim > get_bounding_box(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
virtual void fill_fe_immersed_surface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const NonMatching::ImmersedSurfaceQuadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
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 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
virtual boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
virtual Point< spacedim > get_center(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const bool map_barycenter_of_reference_cell=true) const
virtual std::unique_ptr< InternalDataBase > get_face_data(const UpdateFlags update_flags, const hp::QCollection< dim - 1 > &quadrature) const
Definition point.h:111
unsigned int size() const
Definition collection.h:264
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
Point< 3 > center
Point< 3 > vertices[4]
unsigned int vertex_indices[2]
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
UpdateFlags
@ update_default
No update.
#define DEAL_II_ASSERT_UNREACHABLE()
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition mapping.cc:294
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
Point< 1 > transform_real_to_unit_cell(const std::array< Point< spacedim >, GeometryInfo< 1 >::vertices_per_cell > &vertices, const Point< spacedim > &p)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation