]> https://gitweb.dealii.org/ - dealii.git/commitdiff
RCM: Update fe: fe, fe_tools
authorMatthias Maier <tamiko@43-1.org>
Fri, 25 May 2018 21:11:51 +0000 (16:11 -0500)
committerMatthias Maier <tamiko@43-1.org>
Wed, 6 Jun 2018 15:19:38 +0000 (10:19 -0500)
include/deal.II/fe/fe.h
include/deal.II/fe/fe_tools.h
include/deal.II/fe/fe_tools.templates.h
include/deal.II/fe/fe_tools_extrapolate.templates.h
include/deal.II/fe/fe_tools_interpolate.templates.h
source/fe/fe_tools_extrapolate.inst.in
source/fe/fe_tools_interpolate.inst.in

index 24fc67c069c3a00da668282b8a1955d3f5180344..0b242bb903e5275a11de77357ed64fe9b879ca83 100644 (file)
@@ -412,7 +412,7 @@ class FESystem;
  * with some weights, then the weights have to be exactly the same as those
  * for constrained nodes on the three other edges with respect to the
  * corresponding nodes on these edges. If this isn't the case, you will get
- * into trouble with the ConstraintMatrix class that is the primary consumer
+ * into trouble with the AffineConstraints class that is the primary consumer
  * of the constraint information: while that class is able to handle
  * constraints that are entered more than once (as is necessary for the case
  * above), it insists that the weights are exactly the same.
@@ -421,9 +421,9 @@ class FESystem;
  * parent face degrees of freedom that contain those on the edges of the
  * parent face; it is possible that some of them are in turn constrained
  * themselves, leading to longer chains of constraints that the
- * ConstraintMatrix class will eventually have to sort out. (The constraints
+ * AffineConstraints class will eventually have to sort out. (The constraints
  * described above are used by the DoFTools::make_hanging_node_constraints()
- * function that constructs a ConstraintMatrix object.) However, this is of no
+ * function that constructs an AffineConstraints object.) However, this is of no
  * concern for the FiniteElement and derived classes since they only act
  * locally on one cell and its immediate neighbor, and do not see the bigger
  * picture. The
index a99b256ef8d590cb8abb8b499e4320877be0bb46..ae6c8e4d54d0be1ec5753afff2c51c777ae2e5ba 100644 (file)
@@ -49,7 +49,8 @@ template <int dim, int spacedim>
 class DoFHandler;
 template <int dim>
 class FiniteElementData;
-class ConstraintMatrix;
+template <typename number>
+class AffineConstraints;
 
 
 
@@ -643,9 +644,10 @@ namespace FETools
    * course Q1 on each cell.
    *
    * For this case (continuous elements on grids with hanging nodes), please
-   * use the @p interpolate() function with an additional ConstraintMatrix
-   * argument, see below, or make the field conforming yourself by calling the
-   * @p distribute function of your hanging node constraints object.
+   * use the @p interpolate() function with an additional AffineConstraints
+   * object as argument, see below, or make the field conforming yourself
+   * by calling the @p distribute function of your hanging node constraints
+   * object.
    */
   template <int dim,
             int spacedim,
@@ -682,11 +684,12 @@ namespace FETools
             class InVector,
             class OutVector>
   void
-  interpolate(const DoFHandlerType1<dim, spacedim> &dof1,
-              const InVector &                      u1,
-              const DoFHandlerType2<dim, spacedim> &dof2,
-              const ConstraintMatrix &              constraints,
-              OutVector &                           u2);
+  interpolate(
+    const DoFHandlerType1<dim, spacedim> &                   dof1,
+    const InVector &                                         u1,
+    const DoFHandlerType2<dim, spacedim> &                   dof2,
+    const AffineConstraints<typename OutVector::value_type> &constraints,
+    OutVector &                                              u2);
 
   /**
    * Compute the interpolation of the @p fe1-function @p u1 to a @p
@@ -695,7 +698,7 @@ namespace FETools
    *
    * Note, that this function does not work on continuous elements at hanging
    * nodes. For that case use the @p back_interpolate function, below, that
-   * takes an additional @p ConstraintMatrix object.
+   * takes an additional @p AffineConstraints object.
    *
    * @p dof1 might be a DoFHandler or a hp::DoFHandler onject.
    *
@@ -728,12 +731,13 @@ namespace FETools
    */
   template <int dim, class InVector, class OutVector, int spacedim>
   void
