deal.II version GIT relicensing-2392-g6e04ac39aa 2025-01-20 09:10:00+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>,
37# ifndef _MSC_VER
38 ReferenceCells::max_n_vertices<dim>()
39# else
41# endif
42 >
44 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
45{
46 boost::container::small_vector<Point<spacedim>,
47# ifndef _MSC_VER
48 ReferenceCells::max_n_vertices<dim>()
49# else
51# endif
52 >
53 vertices;
54 for (const unsigned int i : cell->vertex_indices())
55 vertices.push_back(cell->vertex(i));
56
57 return vertices;
58}
59
60
61
62template <int dim, int spacedim>
63boost::container::small_vector<Point<spacedim>,
64# ifndef _MSC_VER
66# else
68# endif
69 >
72 const unsigned int face_no) const
73{
74 boost::container::small_vector<Point<spacedim>,
75# ifndef _MSC_VER
77# else
79# endif
80 >
81 face_vertices;
82
83 const auto &cell_vertices = get_vertices(cell);
84 const auto &reference_cell = cell->reference_cell();
85 const auto face_orientation = cell->combined_face_orientation(face_no);
86
87 for (const unsigned int v :
88 reference_cell.face_reference_cell(face_no).vertex_indices())
89 {
90 face_vertices.push_back(
91 cell_vertices[reference_cell.face_to_cell_vertices(
92 face_no, v, face_orientation)]);
93 }
94
95 return face_vertices;
96}
97
98
99
100template <int dim, int spacedim>
104 const bool map_barycenter_of_reference_cell) const
105{
106 if (map_barycenter_of_reference_cell)
107 {
108 return transform_unit_to_real_cell(
109 cell, cell->reference_cell().template barycenter<dim>());
110 }
111 else
112 {
113 const auto vertices = get_vertices(cell);
114 Point<spacedim> center;
115 for (const auto &v : vertices)
116 center += v;
117 return center / cell->n_vertices();
118 }
119}
120
121
122
123template <int dim, int spacedim>
126 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
127{
128 if (preserves_vertex_locations())
129 return cell->bounding_box();
130 else
131 return BoundingBox<spacedim>(get_vertices(cell));
132}
133
134
135
136template <int dim, int spacedim>
137void
143{
145}
146
147
148
149template <int dim, int spacedim>
150void
153 const ArrayView<const Point<spacedim>> &real_points,
154 const ArrayView<Point<dim>> &unit_points) const
155{
156 AssertDimension(real_points.size(), unit_points.size());
157 for (unsigned int i = 0; i < real_points.size(); ++i)
158 {
159 try
160 {
161 unit_points[i] = transform_real_to_unit_cell(cell, real_points[i]);
162 }
164 {
165 // If the transformation for this one point failed, mark it
166 // as invalid as described in the documentation.
167 unit_points[i] = Point<dim>();
168 unit_points[i][0] = std::numeric_limits<double>::lowest();
169 }
170 }
171}
172
173
174
175template <int dim, int spacedim>
176Point<dim - 1>
179 const unsigned int face_no,
180 const Point<spacedim> &p) const
181{
182 // The function doesn't make physical sense for dim=1
183 Assert(dim > 1, ExcNotImplemented());
184 // Not implemented for higher dimensions
185 Assert(dim <= 3, ExcNotImplemented());
186
187 Point<dim> unit_cell_pt = transform_real_to_unit_cell(cell, p);
188
189 const unsigned int unit_normal_direction =
191
192 if (dim == 2)
193 {
194 if (unit_normal_direction == 0)
195 return Point<dim - 1>{unit_cell_pt[1]};
196 else if (unit_normal_direction == 1)
197 return Point<dim - 1>{unit_cell_pt[0]};
198 }
199 else if (dim == 3)
200 {
201 if (unit_normal_direction == 0)
202 return Point<dim - 1>{unit_cell_pt[1], unit_cell_pt[2]};
203 else if (unit_normal_direction == 1)
204 return Point<dim - 1>{unit_cell_pt[0], unit_cell_pt[2]};
205 else if (unit_normal_direction == 2)
206 return Point<dim - 1>{unit_cell_pt[0], unit_cell_pt[1]};
207 }
208
209 // We should never get here
211 return {};
212}
213
214
215
216# ifndef DOXYGEN
217template <int dim, int spacedim>
218void
221 const unsigned int face_no,
222 const hp::QCollection<dim - 1> &quadrature,
223 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
225 &output_data) const
226{
227 // base class version, implement overridden function in derived classes
228 AssertDimension(quadrature.size(), 1);
229 fill_fe_face_values(cell, face_no, quadrature[0], internal_data, output_data);
230}
231
232
233
234template <int dim, int spacedim>
235void
238 const unsigned int face_no,
239 const Quadrature<dim - 1> &quadrature,
240 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
242 &output_data) const
243{
244 Assert(false,
245 ExcMessage("Use of a deprecated interface, please implement "
246 "fill_fe_face_values taking a hp::QCollection argument"));
247 (void)cell;
248 (void)face_no;
249 (void)quadrature;
250 (void)internal_data;
251 (void)output_data;
252}
253
254
255
256template <int dim, int spacedim>
257std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
259 const UpdateFlags update_flags,
260 const hp::QCollection<dim - 1> &quadrature) const
261{
262 // base class version, implement overridden function in derived classes
263 return get_face_data(update_flags, quadrature[0]);
264}
265
266
267
268template <int dim, int spacedim>
269std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
271 const UpdateFlags update_flags,
272 const Quadrature<dim - 1> &quadrature) const
273{
274 Assert(false,
275 ExcMessage("Use of a deprecated interface, please implement "
276 "fill_fe_face_values taking a hp::QCollection argument"));
277 (void)update_flags;
278 (void)quadrature;
279
280 return std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>();
281}
282# endif
283
284/* ---------------------------- InternalDataBase --------------------------- */
285
286
287template <int dim, int spacedim>
289 : update_each(update_default)
290{}
291
292
293
294template <int dim, int spacedim>
295void
297 const Quadrature<dim> &)
298{
300}
301
302
303
304template <int dim, int spacedim>
305std::size_t
307{
308 return sizeof(*this);
309}
310#endif
311
312/* ------------------------------ Global functions ------------------------- */
313
314template <int dim, int spacedim>
317{
318 const auto &reference_cells = triangulation.get_reference_cells();
319 Assert(reference_cells.size() == 1,
321 "This function can only work for triangulations that "
322 "use only a single cell type -- for example, only triangles "
323 "or only quadrilaterals. For mixed meshes, there is no "
324 "single linear mapping object that can be used for all "
325 "cells of the triangulation. The triangulation you are "
326 "passing to this function uses multiple cell types."));
327
328 return reference_cells.front()
329 .template get_default_linear_mapping<dim, spacedim>();
330}
331
332
333
334// explicit instantiations
335#include "mapping.inst"
336
337
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 >, ReferenceCells::max_n_vertices< dim >() > 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:315
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
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:316
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
constexpr unsigned int max_n_vertices()
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