]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Move transform_unit_to_real_cell from MappingQ1 to MappingQGeneric.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Tue, 8 Sep 2015 15:24:18 +0000 (10:24 -0500)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Thu, 10 Sep 2015 19:36:45 +0000 (14:36 -0500)
include/deal.II/fe/mapping_q.h
include/deal.II/fe/mapping_q1.h
include/deal.II/fe/mapping_q_generic.h
source/fe/mapping_q.cc
source/fe/mapping_q1.cc
source/fe/mapping_q_generic.cc

index afe0414df2178428a715ea5c7015a8b394c75ea7..96bce39faff817c84c518c39ede1041c41b8d332 100644 (file)
@@ -450,9 +450,19 @@ protected:
   const FE_Q<dim> feq;
 
   /**
-   * Pointer to a Q1 mapping to be used on interior cells unless
+   * Pointer to a Q1 mapping. This mapping is used on interior cells unless
    * use_mapping_q_on_all_cells was set in the call to the
-   * constructor.
+   * constructor. The mapping is also used on any cell in the
+   * transform_real_to_unit_cell() to compute a cheap initial
+   * guess for the position of the point before we employ the
+   * more expensive Newton iteration using the full mapping.
+   *
+   * @note MappingQEulerian resets this pointer to an object of type
+   *   MappingQ1Eulerian to ensure that the Q1 mapping also knows
+   *   about the proper shifts and transformations of the Eulerian
+   *   displacements. This also means that we really need to store
+   *   our own Q1 mapping here, rather than simply resorting to
+   *   StaticMappingQ1::mapping.
    */
   std_cxx11::unique_ptr<const MappingQ1<dim,spacedim> > q1_mapping;
 
