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;
* @{
*/
- // 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>
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.
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
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.
||
(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),
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),
// 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);
}
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 &)
{
-
-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 */
+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