-  back_interpolate(const DoFHandler<dim, spacedim> &dof1,
-                   const ConstraintMatrix &         constraints1,
-                   const InVector &                 u1,
-                   const DoFHandler<dim, spacedim> &dof2,
-                   const ConstraintMatrix &         constraints2,
-                   OutVector &                      u1_interpolated);
+  back_interpolate(
+    const DoFHandler<dim, spacedim> &                        dof1,
+    const AffineConstraints<typename OutVector::value_type> &constraints1,
+    const InVector &                                         u1,
+    const DoFHandler<dim, spacedim> &                        dof2,
+    const AffineConstraints<typename OutVector::value_type> &constraints2,
+    OutVector &                                              u1_interpolated);
 
   /**
    * Compute $(Id-I_h)z_1$ for a given @p dof1-function $z_1$, where $I_h$ is
@@ -742,7 +746,7 @@ namespace FETools
    *
    * Note, that this function does not work for continuous elements at hanging
    * nodes. For that case use the @p interpolation_difference function, below,
-   * that takes an additional @p ConstraintMatrix object.
+   * that takes an additional @p AffineConstraints object.
    */
   template <int dim, class InVector, class OutVector, int spacedim>
   void
@@ -765,12 +769,13 @@ namespace FETools
    */
   template <int dim, class InVector, class OutVector, int spacedim>
   void
-  interpolation_difference(const DoFHandler<dim, spacedim> &dof1,
-                           const ConstraintMatrix &         constraints1,
-                           const InVector &                 z1,
-                           const DoFHandler<dim, spacedim> &dof2,
-                           const ConstraintMatrix &         constraints2,
-                           OutVector &                      z1_difference);
+  interpolation_difference(
+    const DoFHandler<dim, spacedim> &                        dof1,
+    const AffineConstraints<typename OutVector::value_type> &constraints1,
+    const InVector &                                         z1,
+    const DoFHandler<dim, spacedim> &                        dof2,
+    const AffineConstraints<typename OutVector::value_type> &constraints2,
+    OutVector &                                              z1_difference);
 
 
 
@@ -842,7 +847,7 @@ namespace FETools
    * @note The resulting field does not satisfy continuity requirements of the
    * given finite elements if the algorithm outlined above is used. When you
    * use continuous elements on grids with hanging nodes, please use the @p
-   * extrapolate function with an additional ConstraintMatrix argument, see
+   * extrapolate function with an additional AffineConstraints argument, see
    * below.
    *
    * @note Since this function operates on patches of cells, it requires that
@@ -870,11 +875,12 @@ namespace FETools
    */
   template <int dim, class InVector, class OutVector, int spacedim>
   void