index 4ec4d7b928fe1e705af726dd90a926bdb8e990f7..fae3b12478cd3b15bcdfc539bbe18abef55d4d94 100644 (file)
@@ -69,12 +69,6 @@ public:
    * @{
    */
 
-  // for documentation, see the Mapping base class
-  virtual
-  Point<spacedim>
-  transform_unit_to_real_cell (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
-                               const Point<dim>                                 &p) const;
-
   // for documentation, see the Mapping base class
   virtual
   Point<dim>
@@ -146,25 +140,6 @@ protected:
   transform_real_to_unit_cell_initial_guess (const std::vector<Point<spacedim> > &vertex,
                                              const Point<spacedim>                            &p) const;
 
-
-  /**
-   * Transforms a point @p p on the unit cell to the point @p p_real on the
-   * real cell @p cell and returns @p p_real.
-   *
-   * This function is called by @p transform_unit_to_real_cell and multiple
-   * times (through the Newton iteration) by @p
-   * transform_real_to_unit_cell_internal.
-   *
-   * Takes a reference to an @p InternalData that must already include the
-   * shape values at point @p p and the mapping support points of the cell.
-   *
-   * This @p InternalData argument avoids multiple computations of the shape
-   * values at point @p p and especially multiple computations of the mapping
-   * support points.
-   */
-  Point<spacedim>
-  transform_unit_to_real_cell_internal (const InternalData &mdata) const;
-
   /**
    * Transforms the point @p p on the real cell to the corresponding point on
    * the unit cell @p cell by a Newton iteration.
index cb5df7fd6012007735ab8b0a8cad5b63c168c182..5eae7b0f9e726e8d8756401bf2ad26672ffffb6f 100644 (file)
@@ -78,6 +78,20 @@ public:
   virtual
   bool preserves_vertex_locations () const;
 
+  /**
+   * @name Mapping points between reference and real cells
+   * @{
+   */
+
+  // for documentation, see the Mapping base class
+  virtual
+  Point<spacedim>
+  transform_unit_to_real_cell (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
+                               const Point<dim>                                 &p) const;
+
+  /**
+   * @}
+   */
 
   /**
    * @name Functions to transform tensors from reference to real coordinates
@@ -454,6 +468,24 @@ protected:
   compute_mapping_support_points (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
                                   std::vector<Point<spacedim> > &a) const = 0;
 
+  /**
+   * Transforms a point @p p on the unit cell to the point @p p_real on the
+   * real cell @p cell and returns @p p_real.
+   *
+   * This function is called by @p transform_unit_to_real_cell and multiple
+   * times (through the Newton iteration) by @p
+   * transform_real_to_unit_cell_internal.
+   *
+   * Takes a reference to an @p InternalData that must already include the
+   * shape values at point @p p and the mapping support points of the cell.
+   *
+   * This @p InternalData argument avoids multiple computations of the shape
+   * values at point @p p and especially multiple computations of the mapping
+   * support points.
+   */
+  Point<spacedim>
+  transform_unit_to_real_cell_internal (const InternalData &mdata) const;
+
   /**
    * Make MappingQ a friend since it needs to call the
    * fill_fe_values() functions on its MappingQ1 sub-object.
index 4447272c2b64d4e01733f181f2c6016e40c20b2a..13443f4338c871fb3cdffc6b2633273a142fbf51 100644 (file)
@@ -75,11 +75,9 @@ MappingQ<dim,spacedim>::MappingQ (const unsigned int degree,
                               ||
                               (dim != spacedim)),
   feq(degree),
-  q1_mapping (this->use_mapping_q_on_all_cells
-              ?
-              0
-              :
-              new MappingQ1<dim,spacedim>()),
+  // create a Q1 mapping for use on interior cells (if necessary)
+  // or to create a good initial guess in transform_real_to_unit_cell()
+  q1_mapping (new MappingQ1<dim,spacedim>()),
   line_support_points(degree+1)
 {
   Assert(n_inner+n_outer==Utilities::fixed_power<dim>(degree+1),
@@ -105,11 +103,9 @@ MappingQ<dim,spacedim>::MappingQ (const MappingQ<dim,spacedim> &mapping)
   n_outer(mapping.n_outer),
   use_mapping_q_on_all_cells (mapping.use_mapping_q_on_all_cells),
   feq(mapping.get_degree()),
-  q1_mapping (mapping.q1_mapping != 0
-              ?
-              dynamic_cast<MappingQ1<dim,spacedim>*>(mapping.q1_mapping->clone())
-              :
-              0),
+  // clone the Q1 mapping for use on interior cells (if necessary)
+  // or to create a good initial guess in transform_real_to_unit_cell()
+  q1_mapping (dynamic_cast<MappingQ1<dim,spacedim>*>(mapping.q1_mapping->clone())),
   line_support_points(mapping.line_support_points)
 {
   Assert(n_inner+n_outer==Utilities::fixed_power<dim>(this->polynomial_degree+1),
@@ -1001,19 +997,9 @@ transform_unit_to_real_cell (const typename Triangulation<dim,spacedim>::cell_it
   // mapping, then either use our own facilities or that of the Q1
   // mapping we store
   if (use_mapping_q_on_all_cells || cell->has_boundary_lines())
-    {
-      const Quadrature<dim> point_quadrature(p);
-      std_cxx11::unique_ptr<InternalData> mdata (get_data(update_quadrature_points,
-                                                          point_quadrature));
-
-      compute_mapping_support_points(cell, mdata->mapping_support_points);
-
-      return this->transform_unit_to_real_cell_internal(*mdata);
-    }
+    return this->MappingQGeneric<dim,spacedim>::transform_unit_to_real_cell (cell, p);
   else
-    {
-      return q1_mapping->transform_unit_to_real_cell (cell, p);
-    }
+    return q1_mapping->transform_unit_to_real_cell (cell, p);
 }
 
 
@@ -1036,8 +1022,7 @@ transform_real_to_unit_cell (const typename Triangulation<dim,spacedim>::cell_it
   Point<dim> initial_p_unit;
   try
     {
-      initial_p_unit
-        = StaticMappingQ1<dim,spacedim>::mapping.transform_real_to_unit_cell(cell, p);
+      initial_p_unit = q1_mapping->transform_real_to_unit_cell(cell, p);
     }
   catch (const typename Mapping<dim,spacedim>::ExcTransformationFailed &)
     {
index e3e1777a16be9e26ba71edc0497d5ba73c53d625..f52b098027191c2b56280a05ae65a7b354ab1df4 100644 (file)
@@ -205,60 +205,6 @@ MappingQ1<dim,spacedim>::compute_mapping_support_points(
 
 
 
-
-template<int dim, int spacedim>
-Point<spacedim>
-MappingQ1<dim,spacedim>::transform_unit_to_real_cell (
-  const typename Triangulation<dim,spacedim>::cell_iterator &cell,
-  const Point<dim> &p) const
-{
-  const Quadrature<dim> point_quadrature(p);
-
-  //TODO: Use get_data() here once MappingQ is no longer derived from
-  //MappingQ1. this doesn't currently work because we here really need
-  //a Q1 InternalData, but MappingQGeneric produces one with the
-  //polynomial degree of the MappingQ
-  std_cxx11::unique_ptr<InternalData> mdata (new InternalData(1));
-  mdata->initialize (this->requires_update_flags (update_quadrature_points),
-                     point_quadrature, 1);
-
-  // compute the mapping support
-  // points
-  compute_mapping_support_points(cell, mdata->mapping_support_points);
-
-  // Mapping support points can be
-  // bigger than necessary. If this
-  // is the case, force them to be
-  // Q1.
-  if (mdata->mapping_support_points.size() > mdata->shape_values.size())
-    mdata->mapping_support_points.resize
-    (GeometryInfo<dim>::vertices_per_cell);
-
-
-  return transform_unit_to_real_cell_internal(*mdata);
-}
-
-
-
-template<int dim, int spacedim>
-Point<spacedim>
-MappingQ1<dim,spacedim>::
-transform_unit_to_real_cell_internal (const InternalData &data) const
-{
-  AssertDimension (data.shape_values.size(),
-                   data.mapping_support_points.size());
-
-  // use now the InternalData to
-  // compute the point in real space.
-  Point<spacedim> p_real;
-  for (unsigned int i=0; i<data.mapping_support_points.size(); ++i)
-    p_real += data.mapping_support_points[i] * data.shape(0,i);
-
-  return p_real;
-}
-
-
-
 /* For an explanation of the  KA and Kb
    arrays see the comments in the declaration of
    transform_real_to_unit_cell_initial_guess */
index 06a42dd06fa3de7043e89c8aca751e967ebfec5f..751fec4621abe5fb05d3ae07fde63b9232f9dc3d 100644 (file)
@@ -785,6 +785,49 @@ MappingQGeneric<dim,spacedim>::get_degree() const
 
 
 
+template<int dim, int spacedim>
+Point<spacedim>
+MappingQGeneric<dim,spacedim>::
+transform_unit_to_real_cell (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
+                             const Point<dim> &p) const
+{
+  //TODO: This function is inefficient -- it might as well evaluate
+  //the shape functions directly, rather than first building the
+  //InternalData object and then working with it
+
+  const Quadrature<dim> point_quadrature(p);
+  std_cxx11::unique_ptr<InternalData> mdata (new InternalData(polynomial_degree));
+  mdata->initialize (this->requires_update_flags(update_quadrature_points),
+                     point_quadrature,
+                     1);
+
+  // compute the mapping support
+  // points
+  compute_mapping_support_points(cell, mdata->mapping_support_points);
+  return transform_unit_to_real_cell_internal(*mdata);
+}
+
+
+
+template<int dim, int spacedim>
+Point<spacedim>
+MappingQGeneric<dim,spacedim>::
+transform_unit_to_real_cell_internal (const InternalData &data) const
+{
+  AssertDimension (data.shape_values.size(),
+                   data.mapping_support_points.size());
+
+  // use now the InternalData to
+  // compute the point in real space.
+  Point<spacedim> p_real;
+  for (unsigned int i=0; i<data.mapping_support_points.size(); ++i)
+    p_real += data.mapping_support_points[i] * data.shape(0,i);
+
+  return p_real;
+}
+
+
+
 template<int dim, int spacedim>
 UpdateFlags
 MappingQGeneric<dim,spacedim>::requires_update_flags (const UpdateFlags in) const

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.