Reference documentation for deal.II version 9.0.0
mapping_c1.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 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 #include <deal.II/base/std_cxx14/memory.h>
17 
18 #include <deal.II/grid/tria_accessor.h>
19 #include <deal.II/grid/tria_iterator.h>
20 #include <deal.II/grid/tria_boundary.h>
21 
22 #include <deal.II/fe/mapping_c1.h>
23 
24 #include <cmath>
25 
26 DEAL_II_NAMESPACE_OPEN
27 
28 
29 
30 template <int dim, int spacedim>
32  :
33  MappingQGeneric<dim,spacedim> (3)
34 {}
35 
36 
37 
38 template <int dim, int spacedim>
40  :
41  MappingQ<dim,spacedim> (3)
42 {
43  Assert (dim > 1, ExcImpossibleInDim(dim));
44 
45  // replace the mapping_qp objects of the base class by something
46  // that knows about generating data points based on the geometry
47  //
48  // we only need to replace the Qp mapping because that's the one that's
49  // used on boundary cells where it matters
50  this->qp_mapping = std::make_shared<MappingC1<dim,spacedim>::MappingC1Generic>();
51 }
52 
53 
54 
55 template <>
56 void
58  std::vector<Point<1> > &) const
59 {
60  const unsigned int dim = 1;
61  (void)dim;
62  Assert (dim > 1, ExcImpossibleInDim(dim));
63 }
64 
65 
66 
67 template <>
68 void
70  std::vector<Point<2> > &a) const
71 {
72  const unsigned int dim = 2;
73  const std::array<double, 2> interior_gl_points
74  {{0.5 - 0.5*std::sqrt(1.0/5.0), 0.5 + 0.5*std::sqrt(1.0/5.0)}};
75 
76  // loop over each of the lines, and if it is at the boundary, then first get
77  // the boundary description and second compute the points on it. if not at
78  // the boundary, get the respective points from another function
79  for (unsigned int line_no=0; line_no<GeometryInfo<dim>::lines_per_cell; ++line_no)
80  {
81  const Triangulation<dim>::line_iterator line = cell->line(line_no);
82 
83  if (line->at_boundary())
84  {
85  // first get the normal vectors at the two vertices of this line
86  // from the boundary description
87  const Manifold<dim> &manifold = line->get_manifold();
88 
89  Manifold<dim>::FaceVertexNormals face_vertex_normals;
90  manifold.get_normals_at_vertices (line, face_vertex_normals);
91 
92  // then transform them into interpolation points for a cubic
93  // polynomial
94  //
95  // for this, note that if we describe the boundary curve as a
96  // polynomial in tangential coordinate @p{t=0..1} (along the line)
97  // and @p{s} in normal direction, then the cubic mapping is such
98  // that @p{s = a*t**3 + b*t**2 + c*t + d}, and we want to determine
99  // the interpolation points at @p{t=0.276} and @p{t=0.724}
100  // (Gauss-Lobatto points). Since at @p{t=0,1} we want a vertex which
101  // is actually at the boundary, we know that @p{d=0} and @p{a=-b-c},
102  // which gives @p{s(0.276)} and @p{s(0.726)} in terms of @p{b,c}. As
103  // side-conditions, we want that the derivatives at @p{t=0} and
104  // @p{t=1}, i.e. at the vertices match those returned by the
105  // boundary.
106  //
107  // The task is then first to determine the coefficients from the
108  // tangentials. for that, first rotate the tangents of @p{s(t)} into
109  // the global coordinate system. they are @p{A (1,c)} and @p{A
110  // (1,-b-2c)} with @p{A} the rotation matrix, since the tangentials
111  // in the coordinate system relative to the line are @p{(1,c)} and
112  // @p{(1,-b-2c)} at the two vertices, respectively. We then have to
113  // make sure by matching @p{b,c} that these tangentials are
114  // orthogonal to the normals returned by the boundary object
115  const Tensor<1,2> coordinate_vector = line->vertex(1) - line->vertex(0);
116  const double h = std::sqrt(coordinate_vector * coordinate_vector);
117  Tensor<1,2> coordinate_axis = coordinate_vector;
118  coordinate_axis /= h;
119 
120  const double alpha = std::atan2(coordinate_axis[1], coordinate_axis[0]);
121  const double c = -((face_vertex_normals[0][1] * std::sin(alpha)
122  +face_vertex_normals[0][0] * std::cos(alpha)) /
123  (face_vertex_normals[0][1] * std::cos(alpha)
124  -face_vertex_normals[0][0] * std::sin(alpha)));
125  const double b = ((face_vertex_normals[1][1] * std::sin(alpha)
126  +face_vertex_normals[1][0] * std::cos(alpha)) /
127  (face_vertex_normals[1][1] * std::cos(alpha)
128  -face_vertex_normals[1][0] * std::sin(alpha)))
129  -2*c;
130 
131  const double t1 = interior_gl_points[0];
132  const double t2 = interior_gl_points[1];
133  const double s_t1 = (((-b-c)*t1+b)*t1+c)*t1;
134  const double s_t2 = (((-b-c)*t2+b)*t2+c)*t2;
135 
136  // next evaluate the so determined cubic polynomial at the points
137  // 1/3 and 2/3, first in unit coordinates
138  const Point<2> new_unit_points[2] = { Point<2>(t1, s_t1),
139  Point<2>(t2, s_t2)
140  };
141  // then transform these points to real coordinates by rotating,
142  // scaling and shifting
143  for (unsigned int i=0; i<2; ++i)
144  {
145  Point<2> real_point (std::cos(alpha) * new_unit_points[i][0]
146  - std::sin(alpha) * new_unit_points[i][1],
147  std::sin(alpha) * new_unit_points[i][0]
148  + std::cos(alpha) * new_unit_points[i][1]);
149  real_point *= h;
150  real_point += line->vertex(0);
151  a.push_back (real_point);
152  };
153  }
154  else
155  // not at boundary, so just use scaled Gauss-Lobatto points (i.e., use
156  // plain straight lines).
157  {
158  // Note that the zeroth Gauss-Lobatto point is a boundary point, so
159  // we push back mapped versions of the first and second.
160  a.push_back((1.0 - interior_gl_points[0])*line->vertex(0) +
161  (interior_gl_points[0]*line->vertex(1)));
162  a.push_back((1.0 - interior_gl_points[1])*line->vertex(0) +
163  (interior_gl_points[1]*line->vertex(1)));
164  }
165  }
166 }
167 
168 
169 
170 template <int dim, int spacedim>
171 void
173  std::vector<Point<dim> > &) const
174 {
175  Assert (false, ExcNotImplemented());
176 }
177 
178 
179 
180 template <>
181 void
183  std::vector<Point<1> > &) const
184 {
185  const unsigned int dim = 1;
186  (void)dim;
187  Assert (dim > 2, ExcImpossibleInDim(dim));
188 }
189 
190 
191 
192 template <>
193 void
195  std::vector<Point<2> > &) const
196 {
197  const unsigned int dim = 2;
198  (void)dim;
199  Assert (dim > 2, ExcImpossibleInDim(dim));
200 }
201 
202 
203 
204 template <int dim, int spacedim>
205 void
207  std::vector<Point<dim> > &) const
208 {
209  Assert (false, ExcNotImplemented());
210 }
211 
212 
213 
214 template <int dim, int spacedim>
215 std::unique_ptr<Mapping<dim,spacedim> >
217 {
218  return std_cxx14::make_unique<MappingC1<dim,spacedim>>();
219 }
220 
221 
222 
223 
224 // explicit instantiations
225 #include "mapping_c1.inst"
226 
227 
228 DEAL_II_NAMESPACE_CLOSE
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
Definition: mapping_c1.cc:216
Definition: point.h:104
virtual void add_quad_support_points(const typename Triangulation< dim >::cell_iterator &cell, std::vector< Point< dim > > &a) const override
Definition: mapping_c1.cc:206
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1142
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
Definition: mpi.h:53
std::shared_ptr< const MappingQGeneric< dim, spacedim > > qp_mapping
Definition: mapping_q.h:375
virtual void add_line_support_points(const typename Triangulation< dim >::cell_iterator &cell, std::vector< Point< dim > > &a) const override
Definition: mapping_c1.cc:172
static ::ExceptionBase & ExcNotImplemented()
virtual void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, FaceVertexNormals &face_vertex_normals) const
Definition: manifold.cc:291
const Manifold< dim, spacedim > & get_manifold(const types::manifold_id number) const
Definition: tria.cc:9083