]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Bugfix: Use correct number type in dof_tools_constraints.(cc|inst.in)
authorMatthias Maier <tamiko@43-1.org>
Mon, 28 May 2018 19:58:42 +0000 (14:58 -0500)
committerMatthias Maier <tamiko@43-1.org>
Wed, 6 Jun 2018 15:19:38 +0000 (10:19 -0500)
DoFTools::make_hanging_node_constraints always uses double internally,
but is also instantiated for std::complex<double>

source/dofs/dof_tools_constraints.cc
source/dofs/dof_tools_constraints.inst.in

index d9d06f5eabe487c850ce3033d255b7c0a092c0ec..1b055b046efb30e79d0375dea3e495ccc78cf3f6 100644 (file)
@@ -480,13 +480,13 @@ namespace DoFTools
        * It also suppresses very small entries in the constraint matrix to
        * avoid making the sparsity pattern fuller than necessary.
        */
-      template <typename number>
+      template <typename number1, typename number2>
       void
       filter_constraints(
         const std::vector<types::global_dof_index> &master_dofs,
         const std::vector<types::global_dof_index> &slave_dofs,
-        const FullMatrix<number> &                  face_constraints,
-        AffineConstraints<number> &                 constraints)
+        const FullMatrix<number1> &                 face_constraints,
+        AffineConstraints<number2> &                constraints)
       {
         Assert(face_constraints.n() == master_dofs.size(),
                ExcDimensionMismatch(master_dofs.size(), face_constraints.n()));
@@ -526,7 +526,7 @@ namespace DoFTools
                 {
                   // add up the absolute values of all constraints in this
                   // line to get a measure of their absolute size
-                  number abs_sum = 0;
+                  number1 abs_sum = 0;
                   for (unsigned int i = 0; i < n_master_dofs; ++i)
                     abs_sum += std::abs(face_constraints(row, i));
 
@@ -1050,7 +1050,7 @@ namespace DoFTools
       // a matrix to be used for constraints below. declared here and
       // simply resized down below to avoid permanent re-allocation of
       // memory
-      FullMatrix<number> constraint_matrix;
+      FullMatrix<double> constraint_matrix;
 
       // similarly have arrays that will hold master and slave dof numbers,
       // as well as a scratch array needed for the complicated case below
@@ -1062,9 +1062,9 @@ namespace DoFTools
       // different (or the same) finite elements. we compute them only
       // once, namely the first time they are needed, and then just reuse
       // them
-      Table<2, std::unique_ptr<FullMatrix<number>>> face_interpolation_matrices(
+      Table<2, std::unique_ptr<FullMatrix<double>>> face_interpolation_matrices(
         n_finite_elements(dof_handler), n_finite_elements(dof_handler));
-      Table<3, std::unique_ptr<FullMatrix<number>>>
+      Table<3, std::unique_ptr<FullMatrix<double>>>
         subface_interpolation_matrices(
           n_finite_elements(dof_handler),
           n_finite_elements(dof_handler),
@@ -1075,7 +1075,7 @@ namespace DoFTools
       // these two matrices are derived from the face interpolation matrix
       // as described in the @ref hp_paper "hp paper"
       Table<2,
-            std::unique_ptr<std::pair<FullMatrix<number>, FullMatrix<number>>>>
+            std::unique_ptr<std::pair<FullMatrix<double>, FullMatrix<double>>>>
         split_face_interpolation_matrices(n_finite_elements(dof_handler),
                                           n_finite_elements(dof_handler));
 
@@ -1349,13 +1349,13 @@ namespace DoFTools
                           split_face_interpolation_matrices
                             [dominating_fe_index][cell->active_fe_index()]);
 
-                        const FullMatrix<number>
+                        const FullMatrix<double>
                           &restrict_mother_to_virtual_master_inv =
                             (split_face_interpolation_matrices
                                [dominating_fe_index][cell->active_fe_index()]
                                  ->first);
 
-                        const FullMatrix<number>
+                        const FullMatrix<double>
                           &restrict_mother_to_virtual_slave =
                             (split_face_interpolation_matrices
                                [dominating_fe_index][cell->active_fe_index()]
@@ -1440,7 +1440,7 @@ namespace DoFTools
                               subface_interpolation_matrices
                                 [dominating_fe_index][subface_fe_index][sf]);
 
-                            const FullMatrix<number>
+                            const FullMatrix<double>
                               &restrict_subface_to_virtual = *(
                                 subface_interpolation_matrices
                                   [dominating_fe_index][subface_fe_index][sf]);
@@ -1659,13 +1659,13 @@ namespace DoFTools
                                 [dominating_fe_index][cell->active_fe_index()]);
 
                             const FullMatrix<
-                              number> &restrict_mother_to_virtual_master_inv =
+                              double> &restrict_mother_to_virtual_master_inv =
                               (split_face_interpolation_matrices
                                  [dominating_fe_index][cell->active_fe_index()]
                                    ->first);
 
                             const FullMatrix<
-                              number> &restrict_mother_to_virtual_slave =
+                              double> &restrict_mother_to_virtual_slave =
                               (split_face_interpolation_matrices
                                  [dominating_fe_index][cell->active_fe_index()]
                                    ->second);
@@ -1724,7 +1724,7 @@ namespace DoFTools
                                 [dominating_fe_index]
                                 [neighbor->active_fe_index()]);
 
-                            const FullMatrix<number>
+                            const FullMatrix<double>
                               &restrict_secondface_to_virtual =
                                 *(face_interpolation_matrices
                                     [dominating_fe_index]
index 2a6443182c7cd33d31f48d05e1968885d9ac18af..a9aae2a9629be0839073a20a4e7426f1a1f12665 100644 (file)
@@ -20,6 +20,10 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : DIMENSIONS)
     template void DoFTools::make_hanging_node_constraints(
       const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       AffineConstraints<double> &);
+
+    template void DoFTools::make_hanging_node_constraints(
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
+      AffineConstraints<std::complex<double>> &);
 #endif
   }
 
@@ -29,6 +33,10 @@ for (DH : DOFHANDLER_TEMPLATES; deal_II_dimension : DIMENSIONS)
       const DH<deal_II_dimension, deal_II_dimension> &,
       AffineConstraints<double> &);
 
+    template void DoFTools::make_hanging_node_constraints(
+      const DH<deal_II_dimension, deal_II_dimension> &,
+      AffineConstraints<std::complex<double>> &);
+
 #if deal_II_dimension != 1
     template void DoFTools::make_periodicity_constraints(
       const DH<deal_II_dimension>::face_iterator &,

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.