-  extrapolate(const DoFHandler<dim, spacedim> &dof1,
-              const InVector &                 z1,
-              const DoFHandler<dim, spacedim> &dof2,
-              const ConstraintMatrix &         constraints,
-              OutVector &                      z2);
+  extrapolate(
+    const DoFHandler<dim, spacedim> &                        dof1,
+    const InVector &                                         z1,
+    const DoFHandler<dim, spacedim> &                        dof2,
+    const AffineConstraints<typename OutVector::value_type> &constraints,
+    OutVector &                                              z2);
 
   //@}
   /**
@@ -1405,7 +1411,7 @@ namespace FETools
                  << "You are using continuous elements on a grid with "
                  << "hanging nodes but without providing hanging node "
                  << "constraints. Use the respective function with "
-                 << "additional ConstraintMatrix argument(s), instead."
+                 << "additional AffineConstraints argument(s), instead."
                  << (arg1 ? "" : ""));
   /**
    * You need at least two grid levels.
index 7370518a3a9fe052ec754d0bd2ae9cef1095b3e5..132c4aa5ff77060461295cf95a108be93d2cbfef 100644 (file)
@@ -60,7 +60,6 @@
 
 #include <deal.II/hp/dof_handler.h>
 
-#include <deal.II/lac/constraint_matrix.h>
 #include <deal.II/lac/full_matrix.h>
 #include <deal.II/lac/householder.h>
 
index 1dd3c2a448a54d3bde7ab2aceeb03c7bbf41e5ec..c9a711203b98022f299a2f53aea4150641f98bfa 100644 (file)
@@ -28,8 +28,8 @@
 
 #include <deal.II/grid/tria.h>
 
+#include <deal.II/lac/affine_constraints.h>
 #include <deal.II/lac/block_vector.h>
-#include <deal.II/lac/constraint_matrix.h>
 #include <deal.II/lac/la_parallel_block_vector.h>
 #include <deal.II/lac/la_parallel_vector.h>
 #include <deal.II/lac/la_vector.h>
@@ -1738,7 +1738,7 @@ namespace FETools
               const DoFHandler<dim, spacedim> &dof2,
               OutVector &                      u2)
   {
-    ConstraintMatrix dummy;
+    AffineConstraints<typename OutVector::value_type> dummy;
     dummy.close();
     extrapolate(dof1, u1, dof2, dummy, u2);
   }
@@ -1747,11 +1747,12 @@ namespace FETools
 
   template <int dim, class InVector, class OutVector, int spacedim>
   void
-  extrapolate(const DoFHandler<dim, spacedim> &dof1,
-              const InVector &                 u1,
-              const DoFHandler<dim, spacedim> &dof2,
-              const ConstraintMatrix &         constraints,
-              OutVector &                      u2)
+  extrapolate(
+    const DoFHandler<dim, spacedim> &                        dof1,
+    const InVector &                                         u1,
+    const DoFHandler<dim, spacedim> &                        dof2,
+    const AffineConstraints<typename OutVector::value_type> &constraints,
+    OutVector &                                              u2)
   {
     Assert(dof1.get_fe(0).n_components() == dof2.get_fe(0).n_components(),
            ExcDimensionMismatch(dof1.get_fe(0).n_components(),
index a576829a61b79fc0e7f76042da628060249e3ea0..0b31e77723a6c0a35fbcc54b5e88c3a1166a5188 100644 (file)
@@ -38,8 +38,8 @@
 
 #include <deal.II/hp/dof_handler.h>
 
+#include <deal.II/lac/affine_constraints.h>
 #include <deal.II/lac/block_vector.h>
-#include <deal.II/lac/constraint_matrix.h>
 #include <deal.II/lac/la_parallel_block_vector.h>
 #include <deal.II/lac/la_parallel_vector.h>
 #include <deal.II/lac/la_vector.h>
@@ -69,7 +69,7 @@ namespace FETools
               const DoFHandlerType2<dim, spacedim> &dof2,
               OutVector &                           u2)
   {
-    ConstraintMatrix dummy;
+    AffineConstraints<typename OutVector::value_type> dummy;
     dummy.close();
     interpolate(dof1, u1, dof2, dummy, u2);
   }
@@ -83,11 +83,12 @@ namespace FETools
             class InVector,
             class OutVector>
   void
-  interpolate(const DoFHandlerType1<dim, spacedim> &dof1,
-              const InVector &                      u1,
-              const DoFHandlerType2<dim, spacedim> &dof2,
-              const ConstraintMatrix &              constraints,
-              OutVector &                           u2)
+  interpolate(
+    const DoFHandlerType1<dim, spacedim> &                   dof1,
+    const InVector &                                         u1,
+    const DoFHandlerType2<dim, spacedim> &                   dof2,
+    const AffineConstraints<typename OutVector::value_type> &constraints,
+    OutVector &                                              u2)
   {
     Assert(&dof1.get_triangulation() == &dof2.get_triangulation(),
            ExcTriangulationMismatch());
@@ -321,13 +322,9 @@ namespace FETools
       if ((cell->subdomain_id() == subdomain_id) ||
           (subdomain_id == numbers::invalid_subdomain_id))
         {
-          // For continuous elements on
-          // grids with hanging nodes we
-          // need hanging node
-          // constraints. Consequently,
-          // when the elements are
-          // continuous no hanging node
-          // constraints are allowed.
+          // For continuous elements on grids with hanging nodes we need
+          // hanging node constraints. Consequently, when the elements are
+          // continuous no hanging node constraints are allowed.
           const bool hanging_nodes_not_allowed =
             (cell->get_fe().dofs_per_vertex != 0) || (fe2.dofs_per_vertex != 0);
 
@@ -341,8 +338,7 @@ namespace FETools
 
           const unsigned int dofs_per_cell1 = cell->get_fe().dofs_per_cell;
 
-          // make sure back_interpolation
-          // matrix is available
+          // make sure back_interpolation matrix is available
           if (interpolation_matrices[&cell->get_fe()] == nullptr)
             {
               interpolation_matrices[&cell->get_fe()] =
@@ -373,12 +369,13 @@ namespace FETools
     {
       template <int dim, int spacedim, class InVector>
       void
-      back_interpolate(const DoFHandler<dim, spacedim> &dof1,
-                       const ConstraintMatrix &         constraints1,
-                       const InVector &                 u1,
-                       const DoFHandler<dim, spacedim> &dof2,
-                       const ConstraintMatrix &         constraints2,
-                       InVector &                       u1_interpolated)
+      back_interpolate(
+        const DoFHandler<dim, spacedim> &                       dof1,
+        const AffineConstraints<typename InVector::value_type> &constraints1,
+        const InVector &                                        u1,
+        const DoFHandler<dim, spacedim> &                       dof2,
+        const AffineConstraints<typename InVector::value_type> &constraints2,
+        InVector &                                              u1_interpolated)
       {
         Vector<typename InVector::value_type> u2(dof2.n_dofs());
         interpolate(dof1, u1, dof2, constraints2, u2);
@@ -389,12 +386,15 @@ namespace FETools
 #ifdef DEAL_II_WITH_PETSC
       template <int dim, int spacedim>
       void
-      back_interpolate(const DoFHandler<dim, spacedim> & dof1,
-                       const ConstraintMatrix &          constraints1,
-                       const PETScWrappers::MPI::Vector &u1,
-                       const DoFHandler<dim, spacedim> & dof2,
-                       const ConstraintMatrix &          constraints2,
-                       PETScWrappers::MPI::Vector &      u1_interpolated)
+      back_interpolate(
+        const DoFHandler<dim, spacedim> &dof1,
+        const AffineConstraints<PETScWrappers::MPI::Vector::value_type>
+          &                               constraints1,
+        const PETScWrappers::MPI::Vector &u1,
+        const DoFHandler<dim, spacedim> & dof2,
+        const AffineConstraints<PETScWrappers::MPI::Vector::value_type>
+          &                         constraints2,
+        PETScWrappers::MPI::Vector &u1_interpolated)
       {
         // if u1 is a parallel distributed PETSc vector, we create a
         // vector u2 with based on the sets of locally owned and relevant
@@ -419,12 +419,15 @@ namespace FETools
 #ifdef DEAL_II_WITH_TRILINOS
       template <int dim, int spacedim>
       void
-      back_interpolate(const DoFHandler<dim, spacedim> &    dof1,
-                       const ConstraintMatrix &             constraints1,
-                       const TrilinosWrappers::MPI::Vector &u1,
-                       const DoFHandler<dim, spacedim> &    dof2,
-                       const ConstraintMatrix &             constraints2,
-                       TrilinosWrappers::MPI::Vector &      u1_interpolated)
+      back_interpolate(
+        const DoFHandler<dim, spacedim> &dof1,
+        const AffineConstraints<
+          typename TrilinosWrappers::MPI::Vector::value_type> &constraints1,
+        const TrilinosWrappers::MPI::Vector &                  u1,
+        const DoFHandler<dim, spacedim> &                      dof2,
+        const AffineConstraints<
+          typename TrilinosWrappers::MPI::Vector::value_type> &constraints2,
+        TrilinosWrappers::MPI::Vector &                        u1_interpolated)
       {
         // if u1 is a parallel distributed Trilinos vector, we create a
         // vector u2 with based on the sets of locally owned and relevant
@@ -450,10 +453,10 @@ namespace FETools
       void
       back_interpolate(
         const DoFHandler<dim, spacedim> &                 dof1,
-        const ConstraintMatrix &                          constraints1,
+        const AffineConstraints<Number> &                 constraints1,
         const LinearAlgebra::distributed::Vector<Number> &u1,
         const DoFHandler<dim, spacedim> &                 dof2,
-        const ConstraintMatrix &                          constraints2,
+        const AffineConstraints<Number> &                 constraints2,
         LinearAlgebra::distributed::Vector<Number> &      u1_interpolated)
       {
         const IndexSet &dof2_locally_owned_dofs = dof2.locally_owned_dofs();
@@ -477,12 +480,13 @@ namespace FETools
 
   template <int dim, class InVector, class OutVector, int spacedim>
   void
-  back_interpolate(const DoFHandler<dim, spacedim> &dof1,
-                   const ConstraintMatrix &         constraints1,
-                   const InVector &                 u1,
-                   const DoFHandler<dim, spacedim> &dof2,
-                   const ConstraintMatrix &         constraints2,
-                   OutVector &                      u1_interpolated)
+  back_interpolate(
+    const DoFHandler<dim, spacedim> &                        dof1,
+    const AffineConstraints<typename OutVector::value_type> &constraints1,
+    const InVector &                                         u1,
+    const DoFHandler<dim, spacedim> &                        dof2,
+    const AffineConstraints<typename OutVector::value_type> &constraints2,
+    OutVector &                                              u1_interpolated)
   {
     // For discontinuous elements without constraints take the simpler version
     // of the back_interpolate function.
@@ -594,12 +598,13 @@ namespace FETools
     {
       template <int dim, class InVector, class OutVector, int spacedim>
       void
-      interpolation_difference(const DoFHandler<dim, spacedim> &dof1,
-                               const ConstraintMatrix &         constraints1,
-                               const InVector &                 u1,
-                               const DoFHandler<dim, spacedim> &dof2,
-                               const ConstraintMatrix &         constraints2,
-                               OutVector &                      u1_difference)
+      interpolation_difference(
+        const DoFHandler<dim, spacedim> &                        dof1,
+        const AffineConstraints<typename OutVector::value_type> &constraints1,
+        const InVector &                                         u1,
+        const DoFHandler<dim, spacedim> &                        dof2,
+        const AffineConstraints<typename OutVector::value_type> &constraints2,
+        OutVector &                                              u1_difference)
       {
         back_interpolate(
           dof1, constraints1, u1, dof2, constraints2, u1_difference);
@@ -610,12 +615,15 @@ namespace FETools
 #ifdef DEAL_II_WITH_TRILINOS
       template <int dim, int spacedim>
       void
-      interpolation_difference(const DoFHandler<dim, spacedim> &dof1,
-                               const ConstraintMatrix &         constraints1,
-                               const TrilinosWrappers::MPI::Vector &u1,
-                               const DoFHandler<dim, spacedim> &    dof2,
-                               const ConstraintMatrix &       constraints2,
-                               TrilinosWrappers::MPI::Vector &u1_difference)
+      interpolation_difference(
+        const DoFHandler<dim, spacedim> &dof1,
+        const AffineConstraints<TrilinosWrappers::MPI::Vector::value_type>
+          &                                  constraints1,
+        const TrilinosWrappers::MPI::Vector &u1,
+        const DoFHandler<dim, spacedim> &    dof2,
+        const AffineConstraints<TrilinosWrappers::MPI::Vector::value_type>
+          &                            constraints2,
+        TrilinosWrappers::MPI::Vector &u1_difference)
       {
         back_interpolate(
           dof1, constraints1, u1, dof2, constraints2, u1_difference);
@@ -638,12 +646,13 @@ namespace FETools
 
   template <int dim, class InVector, class OutVector, int spacedim>
   void
-  interpolation_difference(const DoFHandler<dim, spacedim> &dof1,
-                           const ConstraintMatrix &         constraints1,
-                           const InVector &                 u1,
-                           const DoFHandler<dim, spacedim> &dof2,
-                           const ConstraintMatrix &         constraints2,
-                           OutVector &                      u1_difference)
+  interpolation_difference(
+    const DoFHandler<dim, spacedim> &                        dof1,
+    const AffineConstraints<typename OutVector::value_type> &constraints1,
+    const InVector &                                         u1,
+    const DoFHandler<dim, spacedim> &                        dof2,
+    const AffineConstraints<typename OutVector::value_type> &constraints2,
+    OutVector &                                              u1_difference)
   {
     // For discontinuous elements
     // without constraints take the
index 05c06edc430432d78639733e44fb3e3097ce45ab..3f09588d6ff1d3c66e428203f32f315ab432ef2d 100644 (file)
@@ -31,7 +31,7 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS;
       extrapolate<deal_II_dimension>(const DoFHandler<deal_II_dimension> &,
                                      const VEC &,
                                      const DoFHandler<deal_II_dimension> &,
-                                     const ConstraintMatrix &,
+                                     const AffineConstraints<VEC::value_type> &,
                                      VEC &);
 #endif
     \}
index 9e62e5b3e9e45cbba79d0996060a56b580b8b204..a301dac43ea5fb9ba5a1fd3590fb07982ceac749 100644 (file)
@@ -20,6 +20,7 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS;
   {
     namespace FETools
     \{
+
 #if deal_II_dimension <= deal_II_space_dimension
       template void
       interpolate<deal_II_dimension, deal_II_space_dimension>(
@@ -33,7 +34,7 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS;
         const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
         const Vector &,
         const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
-        const ConstraintMatrix &,
+        const AffineConstraints<Vector::value_type> &,
         Vector &);
 #endif
     \}
@@ -54,7 +55,7 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
       interpolate<deal_II_dimension>(const hp::DoFHandler<deal_II_dimension> &,
                                      const Vector<double> &,
                                      const hp::DoFHandler<deal_II_dimension> &,
-                                     const ConstraintMatrix &,
+                                     const AffineConstraints<double> &,
                                      Vector<double> &);
 
       template void
@@ -67,7 +68,7 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
       interpolate<deal_II_dimension>(const hp::DoFHandler<deal_II_dimension> &,
                                      const Vector<float> &,
                                      const hp::DoFHandler<deal_II_dimension> &,
-                                     const ConstraintMatrix &,
+                                     const AffineConstraints<float> &,
                                      Vector<float> &);
 #endif
     \}
@@ -94,12 +95,13 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS;
         VEC &);
 
       template void
-      back_interpolate<deal_II_dimension>(const DoFHandler<deal_II_dimension> &,
-                                          const ConstraintMatrix &,
-                                          const VEC &,
-                                          const DoFHandler<deal_II_dimension> &,
-                                          const ConstraintMatrix &,
-                                          VEC &);
+      back_interpolate<deal_II_dimension>(
+        const DoFHandler<deal_II_dimension> &,
+        const AffineConstraints<VEC::value_type> &,
+        const VEC &,
+        const DoFHandler<deal_II_dimension> &,
+        const AffineConstraints<VEC::value_type> &,
+        VEC &);
 
       template void
       interpolation_difference<deal_II_dimension>(
@@ -111,10 +113,10 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS;
       template void
       interpolation_difference<deal_II_dimension>(
         const DoFHandler<deal_II_dimension> &,
-        const ConstraintMatrix &,
+        const AffineConstraints<VEC::value_type> &,
         const VEC &,
         const DoFHandler<deal_II_dimension> &,
-        const ConstraintMatrix &,
+        const AffineConstraints<VEC::value_type> &,
         VEC &);
 
       template void

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.