]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Deprecate constructor and replace it by another one.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Sat, 25 Jun 2016 05:11:25 +0000 (00:11 -0500)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Sat, 25 Jun 2016 05:11:48 +0000 (00:11 -0500)
Also update the documentation in a number of places.

doc/news/changes.h
include/deal.II/fe/mapping_q1_eulerian.h
source/fe/mapping_q1_eulerian.cc

index cef4132d1ad8e0707a5a0a1686ba5faf4f85b044..f13477b36d7e6edb4c586c4eff5066f4f3e35ab9 100644 (file)
@@ -139,6 +139,13 @@ inconvenience this causes.
 <h3>General</h3>
 
 <ol>
+ <li> Changed: As MappingQEulerian before, MappingQ1Eulerian has gained
+ a second constructor that reverts the order of the arguments to indicate
+ which DoFHandler a vector is based on. The old constructor is now
+ deprecated and will be removed in a future version.
+ <br>
+ (Wolfgang Bangerth, 2016/06/25)
+ </li>
 
  <li> New: Add new classes to expand a scalar finite element solution into
  the orthogonal bases FESeries::Fourier and FESeries::Legendre. Also
index efc3f6a65424529a931298b41e7e016c807c17dc..69f7302f828eff79faada24a55fa2edb4c74638c 100644 (file)
@@ -31,7 +31,9 @@ template <typename> class Vector;
 /*@{*/
 
 /**
- * Eulerian mapping of general unit cells by $d$-linear shape functions. Each
+ * This class provides a mapping that adds to the location of each cell
+ * a $d$-linear displacement field. (The generalization to higher order
+ * polynomials is provided in the MappingQEulerian class.) Each
  * cell is thus shifted in space by values given to the mapping through a
  * finite element field.
  *
@@ -90,19 +92,29 @@ template <int dim, typename VectorType = Vector<double>, int spacedim=dim >
 class MappingQ1Eulerian : public MappingQGeneric<dim,spacedim>
 {
 public:
+  /**
+   * Constructor.
+   *
+   * @param[in] euler_dof_handler A DoFHandler object that defines a finite
+   * element space. This space needs to have exactly dim components
+   * and these will be considered displacements
+   * relative to the original positions of the cells of the triangulation.
+   * This DoFHandler must be based on a <code>FESystem(FE_Q(1),dim)</code>
+   * finite element.
+   * @param[in] euler_vector A finite element function in the space defined by
+   * the first argument. The dim components of this function will be
+   * interpreted as the displacement we use in defining the mapping, relative
+   * to the location of cells of the underlying triangulation.
+   */
+  MappingQ1Eulerian (const DoFHandler<dim,spacedim> &euler_dof_handler,
+                     const VectorType               &euler_vector);
 
   /**
-   * Constructor. It takes a <tt>Vector<double> &</tt> as its first argument
-   * to specify the transformation of the whole problem from the reference to
-   * the current configuration. The organization of the elements in the @p
-   * Vector must follow the concept how deal.II stores solutions that are
-   * associated to a triangulation.  This is automatically the case if the @p
-   * Vector represents the solution of the previous step of a nonlinear
-   * problem. Alternatively, the @p Vector can be initialized by
-   * <tt>DoFAccessor::set_dof_values()</tt>.
+   * @deprecated Use the constructor with the reverse order of first and
+   * second argument.
    */
-  MappingQ1Eulerian (const VectorType  &euler_transform_vectors,
-                     const DoFHandler<dim,spacedim> &shiftmap_dof_handler);
+  MappingQ1Eulerian (const VectorType               &euler_vector,
+                     const DoFHandler<dim,spacedim> &euler_dof_handler) DEAL_II_DEPRECATED;
 
   /**
    * Return the mapped vertices of the cell. For the current class, this
index 35919e3493a856b204bbfb100fbd3270159cf3ae..0b132b8baf0f7e5f4298b1f732f4d9219140ddaf 100644 (file)
 DEAL_II_NAMESPACE_OPEN
 
 
+template <int dim, class EulerVectorType, int spacedim>
+MappingQ1Eulerian<dim, EulerVectorType, spacedim>::
+MappingQ1Eulerian (const DoFHandler<dim,spacedim> &shiftmap_dof_handler,
+                   const EulerVectorType  &euler_transform_vectors)
+  :
+  MappingQGeneric<dim,spacedim>(1),
+  euler_transform_vectors(&euler_transform_vectors),
+  shiftmap_dof_handler(&shiftmap_dof_handler)
+{}
+
+
 template <int dim, class EulerVectorType, int spacedim>
 MappingQ1Eulerian<dim, EulerVectorType, spacedim>::
 MappingQ1Eulerian (const EulerVectorType  &euler_transform_vectors,

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.